Loading...
Searching...
No Matches
uniformFixedValueFvPatchField.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const fvPatch& p,
38)
40 fixedValueFvPatchField<Type>(p, iF),
41 refValueFunc_(nullptr)
42{}
43
44
45template<class Type>
47(
48 const fvPatch& p,
50 const Field<Type>& fld
51)
53 fixedValueFvPatchField<Type>(p, iF, fld),
54 refValueFunc_(nullptr)
55{}
56
57
58template<class Type>
60(
61 const fvPatch& p,
63 const dictionary& dict
64)
65:
66 fixedValueFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
67 refValueFunc_(PatchFunction1<Type>::New(p.patch(), "uniformValue", dict))
68{
69 if (!this->readValueEntry(dict))
70 {
71 // Ensure field has reasonable initial values
72 this->extrapolateInternal();
73
74 // Evaluate to assign a value
75 this->evaluate();
76 }
77}
78
79
80template<class Type>
82(
83 const uniformFixedValueFvPatchField<Type>& ptf,
84 const fvPatch& p,
85 const DimensionedField<Type, volMesh>& iF,
86 const fvPatchFieldMapper& mapper
87)
88:
89 fixedValueFvPatchField<Type>(p, iF), // Don't map
90 refValueFunc_(ptf.refValueFunc_.clone(p.patch()))
91{
92 if (mapper.direct() && !mapper.hasUnmapped())
93 {
94 // Use mapping instead of re-evaluation
95 this->map(ptf, mapper);
96 }
97 else
98 {
99 // Evaluate since value not mapped
100 this->evaluate();
101 }
102}
103
104
105template<class Type>
107(
109)
111 fixedValueFvPatchField<Type>(ptf),
112 refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch()))
113{}
114
115
116template<class Type>
118(
121)
122:
123 fixedValueFvPatchField<Type>(ptf, iF),
124 refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch()))
125{}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
130template<class Type>
132(
133 const fvPatchFieldMapper& mapper
134)
135{
137
138 if (refValueFunc_)
139 {
140 refValueFunc_().autoMap(mapper);
141
142 if (refValueFunc_().constant())
143 {
144 // If mapper is not dependent on time we're ok to evaluate
145 this->evaluate();
146 }
147 }
148}
149
150
151template<class Type>
153(
154 const fvPatchField<Type>& ptf,
155 const labelList& addr
156)
157{
159
160 const uniformFixedValueFvPatchField& tiptf =
162
163 if (refValueFunc_ && tiptf.refValueFunc_)
165 refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
166 }
167}
168
169
170template<class Type>
172{
173 if (this->updated())
174 {
175 return;
176 }
177
178 const scalar t = this->db().time().timeOutputValue();
180 fvPatchField<Type>::operator==(refValueFunc_->value(t));
182}
183
184
185template<class Type>
187{
189 if (refValueFunc_)
190 {
191 refValueFunc_->writeData(os);
192 }
194}
195
196
197// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
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?
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void write(Ostream &) const
Write.
virtual void operator==(const fvPatchField< Type > &)
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
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, sets updated() to false.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a uniform fixed value condition.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
uniformFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
Different types of constants.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
dictionary dict