patchExprDriverTemplates.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 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 
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Type>
35 (
36  const word& name
37 ) const
38 {
39  bool isPointVal = false;
40  bool isUniformVal = false;
41 
42  tmp<Field<Type>> tfield;
43 
44  if (hasVariable(name) && variable(name).isType<Type>())
45  {
46  const expressions::exprResult& var = variable(name);
47 
48  isPointVal = var.isPointValue();
49  isUniformVal = var.isUniform();
50 
51  tfield = var.cref<Type>().clone();
52  }
53  else if (isGlobalVariable<Type>(name, false))
54  {
55  const expressions::exprResult& var = lookupGlobal(name);
56 
57  isUniformVal = var.isUniform();
58 
59  tfield = var.cref<Type>().clone();
60  }
61 
62  if (tfield.valid())
63  {
64  const label fldLen = tfield().size();
65  const label len = (isPointVal ? this->pointSize() : this->size());
66 
67  if (returnReduce((fldLen == len), andOp<bool>()))
68  {
69  return tfield;
70  }
71 
72  if (!isUniformVal)
73  {
75  << "Variable " << name
76  << " does not fit the size and is not a uniform value." << nl
77  << "Using average value" << endl;
78  }
79 
80  return tmp<Field<Type>>::New(this->size(), gAverage(tfield));
81  }
82 
83  return tfield;
84 }
85 
86 
87 template<class Type>
90 {
91  return getField<Type>(name);
92 }
93 
94 
95 template<class Type>
98 {
99  return getField<Type>(name);
100 }
101 
102 
103 template<class Type>
106 {
107  return getField<Type>(name);
108 }
109 
110 
111 template<class Type>
114 {
115  tmp<Field<Type>> tfield = getVariableIfAvailable<Type>(name);
116 
117  if (tfield.valid())
118  {
119  return tfield;
120  }
121 
122 
126 
127  const objectRegistry& obr = this->mesh().thisDb();
128 
129  const vfieldType* vfield = obr.findObject<vfieldType>(name);
130  const sfieldType* sfield = obr.findObject<sfieldType>(name);
131  const pfieldType* pfield = obr.findObject<pfieldType>(name);
132 
133  // Local, temporary storage
134  tmp<vfieldType> t_vfield;
135  tmp<sfieldType> t_sfield;
136  tmp<pfieldType> t_pfield;
137 
138  if (searchFiles() && !vfield && !sfield && !pfield)
139  {
140  const word fldType = this->getTypeOfField(name);
141 
142  if (fldType == vfieldType::typeName)
143  {
144  t_vfield = this->readAndRegister<vfieldType>(name, mesh());
145  vfield = t_vfield.get();
146  }
147  else if (fldType == sfieldType::typeName)
148  {
149  t_sfield = this->readAndRegister<sfieldType>(name, mesh());
150  sfield = t_sfield.get();
151  }
152  else if (fldType == pfieldType::typeName)
153  {
154  t_pfield = this->readAndRegister<pfieldType>
155  (
156  name,
158  );
159  pfield = t_pfield.get();
160  }
161  }
162 
163  const label patchIndex = patch_.index();
164 
165  if (vfield)
166  {
167  return tmp<Field<Type>>::New
168  (
169  vfield->boundaryField()[patchIndex]
170  );
171  }
172 
173  if (sfield)
174  {
175  return tmp<Field<Type>>::New
176  (
177  sfield->boundaryField()[patchIndex]
178  );
179  }
180 
181  if (pfield)
182  {
183  return pfield->boundaryField()[patchIndex].patchInternalField();
184  }
185 
186 
188  << "No field '" << name << "' of type "
189  << pTraits<Type>::typeName << nl << nl
190  << vfieldType::typeName << " Fields: "
191  << flatOutput(obr.sortedNames<vfieldType>()) << nl
192  << sfieldType::typeName << " Fields: "
193  << flatOutput(obr.sortedNames<sfieldType>()) << nl
194  << pfieldType::typeName << " Fields: "
195  << flatOutput(obr.sortedNames<pfieldType>()) << nl
196  << exit(FatalError);
197 
198  return tmp<Field<Type>>::New();
199 }
200 
201 
202 template<class Type>
205 (
206  const Field<Type>& field
207 ) const
208 {
209  primitivePatchInterpolation interp(patch_.patch());
210 
211  return interp.pointToFaceInterpolate(field);
212 }
213 
214 
215 template<class Type>
218 (
219  const Field<Type>& field
220 ) const
221 {
222  primitivePatchInterpolation interp(patch_.patch());
223 
224  return interp.faceToPointInterpolate(field);
225 }
226 
227 
228 // ************************************************************************* //
Foam::objectRegistry::sortedNames
wordList sortedNames() const
The sorted names of all objects.
Definition: objectRegistry.C:170
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::expressions::patchExpr::parseDriver::getVariableIfAvailable
tmp< Field< Type > > getVariableIfAvailable(const word &fldName) const
Retrieve variable as field if possible.
Foam::PrimitivePatchInterpolation::pointToFaceInterpolate
tmp< Field< Type > > pointToFaceInterpolate(const Field< Type > &pf) const
Interpolate from points to faces.
Definition: PrimitivePatchInterpolation.C:233
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::expressions::patchExpr::parseDriver::getField
tmp< Field< Type > > getField(const word &fldName)
Return named field.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
primitivePatchInterpolation.H
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::expressions::patchExpr::parseDriver::getSurfaceField
tmp< Field< Type > > getSurfaceField(const word &fldName)
Retrieve field (surface field)
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:62
Foam::expressions::exprResult
A polymorphic field/result from evaluating an expression.
Definition: exprResult.H:128
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::tmp::get
T * get() noexcept
Return pointer without nullptr checking.
Definition: tmpI.H:227
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::andOp
Definition: ops.H:233
field
rDeltaTY field()
Foam::expressions::exprResult::isUniform
bool isUniform() const
True if single, uniform value.
Definition: exprResultI.H:265
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::expressions::exprResult::cref
const Field< Type > & cref() const
Return const reference to the field.
Foam::PrimitivePatchInterpolation::faceToPointInterpolate
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
Definition: PrimitivePatchInterpolation.C:176
Foam::expressions::patchExpr::parseDriver::pointToFace
tmp< Field< Type > > pointToFace(const Field< Type > &field) const
Interpolate point to face values.
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
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
Foam::PrimitivePatchInterpolation
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
Definition: PrimitivePatchInterpolation.H:53
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::flatOutput
FlatOutput< Container > flatOutput(const Container &obj, label len=0)
Global flatOutput function.
Definition: FlatOutput.H:85
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::expressions::exprResult::isPointValue
bool isPointValue(const bool isPointVal=true) const
True if representing point values, or test if same as isPointVal.
Definition: exprResultI.H:257
Foam::expressions::patchExpr::parseDriver::getVolField
tmp< Field< Type > > getVolField(const word &fldName)
Retrieve field (vol field)
Foam::expressions::patchExpr::parseDriver::faceToPoint
tmp< Field< Type > > faceToPoint(const Field< Type > &field) const
Interpolate face to point.
Foam::tmp::valid
bool valid() const noexcept
Definition: tmpI.H:206
Foam::expressions::patchExpr::parseDriver::getPointField
tmp< Field< Type > > getPointField(const word &fldName)
Retrieve field (point field)
Foam::GeometricField< Type, fvPatchField, volMesh >
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294