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-2021 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,
270  const word& fieldName
271 )
272 :
273  patch_(p),
274  fieldTableName_(word::null),
275  dataDir_(time().constant()/"boundaryData"/p.name()),
276  pointsName_("points"),
277  sampleName_(word::null),
278  setAverage_(false),
279  perturb_(0),
280  mapperPtr_(nullptr),
281  sampleTimes_(0),
282  startSampleTime_(-1),
283  startSampledValues_(0),
284  startAverage_(Zero),
285  endSampleTime_(-1),
286  endSampledValues_(0),
287  endAverage_(Zero),
288  offset_()
289 {}
290 
291 
292 template<class Type>
294 (
295  const fvPatch& p,
296  const dictionary& dict,
297  const word& fieldName
298 )
299 :
300  patch_(p),
301  fieldTableName_(dict.lookupOrDefault("fieldTable", fieldName)),
302  dataDir_
303  (
304  dict.lookupOrDefault
305  (
306  "dataDir",
307  time().constant()/"boundaryData"/p.name()
308  )
309  ),
310  pointsName_(dict.lookupOrDefault<fileName>("points", "points")),
311  sampleName_(dict.lookupOrDefault("sample", word::null)),
312  setAverage_(dict.lookupOrDefault("setAverage", false)),
313  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
314  mapMethod_
315  (
316  dict.lookupOrDefault<word>
317  (
318  "mapMethod",
319  "planarInterpolation"
320  )
321  ),
322  mapperPtr_(nullptr),
323  sampleTimes_(0),
324  startSampleTime_(-1),
325  startSampledValues_(0),
326  startAverage_(Zero),
327  endSampleTime_(-1),
328  endSampledValues_(0),
329  endAverage_(Zero),
330  offset_()
331 {
332  dataDir_.expand();
333  pointsName_.expand();
334  sampleName_.expand();
335 
336  if (dict.found("offset"))
337  {
338  offset_ = Function1<Type>::New("offset", dict);
339  }
340 
341  if
342  (
343  mapMethod_ != "planarInterpolation"
344  && mapMethod_ != "nearest"
345  )
346  {
348  (
349  dict
350  ) << "mapMethod should be one of 'planarInterpolation'"
351  << ", 'nearest'" << exit(FatalIOError);
352  }
353 }
354 
355 
356 template<class Type>
359 (
361 )
362 :
363  patch_(ptf.patch_),
364  fieldTableName_(ptf.fieldTableName_),
365  dataDir_(ptf.dataDir_),
366  pointsName_(ptf.pointsName_),
367  sampleName_(ptf.sampleName_),
368  setAverage_(ptf.setAverage_),
369  perturb_(ptf.perturb_),
370  mapMethod_(ptf.mapMethod_),
371  mapperPtr_(nullptr),
372  sampleTimes_(ptf.sampleTimes_),
373  startSampleTime_(ptf.startSampleTime_),
374  startSampledValues_(ptf.startSampledValues_),
375  startAverage_(ptf.startAverage_),
376  endSampleTime_(ptf.endSampleTime_),
377  endSampledValues_(ptf.endSampledValues_),
378  endAverage_(ptf.endAverage_),
379  offset_(ptf.offset_, false)
380 {}
381 
382 
383 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
384 
385 template<class Type>
387 (
388  const fvPatchFieldMapper& m
389 )
390 {
391  if (startSampledValues_.size())
392  {
393  m(startSampledValues_, startSampledValues_);
394  m(endSampledValues_, endSampledValues_);
395  }
396  // Clear interpolator
397  mapperPtr_.clear();
398  startSampleTime_ = -1;
399  endSampleTime_ = -1;
400 }
401 
402 
403 template<class Type>
405 (
407  const labelList& addr
408 )
409 {
410  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
411  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
412 
413  // Clear interpolator
414  mapperPtr_.clear();
415  startSampleTime_ = -1;
416  endSampleTime_ = -1;
417 }
418 
419 
420 template<class Type>
422 {
423  checkTable();
424 
425  // Interpolate between the sampled data
426 
427  tmp<Field<Type>> tfld(new Field<Type>(patch_.size()));
428  Field<Type>& fld = tfld.ref();
429 
430  Type wantedAverage;
431 
432  if (endSampleTime_ == -1)
433  {
434  // Only start value
435  if (debug)
436  {
437  Pout<< "updateCoeffs : Sampled, non-interpolated values"
438  << " from start time:"
439  << sampleTimes_[startSampleTime_].name() << nl;
440  }
441 
442  fld = startSampledValues_;
443  wantedAverage = startAverage_;
444  }
445  else
446  {
447  const scalar start = sampleTimes_[startSampleTime_].value();
448  const scalar end = sampleTimes_[endSampleTime_].value();
449 
450  const scalar s = (time().value() - start)/(end - start);
451 
452  if (debug)
453  {
454  Pout<< "updateCoeffs : Sampled, interpolated values"
455  << " between start time:"
456  << sampleTimes_[startSampleTime_].name()
457  << " and end time:" << sampleTimes_[endSampleTime_].name()
458  << " with weight:" << s << endl;
459  }
460 
461  fld = (1 - s)*startSampledValues_ + s*endSampledValues_;
462  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
463  }
464 
465  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
466  // offsetting.
467  if (setAverage_)
468  {
469  Type averagePsi = gSum(patch_.magSf()*fld)/gSum(patch_.magSf());
470 
471  if (debug)
472  {
473  Pout<< "updateCoeffs :"
474  << " actual average:" << averagePsi
475  << " wanted average:" << wantedAverage
476  << endl;
477  }
478 
479  if (mag(averagePsi) < vSmall)
480  {
481  // Field too small to scale. Offset instead.
482  const Type offset = wantedAverage - averagePsi;
483  if (debug)
484  {
485  Pout<< "updateCoeffs :"
486  << " offsetting with:" << offset << endl;
487  }
488  fld += offset;
489  }
490  else
491  {
492  const scalar scale = mag(wantedAverage)/mag(averagePsi);
493 
494  if (debug)
495  {
496  Pout<< "updateCoeffs :"
497  << " scaling with:" << scale << endl;
498  }
499  fld *= scale;
500  }
501  }
502 
503  // Apply offset to mapped values
504  if (offset_.valid())
505  {
506  const scalar t = time().timeOutputValue();
507  fld += offset_->value(t);
508  }
509 
510  if (debug)
511  {
512  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(fld)
513  << " max:" << gMax(fld)
514  << " avg:" << gAverage(fld) << endl;
515  }
516 
517  return tfld;
518 }
519 
520 
521 template<class Type>
523 (
524  Ostream& os
525 ) const
526 {
528  (
529  os,
530  "dataDir",
531  time().constant()/"boundaryData"/patch_.name(),
532  dataDir_
533  );
534 
535  writeEntryIfDifferent(os, "points", fileName("points"), pointsName_);
536  writeEntryIfDifferent(os, "sample", fileName::null, sampleName_);
537  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
538  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
539 
540  writeEntry(os, "fieldTable", fieldTableName_);
541 
543  (
544  os,
545  "mapMethod",
546  word("planarInterpolation"),
547  mapMethod_
548  );
549 
550  if (offset_.valid())
551  {
552  writeEntry(os, offset_());
553  }
554 }
555 
556 
557 // ************************************************************************* //
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
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
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
A class for handling file names.
Definition: fileName.H:79
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Type gMin(const FieldField< Field, Type > &f)
tmp< Field< Type > > map()
Return the current mapped patch field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Patch field mapper which interpolates the values from a set of supplied points in space and time...
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
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))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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))
Type gSum(const FieldField< Field, Type > &f)
void rmap(const timeVaryingMappedFvPatchField< Type > &, const labelList &)
Reverse map the given timeVaryingMappedFvPatchField.
Pre-declare SubField and related Field type.
Definition: Field.H:56
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
Foam::fvPatchFieldMapper.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
timeVaryingMappedFvPatchField(const fvPatch &, const word &fieldName)
Construct from patch and field name.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const word & name() const
Return name.
Definition: fvPatch.H:144
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
IOerror FatalIOError