interpolation.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-2022 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 
26 #include "interpolation.H"
27 #include "volFields.H"
28 #include "polyMesh.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const VolField<Type>& psi
37 )
38 :
39  psi_(psi),
40  mesh_(psi.mesh())
41 {}
42 
43 
44 template<class Type>
46 :
47  psi_(i.psi_),
48  mesh_(i.mesh_)
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
53 
54 template<class Type>
56 (
57  const word& interpolationType,
58  const VolField<Type>& psi
59 )
60 {
61  typename dictionaryConstructorTable::iterator cstrIter =
62  dictionaryConstructorTablePtr_->find(interpolationType);
63 
64  if (cstrIter == dictionaryConstructorTablePtr_->end())
65  {
67  << "Unknown interpolation type " << interpolationType
68  << " for field " << psi.name() << nl << nl
69  << "Valid interpolation types : " << endl
70  << dictionaryConstructorTablePtr_->sortedToc()
71  << exit(FatalError);
72  }
73 
74  return autoPtr<interpolation<Type>>(cstrIter()(psi));
75 }
76 
77 
78 template<class Type>
80 (
81  const dictionary& interpolationSchemes,
82  const VolField<Type>& psi
83 )
84 {
85  return New(word(interpolationSchemes.lookup(psi.name())), psi);
86 }
87 
88 
89 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90 
91 template<class Type>
93 (
94  const barycentric& coordinates,
95  const tetIndices& tetIs,
96  const label facei
97 ) const
98 {
99  return
101  (
102  tetIs.tet(mesh_).barycentricToPoint(coordinates),
103  tetIs.cell(),
104  facei
105  );
106 }
107 
108 
109 template<class Type, class InterpolationType>
112 (
114  const labelList& celli,
115  const labelList& facei
116 ) const
117 {
118  tmp<Field<Type>> tField(new Field<Type>(position.size()));
119 
120  Field<Type>& field = tField.ref();
121 
122  forAll(field, i)
123  {
124  field[i] =
125  static_cast<const InterpolationType&>(*this).interpolate
126  (
127  position[i],
128  celli[i],
129  isNull(facei) ? -1 : facei[i]
130  );
131  }
132 
133  return tField;
134 }
135 
136 
137 template<class Type, class InterpolationType>
140 (
142  const labelList& celli,
143  const labelList& tetFacei,
144  const labelList& tetPti,
145  const labelList& facei
146 ) const
147 {
149 
150  Field<Type>& field = tField.ref();
151 
152  forAll(field, i)
153  {
154  field[i] =
155  static_cast<const InterpolationType&>(*this).interpolate
156  (
157  coordinates[i],
158  tetIndices(celli[i], tetFacei[i], tetPti[i]),
159  isNull(facei) ? -1 : facei[i]
160  );
161  }
162 
163  return tField;
164 }
165 
166 
167 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:307
static direction size()
Return the number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpaceI.H:83
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
virtual tmp< Field< Type > > interpolate(const vectorField &position, const labelList &celli, const labelList &facei=NullObjectRef< labelList >()) const
Interpolate field to the given points in the given cells.
Abstract base class for interpolation.
Definition: interpolation.H:55
static autoPtr< interpolation< Type > > New(const word &interpolationType, const VolField< Type > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:56
interpolation(const VolField< Type > &psi)
Construct from components.
Definition: interpolation.C:35
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
label cell() const
Return the cell.
Definition: tetIndicesI.H:31
tetPointRef tet(const polyMesh &mesh, const pointField &meshPoints, const pointField &cellCentres) const
Return the geometry corresponding to this tet and the given.
Definition: tetIndicesI.H:109
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
Definition: tetrahedronI.H:228
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const volScalarField & psi
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
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:258
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267