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-2016 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 "AverageIOField.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  fieldTableName_(iF.name()),
42  setAverage_(false),
43  perturb_(0),
44  mapperPtr_(NULL),
45  sampleTimes_(0),
46  startSampleTime_(-1),
47  startSampledValues_(0),
48  startAverage_(Zero),
49  endSampleTime_(-1),
50  endSampledValues_(0),
51  endAverage_(Zero),
52  offset_()
53 {}
54 
55 
56 template<class Type>
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
66  fieldTableName_(iF.name()),
67  setAverage_(readBool(dict.lookup("setAverage"))),
68  perturb_(dict.lookupOrDefault("perturb", 1e-5)),
69  mapMethod_
70  (
71  dict.lookupOrDefault<word>
72  (
73  "mapMethod",
74  "planarInterpolation"
75  )
76  ),
77  mapperPtr_(NULL),
78  sampleTimes_(0),
79  startSampleTime_(-1),
80  startSampledValues_(0),
81  startAverage_(Zero),
82  endSampleTime_(-1),
83  endSampledValues_(0),
84  endAverage_(Zero),
85  offset_(Function1<Type>::New("offset", dict))
86 {
87  if
88  (
89  mapMethod_ != "planarInterpolation"
90  && mapMethod_ != "nearest"
91  )
92  {
94  (
95  dict
96  ) << "mapMethod should be one of 'planarInterpolation'"
97  << ", 'nearest'" << exit(FatalIOError);
98  }
99 
100 
101  dict.readIfPresent("fieldTable", fieldTableName_);
102 
103  if (dict.found("value"))
104  {
105  fvPatchField<Type>::operator==(Field<Type>("value", dict, p.size()));
106  }
107  else
108  {
109  // Note: we use evaluate() here to trigger updateCoeffs followed
110  // by re-setting of fvatchfield::updated_ flag. This is
111  // so if first use is in the next time step it retriggers
112  // a new update.
113  this->evaluate(Pstream::blocking);
114  }
115 }
116 
117 
118 template<class Type>
121 (
123  const fvPatch& p,
125  const fvPatchFieldMapper& mapper
126 )
127 :
128  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
129  fieldTableName_(ptf.fieldTableName_),
130  setAverage_(ptf.setAverage_),
131  perturb_(ptf.perturb_),
132  mapMethod_(ptf.mapMethod_),
133  mapperPtr_(NULL),
134  sampleTimes_(0),
135  startSampleTime_(-1),
136  startSampledValues_(0),
137  startAverage_(Zero),
138  endSampleTime_(-1),
139  endSampledValues_(0),
140  endAverage_(Zero),
141  offset_(ptf.offset_, false)
142 {}
143 
144 
145 template<class Type>
148 (
150 )
151 :
153  fieldTableName_(ptf.fieldTableName_),
154  setAverage_(ptf.setAverage_),
155  perturb_(ptf.perturb_),
156  mapMethod_(ptf.mapMethod_),
157  mapperPtr_(NULL),
158  sampleTimes_(ptf.sampleTimes_),
159  startSampleTime_(ptf.startSampleTime_),
160  startSampledValues_(ptf.startSampledValues_),
161  startAverage_(ptf.startAverage_),
162  endSampleTime_(ptf.endSampleTime_),
163  endSampledValues_(ptf.endSampledValues_),
164  endAverage_(ptf.endAverage_),
165  offset_(ptf.offset_, false)
166 {}
167 
168 
169 template<class Type>
172 (
175 )
176 :
178  fieldTableName_(ptf.fieldTableName_),
179  setAverage_(ptf.setAverage_),
180  perturb_(ptf.perturb_),
181  mapMethod_(ptf.mapMethod_),
182  mapperPtr_(NULL),
183  sampleTimes_(ptf.sampleTimes_),
184  startSampleTime_(ptf.startSampleTime_),
185  startSampledValues_(ptf.startSampledValues_),
186  startAverage_(ptf.startAverage_),
187  endSampleTime_(ptf.endSampleTime_),
188  endSampledValues_(ptf.endSampledValues_),
189  endAverage_(ptf.endAverage_),
190  offset_(ptf.offset_, false)
191 {}
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
196 template<class Type>
198 (
199  const fvPatchFieldMapper& m
200 )
201 {
203  if (startSampledValues_.size())
204  {
205  startSampledValues_.autoMap(m);
206  endSampledValues_.autoMap(m);
207  }
208  // Clear interpolator
209  mapperPtr_.clear();
210  startSampleTime_ = -1;
211  endSampleTime_ = -1;
212 }
213 
214 
215 template<class Type>
217 (
218  const fvPatchField<Type>& ptf,
219  const labelList& addr
220 )
221 {
223 
225  refCast<const timeVaryingMappedFixedValueFvPatchField<Type>>(ptf);
226 
227  startSampledValues_.rmap(tiptf.startSampledValues_, addr);
228  endSampledValues_.rmap(tiptf.endSampledValues_, addr);
229 
230  // Clear interpolator
231  mapperPtr_.clear();
232  startSampleTime_ = -1;
233  endSampleTime_ = -1;
234 }
235 
236 
237 template<class Type>
239 {
240  // Initialise
241  if (mapperPtr_.empty())
242  {
243  pointIOField samplePoints
244  (
245  IOobject
246  (
247  "points",
248  this->db().time().constant(),
249  "boundaryData"/this->patch().name(),
250  this->db(),
251  IOobject::MUST_READ,
252  IOobject::AUTO_WRITE,
253  false
254  )
255  );
256 
257  const fileName samplePointsFile = samplePoints.filePath();
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
362  (
363  IOobject
364  (
365  fieldTableName_,
366  this->db().time().constant(),
367  "boundaryData"
368  /this->patch().name()
369  /sampleTimes_[startSampleTime_].name(),
370  this->db(),
371  IOobject::MUST_READ,
372  IOobject::AUTO_WRITE,
373  false
374  )
375  );
376 
377  if (vals.size() != mapperPtr_().sourceSize())
378  {
380  << "Number of values (" << vals.size()
381  << ") differs from the number of points ("
382  << mapperPtr_().sourceSize()
383  << ") in file " << vals.objectPath() << exit(FatalError);
384  }
385 
386  startAverage_ = vals.average();
387  startSampledValues_ = mapperPtr_().interpolate(vals);
388  }
389  }
390 
391  if (hi != endSampleTime_)
392  {
393  endSampleTime_ = hi;
394 
395  if (endSampleTime_ == -1)
396  {
397  // endTime no longer valid. Might as well clear endValues.
398  if (debug)
399  {
400  Pout<< "checkTable : Clearing endValues" << endl;
401  }
402  endSampledValues_.clear();
403  }
404  else
405  {
406  if (debug)
407  {
408  Pout<< "checkTable : Reading endValues from "
409  << "boundaryData"
410  /this->patch().name()
411  /sampleTimes_[endSampleTime_].name()
412  << endl;
413  }
414 
415  // Reread values and interpolate
417  (
418  IOobject
419  (
420  fieldTableName_,
421  this->db().time().constant(),
422  "boundaryData"
423  /this->patch().name()
424  /sampleTimes_[endSampleTime_].name(),
425  this->db(),
426  IOobject::MUST_READ,
427  IOobject::AUTO_WRITE,
428  false
429  )
430  );
431 
432  if (vals.size() != mapperPtr_().sourceSize())
433  {
435  << "Number of values (" << vals.size()
436  << ") differs from the number of points ("
437  << mapperPtr_().sourceSize()
438  << ") in file " << vals.objectPath() << exit(FatalError);
439  }
440 
441  endAverage_ = vals.average();
442  endSampledValues_ = mapperPtr_().interpolate(vals);
443  }
444  }
445 }
446 
447 
448 template<class Type>
450 {
451  if (this->updated())
452  {
453  return;
454  }
455 
456 
457  checkTable();
458 
459  // Interpolate between the sampled data
460 
461  Type wantedAverage;
462 
463  if (endSampleTime_ == -1)
464  {
465  // only start value
466  if (debug)
467  {
468  Pout<< "updateCoeffs : Sampled, non-interpolated values"
469  << " from start time:"
470  << sampleTimes_[startSampleTime_].name() << nl;
471  }
472 
473  this->operator==(startSampledValues_);
474  wantedAverage = startAverage_;
475  }
476  else
477  {
478  scalar start = sampleTimes_[startSampleTime_].value();
479  scalar end = sampleTimes_[endSampleTime_].value();
480 
481  scalar s = (this->db().time().value() - start)/(end - start);
482 
483  if (debug)
484  {
485  Pout<< "updateCoeffs : Sampled, interpolated values"
486  << " between start time:"
487  << sampleTimes_[startSampleTime_].name()
488  << " and end time:" << sampleTimes_[endSampleTime_].name()
489  << " with weight:" << s << endl;
490  }
491 
492  this->operator==((1 - s)*startSampledValues_ + s*endSampledValues_);
493  wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
494  }
495 
496  // Enforce average. Either by scaling (if scaling factor > 0.5) or by
497  // offsetting.
498  if (setAverage_)
499  {
500  const Field<Type>& fld = *this;
501 
502  Type averagePsi =
503  gSum(this->patch().magSf()*fld)
504  /gSum(this->patch().magSf());
505 
506  if (debug)
507  {
508  Pout<< "updateCoeffs :"
509  << " actual average:" << averagePsi
510  << " wanted average:" << wantedAverage
511  << endl;
512  }
513 
514  if (mag(averagePsi) < VSMALL)
515  {
516  // Field too small to scale. Offset instead.
517  const Type offset = wantedAverage - averagePsi;
518  if (debug)
519  {
520  Pout<< "updateCoeffs :"
521  << " offsetting with:" << offset << endl;
522  }
523  this->operator==(fld + offset);
524  }
525  else
526  {
527  const scalar scale = mag(wantedAverage)/mag(averagePsi);
528 
529  if (debug)
530  {
531  Pout<< "updateCoeffs :"
532  << " scaling with:" << scale << endl;
533  }
534  this->operator==(scale*fld);
535  }
536  }
537 
538  // apply offset to mapped values
539  const scalar t = this->db().time().timeOutputValue();
540  this->operator==(*this + offset_->value(t));
541 
542  if (debug)
543  {
544  Pout<< "updateCoeffs : set fixedValue to min:" << gMin(*this)
545  << " max:" << gMax(*this)
546  << " avg:" << gAverage(*this) << endl;
547  }
548 
550 }
551 
552 
553 template<class Type>
555 (
556  Ostream& os
557 ) const
558 {
560  os.writeKeyword("setAverage") << setAverage_ << token::END_STATEMENT << nl;
561  if (perturb_ != 1e-5)
562  {
563  os.writeKeyword("perturb") << perturb_ << token::END_STATEMENT << nl;
564  }
565 
566  if (fieldTableName_ != this->internalField().name())
567  {
568  os.writeKeyword("fieldTable") << fieldTableName_
569  << token::END_STATEMENT << nl;
570  }
571 
572  if
573  (
574  (
575  !mapMethod_.empty()
576  && mapMethod_ != "planarInterpolation"
577  )
578  )
579  {
580  os.writeKeyword("mapMethod") << mapMethod_
581  << token::END_STATEMENT << nl;
582  }
583 
584  offset_->writeData(os);
585 
586  this->writeEntry("value", os);
587 }
588 
589 
590 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
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
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.
A primitive field + average with IO.
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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
bool readBool(Istream &)
Definition: boolIO.C:60
const Type & average() const
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:179
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))
fileName filePath() const
Return complete path + object name if the file exists.
Definition: IOobject.C:299
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
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
Foam::fvPatchFieldMapper.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
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
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
Type gMax(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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)
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:249
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
runTime write()
const word & name() const
Return name.
Definition: IOobject.H:260
Interpolates between two sets of unstructured points using 2D Delaunay triangulation. Used in e.g. timeVaryingMapped bcs.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
IOerror FatalIOError