timeVaryingMappedFixedValuePointPatchField.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) 2012-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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const pointPatch& p,
39 )
40 :
42  fieldTableName_(iF.name()),
43  setAverage_(false),
44  perturb_(0),
45  mapperPtr_(nullptr),
46  sampleTimes_(0),
47  startSampleTime_(-1),
48  startSampledValues_(0),
49  startAverage_(Zero),
50  endSampleTime_(-1),
51  endSampledValues_(0),
52  endAverage_(Zero),
53  offset_()
54 {}
55 
56 
57 template<class Type>
60 (
61  const pointPatch& p,
63  const dictionary& dict
64 )
65 :
67  fieldTableName_(iF.name()),
68  setAverage_(dict.lookupOrDefault("setAverage", false)),
69  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
70  mapMethod_
71  (
72  dict.lookupOrDefault<word>
73  (
74  "mapMethod",
75  "planarInterpolation"
76  )
77  ),
78  mapperPtr_(nullptr),
79  sampleTimes_(0),
80  startSampleTime_(-1),
81  startSampledValues_(0),
82  startAverage_(Zero),
83  endSampleTime_(-1),
84  endSampledValues_(0),
85  endAverage_(Zero),
86  offset_()
87 {
88  if (dict.found("offset"))
89  {
90  offset_ = Function1<Type>::New("offset", dict);
91  }
92 
93  if
94  (
95  mapMethod_ != "planarInterpolation"
96  && mapMethod_ != "nearest"
97  )
98  {
100  (
101  dict
102  ) << "mapMethod should be one of 'planarInterpolation'"
103  << ", 'nearest'" << exit(FatalIOError);
104  }
105 
106  dict.readIfPresent("fieldTableName", fieldTableName_);
107 
108  if (dict.found("value"))
109  {
111  (
112  Field<Type>("value", dict, p.size())
113  );
114  }
115  else
116  {
117  // Note: use evaluate to do updateCoeffs followed by a reset
118  // of the pointPatchField::updated_ flag. This is
119  // so if first use is in the next time step it retriggers
120  // a new update.
121  pointPatchField<Type>::evaluate(Pstream::commsTypes::blocking);
122  }
123 }
124 
125 
126 template<class Type>
129 (
131  const pointPatch& p,
133  const pointPatchFieldMapper& mapper
134 )
135 :
136  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
137  fieldTableName_(ptf.fieldTableName_),
138  setAverage_(ptf.setAverage_),
139  perturb_(ptf.perturb_),
140  mapMethod_(ptf.mapMethod_),
141  mapperPtr_(nullptr),
142  sampleTimes_(0),
143  startSampleTime_(-1),
144  startSampledValues_(0),
145  startAverage_(Zero),
146  endSampleTime_(-1),
147  endSampledValues_(0),
148  endAverage_(Zero),
149  offset_(ptf.offset_, false)
150 {}
151 
152 
153 template<class Type>
156 (
159 )
160 :
162  fieldTableName_(ptf.fieldTableName_),
163  setAverage_(ptf.setAverage_),
164  perturb_(ptf.perturb_),
165  mapMethod_(ptf.mapMethod_),
166  mapperPtr_(ptf.mapperPtr_),
167  sampleTimes_(ptf.sampleTimes_),
168  startSampleTime_(ptf.startSampleTime_),
169  startSampledValues_(ptf.startSampledValues_),
170  startAverage_(ptf.startAverage_),
171  endSampleTime_(ptf.endSampleTime_),
172  endSampledValues_(ptf.endSampledValues_),
173  endAverage_(ptf.endAverage_),
174  offset_(ptf.offset_, false)
175 {}
176 
177 
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179 
180 template<class Type>
182 (
183  const pointPatchFieldMapper& m
184 )
185 {
187  if (startSampledValues_.size())
188  {
189  m(startSampledValues_, startSampledValues_);
190  m(endSampledValues_, endSampledValues_);
191  }
192  // Clear interpolator
193  mapperPtr_.clear();
194  startSampleTime_ = -1;
195  endSampleTime_ = -1;
196 }
197 
198 
199 template<class Type>
201 (
202  const pointPatchField<Type>& ptf,
203  const labelList& addr
204 )
205 {
207 
209  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
210 
211  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
212  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
213 
214  // Clear interpolator
215  mapperPtr_.clear();
216  startSampleTime_ = -1;
217  endSampleTime_ = -1;
218 }
219 
220 
221 template<class Type>
223 {
224  // Initialise
225  if (startSampleTime_ == -1 && endSampleTime_ == -1)
226  {
227  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
228 
229  // Read the initial point position
230  pointField meshPts;
231 
232  if (pMesh.pointsInstance() == pMesh.facesInstance())
233  {
234  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
235  }
236  else
237  {
238  // Load points from facesInstance
239  if (debug)
240  {
241  Info<< "Reloading points0 from " << pMesh.facesInstance()
242  << endl;
243  }
244 
245  pointIOField points0
246  (
247  IOobject
248  (
249  "points",
250  pMesh.facesInstance(),
251  polyMesh::meshSubDir,
252  pMesh,
253  IOobject::MUST_READ,
254  IOobject::NO_WRITE,
255  false
256  )
257  );
258  meshPts = pointField(points0, this->patch().meshPoints());
259  }
260 
261  // Reread values and interpolate
262  fileName samplePointsFile
263  (
264  this->db().time().constant()
265  /"boundaryData"
266  /this->patch().name()
267  /"points"
268  );
269 
270  pointField samplePoints((IFstream(samplePointsFile)()));
271 
272  // tbd: run-time selection
273  bool nearestOnly =
274  (
275  !mapMethod_.empty()
276  && mapMethod_ != "planarInterpolation"
277  );
278 
279  // Allocate the interpolator
280  mapperPtr_.reset
281  (
283  (
284  samplePoints,
285  meshPts,
286  perturb_,
287  nearestOnly
288  )
289  );
290 
291  // Read the times for which data is available
292 
293  const fileName samplePointsDir = samplePointsFile.path();
294  sampleTimes_ = Time::findTimes(samplePointsDir);
295 
296  if (debug)
297  {
298  Info<< "timeVaryingMappedFixedValuePointPatchField : In directory "
299  << samplePointsDir << " found times "
300  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
301  << endl;
302  }
303  }
304 
305  // Find current time in sampleTimes
306  label lo = -1;
307  label hi = -1;
308 
309  bool foundTime = mapperPtr_().findTime
310  (
311  sampleTimes_,
312  startSampleTime_,
313  this->db().time().value(),
314  lo,
315  hi
316  );
317 
318  if (!foundTime)
319  {
321  << "Cannot find starting sampling values for current time "
322  << this->db().time().value() << nl
323  << "Have sampling values for times "
324  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
325  << "In directory "
326  << this->db().time().constant()/"boundaryData"/this->patch().name()
327  << "\n on patch " << this->patch().name()
328  << " of field " << fieldTableName_
329  << exit(FatalError);
330  }
331 
332 
333  // Update sampled data fields.
334 
335  if (lo != startSampleTime_)
336  {
337  startSampleTime_ = lo;
338 
339  if (startSampleTime_ == endSampleTime_)
340  {
341  // No need to reread since are end values
342  if (debug)
343  {
344  Pout<< "checkTable : Setting startValues to (already read) "
345  << "boundaryData"
346  /this->patch().name()
347  /sampleTimes_[startSampleTime_].name()
348  << endl;
349  }
350  startSampledValues_ = endSampledValues_;
351  startAverage_ = endAverage_;
352  }
353  else
354  {
355  if (debug)
356  {
357  Pout<< "checkTable : Reading startValues from "
358  << "boundaryData"
359  /this->patch().name()
360  /sampleTimes_[lo].name()
361  << endl;
362  }
363 
364  // Reread values and interpolate
365  fileName valsFile
366  (
367  this->db().time().constant()
368  /"boundaryData"
369  /this->patch().name()
370  /sampleTimes_[startSampleTime_].name()
371  /fieldTableName_
372  );
373 
374  Field<Type> vals;
375 
376  if (setAverage_)
377  {
378  AverageField<Type> avals((IFstream(valsFile)()));
379  vals = avals;
380  startAverage_ = avals.average();
381  }
382  else
383  {
384  (IFstream(valsFile)()) >> vals;
385  }
386 
387  if (vals.size() != mapperPtr_().sourceSize())
388  {
390  << "Number of values (" << vals.size()
391  << ") differs from the number of points ("
392  << mapperPtr_().sourceSize()
393  << ") in file " << valsFile << exit(FatalError);
394  }
395 
396  startSampledValues_ = mapperPtr_().interpolate(vals);
397  }
398  }
399 
400  if (hi != endSampleTime_)
401  {
402  endSampleTime_ = hi;
403 
404  if (endSampleTime_ == -1)
405  {
406  // endTime no longer valid. Might as well clear endValues.
407  if (debug)
408  {
409  Pout<< "checkTable : Clearing endValues" << endl;
410  }
411  endSampledValues_.clear();
412  }
413  else
414  {
415  if (debug)
416  {
417  Pout<< "checkTable : Reading endValues from "
418  << "boundaryData"
419  /this->patch().name()
420  /sampleTimes_[endSampleTime_].name()
421  << endl;
422  }
423 
424  // Reread values and interpolate
425  fileName valsFile
426  (
427  this->db().time().constant()
428  /"boundaryData"
429  /this->patch().name()
430  /sampleTimes_[endSampleTime_].name()
431  /fieldTableName_
432  );
433 
434  Field<Type> vals;
435 
436  if (setAverage_)
437  {
438  AverageField<Type> avals((IFstream(valsFile)()));
439  vals = avals;
440  endAverage_ = avals.average();
441  }
442  else
443  {
444  (IFstream(valsFile)()) >> vals;
445  }
446 
447  if (vals.size() != mapperPtr_().sourceSize())
448  {
450  << "Number of values (" << vals.size()
451  << ") differs from the number of points ("
452  << mapperPtr_().sourceSize()
453  << ") in file " << valsFile << exit(FatalError);
454  }
455 
456  endSampledValues_ = mapperPtr_().interpolate(vals);
457  }
458  }
459 }
460 
461 
462 template<class Type>
464 {
465  if (this->updated())
466  {
467  return;
468  }
469 
470  checkTable();
471 
472  // Interpolate between the sampled data
473 
474  Type wantedAverage;
475 
476  if (endSampleTime_ == -1)
477  {
478  // only start value
479  if (debug)
480  {
481  Pout<< "updateCoeffs : Sampled, non-interpolated values"
482  << " from start time:"
483  << sampleTimes_[startSampleTime_].name() << nl;
484  }
485 
486  this->operator==(startSampledValues_);
487  wantedAverage = startAverage_;
488  }
489  else
490  {
491  scalar start = sampleTimes_[startSampleTime_].value();
492  scalar end = sampleTimes_[endSampleTime_].value();
493 
494  scalar s = (this->db().time().value()-start)/(end-start);
495 
496  if (debug)
497  {
498  Pout<< "updateCoeffs : 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  this->operator==((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  const Field<Type>& fld = *this;
514 
515  Type averagePsi = gAverage(fld);
516 
517  if (debug)
518  {
519  Pout<< "updateCoeffs :"
520  << " actual average:" << averagePsi
521  << " wanted average:" << wantedAverage
522  << endl;
523  }
524 
525  if (mag(averagePsi) < vSmall)
526  {
527  // Field too small to scale. Offset instead.
528  const Type offset = wantedAverage - averagePsi;
529  if (debug)
530  {
531  Pout<< "updateCoeffs :"
532  << " offsetting with:" << offset << endl;
533  }
534  this->operator==(fld+offset);
535  }
536  else
537  {
538  const scalar scale = mag(wantedAverage)/mag(averagePsi);
539 
540  if (debug)
541  {
542  Pout<< "updateCoeffs :"
543  << " scaling with:" << scale << endl;
544  }
545  this->operator==(scale*fld);
546  }
547  }
548 
549  // Apply offset to mapped values
550  if (offset_.valid())
551  {
552  const scalar t = this->db().time().timeOutputValue();
553  this->operator==(*this + offset_->value(t));
554  }
555 
556  if (debug)
557  {
558  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
559  << " max:" << gMax(*this)
560  << " avg:" << gAverage(*this) << endl;
561  }
562 
564 }
565 
566 
567 template<class Type>
569 (
570  Ostream& os
571 ) const
572 {
574 
575  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
576 
577  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
578 
580  (
581  os,
582  "fieldTable",
583  this->internalField().name(),
584  fieldTableName_
585  );
586 
588  (
589  os,
590  "mapMethod",
591  word("planarInterpolation"),
592  mapMethod_
593  );
594 
595  if (offset_.valid())
596  {
597  writeEntry(os, offset_());
598  }
599 }
600 
601 
602 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
dictionary dict
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
const word & name() const
Return name.
Definition: IOobject.H:303
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
A class for handling file names.
Definition: fileName.H:79
void checkTable()
Find boundary data in between current time and interpolate.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:808
A time-varying form of a mapped fixed value boundary condition.
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A primitive field with a separate average value.
Definition: AverageField.H:49
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Type gMin(const FieldField< Field, Type > &f)
A FixedValue boundary condition for pointField.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
Foam::pointPatchFieldMapper.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
label size() const
Return size.
timeVaryingMappedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
Abstract base class for point-mesh patch fields.
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
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
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))
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)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
const polyMesh & mesh() const
Return the mesh reference.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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)
Input from file stream.
Definition: IFstream.H:81
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const Type & average() const
Definition: AverageField.C:61
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Interpolates between two sets of unstructured points using 2D Delaunay triangulation. Used in e.g. timeVaryingMapped bcs.
IOerror FatalIOError