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-2022 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  const pointPatchField<Type>& ptf
225 )
226 {
228 
230  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
231 
232  startSampledValues_.reset(tiptf.startSampledValues_);
233  endSampledValues_.reset(tiptf.endSampledValues_);
234 
235  // Clear interpolator
236  mapperPtr_.clear();
237  startSampleTime_ = -1;
238  endSampleTime_ = -1;
239 }
240 
241 
242 template<class Type>
244 {
245  // Initialise
246  if (startSampleTime_ == -1 && endSampleTime_ == -1)
247  {
248  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
249 
250  // Read the initial point position
251  pointField meshPts;
252 
253  if (pMesh.pointsInstance() == pMesh.facesInstance())
254  {
255  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
256  }
257  else
258  {
259  // Load points from facesInstance
260  if (debug)
261  {
262  Info<< "Reloading points0 from " << pMesh.facesInstance()
263  << endl;
264  }
265 
266  pointIOField points0
267  (
268  IOobject
269  (
270  "points",
271  pMesh.facesInstance(),
272  polyMesh::meshSubDir,
273  pMesh,
274  IOobject::MUST_READ,
275  IOobject::NO_WRITE,
276  false
277  )
278  );
279  meshPts = pointField(points0, this->patch().meshPoints());
280  }
281 
282  // Reread values and interpolate
283  fileName samplePointsFile
284  (
285  this->db().time().constant()
286  /"boundaryData"
287  /this->patch().name()
288  /"points"
289  );
290 
291  pointField samplePoints((IFstream(samplePointsFile)()));
292 
293  // tbd: run-time selection
294  bool nearestOnly =
295  (
296  !mapMethod_.empty()
297  && mapMethod_ != "planarInterpolation"
298  );
299 
300  // Allocate the interpolator
301  mapperPtr_.reset
302  (
304  (
305  samplePoints,
306  meshPts,
307  perturb_,
308  nearestOnly
309  )
310  );
311 
312  // Read the times for which data is available
313 
314  const fileName samplePointsDir = samplePointsFile.path();
315  sampleTimes_ = Time::findTimes(samplePointsDir);
316 
317  if (debug)
318  {
319  Info<< "timeVaryingMappedFixedValuePointPatchField : In directory "
320  << samplePointsDir << " found times "
321  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
322  << endl;
323  }
324  }
325 
326  // Find current time in sampleTimes
327  label lo = -1;
328  label hi = -1;
329 
330  bool foundTime = mapperPtr_().findTime
331  (
332  sampleTimes_,
333  startSampleTime_,
334  this->db().time().value(),
335  lo,
336  hi
337  );
338 
339  if (!foundTime)
340  {
342  << "Cannot find starting sampling values for current time "
343  << this->db().time().value() << nl
344  << "Have sampling values for times "
345  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
346  << "In directory "
347  << this->db().time().constant()/"boundaryData"/this->patch().name()
348  << "\n on patch " << this->patch().name()
349  << " of field " << fieldTableName_
350  << exit(FatalError);
351  }
352 
353 
354  // Update sampled data fields.
355 
356  if (lo != startSampleTime_)
357  {
358  startSampleTime_ = lo;
359 
360  if (startSampleTime_ == endSampleTime_)
361  {
362  // No need to reread since are end values
363  if (debug)
364  {
365  Pout<< "checkTable : Setting startValues to (already read) "
366  << "boundaryData"
367  /this->patch().name()
368  /sampleTimes_[startSampleTime_].name()
369  << endl;
370  }
371  startSampledValues_ = endSampledValues_;
372  startAverage_ = endAverage_;
373  }
374  else
375  {
376  if (debug)
377  {
378  Pout<< "checkTable : Reading startValues from "
379  << "boundaryData"
380  /this->patch().name()
381  /sampleTimes_[lo].name()
382  << endl;
383  }
384 
385  // Reread values and interpolate
386  fileName valsFile
387  (
388  this->db().time().constant()
389  /"boundaryData"
390  /this->patch().name()
391  /sampleTimes_[startSampleTime_].name()
392  /fieldTableName_
393  );
394 
395  Field<Type> vals;
396 
397  if (setAverage_)
398  {
399  AverageField<Type> avals((IFstream(valsFile)()));
400  vals = avals;
401  startAverage_ = avals.average();
402  }
403  else
404  {
405  (IFstream(valsFile)()) >> vals;
406  }
407 
408  if (vals.size() != mapperPtr_().sourceSize())
409  {
411  << "Number of values (" << vals.size()
412  << ") differs from the number of points ("
413  << mapperPtr_().sourceSize()
414  << ") in file " << valsFile << exit(FatalError);
415  }
416 
417  startSampledValues_ = mapperPtr_().interpolate(vals);
418  }
419  }
420 
421  if (hi != endSampleTime_)
422  {
423  endSampleTime_ = hi;
424 
425  if (endSampleTime_ == -1)
426  {
427  // endTime no longer valid. Might as well clear endValues.
428  if (debug)
429  {
430  Pout<< "checkTable : Clearing endValues" << endl;
431  }
432  endSampledValues_.clear();
433  }
434  else
435  {
436  if (debug)
437  {
438  Pout<< "checkTable : Reading endValues from "
439  << "boundaryData"
440  /this->patch().name()
441  /sampleTimes_[endSampleTime_].name()
442  << endl;
443  }
444 
445  // Reread values and interpolate
446  fileName valsFile
447  (
448  this->db().time().constant()
449  /"boundaryData"
450  /this->patch().name()
451  /sampleTimes_[endSampleTime_].name()
452  /fieldTableName_
453  );
454 
455  Field<Type> vals;
456 
457  if (setAverage_)
458  {
459  AverageField<Type> avals((IFstream(valsFile)()));
460  vals = avals;
461  endAverage_ = avals.average();
462  }
463  else
464  {
465  (IFstream(valsFile)()) >> vals;
466  }
467 
468  if (vals.size() != mapperPtr_().sourceSize())
469  {
471  << "Number of values (" << vals.size()
472  << ") differs from the number of points ("
473  << mapperPtr_().sourceSize()
474  << ") in file " << valsFile << exit(FatalError);
475  }
476 
477  endSampledValues_ = mapperPtr_().interpolate(vals);
478  }
479  }
480 }
481 
482 
483 template<class Type>
485 {
486  if (this->updated())
487  {
488  return;
489  }
490 
491  checkTable();
492 
493  // Interpolate between the sampled data
494 
495  Type wantedAverage;
496 
497  if (endSampleTime_ == -1)
498  {
499  // only start value
500  if (debug)
501  {
502  Pout<< "updateCoeffs : Sampled, non-interpolated values"
503  << " from start time:"
504  << sampleTimes_[startSampleTime_].name() << nl;
505  }
506 
507  this->operator==(startSampledValues_);
508  wantedAverage = startAverage_;
509  }
510  else
511  {
512  scalar start = sampleTimes_[startSampleTime_].value();
513  scalar end = sampleTimes_[endSampleTime_].value();
514 
515  scalar s = (this->db().time().value()-start)/(end-start);
516 
517  if (debug)
518  {
519  Pout<< "updateCoeffs : Sampled, interpolated values"
520  << " between start time:"
521  << sampleTimes_[startSampleTime_].name()
522  << " and end time:" << sampleTimes_[endSampleTime_].name()
523  << " with weight:" << s << endl;
524  }
525 
526  this->operator==((1-s)*startSampledValues_ + s*endSampledValues_);
527  wantedAverage = (1-s)*startAverage_ + s*endAverage_;
528  }
529 
530  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
531  // offsetting.
532  if (setAverage_)
533  {
534  const Field<Type>& fld = *this;
535 
536  Type averagePsi = gAverage(fld);
537 
538  if (debug)
539  {
540  Pout<< "updateCoeffs :"
541  << " actual average:" << averagePsi
542  << " wanted average:" << wantedAverage
543  << endl;
544  }
545 
546  if (mag(averagePsi) < vSmall)
547  {
548  // Field too small to scale. Offset instead.
549  const Type offset = wantedAverage - averagePsi;
550  if (debug)
551  {
552  Pout<< "updateCoeffs :"
553  << " offsetting with:" << offset << endl;
554  }
555  this->operator==(fld+offset);
556  }
557  else
558  {
559  const scalar scale = mag(wantedAverage)/mag(averagePsi);
560 
561  if (debug)
562  {
563  Pout<< "updateCoeffs :"
564  << " scaling with:" << scale << endl;
565  }
566  this->operator==(scale*fld);
567  }
568  }
569 
570  // Apply offset to mapped values
571  if (offset_.valid())
572  {
573  const scalar t = this->db().time().userTimeValue();
574  this->operator==(*this + offset_->value(t));
575  }
576 
577  if (debug)
578  {
579  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
580  << " max:" << gMax(*this)
581  << " avg:" << gAverage(*this) << endl;
582  }
583 
585 }
586 
587 
588 template<class Type>
590 (
591  Ostream& os
592 ) const
593 {
595 
596  writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
597 
598  writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
599 
601  (
602  os,
603  "fieldTable",
604  this->internalField().name(),
605  fieldTableName_
606  );
607 
609  (
610  os,
611  "mapMethod",
612  word("planarInterpolation"),
613  mapMethod_
614  );
615 
616  if (offset_.valid())
617  {
618  writeEntry(os, offset_());
619  }
620 }
621 
622 
623 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
dictionary dict
const word & name() const
Return name.
Definition: IOobject.H:315
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
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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:888
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:306
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:1211
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:882
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
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:318
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:76
volScalarField & p
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
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.
virtual void reset(const pointPatchField< Type > &)
Reset the pointPatchField to the given pointPatchField.
IOerror FatalIOError