43 fieldTableName_(iF.
name()),
44 setAverage_(
dict.lookupOrDefault(
"setAverage", false)),
45 perturb_(
dict.lookupOrDefault(
"perturb", 1
e-5)),
57 startSampledValues_(0),
64 if (
dict.found(
"offset"))
70 this->
db().time().userUnits(),
78 mapMethod_ !=
"planarInterpolation"
79 && mapMethod_ !=
"nearest"
85 ) <<
"mapMethod should be one of 'planarInterpolation'"
89 dict.readIfPresent(
"fieldTableName", fieldTableName_);
91 if (
dict.found(
"value"))
120 fieldTableName_(ptf.fieldTableName_),
121 setAverage_(ptf.setAverage_),
122 perturb_(ptf.perturb_),
123 mapMethod_(ptf.mapMethod_),
126 startSampleTime_(-1),
127 startSampledValues_(0),
130 endSampledValues_(0),
132 offset_(ptf.offset_, false)
145 fieldTableName_(ptf.fieldTableName_),
146 setAverage_(ptf.setAverage_),
147 perturb_(ptf.perturb_),
148 mapMethod_(ptf.mapMethod_),
149 mapperPtr_(ptf.mapperPtr_),
150 sampleTimes_(ptf.sampleTimes_),
151 startSampleTime_(ptf.startSampleTime_),
152 startSampledValues_(ptf.startSampledValues_),
153 startAverage_(ptf.startAverage_),
154 endSampleTime_(ptf.endSampleTime_),
155 endSampledValues_(ptf.endSampledValues_),
156 endAverage_(ptf.endAverage_),
157 offset_(ptf.offset_, false)
173 refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
175 mapper(startSampledValues_, tiptf.startSampledValues_);
176 mapper(endSampledValues_, tiptf.endSampledValues_);
180 startSampleTime_ = -1;
194 refCast<const timeVaryingMappedFixedValuePointPatchField<Type>>(ptf);
196 startSampledValues_.
reset(tiptf.startSampledValues_);
197 endSampledValues_.reset(tiptf.endSampledValues_);
201 startSampleTime_ = -1;
210 if (startSampleTime_ == -1 && endSampleTime_ == -1)
243 meshPts =
pointField(points0, this->patch().meshPoints());
249 this->db().time().constant()
251 /this->patch().
name()
261 && mapMethod_ !=
"planarInterpolation"
278 const fileName samplePointsDir = samplePointsFile.
path();
283 Info<<
"timeVaryingMappedFixedValuePointPatchField : In directory "
284 << samplePointsDir <<
" found times "
294 bool foundTime = mapperPtr_().findTime
298 this->db().time().value(),
306 <<
"Cannot find starting sampling values for current time "
307 << this->db().time().value() <<
nl
308 <<
"Have sampling values for times "
311 << this->db().time().constant()/
"boundaryData"/this->patch().name()
312 <<
"\n on patch " << this->patch().name()
313 <<
" of field " << fieldTableName_
320 if (lo != startSampleTime_)
322 startSampleTime_ = lo;
324 if (startSampleTime_ == endSampleTime_)
329 Pout<<
"checkTable : Setting startValues to (already read) "
331 /this->patch().
name()
332 /sampleTimes_[startSampleTime_].
name()
335 startSampledValues_ = endSampledValues_;
336 startAverage_ = endAverage_;
342 Pout<<
"checkTable : Reading startValues from "
344 /this->patch().
name()
345 /sampleTimes_[lo].
name()
352 this->db().time().constant()
354 /this->patch().
name()
355 /sampleTimes_[startSampleTime_].
name()
365 startAverage_ = avals.
average();
372 if (vals.
size() != mapperPtr_().sourceSize())
375 <<
"Number of values (" << vals.
size()
376 <<
") differs from the number of points ("
377 << mapperPtr_().sourceSize()
381 startSampledValues_ = mapperPtr_().interpolate(vals);
385 if (hi != endSampleTime_)
389 if (endSampleTime_ == -1)
394 Pout<<
"checkTable : Clearing endValues" <<
endl;
396 endSampledValues_.clear();
402 Pout<<
"checkTable : Reading endValues from "
404 /this->patch().
name()
405 /sampleTimes_[endSampleTime_].
name()
412 this->db().time().constant()
414 /this->patch().
name()
415 /sampleTimes_[endSampleTime_].
name()
432 if (vals.
size() != mapperPtr_().sourceSize())
435 <<
"Number of values (" << vals.
size()
436 <<
") differs from the number of points ("
437 << mapperPtr_().sourceSize()
441 endSampledValues_ = mapperPtr_().interpolate(vals);
461 if (endSampleTime_ == -1)
466 Pout<<
"updateCoeffs : Sampled, non-interpolated values"
467 <<
" from start time:"
468 << sampleTimes_[startSampleTime_].
name() <<
nl;
472 wantedAverage = startAverage_;
476 scalar start = sampleTimes_[startSampleTime_].value();
477 scalar end = sampleTimes_[endSampleTime_].value();
479 scalar
s = (this->db().time().value()-start)/(end-start);
483 Pout<<
"updateCoeffs : Sampled, interpolated values"
484 <<
" between start time:"
485 << sampleTimes_[startSampleTime_].
name()
486 <<
" and end time:" << sampleTimes_[endSampleTime_].
name()
487 <<
" with weight:" <<
s <<
endl;
490 this->
operator==((1-s)*startSampledValues_ +
s*endSampledValues_);
491 wantedAverage = (1-
s)*startAverage_ +
s*endAverage_;
504 Pout<<
"updateCoeffs :"
505 <<
" actual average:" << averagePsi
506 <<
" wanted average:" << wantedAverage
510 if (
mag(averagePsi) < vSmall)
513 const Type
offset = wantedAverage - averagePsi;
516 Pout<<
"updateCoeffs :"
523 const scalar scale =
mag(wantedAverage)/
mag(averagePsi);
527 Pout<<
"updateCoeffs :"
528 <<
" scaling with:" << scale <<
endl;
537 this->
operator==(*
this + offset_->value(this->db().time().value()));
542 Pout<<
"updateCoeffs : set fixedValue to min:" <<
gMin(*
this)
543 <<
" max:" <<
gMax(*
this)
567 this->internalField().
name(),
575 word(
"planarInterpolation"),
A primitive field with a separate average value.
const Type & average() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
A primitive field of type <Type> with automated input and output.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
void size(const label)
Override size to be inconsistent with allocated storage.
virtual const fileName & name() const
Return the name of the stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
static instantList findTimes(const fileName &, const word &constantName="constant")
Search a given directory for valid time directories.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
A class for handling file names.
word name() const
Return file name (part beyond last /)
fileName path() const
Return directory path name (part before last /)
A FixedValue boundary condition for pointField.
Abstract base class for point-mesh patch fields.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
const objectRegistry & db() const
Return local objectRegistry.
Basic pointPatch represents a set of points from the mesh.
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
static wordList timeNames(const instantList &)
Helper: extract words of times.
const polyMesh & mesh() const
Return the mesh reference.
Mesh consisting of general polyhedral cells.
const fileName & facesInstance() const
Return the current instance directory for faces.
const fileName & pointsInstance() const
Return the current instance directory for points.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual const pointField & points() const
Return raw points.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
A time-varying form of a mapped fixed value boundary condition.
timeVaryingMappedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void checkTable()
Find boundary data in between current time and interpolate.
virtual void reset(const pointPatchField< Type > &)
Reset the pointPatchField to the given pointPatchField.
virtual void map(const pointPatchField< Type > &, const fieldMapper &)
Map the given pointPatchField onto this pointPatchField.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const pointPatchField< Type > &)
Reset the pointPatchField to the given pointPatchField.
virtual void map(const pointPatchField< Type > &, const fieldMapper &)
Map the given pointPatchField onto this pointPatchField.
A class for handling words, derived from string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar e
Elementary charge.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
vectorField pointField
pointField is a vectorField.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Type gMin(const FieldField< Field, Type > &f)
void offset(label &lst, const label o)
Type gMax(const FieldField< Field, Type > &f)