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-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
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "geometricOneField.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 }
73
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);
118}
119
120
121// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123template<class Type>
126(
128 const tmp<surfaceScalarField>& tlambdas,
130)
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();
170
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)
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;
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();
264
265 for (label fi=0; fi<P.size(); fi++)
266 {
267 sfi[fi] = Sfi[fi] & (lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]]);
268 }
269
270 // Interpolate across coupled patches using given lambdas
271
273 Boundary& sfbf = sf.boundaryFieldRef();
274
275 forAll(lambdas.boundaryField(), pi)
276 {
277 const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
278 const typename SFType::Patch& pSf = Sf.boundaryField()[pi];
279 fvsPatchField<RetType>& psf = sfbf[pi];
280
281 if (vf.boundaryField()[pi].coupled())
282 {
283 psf =
284 pSf
285 & (
286 pLambda*vf.boundaryField()[pi].patchInternalField()
287 + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField()
288 );
289 }
290 else
291 {
292 psf = pSf & vf.boundaryField()[pi];
293 }
294 }
295
296 tlambdas.clear();
297
298// tsf.ref().oriented() = Sf.oriented();
299
300 return tsf;
301}
302
303
304template<class Type>
307(
309 const tmp<surfaceScalarField>& tlambdas
310)
311{
312 return dotInterpolate(geometricOneField(), vf, tlambdas);
313}
314
315
316template<class Type>
318<
320 <
324 >
325>
327(
328 const surfaceVectorField& Sf,
330) const
331{
332 if (surfaceInterpolation::debug)
333 {
335 << "Interpolating "
336 << vf.type() << " "
337 << vf.name()
338 << " from cells to faces"
339 << endl;
340 }
341
342 tmp
343 <
345 <
349 >
350 > tsf = dotInterpolate(Sf, vf, weights(vf));
351
352 tsf.ref().oriented() = Sf.oriented();
353
354 if (corrected())
355 {
356 tsf.ref() += Sf & correction(vf);
357 }
358
359 return tsf;
360}
361
362
363template<class Type>
365<
367 <
371 >
372>
374(
375 const surfaceVectorField& Sf,
377) const
378{
379 tmp
380 <
382 <
386 >
387 > tSfDotinterpVf = dotInterpolate(Sf, tvf());
388
389 tvf.clear();
390 return tSfDotinterpVf;
391}
392
393
394template<class Type>
397(
399) const
400{
401 if (surfaceInterpolation::debug)
402 {
404 << "Interpolating "
405 << vf.type() << " "
406 << vf.name()
407 << " from cells to faces"
408 << endl;
409 }
410
412 = interpolate(vf, weights(vf));
413
414 if (corrected())
415 {
416 tsf.ref() += correction(vf);
417 }
418
419 return tsf;
420}
421
422
423template<class Type>
426(
428) const
429{
431 = interpolate(tvf());
432 tvf.clear();
433 return tinterpVf;
434}
435
436
437// ************************************************************************* //
scalar y
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
const orientedType & oriented() const noexcept
Return oriented type.
Generic templated field type.
Definition: Field.H:82
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
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
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
bool interpolate() const noexcept
Same as isPointData()
Abstract base class for surface interpolation schemes.
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:52
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
type
Volume classification types.
Definition: volumeType.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
mesh interpolate(rAU)
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define InfoInFunction
Report an information message using Foam::Info.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.
const Vector< label > N(dict.get< Vector< label > >("N"))