surfaceInterpolationScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "geometricOneField.H"
30 #include "coupledFvPatchField.H"
31 
32 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
33 
34 template<class Type>
37 (
38  const fvMesh& mesh,
39  Istream& schemeData
40 )
41 {
42  if (schemeData.eof())
43  {
45  (
46  schemeData
47  ) << "Discretisation scheme not specified"
48  << endl << endl
49  << "Valid schemes are :" << endl
50  << MeshConstructorTablePtr_->sortedToc()
51  << exit(FatalIOError);
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  typename MeshConstructorTable::iterator constructorIter =
62  MeshConstructorTablePtr_->find(schemeName);
63 
64  if (constructorIter == MeshConstructorTablePtr_->end())
65  {
67  (
68  schemeData
69  ) << "Unknown discretisation scheme "
70  << schemeName << nl << nl
71  << "Valid schemes are :" << endl
72  << MeshConstructorTablePtr_->sortedToc()
73  << exit(FatalIOError);
74  }
75 
76  return constructorIter()(mesh, schemeData);
77 }
78 
79 
80 template<class Type>
83 (
84  const fvMesh& mesh,
85  const surfaceScalarField& faceFlux,
86  Istream& schemeData
87 )
88 {
89  if (schemeData.eof())
90  {
92  (
93  schemeData
94  ) << "Discretisation scheme not specified"
95  << endl << endl
96  << "Valid schemes are :" << endl
97  << MeshFluxConstructorTablePtr_->sortedToc()
98  << exit(FatalIOError);
99  }
100 
101  const word schemeName(schemeData);
102 
103  if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
104  {
106  << "Discretisation scheme = " << schemeName << endl;
107  }
108 
109  typename MeshFluxConstructorTable::iterator constructorIter =
110  MeshFluxConstructorTablePtr_->find(schemeName);
111 
112  if (constructorIter == MeshFluxConstructorTablePtr_->end())
113  {
115  (
116  schemeData
117  ) << "Unknown discretisation scheme "
118  << schemeName << nl << nl
119  << "Valid schemes are :" << endl
120  << MeshFluxConstructorTablePtr_->sortedToc()
121  << exit(FatalIOError);
122  }
123 
124  return constructorIter()(mesh, faceFlux, schemeData);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
130 template<class Type>
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
137 template<class Type>
140 (
141  const VolField<Type>& vf,
142  const tmp<surfaceScalarField>& tlambdas,
143  const tmp<surfaceScalarField>& tys
144 )
145 {
146  if (surfaceInterpolation::debug)
147  {
149  << "Interpolating "
150  << vf.type() << " "
151  << vf.name()
152  << " from cells to faces "
153  "without explicit correction"
154  << endl;
155  }
156 
157  const surfaceScalarField& lambdas = tlambdas();
158  const surfaceScalarField& ys = tys();
159 
160  const Field<Type>& vfi = vf;
161  const scalarField& lambda = lambdas;
162  const scalarField& y = ys;
163 
164  const fvMesh& mesh = vf.mesh();
165  const labelUList& P = mesh.owner();
166  const labelUList& N = mesh.neighbour();
167 
169  (
171  (
172  "interpolate("+vf.name()+')',
173  mesh,
174  vf.dimensions()
175  )
176  );
177  SurfaceField<Type>& sf = tsf.ref();
178 
179  Field<Type>& sfi = sf.primitiveFieldRef();
180 
181  for (label fi=0; fi<P.size(); fi++)
182  {
183  sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
184  }
185 
186 
187  // Interpolate across coupled patches using given lambdas and ys
188  typename SurfaceField<Type>::
189  Boundary& sfbf = sf.boundaryFieldRef();
190 
191  forAll(lambdas.boundaryField(), pi)
192  {
193  const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
194  const fvsPatchScalarField& pY = ys.boundaryField()[pi];
195 
196  if (vf.boundaryField()[pi].coupled())
197  {
198  sfbf[pi] =
199  pLambda*vf.boundaryField()[pi].patchInternalField()
200  + pY*vf.boundaryField()[pi].patchNeighbourField();
201  }
202  else
203  {
204  sfbf[pi] = vf.boundaryField()[pi];
205  }
206  }
207 
208  tlambdas.clear();
209  tys.clear();
210 
211  return tsf;
212 }
213 
214 
215 template<class Type>
216 template<class SFType>
217 Foam::tmp
218 <
220  <
222  >
223 >
225 (
226  const SFType& Sf,
227  const VolField<Type>& vf,
228  const tmp<surfaceScalarField>& tlambdas
229 )
230 {
231  if (surfaceInterpolation::debug)
232  {
234  << "Interpolating "
235  << vf.type() << " "
236  << vf.name()
237  << " from cells to faces "
238  "without explicit correction"
239  << endl;
240  }
241 
243  RetType;
244 
245  const surfaceScalarField& lambdas = tlambdas();
246 
247  const Field<Type>& vfi = vf;
248  const scalarField& lambda = lambdas;
249 
250  const fvMesh& mesh = vf.mesh();
251  const labelUList& P = mesh.owner();
252  const labelUList& N = mesh.neighbour();
253 
255  (
257  (
258  "interpolate("+vf.name()+')',
259  mesh,
260  Sf.dimensions()*vf.dimensions()
261  )
262  );
263  SurfaceField<RetType>& sf = tsf.ref();
264 
265  Field<RetType>& sfi = sf.primitiveFieldRef();
266 
267  const typename SFType::Internal& Sfi = Sf();
268 
269  for (label fi=0; fi<P.size(); fi++)
270  {
271  sfi[fi] = Sfi[fi] & (lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]]);
272  }
273 
274  // Interpolate across coupled patches using given lambdas
275 
276  typename SurfaceField<RetType>::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  & (
289  pLambda*vf.boundaryField()[pi].patchInternalField()
290  + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField()
291  );
292  }
293  else
294  {
295  psf = pSf & vf.boundaryField()[pi];
296  }
297  }
298 
299  tlambdas.clear();
300 
301  return tsf;
302 }
303 
304 
305 template<class Type>
308 (
309  const VolField<Type>& vf,
310  const tmp<surfaceScalarField>& tlambdas
311 )
312 {
313  return dotInterpolate(geometricOneField(), vf, tlambdas);
314 }
315 
316 
317 template<class Type>
318 Foam::tmp
319 <
321  <
323  >
324 >
326 (
327  const surfaceVectorField& Sf,
328  const VolField<Type>& vf
329 ) const
330 {
331  if (surfaceInterpolation::debug)
332  {
334  << "Interpolating "
335  << vf.type() << " "
336  << vf.name()
337  << " from cells to faces"
338  << endl;
339  }
340 
341  tmp
342  <
344  > tsf = dotInterpolate(Sf, vf, weights(vf));
345 
346  if (corrected())
347  {
348  tsf.ref() += Sf & correction(vf);
349  }
350 
351  return tsf;
352 }
353 
354 
355 template<class Type>
356 Foam::tmp
357 <
359  <
361  >
362 >
364 (
365  const surfaceVectorField& Sf,
366  const tmp<VolField<Type>>& tvf
367 ) const
368 {
369  tmp
370  <
372  > tSfDotinterpVf = dotInterpolate(Sf, tvf());
373  tvf.clear();
374  return tSfDotinterpVf;
375 }
376 
377 
378 template<class Type>
381 (
382  const VolField<Type>& vf
383 ) const
384 {
385  if (surfaceInterpolation::debug)
386  {
388  << "Interpolating "
389  << vf.type() << " "
390  << vf.name()
391  << " from cells to faces"
392  << endl;
393  }
394 
396  = interpolate(vf, weights(vf));
397 
398  if (corrected())
399  {
400  tsf.ref() += correction(vf);
401  }
402 
403  return tsf;
404 }
405 
406 
407 template<class Type>
410 (
411  const tmp<VolField<Type>>& tvf
412 ) const
413 {
414  tmp<SurfaceField<Type>> tinterpVf
415  = interpolate(tvf());
416  tvf.clear();
417  return tinterpVf;
418 }
419 
420 
421 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:336
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:451
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:457
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
Definition: products.H:115
Abstract base class for surface interpolation schemes.
static tmp< SurfaceField< typename innerProduct< typename SFType::value_type, Type >::type > > dotInterpolate(const SFType &Sf, const VolField< Type > &vf, const tmp< surfaceScalarField > &tlambdas)
Return the face-interpolate of the given cell field.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
volScalarField sf(fieldObject, mesh)
dimensionedScalar lambda(viscosity->lookup("lambda"))
#define InfoInFunction
Report an information message using Foam::Info.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
IOerror FatalIOError
SurfaceField< vector > surfaceVectorField
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Foam::surfaceFields.