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-2020 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 hasPointData = false;
40 
42 
43  if (hasVariable(name) && variable(name).isType<Type>())
44  {
45  tvar.cref(variable(name));
46  hasPointData = tvar().isPointData();
47  }
48  else if (isGlobalVariable<Type>(name))
49  {
50  tvar.cref(lookupGlobal(name));
51  }
52 
53 
54  if (tvar.valid())
55  {
56  const auto& var = tvar.cref();
57  const Field<Type>& vals = var.cref<Type>();
58 
59  const label len = (hasPointData ? this->pointSize() : this->size());
60 
61  if (returnReduce((vals.size() == len), andOp<bool>()))
62  {
63  // Return a copy of the field
64  return tmp<Field<Type>>::New(vals);
65  }
66 
67  if (!var.isUniform())
68  {
70  << "Variable " << name
71  << " is nonuniform and does not fit the size"
72  << ". Using average" << endl;
73  }
74 
75  return tmp<Field<Type>>::New(this->size(), gAverage(vals));
76  }
77 
78  return nullptr;
79 }
80 
81 
82 template<class Type>
85 {
86  return getField<Type>(name);
87 }
88 
89 
90 template<class Type>
93 {
94  return getField<Type>(name);
95 }
96 
97 
98 template<class Type>
101 {
102  return getField<Type>(name);
103 }
104 
105 
106 template<class Type>
109 {
110  tmp<Field<Type>> tfield = getVariableIfAvailable<Type>(name);
111 
112  if (tfield.valid())
113  {
114  return tfield;
115  }
116 
117 
121 
122  const objectRegistry& obr = this->mesh().thisDb();
123 
124  const vfieldType* vfield = obr.findObject<vfieldType>(name);
125  const sfieldType* sfield = obr.findObject<sfieldType>(name);
126  const pfieldType* pfield = obr.findObject<pfieldType>(name);
127 
128  // Local, temporary storage
129  tmp<vfieldType> t_vfield;
130  tmp<sfieldType> t_sfield;
131  tmp<pfieldType> t_pfield;
132 
133  if (searchFiles() && !vfield && !sfield && !pfield)
134  {
135  const word fldType = this->getTypeOfField(name);
136 
137  if (fldType == vfieldType::typeName)
138  {
139  t_vfield = this->readAndRegister<vfieldType>(name, mesh());
140  vfield = t_vfield.get();
141  }
142  else if (fldType == sfieldType::typeName)
143  {
144  t_sfield = this->readAndRegister<sfieldType>(name, mesh());
145  sfield = t_sfield.get();
146  }
147  else if (fldType == pfieldType::typeName)
148  {
149  t_pfield = this->readAndRegister<pfieldType>
150  (
151  name,
153  );
154  pfield = t_pfield.get();
155  }
156  }
157 
158  const label patchIndex = patch_.index();
159 
160  if (vfield)
161  {
162  return tmp<Field<Type>>::New
163  (
164  vfield->boundaryField()[patchIndex]
165  );
166  }
167 
168  if (sfield)
169  {
170  return tmp<Field<Type>>::New
171  (
172  sfield->boundaryField()[patchIndex]
173  );
174  }
175 
176  if (pfield)
177  {
178  return pfield->boundaryField()[patchIndex].patchInternalField();
179  }
180 
181 
183  << "No field '" << name << "' of type "
184  << pTraits<Type>::typeName << nl << nl
185  << vfieldType::typeName << " Fields: "
186  << flatOutput(obr.sortedNames<vfieldType>()) << nl
187  << sfieldType::typeName << " Fields: "
188  << flatOutput(obr.sortedNames<sfieldType>()) << nl
189  << pfieldType::typeName << " Fields: "
190  << flatOutput(obr.sortedNames<pfieldType>()) << nl
191  << exit(FatalError);
192 
193  return tmp<Field<Type>>::New();
194 }
195 
196 
197 template<class Type>
200 (
201  const Field<Type>& field
202 ) const
203 {
204  primitivePatchInterpolation interp(patch_.patch());
205 
206  return interp.pointToFaceInterpolate(field);
207 }
208 
209 
210 template<class Type>
213 (
214  const Field<Type>& field
215 ) const
216 {
217  primitivePatchInterpolation interp(patch_.patch());
218 
219  return interp.faceToPointInterpolate(field);
220 }
221 
222 
223 // ************************************************************************* //
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:61
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:234
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:350
primitivePatchInterpolation.H
Foam::tmp::get
T * get() noexcept
Return pointer without nullptr checking.
Definition: tmp.H:190
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::refPtr::valid
bool valid() const noexcept
True for non-null pointer/reference.
Definition: refPtr.H:160
Foam::expressions::patchExpr::parseDriver::getSurfaceField
tmp< Field< Type > > getSurfaceField(const word &fldName)
Retrieve field (surface field)
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::refPtr::cref
const T & cref() const
Definition: refPtrI.H:187
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::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:217
Foam::PrimitivePatchInterpolation::faceToPointInterpolate
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
Definition: PrimitivePatchInterpolation.C:177
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:381
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:385
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
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
True for non-null pointer/reference.
Definition: tmp.H:175
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:303
Foam::refPtr
A class for managing references or pointers (no reference counting)
Definition: PtrList.H:60