faceMDLimitedGrads.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) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
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 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceMDLimitedGrad.H"
30 #include "cellMDLimitedGrad.H"
31 #include "gaussGrad.H"
32 #include "fvMesh.H"
33 #include "volMesh.H"
34 #include "surfaceMesh.H"
35 #include "volFields.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 makeFvGradScheme(faceMDLimitedGrad)
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<>
47 (
48  const volScalarField& vsf,
49  const word& name
50 ) const
51 {
52  const fvMesh& mesh = vsf.mesh();
53 
54  tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
55 
56  if (k_ < SMALL)
57  {
58  return tGrad;
59  }
60 
61  volVectorField& g = tGrad.ref();
62 
63  const labelUList& owner = mesh.owner();
64  const labelUList& neighbour = mesh.neighbour();
65 
66  const volVectorField& C = mesh.C();
67  const surfaceVectorField& Cf = mesh.Cf();
68 
69  const scalar rk = (1.0/k_ - 1.0);
70 
71  forAll(owner, facei)
72  {
73  const label own = owner[facei];
74  const label nei = neighbour[facei];
75 
76  const scalar vsfOwn = vsf[own];
77  const scalar vsfNei = vsf[nei];
78 
79  scalar maxFace = max(vsfOwn, vsfNei);
80  scalar minFace = min(vsfOwn, vsfNei);
81 
82  if (k_ < 1.0)
83  {
84  const scalar maxMinFace = rk*(maxFace - minFace);
85  maxFace += maxMinFace;
86  minFace -= maxMinFace;
87  }
88 
89  // owner side
90  cellMDLimitedGrad<scalar>::limitFace
91  (
92  g[own],
93  maxFace - vsfOwn,
94  minFace - vsfOwn,
95  Cf[facei] - C[own]
96  );
97 
98  // neighbour side
99  cellMDLimitedGrad<scalar>::limitFace
100  (
101  g[nei],
102  maxFace - vsfNei,
103  minFace - vsfNei,
104  Cf[facei] - C[nei]
105  );
106  }
107 
108  const volScalarField::Boundary& bsf = vsf.boundaryField();
109 
110  forAll(bsf, patchi)
111  {
112  const fvPatchScalarField& psf = bsf[patchi];
113 
114  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
115  const vectorField& pCf = Cf.boundaryField()[patchi];
116 
117  if (psf.coupled())
118  {
119  const scalarField psfNei(psf.patchNeighbourField());
120 
121  forAll(pOwner, pFacei)
122  {
123  const label own = pOwner[pFacei];
124 
125  scalar vsfOwn = vsf[own];
126  scalar vsfNei = psfNei[pFacei];
127 
128  scalar maxFace = max(vsfOwn, vsfNei);
129  scalar minFace = min(vsfOwn, vsfNei);
130 
131  if (k_ < 1.0)
132  {
133  const scalar maxMinFace = rk*(maxFace - minFace);
134  maxFace += maxMinFace;
135  minFace -= maxMinFace;
136  }
137 
138  cellMDLimitedGrad<scalar>::limitFace
139  (
140  g[own],
141  maxFace - vsfOwn,
142  minFace - vsfOwn,
143  pCf[pFacei] - C[own]
144  );
145  }
146  }
147  else if (psf.fixesValue())
148  {
149  forAll(pOwner, pFacei)
150  {
151  const label own = pOwner[pFacei];
152 
153  const scalar vsfOwn = vsf[own];
154  const scalar vsfNei = psf[pFacei];
155 
156  scalar maxFace = max(vsfOwn, vsfNei);
157  scalar minFace = min(vsfOwn, vsfNei);
158 
159  if (k_ < 1.0)
160  {
161  const scalar maxMinFace = rk*(maxFace - minFace);
162  maxFace += maxMinFace;
163  minFace -= maxMinFace;
164  }
165 
166  cellMDLimitedGrad<scalar>::limitFace
167  (
168  g[own],
169  maxFace - vsfOwn,
170  minFace - vsfOwn,
171  pCf[pFacei] - C[own]
172  );
173  }
174  }
175  }
176 
177  g.correctBoundaryConditions();
179 
180  return tGrad;
181 }
182 
183 
184 template<>
187 (
188  const volVectorField& vvf,
189  const word& name
190 ) const
191 {
192  const fvMesh& mesh = vvf.mesh();
193 
194  tmp<volTensorField> tGrad = basicGradScheme_().calcGrad(vvf, name);
195 
196  if (k_ < SMALL)
197  {
198  return tGrad;
199  }
200 
201  volTensorField& g = tGrad.ref();
202 
203  const labelUList& owner = mesh.owner();
204  const labelUList& neighbour = mesh.neighbour();
205 
206  const volVectorField& C = mesh.C();
207  const surfaceVectorField& Cf = mesh.Cf();
208 
209  const scalar rk = (1.0/k_ - 1.0);
210 
211  forAll(owner, facei)
212  {
213  const label own = owner[facei];
214  const label nei = neighbour[facei];
215 
216  const vector& vvfOwn = vvf[own];
217  const vector& vvfNei = vvf[nei];
218 
219  vector maxFace(max(vvfOwn, vvfNei));
220  vector minFace(min(vvfOwn, vvfNei));
221 
222  if (k_ < 1.0)
223  {
224  const vector maxMinFace(rk*(maxFace - minFace));
225  maxFace += maxMinFace;
226  minFace -= maxMinFace;
227  }
228 
229  // owner side
231  (
232  g[own],
233  maxFace - vvfOwn,
234  minFace - vvfOwn,
235  Cf[facei] - C[own]
236  );
237 
238 
239  // neighbour side
241  (
242  g[nei],
243  maxFace - vvfNei,
244  minFace - vvfNei,
245  Cf[facei] - C[nei]
246  );
247  }
248 
249 
250  const volVectorField::Boundary& bvf = vvf.boundaryField();
251 
252  forAll(bvf, patchi)
253  {
254  const fvPatchVectorField& psf = bvf[patchi];
255 
256  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
257  const vectorField& pCf = Cf.boundaryField()[patchi];
258 
259  if (psf.coupled())
260  {
261  const vectorField psfNei(psf.patchNeighbourField());
262 
263  forAll(pOwner, pFacei)
264  {
265  const label own = pOwner[pFacei];
266 
267  const vector& vvfOwn = vvf[own];
268  const vector& vvfNei = psfNei[pFacei];
269 
270  vector maxFace(max(vvfOwn, vvfNei));
271  vector minFace(min(vvfOwn, vvfNei));
272 
273  if (k_ < 1.0)
274  {
275  const vector maxMinFace(rk*(maxFace - minFace));
276  maxFace += maxMinFace;
277  minFace -= maxMinFace;
278  }
279 
281  (
282  g[own],
283  maxFace - vvfOwn, minFace - vvfOwn,
284  pCf[pFacei] - C[own]
285  );
286  }
287  }
288  else if (psf.fixesValue())
289  {
290  forAll(pOwner, pFacei)
291  {
292  const label own = pOwner[pFacei];
293 
294  const vector& vvfOwn = vvf[own];
295  const vector& vvfNei = psf[pFacei];
296 
297  vector maxFace(max(vvfOwn, vvfNei));
298  vector minFace(min(vvfOwn, vvfNei));
299 
300  if (k_ < 1.0)
301  {
302  const vector maxMinFace(rk*(maxFace - minFace));
303  maxFace += maxMinFace;
304  minFace -= maxMinFace;
305  }
306 
308  (
309  g[own],
310  maxFace - vvfOwn,
311  minFace - vvfOwn,
312  pCf[pFacei] - C[own]
313  );
314  }
315  }
316  }
317 
318  g.correctBoundaryConditions();
320 
321  return tGrad;
322 }
323 
324 
325 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
volMesh.H
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:332
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
makeFvGradScheme
makeFvGradScheme(faceMDLimitedGrad)
Definition: faceMDLimitedGrads.C:40
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:345
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
faceMDLimitedGrad.H
Foam::Field< vector >
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
surfaceMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:448
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::Vector< scalar >
Foam::fv::faceMDLimitedGrad::calcGrad
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
fixedValueFvPatchFields.H
Foam::UList< label >
Foam::fv::cellMDLimitedGrad
cellMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
Definition: cellMDLimitedGrad.H:65
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:59
Foam::GeometricField< vector, fvPatchField, volMesh >
cellMDLimitedGrad.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62