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_
89  (
90  Function1<Type>::NewIfPresent("offset", dict, word::null, &this->db())
91  )
92 {
93  if
94  (
95  mapMethod_ != "planarInterpolation"
96  && mapMethod_ != "nearest"
97  )
98  {
100  << "mapMethod should be one of 'planarInterpolation'"
101  << ", 'nearest'" << exit(FatalIOError);
102  }
103 
104  dict.readIfPresent("fieldTableName", fieldTableName_);
105 
106  if (dict.found("value"))
107  {
109  (
110  Field<Type>("value", dict, p.size())
111  );
112  }
113  else
114  {
115  // Note: use evaluate to do updateCoeffs followed by a reset
116  // of the pointPatchField::updated_ flag. This is
117  // so if first use is in the next time step it retriggers
118  // a new update.
119  pointPatchField<Type>::evaluate(Pstream::commsTypes::blocking);
120  }
121 }
122 
123 
124 template<class Type>
127 (
129  const pointPatch& p,
131  const pointPatchFieldMapper& mapper
132 )
133 :
134  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
135  fieldTableName_(ptf.fieldTableName_),
136  setAverage_(ptf.setAverage_),
137  perturb_(ptf.perturb_),
138  mapMethod_(ptf.mapMethod_),
139  mapperPtr_(nullptr),
140  sampleTimes_(0),
141  startSampleTime_(-1),
142  startSampledValues_(0),
143  startAverage_(Zero),
144  endSampleTime_(-1),
145  endSampledValues_(0),
146  endAverage_(Zero),
147  offset_(ptf.offset_.clone())
148 {}
149 
150 
151 template<class Type>
154 (
156 )
157 :
159  fieldTableName_(ptf.fieldTableName_),
160  setAverage_(ptf.setAverage_),
161  perturb_(ptf.perturb_),
162  mapMethod_(ptf.mapMethod_),
163  mapperPtr_(ptf.mapperPtr_),
164  sampleTimes_(ptf.sampleTimes_),
165  startSampleTime_(ptf.startSampleTime_),
166  startSampledValues_(ptf.startSampledValues_),
167  startAverage_(ptf.startAverage_),
168  endSampleTime_(ptf.endSampleTime_),
169  endSampledValues_(ptf.endSampledValues_),
170  endAverage_(ptf.endAverage_),
171  offset_(ptf.offset_.clone())
172 {}
173 
174 
175 template<class Type>
178 (
181 )
182 :
184  fieldTableName_(ptf.fieldTableName_),
185  setAverage_(ptf.setAverage_),
186  perturb_(ptf.perturb_),
187  mapMethod_(ptf.mapMethod_),
188  mapperPtr_(ptf.mapperPtr_),
189  sampleTimes_(ptf.sampleTimes_),
190  startSampleTime_(ptf.startSampleTime_),
191  startSampledValues_(ptf.startSampledValues_),
192  startAverage_(ptf.startAverage_),
193  endSampleTime_(ptf.endSampleTime_),
194  endSampledValues_(ptf.endSampledValues_),
195  endAverage_(ptf.endAverage_),
196  offset_(ptf.offset_.clone())
197 {}
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
202 template<class Type>
204 (
205  const pointPatchFieldMapper& m
206 )
207 {
209  if (startSampledValues_.size())
210  {
211  startSampledValues_.autoMap(m);
212  endSampledValues_.autoMap(m);
213  }
214  // Clear interpolator
215  mapperPtr_.clear();
216  startSampleTime_ = -1;
217  endSampleTime_ = -1;
218 }
219 
220 
221 template<class Type>
223 (
224  const pointPatchField<Type>& ptf,
225  const labelList& addr
226 )
227 {
229 
231  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
232 
233  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
234  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
235 
236  // Clear interpolator
237  mapperPtr_.clear();
238  startSampleTime_ = -1;
239  endSampleTime_ = -1;
240 }
241 
242 
243 template<class Type>
245 {
246  const Time& time = this->db().time();
247 
248  // Initialise
249  if (startSampleTime_ == -1 && endSampleTime_ == -1)
250  {
251  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
252 
253  // Read the initial point position
254  pointField meshPts;
255 
256  if (pMesh.pointsInstance() == pMesh.facesInstance())
257  {
258  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
259  }
260  else
261  {
262  // Load points from facesInstance
263  if (debug)
264  {
265  Info<< "Reloading points0 from " << pMesh.facesInstance()
266  << endl;
267  }
268 
270  (
271  IOobject
272  (
273  "points",
274  pMesh.facesInstance(),
275  polyMesh::meshSubDir,
276  pMesh,
277  IOobject::MUST_READ,
278  IOobject::NO_WRITE,
279  false
280  )
281  );
282  meshPts = pointField(points0, this->patch().meshPoints());
283  }
284 
285  // Reread values and interpolate
286  const fileName samplePointsFile
287  (
288  time.caseConstant()
289  /"boundaryData"
290  /this->patch().name()
291  /"points"
292  );
293 
294  IOobject io
295  (
296  samplePointsFile, // absolute path
297  time,
298  IOobject::MUST_READ,
299  IOobject::NO_WRITE,
300  false, // no need to register
301  true // is global object (currently not used)
302  );
303 
304  // Read data
305  const rawIOField<point> samplePoints(io, false);
306 
307  // tbd: run-time selection
308  bool nearestOnly =
309  (
310  !mapMethod_.empty()
311  && mapMethod_ != "planarInterpolation"
312  );
313 
314  // Allocate the interpolator
315  mapperPtr_.reset
316  (
318  (
319  samplePoints,
320  meshPts,
321  perturb_,
322  nearestOnly
323  )
324  );
325 
326  // Read the times for which data is available
327 
328  const fileName samplePointsDir = samplePointsFile.path();
329  sampleTimes_ = Time::findTimes(samplePointsDir);
330 
331  if (debug)
332  {
333  Info<< "timeVaryingMappedFixedValuePointPatchField : In directory "
334  << samplePointsDir << " found times "
335  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
336  << endl;
337  }
338  }
339 
340  // Find current time in sampleTimes
341  label lo = -1;
342  label hi = -1;
343 
344  bool foundTime = mapperPtr_().findTime
345  (
346  sampleTimes_,
347  startSampleTime_,
348  time.value(),
349  lo,
350  hi
351  );
352 
353  if (!foundTime)
354  {
356  << "Cannot find starting sampling values for current time "
357  << time.value() << nl
358  << "Have sampling values for times "
359  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
360  << "In directory "
361  << time.constant()/"boundaryData"/this->patch().name()
362  << "\n on patch " << this->patch().name()
363  << " of field " << fieldTableName_
364  << exit(FatalError);
365  }
366 
367 
368  // Update sampled data fields.
369 
370  if (lo != startSampleTime_)
371  {
372  startSampleTime_ = lo;
373 
374  if (startSampleTime_ == endSampleTime_)
375  {
376  // No need to reread since are end values
377  if (debug)
378  {
379  Pout<< "checkTable : Setting startValues to (already read) "
380  << "boundaryData"
381  /this->patch().name()
382  /sampleTimes_[startSampleTime_].name()
383  << endl;
384  }
385  startSampledValues_ = endSampledValues_;
386  startAverage_ = endAverage_;
387  }
388  else
389  {
390  if (debug)
391  {
392  Pout<< "checkTable : Reading startValues from "
393  << "boundaryData"
394  /this->patch().name()
395  /sampleTimes_[lo].name()
396  << endl;
397  }
398 
399  // Reread values and interpolate
400  const fileName valsFile
401  (
402  time.caseConstant()
403  /"boundaryData"
404  /this->patch().name()
405  /sampleTimes_[startSampleTime_].name()
406  /fieldTableName_
407  );
408 
409  IOobject io
410  (
411  valsFile, // absolute path
412  time,
413  IOobject::MUST_READ,
414  IOobject::NO_WRITE,
415  false, // no need to register
416  true // is global object (currently not used)
417  );
418 
419  const rawIOField<Type> vals(io, setAverage_);
420  if (setAverage_)
421  {
422  startAverage_ = vals.average();
423  }
424 
425  if (vals.size() != mapperPtr_().sourceSize())
426  {
428  << "Number of values (" << vals.size()
429  << ") differs from the number of points ("
430  << mapperPtr_().sourceSize()
431  << ") in file " << valsFile << exit(FatalError);
432  }
433 
434  startSampledValues_ = mapperPtr_().interpolate(vals);
435  }
436  }
437 
438  if (hi != endSampleTime_)
439  {
440  endSampleTime_ = hi;
441 
442  if (endSampleTime_ == -1)
443  {
444  // endTime no longer valid. Might as well clear endValues.
445  if (debug)
446  {
447  Pout<< "checkTable : Clearing endValues" << endl;
448  }
449  endSampledValues_.clear();
450  }
451  else
452  {
453  if (debug)
454  {
455  Pout<< "checkTable : Reading endValues from "
456  << "boundaryData"
457  /this->patch().name()
458  /sampleTimes_[endSampleTime_].name()
459  << endl;
460  }
461 
462  // Reread values and interpolate
463  const fileName valsFile
464  (
465  time.caseConstant()
466  /"boundaryData"
467  /this->patch().name()
468  /sampleTimes_[endSampleTime_].name()
469  /fieldTableName_
470  );
471 
472  IOobject io
473  (
474  valsFile, // absolute path
475  time,
476  IOobject::MUST_READ,
477  IOobject::NO_WRITE,
478  false, // no need to register
479  true // is global object (currently not used)
480  );
481 
482 
483  const rawIOField<Type> vals(io, setAverage_);
484  if (setAverage_)
485  {
486  endAverage_ = vals.average();
487  }
488 
489  if (vals.size() != mapperPtr_().sourceSize())
490  {
492  << "Number of values (" << vals.size()
493  << ") differs from the number of points ("
494  << mapperPtr_().sourceSize()
495  << ") in file " << valsFile << exit(FatalError);
496  }
497 
498  endSampledValues_ = mapperPtr_().interpolate(vals);
499  }
500  }
501 }
502 
503 
504 template<class Type>
506 {
507  if (this->updated())
508  {
509  return;
510  }
511 
512  checkTable();
513 
514  // Interpolate between the sampled data
515 
516  Type wantedAverage;
517 
518  if (endSampleTime_ == -1)
519  {
520  // only start value
521  if (debug)
522  {
523  Pout<< "updateCoeffs : Sampled, non-interpolated values"
524  << " from start time:"
525  << sampleTimes_[startSampleTime_].name() << nl;
526  }
527 
528  this->operator==(startSampledValues_);
529  wantedAverage = startAverage_;
530  }
531  else
532  {
533  scalar start = sampleTimes_[startSampleTime_].value();
534  scalar end = sampleTimes_[endSampleTime_].value();
535 
536  scalar s = (this->db().time().value()-start)/(end-start);
537 
538  if (debug)
539  {
540  Pout<< "updateCoeffs : Sampled, interpolated values"
541  << " between start time:"
542  << sampleTimes_[startSampleTime_].name()
543  << " and end time:" << sampleTimes_[endSampleTime_].name()
544  << " with weight:" << s << endl;
545  }
546 
547  this->operator==((1-s)*startSampledValues_ + s*endSampledValues_);
548  wantedAverage = (1-s)*startAverage_ + s*endAverage_;
549  }
550 
551  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
552  // offsetting.
553  if (setAverage_)
554  {
555  const Field<Type>& fld = *this;
556 
557  Type averagePsi = gAverage(fld);
558 
559  if (debug)
560  {
561  Pout<< "updateCoeffs :"
562  << " actual average:" << averagePsi
563  << " wanted average:" << wantedAverage
564  << endl;
565  }
566 
567  if (mag(averagePsi) < VSMALL)
568  {
569  // Field too small to scale. Offset instead.
570  const Type offset = wantedAverage - averagePsi;
571  if (debug)
572  {
573  Pout<< "updateCoeffs :"
574  << " offsetting with:" << offset << endl;
575  }
576  this->operator==(fld+offset);
577  }
578  else
579  {
580  const scalar scale = mag(wantedAverage)/mag(averagePsi);
581 
582  if (debug)
583  {
584  Pout<< "updateCoeffs :"
585  << " scaling with:" << scale << endl;
586  }
587  this->operator==(scale*fld);
588  }
589  }
590 
591  // Apply offset to mapped values
592  if (offset_)
593  {
594  const scalar t = this->db().time().timeOutputValue();
595  this->operator==(*this + offset_->value(t));
596  }
597 
598  if (debug)
599  {
600  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
601  << " max:" << gMax(*this)
602  << " avg:" << gAverage(*this) << endl;
603  }
604 
606 }
607 
608 
609 template<class Type>
611 (
612  Ostream& os
613 ) const
614 {
616 
617  os.writeEntryIfDifferent("setAverage", Switch(false), setAverage_);
618  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
619 
620  os.writeEntryIfDifferent
621  (
622  "fieldTable",
623  this->internalField().name(),
624  fieldTableName_
625  );
626 
627  os.writeEntryIfDifferent<word>
628  (
629  "mapMethod",
630  "planarInterpolation",
631  mapMethod_
632  );
633 
634  if (offset_)
635  {
636  offset_->writeData(os);
637  }
638 }
639 
640 
641 // ************************************************************************* //
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::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
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:176
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:505
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
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:369
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: propellerInfo.H:291
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 (stdout output on master, null elsewhere)
Foam::timeVaryingMappedFixedValuePointPatchField
A time-varying form of a mapped fixed value boundary condition.
Definition: timeVaryingMappedFixedValuePointPatchField.H:57
Foam::OSstream::name
virtual const fileName & name() const
Get 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:223
Foam::timeVaryingMappedFixedValuePointPatchField::write
virtual void write(Ostream &) const
Write.
Definition: timeVaryingMappedFixedValuePointPatchField.C:611
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:123
os
OBJstream os(runTime.globalPath()/outputName)
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:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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:108
Foam::timeVaryingMappedFixedValuePointPatchField::autoMap
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: timeVaryingMappedFixedValuePointPatchField.C:204
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:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
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:244
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37