timeVaryingMappedFixedValuePointPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "Time.H"
31 #include "rawIOField.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
38 (
39  const pointPatch& p,
41 )
42 :
43  fixedValuePointPatchField<Type>(p, iF),
44  setAverage_(false),
45  perturb_(0),
46  fieldTableName_(iF.name()),
47  pointsName_("points"),
48  mapMethod_(),
49  mapperPtr_(nullptr),
50  sampleTimes_(),
51  sampleIndex_(-1, -1),
52  sampleAverage_(Zero, Zero),
53  sampleValues_(),
54  offset_(nullptr)
55 {}
56 
57 
58 template<class Type>
61 (
62  const pointPatch& p,
64  const dictionary& dict
65 )
66 :
67  fixedValuePointPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
68  setAverage_(dict.getOrDefault("setAverage", false)),
69  perturb_(dict.getOrDefault("perturb", 1e-5)),
70  fieldTableName_(iF.name()),
71  pointsName_(dict.getOrDefault<word>("points", "points")),
72  mapMethod_(),
73  mapperPtr_(nullptr),
74  sampleTimes_(),
75  sampleIndex_(-1, -1),
76  sampleAverage_(Zero, Zero),
77  sampleValues_(),
78  offset_
79  (
80  Function1<Type>::NewIfPresent("offset", dict, word::null, &this->db())
81  )
82 {
83  if
84  (
85  dict.readIfPresent("mapMethod", mapMethod_)
86  && !mapMethod_.empty()
87  && mapMethod_ != "nearest"
88  && !mapMethod_.starts_with("planar")
89  )
90  {
92  << "Unknown mapMethod type " << mapMethod_
93  << "\n\nValid mapMethod types :\n"
94  << "(nearest planar)" << nl
95  << exit(FatalIOError);
96  }
97 
98  dict.readIfPresentCompat
99  (
100  "fieldTable", {{"fieldTableName", 2206}},
101  fieldTableName_
102  );
103 
104  if (!this->readValueEntry(dict))
105  {
106  // Note: use evaluate to do updateCoeffs followed by a reset
107  // of the pointPatchField::updated_ flag. This is
108  // so if first use is in the next time step it retriggers
109  // a new update.
111  }
112 }
113 
114 
115 template<class Type>
118 (
119  const timeVaryingMappedFixedValuePointPatchField<Type>& ptf,
120  const pointPatch& p,
121  const DimensionedField<Type, pointMesh>& iF,
122  const pointPatchFieldMapper& mapper
123 )
124 :
125  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
126  setAverage_(ptf.setAverage_),
127  perturb_(ptf.perturb_),
128  fieldTableName_(ptf.fieldTableName_),
129  pointsName_(ptf.pointsName_),
130  mapMethod_(ptf.mapMethod_),
131  mapperPtr_(nullptr),
132  sampleTimes_(),
133  sampleIndex_(-1, -1),
134  sampleAverage_(Zero, Zero),
135  sampleValues_(),
136  offset_(ptf.offset_.clone())
137 {}
138 
139 
140 template<class Type>
143 (
145 )
146 :
147  fixedValuePointPatchField<Type>(ptf),
148  setAverage_(ptf.setAverage_),
149  perturb_(ptf.perturb_),
150  fieldTableName_(ptf.fieldTableName_),
151  pointsName_(ptf.pointsName_),
152  mapMethod_(ptf.mapMethod_),
153  mapperPtr_(ptf.mapperPtr_),
154  sampleTimes_(ptf.sampleTimes_),
155  sampleIndex_(ptf.sampleIndex_),
156  sampleAverage_(ptf.sampleAverage_),
157  sampleValues_(ptf.sampleValues_),
158  offset_(ptf.offset_.clone())
159 {}
160 
161 
162 template<class Type>
165 (
168 )
169 :
170  fixedValuePointPatchField<Type>(ptf, iF),
171  setAverage_(ptf.setAverage_),
172  perturb_(ptf.perturb_),
173  fieldTableName_(ptf.fieldTableName_),
174  pointsName_(ptf.pointsName_),
175  mapMethod_(ptf.mapMethod_),
176  mapperPtr_(ptf.mapperPtr_),
177  sampleTimes_(ptf.sampleTimes_),
178  sampleIndex_(ptf.sampleIndex_),
179  sampleAverage_(ptf.sampleAverage_),
180  sampleValues_(ptf.sampleValues_),
181  offset_(ptf.offset_.clone())
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
187 template<class Type>
189 (
190  const pointPatchFieldMapper& m
191 )
192 {
194 
195  if (sampleValues_.first().size())
196  {
197  sampleValues_.first().autoMap(m);
198  }
199 
200  if (sampleValues_.second().size())
201  {
202  sampleValues_.second().autoMap(m);
203  }
204 
205  // Clear interpolator
206  mapperPtr_.reset(nullptr);
207  sampleIndex_ = labelPair(-1, -1);
208 }
209 
210 
211 template<class Type>
213 (
214  const pointPatchField<Type>& ptf,
215  const labelList& addr
216 )
217 {
219 
221  refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
222 
223  sampleValues_.first().rmap(tiptf.sampleValues_.first(), addr);
224  sampleValues_.second().rmap(tiptf.sampleValues_.second(), addr);
225 
226  // Clear interpolator
227  mapperPtr_.reset(nullptr);
228  sampleIndex_ = labelPair(-1, -1);
229 }
230 
231 
232 template<class Type>
234 (
235  const label sampleIndex,
236  Field<Type>& field,
237  Type& avg
238 ) const
239 {
240  tmp<Field<Type>> tvalues;
241 
242  // Update sampled data fields
243  {
244  const Time& time = this->db().time();
245 
246  if (debug)
247  {
248  Pout<< "checkTable : Reading values from "
249  <<
250  (
251  "boundaryData"
252  / this->patch().name()
253  / sampleTimes_[sampleIndex].name()
254  / fieldTableName_
255  ) << endl;
256  }
257 
258  // Reread values and interpolate
259  const fileName valsFile
260  (
261  time.caseConstant()
262  /"boundaryData"
263  /this->patch().name()
264  /sampleTimes_[sampleIndex].name()
265  /fieldTableName_
266  );
267 
268  IOobject io
269  (
270  valsFile, // absolute path
271  time,
275  true // is global object (currently not used)
276  );
277 
278  rawIOField<Type> vals(io, setAverage_);
279 
280  if (vals.hasAverage())
281  {
282  avg = vals.average();
283  }
284 
285  if (vals.size() != mapperPtr_().sourceSize())
286  {
288  << "Number of values (" << vals.size()
289  << ") differs from the number of points ("
290  << mapperPtr_().sourceSize()
291  << ") in file " << valsFile << exit(FatalError);
292  }
293 
294  tvalues = tmp<Field<Type>>::New(std::move(vals.field()));
295  }
296 
297  // From input values to interpolated (sampled) positions
298  field = mapperPtr_().interpolate(tvalues);
299 }
300 
301 
302 template<class Type>
304 (
305  const scalar t
306 )
307 {
308  const Time& time = this->db().time();
309 
310  // Initialise
311  if (sampleIndex_.first() == -1 && sampleIndex_.second() == -1)
312  {
313  const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
314 
315  // Read the initial point position
316  pointField meshPts;
317 
318  if (pMesh.pointsInstance() == pMesh.facesInstance())
319  {
320  meshPts = pointField(pMesh.points(), this->patch().meshPoints());
321  }
322  else
323  {
324  // Load points from facesInstance
326  << "Reloading points0 from " << pMesh.facesInstance()
327  << endl;
328 
330  (
331  IOobject
332  (
333  "points",
334  pMesh.facesInstance(),
336  pMesh,
340  )
341  );
342  meshPts = pointField(points0, this->patch().meshPoints());
343  }
344 
345  // Reread values and interpolate
346  const fileName samplePointsFile
347  (
348  time.caseConstant()
349  /"boundaryData"
350  /this->patch().name()
351  /pointsName_
352  );
353 
354  IOobject io
355  (
356  samplePointsFile, // absolute path
357  time,
361  true // is global object (currently not used)
362  );
363 
364  // Read data (no average value!)
365  const rawIOField<point> samplePoints(io);
366 
367  // tbd: run-time selection
368  const bool nearestOnly =
369  (
370  !mapMethod_.empty() && !mapMethod_.starts_with("planar")
371  );
372 
373  // Allocate the interpolator
374  mapperPtr_.reset
375  (
376  new pointToPointPlanarInterpolation
377  (
378  samplePoints,
379  meshPts,
380  perturb_,
381  nearestOnly
382  )
383  );
384 
385 
386  // Read the times for which data is available
387  const fileName samplePointsDir = samplePointsFile.path();
388  sampleTimes_ = Time::findTimes(samplePointsDir);
389 
391  << "Found times "
393  << endl;
394  }
395 
396  // Find range of current time indices in sampleTimes
397  Pair<label> timeIndices = instant::findRange
398  (
399  sampleTimes_,
400  t, // time.value(),
401  sampleIndex_.first()
402  );
403 
404  if (timeIndices.first() < 0)
405  {
407  << "Cannot find starting sampling values for current time "
408  << time.value() << nl
409  << "Have sampling values for times "
411  << "In directory "
412  << time.constant()/"boundaryData"/this->patch().name()
413  << "\n on patch " << this->patch().name()
414  << " of field " << fieldTableName_
415  << exit(FatalError);
416  }
417 
418 
419  // Update sampled data fields.
420 
421  if (sampleIndex_.first() != timeIndices.first())
422  {
423  sampleIndex_.first() = timeIndices.first();
424 
425  if (sampleIndex_.first() == sampleIndex_.second())
426  {
427  // Can reuse previous end values
428  sampleValues_.first() = sampleValues_.second();
429  sampleAverage_.first() = sampleAverage_.second();
430  }
431  else
432  {
433  // Update first() values
434  this->updateSampledValues
435  (
436  sampleIndex_.first(),
437  sampleValues_.first(),
438  sampleAverage_.first()
439  );
440  }
441  }
442 
443  if (sampleIndex_.second() != timeIndices.second())
444  {
445  sampleIndex_.second() = timeIndices.second();
446 
447  if (sampleIndex_.second() == -1)
448  {
449  // Index no longer valid - clear values accordingly
450  sampleValues_.second().clear();
451  }
452  else
453  {
454  // Update second() values
455  this->updateSampledValues
456  (
457  sampleIndex_.second(),
458  sampleValues_.second(),
459  sampleAverage_.second()
460  );
461  }
462  }
463 }
464 
465 
466 template<class Type>
468 {
469  if (this->updated())
470  {
471  return;
472  }
473 
474  // Current time value
475  const scalar x = this->db().time().value();
476 
477  checkTable(x);
478 
479  // Interpolate between the sampled data
480  auto& fld = static_cast<Field<Type>&>(*this);
481  Type wantedAverage;
482 
483  if (sampleIndex_.second() == -1)
484  {
485  // Only start value
486  fld = sampleValues_.first();
487  wantedAverage = sampleAverage_.first();
488  }
489  else
490  {
491  const scalar beg = sampleTimes_[sampleIndex_.first()].value();
492  const scalar end = sampleTimes_[sampleIndex_.second()].value();
493  const scalar s = (x - beg)/(end - beg);
494 
495  fld = (1 - s)*sampleValues_.first() + s*sampleValues_.second();
496 
497  wantedAverage
498  = (1 - s)*sampleAverage_.first() + s*sampleAverage_.second();
499 
501  << "Sampled, interpolated values"
502  << " between time:"
503  << sampleTimes_[sampleIndex_.first()].name()
504  << " and time:" << sampleTimes_[sampleIndex_.second()].name()
505  << " with weight:" << s << endl;
506  }
507 
508  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
509  // offsetting.
510  if (setAverage_)
511  {
512  Type averagePsi = gAverage(fld);
513 
514  if (debug)
515  {
516  Pout<< "updateCoeffs :"
517  << " actual average:" << averagePsi
518  << " wanted average:" << wantedAverage
519  << endl;
520  }
521 
522  if (mag(averagePsi) < VSMALL)
523  {
524  // Field too small to scale. Offset instead.
525  const Type offset = wantedAverage - averagePsi;
526  if (debug)
527  {
528  Pout<< "updateCoeffs :"
529  << " offsetting with:" << offset << endl;
530  }
531  fld += offset;
532  }
533  else
534  {
535  const scalar scale = mag(wantedAverage)/mag(averagePsi);
536 
537  if (debug)
538  {
539  Pout<< "updateCoeffs :"
540  << " scaling with:" << scale << endl;
541  }
542  fld *= scale;
543  }
544  }
545 
546  // Apply offset to mapped values
547  if (offset_)
548  {
549  const scalar t = this->db().time().timeOutputValue();
550  fld += offset_->value(t);
551  }
552 
553  if (debug)
554  {
555  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
556  << " max:" << gMax(*this)
557  << " avg:" << gAverage(*this) << endl;
558  }
561 }
562 
563 
564 template<class Type>
566 (
567  Ostream& os
568 ) const
569 {
571 
573  (
574  "fieldTable",
575  this->internalField().name(),
576  fieldTableName_
577  );
578 
579  if (!pointsName_.empty())
580  {
581  os.writeEntryIfDifferent<word>("points", "points", pointsName_);
582  }
583 
584  if (!mapMethod_.empty() && !mapMethod_.starts_with("planar"))
585  {
586  os.writeEntry("mapMethod", mapMethod_);
587  }
588 
589  if (setAverage_)
590  {
591  os.writeEntry("setAverage", setAverage_);
592  }
593 
594  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
595 
596  if (offset_)
597  {
598  offset_->writeData(os);
599  }
600 }
601 
602 
603 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
static wordList timeNames(const instantList &times)
Helper: extract words of times.
dictionary dict
"blocking" : (MPI_Bsend, MPI_Recv)
rDeltaTY field()
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A time-varying form of a mapped fixed value boundary condition.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
static Pair< label > findRange(const UList< instant > &times, const scalar timeVal, const label start=-1)
Find lower/upper indices for given time value in list of instances (linear search) continuing after t...
Definition: instant.C:54
A FixedValue boundary condition for pointField.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Foam::pointPatchFieldMapper.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
Ignore writing from objectRegistry::writeObject()
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
timeVaryingMappedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
Abstract base class for point-mesh patch fields.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Type gMax(const FieldField< Field, Type > &f)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
OBJstream os(runTime.globalPath()/outputName)
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:109
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
bool starts_with(char c) const
True if string starts with given character (cf. C++20)
Definition: string.H:435
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
Type gAverage(const FieldField< Field, Type > &f)
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127