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 -------------------------------------------------------------------------------
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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  refValueFunc_(nullptr)
42 {}
43 
44 
45 template<class Type>
47 (
48  const fvPatch& p,
50  const Field<Type>& fld
51 )
52 :
54  refValueFunc_(nullptr)
55 {}
56 
57 
58 template<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 
80 template<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 
105 template<class Type>
107 (
108  const uniformFixedValueFvPatchField<Type>& ptf
109 )
110 :
112  refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch()))
113 {}
114 
115 
116 template<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 
130 template<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 
151 template<class Type>
153 (
154  const fvPatchField<Type>& ptf,
155  const labelList& addr
156 )
157 {
159 
160  const uniformFixedValueFvPatchField& tiptf =
161  refCast<const uniformFixedValueFvPatchField>(ptf);
162 
163  if (refValueFunc_ && tiptf.refValueFunc_)
164  {
165  refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
166  }
167 }
168 
169 
170 template<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 
185 template<class Type>
187 {
189  if (refValueFunc_)
190  {
191  refValueFunc_->writeData(os);
192  }
194 }
195 
196 
197 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:236
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
This boundary condition provides a uniform fixed value condition.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:546
Generic templated field type.
Definition: Field.H:62
A FieldMapper for finite-volume patch fields.
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 rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:299
uniformFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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...
virtual void write(Ostream &os) const
Write.
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:329
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
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.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.