MappedFile.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) 2018-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 
28 #include "polyMesh.H"
29 #include "rawIOField.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const polyPatch& pp,
37  const word& type,
38  const word& entryName,
39  const dictionary& dict,
40  const bool faceValues
41 )
42 :
43  PatchFunction1<Type>(pp, entryName, dict, faceValues),
44  dictConstructed_(true),
45  fieldTableName_(entryName),
46  setAverage_(dict.getOrDefault("setAverage", false)),
47  perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
48  pointsName_(dict.getOrDefault<word>("points", "points")),
49  mapMethod_
50  (
51  dict.getOrDefault<word>
52  (
53  "mapMethod",
54  "planarInterpolation"
55  )
56  ),
57  mapperPtr_(nullptr),
58  sampleTimes_(0),
59  startSampleTime_(-1),
60  startSampledValues_(0),
61  startAverage_(Zero),
62  endSampleTime_(-1),
63  endSampledValues_(0),
64  endAverage_(Zero),
65  offset_(nullptr)
66 {
67  if (dict.found("offset"))
68  {
69  offset_ = Function1<Type>::New("offset", dict);
70  }
71 
72  if
73  (
74  mapMethod_ != "planarInterpolation"
75  && mapMethod_ != "nearest"
76  )
77  {
79  << "mapMethod should be one of 'planarInterpolation'"
80  << ", 'nearest'" << exit(FatalIOError);
81  }
82 
83  dict.readIfPresent("fieldTable", fieldTableName_);
84 }
85 
86 
87 template<class Type>
89 (
90  const polyPatch& pp,
91  const word& entryName,
92  const dictionary& dict,
93  const word& fieldTableName,
94  const bool faceValues
95 )
96 :
97  PatchFunction1<Type>(pp, entryName, dict, faceValues),
98  dictConstructed_(false),
99  fieldTableName_(fieldTableName),
100  setAverage_(dict.getOrDefault("setAverage", false)),
101  perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
102  pointsName_(dict.getOrDefault<word>("points", "points")),
103  mapMethod_
104  (
105  dict.getOrDefault<word>
106  (
107  "mapMethod",
108  "planarInterpolation"
109  )
110  ),
111  mapperPtr_(nullptr),
112  sampleTimes_(0),
113  startSampleTime_(-1),
114  startSampledValues_(0),
115  startAverage_(Zero),
116  endSampleTime_(-1),
117  endSampledValues_(0),
118  endAverage_(Zero),
119  offset_(nullptr)
120 {
121  if (dict.found("offset"))
122  {
123  offset_ = Function1<Type>::New("offset", dict);
124  }
125 
126  if
127  (
128  mapMethod_ != "planarInterpolation"
129  && mapMethod_ != "nearest"
130  )
131  {
133  << "mapMethod should be one of 'planarInterpolation'"
134  << ", 'nearest'" << exit(FatalIOError);
135  }
136 }
137 
138 
139 template<class Type>
141 (
142  const MappedFile<Type>& ut
143 )
144 :
146  dictConstructed_(ut.dictConstructed_),
147  fieldTableName_(ut.fieldTableName_),
148  setAverage_(ut.setAverage_),
149  perturb_(ut.perturb_),
150  pointsName_(ut.pointsName_),
151  mapMethod_(ut.mapMethod_),
152  mapperPtr_(ut.mapperPtr_.clone()),
153  sampleTimes_(ut.sampleTimes_),
154  startSampleTime_(ut.startSampleTime_),
155  startSampledValues_(ut.startSampledValues_),
156  startAverage_(ut.startAverage_),
157  endSampleTime_(ut.endSampleTime_),
158  endSampledValues_(ut.endSampledValues_),
159  endAverage_(ut.endAverage_),
160  offset_(ut.offset_.clone())
161 {}
162 
163 
164 template<class Type>
166 (
167  const MappedFile<Type>& ut,
168  const polyPatch& pp
169 )
170 :
171  PatchFunction1<Type>(ut, pp),
172  dictConstructed_(ut.dictConstructed_),
173  fieldTableName_(ut.fieldTableName_),
174  setAverage_(ut.setAverage_),
175  perturb_(ut.perturb_),
176  pointsName_(ut.pointsName_),
177  mapMethod_(ut.mapMethod_),
178  mapperPtr_(ut.mapperPtr_.clone()),
179  sampleTimes_(ut.sampleTimes_),
180  startSampleTime_(ut.startSampleTime_),
181  startSampledValues_(ut.startSampledValues_),
182  startAverage_(ut.startAverage_),
183  endSampleTime_(ut.endSampleTime_),
184  endSampledValues_(ut.endSampledValues_),
185  endAverage_(ut.endAverage_),
186  offset_(ut.offset_.clone())
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
192 template<class Type>
194 (
195  const FieldMapper& mapper
196 )
197 {
199 
200  if (startSampledValues_.size())
201  {
202  startSampledValues_.autoMap(mapper);
203  endSampledValues_.autoMap(mapper);
204  }
205  // Clear interpolator
206  mapperPtr_.clear();
207  startSampleTime_ = -1;
208  endSampleTime_ = -1;
209 }
210 
211 
212 template<class Type>
214 (
215  const PatchFunction1<Type>& pf1,
216  const labelList& addr
217 )
218 {
219  PatchFunction1<Type>::rmap(pf1, addr);
220 
222  refCast<const PatchFunction1Types::MappedFile<Type>>(pf1);
223 
224  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
225  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
226 
227  // Clear interpolator
228  mapperPtr_.clear();
229  startSampleTime_ = -1;
230  endSampleTime_ = -1;
231 }
232 
233 
234 template<class Type>
236 (
237  const scalar t
238 ) const
239 {
240  const polyMesh& mesh = this->patch_.boundaryMesh().mesh();
241  const Time& time = mesh.time();
242 
243  // Initialise
244  if (!mapperPtr_)
245  {
246  // Reread values and interpolate
247  const fileName samplePointsFile
248  (
249  time.globalPath()
250  /time.constant()
251  /"boundaryData"
252  /this->patch_.name()
253  /pointsName_
254  );
255 
256  IOobject io
257  (
258  samplePointsFile, // absolute path
259  time,
260  IOobject::MUST_READ,
261  IOobject::NO_WRITE,
262  false, // no need to register
263  true // is global object (currently not used)
264  );
265 
266  // Read data
267  const rawIOField<point> samplePoints(io, false);
268 
269  DebugInfo
270  << "Read " << samplePoints.size() << " sample points from "
271  << samplePointsFile << endl;
272 
273 
274  // tbd: run-time selection
275  bool nearestOnly =
276  (
277  !mapMethod_.empty()
278  && mapMethod_ != "planarInterpolation"
279  );
280 
281  // Allocate the interpolator
282  if (this->faceValues_)
283  {
284  mapperPtr_.reset
285  (
287  (
288  samplePoints,
289  this->localPosition(this->patch_.faceCentres()),
290  perturb_,
291  nearestOnly
292  )
293  );
294  }
295  else
296  {
297  mapperPtr_.reset
298  (
300  (
301  samplePoints,
302  this->localPosition(this->patch_.localPoints()),
303  perturb_,
304  nearestOnly
305  )
306  );
307  }
308 
309 
310  // Read the times for which data is available
311  const fileName samplePointsDir = samplePointsFile.path();
312  sampleTimes_ = Time::findTimes(samplePointsDir);
313 
314  DebugInfo
315  << "In directory "
316  << samplePointsDir << " found times "
317  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
318  << endl;
319  }
320 
321 
322  // Find current time in sampleTimes
323  label lo = -1;
324  label hi = -1;
325 
326  bool foundTime = mapperPtr_().findTime
327  (
328  sampleTimes_,
329  startSampleTime_,
330  t, //mesh.time().value(),
331  lo,
332  hi
333  );
334 
335  if (!foundTime)
336  {
338  << "Cannot find starting sampling values for index "
339  << t << nl
340  << "Have sampling values for "
341  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
342  << "In directory "
343  << time.constant()/"boundaryData"/this->patch_.name()
344  << "\n on patch " << this->patch_.name()
345  << " of field " << fieldTableName_
346  << exit(FatalError);
347  }
348 
349 
350  // Update sampled data fields.
351 
352  if (lo != startSampleTime_)
353  {
354  startSampleTime_ = lo;
355 
356  if (startSampleTime_ == endSampleTime_)
357  {
358  // No need to reread since are end values
359  if (debug)
360  {
361  Pout<< "checkTable : Setting startValues to (already read) "
362  << "boundaryData"
363  /this->patch_.name()
364  /sampleTimes_[startSampleTime_].name()
365  << endl;
366  }
367  startSampledValues_ = endSampledValues_;
368  startAverage_ = endAverage_;
369  }
370  else
371  {
372  if (debug)
373  {
374  Pout<< "checkTable : Reading startValues from "
375  << "boundaryData"
376  /this->patch_.name()
377  /sampleTimes_[lo].name()
378  << endl;
379  }
380 
381 
382  // Reread values and interpolate
383  const fileName valsFile
384  (
385  time.globalPath()
386  /time.constant()
387  /"boundaryData"
388  /this->patch_.name()
389  /sampleTimes_[startSampleTime_].name()
390  /fieldTableName_
391  );
392 
393  IOobject io
394  (
395  valsFile, // absolute path
396  time,
397  IOobject::MUST_READ,
398  IOobject::NO_WRITE,
399  false, // no need to register
400  true // is global object (currently not used)
401  );
402 
403  const rawIOField<Type> vals(io, setAverage_);
404  if (setAverage_)
405  {
406  startAverage_ = vals.average();
407  }
408 
409  if (vals.size() != mapperPtr_().sourceSize())
410  {
412  << "Number of values (" << vals.size()
413  << ") differs from the number of points ("
414  << mapperPtr_().sourceSize()
415  << ") in file " << valsFile << exit(FatalError);
416  }
417 
418  startSampledValues_ = mapperPtr_().interpolate(vals);
419  }
420  }
421 
422  if (hi != endSampleTime_)
423  {
424  endSampleTime_ = hi;
425 
426  if (endSampleTime_ == -1)
427  {
428  // endTime no longer valid. Might as well clear endValues.
429  if (debug)
430  {
431  Pout<< "checkTable : Clearing endValues" << endl;
432  }
433  endSampledValues_.clear();
434  }
435  else
436  {
437  if (debug)
438  {
439  Pout<< "checkTable : Reading endValues from "
440  << "boundaryData"
441  /this->patch_.name()
442  /sampleTimes_[endSampleTime_].name()
443  << endl;
444  }
445 
446  // Reread values and interpolate
447  fileName valsFile
448  (
449  time.globalPath()
450  /time.constant()
451  /"boundaryData"
452  /this->patch_.name()
453  /sampleTimes_[endSampleTime_].name()
454  /fieldTableName_
455  );
456 
457  IOobject io
458  (
459  valsFile, // absolute path
460  time,
461  IOobject::MUST_READ,
462  IOobject::NO_WRITE,
463  false, // no need to register
464  true // is global object (currently not used)
465  );
466 
467  const rawIOField<Type> vals(io, setAverage_);
468  if (setAverage_)
469  {
470  endAverage_ = vals.average();
471  }
472 
473  if (vals.size() != mapperPtr_().sourceSize())
474  {
476  << "Number of values (" << vals.size()
477  << ") differs from the number of points ("
478  << mapperPtr_().sourceSize()
479  << ") in file " << valsFile << exit(FatalError);
480  }
481 
482  endSampledValues_ = mapperPtr_().interpolate(vals);
483  }
484  }
485 }
486 
487 
488 template<class Type>
491 (
492  const scalar x
493 ) const
494 {
495  checkTable(x);
496 
497  auto tfld = tmp<Field<Type>>::New(startSampledValues_.size());
498  auto& fld = tfld.ref();
499  Type wantedAverage;
500 
501  if (endSampleTime_ == -1)
502  {
503  // Only start value
504  if (debug)
505  {
506  Pout<< "MappedFile<Type>::value : Sampled, non-interpolated values"
507  << " from start time:"
508  << sampleTimes_[startSampleTime_].name() << nl;
509  }
510 
511  fld = startSampledValues_;
512  wantedAverage = startAverage_;
513  }
514  else
515  {
516  scalar start = sampleTimes_[startSampleTime_].value();
517  scalar end = sampleTimes_[endSampleTime_].value();
518 
519  scalar s = (x - start)/(end - start);
520 
521  if (debug)
522  {
523  Pout<< "MappedFile<Type>::value : Sampled, interpolated values"
524  << " between start time:"
525  << sampleTimes_[startSampleTime_].name()
526  << " and end time:" << sampleTimes_[endSampleTime_].name()
527  << " with weight:" << s << endl;
528  }
529 
530  fld = ((1 - s)*startSampledValues_ + s*endSampledValues_);
531  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
532  }
533 
534  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
535  // offsetting.
536  if (setAverage_)
537  {
538  Type averagePsi;
539  if (this->faceValues_)
540  {
541  const scalarField magSf(mag(this->patch_.faceAreas()));
542  averagePsi = gSum(magSf*fld)/gSum(magSf);
543  }
544  else
545  {
546  averagePsi = gAverage(fld);
547  }
548 
549  if (debug)
550  {
551  Pout<< "MappedFile<Type>::value :"
552  << " actual average:" << averagePsi
553  << " wanted average:" << wantedAverage
554  << endl;
555  }
556 
557  if (mag(averagePsi) < VSMALL)
558  {
559  // Field too small to scale. Offset instead.
560  const Type offset = wantedAverage - averagePsi;
561  if (debug)
562  {
563  Pout<< "MappedFile<Type>::value :"
564  << " offsetting with:" << offset << endl;
565  }
566  fld += offset;
567  }
568  else
569  {
570  const scalar scale = mag(wantedAverage)/mag(averagePsi);
571 
572  if (debug)
573  {
574  Pout<< "MappedFile<Type>::value :"
575  << " scaling with:" << scale << endl;
576  }
577  fld *= scale;
578  }
579  }
580 
581  // Apply offset to mapped values
582  if (offset_)
583  {
584  fld += offset_->value(x);
585  }
586 
587  if (debug)
588  {
589  Pout<< "MappedFile<Type>::value : set fixedValue to min:" << gMin(fld)
590  << " max:" << gMax(fld)
591  << " avg:" << gAverage(fld) << endl;
592  }
593 
594  return this->transform(tfld);
595 }
596 
597 
598 template<class Type>
601 (
602  const scalar x1,
603  const scalar x2
604 ) const
605 {
607  return nullptr;
608 }
609 
610 
611 template<class Type>
613 (
614  Ostream& os
615 ) const
616 {
618 
619  // Check if field name explicitly provided
620  // (e.g. through timeVaryingMapped bc)
621  if (dictConstructed_)
622  {
623  os.writeEntry(this->name(), type());
624 
626  (
627  "fieldTable",
628  this->name(),
629  fieldTableName_
630  );
631  }
632 
633  if (setAverage_)
634  {
635  os.writeEntry("setAverage", setAverage_);
636  }
637 
638  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
639 
640  os.writeEntryIfDifferent<word>("points", "points", pointsName_);
641 
643  (
644  "mapMethod",
645  "planarInterpolation",
646  mapMethod_
647  );
648 
649  if (offset_)
650  {
651  offset_->writeData(os);
652  }
653 }
654 
655 
656 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PatchFunction1Types::MappedFile::autoMap
virtual void autoMap(const FieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
Definition: MappedFile.C:194
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
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
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
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::PatchFunction1Types::MappedFile::writeData
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: MappedFile.C:613
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:209
writeData
const bool writeData(pdfDictionary.get< bool >("writeData"))
Foam::FieldMapper
Abstract base class to hold the Field mapping addressing and weights.
Definition: FieldMapper.H:49
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::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:436
Foam::Field< scalar >
Foam::rawIOField
Like IOField but falls back to raw IFstream if no header found. Optionally reads average value....
Definition: rawIOField.H:52
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:67
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
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::PatchFunction1Types::MappedFile::integrate
virtual tmp< Field< Type > > integrate(const scalar x1, const scalar x2) const
Integrate between two values.
Definition: MappedFile.C:601
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:63
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PatchFunction1Types::MappedFile
Patch value mapping from file.
Definition: MappedFile.H:101
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:372
Foam::PatchFunction1Types::MappedFile::MappedFile
MappedFile(const polyPatch &pp, const word &type, const word &entryName, const dictionary &dict, const bool faceValues=true)
Construct from entry name and dictionary.
Definition: MappedFile.C:35
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:354
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::List< label >
Foam::PatchFunction1Types::MappedFile::rmap
virtual void rmap(const PatchFunction1< Type > &pf1, const labelList &addr)
Reverse map the given PatchFunction1 onto this PatchFunction1.
Definition: MappedFile.C:214
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
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::TimePaths::globalPath
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:72
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:392
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::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::PatchFunction1Types::MappedFile::value
virtual tmp< Field< Type > > value(const scalar) const
Return MappedFile value.
Definition: MappedFile.C:491