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-2022 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 (returnReduceAnd(vals.size() == len))
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  const objectRegistry& obr = this->mesh().thisDb();
118  const label patchIndex = patch_.index();
119 
120 
121  // Local, temporary storage and/or lookup values
122  bool found = false;
123  tmp<VolumeField<Type>> vfield;
124  tmp<SurfaceField<Type>> sfield;
125  tmp<PointField<Type>> pfield;
126 
127  for (int checki = 0; !found && checki < 2; ++checki)
128  {
129  // Check 0: object context (first)
130  // Check 1: regular objectRegistry
131  const regIOobject* ioptr =
132  (
133  (checki == 0)
135  : obr.cfindIOobject(name)
136  );
137  if (!ioptr) continue;
138 
139  if (!found)
140  {
141  vfield.cref(dynamic_cast<const VolumeField<Type>*>(ioptr));
142  found = vfield.valid();
143  }
144  if (!found)
145  {
146  sfield.cref(dynamic_cast<const SurfaceField<Type>*>(ioptr));
147  found = sfield.valid();
148  }
149  if (!found)
150  {
151  pfield.cref(dynamic_cast<const PointField<Type>*>(ioptr));
152  found = pfield.valid();
153  }
154  }
155 
156 
157  // Finally, search files if necessary (and permitted)
158  if (!found && searchFiles())
159  {
160  const word fldType = this->getTypeOfField(name);
161 
162  if (fldType == VolumeField<Type>::typeName)
163  {
164  vfield = this->readAndRegister<VolumeField<Type>>(name, mesh());
165  }
166  else if (fldType == SurfaceField<Type>::typeName)
167  {
168  sfield = this->readAndRegister<SurfaceField<Type>>(name, mesh());
169  }
170  else if (fldType == PointField<Type>::typeName)
171  {
172  pfield = this->readAndRegister<PointField<Type>>
173  (
174  name,
176  );
177  }
178  }
179 
180 
181  if (vfield.valid())
182  {
183  return tmp<Field<Type>>::New
184  (
185  vfield().boundaryField()[patchIndex]
186  );
187  }
188  if (sfield.valid())
189  {
190  return tmp<Field<Type>>::New
191  (
192  sfield().boundaryField()[patchIndex]
193  );
194  }
195  if (pfield.valid())
196  {
197  return pfield().boundaryField()[patchIndex].patchInternalField();
198  }
199 
200 
202  << "No field '" << name << "' of type "
203  << pTraits<Type>::typeName << nl << nl;
204 
205  FatalError
206  << VolumeField<Type>::typeName << " Fields: "
207  << flatOutput(obr.sortedNames<VolumeField<Type>>()) << nl;
208 
209  FatalError
210  << SurfaceField<Type>::typeName << " Fields: "
211  << flatOutput(obr.sortedNames<SurfaceField<Type>>()) << nl;
212 
213  FatalError
214  << PointField<Type>::typeName << " Fields: "
215  << flatOutput(obr.sortedNames<PointField<Type>>()) << nl;
216 
217  FatalError
218  << exit(FatalError);
219 
221 }
222 
223 
224 template<class Type>
227 (
228  const word& name
229 )
230 {
231  tmp<Field<Type>> tfield = getVariableIfAvailable<Type>(name);
232 
233  if (tfield.valid())
234  {
235  return tfield;
236  }
237 
238  const objectRegistry& obr = this->mesh().thisDb();
239  const label patchIndex = patch_.index();
240 
241 
242  // Local, temporary storage and/or lookup values
243  bool found = false;
244  tmp<VolumeField<Type>> vfield;
245  tmp<PointField<Type>> pfield;
246 
247  for (int checki = 0; !found && checki < 2; ++checki)
248  {
249  // Check 0: object context (first)
250  // Check 1: regular objectRegistry
251  const regIOobject* ioptr =
252  (
253  (checki == 0)
255  : obr.cfindIOobject(name)
256  );
257  if (!ioptr) continue;
258 
259  if (!found)
260  {
261  vfield.cref(dynamic_cast<const VolumeField<Type>*>(ioptr));
262  found = vfield.valid();
263  }
264  if (!found)
265  {
266  pfield.cref(dynamic_cast<const PointField<Type>*>(ioptr));
267  found = pfield.valid();
268  }
269  }
270 
271 
272  // Finally, search files if necessary (and permitted)
273  if (!found && searchFiles())
274  {
275  const word fldType = this->getTypeOfField(name);
276 
277  if (fldType == VolumeField<Type>::typeName)
278  {
279  vfield = this->readAndRegister<VolumeField<Type>>(name, mesh());
280  }
281  else if (fldType == PointField<Type>::typeName)
282  {
283  pfield = this->readAndRegister<PointField<Type>>
284  (
285  name,
287  );
288  }
289  }
290 
291 
292  if (vfield.valid())
293  {
294  return vfield().boundaryField()[patchIndex].patchInternalField();
295  }
296  if (pfield.valid())
297  {
298  return pfield().boundaryField()[patchIndex].patchInternalField();
299  }
300 
301 
303  << "No field '" << name << "' of type "
304  << pTraits<Type>::typeName << nl << nl;
305 
306  FatalError
307  << VolumeField<Type>::typeName << " Fields: "
308  << flatOutput(obr.sortedNames<VolumeField<Type>>()) << nl;
309 
310  FatalError
311  << PointField<Type>::typeName << " Fields: "
312  << flatOutput(obr.sortedNames<PointField<Type>>()) << nl;
313 
314  FatalError
315  << exit(FatalError);
316 
318 }
319 
320 
321 template<class Type>
324 (
325  const word& name
326 )
327 {
328  tmp<Field<Type>> tfield = getVariableIfAvailable<Type>(name);
329 
330  if (tfield.valid())
331  {
332  return tfield;
333  }
334 
335  const objectRegistry& obr = this->mesh().thisDb();
336  const label patchIndex = patch_.index();
337 
338 
339  // Local, temporary storage and/or lookup values
340  bool found = false;
341  tmp<VolumeField<Type>> vfield;
342 
343  for (int checki = 0; !found && checki < 2; ++checki)
344  {
345  // Check 0: object context (first)
346  // Check 1: regular objectRegistry
347  const regIOobject* ioptr =
348  (
349  (checki == 0)
351  : obr.cfindIOobject(name)
352  );
353  if (!ioptr) continue;
354 
355  if (!found)
356  {
357  vfield.cref(dynamic_cast<const VolumeField<Type>*>(ioptr));
358  found = vfield.valid();
359  }
360  }
361 
362 
363  // Finally, search files if necessary (and permitted)
364  if (!found && searchFiles())
365  {
366  const word fldType = this->getTypeOfField(name);
367 
368  if (fldType == VolumeField<Type>::typeName)
369  {
370  vfield = this->readAndRegister<VolumeField<Type>>(name, mesh());
371  }
372  }
373 
374 
375  if (vfield.valid())
376  {
377  return vfield().boundaryField()[patchIndex].patchNeighbourField();
378  }
379 
380 
382  << "No field '" << name << "' of type "
383  << pTraits<Type>::typeName << nl << nl;
384 
385  FatalError
386  << VolumeField<Type>::typeName << " Fields: "
387  << flatOutput(obr.sortedNames<VolumeField<Type>>()) << nl;
388 
389  FatalError
390  << exit(FatalError);
391 
393 }
394 
395 
396 template<class Type>
399 (
400  const word& name
401 )
402 {
403  tmp<Field<Type>> tfield = getVariableIfAvailable<Type>(name);
404 
405  if (tfield.valid())
406  {
407  return tfield;
408  }
409 
410  const objectRegistry& obr = this->mesh().thisDb();
411  const label patchIndex = patch_.index();
412 
413 
414  // Local, temporary storage and/or lookup values
415  bool found = false;
416  tmp<VolumeField<Type>> vfield;
417 
418  for (int checki = 0; !found && checki < 2; ++checki)
419  {
420  // Check 0: object context (first)
421  // Check 1: regular objectRegistry
422  const regIOobject* ioptr =
423  (
424  (checki == 0)
426  : obr.cfindIOobject(name)
427  );
428  if (!ioptr) continue;
429 
430  if (!found)
431  {
432  vfield.cref(dynamic_cast<const VolumeField<Type>*>(ioptr));
433  found = vfield.valid();
434  }
435  }
436 
437 
438  // Finally, search files if necessary (and permitted)
439  if (!found && searchFiles())
440  {
441  const word fldType = this->getTypeOfField(name);
442 
443  if (fldType == VolumeField<Type>::typeName)
444  {
445  vfield = this->readAndRegister<VolumeField<Type>>(name, mesh());
446  }
447  }
448 
449 
450  if (vfield.valid())
451  {
452  return vfield().boundaryField()[patchIndex].snGrad();
453  }
454 
455 
457  << "No field '" << name << "' of type "
458  << pTraits<Type>::typeName << nl << nl;
459 
460  FatalError
461  << VolumeField<Type>::typeName << " Fields: "
462  << flatOutput(obr.sortedNames<VolumeField<Type>>()) << nl;
463 
464  FatalError
465  << exit(FatalError);
466 
468 }
469 
470 
471 template<class Type>
474 (
475  const Field<Type>& field
476 ) const
477 {
478  primitivePatchInterpolation interp(patch_.patch());
479 
481 }
482 
483 
484 template<class Type>
487 (
488  const Field<Type>& field
489 ) const
490 {
491  primitivePatchInterpolation interp(patch_.patch());
492 
493  return interp.pointToFaceInterpolate(field);
494 }
495 
496 
497 // ************************************************************************* //
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:472
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
tmp< Field< Type > > getField(const word &fldName)
Return named field.
rDeltaTY field()
tmp< Field< Type > > getPointField(const word &fldName)
Retrieve field (point field)
tmp< Field< Type > > pointToFaceInterpolate(const Field< Type > &pf) const
Interpolate from points to faces.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
tmp< Field< Type > > faceToPoint(const Field< Type > &field) const
Interpolate face to point.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
tmp< Field< Type > > getSurfaceField(const word &fldName)
Retrieve field (surface field)
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool hasVariable(const word &name) const
True if named variable exists.
Definition: fvExprDriverI.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: refPtr.H:512
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: tmpI.H:221
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const regIOobject * cfindContextIOobject(const word &name) const
Find named context field, if it exists.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
static autoPtr< fvExprDriver > New(const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected value driver.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: refPtrI.H:216
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
tmp< Field< Type > > getVariableIfAvailable(const word &fldName) const
Retrieve variable as field if possible.
tmp< Field< Type > > patchNeighbourField(const word &fldName)
Return patchField on the opposite patch of a coupled patch.
tmp< Field< Type > > patchInternalField(const word &fldName)
Return internal field next to patch.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
tmp< Field< Type > > pointToFace(const Field< Type > &field) const
Interpolate point to face values.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
wordList sortedNames() const
The sorted names of all objects.
virtual exprResult & variable(const word &name)
Non-const access to the named variable (sub-classes only)
Definition: fvExprDriverI.H:52
virtual label size() const
The natural field size for the expression.
tmp< Field< Type > > getVolField(const word &fldName)
Retrieve field (vol field)
const word typeName("volScalarField::Internal")
virtual label pointSize() const
The point field size for the expression.
tmp< Field< Type > > patchNormalField(const word &fldName)
Return surface normal field (snGrad)
#define WarningInFunction
Report a warning using Foam::Warning.
Type gAverage(const FieldField< Field, Type > &f)
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
const exprResult & lookupGlobal(const word &name) const
Return the global variable if available or a null result.
Definition: fvExprDriver.C:695
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
bool found
const regIOobject * cfindIOobject(const word &name, const bool recursive=false) const
Return const pointer to the regIOobject.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225