All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
leastSquaresStencilOpt.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 Class
25  Foam::fvsc::leastSquaresOpt::leastSquaresStencilOpt
26 Description
27  Methods for optimal calculating of directional derivative.
28  With parallel realisation.
29 Source Files
30  leastSquaresStencilOpt.C
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef leastSquaresOptStencil_H
34 #define leastSquaresOptStencil_H
35 
36 #include "fvscStencil.H"
37 #include "regIOobject.H"
38 #include "labelList.H"
39 
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "surfaceMesh.H"
43 
44 #include "fvCFD.H"
45 #include "zeroGradientFvPatchFields.H"
46 #include "vector.H"
47 #include "List.H"
48 #include "leastSquaresBase.H"
49 
50 namespace Foam
51 {
52 
53 namespace fvsc
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class leastSquaresOpt Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class leastSquaresOpt
61 :
62  public fvscStencil, public leastSquaresBase
63 {
64 
65 protected:
66 
67  void faceScalarDer(const Field<scalar>& iF,const Field<scalar>& sF,int com, surfaceScalarField& rField);
68  void faceScalarDer(const tmp<Field<scalar>>& tiF,const tmp<Field<scalar>>& tsF, int com, tmp<surfaceScalarField>& trField);
69 
70  void faceScalarDer
71  (
72  const List3<scalar>& procVfValues,
73  const surfaceScalarField& sF,
74  int derComp,
75  surfaceScalarField& rField
76  );
77 
78  void faceScalarDer
79  (
80  const tmp<List3<scalar>>& tprocVfValues,
81  const tmp<surfaceScalarField>& tsF,
82  int derComp,
83  tmp<surfaceScalarField>& trField
84  );
85 
87 
88 public:
89 
90  TypeName ("leastSquaresOpt");
91 
92 // Constructors
93 
94  //- Construct from IOobject. Optional flag for if IOobject is the
95  // top level regIOobject.
96  leastSquaresOpt(const IOobject&);
97 
98  //- Destructor
100 
101  tmp<surfaceVectorField> Grad(const volScalarField& iF);
102  tmp<surfaceVectorField> Grad(const volScalarField& iF, const surfaceScalarField&);
103 
104  tmp<surfaceTensorField> Grad(const volVectorField& iVF);
105  tmp<surfaceTensorField> Grad(const volVectorField& iVF, const surfaceVectorField&);
106 
107  tmp<surfaceScalarField> Div(const volVectorField& iVF);
108 
109  tmp<surfaceVectorField> Div(const volTensorField& iTF);
110 
111 };
112 
113 } //fvsc
114 
115 } //Foam
116 
117 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118 
119 #endif
120 
tmp< surfaceScalarField > Div(const volVectorField &iVF)
Calculate divergence of volume vector field on the faces.
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
void faceScalarDer(const Field< scalar > &iF, const Field< scalar > &sF, int com, surfaceScalarField &rField)
TypeName("leastSquaresOpt")
fvscStencil(const IOobject &io)
Construct from components.
Methods calculating of differential operators.
leastSquaresOpt(const IOobject &)
Construct from IOobject. Optional flag for if IOobject is the.