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