timeVaryingMappedFixedValuePointPatchField.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "Time.H"
31 #include "rawIOField.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
38 (
39  const pointPatch& p,
41 )
42 :
44  fieldTableName_(iF.name()),
45  setAverage_(false),
46  perturb_(0),
47  mapperPtr_(nullptr),
48  sampleTimes_(0),
49  startSampleTime_(-1),
50  startSampledValues_(0),
51  startAverage_(Zero),
52  endSampleTime_(-1),
53  endSampledValues_(0),
54  endAverage_(Zero),
55  offset_()
56 {}
57 
58 
59 template<class Type>
62 (
63  const pointPatch& p,
65  const dictionary& dict
66 )
67 :
69  fieldTableName_(iF.name()),
70  setAverage_(dict.getOrDefault("setAverage", false)),
71  perturb_(dict.getOrDefault("perturb", 1e-5)),
72  mapMethod_
73  (
74  dict.getOrDefault<word>
75  (
76  "mapMethod",
77  "planarInterpolation"
78  )
79  ),
80  mapperPtr_(nullptr),
81  sampleTimes_(0),
82  startSampleTime_(-1),
83  startSampledValues_(0),
84  startAverage_(Zero),
85  endSampleTime_(-1),
86  endSampledValues_(0),
87  endAverage_(Zero),
88  offset_(Function1<Type>::NewIfPresent("offset", dict))
89 {
90  if
91  (
92  mapMethod_ != "planarInterpolation"
93  && mapMethod_ != "nearest"
94  )
95  {
97  << "mapMethod should be one of 'planarInterpolation'"
98  << ", 'nearest'" << exit(FatalIOError);
99  }
100 
101  dict.readIfPresent("fieldTableName", fieldTableName_);
102 
103  if (dict.found("value"))
104  {
106  (
107  Field<Type>("value", dict, p.size())
108  );
109  }
110  else
111  {
112  // Note: use evaluate to do updateCoeffs followed by a reset
113  // of the pointPatchField::updated_ flag. This is
114  // so if first use is in the next time step it retriggers
115  // a new update.
116  pointPatchField<Type>::evaluate(Pstream::commsTypes::blocking);
117  }
118 }
119 
120 
121 template<class Type>
124 (
126  const pointPatch& p,
128  const pointPatchFieldMapper& mapper
129 )
130 :
131  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
132  fieldTableName_(ptf.fieldTableName_),
133  setAverage_(ptf.setAverage_),
134  perturb_(ptf.perturb_),
135  mapMethod_(ptf.mapMethod_),
136  mapperPtr_(nullptr),
137  sampleTimes_(0),
138  startSampleTime_(-1),
139  startSampledValues_(0),
140  startAverage_(Zero),
141  endSampleTime_(-1),
142  endSampledValues_(0),
143  endAverage_(Zero),
144  offset_(ptf.offset_.clone())
145 {}
146 
147 
148 template<class Type>
151 (
153 )
154 :
156  fieldTableName_(ptf.fieldTableName_),
157  setAverage_(ptf.setAverage_),
158  perturb_(ptf.perturb_),
159  mapMethod_(ptf.mapMethod_),
160  mapperPtr_(ptf.mapperPtr_),
161  sampleTimes_(ptf.sampleTimes_),
162  startSampleTime_(ptf.startSampleTime_),
163  startSampledValues_(ptf.startSampledValues_),
164  startAverage_(ptf.startAverage_),
165  endSampleTime_(ptf.endSampleTime_),
166  endSampledValues_(ptf.endSampledValues_),
167  endAverage_(ptf.endAverage_),
168  offset_(ptf.offset_.clone())
169 {}
170 
171 
172 template<class Type>
175 (
178 )
179 :
181  fieldTableName_(ptf.fieldTableName_),
182  setAverage_(ptf.setAverage_),
183  perturb_(ptf.perturb_),
184  mapMethod_(ptf.mapMethod_),
185  mapperPtr_(ptf.mapperPtr_),
186  sampleTimes_(ptf.sampleTimes_),
187  startSampleTime_(ptf.startSampleTime_),
188  startSampledValues_(ptf.startSampledValues_),
189  startAverage_(ptf.startAverage_),
190  endSampleTime_(ptf.endSampleTime_),
191  endSampledValues_(ptf.endSampledValues_),
192  endAverage_(ptf.endAverage_),
193  offset_(ptf.offset_.clone())
194 {}
195 
196 
197 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198 
199 template<class Type>
201 (
202  const pointPatchFieldMapper& m
203 )
204 {
206  if (startSampledValues_.size())
207  {
208  startSampledValues_.autoMap(m);
209  endSampledValues_.autoMap(m);
210  }
211  // Clear interpolator
212  mapperPtr_.clear();
213  startSampleTime_ = -1;
214  endSampleTime_ = -1;
215 }
216 
217 
218 template<class Type>
220 (
221  const pointPatchField<Type>& ptf,
222  const labelList& addr
223 )
224 {
226 
228  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
229 
230  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
231  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
232 
233  // Clear interpolator
234  mapperPtr_.clear();
235  startSampleTime_ = -1;
236  endSampleTime_ = -1;
237 }
238 
239 
240 template<class Type>
242 {
243  const Time& time = this->db().time();
244 
245  // Initialise
246  if (startSampleTime_ == -1 && endSampleTime_ == -1)
247  {
248  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
249 
250  // Read the initial point position
251  pointField meshPts;
252 
253  if (pMesh.pointsInstance() == pMesh.facesInstance())
254  {
255  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
256  }
257  else
258  {
259  // Load points from facesInstance
260  if (debug)
261  {
262  Info<< "Reloading points0 from " << pMesh.facesInstance()
263  << endl;
264  }
265 
267  (
268  IOobject
269  (
270  "points",
271  pMesh.facesInstance(),
272  polyMesh::meshSubDir,
273  pMesh,
274  IOobject::MUST_READ,
275  IOobject::NO_WRITE,
276  false
277  )
278  );
279  meshPts = pointField(points0, this->patch().meshPoints());
280  }
281 
282  // Reread values and interpolate
283  const fileName samplePointsFile
284  (
285  time.caseConstant()
286  /"boundaryData"
287  /this->patch().name()
288  /"points"
289  );
290 
291  IOobject io
292  (
293  samplePointsFile, // absolute path
294  time,
295  IOobject::MUST_READ,
296  IOobject::NO_WRITE,
297  false, // no need to register
298  true // is global object (currently not used)
299  );
300 
301  // Read data
302  const rawIOField<point> samplePoints(io, false);
303 
304  // tbd: run-time selection
305  bool nearestOnly =
306  (
307  !mapMethod_.empty()
308  && mapMethod_ != "planarInterpolation"
309  );
310 
311  // Allocate the interpolator
312  mapperPtr_.reset
313  (
315  (
316  samplePoints,
317  meshPts,
318  perturb_,
319  nearestOnly
320  )
321  );
322 
323  // Read the times for which data is available
324 
325  const fileName samplePointsDir = samplePointsFile.path();
326  sampleTimes_ = Time::findTimes(samplePointsDir);
327 
328  if (debug)
329  {
330  Info<< "timeVaryingMappedFixedValuePointPatchField : In directory "
331  << samplePointsDir << " found times "
332  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
333  << endl;
334  }
335  }
336 
337  // Find current time in sampleTimes
338  label lo = -1;
339  label hi = -1;
340 
341  bool foundTime = mapperPtr_().findTime
342  (
343  sampleTimes_,
344  startSampleTime_,
345  time.value(),
346  lo,
347  hi
348  );
349 
350  if (!foundTime)
351  {
353  << "Cannot find starting sampling values for current time "
354  << time.value() << nl
355  << "Have sampling values for times "
356  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
357  << "In directory "
358  << time.constant()/"boundaryData"/this->patch().name()
359  << "\n on patch " << this->patch().name()
360  << " of field " << fieldTableName_
361  << exit(FatalError);
362  }
363 
364 
365  // Update sampled data fields.
366 
367  if (lo != startSampleTime_)
368  {
369  startSampleTime_ = lo;
370 
371  if (startSampleTime_ == endSampleTime_)
372  {
373  // No need to reread since are end values
374  if (debug)
375  {
376  Pout<< "checkTable : Setting startValues to (already read) "
377  << "boundaryData"
378  /this->patch().name()
379  /sampleTimes_[startSampleTime_].name()
380  << endl;
381  }
382  startSampledValues_ = endSampledValues_;
383  startAverage_ = endAverage_;
384  }
385  else
386  {
387  if (debug)
388  {
389  Pout<< "checkTable : Reading startValues from "
390  << "boundaryData"
391  /this->patch().name()
392  /sampleTimes_[lo].name()
393  << endl;
394  }
395 
396  // Reread values and interpolate
397  const fileName valsFile
398  (
399  time.caseConstant()
400  /"boundaryData"
401  /this->patch().name()
402  /sampleTimes_[startSampleTime_].name()
403  /fieldTableName_
404  );
405 
406  IOobject io
407  (
408  valsFile, // absolute path
409  time,
410  IOobject::MUST_READ,
411  IOobject::NO_WRITE,
412  false, // no need to register
413  true // is global object (currently not used)
414  );
415 
416  const rawIOField<Type> vals(io, setAverage_);
417  if (setAverage_)
418  {
419  startAverage_ = vals.average();
420  }
421 
422  if (vals.size() != mapperPtr_().sourceSize())
423  {
425  << "Number of values (" << vals.size()
426  << ") differs from the number of points ("
427  << mapperPtr_().sourceSize()
428  << ") in file " << valsFile << exit(FatalError);
429  }
430 
431  startSampledValues_ = mapperPtr_().interpolate(vals);
432  }
433  }
434 
435  if (hi != endSampleTime_)
436  {
437  endSampleTime_ = hi;
438 
439  if (endSampleTime_ == -1)
440  {
441  // endTime no longer valid. Might as well clear endValues.
442  if (debug)
443  {
444  Pout<< "checkTable : Clearing endValues" << endl;
445  }
446  endSampledValues_.clear();
447  }
448  else
449  {
450  if (debug)
451  {
452  Pout<< "checkTable : Reading endValues from "
453  << "boundaryData"
454  /this->patch().name()
455  /sampleTimes_[endSampleTime_].name()
456  << endl;
457  }
458 
459  // Reread values and interpolate
460  const fileName valsFile
461  (
462  time.caseConstant()
463  /"boundaryData"
464  /this->patch().name()
465  /sampleTimes_[endSampleTime_].name()
466  /fieldTableName_
467  );
468 
469  IOobject io
470  (
471  valsFile, // absolute path
472  time,
473  IOobject::MUST_READ,
474  IOobject::NO_WRITE,
475  false, // no need to register
476  true // is global object (currently not used)
477  );
478 
479 
480  const rawIOField<Type> vals(io, setAverage_);
481  if (setAverage_)
482  {
483  endAverage_ = vals.average();
484  }
485 
486  if (vals.size() != mapperPtr_().sourceSize())
487  {
489  << "Number of values (" << vals.size()
490  << ") differs from the number of points ("
491  << mapperPtr_().sourceSize()
492  << ") in file " << valsFile << exit(FatalError);
493  }
494 
495  endSampledValues_ = mapperPtr_().interpolate(vals);
496  }
497  }
498 }
499 
500 
501 template<class Type>
503 {
504  if (this->updated())
505  {
506  return;
507  }
508 
509  checkTable();
510 
511  // Interpolate between the sampled data
512 
513  Type wantedAverage;
514 
515  if (endSampleTime_ == -1)
516  {
517  // only start value
518  if (debug)
519  {
520  Pout<< "updateCoeffs : Sampled, non-interpolated values"
521  << " from start time:"
522  << sampleTimes_[startSampleTime_].name() << nl;
523  }
524 
525  this->operator==(startSampledValues_);
526  wantedAverage = startAverage_;
527  }
528  else
529  {
530  scalar start = sampleTimes_[startSampleTime_].value();
531  scalar end = sampleTimes_[endSampleTime_].value();
532 
533  scalar s = (this->db().time().value()-start)/(end-start);
534 
535  if (debug)
536  {
537  Pout<< "updateCoeffs : Sampled, interpolated values"
538  << " between start time:"
539  << sampleTimes_[startSampleTime_].name()
540  << " and end time:" << sampleTimes_[endSampleTime_].name()
541  << " with weight:" << s << endl;
542  }
543 
544  this->operator==((1-s)*startSampledValues_ + s*endSampledValues_);
545  wantedAverage = (1-s)*startAverage_ + s*endAverage_;
546  }
547 
548  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
549  // offsetting.
550  if (setAverage_)
551  {
552  const Field<Type>& fld = *this;
553 
554  Type averagePsi = gAverage(fld);
555 
556  if (debug)
557  {
558  Pout<< "updateCoeffs :"
559  << " actual average:" << averagePsi
560  << " wanted average:" << wantedAverage
561  << endl;
562  }
563 
564  if (mag(averagePsi) < VSMALL)
565  {
566  // Field too small to scale. Offset instead.
567  const Type offset = wantedAverage - averagePsi;
568  if (debug)
569  {
570  Pout<< "updateCoeffs :"
571  << " offsetting with:" << offset << endl;
572  }
573  this->operator==(fld+offset);
574  }
575  else
576  {
577  const scalar scale = mag(wantedAverage)/mag(averagePsi);
578 
579  if (debug)
580  {
581  Pout<< "updateCoeffs :"
582  << " scaling with:" << scale << endl;
583  }
584  this->operator==(scale*fld);
585  }
586  }
587 
588  // Apply offset to mapped values
589  if (offset_)
590  {
591  const scalar t = this->db().time().timeOutputValue();
592  this->operator==(*this + offset_->value(t));
593  }
594 
595  if (debug)
596  {
597  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
598  << " max:" << gMax(*this)
599  << " avg:" << gAverage(*this) << endl;
600  }
601 
603 }
604 
605 
606 template<class Type>
608 (
609  Ostream& os
610 ) const
611 {
613 
614  os.writeEntryIfDifferent("setAverage", Switch(false), setAverage_);
615  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
616 
618  (
619  "fieldTable",
620  this->internalField().name(),
621  fieldTableName_
622  );
623 
625  (
626  "mapMethod",
627  "planarInterpolation",
628  mapMethod_
629  );
630 
631  if (offset_)
632  {
633  offset_->writeData(os);
634  }
635 }
636 
637 
638 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:244
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
rawIOField.H
s
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:186
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::timeVaryingMappedFixedValuePointPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: timeVaryingMappedFixedValuePointPatchField.C:502
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:209
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:852
Foam::pointToPointPlanarInterpolation
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
Definition: pointToPointPlanarInterpolation.H:53
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:58
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:60
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:48
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:846
Foam::fixedValuePointPatchField
A FixedValue boundary condition for pointField.
Definition: fixedValuePointPatchField.H:51
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::rawIOField
Like IOField but falls back to raw IFstream if no header found. Optionally reads average value....
Definition: rawIOField.H:52
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::timeVaryingMappedFixedValuePointPatchField
A time-varying form of a mapped fixed value boundary condition.
Definition: timeVaryingMappedFixedValuePointPatchField.H:57
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:107
Foam::timeVaryingMappedFixedValuePointPatchField::rmap
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
Definition: timeVaryingMappedFixedValuePointPatchField.C:220
Foam::timeVaryingMappedFixedValuePointPatchField::write
virtual void write(Ostream &) const
Write.
Definition: timeVaryingMappedFixedValuePointPatchField.C:608
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
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::timeVaryingMappedFixedValuePointPatchField::timeVaryingMappedFixedValuePointPatchField
timeVaryingMappedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
Definition: timeVaryingMappedFixedValuePointPatchField.C:38
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Time.H
timeVaryingMappedFixedValuePointPatchField.H
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
Foam::rawIOField::average
const Type & average() const
Definition: rawIOField.H:83
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::TimePaths::caseConstant
fileName caseConstant() const
Definition: TimePathsI.H:100
Foam::timeVaryingMappedFixedValuePointPatchField::autoMap
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: timeVaryingMappedFixedValuePointPatchField.C:201
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::timeVaryingMappedFixedValuePointPatchField::checkTable
void checkTable()
Find boundary data inbetween current time and interpolate.
Definition: timeVaryingMappedFixedValuePointPatchField.C:241
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54