Loading...
Searching...
No Matches
fusedGaussGrad.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 Copyright (C) 2018-2021 OpenCFD Ltd.
10 Copyright (C) 2024 M. Janssens
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "fusedGaussGrad.H"
32#include "fvcSurfaceOps.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36template<class Type>
38<
40 <
44 >
45>
47(
49 const word& name
50)
51{
52 typedef typename outerProduct<vector, Type>::type GradType;
54
55 const fvMesh& mesh = ssf.mesh();
56
57 DebugPout<< "fusedGaussGrad<Type>::gradf on " << ssf.name()
58 << " with name " << name << endl;
59
61 (
62 new GradFieldType
63 (
65 (
66 name,
67 ssf.instance(),
68 mesh,
71 ),
72 mesh,
75 )
76 );
77 GradFieldType& gGrad = tgGrad.ref();
78
79 const labelUList& owner = mesh.owner();
80 const labelUList& neighbour = mesh.neighbour();
81 const vectorField& Sf = mesh.Sf();
82
83 Field<GradType>& igGrad = gGrad;
84 const Field<Type>& issf = ssf;
85
86 forAll(owner, facei)
87 {
88 const GradType Sfssf = Sf[facei]*issf[facei];
89
90 igGrad[owner[facei]] += Sfssf;
91 igGrad[neighbour[facei]] -= Sfssf;
92 }
93
94 forAll(mesh.boundary(), patchi)
95 {
96 const labelUList& pFaceCells =
97 mesh.boundary()[patchi].faceCells();
98
99 const vectorField& pSf = mesh.Sf().boundaryField()[patchi];
100
101 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
102
103 forAll(mesh.boundary()[patchi], facei)
104 {
105 igGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
106 }
107 }
108
109 igGrad /= mesh.V();
110
111 gGrad.correctBoundaryConditions();
112
113 return tgGrad;
114}
115
116
117template<class Type>
119<
121 <
125 >
126>
128(
130 const word& name
131) const
132{
133 typedef typename outerProduct<vector, Type>::type GradType;
135
136 const fvMesh& mesh = vf.mesh();
137
138 DebugPout<< "fusedGaussGrad<Type>::calcGrad on " << vf.name()
139 << " with name " << name << endl;
140
141 tmp<GradFieldType> tgGrad
142 (
143 new GradFieldType
144 (
146 (
147 name,
148 vf.instance(),
149 mesh,
152 ),
153 mesh,
156 )
157 );
158 GradFieldType& gGrad = tgGrad.ref();
159
160 if (this->tinterpScheme_().corrected())
161 {
162 const auto tfaceCorr(this->tinterpScheme_().correction(vf));
163 auto& faceCorr = tfaceCorr();
164
165 DebugPout<< "fusedGaussGrad<Type>::calcGrad corrected interpScheme "
166 << this->tinterpScheme_().type() << endl;
167
168 const auto interpolate = [&]
169 (
170 const vector& area,
171 const scalar lambda,
172
173 const Type& ownVal,
174 const Type& neiVal,
175
176 const Type& correction
177
178 ) -> GradType
179 {
180 return area*((lambda*(ownVal - neiVal) + neiVal) + correction);
181 };
182
184 (
185 this->tinterpScheme_().weights(vf),
186 vf,
187 faceCorr,
189 gGrad,
190 false
191 );
192 }
193 else
194 {
195 DebugPout<< "fusedGaussGrad<Type>::calcGrad uncorrected interpScheme "
196 << this->tinterpScheme_().type() << endl;
197
198 const auto interpolate = [&]
199 (
200 const vector& area,
201 const scalar lambda,
202 const Type& ownVal,
203 const Type& neiVal
204 ) -> GradType
205 {
206 return area*(lambda*(ownVal - neiVal) + neiVal);
207 };
208
210 (
211 tinterpScheme_().weights(vf),
212 vf,
214 gGrad,
215 false
216 );
217 }
218
219 gGrad.primitiveFieldRef() /= mesh.V();
220
221 gGrad.correctBoundaryConditions();
222
224
225 return tgGrad;
226}
227
228
229template<class Type>
231(
233 <
236 volMesh
237 >& gGrad,
239) const
240{
241 typedef typename outerProduct<vector, Type>::type GradType;
242
243 const fvMesh& mesh = vf.mesh();
244
245 DebugPout<< "fusedGaussGrad<Type>::calcGrad on " << vf.name()
246 << " into " << gGrad.name() << endl;
247
249
250 if (this->tinterpScheme_().corrected())
251 {
252 const auto tfaceCorr(this->tinterpScheme_().correction(vf));
253 auto& faceCorr = tfaceCorr();
254
255 DebugPout<< "fusedGaussGrad<Type>::calcGrad corrected interpScheme "
256 << this->tinterpScheme_().type() << endl;
257
258 const auto interpolate = [&]
259 (
260 const vector& area,
261 const scalar lambda,
262
263 const Type& ownVal,
264 const Type& neiVal,
265
266 const Type& correction
267
268 ) -> GradType
269 {
270 return area*((lambda*(ownVal - neiVal) + neiVal) + correction);
271 };
272
274 (
275 this->tinterpScheme_().weights(vf),
276 vf,
277 faceCorr,
279 gGrad,
280 false
281 );
282 }
283 else
284 {
285 DebugPout<< "fusedGaussGrad<Type>::calcGrad uncorrected interpScheme "
286 << this->tinterpScheme_().type() << endl;
287
288 const auto interpolate = [&]
289 (
290 const vector& area,
291 const scalar lambda,
292 const Type& ownVal,
293 const Type& neiVal
294 ) -> GradType
295 {
296 return area*(lambda*(ownVal - neiVal) + neiVal);
297 };
298
299 fvc::surfaceSum
300 (
301 tinterpScheme_().weights(vf),
302 vf,
304 gGrad,
305 false
306 );
307 }
308
309 gGrad.primitiveFieldRef() /= mesh.V();
310
311 gGrad.correctBoundaryConditions();
313 correctBoundaryConditions(vf, gGrad);
314}
315
316
317template<class Type>
318template<class GradType>
320(
323)
324{
325 DebugPout<< "fusedGaussGrad<Type>::correctBoundaryConditions on "
326 << vf.name() << " with gGrad " << gGrad.name() << endl;
327
328 const fvMesh& mesh = vf.mesh();
329 auto& gGradbf = gGrad.boundaryFieldRef();
330
331 forAll(vf.boundaryField(), patchi)
332 {
333 if (!vf.boundaryField()[patchi].coupled())
334 {
335 const auto& pSf = mesh.Sf().boundaryField()[patchi];
336 const auto tsnGrad(vf.boundaryField()[patchi].snGrad());
337 const auto& snGrad = tsnGrad();
338 auto& pgrad = gGradbf[patchi];
339
340 forAll(pgrad, facei)
341 {
342 const vector n(pSf[facei]/mag(pSf[facei]));
343 const Type uncorrectSnGrad(n & pgrad[facei]);
344 pgrad[facei] += n*(snGrad[facei] - uncorrectSnGrad);
345 }
346 }
347 }
348}
349
350
351// ************************************************************************* //
label n
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.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
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
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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< 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 for optional caching.
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< GradType, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
static tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > gradf(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const word &name)
Return the gradient of the given field calculated using Gauss' theorem on the given surface field.
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
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
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
mesh interpolate(rAU)
dynamicFvMesh & mesh
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
#define DebugPout
Report an information message using Foam::Pout.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
cellMask correctBoundaryConditions()
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299