Loading...
Searching...
No Matches
surfaceInterpolationScheme.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) 2019-2023 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
30#include "volFields.H"
31#include "surfaceFields.H"
33#include "coupledFvPatchField.H"
34
35// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
36
37template<class Type>
40(
41 const fvMesh& mesh,
42 Istream& schemeData
43)
44{
45 if (schemeData.eof())
46 {
47 FatalIOErrorInFunction(schemeData)
48 << "Discretisation scheme not specified\n\n"
49 << "Valid schemes:\n"
50 << MeshConstructorTablePtr_->sortedToc()
52 }
53
54 const word schemeName(schemeData);
55
56 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
57 {
58 InfoInFunction << "Discretisation scheme = " << schemeName << endl;
59 }
60
61 auto* ctorPtr = MeshConstructorTable(schemeName);
62
63 if (!ctorPtr)
64 {
66 (
67 schemeData,
68 "discretisation",
69 schemeName,
70 *MeshConstructorTablePtr_
71 ) << exit(FatalIOError);
72 }
74 return ctorPtr(mesh, schemeData);
75}
76
77
78template<class Type>
81(
82 const fvMesh& mesh,
83 const surfaceScalarField& faceFlux,
84 Istream& schemeData
85)
86{
87 if (schemeData.eof())
88 {
89 FatalIOErrorInFunction(schemeData)
90 << "Discretisation scheme not specified"
91 << endl << endl
92 << "Valid schemes are :" << endl
93 << MeshConstructorTablePtr_->sortedToc()
95 }
96
97 const word schemeName(schemeData);
98
99 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
100 {
101 InfoInFunction << "Discretisation scheme = " << schemeName << endl;
102 }
103
104 auto* ctorPtr = MeshFluxConstructorTable(schemeName);
105
106 if (!ctorPtr)
107 {
109 (
110 schemeData,
111 "discretisation",
112 schemeName,
113 *MeshFluxConstructorTablePtr_
114 ) << exit(FatalIOError);
115 }
116
117 return ctorPtr(mesh, faceFlux, schemeData);
119
120
121// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122
123template<class Type>
126(
128 const tmp<surfaceScalarField>& tlambdas,
129 const tmp<surfaceScalarField>& tys
131{
132 if (surfaceInterpolation::debug)
133 {
135 << "Interpolating "
136 << vf.type() << " "
137 << vf.name()
138 << " from cells to faces without explicit correction"
139 << endl;
140 }
141
142 const surfaceScalarField& lambdas = tlambdas();
143 const surfaceScalarField& ys = tys();
144
145 const Field<Type>& vfi = vf;
146 const scalarField& lambda = lambdas;
147 const scalarField& y = ys;
148
149 const fvMesh& mesh = vf.mesh();
150 const labelUList& P = mesh.owner();
151 const labelUList& N = mesh.neighbour();
152
154 (
156 (
158 (
159 "interpolate("+vf.name()+')',
160 vf.instance(),
161 vf.db()
162 ),
163 mesh,
164 vf.dimensions()
165 )
166 );
168
169 Field<Type>& sfi = sf.primitiveFieldRef();
171 for (label fi=0; fi<P.size(); fi++)
172 {
173 sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
174 }
175
176
177 // Interpolate across coupled patches using given lambdas and ys
179 Boundary& sfbf = sf.boundaryFieldRef();
180
181 forAll(lambdas.boundaryField(), pi)
182 {
183 const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
184 const fvsPatchScalarField& pY = ys.boundaryField()[pi];
185
186 if (vf.boundaryField()[pi].coupled())
187 {
188 sfbf[pi] =
189 pLambda*vf.boundaryField()[pi].patchInternalField()
190 + pY*vf.boundaryField()[pi].patchNeighbourField();
191 }
192 else
193 {
194 sfbf[pi] = vf.boundaryField()[pi];
195 }
196 }
197
198 tlambdas.clear();
199 tys.clear();
200
201 return tsf;
202}
203
204
205template<class Type>
206template<class SFType>
208<
210 <
214 >
215>
217(
218 const SFType& Sf,
220 const tmp<surfaceScalarField>& tlambdas
221)
222{
223 if (surfaceInterpolation::debug)
224 {
226 << "Interpolating "
227 << vf.type() << " "
228 << vf.name()
229 << " from cells to faces without explicit correction"
230 << endl;
231 }
232
234 RetType;
235
236 const surfaceScalarField& lambdas = tlambdas();
237
238 const Field<Type>& vfi = vf;
239 const scalarField& lambda = lambdas;
240
241 const fvMesh& mesh = vf.mesh();
242 const labelUList& P = mesh.owner();
243 const labelUList& N = mesh.neighbour();
244
246 (
248 (
250 (
251 "interpolate("+vf.name()+')',
252 vf.instance(),
253 vf.db()
254 ),
255 mesh,
256 Sf.dimensions()*vf.dimensions()
257 )
258 );
260
262
263 const typename SFType::Internal& Sfi = Sf.internalField();
264
265 for (label fi=0; fi<P.size(); fi++)
266 {
267 // Same as:
268 // sfi[fi] = Sfi[fi] & lerp(vfi[N[fi]], vfi[P[fi]], lambda[fi]);
269 // but maybe the compiler notices the fused multiply add form
270 sfi[fi] = Sfi[fi] & (lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]]);
271 }
272
273 // Interpolate across coupled patches using given lambdas
274
276 Boundary& sfbf = sf.boundaryFieldRef();
277
278 forAll(lambdas.boundaryField(), pi)
279 {
280 const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
281 const typename SFType::Patch& pSf = Sf.boundaryField()[pi];
282 fvsPatchField<RetType>& psf = sfbf[pi];
283
284 if (vf.boundaryField()[pi].coupled())
285 {
286 psf =
287 pSf
288 & lerp
289 (
290 vf.boundaryField()[pi].patchNeighbourField(),
291 vf.boundaryField()[pi].patchInternalField(),
292 pLambda
293 );
294 }
295 else
296 {
297 psf = pSf & vf.boundaryField()[pi];
298 }
299 }
300
301 tlambdas.clear();
302
303// tsf.ref().oriented() = Sf.oriented();
312(
314 const tmp<surfaceScalarField>& tlambdas
315)
316{
317 return dotInterpolate(geometricOneField(), vf, tlambdas);
318}
319
320
321template<class Type>
323<
325 <
329 >
330>
332(
333 const surfaceVectorField& Sf,
335) const
336{
337 if (surfaceInterpolation::debug)
338 {
340 << "Interpolating "
341 << vf.type() << " "
342 << vf.name()
343 << " from cells to faces"
344 << endl;
345 }
346
347 tmp
348 <
350 <
354 >
355 > tsf = dotInterpolate(Sf, vf, weights(vf));
356
357 tsf.ref().oriented() = Sf.oriented();
358
359 if (corrected())
360 {
361 tsf.ref() += Sf & correction(vf);
362 }
363
364 return tsf;
365}
366
367
368template<class Type>
370<
391 >
392 > tSfDotinterpVf = dotInterpolate(Sf, tvf());
393
394 tvf.clear();
395 return tSfDotinterpVf;
396}
397
398
399template<class Type>
402(
404) const
405{
406 if (surfaceInterpolation::debug)
407 {
409 << "Interpolating "
410 << vf.type() << " "
411 << vf.name()
412 << " from cells to faces"
413 << endl;
414 }
415
416 tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf
417 = interpolate(vf, weights(vf));
418
419 if (corrected())
420 {
421 tsf.ref() += correction(vf);
422 }
433) const
434{
436 = interpolate(tvf());
437 tvf.clear();
438 return tinterpVf;
439}
440
441
442// ************************************************************************* //
scalar y
constexpr scalar pi(M_PI)
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
orientedType oriented() const noexcept
Return oriented type.
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.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
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.
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 objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
bool eof() const noexcept
True if end of input seen.
Definition IOstream.H:289
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
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
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) -2 >::type type
Definition products.H:155
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
Return the face-interpolate of the given cell field.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Mesh data needed to do the Finite Volume discretisation.
Definition surfaceMesh.H:47
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh interpolate(rAU)
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define InfoInFunction
Report an information message using Foam::Info.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
const Vector< label > N(dict.get< Vector< label > >("N"))