timeVaryingMappedFixedValueFvPatchField.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) 2011-2018 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const word& timeName
38 ) const
39 {
40  const fileName fieldFileName
41  (
42  dataDir_/timeName/sampleName_/fieldTableName_
43  );
44 
45  const fileName typeFieldFileName
46  (
47  dataDir_/timeName/sampleName_
48  /pTraits<Type>::typeName + Field<Type>::typeName
49  /fieldTableName_
50  );
51 
52  if (exists(fieldFileName))
53  {
54  return fieldFileName;
55  }
56  else if (exists(typeFieldFileName))
57  {
58  return typeFieldFileName;
59  }
60  else
61  {
63  << "Cannot find field file "
64  << fieldFileName << " " << typeFieldFileName
65  << exit(FatalError);
66 
67  return fileName::null;
68  }
69 }
70 
71 
72 template<class Type>
74 {
75  // Initialise
76  if (mapperPtr_.empty())
77  {
78  // Reread values and interpolate
79  const fileName samplePointsFile(dataDir_/pointsName_);
80 
81  pointField samplePoints((IFstream(samplePointsFile)()));
82 
83  if (debug)
84  {
85  Info<< "timeVaryingMappedFixedValueFvPatchField :"
86  << " Read " << samplePoints.size() << " sample points from "
87  << samplePointsFile << endl;
88  }
89 
90 
91  // tbd: run-time selection
92  bool nearestOnly
93  (
94  !mapMethod_.empty()
95  && mapMethod_ != "planarInterpolation"
96  );
97 
98  // Allocate the interpolator
99  mapperPtr_.reset
100  (
101  new pointToPointPlanarInterpolation
102  (
103  samplePoints,
104  this->patch().patch().faceCentres(),
105  perturb_,
106  nearestOnly
107  )
108  );
109 
110  // Read the times for which data is available
111  sampleTimes_ = Time::findTimes(dataDir_);
112 
113  if (debug)
114  {
115  Info<< "timeVaryingMappedFixedValueFvPatchField : In directory "
116  << dataDir_ << " found times "
117  << pointToPointPlanarInterpolation::timeNames(sampleTimes_)
118  << endl;
119  }
120  }
121 
122 
123  // Find current time in sampleTimes
124  label lo = -1;
125  label hi = -1;
126 
127  bool foundTime = mapperPtr_().findTime
128  (
129  sampleTimes_,
130  startSampleTime_,
131  this->db().time().value(),
132  lo,
133  hi
134  );
135 
136  if (!foundTime)
137  {
139  << "Cannot find starting sampling values for current time "
140  << this->db().time().value() << nl
141  << "Have sampling values for times "
142  << pointToPointPlanarInterpolation::timeNames(sampleTimes_) << nl
143  << "In directory " << dataDir_ << " of field " << fieldTableName_
144  << exit(FatalError);
145  }
146 
147 
148  // Update sampled data fields.
149 
150  if (lo != startSampleTime_)
151  {
152  startSampleTime_ = lo;
153 
154  if (startSampleTime_ == endSampleTime_)
155  {
156  // No need to reread since are end values
157  if (debug)
158  {
159  Pout<< "checkTable : Setting startValues to (already read) "
160  << dataDir_/sampleTimes_[startSampleTime_].name()
161  << endl;
162  }
163  startSampledValues_ = endSampledValues_;
164  startAverage_ = endAverage_;
165  }
166  else
167  {
168  if (debug)
169  {
170  Pout<< "checkTable : Reading startValues from "
171  << dataDir_/sampleTimes_[lo].name()
172  << endl;
173  }
174 
175  // Reread values and interpolate
176  const fileName valsFile
177  (
178  findFieldFile(sampleTimes_[startSampleTime_].name())
179  );
180 
181  Field<Type> vals;
182 
183  if (setAverage_)
184  {
185  AverageField<Type> avals((IFstream(valsFile)()));
186  vals = avals;
187  startAverage_ = avals.average();
188  }
189  else
190  {
191  IFstream(valsFile)() >> vals;
192  }
193 
194  if (vals.size() != mapperPtr_().sourceSize())
195  {
197  << "Number of values (" << vals.size()
198  << ") differs from the number of points ("
199  << mapperPtr_().sourceSize()
200  << ") in file " << valsFile << exit(FatalError);
201  }
202 
203  startSampledValues_ = mapperPtr_().interpolate(vals);
204  }
205  }
206 
207  if (hi != endSampleTime_)
208  {
209  endSampleTime_ = hi;
210 
211  if (endSampleTime_ == -1)
212  {
213  // endTime no longer valid. Might as well clear endValues.
214  if (debug)
215  {
216  Pout<< "checkTable : Clearing endValues" << endl;
217  }
218  endSampledValues_.clear();
219  }
220  else
221  {
222  if (debug)
223  {
224  Pout<< "checkTable : Reading endValues from "
225  << dataDir_/sampleTimes_[endSampleTime_].name()
226  << endl;
227  }
228 
229  // Reread values and interpolate
230  const fileName valsFile
231  (
232  findFieldFile(sampleTimes_[endSampleTime_].name())
233  );
234 
235  Field<Type> vals;
236 
237  if (setAverage_)
238  {
239  AverageField<Type> avals((IFstream(valsFile)()));
240  vals = avals;
241  endAverage_ = avals.average();
242  }
243  else
244  {
245  IFstream(valsFile)() >> vals;
246  }
247 
248  if (vals.size() != mapperPtr_().sourceSize())
249  {
251  << "Number of values (" << vals.size()
252  << ") differs from the number of points ("
253  << mapperPtr_().sourceSize()
254  << ") in file " << valsFile << exit(FatalError);
255  }
256 
257  endSampledValues_ = mapperPtr_().interpolate(vals);
258  }
259  }
260 }
261 
262 
263 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
264 
265 template<class Type>
268 (
269  const fvPatch& p,
271 )
272 :
274  fieldTableName_(iF.name()),
275  dataDir_(this->db().time().constant()/"boundaryData"/this->patch().name()),
276  pointsName_("points"),
277  sampleName_(word::null),
278  setAverage_(false),
279  perturb_(0),
280  mapperPtr_(nullptr),
281  sampleTimes_(0),
282  startSampleTime_(-1),
283  startSampledValues_(0),
284  startAverage_(Zero),
285  endSampleTime_(-1),
286  endSampledValues_(0),
287  endAverage_(Zero),
288  offset_()
289 {}
290 
291 
292 template<class Type>
295 (
296  const fvPatch& p,
298  const dictionary& dict
299 )
300 :
301  fixedValueFvPatchField<Type>(p, iF, dict, false),
302  fieldTableName_(dict.lookupOrDefault("fieldTable", iF.name())),
303  dataDir_
304  (
305  dict.lookupOrDefault
306  (
307  "dataDir",
308  this->db().time().constant()/"boundaryData"/this->patch().name()
309  )
310  ),
311  pointsName_(dict.lookupOrDefault<fileName>("points", "points")),
312  sampleName_(dict.lookupOrDefault("sample", word::null)),
313  setAverage_(dict.lookupOrDefault("setAverage", false)),
314  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
315  mapMethod_
316  (
317  dict.lookupOrDefault<word>
318  (
319  "mapMethod",
320  "planarInterpolation"
321  )
322  ),
323  mapperPtr_(nullptr),
324  sampleTimes_(0),
325  startSampleTime_(-1),
326  startSampledValues_(0),
327  startAverage_(Zero),
328  endSampleTime_(-1),
329  endSampledValues_(0),
330  endAverage_(Zero),
331  offset_()
332 {
333  dataDir_.expand();
334  pointsName_.expand();
335  sampleName_.expand();
336 
337  if (dict.found("offset"))
338  {
339  offset_ = Function1<Type>::New("offset", dict);
340  }
341 
342  if
343  (
344  mapMethod_ != "planarInterpolation"
345  && mapMethod_ != "nearest"
346  )
347  {
349  (
350  dict
351  ) << "mapMethod should be one of 'planarInterpolation'"
352  << ", 'nearest'" << exit(FatalIOError);
353  }
354 
355  if (dict.found("value"))
356  {
357  fvPatchField<Type>::operator==(Field<Type>("value", dict, p.size()));
358  }
359  else
360  {
361  // Note: we use evaluate() here to trigger updateCoeffs followed
362  // by re-setting of fvatchfield::updated_ flag. This is
363  // so if first use is in the next time step it retriggers
364  // a new update.
365  this->evaluate(Pstream::commsTypes::blocking);
366  }
367 }
368 
369 
370 template<class Type>
373 (
375  const fvPatch& p,
377  const fvPatchFieldMapper& mapper
378 )
379 :
380  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
381  fieldTableName_(ptf.fieldTableName_),
382  dataDir_(ptf.dataDir_),
383  pointsName_(ptf.pointsName_),
384  sampleName_(ptf.sampleName_),
385  setAverage_(ptf.setAverage_),
386  perturb_(ptf.perturb_),
387  mapMethod_(ptf.mapMethod_),
388  mapperPtr_(nullptr),
389  sampleTimes_(0),
390  startSampleTime_(-1),
391  startSampledValues_(0),
392  startAverage_(Zero),
393  endSampleTime_(-1),
394  endSampledValues_(0),
395  endAverage_(Zero),
396  offset_(ptf.offset_, false)
397 {}
398 
399 
400 template<class Type>
403 (
405 )
406 :
408  fieldTableName_(ptf.fieldTableName_),
409  dataDir_(ptf.dataDir_),
410  pointsName_(ptf.pointsName_),
411  sampleName_(ptf.sampleName_),
412  setAverage_(ptf.setAverage_),
413  perturb_(ptf.perturb_),
414  mapMethod_(ptf.mapMethod_),
415  mapperPtr_(nullptr),
416  sampleTimes_(ptf.sampleTimes_),
417  startSampleTime_(ptf.startSampleTime_),
418  startSampledValues_(ptf.startSampledValues_),
419  startAverage_(ptf.startAverage_),
420  endSampleTime_(ptf.endSampleTime_),
421  endSampledValues_(ptf.endSampledValues_),
422  endAverage_(ptf.endAverage_),
423  offset_(ptf.offset_, false)
424 {}
425 
426 
427 template<class Type>
430 (
433 )
434 :
436  fieldTableName_(ptf.fieldTableName_),
437  dataDir_(ptf.dataDir_),
438  pointsName_(ptf.pointsName_),
439  sampleName_(ptf.sampleName_),
440  setAverage_(ptf.setAverage_),
441  perturb_(ptf.perturb_),
442  mapMethod_(ptf.mapMethod_),
443  mapperPtr_(nullptr),
444  sampleTimes_(ptf.sampleTimes_),
445  startSampleTime_(ptf.startSampleTime_),
446  startSampledValues_(ptf.startSampledValues_),
447  startAverage_(ptf.startAverage_),
448  endSampleTime_(ptf.endSampleTime_),
449  endSampledValues_(ptf.endSampledValues_),
450  endAverage_(ptf.endAverage_),
451  offset_(ptf.offset_, false)
452 {}
453 
454 
455 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
456 
457 template<class Type>
459 (
460  const fvPatchFieldMapper& m
461 )
462 {
464  if (startSampledValues_.size())
465  {
466  startSampledValues_.autoMap(m);
467  endSampledValues_.autoMap(m);
468  }
469  // Clear interpolator
470  mapperPtr_.clear();
471  startSampleTime_ = -1;
472  endSampleTime_ = -1;
473 }
474 
475 
476 template<class Type>
478 (
479  const fvPatchField<Type>& ptf,
480  const labelList& addr
481 )
482 {
484 
486  refCast<const timeVaryingMappedFixedValueFvPatchField<Type>>(ptf);
487 
488  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
489  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
490 
491  // Clear interpolator
492  mapperPtr_.clear();
493  startSampleTime_ = -1;
494  endSampleTime_ = -1;
495 }
496 
497 
498 template<class Type>
500 {
501  if (this->updated())
502  {
503  return;
504  }
505 
506 
507  checkTable();
508 
509  // Interpolate between the sampled data
510 
511  Type wantedAverage;
512 
513  if (endSampleTime_ == -1)
514  {
515  // Only start value
516  if (debug)
517  {
518  Pout<< "updateCoeffs : Sampled, non-interpolated values"
519  << " from start time:"
520  << sampleTimes_[startSampleTime_].name() << nl;
521  }
522 
523  this->operator==(startSampledValues_);
524  wantedAverage = startAverage_;
525  }
526  else
527  {
528  scalar start = sampleTimes_[startSampleTime_].value();
529  scalar end = sampleTimes_[endSampleTime_].value();
530 
531  scalar s = (this->db().time().value() - start)/(end - start);
532 
533  if (debug)
534  {
535  Pout<< "updateCoeffs : Sampled, interpolated values"
536  << " between start time:"
537  << sampleTimes_[startSampleTime_].name()
538  << " and end time:" << sampleTimes_[endSampleTime_].name()
539  << " with weight:" << s << endl;
540  }
541 
542  this->operator==((1 - s)*startSampledValues_ + s*endSampledValues_);
543  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
544  }
545 
546  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
547  // offsetting.
548  if (setAverage_)
549  {
550  const Field<Type>& fld = *this;
551 
552  Type averagePsi =
553  gSum(this->patch().magSf()*fld)
554  /gSum(this->patch().magSf());
555 
556  if (debug)
557  {
558  Pout<< "updateCoeffs :"
559  << " actual average:" << averagePsi
560  << " wanted average:" << wantedAverage
561  << endl;
562  }
563 
564  if (mag(averagePsi) < vSmall)
565  {
566  // Field too small to scale. Offset instead.
567  const Type offset = wantedAverage - averagePsi;
568  if (debug)
569  {
570  Pout<< "updateCoeffs :"
571  << " offsetting with:" << offset << endl;
572  }
573  this->operator==(fld + offset);
574  }
575  else
576  {
577  const scalar scale = mag(wantedAverage)/mag(averagePsi);
578 
579  if (debug)
580  {
581  Pout<< "updateCoeffs :"
582  << " scaling with:" << scale << endl;
583  }
584  this->operator==(scale*fld);
585  }
586  }
587 
588  // Apply offset to mapped values
589  if (offset_.valid())
590  {
591  const scalar t = this->db().time().timeOutputValue();
592  this->operator==(*this + offset_->value(t));
593  }
594 
595  if (debug)
596  {
597  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
598  << " max:" << gMax(*this)
599  << " avg:" << gAverage(*this) << endl;
600  }
601 
603 }
604 
605 
606 template<class Type>
608 (
609  Ostream& os
610 ) const
611 {
613 
614  this->writeEntryIfDifferent
615  (
616  os,
617  "dataDir",
618  this->db().time().constant()/"boundaryData"/this->patch().name(),
619  dataDir_
620  );
621 
622  this->writeEntryIfDifferent(os, "points", fileName("points"), pointsName_);
623 
624  this->writeEntryIfDifferent(os, "sample", fileName::null, sampleName_);
625 
626  this->writeEntryIfDifferent(os, "setAverage", Switch(false), setAverage_);
627 
628  this->writeEntryIfDifferent(os, "perturb", scalar(1e-5), perturb_);
629 
630  this->writeEntryIfDifferent
631  (
632  os,
633  "fieldTable",
634  this->internalField().name(),
635  fieldTableName_
636  );
637 
638  this->writeEntryIfDifferent
639  (
640  os,
641  "mapMethod",
642  word("planarInterpolation"),
643  mapMethod_
644  );
645 
646  if (offset_.valid())
647  {
648  offset_->writeData(os);
649  }
650 
651  this->writeEntry("value", os);
652 }
653 
654 
655 // ************************************************************************* //
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:297
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
#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:256
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))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:98
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:179
static const zero Zero
Definition: zero.H:97
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:265
Type gMax(const FieldField< Field, Type > &f)
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.
#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 > &)
volScalarField & p
bool exists(const fileName &, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:509
IOerror FatalIOError