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>
36(
37 const Field<scalar>& limiter,
38 Field<vector>& gIf
39) const
40{
41 gIf *= limiter;
42}
43
44
45template<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
64template<class Type, class Limiter>
66<
68 <
72 >
73>
75(
77 const word& name
78) const
79{
80 const fvMesh& mesh = vsf.mesh();
81
82 tmp
83 <
86 > tGrad = basicGradScheme_().calcGrad(vsf, name);
87
88 if (k_ < SMALL)
89 {
90 return tGrad;
91 }
92
94 <
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 auto& bsf = vsf.boundaryField();
126
127 forAll(bsf, patchi)
128 {
129 const fvPatchField<Type>& psf = bsf[patchi];
130 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
131
132 if (psf.coupled())
133 {
134 const Field<Type> psfNei(psf.patchNeighbourField());
135
136 forAll(pOwner, pFacei)
137 {
138 const label own = pOwner[pFacei];
139 const Type& vsfNei = psfNei[pFacei];
140
141 maxVsf[own] = max(maxVsf[own], vsfNei);
142 minVsf[own] = min(minVsf[own], vsfNei);
143 }
144 }
145 else
146 {
147 forAll(pOwner, pFacei)
148 {
149 const label own = pOwner[pFacei];
150 const Type& vsfNei = psf[pFacei];
151
152 maxVsf[own] = max(maxVsf[own], vsfNei);
153 minVsf[own] = min(minVsf[own], vsfNei);
154 }
155 }
156 }
157
158 maxVsf -= vsf;
159 minVsf -= vsf;
160
161 if (k_ < 1.0)
162 {
163 const Field<Type> maxMinVsf((1.0/k_ - 1.0)*(maxVsf - minVsf));
164 maxVsf += maxMinVsf;
165 minVsf -= maxMinVsf;
166 }
167
168
169 // Create limiter initialized to 1
170 // Note: the limiter is not permitted to be > 1
172
173 forAll(owner, facei)
174 {
175 const label own = owner[facei];
176 const label nei = neighbour[facei];
177
178 // owner side
179 limitFace
180 (
181 limiter[own],
182 maxVsf[own],
183 minVsf[own],
184 (Cf[facei] - C[own]) & g[own]
185 );
186
187 // neighbour side
188 limitFace
189 (
190 limiter[nei],
191 maxVsf[nei],
192 minVsf[nei],
193 (Cf[facei] - C[nei]) & g[nei]
194 );
195 }
196
197 forAll(bsf, patchi)
198 {
199 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
200 const vectorField& pCf = Cf.boundaryField()[patchi];
201
202 forAll(pOwner, pFacei)
203 {
204 const label own = pOwner[pFacei];
205
206 limitFace
207 (
208 limiter[own],
209 maxVsf[own],
210 minVsf[own],
211 ((pCf[pFacei] - C[own]) & g[own])
212 );
213 }
214 }
215
216 if (fv::debug)
217 {
218 Info<< "gradient limiter for: " << vsf.name()
219 << " max = " << gMax(limiter)
220 << " min = " << gMin(limiter)
221 << " average: " << gAverage(limiter) << endl;
222 }
223
224 limitGradient(limiter, g);
225 g.correctBoundaryConditions();
227
228 return tGrad;
229}
230
231
232// ************************************************************************* //
scalar y
Y[inertIndex] max(0.0)
const uniformDimensionedVectorField & g
Graphite solid properties.
Definition: C.H:53
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:350
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
Basic second-order gradient scheme using face-interpolation and Gauss' theorem.
Definition: gaussGrad.H:66
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
Definition: tmp.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:54
type
Volume classification types.
Definition: volumeType.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Type gAverage(const FieldField< Field, Type > &f)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Type gMax(const FieldField< Field, Type > &f)
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333