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