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-2019 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 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const word& timeName
38 ) const
39 {
40  const fileName fieldFileName
41  (
42  dataDir_/timeName/sampleName_/fieldTableName_
43  );
44 
45  const fileName typeFieldFileName
46  (
47  dataDir_/timeName/sampleName_
48  /pTraits<Type>::typeName + Field<Type>::typeName
49  /fieldTableName_
50  );
51 
52  if (exists(fieldFileName))
53  {
54  return fieldFileName;
55  }
56  else if (exists(typeFieldFileName))
57  {
58  return typeFieldFileName;
59  }
60  else
61  {
63  << "Cannot find field file "
64  << fieldFileName << " " << typeFieldFileName
65  << exit(FatalError);
66 
67  return fileName::null;
68  }
69 }
70 
71 
72 template<class Type>
74 {
75  // Initialise
76  if (mapperPtr_.empty())
77  {
78  // Reread values and interpolate
79  const fileName samplePointsFile(dataDir_/pointsName_);
80 
81  pointField samplePoints((IFstream(samplePointsFile)()));
82 
83  if (debug)
84  {
85  Info<< "timeVaryingMappedFvPatchField :"
86  << " Read " << samplePoints.size() << " sample points from "
87  << samplePointsFile << endl;
88  }
89 
90 
91  // tbd: run-time selection
92  const bool nearestOnly
93  (
94  !mapMethod_.empty()
95  && mapMethod_ != "planarInterpolation"
96  );
97 
98  // Allocate the interpolator
99  mapperPtr_.reset
100  (
101  new pointToPointPlanarInterpolation
102  (
103  samplePoints,
104  patch_.patch().faceCentres(),
105  perturb_,
106  nearestOnly
107  )
108  );
109 
110  // Read the times for which data is available
111  sampleTimes_ = Time::findTimes(dataDir_);
112 
113  if (debug)
114  {
115  Info<< "timeVaryingMappedFvPatchField : In directory "
116  << dataDir_ << " found times "
117  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
118  << endl;
119  }
120  }
121 
122 
123  // Find current time in sampleTimes
124  label lo = -1;
125  label hi = -1;
126 
127  bool foundTime = mapperPtr_().findTime
128  (
129  sampleTimes_,
130  startSampleTime_,
131  time().value(),
132  lo,
133  hi
134  );
135 
136  if (!foundTime)
137  {
139  << "Cannot find starting sampling values for current time "
140  << time().value() << nl
141  << "Have sampling values for times "
142  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
143  << "In directory " << dataDir_ << " of field " << fieldTableName_
144  << exit(FatalError);
145  }
146 
147 
148  // Update sampled data fields.
149 
150  if (lo != startSampleTime_)
151  {
152  startSampleTime_ = lo;
153 
154  if (startSampleTime_ == endSampleTime_)
155  {
156  // No need to reread since are end values
157  if (debug)
158  {
159  Pout<< "checkTable : Setting startValues to (already read) "
160  << dataDir_/sampleTimes_[startSampleTime_].name()
161  << endl;
162  }
163  startSampledValues_ = endSampledValues_;
164  startAverage_ = endAverage_;
165  }
166  else
167  {
168  if (debug)
169  {
170  Pout<< "checkTable : Reading startValues from "
171  << dataDir_/sampleTimes_[lo].name()
172  << endl;
173  }
174 
175  // Reread values and interpolate
176  const fileName valsFile
177  (
178  findFieldFile(sampleTimes_[startSampleTime_].name())
179  );
180 
181  Field<Type> vals;
182 
183  if (setAverage_)
184  {
185  AverageField<Type> avals((IFstream(valsFile)()));
186  vals = avals;
187  startAverage_ = avals.average();
188  }
189  else
190  {
191  IFstream(valsFile)() >> vals;
192  }
193 
194  if (vals.size() != mapperPtr_().sourceSize())
195  {
197  << "Number of values (" << vals.size()
198  << ") differs from the number of points ("
199  << mapperPtr_().sourceSize()
200  << ") in file " << valsFile << exit(FatalError);
201  }
202 
203  startSampledValues_ = mapperPtr_().interpolate(vals);
204  }
205  }
206 
207  if (hi != endSampleTime_)
208  {
209  endSampleTime_ = hi;
210 
211  if (endSampleTime_ == -1)
212  {
213  // endTime no longer valid. Might as well clear endValues.
214  if (debug)
215  {
216  Pout<< "checkTable : Clearing endValues" << endl;
217  }
218  endSampledValues_.clear();
219  }
220  else
221  {
222  if (debug)
223  {
224  Pout<< "checkTable : Reading endValues from "
225  << dataDir_/sampleTimes_[endSampleTime_].name()
226  << endl;
227  }
228 
229  // Reread values and interpolate
230  const fileName valsFile
231  (
232  findFieldFile(sampleTimes_[endSampleTime_].name())
233  );
234 
235  Field<Type> vals;
236 
237  if (setAverage_)
238  {
239  AverageField<Type> avals((IFstream(valsFile)()));
240  vals = avals;
241  endAverage_ = avals.average();
242  }
243  else
244  {
245  IFstream(valsFile)() >> vals;
246  }
247 
248  if (vals.size() != mapperPtr_().sourceSize())
249  {
251  << "Number of values (" << vals.size()
252  << ") differs from the number of points ("
253  << mapperPtr_().sourceSize()
254  << ") in file " << valsFile << exit(FatalError);
255  }
256 
257  endSampledValues_ = mapperPtr_().interpolate(vals);
258  }
259  }
260 }
261 
262 
263 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
264 
265 template<class Type>
267 (
268  const fvPatch& p,
269  const word& fieldName
270 )
271 :
272  patch_(p),
273  fieldTableName_(word::null),
274  dataDir_(time().constant()/"boundaryData"/p.name()),
275  pointsName_("points"),
276  sampleName_(word::null),
277  setAverage_(false),
278  perturb_(0),
279  mapperPtr_(nullptr),
280  sampleTimes_(0),
281  startSampleTime_(-1),
282  startSampledValues_(0),
283  startAverage_(Zero),
284  endSampleTime_(-1),
285  endSampledValues_(0),
286  endAverage_(Zero),
287  offset_()
288 {}
289 
290 
291 template<class Type>
293 (
294  const fvPatch& p,
295  const dictionary& dict,
296  const word& fieldName
297 )
298 :
299  patch_(p),
300  fieldTableName_(dict.lookupOrDefault("fieldTable", fieldName)),
301  dataDir_
302  (
303  dict.lookupOrDefault
304  (
305  "dataDir",
306  time().constant()/"boundaryData"/p.name()
307  )
308  ),
309  pointsName_(dict.lookupOrDefault<fileName>("points", "points")),
310  sampleName_(dict.lookupOrDefault("sample", word::null)),
311  setAverage_(dict.lookupOrDefault("setAverage", false)),
312  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
313  mapMethod_
314  (
315  dict.lookupOrDefault<word>
316  (
317  "mapMethod",
318  "planarInterpolation"
319  )
320  ),
321  mapperPtr_(nullptr),
322  sampleTimes_(0),
323  startSampleTime_(-1),
324  startSampledValues_(0),
325  startAverage_(Zero),
326  endSampleTime_(-1),
327  endSampledValues_(0),
328  endAverage_(Zero),
329  offset_()
330 {
331  dataDir_.expand();
332  pointsName_.expand();
333  sampleName_.expand();
334 
335  if (dict.found("offset"))
336  {
337  offset_ = Function1<Type>::New("offset", dict);
338  }
339 
340  if
341  (
342  mapMethod_ != "planarInterpolation"
343  && mapMethod_ != "nearest"
344  )
345  {
347  (
348  dict
349  ) << "mapMethod should be one of 'planarInterpolation'"
350  << ", 'nearest'" << exit(FatalIOError);
351  }
352 }
353 
354 
355 template<class Type>
358 (
360 )
361 :
362  patch_(ptf.patch_),
363  fieldTableName_(ptf.fieldTableName_),
364  dataDir_(ptf.dataDir_),
365  pointsName_(ptf.pointsName_),
366  sampleName_(ptf.sampleName_),
367  setAverage_(ptf.setAverage_),
368  perturb_(ptf.perturb_),
369  mapMethod_(ptf.mapMethod_),
370  mapperPtr_(nullptr),
371  sampleTimes_(ptf.sampleTimes_),
372  startSampleTime_(ptf.startSampleTime_),
373  startSampledValues_(ptf.startSampledValues_),
374  startAverage_(ptf.startAverage_),
375  endSampleTime_(ptf.endSampleTime_),
376  endSampledValues_(ptf.endSampledValues_),
377  endAverage_(ptf.endAverage_),
378  offset_(ptf.offset_, false)
379 {}
380 
381 
382 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
383 
384 template<class Type>
386 (
387  const fvPatchFieldMapper& m
388 )
389 {
390  if (startSampledValues_.size())
391  {
392  m(startSampledValues_, startSampledValues_);
393  m(endSampledValues_, endSampledValues_);
394  }
395  // Clear interpolator
396  mapperPtr_.clear();
397  startSampleTime_ = -1;
398  endSampleTime_ = -1;
399 }
400 
401 
402 template<class Type>
404 (
406  const labelList& addr
407 )
408 {
409  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
410  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
411 
412  // Clear interpolator
413  mapperPtr_.clear();
414  startSampleTime_ = -1;
415  endSampleTime_ = -1;
416 }
417 
418 
419 template<class Type>
421 {
422  checkTable();
423 
424  // Interpolate between the sampled data
425 
426  tmp<Field<Type>> tfld(new Field<Type>(patch_.size()));
427  Field<Type>& fld = tfld.ref();
428 
429  Type wantedAverage;
430 
431  if (endSampleTime_ == -1)
432  {
433  // Only start value
434  if (debug)
435  {
436  Pout<< "updateCoeffs : Sampled, non-interpolated values"
437  << " from start time:"
438  << sampleTimes_[startSampleTime_].name() << nl;
439  }
440 
441  fld = startSampledValues_;
442  wantedAverage = startAverage_;
443  }
444  else
445  {
446  const scalar start = sampleTimes_[startSampleTime_].value();
447  const scalar end = sampleTimes_[endSampleTime_].value();
448 
449  const scalar s = (time().value() - start)/(end - start);
450 
451  if (debug)
452  {
453  Pout<< "updateCoeffs : Sampled, interpolated values"
454  << " between start time:"
455  << sampleTimes_[startSampleTime_].name()
456  << " and end time:" << sampleTimes_[endSampleTime_].name()
457  << " with weight:" << s << endl;
458  }
459 
460  fld = (1 - s)*startSampledValues_ + s*endSampledValues_;
461  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
462  }
463 
464  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
465  // offsetting.
466  if (setAverage_)
467  {
468  Type averagePsi = gSum(patch_.magSf()*fld)/gSum(patch_.magSf());
469 
470  if (debug)
471  {
472  Pout<< "updateCoeffs :"
473  << " actual average:" << averagePsi
474  << " wanted average:" << wantedAverage
475  << endl;
476  }
477 
478  if (mag(averagePsi) < vSmall)
479  {
480  // Field too small to scale. Offset instead.
481  const Type offset = wantedAverage - averagePsi;
482  if (debug)
483  {
484  Pout<< "updateCoeffs :"
485  << " offsetting with:" << offset << endl;
486  }
487  fld += offset;
488  }
489  else
490  {
491  const scalar scale = mag(wantedAverage)/mag(averagePsi);
492 
493  if (debug)
494  {
495  Pout<< "updateCoeffs :"
496  << " scaling with:" << scale << endl;
497  }
498  fld *= scale;
499  }
500  }
501 
502  // Apply offset to mapped values
503  if (offset_.valid())
504  {
505  const scalar t = time().timeOutputValue();
506  fld += offset_->value(t);
507  }
508 
509  if (debug)
510  {
511  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(fld)
512  << " max:" << gMax(fld)
513  << " avg:" << gAverage(fld) << endl;
514  }
515 
516  return tfld;
517 }
518 
519 
520 template<class Type>
522 (
523  Ostream& os
524 ) const
525 {
527  (
528  os,
529  "dataDir",
530  time().constant()/"boundaryData"/patch_.name(),
531  dataDir_
532  );
533 
534  writeEntryIfDifferent(os, "points", fileName("points"), pointsName_);
535  writeEntryIfDifferent(os, "sample", fileName::null, sampleName_);
536  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
537  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
538 
539  writeEntry(os, "fieldTable", fieldTableName_);
540 
542  (
543  os,
544  "mapMethod",
545  word("planarInterpolation"),
546  mapMethod_
547  );
548 
549  if (offset_.valid())
550  {
551  writeEntry(os, offset_());
552  }
553 }
554 
555 
556 // ************************************************************************* //
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:667
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:174
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:61
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...
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
const dimensionedScalar & e
Elementary charge.
Definition: doubleScalar.H:105
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:297
Foam::fvPatchFieldMapper.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:183
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:331
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:143
A class for managing temporary objects.
Definition: PtrList.H:53
IOerror FatalIOError