cellMDLimitedGrads.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 "cellMDLimitedGrad.H"
30 #include "gaussGrad.H"
31 #include "fvMesh.H"
32 #include "volMesh.H"
33 #include "surfaceMesh.H"
34 #include "volFields.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 makeFvGradScheme(cellMDLimitedGrad)
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<>
46 (
47  const volScalarField& vsf,
48  const word& name
49 ) const
50 {
51  const fvMesh& mesh = vsf.mesh();
52 
53  tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
54 
55  if (k_ < SMALL)
56  {
57  return tGrad;
58  }
59 
60  volVectorField& g = tGrad.ref();
61 
62  const labelUList& owner = mesh.owner();
63  const labelUList& neighbour = mesh.neighbour();
64 
65  const volVectorField& C = mesh.C();
66  const surfaceVectorField& Cf = mesh.Cf();
67 
68  scalarField maxVsf(vsf.primitiveField());
69  scalarField minVsf(vsf.primitiveField());
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  maxVsf[own] = max(maxVsf[own], vsfNei);
80  minVsf[own] = min(minVsf[own], vsfNei);
81 
82  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
83  minVsf[nei] = min(minVsf[nei], vsfOwn);
84  }
85 
86 
87  const volScalarField::Boundary& bsf = vsf.boundaryField();
88 
89  forAll(bsf, patchi)
90  {
91  const fvPatchScalarField& psf = bsf[patchi];
92 
93  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
94 
95  if (psf.coupled())
96  {
97  const scalarField psfNei(psf.patchNeighbourField());
98 
99  forAll(pOwner, pFacei)
100  {
101  const label own = pOwner[pFacei];
102  const scalar vsfNei = psfNei[pFacei];
103 
104  maxVsf[own] = max(maxVsf[own], vsfNei);
105  minVsf[own] = min(minVsf[own], vsfNei);
106  }
107  }
108  else
109  {
110  forAll(pOwner, pFacei)
111  {
112  const label own = pOwner[pFacei];
113  const scalar vsfNei = psf[pFacei];
114 
115  maxVsf[own] = max(maxVsf[own], vsfNei);
116  minVsf[own] = min(minVsf[own], vsfNei);
117  }
118  }
119  }
120 
121  maxVsf -= vsf;
122  minVsf -= vsf;
123 
124  if (k_ < 1.0)
125  {
126  const scalarField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
127  maxVsf += maxMinVsf;
128  minVsf -= maxMinVsf;
129 
130  //maxVsf *= 1.0/k_;
131  //minVsf *= 1.0/k_;
132  }
133 
134 
135  forAll(owner, facei)
136  {
137  const label own = owner[facei];
138  const label nei = neighbour[facei];
139 
140  // owner side
141  limitFace
142  (
143  g[own],
144  maxVsf[own],
145  minVsf[own],
146  Cf[facei] - C[own]
147  );
148 
149  // neighbour side
150  limitFace
151  (
152  g[nei],
153  maxVsf[nei],
154  minVsf[nei],
155  Cf[facei] - C[nei]
156  );
157  }
158 
159 
160  forAll(bsf, patchi)
161  {
162  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
163  const vectorField& pCf = Cf.boundaryField()[patchi];
164 
165  forAll(pOwner, pFacei)
166  {
167  const label own = pOwner[pFacei];
168 
169  limitFace
170  (
171  g[own],
172  maxVsf[own],
173  minVsf[own],
174  pCf[pFacei] - C[own]
175  );
176  }
177  }
178 
179  g.correctBoundaryConditions();
181 
182  return tGrad;
183 }
184 
185 
186 template<>
189 (
190  const volVectorField& vsf,
191  const word& name
192 ) const
193 {
194  const fvMesh& mesh = vsf.mesh();
195 
196  tmp<volTensorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
197 
198  if (k_ < SMALL)
199  {
200  return tGrad;
201  }
202 
203  volTensorField& g = tGrad.ref();
204 
205  const labelUList& owner = mesh.owner();
206  const labelUList& neighbour = mesh.neighbour();
207 
208  const volVectorField& C = mesh.C();
209  const surfaceVectorField& Cf = mesh.Cf();
210 
211  vectorField maxVsf(vsf.primitiveField());
212  vectorField minVsf(vsf.primitiveField());
213 
214  forAll(owner, facei)
215  {
216  const label own = owner[facei];
217  const label nei = neighbour[facei];
218 
219  const vector& vsfOwn = vsf[own];
220  const vector& vsfNei = vsf[nei];
221 
222  maxVsf[own] = max(maxVsf[own], vsfNei);
223  minVsf[own] = min(minVsf[own], vsfNei);
224 
225  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
226  minVsf[nei] = min(minVsf[nei], vsfOwn);
227  }
228 
229 
230  const volVectorField::Boundary& bsf = vsf.boundaryField();
231 
232  forAll(bsf, patchi)
233  {
234  const fvPatchVectorField& psf = bsf[patchi];
235  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
236 
237  if (psf.coupled())
238  {
239  const vectorField psfNei(psf.patchNeighbourField());
240 
241  forAll(pOwner, pFacei)
242  {
243  const label own = pOwner[pFacei];
244  const vector& vsfNei = psfNei[pFacei];
245 
246  maxVsf[own] = max(maxVsf[own], vsfNei);
247  minVsf[own] = min(minVsf[own], vsfNei);
248  }
249  }
250  else
251  {
252  forAll(pOwner, pFacei)
253  {
254  const label own = pOwner[pFacei];
255  const vector& vsfNei = psf[pFacei];
256 
257  maxVsf[own] = max(maxVsf[own], vsfNei);
258  minVsf[own] = min(minVsf[own], vsfNei);
259  }
260  }
261  }
262 
263  maxVsf -= vsf;
264  minVsf -= vsf;
265 
266  if (k_ < 1.0)
267  {
268  const vectorField maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
269  maxVsf += maxMinVsf;
270  minVsf -= maxMinVsf;
271 
272  //maxVsf *= 1.0/k_;
273  //minVsf *= 1.0/k_;
274  }
275 
276 
277  forAll(owner, facei)
278  {
279  const label own = owner[facei];
280  const label nei = neighbour[facei];
281 
282  // owner side
283  limitFace
284  (
285  g[own],
286  maxVsf[own],
287  minVsf[own],
288  Cf[facei] - C[own]
289  );
290 
291  // neighbour side
292  limitFace
293  (
294  g[nei],
295  maxVsf[nei],
296  minVsf[nei],
297  Cf[facei] - C[nei]
298  );
299  }
300 
301 
302  forAll(bsf, patchi)
303  {
304  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
305  const vectorField& pCf = Cf.boundaryField()[patchi];
306 
307  forAll(pOwner, pFacei)
308  {
309  const label own = pOwner[pFacei];
310 
311  limitFace
312  (
313  g[own],
314  maxVsf[own],
315  minVsf[own],
316  pCf[pFacei] - C[own]
317  );
318  }
319  }
320 
321  g.correctBoundaryConditions();
323 
324  return tGrad;
325 }
326 
327 
328 // ************************************************************************* //
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
makeFvGradScheme
makeFvGradScheme(cellMDLimitedGrad)
Definition: cellMDLimitedGrads.C:39
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::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
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
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
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::fv::cellMDLimitedGrad::calcGrad
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Foam::Vector< scalar >
fixedValueFvPatchFields.H
Foam::UList< label >
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