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