All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
GaussVolPointStencil.C
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
25  grpGaussVolPoint
26 Class
27  Foam::fvsc::GaussVolPoint::GaussVolPointStencil
28 Description
29  Methods calculation of differential operators without tangential derivatives
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvc.H"
34 #include "GaussVolPointStencil.H"
35 #include "addToRunTimeSelectionTable.H"
36 #include "volPointInterpolation.H"
37 #include "processorFvPatchField.H"
38 
39 namespace Foam
40 {
41 namespace fvsc
42 {
43  defineTypeNameAndDebug(GaussVolPoint,0);
45  (
46  fvscStencil,
48  components
49  );
50 }
51 }
52 
53 // constructors
55 :
56  fvscStencil(io),
57  GaussVolPointBase(refCast<const fvMesh>(io.db()))
58 {
59 }
60 
62 {
63 }
64 
65 //- Calculate gradient of volume scalar function on the faces
66 //
67 // \param iF Internal scalar field.
68 // Allowable values: constant reference to the volScalarField.
69 //
70 // \return Gradient of iF (vector field) which was computed on the faces of mesh.
71 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::GaussVolPoint::Grad(const volScalarField& vF)
72 {
73  const_cast<volScalarField&>(vF).correctBoundaryConditions();
74  tmp<surfaceVectorField> tgradIF(vector::zero * fvc::snGrad(vF));
75  //tmp<surfaceVectorField> tgradIF(nf_ * fvc::snGrad(vF));
76  surfaceVectorField& gradIf = tgradIF.ref();
77 
78  this->faceGrad(vF, gradIf);
79 
80  return tgradIF;
81 };
82 
83 //- Calculate gradient of volume vector field on the faces.
84 //
85 // \param iVF Internal vector field.
86 // Allowable values: constant reference to the volVectorField.
87 //
88 // \return Gradient of iVF (tensor field) which was computed on the faces of mesh.
89 Foam::tmp<Foam::surfaceTensorField> Foam::fvsc::GaussVolPoint::Grad(const volVectorField& iVF)
90 {
91  const_cast<volVectorField&>(iVF).correctBoundaryConditions();
92  tmp<surfaceTensorField> tgradIVF(tensor::zero * fvc::snGrad(iVF.component(0)));
93  //tmp<surfaceTensorField> tgradIVF(nf_ * fvc::snGrad(iVF));
94  surfaceTensorField& gradIVF = tgradIVF.ref();
95 
96  this->faceGrad(iVF, gradIVF);
97 
98  return tgradIVF;
99 };
101 Foam::tmp<Foam::surfaceScalarField> Foam::fvsc::GaussVolPoint::Div(const volVectorField& iVF)
102 {
103  const_cast<volVectorField&>(iVF).correctBoundaryConditions();
104  tmp<surfaceScalarField> tdivIVF(0.0 * fvc::snGrad(iVF.component(0)));
105  //tmp<surfaceScalarField> tdivIVF(nf_ & fvc::snGrad(iVF));
106  surfaceScalarField& divIVF = tdivIVF.ref();
107 
108  this->faceDiv(iVF, divIVF);
109 
110  return tdivIVF;
111 };
112 
113 //- Calculate divergence of volume tensor field on the faces.
114 //
115 // \param iTF Internal tensor field.
116 // Allowable values: constant reference to the volTensorField.
117 //
118 // \return Divergence of iTF (vector field) which was computed on the faces of mesh.
119 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::GaussVolPoint::Div(const volTensorField& iTF)
120 {
121  const_cast<volTensorField&>(iTF).correctBoundaryConditions();
122  tmp<surfaceVectorField> tdivITF(vector::zero*fvc::snGrad(iTF.component(0)));
123  //tmp<surfaceVectorField> tdivITF(nf_*fvc::snGrad(iTF.component(0)));
124  surfaceVectorField& divITF = tdivITF.ref();
125 
126  this->faceDiv(iTF, divITF);
127 
128  return tdivITF;
129 }
130 
131 
132 //
133 //END-OF-FILE
134 //
135 
136 
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
GaussVolPoint(const IOobject &)
Construct from IOobject.
This is a method for approximating derivatives of tangents to a face. They are further used in the ca...
addToRunTimeSelectionTable(fvscStencil, GaussVolPoint, components)
e correctBoundaryConditions()
tmp< surfaceScalarField > Div(const volVectorField &iVF)
defineTypeNameAndDebug(fvscStencil, 0)
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.