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& redirectType,
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  setAverage_(dict.getOrDefault("setAverage", false)),
46  fieldTableName_(entryName),
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_(Function1<Type>::NewIfPresent("offset", dict))
66 {
67  if
68  (
69  mapMethod_ != "planarInterpolation"
70  && mapMethod_ != "nearest"
71  )
72  {
74  << "mapMethod should be one of 'planarInterpolation'"
75  << ", 'nearest'" << exit(FatalIOError);
76  }
77 
78  dict.readIfPresent("fieldTable", fieldTableName_);
79 }
80 
81 
82 template<class Type>
84 (
85  const polyPatch& pp,
86  const word& entryName,
87  const dictionary& dict,
88  const word& fieldTableName,
89  const bool faceValues
90 )
91 :
92  PatchFunction1<Type>(pp, entryName, dict, faceValues),
93  dictConstructed_(false),
94  setAverage_(dict.getOrDefault("setAverage", false)),
95  fieldTableName_(fieldTableName),
96  perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
97  pointsName_(dict.getOrDefault<word>("points", "points")),
98  mapMethod_
99  (
100  dict.getOrDefault<word>
101  (
102  "mapMethod",
103  "planarInterpolation"
104  )
105  ),
106  mapperPtr_(nullptr),
107  sampleTimes_(0),
108  startSampleTime_(-1),
109  startSampledValues_(0),
110  startAverage_(Zero),
111  endSampleTime_(-1),
112  endSampledValues_(0),
113  endAverage_(Zero),
114  offset_(Function1<Type>::NewIfPresent("offset", dict))
115 {
116  if
117  (
118  mapMethod_ != "planarInterpolation"
119  && mapMethod_ != "nearest"
120  )
121  {
123  << "mapMethod should be one of 'planarInterpolation'"
124  << ", 'nearest'" << exit(FatalIOError);
125  }
126 }
127 
128 
129 template<class Type>
131 (
132  const MappedFile<Type>& rhs
133 )
134 :
135  MappedFile<Type>(rhs, rhs.patch())
136 {}
137 
138 
139 template<class Type>
141 (
142  const MappedFile<Type>& rhs,
143  const polyPatch& pp
144 )
145 :
146  PatchFunction1<Type>(rhs, pp),
147  dictConstructed_(rhs.dictConstructed_),
148  setAverage_(rhs.setAverage_),
149  fieldTableName_(rhs.fieldTableName_),
150  perturb_(rhs.perturb_),
151  pointsName_(rhs.pointsName_),
152  mapMethod_(rhs.mapMethod_),
153  mapperPtr_(rhs.mapperPtr_.clone()),
154  sampleTimes_(rhs.sampleTimes_),
155  startSampleTime_(rhs.startSampleTime_),
156  startSampledValues_(rhs.startSampledValues_),
157  startAverage_(rhs.startAverage_),
158  endSampleTime_(rhs.endSampleTime_),
159  endSampledValues_(rhs.endSampledValues_),
160  endAverage_(rhs.endAverage_),
161  offset_(rhs.offset_.clone())
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
167 template<class Type>
169 (
170  const FieldMapper& mapper
171 )
172 {
174 
175  if (startSampledValues_.size())
176  {
177  startSampledValues_.autoMap(mapper);
178  endSampledValues_.autoMap(mapper);
179  }
180  // Clear interpolator
181  mapperPtr_.clear();
182  startSampleTime_ = -1;
183  endSampleTime_ = -1;
184 }
185 
186 
187 template<class Type>
189 (
190  const PatchFunction1<Type>& pf1,
191  const labelList& addr
192 )
193 {
194  PatchFunction1<Type>::rmap(pf1, addr);
195 
197  refCast<const PatchFunction1Types::MappedFile<Type>>(pf1);
198 
199  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
200  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
201 
202  // Clear interpolator
203  mapperPtr_.clear();
204  startSampleTime_ = -1;
205  endSampleTime_ = -1;
206 }
207 
208 
209 template<class Type>
211 (
212  const scalar t
213 ) const
214 {
215  const polyMesh& mesh = this->patch_.boundaryMesh().mesh();
216  const Time& time = mesh.time();
217 
218  // Initialise
219  if (!mapperPtr_)
220  {
221  // Reread values and interpolate
222  const fileName samplePointsFile
223  (
224  time.globalPath()
225  /time.constant()
226  /"boundaryData"
227  /this->patch_.name()
228  /pointsName_
229  );
230 
231  IOobject io
232  (
233  samplePointsFile, // absolute path
234  time,
235  IOobject::MUST_READ,
236  IOobject::NO_WRITE,
237  false, // no need to register
238  true // is global object (currently not used)
239  );
240 
241  // Read data
242  const rawIOField<point> samplePoints(io, false);
243 
244  DebugInfo
245  << "Read " << samplePoints.size() << " sample points from "
246  << samplePointsFile << endl;
247 
248 
249  // tbd: run-time selection
250  bool nearestOnly =
251  (
252  !mapMethod_.empty()
253  && mapMethod_ != "planarInterpolation"
254  );
255 
256  // Allocate the interpolator
257  if (this->faceValues())
258  {
259  mapperPtr_.reset
260  (
262  (
263  samplePoints,
264  this->localPosition(this->patch_.faceCentres()),
265  perturb_,
266  nearestOnly
267  )
268  );
269  }
270  else
271  {
272  mapperPtr_.reset
273  (
275  (
276  samplePoints,
277  this->localPosition(this->patch_.localPoints()),
278  perturb_,
279  nearestOnly
280  )
281  );
282  }
283 
284 
285  // Read the times for which data is available
286  const fileName samplePointsDir = samplePointsFile.path();
287  sampleTimes_ = Time::findTimes(samplePointsDir);
288 
289  DebugInfo
290  << "In directory "
291  << samplePointsDir << " found times "
292  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
293  << endl;
294  }
295 
296 
297  // Find current time in sampleTimes
298  label lo = -1;
299  label hi = -1;
300 
301  bool foundTime = mapperPtr_().findTime
302  (
303  sampleTimes_,
304  startSampleTime_,
305  t, //mesh.time().value(),
306  lo,
307  hi
308  );
309 
310  if (!foundTime)
311  {
313  << "Cannot find starting sampling values for index "
314  << t << nl
315  << "Have sampling values for "
316  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
317  << "In directory "
318  << time.constant()/"boundaryData"/this->patch_.name()
319  << "\n on patch " << this->patch_.name()
320  << " of field " << fieldTableName_
321  << exit(FatalError);
322  }
323 
324 
325  // Update sampled data fields.
326 
327  if (lo != startSampleTime_)
328  {
329  startSampleTime_ = lo;
330 
331  if (startSampleTime_ == endSampleTime_)
332  {
333  // No need to reread since are end values
334  if (debug)
335  {
336  Pout<< "checkTable : Setting startValues to (already read) "
337  << "boundaryData"
338  /this->patch_.name()
339  /sampleTimes_[startSampleTime_].name()
340  << endl;
341  }
342  startSampledValues_ = endSampledValues_;
343  startAverage_ = endAverage_;
344  }
345  else
346  {
347  if (debug)
348  {
349  Pout<< "checkTable : Reading startValues from "
350  << "boundaryData"
351  /this->patch_.name()
352  /sampleTimes_[lo].name()
353  << endl;
354  }
355 
356 
357  // Reread values and interpolate
358  const fileName valsFile
359  (
360  time.globalPath()
361  /time.constant()
362  /"boundaryData"
363  /this->patch_.name()
364  /sampleTimes_[startSampleTime_].name()
365  /fieldTableName_
366  );
367 
368  IOobject io
369  (
370  valsFile, // absolute path
371  time,
372  IOobject::MUST_READ,
373  IOobject::NO_WRITE,
374  false, // no need to register
375  true // is global object (currently not used)
376  );
377 
378  const rawIOField<Type> vals(io, setAverage_);
379  if (setAverage_)
380  {
381  startAverage_ = vals.average();
382  }
383 
384  if (vals.size() != mapperPtr_().sourceSize())
385  {
387  << "Number of values (" << vals.size()
388  << ") differs from the number of points ("
389  << mapperPtr_().sourceSize()
390  << ") in file " << valsFile << exit(FatalError);
391  }
392 
393  startSampledValues_ = mapperPtr_().interpolate(vals);
394  }
395  }
396 
397  if (hi != endSampleTime_)
398  {
399  endSampleTime_ = hi;
400 
401  if (endSampleTime_ == -1)
402  {
403  // endTime no longer valid. Might as well clear endValues.
404  if (debug)
405  {
406  Pout<< "checkTable : Clearing endValues" << endl;
407  }
408  endSampledValues_.clear();
409  }
410  else
411  {
412  if (debug)
413  {
414  Pout<< "checkTable : Reading endValues from "
415  << "boundaryData"
416  /this->patch_.name()
417  /sampleTimes_[endSampleTime_].name()
418  << endl;
419  }
420 
421  // Reread values and interpolate
422  fileName valsFile
423  (
424  time.globalPath()
425  /time.constant()
426  /"boundaryData"
427  /this->patch_.name()
428  /sampleTimes_[endSampleTime_].name()
429  /fieldTableName_
430  );
431 
432  IOobject io
433  (
434  valsFile, // absolute path
435  time,
436  IOobject::MUST_READ,
437  IOobject::NO_WRITE,
438  false, // no need to register
439  true // is global object (currently not used)
440  );
441 
442  const rawIOField<Type> vals(io, setAverage_);
443  if (setAverage_)
444  {
445  endAverage_ = vals.average();
446  }
447 
448  if (vals.size() != mapperPtr_().sourceSize())
449  {
451  << "Number of values (" << vals.size()
452  << ") differs from the number of points ("
453  << mapperPtr_().sourceSize()
454  << ") in file " << valsFile << exit(FatalError);
455  }
456 
457  endSampledValues_ = mapperPtr_().interpolate(vals);
458  }
459  }
460 }
461 
462 
463 template<class Type>
466 (
467  const scalar x
468 ) const
469 {
470  checkTable(x);
471 
472  auto tfld = tmp<Field<Type>>::New(startSampledValues_.size());
473  auto& fld = tfld.ref();
474  Type wantedAverage;
475 
476  if (endSampleTime_ == -1)
477  {
478  // Only start value
479  if (debug)
480  {
481  Pout<< "MappedFile<Type>::value : Sampled, non-interpolated values"
482  << " from start time:"
483  << sampleTimes_[startSampleTime_].name() << nl;
484  }
485 
486  fld = startSampledValues_;
487  wantedAverage = startAverage_;
488  }
489  else
490  {
491  scalar start = sampleTimes_[startSampleTime_].value();
492  scalar end = sampleTimes_[endSampleTime_].value();
493 
494  scalar s = (x - start)/(end - start);
495 
496  if (debug)
497  {
498  Pout<< "MappedFile<Type>::value : Sampled, interpolated values"
499  << " between start time:"
500  << sampleTimes_[startSampleTime_].name()
501  << " and end time:" << sampleTimes_[endSampleTime_].name()
502  << " with weight:" << s << endl;
503  }
504 
505  fld = ((1 - s)*startSampledValues_ + s*endSampledValues_);
506  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
507  }
508 
509  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
510  // offsetting.
511  if (setAverage_)
512  {
513  Type averagePsi;
514  if (this->faceValues())
515  {
516  const scalarField magSf(mag(this->patch_.faceAreas()));
517  averagePsi = gSum(magSf*fld)/gSum(magSf);
518  }
519  else
520  {
521  averagePsi = gAverage(fld);
522  }
523 
524  if (debug)
525  {
526  Pout<< "MappedFile<Type>::value :"
527  << " actual average:" << averagePsi
528  << " wanted average:" << wantedAverage
529  << endl;
530  }
531 
532  if (mag(averagePsi) < VSMALL)
533  {
534  // Field too small to scale. Offset instead.
535  const Type offset = wantedAverage - averagePsi;
536  if (debug)
537  {
538  Pout<< "MappedFile<Type>::value :"
539  << " offsetting with:" << offset << endl;
540  }
541  fld += offset;
542  }
543  else
544  {
545  const scalar scale = mag(wantedAverage)/mag(averagePsi);
546 
547  if (debug)
548  {
549  Pout<< "MappedFile<Type>::value :"
550  << " scaling with:" << scale << endl;
551  }
552  fld *= scale;
553  }
554  }
555 
556  // Apply offset to mapped values
557  if (offset_)
558  {
559  fld += offset_->value(x);
560  }
561 
562  if (debug)
563  {
564  Pout<< "MappedFile<Type>::value : set fixedValue to min:" << gMin(fld)
565  << " max:" << gMax(fld)
566  << " avg:" << gAverage(fld) << endl;
567  }
568 
569  return this->transform(tfld);
570 }
571 
572 
573 template<class Type>
576 (
577  const scalar x1,
578  const scalar x2
579 ) const
580 {
582  return nullptr;
583 }
584 
585 
586 template<class Type>
588 (
589  Ostream& os
590 ) const
591 {
593 
594  // Check if field name explicitly provided
595  // (e.g. through timeVaryingMapped bc)
596  if (dictConstructed_)
597  {
598  os.writeEntry(this->name(), type());
599 
601  (
602  "fieldTable",
603  this->name(),
604  fieldTableName_
605  );
606  }
607 
608  if (setAverage_)
609  {
610  os.writeEntry("setAverage", setAverage_);
611  }
612 
613  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
614 
615  os.writeEntryIfDifferent<word>("points", "points", pointsName_);
616 
618  (
619  "mapMethod",
620  "planarInterpolation",
621  mapMethod_
622  );
623 
624  if (offset_)
625  {
626  offset_->writeData(os);
627  }
628 }
629 
630 
631 // ************************************************************************* //
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:169
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:61
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:588
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
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
polyMesh.H
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
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:445
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:68
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:576
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:59
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:381
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
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:189
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
Foam::PatchFunction1Types::MappedFile::MappedFile
MappedFile(const polyPatch &pp, const word &redirectType, const word &entryName, const dictionary &dict, const bool faceValues=true)
Construct from entry name and dictionary.
Definition: MappedFile.C:35
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
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:466