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