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-2022 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 #include "clockTime.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const bool dictConstructed,
38  const polyPatch& pp,
39  const word& entryName,
40  const dictionary& dict,
41  const word& fieldTableName,
42  const bool faceValues
43 )
44 :
45  PatchFunction1<Type>(pp, entryName, dict, faceValues),
46  dictConstructed_(dictConstructed),
47  setAverage_(dict.getOrDefault("setAverage", false)),
48  perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
49  fieldTableName_(fieldTableName),
50  pointsName_(dict.getOrDefault<word>("points", "points")),
51  mapMethod_(),
52  filterRadius_(dict.getOrDefault<scalar>("filterRadius", 0)),
53  filterSweeps_(dict.getOrDefault<label>("filterSweeps", 0)),
54  filterFieldPtr_(nullptr),
55  readerFormat_(),
56  readerFile_(),
57  readerPtr_(nullptr),
58  mapperPtr_(nullptr),
59  sampleTimes_(),
60  sampleIndex_(-1, -1),
61  sampleAverage_(Zero, Zero),
62  sampleValues_(),
63  offset_(Function1<Type>::NewIfPresent("offset", dict))
64 {
65  if (fieldTableName_.empty())
66  {
67  fieldTableName_ = entryName;
68  }
69 
70  // Simple sanity check
71  if ((filterSweeps_ < 1) || (filterRadius_ <= VSMALL))
72  {
73  filterRadius_ = 0;
74  filterSweeps_ = 0;
75  }
76 
77  if (dict.readIfPresent("sampleFormat", readerFormat_))
78  {
79  dict.readEntry("sampleFile", readerFile_);
80 
81  fileName fName(readerFile_);
82  fName.expand();
83 
84  readerPtr_ = surfaceReader::New
85  (
86  readerFormat_,
87  fName,
88  surfaceReader::formatOptions(dict, readerFormat_, "readOptions")
89  );
90  }
91 
92  if (debug)
93  {
94  Info<< "mappedFile:" << nl;
95  if (readerFormat_.empty())
96  {
97  Info<< " boundary format" << nl;
98  }
99  else
100  {
101  Info<< " format:" << readerFormat_
102  << " file:" << readerFile_ << nl;
103  }
104 
105  Info<< " filter radius=" << filterRadius_
106  << " sweeps=" << filterSweeps_ << endl;
107  }
108 
109  if
110  (
111  dict.readIfPresent("mapMethod", mapMethod_)
112  && !mapMethod_.empty()
113  && mapMethod_ != "nearest"
114  && !mapMethod_.starts_with("planar")
115  )
116  {
118  << "Unknown mapMethod type " << mapMethod_
119  << "\n\nValid mapMethod types :\n"
120  << "(nearest planar)" << nl
122  }
123 }
124 
125 
126 template<class Type>
128 (
129  const polyPatch& pp,
130  const word& redirectType,
131  const word& entryName,
132  const dictionary& dict,
133  const bool faceValues
134 )
135 :
136  MappedFile<Type>
137  (
138  true, // dictConstructed = true
139  pp,
140  entryName,
141  dict,
142  dict.getOrDefault<word>("fieldTable", entryName),
143  faceValues
144  )
145 {}
146 
147 
148 template<class Type>
150 (
151  const polyPatch& pp,
152  const word& entryName,
153  const dictionary& dict,
154  const word& fieldTableName,
155  const bool faceValues
156 )
157 :
158  MappedFile<Type>
159  (
160  false, // dictConstructed = false
161  pp,
162  entryName,
163  dict,
164  fieldTableName,
165  faceValues
166  )
167 {}
168 
169 
170 template<class Type>
172 (
173  const MappedFile<Type>& rhs
174 )
175 :
176  MappedFile<Type>(rhs, rhs.patch())
177 {}
178 
179 
180 template<class Type>
182 (
183  const MappedFile<Type>& rhs,
184  const polyPatch& pp
185 )
186 :
187  PatchFunction1<Type>(rhs, pp),
188  dictConstructed_(rhs.dictConstructed_),
189  setAverage_(rhs.setAverage_),
190  perturb_(rhs.perturb_),
191  fieldTableName_(rhs.fieldTableName_),
192  pointsName_(rhs.pointsName_),
193  mapMethod_(rhs.mapMethod_),
194  filterRadius_(rhs.filterRadius_),
195  filterSweeps_(rhs.filterSweeps_),
196  filterFieldPtr_(nullptr),
197  readerFormat_(rhs.readerFormat_),
198  readerFile_(rhs.readerFile_),
199  readerPtr_(nullptr),
200  mapperPtr_(rhs.mapperPtr_.clone()),
201  sampleTimes_(rhs.sampleTimes_),
202  sampleIndex_(rhs.sampleIndex_),
203  sampleAverage_(rhs.sampleAverage_),
204  sampleValues_(rhs.sampleValues_),
205  offset_(rhs.offset_.clone())
206 {
207  if (!readerFormat_.empty() && !readerFile_.empty())
208  {
209  fileName fName(readerFile_);
210  fName.expand();
211 
212  readerPtr_ = surfaceReader::New(readerFormat_, fName);
213  }
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
219 template<class Type>
221 (
222  const FieldMapper& mapper
223 )
224 {
226 
227  if (sampleValues_.first().size())
228  {
229  sampleValues_.first().autoMap(mapper);
230  }
231  if (sampleValues_.second().size())
232  {
233  sampleValues_.second().autoMap(mapper);
234  }
235 
236  // Clear interpolator
237  filterFieldPtr_.reset(nullptr);
238  mapperPtr_.reset(nullptr);
239  sampleIndex_ = labelPair(-1, -1);
240 }
241 
242 
243 template<class Type>
245 (
246  const PatchFunction1<Type>& pf1,
247  const labelList& addr
248 )
249 {
250  PatchFunction1<Type>::rmap(pf1, addr);
251 
252  const auto& tiptf =
253  refCast<const PatchFunction1Types::MappedFile<Type>>(pf1);
254 
255  if (tiptf.sampleValues_.first().size())
256  {
257  sampleValues_.first().resize(this->size());
258  sampleValues_.first().rmap(tiptf.sampleValues_.first(), addr);
259  }
260 
261  if (tiptf.sampleValues_.second().size())
262  {
263  sampleValues_.second().resize(this->size());
264  sampleValues_.second().rmap(tiptf.sampleValues_.second(), addr);
265  }
266 
267  // Clear interpolator
268  filterFieldPtr_.reset(nullptr);
269  mapperPtr_.reset(nullptr);
270  sampleIndex_ = labelPair(-1, -1);
271 }
272 
273 
274 template<class Type>
276 (
277  const label sampleIndex,
278  Field<Type>& field,
279  Type& avg
280 ) const
281 {
282  tmp<Field<Type>> tvalues;
283 
284  // Update sampled data fields
285  if (readerPtr_)
286  {
287  wordList fieldNames = readerPtr_->fieldNames(sampleIndex);
288 
289  label fieldIndex = fieldNames.find(fieldTableName_);
290 
291  if (fieldIndex < 0)
292  {
294  << "Sample field='" << fieldTableName_
295  << "' not found. Known field names: "
296  << flatOutput(fieldNames) << nl
297  << exit(FatalError);
298  }
299 
300  if (debug)
301  {
302  Pout<< "checkTable : Update index=" << sampleIndex
303  << " field=" << fieldNames[fieldIndex] << endl;
304  }
305 
306  tvalues = readerPtr_->field
307  (
308  sampleIndex,
309  fieldIndex,
310  pTraits<Type>::zero
311  );
312 
313  if (tvalues().size() != mapperPtr_().sourceSize())
314  {
316  << "Number of values (" << tvalues().size()
317  << ") differs from the number of points ("
318  << mapperPtr_().sourceSize() << ")"
319  << exit(FatalError);
320  }
321  }
322  else
323  {
324  const polyMesh& mesh = this->patch_.boundaryMesh().mesh();
325  const Time& time = mesh.time();
326 
327  if (debug)
328  {
329  Pout<< "checkTable : Update index=" << sampleIndex
330  << " Reading values from "
331  <<
332  (
333  "boundaryData"
334  / this->patch_.name()
335  / sampleTimes_[sampleIndex].name()
336  / fieldTableName_
337  ) << endl;
338  }
339 
340  // Reread values and interpolate
341  const fileName valsFile
342  (
343  time.globalPath()
344  /time.constant()
345  /mesh.dbDir() // region
346  /"boundaryData"
347  /this->patch_.name()
348  /sampleTimes_[sampleIndex].name()
349  /fieldTableName_
350  );
351 
352  IOobject io
353  (
354  valsFile, // absolute path
355  time,
359  true // is global object (currently not used)
360  );
361 
362  rawIOField<Type> vals(io, setAverage_);
363 
364  if (vals.hasAverage())
365  {
366  avg = vals.average();
367  }
368 
369  if (vals.size() != mapperPtr_().sourceSize())
370  {
372  << "Number of values (" << vals.size()
373  << ") differs from the number of points ("
374  << mapperPtr_().sourceSize()
375  << ") in file " << valsFile
376  << exit(FatalError);
377  }
378 
379  tvalues = tmp<Field<Type>>::New(std::move(vals.field()));
380  }
381 
382  if (filterFieldPtr_)
383  {
384  DebugInfo
385  << "apply " << filterSweeps_ << " filter sweeps" << endl;
386 
387  tvalues = filterFieldPtr_().evaluate(tvalues, filterSweeps_);
388  }
389 
390  // From input values to interpolated (sampled) positions
391  field = mapperPtr_().interpolate(tvalues);
392 }
393 
394 
395 template<class Type>
397 (
398  const scalar t
399 ) const
400 {
401  const polyMesh& mesh = this->patch_.boundaryMesh().mesh();
402  const Time& time = mesh.time();
403 
404  // Initialise
405  if (!mapperPtr_ && readerPtr_)
406  {
407  clockTime timing;
408 
409  auto& reader = readerPtr_();
410 
411  const meshedSurface& geom = reader.geometry(0);
412 
413  sampleTimes_ = reader.times();
414 
415  // We may have been passed in geometry without any faces
416  // eg, boundaryData
417 
418  const pointField& samplePoints =
419  (
420  geom.nFaces() ? geom.faceCentres() : geom.points()
421  );
422 
423  DebugInfo
424  << "Read " << samplePoints.size() << " sample points from "
425  << readerFile_ << endl
426  << "Found times "
428  << "... in " << timing.timeIncrement() << 's' << endl;
429 
430  // tbd: run-time selection
431  const bool nearestOnly =
432  (
433  !mapMethod_.empty() && !mapMethod_.starts_with("planar")
434  );
435 
436  // Allocate the interpolator
437  if (this->faceValues())
438  {
439  mapperPtr_.reset
440  (
441  new pointToPointPlanarInterpolation
442  (
443  samplePoints,
444  this->localPosition(this->patch_.faceCentres()),
445  perturb_,
446  nearestOnly
447  )
448  );
449  }
450  else
451  {
452  mapperPtr_.reset
453  (
454  new pointToPointPlanarInterpolation
455  (
456  samplePoints,
457  this->localPosition(this->patch_.localPoints()),
458  perturb_,
459  nearestOnly
460  )
461  );
462  }
463 
464  DebugInfo
465  << "Created point/point planar interpolation"
466  << " - in " << timing.timeIncrement() << 's' << endl;
467 
468 
469  // Setup median filter (if any)
470  if (filterSweeps_ > 0)
471  {
472  filterFieldPtr_.reset(new FilterField(geom, filterRadius_));
473 
474  DebugInfo
475  << "Calculated field-filter"
476  << " - in " << timing.timeIncrement() << 's' << endl;
477  }
478  else
479  {
480  filterFieldPtr_.reset(nullptr);
481  }
482  }
483  else if (!mapperPtr_)
484  {
485  clockTime timing;
486 
487  // Reread values and interpolate
488  const fileName samplePointsFile
489  (
490  time.globalPath()
491  /time.constant() // instance
492  /mesh.dbDir() // region
493  /"boundaryData"
494  /this->patch_.name()
495  /pointsName_
496  );
497 
498  IOobject io
499  (
500  samplePointsFile, // absolute path
501  time,
505  true // is global object (currently not used)
506  );
507 
508  // Read data (no average value!)
509  const rawIOField<point> samplePoints(io);
510 
511  // Read the times for which data is available
512  sampleTimes_ = Time::findTimes(samplePointsFile.path());
513 
514  DebugInfo
515  << "Read " << samplePoints.size() << " sample points from "
516  << samplePointsFile << endl
517  << "Found times "
519  << "... in " << timing.timeIncrement() << 's' << endl;
520 
521 
522  // tbd: run-time selection
523  const bool nearestOnly =
524  (
525  !mapMethod_.empty() && !mapMethod_.starts_with("planar")
526  );
527 
528  // Allocate the interpolator
529  if (this->faceValues())
530  {
531  mapperPtr_.reset
532  (
533  new pointToPointPlanarInterpolation
534  (
535  samplePoints,
536  this->localPosition(this->patch_.faceCentres()),
537  perturb_,
538  nearestOnly
539  )
540  );
541  }
542  else
543  {
544  mapperPtr_.reset
545  (
546  new pointToPointPlanarInterpolation
547  (
548  samplePoints,
549  this->localPosition(this->patch_.localPoints()),
550  perturb_,
551  nearestOnly
552  )
553  );
554  }
555 
556  DebugInfo
557  << "Created point/point planar interpolation"
558  << " - in " << timing.timeIncrement() << 's' << endl;
559 
560 
561  // Setup median filter (if any)
562  if (filterSweeps_ > 0)
563  {
564  filterFieldPtr_.reset(new FilterField(samplePoints, filterRadius_));
565 
566  DebugInfo
567  << "Calculated field-filter"
568  << " - in " << timing.timeIncrement() << 's' << endl;
569  }
570  else
571  {
572  filterFieldPtr_.reset(nullptr);
573  }
574  }
575 
576 
577  // Find range of current time indices in sampleTimes
578  labelPair timeIndices = instant::findRange
579  (
580  sampleTimes_,
581  t, //mesh.time().value(),
582  sampleIndex_.first()
583  );
584 
585  if (timeIndices.first() < 0)
586  {
588  << "Cannot find starting sampling values for index "
589  << t << nl
590  << "Have sampling values for "
592  << "In directory "
593  << time.constant()/mesh.dbDir()/"boundaryData"/this->patch_.name()
594  << "\n on patch " << this->patch_.name()
595  << " of field " << fieldTableName_
596  << exit(FatalError);
597  }
598 
599 
600  // Update sampled data fields.
601 
602  if (sampleIndex_.first() != timeIndices.first())
603  {
604  sampleIndex_.first() = timeIndices.first();
605 
606  if (sampleIndex_.first() == sampleIndex_.second())
607  {
608  // Can reuse previous end values
609  sampleValues_.first() = sampleValues_.second();
610  sampleAverage_.first() = sampleAverage_.second();
611  }
612  else
613  {
614  // Update first() values
615  this->updateSampledValues
616  (
617  sampleIndex_.first(),
618  sampleValues_.first(),
619  sampleAverage_.first()
620  );
621  }
622  }
623 
624  if (sampleIndex_.second() != timeIndices.second())
625  {
626  sampleIndex_.second() = timeIndices.second();
627 
628  if (sampleIndex_.second() == -1)
629  {
630  // Index no longer valid - clear values accordingly
631  sampleValues_.second().clear();
632  }
633  else
634  {
635  // Update second() values
636  this->updateSampledValues
637  (
638  sampleIndex_.second(),
639  sampleValues_.second(),
640  sampleAverage_.second()
641  );
642  }
643  }
644 }
645 
646 
647 template<class Type>
650 (
651  const scalar x
652 ) const
653 {
654  checkTable(x);
655 
656  // Interpolate between the sampled data
657  auto tfld = tmp<Field<Type>>::New();
658  auto& fld = tfld.ref();
659  Type wantedAverage;
660 
661  if (sampleIndex_.second() == -1)
662  {
663  // Only start value
664  fld = sampleValues_.first();
665  wantedAverage = sampleAverage_.first();
666  }
667  else
668  {
669  const scalar beg = sampleTimes_[sampleIndex_.first()].value();
670  const scalar end = sampleTimes_[sampleIndex_.second()].value();
671  const scalar s = (x - beg)/(end - beg);
672 
673  fld = (1 - s)*sampleValues_.first() + s*sampleValues_.second();
674 
675  wantedAverage
676  = (1 - s)*sampleAverage_.first() + s*sampleAverage_.second();
677 
678  DebugInfo
679  << "MappedFile<Type>::value : "
680  << "Sampled, interpolated values"
681  << " between time:"
682  << sampleTimes_[sampleIndex_.first()].name()
683  << " and time:" << sampleTimes_[sampleIndex_.second()].name()
684  << " with weight:" << s << endl;
685  }
686 
687  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
688  // offsetting.
689  if (setAverage_)
690  {
691  Type averagePsi;
692 
693  if (this->faceValues())
694  {
695  const scalarField magSf(mag(this->patch_.faceAreas()));
696  averagePsi = gSum(magSf*fld)/gSum(magSf);
697  }
698  else
699  {
700  averagePsi = gAverage(fld);
701  }
702 
703  if (debug)
704  {
705  Pout<< "MappedFile<Type>::value :"
706  << " actual average:" << averagePsi
707  << " wanted average:" << wantedAverage
708  << endl;
709  }
710 
711  if (mag(averagePsi) < VSMALL)
712  {
713  // Field too small to scale. Offset instead.
714  const Type offset = wantedAverage - averagePsi;
715  if (debug)
716  {
717  Pout<< "MappedFile<Type>::value :"
718  << " offsetting with:" << offset << endl;
719  }
720  fld += offset;
721  }
722  else
723  {
724  const scalar scale = mag(wantedAverage)/mag(averagePsi);
725 
726  if (debug)
727  {
728  Pout<< "MappedFile<Type>::value :"
729  << " scaling with:" << scale << endl;
730  }
731  fld *= scale;
732  }
733  }
734 
735  // Apply offset to mapped values
736  if (offset_)
737  {
738  fld += offset_->value(x);
739  }
740 
741  if (debug)
742  {
743  Pout<< "MappedFile<Type>::value : set fixedValue to min:" << gMin(fld)
744  << " max:" << gMax(fld)
745  << " avg:" << gAverage(fld) << endl;
746  }
747 
748  return this->transform(tfld);
749 }
750 
751 
752 template<class Type>
755 (
756  const scalar x1,
757  const scalar x2
758 ) const
759 {
761  return nullptr;
762 }
763 
764 
765 template<class Type>
767 (
768  Ostream& os
769 ) const
770 {
771  if (!readerFormat_.empty() && !readerFile_.empty())
772  {
773  os.writeEntry("readerFormat", readerFormat_);
774  os.writeEntry("readerFile", readerFile_);
775  }
776 
777  os.writeEntryIfDifferent
778  (
779  "fieldTable",
780  this->name(),
781  fieldTableName_
782  );
783 
784  if (!pointsName_.empty())
785  {
786  os.writeEntryIfDifferent<word>("points", "points", pointsName_);
787  }
788 
789  if (!mapMethod_.empty() && !mapMethod_.starts_with("planar"))
790  {
791  os.writeEntry("mapMethod", mapMethod_);
792  }
793 
794  if (setAverage_)
795  {
796  os.writeEntry("setAverage", setAverage_);
797  }
798 
799  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
800 
801  if (filterSweeps_ >= 1)
802  {
803  os.writeEntry("filterRadius", filterRadius_);
804  os.writeEntry("filterSweeps", filterSweeps_);
805  }
806 
807  if (offset_)
808  {
809  offset_->writeData(os);
810  }
811 }
812 
813 
814 template<class Type>
816 (
817  Ostream& os
818 ) const
819 {
821 
822  // Check if field name explicitly provided
823  // (e.g. through timeVaryingMapped bc)
824  if (dictConstructed_)
825  {
826  os.writeEntry(this->name(), type());
827 
828  os.beginBlock(word(this->name() + "Coeffs"));
829  writeEntries(os);
830  os.endBlock();
831  }
832  else
833  {
834  // Note that usually dictConstructed = true. The
835  // construct-from-dictionary (dictConstructed_ = false)
836  // should only be used if there is only
837  // a single potential MappedFile in the local scope.
838  writeEntries(os);
839  }
840 }
841 
842 
843 // ************************************************************************* //
Patch value mapping from a set of values stored in a file and a set of unstructured points using the ...
Definition: MappedFile.H:170
static wordList timeNames(const instantList &times)
Helper: extract words of times.
dictionary dict
rDeltaTY field()
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
static Pair< label > findRange(const UList< instant > &times, const scalar timeVal, const label start=-1)
Find lower/upper indices for given time value in list of instances (linear search) continuing after t...
Definition: instant.C:54
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual tmp< Field< Type > > value(const scalar) const
Return MappedFile value.
Definition: MappedFile.C:643
Ignore writing from objectRegistry::writeObject()
const polyPatch const word const word const dictionary & dict
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual tmp< Field< Type > > integrate(const scalar x1, const scalar x2) const
Integrate between two values.
Definition: MappedFile.C:748
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Abstract base class to hold the Field mapping addressing and weights.
Definition: FieldMapper.H:43
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
MeshedSurface< face > meshedSurface
virtual void rmap(const PatchFunction1< Type > &pf1, const labelList &addr)
Reverse map the given PatchFunction1 onto this PatchFunction1.
Definition: MappedFile.C:238
virtual void autoMap(const FieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
const polyPatch const word const word & entryName
#define DebugInfo
Report an information message using Foam::Info.
virtual void writeData(Ostream &os) const
Write in dictionary format.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Type gMax(const FieldField< Field, Type > &f)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
OBJstream os(runTime.globalPath()/outputName)
virtual void writeEntries(Ostream &os) const
Write coefficient entries in dictionary format.
Definition: MappedFile.C:760
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
Definition: surfaceReader.C:36
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:109
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))
List< word > wordList
List of word.
Definition: fileName.H:59
bool starts_with(char c) const
True if string starts with given character (cf. C++20)
Definition: string.H:435
static autoPtr< surfaceReader > New(const word &readType, const fileName &fName, const dictionary &options=dictionary())
Return a reference to the selected surfaceReader.
Definition: surfaceReader.C:71
virtual void autoMap(const FieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
Definition: MappedFile.C:214
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Type gAverage(const FieldField< Field, Type > &f)
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:90
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
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))
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
virtual void rmap(const PatchFunction1< Type > &rhs, const labelList &addr)
Reverse map the given PatchFunction1 onto this PatchFunction1.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: MappedFile.C:809