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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
35 // Alternative
36 // refCast<const facePointPatch>(p).patch()
37 
38 template<class Type>
39 const Foam::polyPatch&
41 {
42  const polyMesh& mesh = p.boundaryMesh().mesh()();
43  label patchi = mesh.boundaryMesh().findPatchID(p.name());
44 
45  if (patchi == -1)
46  {
48  << "Cannot use uniformFixedValue on patch " << p.name()
49  << " since there is no underlying mesh patch" << exit(FatalError);
50  }
51  return mesh.boundaryMesh()[patchi];
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
56 
57 template<class Type>
60 (
61  const pointPatch& p,
63 )
64 :
65  fixedValuePointPatchField<Type>(p, iF),
66  refValueFunc_(nullptr)
67 {}
68 
69 
70 template<class Type>
73 (
74  const pointPatch& p,
76  const dictionary& dict
77 )
78 :
79  fixedValuePointPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
80  refValueFunc_
81  (
82  PatchFunction1<Type>::New
83  (
84  this->getPatch(p),
85  "uniformValue",
86  dict,
87  false // generate point values
88  )
89  )
90 {
91  if (!this->readValueEntry(dict))
92  {
93  // Ensure field has reasonable initial values
94  this->extrapolateInternal();
95 
96  // Evaluate to assign a value
97  this->evaluate();
98  }
99 }
100 
101 
102 template<class Type>
105 (
106  const uniformFixedValuePointPatchField<Type>& ptf,
107  const pointPatch& p,
108  const DimensionedField<Type, pointMesh>& iF,
109  const pointPatchFieldMapper& mapper
110 )
111 :
112  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
113  refValueFunc_(ptf.refValueFunc_.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 
128 template<class Type>
131 (
132  const uniformFixedValuePointPatchField<Type>& ptf
133 )
134 :
135  fixedValuePointPatchField<Type>(ptf),
136  refValueFunc_(ptf.refValueFunc_.clone(this->getPatch(this->patch())))
137 {}
138 
139 
140 template<class Type>
143 (
146 )
147 :
148  fixedValuePointPatchField<Type>(ptf, iF),
149  refValueFunc_(ptf.refValueFunc_.clone(this->getPatch(this->patch())))
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
155 template<class Type>
157 (
158  const pointPatchFieldMapper& mapper
159 )
160 {
162 
163  if (refValueFunc_)
164  {
165  refValueFunc_().autoMap(mapper);
166 
167  if (refValueFunc_().constant())
168  {
169  // If mapper is not dependent on time we're ok to evaluate
170  this->evaluate();
171  }
172  }
173 }
174 
175 
176 template<class Type>
178 (
179  const pointPatchField<Type>& ptf,
180  const labelList& addr
181 )
182 {
184 
185  const uniformFixedValuePointPatchField& tiptf =
186  refCast<const uniformFixedValuePointPatchField>(ptf);
187 
188  if (refValueFunc_ && tiptf.refValueFunc_)
189  {
190  refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
191  }
192 }
193 
194 
195 template<class Type>
197 {
198  if (this->updated())
199  {
200  return;
201  }
202  const scalar t = this->db().time().timeOutputValue();
203 
204  valuePointPatchField<Type>::operator=(refValueFunc_->value(t));
206 }
207 
208 
209 template<class Type>
211 write(Ostream& os) const
212 {
213  // Note: write value
215  if (refValueFunc_)
216  {
217  refValueFunc_->writeData(os);
218  }
219 }
220 
221 
222 // ************************************************************************* //
void extrapolateInternal()
Assign the patch field from the internal field.
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A FixedValue boundary condition for pointField.
Foam::pointPatchFieldMapper.
virtual void operator=(const valuePointPatchField< Type > &)
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
dynamicFvMesh & mesh
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:303
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
OBJstream os(runTime.globalPath()/outputName)
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Enables the specification of a uniform fixed value condition.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
uniformFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.