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 | Copyright (C) 2011-2015 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 Description
25  Abstract base class for surface interpolation schemes.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "coupledFvPatchField.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
40 
41 // Return weighting factors for scheme given by name in dictionary
42 template<class Type>
43 tmp<surfaceInterpolationScheme<Type> > surfaceInterpolationScheme<Type>::New
44 (
45  const fvMesh& mesh,
46  Istream& schemeData
47 )
48 {
49  if (schemeData.eof())
50  {
52  (
53  schemeData
54  ) << "Discretisation scheme not specified"
55  << endl << endl
56  << "Valid schemes are :" << endl
57  << MeshConstructorTablePtr_->sortedToc()
58  << exit(FatalIOError);
59  }
60 
61  const word schemeName(schemeData);
62 
63  if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
64  {
65  Info<< "surfaceInterpolationScheme<Type>::New"
66  "(const fvMesh&, Istream&)"
67  " : discretisation scheme = "
68  << schemeName
69  << endl;
70  }
71 
72  typename MeshConstructorTable::iterator constructorIter =
73  MeshConstructorTablePtr_->find(schemeName);
74 
75  if (constructorIter == MeshConstructorTablePtr_->end())
76  {
78  (
79  schemeData
80  ) << "Unknown discretisation scheme "
81  << schemeName << nl << nl
82  << "Valid schemes are :" << endl
83  << MeshConstructorTablePtr_->sortedToc()
84  << exit(FatalIOError);
85  }
86 
87  return constructorIter()(mesh, schemeData);
88 }
89 
90 
91 // Return weighting factors for scheme given by name in dictionary
92 template<class Type>
94 (
95  const fvMesh& mesh,
96  const surfaceScalarField& faceFlux,
97  Istream& schemeData
98 )
99 {
100  if (schemeData.eof())
101  {
103  (
104  schemeData
105  ) << "Discretisation scheme not specified"
106  << endl << endl
107  << "Valid schemes are :" << endl
108  << MeshConstructorTablePtr_->sortedToc()
109  << exit(FatalIOError);
110  }
111 
112  const word schemeName(schemeData);
113 
114  if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
115  {
116  Info<< "surfaceInterpolationScheme<Type>::New"
117  "(const fvMesh&, const surfaceScalarField&, Istream&)"
118  " : discretisation scheme = "
119  << schemeName
120  << endl;
121  }
122 
123  typename MeshFluxConstructorTable::iterator constructorIter =
124  MeshFluxConstructorTablePtr_->find(schemeName);
125 
126  if (constructorIter == MeshFluxConstructorTablePtr_->end())
127  {
129  (
130  schemeData
131  ) << "Unknown discretisation scheme "
132  << schemeName << nl << nl
133  << "Valid schemes are :" << endl
134  << MeshFluxConstructorTablePtr_->sortedToc()
135  << exit(FatalIOError);
136  }
137 
138  return constructorIter()(mesh, faceFlux, schemeData);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
143 
144 template<class Type>
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 //- Return the face-interpolate of the given cell field
152 // with the given owner and neighbour weighting factors
153 template<class Type>
156 (
158  const tmp<surfaceScalarField>& tlambdas,
159  const tmp<surfaceScalarField>& tys
160 )
161 {
162  if (surfaceInterpolation::debug)
163  {
164  Info<< "surfaceInterpolationScheme<Type>::uncorrectedInterpolate"
165  "(const GeometricField<Type, fvPatchField, volMesh>&, "
166  "const tmp<surfaceScalarField>&, "
167  "const tmp<surfaceScalarField>&) : "
168  "interpolating "
169  << vf.type() << " "
170  << vf.name()
171  << " from cells to faces "
172  "without explicit correction"
173  << endl;
174  }
175 
176  const surfaceScalarField& lambdas = tlambdas();
177  const surfaceScalarField& ys = tys();
178 
179  const Field<Type>& vfi = vf.internalField();
180  const scalarField& lambda = lambdas.internalField();
181  const scalarField& y = ys.internalField();
182 
183  const fvMesh& mesh = vf.mesh();
184  const labelUList& P = mesh.owner();
185  const labelUList& N = mesh.neighbour();
186 
188  (
190  (
191  IOobject
192  (
193  "interpolate("+vf.name()+')',
194  vf.instance(),
195  vf.db()
196  ),
197  mesh,
198  vf.dimensions()
199  )
200  );
202 
203  Field<Type>& sfi = sf.internalField();
204 
205  for (label fi=0; fi<P.size(); fi++)
206  {
207  sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
208  }
209 
210 
211  // Interpolate across coupled patches using given lambdas and ys
212 
213  forAll(lambdas.boundaryField(), pi)
214  {
215  const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
216  const fvsPatchScalarField& pY = ys.boundaryField()[pi];
217 
218  if (vf.boundaryField()[pi].coupled())
219  {
220  sf.boundaryField()[pi] =
221  pLambda*vf.boundaryField()[pi].patchInternalField()
222  + pY*vf.boundaryField()[pi].patchNeighbourField();
223  }
224  else
225  {
226  sf.boundaryField()[pi] = vf.boundaryField()[pi];
227  }
228  }
229 
230  tlambdas.clear();
231  tys.clear();
232 
233  return tsf;
234 }
235 
236 
237 //- Return the face-interpolate of the given cell field
238 // with the given weighting factors
239 template<class Type>
242 (
244  const tmp<surfaceScalarField>& tlambdas
245 )
246 {
247  if (surfaceInterpolation::debug)
248  {
249  Info<< "surfaceInterpolationScheme<Type>::interpolate"
250  "(const GeometricField<Type, fvPatchField, volMesh>&, "
251  "const tmp<surfaceScalarField>&) : "
252  "interpolating "
253  << vf.type() << " "
254  << vf.name()
255  << " from cells to faces "
256  "without explicit correction"
257  << endl;
258  }
259 
260  const surfaceScalarField& lambdas = tlambdas();
261 
262  const Field<Type>& vfi = vf.internalField();
263  const scalarField& lambda = lambdas.internalField();
264 
265  const fvMesh& mesh = vf.mesh();
266  const labelUList& P = mesh.owner();
267  const labelUList& N = mesh.neighbour();
268 
270  (
272  (
273  IOobject
274  (
275  "interpolate("+vf.name()+')',
276  vf.instance(),
277  vf.db()
278  ),
279  mesh,
280  vf.dimensions()
281  )
282  );
284 
285  Field<Type>& sfi = sf.internalField();
286 
287  for (label fi=0; fi<P.size(); fi++)
288  {
289  sfi[fi] = lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]];
290  }
291 
292  // Interpolate across coupled patches using given lambdas
293 
294  forAll(lambdas.boundaryField(), pi)
295  {
296  const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
297 
298  if (vf.boundaryField()[pi].coupled())
299  {
300  tsf().boundaryField()[pi] =
301  pLambda*vf.boundaryField()[pi].patchInternalField()
302  + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField();
303  }
304  else
305  {
306  sf.boundaryField()[pi] = vf.boundaryField()[pi];
307  }
308  }
309 
310  tlambdas.clear();
311 
312  return tsf;
313 }
314 
315 
316 //- Return the face-interpolate of the given cell field
317 // with explicit correction
318 template<class Type>
321 (
323 ) const
324 {
325  if (surfaceInterpolation::debug)
326  {
327  Info<< "surfaceInterpolationScheme<Type>::interpolate"
328  "(const GeometricField<Type, fvPatchField, volMesh>&) : "
329  "interpolating "
330  << vf.type() << " "
331  << vf.name()
332  << " from cells to faces"
333  << endl;
334  }
335 
337  = interpolate(vf, weights(vf));
338 
339  if (corrected())
340  {
341  tsf() += correction(vf);
342  }
343 
344  return tsf;
345 }
346 
347 
348 //- Return the face-interpolate of the given cell field
349 // with explicit correction
350 template<class Type>
353 (
355 ) const
356 {
358  = interpolate(tvf());
359  tvf.clear();
360  return tinterpVf;
361 }
362 
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 } // End namespace Foam
367 
368 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::IOstream::eof
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
Foam::surfaceInterpolationScheme::interpolate
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.
Definition: surfaceInterpolationScheme.C:156
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Foam::surfaceInterpolationScheme::~surfaceInterpolationScheme
virtual ~surfaceInterpolationScheme()
Destructor.
Definition: surfaceInterpolationScheme.C:145
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::label
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
Foam::Field< Type >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
sf
volScalarField sf(fieldObject, mesh)
Foam::fvMesh::owner
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::surfaceInterpolationScheme::New
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: surfaceInterpolationScheme.C:44
Foam::UList
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
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:55
Foam::constant::mathematical::pi
const scalar pi(M_PI)
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
coupledFvPatchField.H
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
surfaceInterpolationScheme.H
y
scalar y
Definition: LISASMDCalcMethod1.H:14
N
label N
Definition: createFields.H:22