Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
34template<class Type, class Limiter>
35void Foam::fv::cellLimitedGrad<Type, Limiter>::limitGradient
36(
37 const Field<scalar>& limiter,
38 Field<vector>& gIf
39) const
40{
41 gIf *= limiter;
42}
43
44
45template<class Type, class Limiter>
46void Foam::fv::cellLimitedGrad<Type, Limiter>::limitGradient
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
64template<class Type, class Limiter>
66(
68 <
72 >& g,
74) const
75{
76 const fvMesh& mesh = vsf.mesh();
77
78 basicGradScheme_().calcGrad(g, vsf);
79
80 if (k_ < SMALL)
81 {
82 return;
83 }
84
85 const labelUList& owner = mesh.owner();
86 const labelUList& neighbour = mesh.neighbour();
87
88 const volVectorField& C = mesh.C();
89 const surfaceVectorField& Cf = mesh.Cf();
90
91 Field<Type> maxVsf(vsf.primitiveField());
92 Field<Type> minVsf(vsf.primitiveField());
93
94 forAll(owner, facei)
95 {
96 const label own = owner[facei];
97 const label nei = neighbour[facei];
98
99 const Type& vsfOwn = vsf[own];
100 const Type& vsfNei = vsf[nei];
101
102 maxVsf[own] = max(maxVsf[own], vsfNei);
103 minVsf[own] = min(minVsf[own], vsfNei);
104
105 maxVsf[nei] = max(maxVsf[nei], vsfOwn);
106 minVsf[nei] = min(minVsf[nei], vsfOwn);
107 }
108
109
110 const auto& bsf = vsf.boundaryField();
111
112 forAll(bsf, patchi)
113 {
114 const fvPatchField<Type>& psf = bsf[patchi];
115 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
116
117 if (psf.coupled())
118 {
119 const Field<Type> psfNei(psf.patchNeighbourField());
120
121 forAll(pOwner, pFacei)
122 {
123 const label own = pOwner[pFacei];
124 const Type& vsfNei = psfNei[pFacei];
125
126 maxVsf[own] = max(maxVsf[own], vsfNei);
127 minVsf[own] = min(minVsf[own], vsfNei);
128 }
129 }
130 else
131 {
132 forAll(pOwner, pFacei)
133 {
134 const label own = pOwner[pFacei];
135 const Type& vsfNei = psf[pFacei];
136
137 maxVsf[own] = max(maxVsf[own], vsfNei);
138 minVsf[own] = min(minVsf[own], vsfNei);
139 }
140 }
141 }
142
143 maxVsf -= vsf;
144 minVsf -= vsf;
145
146 if (k_ < 1.0)
147 {
148 const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
149 maxVsf += maxMinVsf;
150 minVsf -= maxMinVsf;
151 }
152
153
154 // Create limiter initialized to 1
155 // Note: the limiter is not permitted to be > 1
156 Field<Type> limiter(vsf.primitiveField().size(), pTraits<Type>::one);
157
158 forAll(owner, facei)
159 {
160 const label own = owner[facei];
161 const label nei = neighbour[facei];
162
163 // owner side
164 limitFace
165 (
166 limiter[own],
167 maxVsf[own],
168 minVsf[own],
169 (Cf[facei] - C[own]) & g[own]
170 );
171
172 // neighbour side
173 limitFace
174 (
175 limiter[nei],
176 maxVsf[nei],
177 minVsf[nei],
178 (Cf[facei] - C[nei]) & g[nei]
179 );
180 }
181
182 forAll(bsf, patchi)
183 {
184 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
185 const vectorField& pCf = Cf.boundaryField()[patchi];
186
187 forAll(pOwner, pFacei)
188 {
189 const label own = pOwner[pFacei];
190
191 limitFace
192 (
193 limiter[own],
194 maxVsf[own],
195 minVsf[own],
196 ((pCf[pFacei] - C[own]) & g[own])
197 );
198 }
199 }
200
201 if (fv::debug)
202 {
203 auto limits = gMinMax(limiter);
204 auto avg = gAverage(limiter);
205
206 Info<< "gradient limiter for: " << vsf.name()
207 << " min = " << limits.min()
208 << " max = " << limits.max()
209 << " average: " << avg << endl;
210 }
211
212 limitGradient(limiter, g);
213 g.correctBoundaryConditions();
214 gaussGrad<Type>::correctBoundaryConditions(vsf, g);
215}
216
217
218template<class Type, class Limiter>
220<
222 <
226 >
227>
229(
231 const word& name
232) const
233{
234 typedef typename outerProduct<vector, Type>::type GradType;
236
237 const fvMesh& mesh = vsf.mesh();
238
240 (
241 new GradFieldType
242 (
244 (
245 name,
246 vsf.instance(),
247 mesh,
250 ),
251 mesh,
252 vsf.dimensions()/dimLength,
254 )
255 );
256
257 calcGrad(tGrad.ref(), vsf);
258
259 return tGrad;
260}
261
262
263// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar y
const uniformDimensionedVectorField & g
Graphite solid properties.
Definition C.H:49
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual bool coupled() const
True if the patch field is coupled.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
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.
void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate) const
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
Definition gaussGrad.C:199
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition faNVDscheme.C:31
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299