timeVaryingMappedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
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 "AverageField.H"
29 #include "IFstream.H"
30 #include "OSspecific.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
37 (
38  const word& timeName
39 ) const
40 {
41  const fileName fieldFileName
42  (
43  dataDir_/timeName/sampleName_/fieldTableName_
44  );
45 
46  const fileName typeFieldFileName
47  (
48  dataDir_/timeName/sampleName_
49  /pTraits<Type>::typeName + Field<Type>::typeName
50  /fieldTableName_
51  );
52 
53  if (exists(fieldFileName))
54  {
55  return fieldFileName;
56  }
57  else if (exists(typeFieldFileName))
58  {
59  return typeFieldFileName;
60  }
61  else
62  {
64  << "Cannot find field file "
65  << fieldFileName << " " << typeFieldFileName
66  << exit(FatalError);
67 
68  return fileName::null;
69  }
70 }
71 
72 
73 template<class Type>
75 {
76  // Initialise
77  if (mapperPtr_.empty())
78  {
79  // Reread values and interpolate
80  const fileName samplePointsFile(dataDir_/pointsName_);
81 
82  pointField samplePoints((IFstream(samplePointsFile)()));
83 
84  if (debug)
85  {
86  Info<< "timeVaryingMappedFvPatchField :"
87  << " Read " << samplePoints.size() << " sample points from "
88  << samplePointsFile << endl;
89  }
90 
91 
92  // tbd: run-time selection
93  const bool nearestOnly
94  (
95  !mapMethod_.empty()
96  && mapMethod_ != "planarInterpolation"
97  );
98 
99  // Allocate the interpolator
100  mapperPtr_.reset
101  (
102  new pointToPointPlanarInterpolation
103  (
104  samplePoints,
105  patch_.patch().faceCentres(),
106  perturb_,
107  nearestOnly
108  )
109  );
110 
111  // Read the times for which data is available
112  sampleTimes_ = Time::findTimes(dataDir_);
113 
114  if (debug)
115  {
116  Info<< "timeVaryingMappedFvPatchField : In directory "
117  << dataDir_ << " found times "
118  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
119  << endl;
120  }
121  }
122 
123 
124  // Find current time in sampleTimes
125  label lo = -1;
126  label hi = -1;
127 
128  bool foundTime = mapperPtr_().findTime
129  (
130  sampleTimes_,
131  startSampleTime_,
132  time().value(),
133  lo,
134  hi
135  );
136 
137  if (!foundTime)
138  {
140  << "Cannot find starting sampling values for current time "
141  << time().value() << nl
142  << "Have sampling values for times "
143  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
144  << "In directory " << dataDir_ << " of field " << fieldTableName_
145  << exit(FatalError);
146  }
147 
148 
149  // Update sampled data fields.
150 
151  if (lo != startSampleTime_)
152  {
153  startSampleTime_ = lo;
154 
155  if (startSampleTime_ == endSampleTime_)
156  {
157  // No need to reread since are end values
158  if (debug)
159  {
160  Pout<< "checkTable : Setting startValues to (already read) "
161  << dataDir_/sampleTimes_[startSampleTime_].name()
162  << endl;
163  }
164  startSampledValues_ = endSampledValues_;
165  startAverage_ = endAverage_;
166  }
167  else
168  {
169  if (debug)
170  {
171  Pout<< "checkTable : Reading startValues from "
172  << dataDir_/sampleTimes_[lo].name()
173  << endl;
174  }
175 
176  // Reread values and interpolate
177  const fileName valsFile
178  (
179  findFieldFile(sampleTimes_[startSampleTime_].name())
180  );
181 
182  Field<Type> vals;
183 
184  if (setAverage_)
185  {
186  AverageField<Type> avals((IFstream(valsFile)()));
187  vals = avals;
188  startAverage_ = avals.average();
189  }
190  else
191  {
192  (IFstream(valsFile)()) >> vals;
193  }
194 
195  if (vals.size() != mapperPtr_().sourceSize())
196  {
198  << "Number of values (" << vals.size()
199  << ") differs from the number of points ("
200  << mapperPtr_().sourceSize()
201  << ") in file " << valsFile << exit(FatalError);
202  }
203 
204  startSampledValues_ = mapperPtr_().interpolate(vals);
205  }
206  }
207 
208  if (hi != endSampleTime_)
209  {
210  endSampleTime_ = hi;
211 
212  if (endSampleTime_ == -1)
213  {
214  // endTime no longer valid. Might as well clear endValues.
215  if (debug)
216  {
217  Pout<< "checkTable : Clearing endValues" << endl;
218  }
219  endSampledValues_.clear();
220  }
221  else
222  {
223  if (debug)
224  {
225  Pout<< "checkTable : Reading endValues from "
226  << dataDir_/sampleTimes_[endSampleTime_].name()
227  << endl;
228  }
229 
230  // Reread values and interpolate
231  const fileName valsFile
232  (
233  findFieldFile(sampleTimes_[endSampleTime_].name())
234  );
235 
236  Field<Type> vals;
237 
238  if (setAverage_)
239  {
240  AverageField<Type> avals((IFstream(valsFile)()));
241  vals = avals;
242  endAverage_ = avals.average();
243  }
244  else
245  {
246  (IFstream(valsFile)()) >> vals;
247  }
248 
249  if (vals.size() != mapperPtr_().sourceSize())
250  {
252  << "Number of values (" << vals.size()
253  << ") differs from the number of points ("
254  << mapperPtr_().sourceSize()
255  << ") in file " << valsFile << exit(FatalError);
256  }
257 
258  endSampledValues_ = mapperPtr_().interpolate(vals);
259  }
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
266 template<class Type>
268 (
269  const fvPatch& p,
271  const dictionary& dict
272 )
273 :
274  patch_(p),
275  internalField_(iF),
276  fieldTableName_(dict.lookupOrDefault("fieldTable", iF.name())),
277  dataDir_
278  (
279  dict.lookupOrDefault
280  (
281  "dataDir",
282  time().constant()/"boundaryData"/p.name()
283  )
284  ),
285  pointsName_(dict.lookupOrDefault<fileName>("points", "points")),
286  sampleName_(dict.lookupOrDefault("sample", word::null)),
287  setAverage_(dict.lookupOrDefault("setAverage", false)),
288  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
289  mapMethod_
290  (
291  dict.lookupOrDefault<word>
292  (
293  "mapMethod",
294  "planarInterpolation"
295  )
296  ),
297  mapperPtr_(nullptr),
298  sampleTimes_(0),
299  startSampleTime_(-1),
300  startSampledValues_(0),
301  startAverage_(Zero),
302  endSampleTime_(-1),
303  endSampledValues_(0),
304  endAverage_(Zero),
305  offset_()
306 {
307  dataDir_.expand();
308  pointsName_.expand();
309  sampleName_.expand();
310 
311  if (dict.found("offset"))
312  {
313  offset_ =
314  Function1<Type>::New("offset", dimTime, iF.dimensions(), dict);
315  }
316 
317  if
318  (
319  mapMethod_ != "planarInterpolation"
320  && mapMethod_ != "nearest"
321  )
322  {
324  << "mapMethod should be one of 'planarInterpolation'"
325  << ", 'nearest'" << exit(FatalIOError);
326  }
327 }
328 
329 
330 template<class Type>
333 (
335 )
336 :
337  patch_(ptf.patch_),
338  internalField_(ptf.internalField_),
339  fieldTableName_(ptf.fieldTableName_),
340  dataDir_(ptf.dataDir_),
341  pointsName_(ptf.pointsName_),
342  sampleName_(ptf.sampleName_),
343  setAverage_(ptf.setAverage_),
344  perturb_(ptf.perturb_),
345  mapMethod_(ptf.mapMethod_),
346  mapperPtr_(nullptr),
347  sampleTimes_(ptf.sampleTimes_),
348  startSampleTime_(ptf.startSampleTime_),
349  startSampledValues_(ptf.startSampledValues_),
350  startAverage_(ptf.startAverage_),
351  endSampleTime_(ptf.endSampleTime_),
352  endSampledValues_(ptf.endSampledValues_),
353  endAverage_(ptf.endAverage_),
354  offset_(ptf.offset_, false)
355 {}
356 
357 
358 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
359 
360 template<class Type>
362 (
364  const fieldMapper& mapper
365 )
366 {
367  mapper(startSampledValues_, tiptf.startSampledValues_);
368  mapper(endSampledValues_, tiptf.endSampledValues_);
369 
370  // Clear interpolator
371  mapperPtr_.clear();
372  startSampleTime_ = -1;
373  endSampleTime_ = -1;
374 }
375 
376 
377 template<class Type>
379 (
381 )
382 {
383  startSampledValues_.reset(tiptf.startSampledValues_);
384  endSampledValues_.reset(tiptf.endSampledValues_);
385 
386  // Clear interpolator
387  mapperPtr_.clear();
388  startSampleTime_ = -1;
389  endSampleTime_ = -1;
390 }
391 
392 
393 template<class Type>
395 {
396  checkTable();
397 
398  // Interpolate between the sampled data
399 
400  tmp<Field<Type>> tfld(new Field<Type>(patch_.size()));
401  Field<Type>& fld = tfld.ref();
402 
403  Type wantedAverage;
404 
405  if (endSampleTime_ == -1)
406  {
407  // Only start value
408  if (debug)
409  {
410  Pout<< "updateCoeffs : Sampled, non-interpolated values"
411  << " from start time:"
412  << sampleTimes_[startSampleTime_].name() << nl;
413  }
414 
415  fld = startSampledValues_;
416  wantedAverage = startAverage_;
417  }
418  else
419  {
420  const scalar start = sampleTimes_[startSampleTime_].value();
421  const scalar end = sampleTimes_[endSampleTime_].value();
422 
423  const scalar s = (time().value() - start)/(end - start);
424 
425  if (debug)
426  {
427  Pout<< "updateCoeffs : Sampled, interpolated values"
428  << " between start time:"
429  << sampleTimes_[startSampleTime_].name()
430  << " and end time:" << sampleTimes_[endSampleTime_].name()
431  << " with weight:" << s << endl;
432  }
433 
434  fld = (1 - s)*startSampledValues_ + s*endSampledValues_;
435  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
436  }
437 
438  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
439  // offsetting.
440  if (setAverage_)
441  {
442  Type averagePsi = gSum(patch_.magSf()*fld)/gSum(patch_.magSf());
443 
444  if (debug)
445  {
446  Pout<< "updateCoeffs :"
447  << " actual average:" << averagePsi
448  << " wanted average:" << wantedAverage
449  << endl;
450  }
451 
452  if (mag(averagePsi) < vSmall)
453  {
454  // Field too small to scale. Offset instead.
455  const Type offset = wantedAverage - averagePsi;
456  if (debug)
457  {
458  Pout<< "updateCoeffs :"
459  << " offsetting with:" << offset << endl;
460  }
461  fld += offset;
462  }
463  else
464  {
465  const scalar scale = mag(wantedAverage)/mag(averagePsi);
466 
467  if (debug)
468  {
469  Pout<< "updateCoeffs :"
470  << " scaling with:" << scale << endl;
471  }
472  fld *= scale;
473  }
474  }
475 
476  // Apply offset to mapped values
477  if (offset_.valid())
478  {
479  fld += offset_->value(time().value());
480  }
481 
482  if (debug)
483  {
484  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(fld)
485  << " max:" << gMax(fld)
486  << " avg:" << gAverage(fld) << endl;
487  }
488 
489  return tfld;
490 }
491 
492 
493 template<class Type>
495 (
496  Ostream& os
497 ) const
498 {
500  (
501  os,
502  "dataDir",
503  time().constant()/"boundaryData"/patch_.name(),
504  dataDir_
505  );
506 
507  writeEntryIfDifferent(os, "points", fileName("points"), pointsName_);
508  writeEntryIfDifferent(os, "sample", fileName::null, sampleName_);
509  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
510  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
511 
512  writeEntry(os, "fieldTable", fieldTableName_);
513 
515  (
516  os,
517  "mapMethod",
518  word("planarInterpolation"),
519  mapMethod_
520  );
521 
522  if (offset_.valid())
523  {
524  writeEntry
525  (
526  os,
527  time().userUnits(),
528  internalField_.dimensions(),
529  offset_()
530  );
531  }
532 }
533 
534 
535 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
static const fileName null
An empty fileName.
Definition: fileName.H:97
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
string & expand(const bool allowEmpty=false)
Expand initial tildes and all occurrences of environment variables.
Definition: string.C:125
Patch field mapper which interpolates the values from a set of supplied points in space and time.
timeVaryingMappedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
void reset(const timeVaryingMappedFvPatchField &)
Reset the fvPatchField to the given fvPatchField.
tmp< Field< Type > > map()
Return the current mapped patch field.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar e
Elementary charge.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
const doubleScalar e
Definition: doubleScalar.H:106
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
Type gMin(const FieldField< Field, Type > &f)
void offset(label &lst, const label o)
static const char nl
Definition: Ostream.H:266
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p