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-2019 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  << MeshConstructorTablePtr_->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 (
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  );
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
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  <
221  typename Foam::innerProduct<typename SFType::value_type, Type>::type,
222  Foam::fvsPatchField,
224  >
225 >
227 (
228  const SFType& Sf,
230  const tmp<surfaceScalarField>& tlambdas
231 )
232 {
233  if (surfaceInterpolation::debug)
234  {
236  << "Interpolating "
237  << vf.type() << " "
238  << vf.name()
239  << " from cells to faces "
240  "without explicit correction"
241  << endl;
242  }
243 
244  typedef typename Foam::innerProduct<typename SFType::value_type, Type>::type
245  RetType;
246 
247  const surfaceScalarField& lambdas = tlambdas();
248 
249  const Field<Type>& vfi = vf;
250  const scalarField& lambda = lambdas;
251 
252  const fvMesh& mesh = vf.mesh();
253  const labelUList& P = mesh.owner();
254  const labelUList& N = mesh.neighbour();
255 
257  (
259  (
260  "interpolate("+vf.name()+')',
261  mesh,
262  Sf.dimensions()*vf.dimensions()
263  )
264  );
266 
267  Field<RetType>& sfi = sf.primitiveFieldRef();
268 
269  const typename SFType::Internal& Sfi = Sf();
270 
271  for (label fi=0; fi<P.size(); fi++)
272  {
273  sfi[fi] = Sfi[fi] & (lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]]);
274  }
275 
276  // Interpolate across coupled patches using given lambdas
277 
279  Boundary& sfbf = sf.boundaryFieldRef();
280 
281  forAll(lambdas.boundaryField(), pi)
282  {
283  const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
284  const typename SFType::Patch& pSf = Sf.boundaryField()[pi];
285  fvsPatchField<RetType>& psf = sfbf[pi];
286 
287  if (vf.boundaryField()[pi].coupled())
288  {
289  psf =
290  pSf
291  & (
292  pLambda*vf.boundaryField()[pi].patchInternalField()
293  + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField()
294  );
295  }
296  else
297  {
298  psf = pSf & vf.boundaryField()[pi];
299  }
300  }
301 
302  tlambdas.clear();
303 
304  return tsf;
305 }
306 
307 
308 template<class Type>
311 (
313  const tmp<surfaceScalarField>& tlambdas
314 )
315 {
316  return dotInterpolate(geometricOneField(), vf, tlambdas);
317 }
318 
319 
320 template<class Type>
321 Foam::tmp
322 <
324  <
325  typename Foam::innerProduct<Foam::vector, Type>::type,
326  Foam::fvsPatchField,
328  >
329 >
331 (
332  const surfaceVectorField& Sf,
334 ) const
335 {
336  if (surfaceInterpolation::debug)
337  {
339  << "Interpolating "
340  << vf.type() << " "
341  << vf.name()
342  << " from cells to faces"
343  << endl;
344  }
345 
346  tmp
347  <
349  <
350  typename Foam::innerProduct<Foam::vector, Type>::type,
353  >
354  > tsf = dotInterpolate(Sf, vf, weights(vf));
355 
356  if (corrected())
357  {
358  tsf.ref() += Sf & correction(vf);
359  }
360 
361  return tsf;
362 }
363 
364 
365 template<class Type>
366 Foam::tmp
367 <
369  <
370  typename Foam::innerProduct<Foam::vector, Type>::type,
371  Foam::fvsPatchField,
373  >
374 >
376 (
377  const surfaceVectorField& Sf,
379 ) const
380 {
381  tmp
382  <
384  <
385  typename Foam::innerProduct<Foam::vector, Type>::type,
388  >
389  > tSfDotinterpVf = dotInterpolate(Sf, tvf());
390  tvf.clear();
391  return tSfDotinterpVf;
392 }
393 
394 
395 template<class Type>
398 (
400 ) const
401 {
402  if (surfaceInterpolation::debug)
403  {
405  << "Interpolating "
406  << vf.type() << " "
407  << vf.name()
408  << " from cells to faces"
409  << endl;
410  }
411 
413  = interpolate(vf, weights(vf));
414 
415  if (corrected())
416  {
417  tsf.ref() += correction(vf);
418  }
419 
420  return tsf;
421 }
422 
423 
424 template<class Type>
427 (
429 ) const
430 {
432  = interpolate(tvf());
433  tvf.clear();
434  return tinterpVf;
435 }
436 
437 
438 // ************************************************************************* //
Foam::surfaceFields.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
Definition: IOobject.H:303
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:47
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:284
scalar y
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
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.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
volScalarField sf(fieldObject, mesh)
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.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
label N
Definition: createFields.H:22
A class for managing temporary objects.
Definition: PtrList.H:53
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Abstract base class for surface interpolation schemes.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.