faceLimitedGrads.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
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 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "faceLimitedGrad.H"
29 #include "gaussGrad.H"
30 #include "fvMesh.H"
31 #include "volMesh.H"
32 #include "surfaceMesh.H"
33 #include "volFields.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 makeFvGradScheme(faceLimitedGrad)
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<>
45 (
46  const volScalarField& vsf,
47  const word& name
48 ) const
49 {
50  const fvMesh& mesh = vsf.mesh();
51 
52  tmp<volVectorField> tGrad = basicGradScheme_().calcGrad(vsf, name);
53 
54  if (k_ < SMALL)
55  {
56  return tGrad;
57  }
58 
59  volVectorField& g = tGrad.ref();
60 
61  const labelUList& owner = mesh.owner();
62  const labelUList& neighbour = mesh.neighbour();
63 
64  const volVectorField& C = mesh.C();
65  const surfaceVectorField& Cf = mesh.Cf();
66 
67  // create limiter
68  scalarField limiter(vsf.primitiveField().size(), 1.0);
69 
70  scalar rk = (1.0/k_ - 1.0);
71 
72  forAll(owner, facei)
73  {
74  label own = owner[facei];
75  label nei = neighbour[facei];
76 
77  scalar vsfOwn = vsf[own];
78  scalar vsfNei = vsf[nei];
79 
80  scalar maxFace = max(vsfOwn, vsfNei);
81  scalar minFace = min(vsfOwn, vsfNei);
82  scalar maxMinFace = rk*(maxFace - minFace);
83  maxFace += maxMinFace;
84  minFace -= maxMinFace;
85 
86  // owner side
87  limitFace
88  (
89  limiter[own],
90  maxFace - vsfOwn, minFace - vsfOwn,
91  (Cf[facei] - C[own]) & g[own]
92  );
93 
94  // neighbour side
95  limitFace
96  (
97  limiter[nei],
98  maxFace - vsfNei, minFace - vsfNei,
99  (Cf[facei] - C[nei]) & g[nei]
100  );
101  }
102 
103  const volScalarField::Boundary& bsf = vsf.boundaryField();
104 
105  forAll(bsf, patchi)
106  {
107  const fvPatchScalarField& psf = bsf[patchi];
108 
109  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
110  const vectorField& pCf = Cf.boundaryField()[patchi];
111 
112  if (psf.coupled())
113  {
114  const scalarField psfNei(psf.patchNeighbourField());
115 
116  forAll(pOwner, pFacei)
117  {
118  label own = pOwner[pFacei];
119 
120  scalar vsfOwn = vsf[own];
121  scalar vsfNei = psfNei[pFacei];
122 
123  scalar maxFace = max(vsfOwn, vsfNei);
124  scalar minFace = min(vsfOwn, vsfNei);
125  scalar maxMinFace = rk*(maxFace - minFace);
126  maxFace += maxMinFace;
127  minFace -= maxMinFace;
128 
129  limitFace
130  (
131  limiter[own],
132  maxFace - vsfOwn, minFace - vsfOwn,
133  (pCf[pFacei] - C[own]) & g[own]
134  );
135  }
136  }
137  else if (psf.fixesValue())
138  {
139  forAll(pOwner, pFacei)
140  {
141  label own = pOwner[pFacei];
142 
143  scalar vsfOwn = vsf[own];
144  scalar vsfNei = psf[pFacei];
145 
146  scalar maxFace = max(vsfOwn, vsfNei);
147  scalar minFace = min(vsfOwn, vsfNei);
148  scalar maxMinFace = rk*(maxFace - minFace);
149  maxFace += maxMinFace;
150  minFace -= maxMinFace;
151 
152  limitFace
153  (
154  limiter[own],
155  maxFace - vsfOwn, minFace - vsfOwn,
156  (pCf[pFacei] - C[own]) & g[own]
157  );
158  }
159  }
160  }
161 
162  if (fv::debug)
163  {
164  Info<< "gradient limiter for: " << vsf.name()
165  << " max = " << gMax(limiter)
166  << " min = " << gMin(limiter)
167  << " average: " << gAverage(limiter) << endl;
168  }
169 
170  g.primitiveFieldRef() *= limiter;
171  g.correctBoundaryConditions();
173 
174  return tGrad;
175 }
176 
177 
178 template<>
181 (
182  const volVectorField& vvf,
183  const word& name
184 ) const
185 {
186  const fvMesh& mesh = vvf.mesh();
187 
188  tmp<volTensorField> tGrad = basicGradScheme_().calcGrad(vvf, name);
189 
190  if (k_ < SMALL)
191  {
192  return tGrad;
193  }
194 
195  volTensorField& g = tGrad.ref();
196 
197  const labelUList& owner = mesh.owner();
198  const labelUList& neighbour = mesh.neighbour();
199 
200  const volVectorField& C = mesh.C();
201  const surfaceVectorField& Cf = mesh.Cf();
202 
203  // create limiter
204  scalarField limiter(vvf.primitiveField().size(), 1.0);
205 
206  scalar rk = (1.0/k_ - 1.0);
207 
208  forAll(owner, facei)
209  {
210  label own = owner[facei];
211  label nei = neighbour[facei];
212 
213  vector vvfOwn = vvf[own];
214  vector vvfNei = vvf[nei];
215 
216  // owner side
217  vector gradf = (Cf[facei] - C[own]) & g[own];
218 
219  scalar vsfOwn = gradf & vvfOwn;
220  scalar vsfNei = gradf & vvfNei;
221 
222  scalar maxFace = max(vsfOwn, vsfNei);
223  scalar minFace = min(vsfOwn, vsfNei);
224  scalar maxMinFace = rk*(maxFace - minFace);
225  maxFace += maxMinFace;
226  minFace -= maxMinFace;
227 
228  limitFace
229  (
230  limiter[own],
231  maxFace - vsfOwn, minFace - vsfOwn,
232  magSqr(gradf)
233  );
234 
235 
236  // neighbour side
237  gradf = (Cf[facei] - C[nei]) & g[nei];
238 
239  vsfOwn = gradf & vvfOwn;
240  vsfNei = gradf & vvfNei;
241 
242  maxFace = max(vsfOwn, vsfNei);
243  minFace = min(vsfOwn, vsfNei);
244 
245  limitFace
246  (
247  limiter[nei],
248  maxFace - vsfNei, minFace - vsfNei,
249  magSqr(gradf)
250  );
251  }
252 
253 
254  const volVectorField::Boundary& bvf = vvf.boundaryField();
255 
256  forAll(bvf, patchi)
257  {
258  const fvPatchVectorField& psf = bvf[patchi];
259 
260  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
261  const vectorField& pCf = Cf.boundaryField()[patchi];
262 
263  if (psf.coupled())
264  {
265  const vectorField psfNei(psf.patchNeighbourField());
266 
267  forAll(pOwner, pFacei)
268  {
269  label own = pOwner[pFacei];
270 
271  vector vvfOwn = vvf[own];
272  vector vvfNei = psfNei[pFacei];
273 
274  vector gradf = (pCf[pFacei] - C[own]) & g[own];
275 
276  scalar vsfOwn = gradf & vvfOwn;
277  scalar vsfNei = gradf & vvfNei;
278 
279  scalar maxFace = max(vsfOwn, vsfNei);
280  scalar minFace = min(vsfOwn, vsfNei);
281  scalar maxMinFace = rk*(maxFace - minFace);
282  maxFace += maxMinFace;
283  minFace -= maxMinFace;
284 
285  limitFace
286  (
287  limiter[own],
288  maxFace - vsfOwn, minFace - vsfOwn,
289  magSqr(gradf)
290  );
291  }
292  }
293  else if (psf.fixesValue())
294  {
295  forAll(pOwner, pFacei)
296  {
297  label own = pOwner[pFacei];
298 
299  vector vvfOwn = vvf[own];
300  vector vvfNei = psf[pFacei];
301 
302  vector gradf = (pCf[pFacei] - C[own]) & g[own];
303 
304  scalar vsfOwn = gradf & vvfOwn;
305  scalar vsfNei = gradf & vvfNei;
306 
307  scalar maxFace = max(vsfOwn, vsfNei);
308  scalar minFace = min(vsfOwn, vsfNei);
309  scalar maxMinFace = rk*(maxFace - minFace);
310  maxFace += maxMinFace;
311  minFace -= maxMinFace;
312 
313  limitFace
314  (
315  limiter[own],
316  maxFace - vsfOwn, minFace - vsfOwn,
317  magSqr(gradf)
318  );
319  }
320  }
321  }
322 
323  if (fv::debug)
324  {
325  Info<< "gradient limiter for: " << vvf.name()
326  << " max = " << gMax(limiter)
327  << " min = " << gMin(limiter)
328  << " average: " << gAverage(limiter) << endl;
329  }
330 
331  g.primitiveFieldRef() *= limiter;
332  g.correctBoundaryConditions();
334 
335  return tGrad;
336 }
337 
338 
339 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
makeFvGradScheme
makeFvGradScheme(faceLimitedGrad)
Definition: faceLimitedGrads.C:38
faceLimitedGrad.H
volMesh.H
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:318
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
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:331
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:83
fvMesh.H
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:434
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::Vector< scalar >
fixedValueFvPatchFields.H
Foam::UList< label >
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
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 >
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::fv::faceLimitedGrad::calcGrad
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
Definition: faceLimitedGrad.H:128