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-2025 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_ = Function1<Type>::New
67  (
68  "offset",
69  this->db().time().userUnits(),
70  iF.dimensions(),
71  dict
72  );
73  }
74 
75  if
76  (
77  mapMethod_ != "planarInterpolation"
78  && mapMethod_ != "nearest"
79  )
80  {
82  (
83  dict
84  ) << "mapMethod should be one of 'planarInterpolation'"
85  << ", 'nearest'" << exit(FatalIOError);
86  }
87 
88  dict.readIfPresent("fieldTableName", fieldTableName_);
89 
90  if (dict.found("value"))
91  {
93  (
94  Field<Type>("value", iF.dimensions(), dict, p.size())
95  );
96  }
97  else
98  {
99  // Note: use evaluate to do updateCoeffs followed by a reset
100  // of the pointPatchField::updated_ flag. This is
101  // so if first use is in the next time step it retriggers
102  // a new update.
104  }
105 }
106 
107 
108 template<class Type>
111 (
113  const pointPatch& p,
115  const fieldMapper& mapper
116 )
117 :
118  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
119  fieldTableName_(ptf.fieldTableName_),
120  setAverage_(ptf.setAverage_),
121  perturb_(ptf.perturb_),
122  mapMethod_(ptf.mapMethod_),
123  mapperPtr_(nullptr),
124  sampleTimes_(0),
125  startSampleTime_(-1),
126  startSampledValues_(0),
127  startAverage_(Zero),
128  endSampleTime_(-1),
129  endSampledValues_(0),
130  endAverage_(Zero),
131  offset_(ptf.offset_, false)
132 {}
133 
134 
135 template<class Type>
138 (
141 )
142 :
143  fixedValuePointPatchField<Type>(ptf, iF),
144  fieldTableName_(ptf.fieldTableName_),
145  setAverage_(ptf.setAverage_),
146  perturb_(ptf.perturb_),
147  mapMethod_(ptf.mapMethod_),
148  mapperPtr_(ptf.mapperPtr_),
149  sampleTimes_(ptf.sampleTimes_),
150  startSampleTime_(ptf.startSampleTime_),
151  startSampledValues_(ptf.startSampledValues_),
152  startAverage_(ptf.startAverage_),
153  endSampleTime_(ptf.endSampleTime_),
154  endSampledValues_(ptf.endSampledValues_),
155  endAverage_(ptf.endAverage_),
156  offset_(ptf.offset_, false)
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
162 template<class Type>
164 (
165  const pointPatchField<Type>& ptf,
166  const fieldMapper& mapper
167 )
168 {
169  reset(ptf);
170 }
171 
172 
173 template<class Type>
175 (
176  const pointPatchField<Type>& ptf
177 )
178 {
180 
182  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
183 
184  startSampledValues_.reset(tiptf.startSampledValues_);
185  endSampledValues_.reset(tiptf.endSampledValues_);
186 
187  // Clear interpolator
188  mapperPtr_.clear();
189  startSampleTime_ = -1;
190  endSampleTime_ = -1;
191 }
192 
193 
194 template<class Type>
196 {
197  // Initialise
198  if (startSampleTime_ == -1 && endSampleTime_ == -1)
199  {
200  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
201 
202  // Read the initial point position
203  pointField meshPts;
204 
205  if (pMesh.pointsInstance() == pMesh.facesInstance())
206  {
207  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
208  }
209  else
210  {
211  // Load points from facesInstance
212  if (debug)
213  {
214  Info<< "Reloading points0 from " << pMesh.facesInstance()
215  << endl;
216  }
217 
218  pointIOField points0
219  (
220  IOobject
221  (
222  "points",
223  pMesh.facesInstance(),
225  pMesh,
228  false
229  )
230  );
231  meshPts = pointField(points0, this->patch().meshPoints());
232  }
233 
234  // Reread values and interpolate
235  fileName samplePointsFile
236  (
237  this->db().time().constant()
238  /"boundaryData"
239  /this->patch().name()
240  /"points"
241  );
242 
243  pointField samplePoints((IFstream(samplePointsFile)()));
244 
245  // tbd: run-time selection
246  bool nearestOnly =
247  (
248  !mapMethod_.empty()
249  && mapMethod_ != "planarInterpolation"
250  );
251 
252  // Allocate the interpolator
253  mapperPtr_.reset
254  (
256  (
257  samplePoints,
258  meshPts,
259  perturb_,
260  nearestOnly
261  )
262  );
263 
264  // Read the times for which data is available
265 
266  const fileName samplePointsDir = samplePointsFile.path();
267  sampleTimes_ = this->db().time().findTimes(samplePointsDir);
268 
269  if (debug)
270  {
271  Info<< "timeVaryingMappedFixedValuePointPatchField : In directory "
272  << samplePointsDir << " found times "
274  << endl;
275  }
276  }
277 
278  // Find current time in sampleTimes
279  label lo = -1;
280  label hi = -1;
281 
282  bool foundTime = mapperPtr_().findTime
283  (
284  sampleTimes_,
285  startSampleTime_,
286  this->db().time().userTimeValue(),
287  lo,
288  hi
289  );
290 
291  if (!foundTime)
292  {
294  << "Cannot find starting sampling values for current time "
295  << this->db().time().userTimeValue() << nl
296  << "Have sampling values for times "
298  << "In directory "
299  << this->db().time().constant()/"boundaryData"/this->patch().name()
300  << "\n on patch " << this->patch().name()
301  << " of field " << fieldTableName_
302  << exit(FatalError);
303  }
304 
305 
306  // Update sampled data fields.
307 
308  if (lo != startSampleTime_)
309  {
310  startSampleTime_ = lo;
311 
312  if (startSampleTime_ == endSampleTime_)
313  {
314  // No need to reread since are end values
315  if (debug)
316  {
317  Pout<< "checkTable : Setting startValues to (already read) "
318  << "boundaryData"
319  /this->patch().name()
320  /sampleTimes_[startSampleTime_].name()
321  << endl;
322  }
323  startSampledValues_ = endSampledValues_;
324  startAverage_ = endAverage_;
325  }
326  else
327  {
328  if (debug)
329  {
330  Pout<< "checkTable : Reading startValues from "
331  << "boundaryData"
332  /this->patch().name()
333  /sampleTimes_[lo].name()
334  << endl;
335  }
336 
337  // Reread values and interpolate
338  fileName valsFile
339  (
340  this->db().time().constant()
341  /"boundaryData"
342  /this->patch().name()
343  /sampleTimes_[startSampleTime_].name()
344  /fieldTableName_
345  );
346 
347  Field<Type> vals;
348 
349  if (setAverage_)
350  {
351  AverageField<Type> avals((IFstream(valsFile)()));
352  vals = avals;
353  startAverage_ = avals.average();
354  }
355  else
356  {
357  (IFstream(valsFile)()) >> vals;
358  }
359 
360  if (vals.size() != mapperPtr_().sourceSize())
361  {
363  << "Number of values (" << vals.size()
364  << ") differs from the number of points ("
365  << mapperPtr_().sourceSize()
366  << ") in file " << valsFile << exit(FatalError);
367  }
368 
369  startSampledValues_ = mapperPtr_().interpolate(vals);
370  }
371  }
372 
373  if (hi != endSampleTime_)
374  {
375  endSampleTime_ = hi;
376 
377  if (endSampleTime_ == -1)
378  {
379  // endTime no longer valid. Might as well clear endValues.
380  if (debug)
381  {
382  Pout<< "checkTable : Clearing endValues" << endl;
383  }
384  endSampledValues_.clear();
385  }
386  else
387  {
388  if (debug)
389  {
390  Pout<< "checkTable : Reading endValues from "
391  << "boundaryData"
392  /this->patch().name()
393  /sampleTimes_[endSampleTime_].name()
394  << endl;
395  }
396 
397  // Reread values and interpolate
398  fileName valsFile
399  (
400  this->db().time().constant()
401  /"boundaryData"
402  /this->patch().name()
403  /sampleTimes_[endSampleTime_].name()
404  /fieldTableName_
405  );
406 
407  Field<Type> vals;
408 
409  if (setAverage_)
410  {
411  AverageField<Type> avals((IFstream(valsFile)()));
412  vals = avals;
413  endAverage_ = avals.average();
414  }
415  else
416  {
417  (IFstream(valsFile)()) >> vals;
418  }
419 
420  if (vals.size() != mapperPtr_().sourceSize())
421  {
423  << "Number of values (" << vals.size()
424  << ") differs from the number of points ("
425  << mapperPtr_().sourceSize()
426  << ") in file " << valsFile << exit(FatalError);
427  }
428 
429  endSampledValues_ = mapperPtr_().interpolate(vals);
430  }
431  }
432 }
433 
434 
435 template<class Type>
437 {
438  if (this->updated())
439  {
440  return;
441  }
442 
443  checkTable();
444 
445  // Interpolate between the sampled data
446 
447  Type wantedAverage;
448 
449  if (endSampleTime_ == -1)
450  {
451  // only start value
452  if (debug)
453  {
454  Pout<< "updateCoeffs : Sampled, non-interpolated values"
455  << " from start time:"
456  << sampleTimes_[startSampleTime_].name() << nl;
457  }
458 
459  this->operator==(startSampledValues_);
460  wantedAverage = startAverage_;
461  }
462  else
463  {
464  scalar start = sampleTimes_[startSampleTime_].value();
465  scalar end = sampleTimes_[endSampleTime_].value();
466 
467  scalar s = (this->db().time().userTimeValue()-start)/(end-start);
468 
469  if (debug)
470  {
471  Pout<< "updateCoeffs : Sampled, interpolated values"
472  << " between start time:"
473  << sampleTimes_[startSampleTime_].name()
474  << " and end time:" << sampleTimes_[endSampleTime_].name()
475  << " with weight:" << s << endl;
476  }
477 
478  this->operator==((1 - s)*startSampledValues_ + s*endSampledValues_);
479  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
480  }
481 
482  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
483  // offsetting.
484  if (setAverage_)
485  {
486  const Field<Type>& fld = *this;
487 
488  Type averagePsi = gAverage(fld);
489 
490  if (debug)
491  {
492  Pout<< "updateCoeffs :"
493  << " actual average:" << averagePsi
494  << " wanted average:" << wantedAverage
495  << endl;
496  }
497 
498  if (mag(averagePsi) < vSmall)
499  {
500  // Field too small to scale. Offset instead.
501  const Type offset = wantedAverage - averagePsi;
502  if (debug)
503  {
504  Pout<< "updateCoeffs :"
505  << " offsetting with:" << offset << endl;
506  }
507  this->operator==(fld + offset);
508  }
509  else
510  {
511  const scalar scale = mag(wantedAverage)/mag(averagePsi);
512 
513  if (debug)
514  {
515  Pout<< "updateCoeffs :"
516  << " scaling with:" << scale << endl;
517  }
518  this->operator==(scale*fld);
519  }
520  }
521 
522  // Apply offset to mapped values
523  if (offset_.valid())
524  {
525  this->operator==
526  (
527  *this + offset_->value(this->db().time().value())
528  );
529  }
530 
531  if (debug)
532  {
533  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
534  << " max:" << gMax(*this)
535  << " avg:" << gAverage(*this) << endl;
536  }
537 
539 }
540 
541 
542 template<class Type>
544 (
545  Ostream& os
546 ) const
547 {
549 
550  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
551 
552  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
553 
555  (
556  os,
557  "fieldTable",
558  this->internalField().name(),
559  fieldTableName_
560  );
561 
563  (
564  os,
565  "mapMethod",
566  word("planarInterpolation"),
567  mapMethod_
568  );
569 
570  if (offset_.valid())
571  {
572  writeEntry(os, offset_());
573  }
574 }
575 
576 
577 // ************************************************************************* //
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
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:284
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:994
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:273
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.
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(lagrangian::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(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::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:258
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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:267
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p