All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
reducedFaceNormalStencil.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 grpreduced
25  This group contains common part of QGD solvers.
26 Class
27  Foam::fvsc::reduced::reducedFaceNormalStencil
28 Description
29  Methods calculating of differential operators without tangential direvatives
30 \*---------------------------------------------------------------------------*/
31 
32 
33 
34 #include "fvc.H"
36 #include "addToRunTimeSelectionTable.H"
37 
38 namespace Foam
39 {
40 namespace fvsc
41 {
42  defineTypeNameAndDebug(reduced,0);
44  (
45  fvscStencil,
47  components
48  );
49 }
50 }
51 
52 // constructors
53 Foam::fvsc::reduced::reduced(const IOobject& io)
54 :
55  fvscStencil(io)
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::reduced::Grad(const volScalarField& vF)
70 {
71  tmp<surfaceVectorField> tgradIF(nf_ * fvc::snGrad(vF));
72 
73  return tgradIF;
74 };
75 
76 //- Calculate gradient of volume vector field on the faces.
77 //
78 // \param iVF Internal vector field.
79 // Allowable values: constant reference to the volVectorField.
80 //
81 // \return Gradient of iVF (tensor field) which was computed on the faces of mesh.
82 Foam::tmp<Foam::surfaceTensorField> Foam::fvsc::reduced::Grad(const volVectorField& iVF)
83 {
84 
85  tmp<surfaceTensorField> tgradIVF(nf_ * fvc::snGrad(iVF));
86 
87  return tgradIVF;
88 };
89 
90 Foam::tmp<Foam::surfaceScalarField> Foam::fvsc::reduced::Div(const volVectorField& iVF)
91 {
92  tmp<surfaceScalarField> tdivIVF(nf_ & fvc::snGrad(iVF));
93 
94  return tdivIVF;
95 };
96 
97 //- Calculate divergence of volume tensor field on the faces.
98 //
99 // \param iTF Internal tensor field.
100 // Allowable values: constant reference to the volTensorField.
101 //
102 // \return Divergence of iTF (vector field) which was computed on the faces of mesh.
103 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::reduced::Div(const volTensorField& iTF)
104 {
105  tmp<surfaceVectorField> tdivITF(nf_ & fvc::snGrad(iTF));
106 
107  return tdivITF;
108 }
109 
110 
111 //
112 //END-OF-FILE
113 //
114 
115 
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)
addToRunTimeSelectionTable(fvscStencil, GaussVolPoint, components)
reduced(const IOobject &)
Construct from IOobject.
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
defineTypeNameAndDebug(fvscStencil, 0)