timeVaryingMappedFixedValueFvPatchField.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "Time.H"
28 #include "AverageIOField.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
40 (
41  const fvPatch& p,
43 )
44 :
46  fieldTableName_(iF.name()),
47  setAverage_(false),
48  perturb_(0),
49  mapperPtr_(NULL),
50  sampleTimes_(0),
51  startSampleTime_(-1),
52  startSampledValues_(0),
53  startAverage_(pTraits<Type>::zero),
54  endSampleTime_(-1),
55  endSampledValues_(0),
56  endAverage_(pTraits<Type>::zero),
57  offset_()
58 {}
59 
60 
61 template<class Type>
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
71  fieldTableName_(iF.name()),
72  setAverage_(readBool(dict.lookup("setAverage"))),
73  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
74  mapMethod_
75  (
77  (
78  "mapMethod",
79  "planarInterpolation"
80  )
81  ),
82  mapperPtr_(NULL),
83  sampleTimes_(0),
84  startSampleTime_(-1),
85  startSampledValues_(0),
86  startAverage_(pTraits<Type>::zero),
87  endSampleTime_(-1),
88  endSampledValues_(0),
89  endAverage_(pTraits<Type>::zero),
90  offset_(DataEntry<Type>::New("offset", dict))
91 {
92  if
93  (
94  mapMethod_ != "planarInterpolation"
95  && mapMethod_ != "nearest"
96  )
97  {
99  (
100  dict
101  ) << "mapMethod should be one of 'planarInterpolation'"
102  << ", 'nearest'" << exit(FatalIOError);
103  }
104 
105 
106  dict.readIfPresent("fieldTableName", fieldTableName_);
107 
108  if (dict.found("value"))
109  {
111  }
112  else
113  {
114  // Note: we use evaluate() here to trigger updateCoeffs followed
115  // by re-setting of fvatchfield::updated_ flag. This is
116  // so if first use is in the next time step it retriggers
117  // a new update.
118  this->evaluate(Pstream::blocking);
119  }
120 }
121 
122 
123 template<class Type>
126 (
128  const fvPatch& p,
130  const fvPatchFieldMapper& mapper
131 )
132 :
133  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
134  fieldTableName_(ptf.fieldTableName_),
135  setAverage_(ptf.setAverage_),
136  perturb_(ptf.perturb_),
137  mapMethod_(ptf.mapMethod_),
138  mapperPtr_(NULL),
139  sampleTimes_(0),
140  startSampleTime_(-1),
141  startSampledValues_(0),
142  startAverage_(pTraits<Type>::zero),
143  endSampleTime_(-1),
144  endSampledValues_(0),
145  endAverage_(pTraits<Type>::zero),
146  offset_(ptf.offset_, false)
147 {}
148 
149 
150 template<class Type>
153 (
155 )
156 :
158  fieldTableName_(ptf.fieldTableName_),
159  setAverage_(ptf.setAverage_),
160  perturb_(ptf.perturb_),
161  mapMethod_(ptf.mapMethod_),
162  mapperPtr_(NULL),
163  sampleTimes_(ptf.sampleTimes_),
164  startSampleTime_(ptf.startSampleTime_),
165  startSampledValues_(ptf.startSampledValues_),
166  startAverage_(ptf.startAverage_),
167  endSampleTime_(ptf.endSampleTime_),
168  endSampledValues_(ptf.endSampledValues_),
169  endAverage_(ptf.endAverage_),
170  offset_(ptf.offset_, false)
171 {}
172 
173 
174 template<class Type>
177 (
180 )
181 :
183  fieldTableName_(ptf.fieldTableName_),
184  setAverage_(ptf.setAverage_),
185  perturb_(ptf.perturb_),
186  mapMethod_(ptf.mapMethod_),
187  mapperPtr_(NULL),
188  sampleTimes_(ptf.sampleTimes_),
189  startSampleTime_(ptf.startSampleTime_),
190  startSampledValues_(ptf.startSampledValues_),
191  startAverage_(ptf.startAverage_),
192  endSampleTime_(ptf.endSampleTime_),
193  endSampledValues_(ptf.endSampledValues_),
194  endAverage_(ptf.endAverage_),
195  offset_(ptf.offset_, false)
196 {}
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class Type>
203 (
204  const fvPatchFieldMapper& m
205 )
206 {
208  if (startSampledValues_.size())
209  {
210  startSampledValues_.autoMap(m);
211  endSampledValues_.autoMap(m);
212  }
213  // Clear interpolator
214  mapperPtr_.clear();
215  startSampleTime_ = -1;
216  endSampleTime_ = -1;
217 }
218 
219 
220 template<class Type>
222 (
223  const fvPatchField<Type>& ptf,
224  const labelList& addr
225 )
226 {
228 
230  refCast<const timeVaryingMappedFixedValueFvPatchField<Type> >(ptf);
231 
232  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
233  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
234 
235  // Clear interpolator
236  mapperPtr_.clear();
237  startSampleTime_ = -1;
238  endSampleTime_ = -1;
239 }
240 
241 
242 template<class Type>
244 {
245  // Initialise
246  if (mapperPtr_.empty())
247  {
248  pointIOField samplePoints
249  (
250  IOobject
251  (
252  "points",
253  this->db().time().caseConstant(),
254  "boundaryData"/this->patch().name(),
255  this->db(),
258  false
259  )
260  );
261 
262  const fileName samplePointsFile = samplePoints.filePath();
263 
264  if (debug)
265  {
266  Info<< "timeVaryingMappedFixedValueFvPatchField :"
267  << " Read " << samplePoints.size() << " sample points from "
268  << samplePointsFile << endl;
269  }
270 
271 
272  // tbd: run-time selection
273  bool nearestOnly =
274  (
275  !mapMethod_.empty()
276  && mapMethod_ != "planarInterpolation"
277  );
278 
279  // Allocate the interpolator
280  mapperPtr_.reset
281  (
283  (
284  samplePoints,
285  this->patch().patch().faceCentres(),
286  perturb_,
287  nearestOnly
288  )
289  );
290 
291  // Read the times for which data is available
292  const fileName samplePointsDir = samplePointsFile.path();
293  sampleTimes_ = Time::findTimes(samplePointsDir);
294 
295  if (debug)
296  {
297  Info<< "timeVaryingMappedFixedValueFvPatchField : In directory "
298  << samplePointsDir << " found times "
300  << endl;
301  }
302  }
303 
304 
305  // Find current time in sampleTimes
306  label lo = -1;
307  label hi = -1;
308 
309  bool foundTime = mapperPtr_().findTime
310  (
311  sampleTimes_,
312  startSampleTime_,
313  this->db().time().value(),
314  lo,
315  hi
316  );
317 
318  if (!foundTime)
319  {
321  << "Cannot find starting sampling values for current time "
322  << this->db().time().value() << nl
323  << "Have sampling values for times "
325  << "In directory "
326  << this->db().time().constant()/"boundaryData"/this->patch().name()
327  << "\n on patch " << this->patch().name()
328  << " of field " << fieldTableName_
329  << exit(FatalError);
330  }
331 
332 
333  // Update sampled data fields.
334 
335  if (lo != startSampleTime_)
336  {
337  startSampleTime_ = lo;
338 
339  if (startSampleTime_ == endSampleTime_)
340  {
341  // No need to reread since are end values
342  if (debug)
343  {
344  Pout<< "checkTable : Setting startValues to (already read) "
345  << "boundaryData"
346  /this->patch().name()
347  /sampleTimes_[startSampleTime_].name()
348  << endl;
349  }
350  startSampledValues_ = endSampledValues_;
351  startAverage_ = endAverage_;
352  }
353  else
354  {
355  if (debug)
356  {
357  Pout<< "checkTable : Reading startValues from "
358  << "boundaryData"
359  /this->patch().name()
360  /sampleTimes_[lo].name()
361  << endl;
362  }
363 
364 
365  // Reread values and interpolate
367  (
368  IOobject
369  (
370  fieldTableName_,
371  this->db().time().caseConstant(),
372  "boundaryData"
373  /this->patch().name()
374  /sampleTimes_[startSampleTime_].name(),
375  this->db(),
378  false
379  )
380  );
381 
382  if (vals.size() != mapperPtr_().sourceSize())
383  {
385  << "Number of values (" << vals.size()
386  << ") differs from the number of points ("
387  << mapperPtr_().sourceSize()
388  << ") in file " << vals.objectPath() << exit(FatalError);
389  }
390 
391  startAverage_ = vals.average();
392  startSampledValues_ = mapperPtr_().interpolate(vals);
393  }
394  }
395 
396  if (hi != endSampleTime_)
397  {
398  endSampleTime_ = hi;
399 
400  if (endSampleTime_ == -1)
401  {
402  // endTime no longer valid. Might as well clear endValues.
403  if (debug)
404  {
405  Pout<< "checkTable : Clearing endValues" << endl;
406  }
407  endSampledValues_.clear();
408  }
409  else
410  {
411  if (debug)
412  {
413  Pout<< "checkTable : Reading endValues from "
414  << "boundaryData"
415  /this->patch().name()
416  /sampleTimes_[endSampleTime_].name()
417  << endl;
418  }
419 
420  // Reread values and interpolate
422  (
423  IOobject
424  (
425  fieldTableName_,
426  this->db().time().caseConstant(),
427  "boundaryData"
428  /this->patch().name()
429  /sampleTimes_[endSampleTime_].name(),
430  this->db(),
433  false
434  )
435  );
436 
437  if (vals.size() != mapperPtr_().sourceSize())
438  {
440  << "Number of values (" << vals.size()
441  << ") differs from the number of points ("
442  << mapperPtr_().sourceSize()
443  << ") in file " << vals.objectPath() << exit(FatalError);
444  }
445 
446  endAverage_ = vals.average();
447  endSampledValues_ = mapperPtr_().interpolate(vals);
448  }
449  }
450 }
451 
452 
453 template<class Type>
455 {
456  if (this->updated())
457  {
458  return;
459  }
460 
461 
462  checkTable();
463 
464  // Interpolate between the sampled data
465 
466  Type wantedAverage;
467 
468  if (endSampleTime_ == -1)
469  {
470  // only start value
471  if (debug)
472  {
473  Pout<< "updateCoeffs : Sampled, non-interpolated values"
474  << " from start time:"
475  << sampleTimes_[startSampleTime_].name() << nl;
476  }
477 
478  this->operator==(startSampledValues_);
479  wantedAverage = startAverage_;
480  }
481  else
482  {
483  scalar start = sampleTimes_[startSampleTime_].value();
484  scalar end = sampleTimes_[endSampleTime_].value();
485 
486  scalar s = (this->db().time().value() - start)/(end - start);
487 
488  if (debug)
489  {
490  Pout<< "updateCoeffs : Sampled, interpolated values"
491  << " between start time:"
492  << sampleTimes_[startSampleTime_].name()
493  << " and end time:" << sampleTimes_[endSampleTime_].name()
494  << " with weight:" << s << endl;
495  }
496 
497  this->operator==((1 - s)*startSampledValues_ + s*endSampledValues_);
498  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
499  }
500 
501  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
502  // offsetting.
503  if (setAverage_)
504  {
505  const Field<Type>& fld = *this;
506 
507  Type averagePsi =
508  gSum(this->patch().magSf()*fld)
509  /gSum(this->patch().magSf());
510 
511  if (debug)
512  {
513  Pout<< "updateCoeffs :"
514  << " actual average:" << averagePsi
515  << " wanted average:" << wantedAverage
516  << endl;
517  }
518 
519  if (mag(averagePsi) < VSMALL)
520  {
521  // Field too small to scale. Offset instead.
522  const Type offset = wantedAverage - averagePsi;
523  if (debug)
524  {
525  Pout<< "updateCoeffs :"
526  << " offsetting with:" << offset << endl;
527  }
528  this->operator==(fld + offset);
529  }
530  else
531  {
532  const scalar scale = mag(wantedAverage)/mag(averagePsi);
533 
534  if (debug)
535  {
536  Pout<< "updateCoeffs :"
537  << " scaling with:" << scale << endl;
538  }
539  this->operator==(scale*fld);
540  }
541  }
542 
543  // apply offset to mapped values
544  const scalar t = this->db().time().timeOutputValue();
545  this->operator==(*this + offset_->value(t));
546 
547  if (debug)
548  {
549  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
550  << " max:" << gMax(*this)
551  << " avg:" << gAverage(*this) << endl;
552  }
553 
555 }
556 
557 
558 template<class Type>
560 {
562  os.writeKeyword("setAverage") << setAverage_ << token::END_STATEMENT << nl;
563  if (perturb_ != 1e-5)
564  {
565  os.writeKeyword("perturb") << perturb_ << token::END_STATEMENT << nl;
566  }
567 
568  if (fieldTableName_ != this->dimensionedInternalField().name())
569  {
570  os.writeKeyword("fieldTableName") << fieldTableName_
571  << token::END_STATEMENT << nl;
572  }
573 
574  if
575  (
576  (
577  !mapMethod_.empty()
578  && mapMethod_ != "planarInterpolation"
579  )
580  )
581  {
582  os.writeKeyword("mapMethod") << mapMethod_
583  << token::END_STATEMENT << nl;
584  }
585 
586  offset_->writeData(os);
587 
588  this->writeEntry("value", os);
589 }
590 
591 
592 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
593 
594 } // End namespace Foam
595 
596 // ************************************************************************* //
Foam::timeVaryingMappedFixedValueFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: timeVaryingMappedFixedValueFvPatchField.C:222
Foam::fvPatchField< Type >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvPatchField::operator==
virtual void operator==(const fvPatchField< Type > &)
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::timeVaryingMappedFixedValueFvPatchField::startAverage_
Type startAverage_
If setAverage: starting average value.
Definition: timeVaryingMappedFixedValueFvPatchField.H:155
Foam::IOobject::filePath
fileName filePath() const
Return complete path + object name if the file exists.
Definition: IOobject.C:336
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField< vector >
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::timeVaryingMappedFixedValueFvPatchField::perturb_
scalar perturb_
Fraction of perturbation (fraction of bounding box) to add.
Definition: timeVaryingMappedFixedValueFvPatchField.H:137
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::timeVaryingMappedFixedValueFvPatchField::endSampleTime_
label endSampleTime_
Current end index in sampleTimes.
Definition: timeVaryingMappedFixedValueFvPatchField.H:158
Foam::timeVaryingMappedFixedValueFvPatchField
This boundary conditions interpolates the values from a set of supplied points in space and time....
Definition: timeVaryingMappedFixedValueFvPatchField.H:124
Foam::timeVaryingMappedFixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: timeVaryingMappedFixedValueFvPatchField.C:559
dimensionedInternalField
rDeltaT dimensionedInternalField()
AverageIOField.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::fileName::path
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:293
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::timeVaryingMappedFixedValueFvPatchField::endAverage_
Type endAverage_
If setAverage: end average value.
Definition: timeVaryingMappedFixedValueFvPatchField.H:164
Foam::pointToPointPlanarInterpolation
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
Definition: pointToPointPlanarInterpolation.H:51
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::FatalIOError
IOerror FatalIOError
Foam::timeVaryingMappedFixedValueFvPatchField::startSampleTime_
label startSampleTime_
Current starting index in sampleTimes.
Definition: timeVaryingMappedFixedValueFvPatchField.H:149
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::timeVaryingMappedFixedValueFvPatchField::fieldTableName_
word fieldTableName_
Name of the field data table, defaults to the name of the field.
Definition: timeVaryingMappedFixedValueFvPatchField.H:131
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:79
Foam::fileName::name
word name() const
Return file name (part beyond last /)
Definition: fileName.C:212
Foam::IOobject::objectPath
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:376
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::timeVaryingMappedFixedValueFvPatchField::checkTable
void checkTable()
Find boundary data inbetween current time and interpolate.
Definition: timeVaryingMappedFixedValueFvPatchField.C:243
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field< Type >
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::timeVaryingMappedFixedValueFvPatchField::timeVaryingMappedFixedValueFvPatchField
timeVaryingMappedFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: timeVaryingMappedFixedValueFvPatchField.C:40
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
Foam::timeVaryingMappedFixedValueFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: timeVaryingMappedFixedValueFvPatchField.C:454
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::timeVaryingMappedFixedValueFvPatchField::offset_
autoPtr< DataEntry< Type > > offset_
Time varying offset values to interpolated data.
Definition: timeVaryingMappedFixedValueFvPatchField.H:167
Foam::timeVaryingMappedFixedValueFvPatchField::setAverage_
bool setAverage_
If true adjust the mapped field to maintain average value.
Definition: timeVaryingMappedFixedValueFvPatchField.H:134
Foam::timeVaryingMappedFixedValueFvPatchField::sampleTimes_
instantList sampleTimes_
List of boundaryData time directories.
Definition: timeVaryingMappedFixedValueFvPatchField.H:146
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fvPatchField< Type >::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:286
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
timeVaryingMappedFixedValueFvPatchField.H
Foam::AverageIOField::average
const Type & average() const
Definition: AverageIOField.H:92
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::timeVaryingMappedFixedValueFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: timeVaryingMappedFixedValueFvPatchField.C:203
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::timeVaryingMappedFixedValueFvPatchField::mapMethod_
word mapMethod_
Interpolation scheme to use.
Definition: timeVaryingMappedFixedValueFvPatchField.H:140
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::Time::findTimes
static instantList findTimes(const fileName &, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: findTimes.C:38
Foam::AverageIOField
A primitive field + average with IO.
Definition: AverageIOField.H:50
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
Foam::fvPatchField< Type >::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:223
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::pointToPointPlanarInterpolation::timeNames
static wordList timeNames(const instantList &)
Helper: extract words of times.
Definition: pointToPointPlanarInterpolation.C:347
Foam::timeVaryingMappedFixedValueFvPatchField::startSampledValues_
Field< Type > startSampledValues_
Interpolated values from startSampleTime.
Definition: timeVaryingMappedFixedValueFvPatchField.H:152
Foam::timeVaryingMappedFixedValueFvPatchField::endSampledValues_
Field< Type > endSampledValues_
Interpolated values from endSampleTime.
Definition: timeVaryingMappedFixedValueFvPatchField.H:161
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51