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