49Foam::uniformFixedValuePointPatchField<Type>::getPolyPatch(
const pointPatch&
p)
51 const polyMesh&
mesh =
p.boundaryMesh().mesh()();
52 return mesh.boundaryMesh().cfindPatch(
p.name());
83 refValueFunc_ = PatchFunction1<Type>::New
92 refPointValueFunc_ = Function1<Type>::New
114 const uniformFixedValuePointPatchField<Type>& ptf,
116 const DimensionedField<Type, pointMesh>& iF,
117 const pointPatchFieldMapper& mapper
120 fixedValuePointPatchField<Type>(ptf,
p, iF, mapper)
122 if (
const polyPatch*
pp = this->getPolyPatch(this->
patch()))
124 refValueFunc_ = ptf.refValueFunc_.clone(*
pp);
127 refPointValueFunc_ = ptf.refPointValueFunc_.clone();
129 if (mapper.direct() && !mapper.hasUnmapped())
132 this->map(ptf, mapper);
154 refValueFunc_ = pfld.refValueFunc_.clone(*
pp);
157 refPointValueFunc_ = pfld.refPointValueFunc_.clone();
171 bool canEvaluate(
false);
175 refValueFunc_().autoMap(mapper);
178 if (refValueFunc_->constant())
183 if (refPointValueFunc_)
186 if (refPointValueFunc_->constant())
202 const pointPatchField<Type>& ptf,
210 if (refValueFunc_ && tiptf.refValueFunc_)
212 refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
224 const scalar t = this->db().time().timeOutputValue();
246 refValueFunc_->writeData(
os);
248 else if (refPointValueFunc_)
250 refPointValueFunc_->writeData(
os);
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual bool hasUnmapped() const =0
Any unmapped values?
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
const pointPatch & patch() const noexcept
Return the patch.
const objectRegistry & db() const
The associated objectRegistry.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
const DimensionedField< Type, pointMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Basic pointPatch represents a set of points from the mesh.
A patch is a list of labels that address the faces in the global face list.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void operator=(const valuePointPatchField< Type > &)
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
OBJstream os(runTime.globalPath()/outputName)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< label > labelList
A List of labels.