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 -------------------------------------------------------------------------------
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 "faceMDLimitedGrad.H"
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(faceMDLimitedGrad)
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  scalar rk = (1.0/k_ - 1.0);
69 
70  forAll(owner, facei)
71  {
72  label own = owner[facei];
73  label nei = neighbour[facei];
74 
75  scalar vsfOwn = vsf[own];
76  scalar vsfNei = vsf[nei];
77 
78  scalar maxFace = max(vsfOwn, vsfNei);
79  scalar minFace = min(vsfOwn, vsfNei);
80 
81  if (k_ < 1.0)
82  {
83  scalar maxMinFace = rk*(maxFace - minFace);
84  maxFace += maxMinFace;
85  minFace -= maxMinFace;
86  }
87 
88  // owner side
89  cellMDLimitedGrad<scalar>::limitFace
90  (
91  g[own],
92  maxFace - vsfOwn,
93  minFace - vsfOwn,
94  Cf[facei] - C[own]
95  );
96 
97  // neighbour side
98  cellMDLimitedGrad<scalar>::limitFace
99  (
100  g[nei],
101  maxFace - vsfNei,
102  minFace - vsfNei,
103  Cf[facei] - C[nei]
104  );
105  }
106 
107  const volScalarField::Boundary& bsf = vsf.boundaryField();
108 
109  forAll(bsf, patchi)
110  {
111  const fvPatchScalarField& psf = bsf[patchi];
112 
113  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
114  const vectorField& pCf = Cf.boundaryField()[patchi];
115 
116  if (psf.coupled())
117  {
118  const scalarField psfNei(psf.patchNeighbourField());
119 
120  forAll(pOwner, pFacei)
121  {
122  label own = pOwner[pFacei];
123 
124  scalar vsfOwn = vsf[own];
125  scalar vsfNei = psfNei[pFacei];
126 
127  scalar maxFace = max(vsfOwn, vsfNei);
128  scalar minFace = min(vsfOwn, vsfNei);
129 
130  if (k_ < 1.0)
131  {
132  scalar maxMinFace = rk*(maxFace - minFace);
133  maxFace += maxMinFace;
134  minFace -= maxMinFace;
135  }
136 
137  cellMDLimitedGrad<scalar>::limitFace
138  (
139  g[own],
140  maxFace - vsfOwn,
141  minFace - vsfOwn,
142  pCf[pFacei] - C[own]
143  );
144  }
145  }
146  else if (psf.fixesValue())
147  {
148  forAll(pOwner, pFacei)
149  {
150  label own = pOwner[pFacei];
151 
152  scalar vsfOwn = vsf[own];
153  scalar vsfNei = psf[pFacei];
154 
155  scalar maxFace = max(vsfOwn, vsfNei);
156  scalar minFace = min(vsfOwn, vsfNei);
157 
158  if (k_ < 1.0)
159  {
160  scalar maxMinFace = rk*(maxFace - minFace);
161  maxFace += maxMinFace;
162  minFace -= maxMinFace;
163  }
164 
165  cellMDLimitedGrad<scalar>::limitFace
166  (
167  g[own],
168  maxFace - vsfOwn,
169  minFace - vsfOwn,
170  pCf[pFacei] - C[own]
171  );
172  }
173  }
174  }
175 
176  g.correctBoundaryConditions();
178 
179  return tGrad;
180 }
181 
182 
183 template<>
186 (
187  const volVectorField& vvf,
188  const word& name
189 ) const
190 {
191  const fvMesh& mesh = vvf.mesh();
192 
193  tmp<volTensorField> tGrad = basicGradScheme_().calcGrad(vvf, name);
194 
195  if (k_ < SMALL)
196  {
197  return tGrad;
198  }
199 
200  volTensorField& g = tGrad.ref();
201 
202  const labelUList& owner = mesh.owner();
203  const labelUList& neighbour = mesh.neighbour();
204 
205  const volVectorField& C = mesh.C();
206  const surfaceVectorField& Cf = mesh.Cf();
207 
208  scalar rk = (1.0/k_ - 1.0);
209 
210  forAll(owner, facei)
211  {
212  label own = owner[facei];
213  label nei = neighbour[facei];
214 
215  vector vvfOwn = vvf[own];
216  vector vvfNei = vvf[nei];
217 
218  vector maxFace = max(vvfOwn, vvfNei);
219  vector minFace = min(vvfOwn, vvfNei);
220 
221  if (k_ < 1.0)
222  {
223  vector maxMinFace = rk*(maxFace - minFace);
224  maxFace += maxMinFace;
225  minFace -= maxMinFace;
226  }
227 
228  // owner side
230  (
231  g[own],
232  maxFace - vvfOwn,
233  minFace - vvfOwn,
234  Cf[facei] - C[own]
235  );
236 
237 
238  // neighbour side
240  (
241  g[nei],
242  maxFace - vvfNei,
243  minFace - vvfNei,
244  Cf[facei] - C[nei]
245  );
246  }
247 
248 
249  const volVectorField::Boundary& bvf = vvf.boundaryField();
250 
251  forAll(bvf, patchi)
252  {
253  const fvPatchVectorField& psf = bvf[patchi];
254 
255  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
256  const vectorField& pCf = Cf.boundaryField()[patchi];
257 
258  if (psf.coupled())
259  {
260  const vectorField psfNei(psf.patchNeighbourField());
261 
262  forAll(pOwner, pFacei)
263  {
264  label own = pOwner[pFacei];
265 
266  vector vvfOwn = vvf[own];
267  vector vvfNei = psfNei[pFacei];
268 
269  vector maxFace = max(vvfOwn, vvfNei);
270  vector minFace = min(vvfOwn, vvfNei);
271 
272  if (k_ < 1.0)
273  {
274  vector maxMinFace = rk*(maxFace - minFace);
275  maxFace += maxMinFace;
276  minFace -= maxMinFace;
277  }
278 
280  (
281  g[own],
282  maxFace - vvfOwn, minFace - vvfOwn,
283  pCf[pFacei] - C[own]
284  );
285  }
286  }
287  else if (psf.fixesValue())
288  {
289  forAll(pOwner, pFacei)
290  {
291  label own = pOwner[pFacei];
292 
293  vector vvfOwn = vvf[own];
294  vector vvfNei = psf[pFacei];
295 
296  vector maxFace = max(vvfOwn, vvfNei);
297  vector minFace = min(vvfOwn, vvfNei);
298 
299  if (k_ < 1.0)
300  {
301  vector maxMinFace = rk*(maxFace - minFace);
302  maxFace += maxMinFace;
303  minFace -= maxMinFace;
304  }
305 
307  (
308  g[own],
309  maxFace - vvfOwn,
310  minFace - vvfOwn,
311  pCf[pFacei] - C[own]
312  );
313  }
314  }
315  }
316 
317  g.correctBoundaryConditions();
319 
320  return tGrad;
321 }
322 
323 
324 // ************************************************************************* //
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
volMesh.H
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
makeFvGradScheme
makeFvGradScheme(faceMDLimitedGrad)
Definition: faceMDLimitedGrads.C:39
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:331
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
faceMDLimitedGrad.H
Foam::Field< vector >
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 >
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
Return the gradient of the given field to the gradScheme::grad.
fixedValueFvPatchFields.H
Foam::UList< label >
Foam::fv::cellMDLimitedGrad
cellMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
Definition: cellMDLimitedGrad.H:64
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 >
cellMDLimitedGrad.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62