volumeExprDriverTemplates.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2019-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "exprOps.H"
29 #include "FieldOps.H"
30 #include "surfaceInterpolate.H"
31 #include "volPointInterpolation.H"
32 #include "interpolatePointToCell.H"
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const Field<Type>& fld
40 )
41 {
42  if (isLogical_)
43  {
44  // Eg, volScalarField -> volLogicalField
45  resultType_.replace("Scalar", "Logical");
46 
47  Field<bool> bools(fld.size());
49 
50  this->result().setResult(std::move(bools), this->isPointData());
51  }
52  else
53  {
54  // Deep copy
55  this->result().setResult(fld, this->isPointData());
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class Type>
64 (
66  bool logical
67 )
68 {
70 
71  resultField_.reset(nullptr);
72 
73  // Characteristics
74  resultType_ = pTraits<fieldType>::typeName;
75  isLogical_ = logical;
76  fieldGeoType_ = VOLUME_DATA;
77 
78  // Assign dimensions
79  if (hasDimensions_ && !logical)
80  {
81  ptr->dimensions().reset(resultDimensions_);
82  }
83 
84  setInternalFieldResult(ptr->primitiveField());
85 
86  // Take ownership
87  resultField_.reset(ptr);
88 }
89 
90 
91 template<class Type>
93 (
95  bool logical
96 )
97 {
99 
100  resultField_.reset(nullptr);
101 
102  // Characteristics
103  resultType_ = pTraits<fieldType>::typeName;
104  isLogical_ = logical;
105  fieldGeoType_ = FACE_DATA;
106 
107  // Assign dimensions
108  if (hasDimensions_ && !logical)
109  {
110  ptr->dimensions().reset(resultDimensions_);
111  }
112 
113  setInternalFieldResult(ptr->primitiveField());
114 
115  // Take ownership
116  resultField_.reset(ptr);
117 }
118 
119 
120 template<class Type>
122 (
124  bool logical
125 )
126 {
128 
129  resultField_.reset(nullptr);
130 
131  // Characteristics
132  resultType_ = pTraits<fieldType>::typeName;
133  isLogical_ = logical;
134  fieldGeoType_ = POINT_DATA;
135 
136  // Assign dimensions
137  if (hasDimensions_ && !logical)
138  {
139  ptr->dimensions().reset(resultDimensions_);
140  }
141 
142  setInternalFieldResult(ptr->primitiveField());
143 
144  // Take ownership
145  resultField_.reset(ptr);
146 }
147 
148 
149 template<class GeomField>
150 const GeomField*
152 {
153  return dynamic_cast<const GeomField*>(resultField_.get());
154 }
155 
156 
157 template<class GeomField>
158 const GeomField*
160 (
161  bool logical,
162  bool dieOnNull
163 ) const
164 {
165  const regIOobject* ptr = resultField_.get();
166 
167  if (dieOnNull && !ptr)
168  {
170  << "No result available. Requested "
172  << exit(FatalError);
173  }
174 
175  if (isLogical_ == logical)
176  {
177  return dynamic_cast<const GeomField*>(ptr);
178  }
179 
180  return nullptr;
181 }
182 
183 
184 template<class Type>
187 (
188  const word& fldName,
189  bool getOldTime
190 )
191 {
193 
194  return this->getOrReadField<fieldType>
195  (
196  fldName,
197  true, // mandatory
198  getOldTime
199  );
200 }
201 
202 
203 template<class Type>
206 (
207  const word& fldName,
208  bool getOldTime
209 )
210 {
212 
213  return this->getOrReadField<fieldType>
214  (
215  fldName,
216  true, // mandatory
217  getOldTime
218  );
219 }
220 
221 
222 template<class Type>
225 (
226  const word& fldName,
227  bool getOldTime
228 )
229 {
231 
232  return this->getOrReadPointField<fieldType>
233  (
234  fldName,
235  true, // mandatory
236  getOldTime
237  );
238 }
239 
240 
241 template<class Type>
244 (
245  const Type& val
246 ) const
247 {
249 
250  return fieldType::New
251  (
252  word("constant.") + word(pTraits<Type>::typeName),
253  mesh(),
254  dimensioned<Type>("", dimless, val)
255  );
256 }
257 
258 
259 template<class Type>
262 (
263  const Type& val
264 ) const
265 {
267 
268  return fieldType::New
269  (
270  word("constant.") + word(pTraits<Type>::typeName),
271  mesh(),
272  dimensioned<Type>("", dimless, val)
273  );
274 }
275 
276 
277 template<class Type>
280 (
281  const Type& val
282 ) const
283 {
285 
286  return fieldType::New
287  (
288  word("constant.") + word(pTraits<Type>::typeName),
289  pointMesh::New(mesh()),
290  dimensioned<Type>("", dimless, val)
291  );
292 }
293 
294 
295 template<class Type>
298 (
300 ) const
301 {
302  return fvc::interpolate(field);
303 }
304 
305 
306 template<class Type>
309 (
311 ) const
312 {
313  volPointInterpolation interp(this->mesh());
314  return interp.interpolate(field);
315 }
316 
317 
318 template<class Type>
321 (
323 ) const
324 {
325  auto tresult = newVolField<Type>();
326  auto& result = tresult.ref();
327 
328  forAll(result,celli)
329  {
330  result[celli] = interpolatePointToCell(field, celli);
331  }
332 
333  return tresult;
334 }
335 
336 
337 // ************************************************************************* //
FieldOps.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::HashSetOps::bools
List< bool > bools(const labelHashSet &locations)
Definition: HashOps.C:81
Foam::expressions::volumeExpr::parseDriver::cellToPoint
tmp< GeometricField< Type, pointPatchField, pointMesh > > cellToPoint(const GeometricField< Type, fvPatchField, volMesh > &field) const
Interpolate cell to point values.
Foam::expressions::volumeExpr::parseDriver::isResultType
const GeoField * isResultType() const
Test if stored result pointer is the specified type.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::interpolatePointToCell
Type interpolatePointToCell(const GeometricField< Type, pointPatchField, pointMesh > &ptf, const label celli)
Definition: interpolatePointToCell.C:34
interpolatePointToCell.H
Interpolates (averages) the vertex values to the cell center.
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Foam::expressions::VOLUME_DATA
Volume data.
Definition: exprFieldAssociation.H:48
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::expressions::FACE_DATA
Face data.
Definition: exprFieldAssociation.H:47
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::expressions::volumeExpr::parseDriver::getPointField
tmp< GeometricField< Type, pointPatchField, pointMesh > > getPointField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
Foam::expressions::volumeExpr::parseDriver::resultField_
autoPtr< regIOobject > resultField_
The results (volume, surface, point)
Definition: volumeExprDriver.H:243
Foam::expressions::volumeExpr::parseDriver::getSurfaceField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > getSurfaceField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
field
rDeltaTY field()
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::FatalError
error FatalError
Foam::expressions::boolOp
Convert [0-1] values (usually scalars) as false/true.
Definition: exprOps.H:56
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::expressions::volumeExpr::parseDriver::setInternalFieldResult
void setInternalFieldResult(const Field< Type > &fld)
Deep-copy the internalField as a result.
Definition: volumeExprDriverTemplates.C:38
Foam::fvc::interpolate
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.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
volPointInterpolation.H
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::expressions::volumeExpr::parseDriver::getVolField
tmp< GeometricField< Type, fvPatchField, volMesh > > getVolField(const word &fldName, bool getOldTime=false)
Retrieve field (vol field)
Foam::expressions::volumeExpr::parseDriver::pointToCell
tmp< GeometricField< Type, fvPatchField, volMesh > > pointToCell(const GeometricField< Type, pointPatchField, pointMesh > &field) const
Interpolate point to cell values.
Foam::expressions::POINT_DATA
Point data.
Definition: exprFieldAssociation.H:46
Foam::expressions::volumeExpr::parseDriver::newVolField
tmp< GeometricField< Type, fvPatchField, volMesh > > newVolField(const Type &val=pTraits< Type >::zero) const
Return a new volume field with the mesh size.
Foam::expressions::volumeExpr::parseDriver::setResult
void setResult(GeometricField< Type, fvPatchField, volMesh > *ptr, bool logical=false)
Set result (vol field)
Definition: volumeExprDriverTemplates.C:64
Foam::FieldOps::assign
void assign(Field< Tout > &result, const Field< T1 > &a, const UnaryOp &op)
Populate a field as the result of a unary operation on an input.
Definition: FieldOps.C:35
Foam::expressions::volumeExpr::parseDriver::newSurfaceField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > newSurfaceField(const Type &val=pTraits< Type >::zero) const
Return a new surface field with the mesh nInternalFaces size.
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::volPointInterpolation
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Definition: volPointInterpolation.H:59
Foam::expressions::volumeExpr::parseDriver::cellToFace
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > cellToFace(const GeometricField< Type, fvPatchField, volMesh > &field) const
Interpolate cell to face values.
Foam::expressions::volumeExpr::parseDriver::newPointField
tmp< GeometricField< Type, pointPatchField, pointMesh > > newPointField(const Type &val=pTraits< Type >::zero) const
Return a new point field with the mesh nPoints size.
exprOps.H
Operations involving expressions.
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189