All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
GaussVolPointBase.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::GaussVolPointBase
28 Description
29  This is a method for approximating derivatives of tangents to a face.
30  They are further used in the calculation of the QGD terms.
31 \*---------------------------------------------------------------------------*/
32 #include "GaussVolPointBase.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "coupledFvPatch.H"
36 #include "processorFvPatch.H"
37 #include "emptyFvPatch.H"
38 #include "wedgeFvPatch.H"
39 #include "fvcSnGrad.H"
40 
42 :
44  GaussVolPointBase2D(mesh),
46 {
47 
48 };
49 
51 {
52 };
53 
54 void Foam::fvsc::GaussVolPointBase::faceGrad(const volScalarField& vF, surfaceVectorField& gradf)
55 {
56  if (vF.mesh().nGeometricD() == 1)
57  {
59  return;
60  }
61  else if (vF.mesh().nGeometricD() == 2)
62  {
64  return;
65  }
66  else
67  {
69  }
70 }
71 
72 void Foam::fvsc::GaussVolPointBase::faceGrad(const volVectorField& vF, surfaceTensorField& gradf)
73 {
74  if (vF.mesh().nGeometricD() == 1)
75  {
77  return;
78  }
79  else if (vF.mesh().nGeometricD() == 2)
80  {
81  surfaceVectorField gradU (vector::zero*fvc::snGrad(vF.component(0)));
82  surfaceVectorField gradV (gradU*0.0);
83  surfaceVectorField gradW (gradU*0.0);
84 
85  GaussVolPointBase2D::faceGrad(vF.component(0), gradU);
86  GaussVolPointBase2D::faceGrad(vF.component(1), gradV);
87  GaussVolPointBase2D::faceGrad(vF.component(2), gradW);
88 
89  //Internal field
90  gradf.primitiveFieldRef().replace(0, gradU.primitiveField().component(0));
91  gradf.primitiveFieldRef().replace(1, gradV.primitiveField().component(0));
92  gradf.primitiveFieldRef().replace(2, gradW.primitiveField().component(0));
93 
94  gradf.primitiveFieldRef().replace(3, gradU.primitiveField().component(1));
95  gradf.primitiveFieldRef().replace(4, gradV.primitiveField().component(1));
96  gradf.primitiveFieldRef().replace(5, gradW.primitiveField().component(1));
97 
98  gradf.primitiveFieldRef().replace(6, gradU.primitiveField().component(2));
99  gradf.primitiveFieldRef().replace(7, gradV.primitiveField().component(2));
100  gradf.primitiveFieldRef().replace(8, gradW.primitiveField().component(2));
101 
102  //Boundary field
103  gradf.boundaryFieldRef().replace(0, gradU.boundaryField().component(0));
104  gradf.boundaryFieldRef().replace(1, gradV.boundaryField().component(0));
105  gradf.boundaryFieldRef().replace(2, gradW.boundaryField().component(0));
106 
107  gradf.boundaryFieldRef().replace(3, gradU.boundaryField().component(1));
108  gradf.boundaryFieldRef().replace(4, gradV.boundaryField().component(1));
109  gradf.boundaryFieldRef().replace(5, gradW.boundaryField().component(1));
110 
111  gradf.boundaryFieldRef().replace(6, gradU.boundaryField().component(2));
112  gradf.boundaryFieldRef().replace(7, gradV.boundaryField().component(2));
113  gradf.boundaryFieldRef().replace(8, gradW.boundaryField().component(2));
114 
115  return;
116  }
117  else
118  {
120  }
121 }
122 
123 void Foam::fvsc::GaussVolPointBase::faceDiv(const volVectorField& vVF, surfaceScalarField& divf)
124 {
125 
126  if (vVF.mesh().nGeometricD() == 1)
127  {
129  return;
130  }
131  else if (vVF.mesh().nGeometricD() == 2)
132  {
134  return;
135  }
136  else
137  {
139  }
140 }
141 
142 void Foam::fvsc::GaussVolPointBase::faceDiv(const volTensorField& vTF, surfaceVectorField& divf)
143 {
144  if (vTF.mesh().nGeometricD() == 1)
145  {
147  return;
148  }
149  else if (vTF.mesh().nGeometricD() == 2)
150  {
152  return;
153  }
154  else
155  {
157  }
158 }
159 
160 //
161 //END-OF-FILE
162 //
163 
164 
void faceDiv(const volVectorField &f, surfaceScalarField &divf)
void faceGrad(const volScalarField &f, surfaceVectorField &gradf)
void faceGrad(const volScalarField &f, surfaceVectorField &gradf)
GaussVolPointBase(const fvMesh &mesh)
void faceDiv(const volVectorField &vVF, surfaceScalarField &divf)
void faceDiv(const volVectorField &f, surfaceScalarField &divf)
void faceGrad(const volScalarField &f, surfaceVectorField &gradf)
void faceDiv(const volVectorField &f, surfaceScalarField &divf)
void faceGrad(const volScalarField &vF, surfaceVectorField &gradf)