All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
GaussVolPointBase1D.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::GaussVolPoint1D
28 Description
29  This is a method for approximating derivatives of tangents to a face (1D case).
30  They are further used in the calculation of the QGD terms.
31 \*---------------------------------------------------------------------------*/
32 #include "GaussVolPointBase1D.H"
33 #include "fvMesh.H"
34 #include "volFields.H"
35 #include "surfaceFields.H"
36 #include "fvcSnGrad.H"
37 
39 :
40  nfRef_(mesh.thisDb().lookupObject<surfaceVectorField>("nf"))
41 {
42 };
43 
44 
46 {
47 }
48 
49 void Foam::fvsc::GaussVolPointBase1D::faceGrad(const volScalarField& f, surfaceVectorField& gradf)
50 {
51  if (f.mesh().nGeometricD() == 1)
52  {
53  gradf = nfRef_() * fvc::snGrad(f);
54  }
55 };
56 
57 void Foam::fvsc::GaussVolPointBase1D::faceGrad(const volVectorField& f, surfaceTensorField& gradf)
58 {
59  if (f.mesh().nGeometricD() == 1)
60  {
61  gradf = nfRef_() * fvc::snGrad(f);
62  }
63 };
64 
65 void Foam::fvsc::GaussVolPointBase1D::faceDiv(const volVectorField& f, surfaceScalarField& divf)
66 {
67  if (f.mesh().nGeometricD() == 1)
68  {
69  divf = nfRef_() & fvc::snGrad(f);
70  }
71 };
72 
73 void Foam::fvsc::GaussVolPointBase1D::faceDiv(const volTensorField& f, surfaceVectorField& divf)
74 {
75  if (f.mesh().nGeometricD() == 1)
76  {
77  divf = nfRef_() & fvc::snGrad(f);
78  }
79 }
80 
81 //
82 //END-OF-FILE
83 //
84 
85 
GaussVolPointBase1D(const fvMesh &mesh)
void faceGrad(const volScalarField &f, surfaceVectorField &gradf)
void faceDiv(const volVectorField &f, surfaceScalarField &divf)