uniformFixedValuePointPatchField.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-2017 OpenFOAM Foundation
9 Copyright (C) 2019 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#include "SubField.H"
31#include "polyPatch.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34
35template<class Type>
36const Foam::polyPatch&
38{
39 const polyMesh& mesh = p.boundaryMesh().mesh()();
40 label patchi = mesh.boundaryMesh().findPatchID(p.name());
41
42 if (patchi == -1)
43 {
45 << "Cannot use uniformFixedValue on patch " << p.name()
46 << " since there is no underlying mesh patch" << exit(FatalError);
47 }
48 return mesh.boundaryMesh()[patchi];
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
53
54template<class Type>
57(
58 const pointPatch& p,
60)
61:
63 uniformValue_()
64{}
65
66
67template<class Type>
70(
71 const pointPatch& p,
73 const dictionary& dict
74)
75:
76 fixedValuePointPatchField<Type>(p, iF, dict, false),
77 uniformValue_
78 (
79 PatchFunction1<Type>::New
80 (
81 this->getPatch(p),
82 "uniformValue",
83 dict,
84 false // generate point values
85 )
86 )
87{
88 if (dict.found("value"))
89 {
91 (
92 Field<Type>("value", dict, p.size())
93 );
94 }
95 else
96 {
97 this->evaluate();
98 }
99}
100
101
102template<class Type>
105(
107 const pointPatch& p,
109 const pointPatchFieldMapper& mapper
110)
111:
112 fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
113 uniformValue_(ptf.uniformValue_.clone(this->getPatch(p)))
114{
115 if (mapper.direct() && !mapper.hasUnmapped())
116 {
117 // Use mapping instead of re-evaluation
118 this->map(ptf, mapper);
119 }
120 else
121 {
122 // Evaluate since value not mapped
123 this->evaluate();
124 }
125}
126
127
128template<class Type>
131(
133)
134:
135 fixedValuePointPatchField<Type>(ptf),
136 uniformValue_(ptf.uniformValue_.clone(this->getPatch(this->patch())))
137{}
138
139
140template<class Type>
143(
146)
147:
148 fixedValuePointPatchField<Type>(ptf, iF),
149 uniformValue_(ptf.uniformValue_.clone(this->getPatch(this->patch())))
150{}
151
152
153// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154
155template<class Type>
157(
158 const pointPatchFieldMapper& mapper
159)
160{
162 uniformValue_().autoMap(mapper);
163
164 if (uniformValue_().constant())
165 {
166 // If mapper is not dependent on time we're ok to evaluate
167 this->evaluate();
168 }
169}
170
171
172template<class Type>
174(
175 const pointPatchField<Type>& ptf,
176 const labelList& addr
177)
178{
180
182 refCast<const uniformFixedValuePointPatchField>(ptf);
183
184 uniformValue_().rmap(tiptf.uniformValue_(), addr);
185}
186
187
188template<class Type>
190{
191 if (this->updated())
192 {
193 return;
194 }
195 const scalar t = this->db().time().timeOutputValue();
196 fixedValuePointPatchField<Type>::operator==(uniformValue_->value(t));
198}
199
200
201template<class Type>
203write(Ostream& os) const
204{
205 // Note: write value
207 uniformValue_->writeData(os);
208}
209
210
211// ************************************************************************* //
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.
Definition: Field.H:82
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:240
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
A FixedValue boundary condition for pointField.
virtual bool write()
Write the output fields.
constant condensation/saturation model.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
Enables the specification of a uniform fixed value boundary condition.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict