cellLimitedGrad.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) 2018 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 "cellLimitedGrad.H"
30 #include "gaussGrad.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class Type, class Limiter>
36 (
37  const Field<scalar>& limiter,
38  Field<vector>& gIf
39 ) const
40 {
41  gIf *= limiter;
42 }
43 
44 
45 template<class Type, class Limiter>
47 (
48  const Field<vector>& limiter,
49  Field<tensor>& gIf
50 ) const
51 {
52  forAll(gIf, celli)
53  {
54  gIf[celli] = tensor
55  (
56  cmptMultiply(limiter[celli], gIf[celli].x()),
57  cmptMultiply(limiter[celli], gIf[celli].y()),
58  cmptMultiply(limiter[celli], gIf[celli].z())
59  );
60  }
61 }
62 
63 
64 template<class Type, class Limiter>
66 <
68  <
72  >
73 >
75 (
76  const GeometricField<Type, fvPatchField, volMesh>& vsf,
77  const word& name
78 ) const
79 {
80  const fvMesh& mesh = vsf.mesh();
81 
82  tmp
83  <
84  GeometricField
85  <typename outerProduct<vector, Type>::type, fvPatchField, volMesh>
86  > tGrad = basicGradScheme_().calcGrad(vsf, name);
87 
88  if (k_ < SMALL)
89  {
90  return tGrad;
91  }
92 
93  GeometricField
94  <
96  fvPatchField,
97  volMesh
98  >& g = tGrad.ref();
99 
100  const labelUList& owner = mesh.owner();
101  const labelUList& neighbour = mesh.neighbour();
102 
103  const volVectorField& C = mesh.C();
104  const surfaceVectorField& Cf = mesh.Cf();
105 
106  Field<Type> maxVsf(vsf.primitiveField());
107  Field<Type> minVsf(vsf.primitiveField());
108 
109  forAll(owner, facei)
110  {
111  const label own = owner[facei];
112  const label nei = neighbour[facei];
113 
114  const Type& vsfOwn = vsf[own];
115  const Type& vsfNei = vsf[nei];
116 
117  maxVsf[own] = max(maxVsf[own], vsfNei);
118  minVsf[own] = min(minVsf[own], vsfNei);
119 
120  maxVsf[nei] = max(maxVsf[nei], vsfOwn);
121  minVsf[nei] = min(minVsf[nei], vsfOwn);
122  }
123 
124 
125  const typename GeometricField<Type, fvPatchField, volMesh>::Boundary& bsf =
126  vsf.boundaryField();
127 
128  forAll(bsf, patchi)
129  {
130  const fvPatchField<Type>& psf = bsf[patchi];
131  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
132 
133  if (psf.coupled())
134  {
135  const Field<Type> psfNei(psf.patchNeighbourField());
136 
137  forAll(pOwner, pFacei)
138  {
139  const label own = pOwner[pFacei];
140  const Type& vsfNei = psfNei[pFacei];
141 
142  maxVsf[own] = max(maxVsf[own], vsfNei);
143  minVsf[own] = min(minVsf[own], vsfNei);
144  }
145  }
146  else
147  {
148  forAll(pOwner, pFacei)
149  {
150  const label own = pOwner[pFacei];
151  const Type& vsfNei = psf[pFacei];
152 
153  maxVsf[own] = max(maxVsf[own], vsfNei);
154  minVsf[own] = min(minVsf[own], vsfNei);
155  }
156  }
157  }
158 
159  maxVsf -= vsf;
160  minVsf -= vsf;
161 
162  if (k_ < 1.0)
163  {
164  const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
165  maxVsf += maxMinVsf;
166  minVsf -= maxMinVsf;
167  }
168 
169 
170  // Create limiter initialized to 1
171  // Note: the limiter is not permitted to be > 1
172  Field<Type> limiter(vsf.primitiveField().size(), pTraits<Type>::one);
173 
174  forAll(owner, facei)
175  {
176  const label own = owner[facei];
177  const label nei = neighbour[facei];
178 
179  // owner side
180  limitFace
181  (
182  limiter[own],
183  maxVsf[own],
184  minVsf[own],
185  (Cf[facei] - C[own]) & g[own]
186  );
187 
188  // neighbour side
189  limitFace
190  (
191  limiter[nei],
192  maxVsf[nei],
193  minVsf[nei],
194  (Cf[facei] - C[nei]) & g[nei]
195  );
196  }
197 
198  forAll(bsf, patchi)
199  {
200  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
201  const vectorField& pCf = Cf.boundaryField()[patchi];
202 
203  forAll(pOwner, pFacei)
204  {
205  const label own = pOwner[pFacei];
206 
207  limitFace
208  (
209  limiter[own],
210  maxVsf[own],
211  minVsf[own],
212  ((pCf[pFacei] - C[own]) & g[own])
213  );
214  }
215  }
216 
217  if (fv::debug)
218  {
219  Info<< "gradient limiter for: " << vsf.name()
220  << " max = " << gMax(limiter)
221  << " min = " << gMin(limiter)
222  << " average: " << gAverage(limiter) << endl;
223  }
224 
225  limitGradient(limiter, g);
226  g.correctBoundaryConditions();
228 
229  return tGrad;
230 }
231 
232 
233 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
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
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:51
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:114
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
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()
gaussGrad.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::fv::cellLimitedGrad::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.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
cellLimitedGrad.H
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:59
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::fv::cellLimitedGrad
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
Definition: cellLimitedGrad.H:63
y
scalar y
Definition: LISASMDCalcMethod1.H:14