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