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-2020 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 // * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
37 )
38 :
39  psi_(psi),
40  mesh_(psi.mesh())
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
45 
46 template<class Type>
48 (
49  const word& interpolationType,
51 )
52 {
53  typename dictionaryConstructorTable::iterator cstrIter =
54  dictionaryConstructorTablePtr_->find(interpolationType);
55 
56  if (cstrIter == dictionaryConstructorTablePtr_->end())
57  {
59  << "Unknown interpolation type " << interpolationType
60  << " for field " << psi.name() << nl << nl
61  << "Valid interpolation types : " << endl
62  << dictionaryConstructorTablePtr_->sortedToc()
63  << exit(FatalError);
64  }
65 
66  return autoPtr<interpolation<Type>>(cstrIter()(psi));
67 }
68 
69 
70 template<class Type>
72 (
73  const dictionary& interpolationSchemes,
75 )
76 {
77  return New(word(interpolationSchemes.lookup(psi.name())), psi);
78 }
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
83 template<class Type>
85 (
86  const barycentric& coordinates,
87  const tetIndices& tetIs,
88  const label facei
89 ) const
90 {
91  return
93  (
94  tetIs.tet(mesh_).barycentricToPoint(coordinates),
95  tetIs.cell(),
96  facei
97  );
98 }
99 
100 
101 template<class Type, class InterpolationType>
104 (
105  const vectorField& position,
106  const labelField& celli,
107  const labelField& facei
108 ) const
109 {
110  tmp<Field<Type>> tField(new Field<Type>(position.size()));
111 
112  Field<Type>& field = tField.ref();
113 
114  forAll(field, i)
115  {
116  field[i] =
117  static_cast<const InterpolationType&>(*this).interpolate
118  (
119  position[i],
120  celli[i],
121  isNull(facei) ? -1 : facei[i]
122  );
123  }
124 
125  return tField;
126 }
127 
128 
129 template<class Type, class InterpolationType>
132 (
133  const Field<barycentric>& coordinates,
134  const labelField& celli,
135  const labelField& tetFacei,
136  const labelField& tetPti,
137  const labelField& facei
138 ) const
139 {
140  tmp<Field<Type>> tField(new Field<Type>(coordinates.size()));
141 
142  Field<Type>& field = tField.ref();
143 
144  forAll(field, i)
145  {
146  field[i] =
147  static_cast<const InterpolationType&>(*this).interpolate
148  (
149  coordinates[i],
150  tetIndices(celli[i], tetFacei[i], tetPti[i]),
151  isNull(facei) ? -1 : facei[i]
152  );
153  }
154 
155  return tField;
156 }
157 
158 
159 // ************************************************************************* //
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:48
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Generic GeometricField class.
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:98
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
Definition: tetrahedronI.H:254
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.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:46
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
label cell() const
Return the cell.
Definition: tetIndicesI.H:28
virtual tmp< Field< Type > > interpolate(const vectorField &position, const labelField &celli, const labelField &facei=NullObjectRef< labelField >()) const
Interpolate field to the given points in the given cells.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
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.
interpolation(const GeometricField< Type, fvPatchField, volMesh > &psi)
Construct from components.
Definition: interpolation.C:35
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const volScalarField & psi
rDeltaTY field()
A class for managing temporary objects.
Definition: PtrList.H:53
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844