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 Group grpleastSquaresOpt
25  This group contains common part of QGD solvers.
26 Class
27  Foam::fvsc::leastSquaresOpt::leastSquaresStencilOpt
28 Description
29  Methods for optimal calculating of directional derivative.
30  With parallel realisation.
31 Source Files
32  leastSquaresStencilOpt.C
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef leastSquaresOptStencil_H
36 #define leastSquaresOptStencil_H
37 
38 #include "fvscStencil.H"
39 #include "regIOobject.H"
40 #include "labelList.H"
41 
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 #include "surfaceMesh.H"
45 
46 #include "fvCFD.H"
47 #include "zeroGradientFvPatchFields.H"
48 #include "vector.H"
49 #include "List.H"
50 #include "leastSquaresBase.H"
51 
52 namespace Foam
53 {
54 
55 namespace fvsc
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class leastSquaresOpt Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class leastSquaresOpt
63 :
64  public fvscStencil, public leastSquaresBase
65 {
66 
67 protected:
68 
69  void faceScalarDer(const Field<scalar>& iF,const Field<scalar>& sF,int com, surfaceScalarField& rField);
70  void faceScalarDer(const tmp<Field<scalar>>& tiF,const tmp<Field<scalar>>& tsF, int com, tmp<surfaceScalarField>& trField);
71 
72  void faceScalarDer
73  (
74  const List3<scalar>& procVfValues,
75  const surfaceScalarField& sF,
76  int derComp,
77  surfaceScalarField& rField
78  );
79 
80  void faceScalarDer
81  (
82  const tmp<List3<scalar>>& tprocVfValues,
83  const tmp<surfaceScalarField>& tsF,
84  int derComp,
85  tmp<surfaceScalarField>& trField
86  );
87 
89 
90 public:
91 
92  TypeName ("leastSquaresOpt");
93 
94 // Constructors
95 
96  //- Construct from IOobject. Optional flag for if IOobject is the
97  // top level regIOobject.
98  leastSquaresOpt(const IOobject&);
99 
100  //- Destructor
102 
103  tmp<surfaceVectorField> Grad(const volScalarField& iF);
104  tmp<surfaceVectorField> Grad(const volScalarField& iF, const surfaceScalarField&);
105 
106  tmp<surfaceTensorField> Grad(const volVectorField& iVF);
107  tmp<surfaceTensorField> Grad(const volVectorField& iVF, const surfaceVectorField&);
108 
109  tmp<surfaceScalarField> Div(const volVectorField& iVF);
110 
111  tmp<surfaceVectorField> Div(const volTensorField& iTF);
112 
113 };
114 
115 } //fvsc
116 
117 } //Foam
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 #endif
122 
This is a method for calculation the differential operators without tangential derivatives. They are further used in the calculation of the QGD terms.
Definition: fvscStencil.H:38
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")
Methods calculating of differential operators.
leastSquaresOpt(const IOobject &)
Construct from IOobject. Optional flag for if IOobject is the.