timeVaryingMappedFixedValueFvPatchField.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 | Copyright (C) 2011-2017 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 fvPatch& 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 fvPatch& p,
63  const dictionary& dict
64 )
65 :
66  fixedValueFvPatchField<Type>(p, iF, dict, false),
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("fieldTable", fieldTableName_);
107 
108  if (dict.found("value"))
109  {
110  fvPatchField<Type>::operator==(Field<Type>("value", dict, p.size()));
111  }
112  else
113  {
114  // Note: we use evaluate() here to trigger updateCoeffs followed
115  // by re-setting of fvatchfield::updated_ flag. This is
116  // so if first use is in the next time step it retriggers
117  // a new update.
118  this->evaluate(Pstream::commsTypes::blocking);
119  }
120 }
121 
122 
123 template<class Type>
126 (
128  const fvPatch& p,
130  const fvPatchFieldMapper& mapper
131 )
132 :
133  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
134  fieldTableName_(ptf.fieldTableName_),
135  setAverage_(ptf.setAverage_),
136  perturb_(ptf.perturb_),
137  mapMethod_(ptf.mapMethod_),
138  mapperPtr_(nullptr),
139  sampleTimes_(0),
140  startSampleTime_(-1),
141  startSampledValues_(0),
142  startAverage_(Zero),
143  endSampleTime_(-1),
144  endSampledValues_(0),
145  endAverage_(Zero),
146  offset_(ptf.offset_, false)
147 {}
148 
149 
150 template<class Type>
153 (
155 )
156 :
158  fieldTableName_(ptf.fieldTableName_),
159  setAverage_(ptf.setAverage_),
160  perturb_(ptf.perturb_),
161  mapMethod_(ptf.mapMethod_),
162  mapperPtr_(nullptr),
163  sampleTimes_(ptf.sampleTimes_),
164  startSampleTime_(ptf.startSampleTime_),
165  startSampledValues_(ptf.startSampledValues_),
166  startAverage_(ptf.startAverage_),
167  endSampleTime_(ptf.endSampleTime_),
168  endSampledValues_(ptf.endSampledValues_),
169  endAverage_(ptf.endAverage_),
170  offset_(ptf.offset_, false)
171 {}
172 
173 
174 template<class Type>
177 (
180 )
181 :
183  fieldTableName_(ptf.fieldTableName_),
184  setAverage_(ptf.setAverage_),
185  perturb_(ptf.perturb_),
186  mapMethod_(ptf.mapMethod_),
187  mapperPtr_(nullptr),
188  sampleTimes_(ptf.sampleTimes_),
189  startSampleTime_(ptf.startSampleTime_),
190  startSampledValues_(ptf.startSampledValues_),
191  startAverage_(ptf.startAverage_),
192  endSampleTime_(ptf.endSampleTime_),
193  endSampledValues_(ptf.endSampledValues_),
194  endAverage_(ptf.endAverage_),
195  offset_(ptf.offset_, false)
196 {}
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class Type>
203 (
204  const fvPatchFieldMapper& m
205 )
206 {
208  if (startSampledValues_.size())
209  {
210  startSampledValues_.autoMap(m);
211  endSampledValues_.autoMap(m);
212  }
213  // Clear interpolator
214  mapperPtr_.clear();
215  startSampleTime_ = -1;
216  endSampleTime_ = -1;
217 }
218 
219 
220 template<class Type>
222 (
223  const fvPatchField<Type>& ptf,
224  const labelList& addr
225 )
226 {
228 
230  refCast<const timeVaryingMappedFixedValueFvPatchField<Type>>(ptf);
231 
232  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
233  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
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 (mapperPtr_.empty())
247  {
248  // Reread values and interpolate
249  fileName samplePointsFile
250  (
251  this->db().time().constant()
252  /"boundaryData"
253  /this->patch().name()
254  /"points"
255  );
256 
257  pointField samplePoints((IFstream(samplePointsFile)()));
258 
259  if (debug)
260  {
261  Info<< "timeVaryingMappedFixedValueFvPatchField :"
262  << " Read " << samplePoints.size() << " sample points from "
263  << samplePointsFile << endl;
264  }
265 
266 
267  // tbd: run-time selection
268  bool nearestOnly =
269  (
270  !mapMethod_.empty()
271  && mapMethod_ != "planarInterpolation"
272  );
273 
274  // Allocate the interpolator
275  mapperPtr_.reset
276  (
278  (
279  samplePoints,
280  this->patch().patch().faceCentres(),
281  perturb_,
282  nearestOnly
283  )
284  );
285 
286  // Read the times for which data is available
287  const fileName samplePointsDir = samplePointsFile.path();
288  sampleTimes_ = Time::findTimes(samplePointsDir);
289 
290  if (debug)
291  {
292  Info<< "timeVaryingMappedFixedValueFvPatchField : In directory "
293  << samplePointsDir << " found times "
294  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
295  << endl;
296  }
297  }
298 
299 
300  // Find current time in sampleTimes
301  label lo = -1;
302  label hi = -1;
303 
304  bool foundTime = mapperPtr_().findTime
305  (
306  sampleTimes_,
307  startSampleTime_,
308  this->db().time().value(),
309  lo,
310  hi
311  );
312 
313  if (!foundTime)
314  {
316  << "Cannot find starting sampling values for current time "
317  << this->db().time().value() << nl
318  << "Have sampling values for times "
319  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
320  << "In directory "
321  << this->db().time().constant()/"boundaryData"/this->patch().name()
322  << "\n on patch " << this->patch().name()
323  << " of field " << fieldTableName_
324  << exit(FatalError);
325  }
326 
327 
328  // Update sampled data fields.
329 
330  if (lo != startSampleTime_)
331  {
332  startSampleTime_ = lo;
333 
334  if (startSampleTime_ == endSampleTime_)
335  {
336  // No need to reread since are end values
337  if (debug)
338  {
339  Pout<< "checkTable : Setting startValues to (already read) "
340  << "boundaryData"
341  /this->patch().name()
342  /sampleTimes_[startSampleTime_].name()
343  << endl;
344  }
345  startSampledValues_ = endSampledValues_;
346  startAverage_ = endAverage_;
347  }
348  else
349  {
350  if (debug)
351  {
352  Pout<< "checkTable : Reading startValues from "
353  << "boundaryData"
354  /this->patch().name()
355  /sampleTimes_[lo].name()
356  << endl;
357  }
358 
359 
360  // Reread values and interpolate
361  fileName valsFile
362  (
363  this->db().time().constant()
364  /"boundaryData"
365  /this->patch().name()
366  /sampleTimes_[startSampleTime_].name()
367  /fieldTableName_
368  );
369 
370  Field<Type> vals;
371 
372  if (setAverage_)
373  {
374  AverageField<Type> avals((IFstream(valsFile)()));
375  vals = avals;
376  startAverage_ = avals.average();
377  }
378  else
379  {
380  IFstream(valsFile)() >> vals;
381  }
382 
383  if (vals.size() != mapperPtr_().sourceSize())
384  {
386  << "Number of values (" << vals.size()
387  << ") differs from the number of points ("
388  << mapperPtr_().sourceSize()
389  << ") in file " << valsFile << exit(FatalError);
390  }
391 
392  startSampledValues_ = mapperPtr_().interpolate(vals);
393  }
394  }
395 
396  if (hi != endSampleTime_)
397  {
398  endSampleTime_ = hi;
399 
400  if (endSampleTime_ == -1)
401  {
402  // endTime no longer valid. Might as well clear endValues.
403  if (debug)
404  {
405  Pout<< "checkTable : Clearing endValues" << endl;
406  }
407  endSampledValues_.clear();
408  }
409  else
410  {
411  if (debug)
412  {
413  Pout<< "checkTable : Reading endValues from "
414  << "boundaryData"
415  /this->patch().name()
416  /sampleTimes_[endSampleTime_].name()
417  << endl;
418  }
419 
420  // Reread values and interpolate
421  fileName valsFile
422  (
423  this->db().time().constant()
424  /"boundaryData"
425  /this->patch().name()
426  /sampleTimes_[endSampleTime_].name()
427  /fieldTableName_
428  );
429 
430  Field<Type> vals;
431 
432  if (setAverage_)
433  {
434  AverageField<Type> avals((IFstream(valsFile)()));
435  vals = avals;
436  endAverage_ = avals.average();
437  }
438  else
439  {
440  IFstream(valsFile)() >> vals;
441  }
442 
443  if (vals.size() != mapperPtr_().sourceSize())
444  {
446  << "Number of values (" << vals.size()
447  << ") differs from the number of points ("
448  << mapperPtr_().sourceSize()
449  << ") in file " << valsFile << exit(FatalError);
450  }
451 
452  endSampledValues_ = mapperPtr_().interpolate(vals);
453  }
454  }
455 }
456 
457 
458 template<class Type>
460 {
461  if (this->updated())
462  {
463  return;
464  }
465 
466 
467  checkTable();
468 
469  // Interpolate between the sampled data
470 
471  Type wantedAverage;
472 
473  if (endSampleTime_ == -1)
474  {
475  // Only start value
476  if (debug)
477  {
478  Pout<< "updateCoeffs : Sampled, non-interpolated values"
479  << " from start time:"
480  << sampleTimes_[startSampleTime_].name() << nl;
481  }
482 
483  this->operator==(startSampledValues_);
484  wantedAverage = startAverage_;
485  }
486  else
487  {
488  scalar start = sampleTimes_[startSampleTime_].value();
489  scalar end = sampleTimes_[endSampleTime_].value();
490 
491  scalar s = (this->db().time().value() - start)/(end - start);
492 
493  if (debug)
494  {
495  Pout<< "updateCoeffs : Sampled, interpolated values"
496  << " between start time:"
497  << sampleTimes_[startSampleTime_].name()
498  << " and end time:" << sampleTimes_[endSampleTime_].name()
499  << " with weight:" << s << endl;
500  }
501 
502  this->operator==((1 - s)*startSampledValues_ + s*endSampledValues_);
503  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
504  }
505 
506  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
507  // offsetting.
508  if (setAverage_)
509  {
510  const Field<Type>& fld = *this;
511 
512  Type averagePsi =
513  gSum(this->patch().magSf()*fld)
514  /gSum(this->patch().magSf());
515 
516  if (debug)
517  {
518  Pout<< "updateCoeffs :"
519  << " actual average:" << averagePsi
520  << " wanted average:" << wantedAverage
521  << endl;
522  }
523 
524  if (mag(averagePsi) < VSMALL)
525  {
526  // Field too small to scale. Offset instead.
527  const Type offset = wantedAverage - averagePsi;
528  if (debug)
529  {
530  Pout<< "updateCoeffs :"
531  << " offsetting with:" << offset << endl;
532  }
533  this->operator==(fld + offset);
534  }
535  else
536  {
537  const scalar scale = mag(wantedAverage)/mag(averagePsi);
538 
539  if (debug)
540  {
541  Pout<< "updateCoeffs :"
542  << " scaling with:" << scale << endl;
543  }
544  this->operator==(scale*fld);
545  }
546  }
547 
548  // Apply offset to mapped values
549  if (offset_.valid())
550  {
551  const scalar t = this->db().time().timeOutputValue();
552  this->operator==(*this + offset_->value(t));
553  }
554 
555  if (debug)
556  {
557  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
558  << " max:" << gMax(*this)
559  << " avg:" << gAverage(*this) << endl;
560  }
561 
563 }
564 
565 
566 template<class Type>
568 (
569  Ostream& os
570 ) const
571 {
573 
574  this->writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
575 
576  this->writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
577 
578  this->writeEntryIfDifferent
579  (
580  os,
581  "fieldTable",
582  this->internalField().name(),
583  fieldTableName_
584  );
585 
586  this->writeEntryIfDifferent
587  (
588  os,
589  "mapMethod",
590  word("planarInterpolation"),
591  mapMethod_
592  );
593 
594  if (offset_.valid())
595  {
596  offset_->writeData(os);
597  }
598 
599  this->writeEntry("value", os);
600 }
601 
602 
603 // ************************************************************************* //
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
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
const word & name() const
Return name.
Definition: IOobject.H:291
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A primitive field with a separate average value.
Definition: AverageField.H:49
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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.
Definition: Switch.H:60
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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))
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
This boundary conditions interpolates the values from a set of supplied points in space and time...
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))
Type gSum(const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
Definition: Field.H:57
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:180
static const zero Zero
Definition: zero.H:91
void checkTable()
Find boundary data inbetween current time and interpolate.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Type gMax(const FieldField< Field, Type > &f)
Input from file stream.
Definition: IFstream.H:81
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
timeVaryingMappedFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Constant dispersed-phase particle diameter model.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
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:250
volScalarField & p
const Type & average() const
Definition: AverageField.C:61
Interpolates between two sets of unstructured points using 2D Delaunay triangulation. Used in e.g. timeVaryingMapped bcs.
IOerror FatalIOError