codedFixedValueFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "dynamicCode.H"
28 #include "dynamicCodeContext.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  {"code", "codeInclude", "localCode"}
37 );
38 
39 template<class Type>
41 (
42  {word::null, word::null, word::null}
43 );
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 template<class Type>
50 (
51  dynamicCode& dynCode,
52  const dynamicCodeContext& context
53 ) const
54 {
55  dynCode.setFilterVariable("typeName", codeName());
56 
57  // Set TemplateType and FieldType filter variables
58  // (for fvPatchField)
59  word fieldType(pTraits<Type>::typeName);
60 
61  // Template type for fvPatchField
62  dynCode.setFilterVariable("TemplateType", fieldType);
63 
64  // Name for fvPatchField - eg, ScalarField, VectorField, ...
65  fieldType[0] = toupper(fieldType[0]);
66  dynCode.setFilterVariable("FieldType", fieldType + "Field");
67 
68  // Compile filtered C template
69  dynCode.addCompileFile(codeTemplateC("codedFixedValueFvPatchField"));
70 
71  // Copy filtered H template
72  dynCode.addCopyFile(codeTemplateH("codedFixedValueFvPatchField"));
73 
74  // Make verbose if debugging
75  dynCode.setFilterVariable("verbose", Foam::name(bool(debug)));
76 
77  if (debug)
78  {
79  Info<<"compile " << codeName() << " sha1: " << context.sha1() << endl;
80  }
81 
82  // Define Make/options
83  dynCode.setMakeOptions
84  (
85  "EXE_INC = -g \\\n"
86  "-I$(LIB_SRC)/finiteVolume/lnInclude \\\n"
87  + context.options()
88  + "\n\nLIB_LIBS = \\\n"
89  + " -lOpenFOAM \\\n"
90  + " -lfiniteVolume \\\n"
91  + context.libs()
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
98 template<class Type>
100 (
101  const fvPatch& p,
103  const dictionary& dict
104 )
105 :
106  fixedValueFvPatchField<Type>(p, iF, dict),
107  codedBase(dict, codeKeys, codeDictVars)
108 {
109  // Compile the library containing user-defined fvPatchField
111 }
112 
113 
114 template<class Type>
116 (
118  const fvPatch& p,
120  const fieldMapper& mapper
121 )
122 :
123  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
124  codedBase(ptf)
125 {}
126 
127 
128 template<class Type>
130 (
133 )
134 :
135  fixedValueFvPatchField<Type>(ptf, iF),
136  codedBase(ptf)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
142 template<class Type>
145 {
146  if (!redirectPatchFieldPtr_.valid())
147  {
148  OStringStream os;
149  writeEntry(os, "type", codeName());
150 
151  if (this->overridesConstraint())
152  {
153  writeEntry(os, "patchType", this->patch().type());
154  }
155 
156  writeEntry(os, "value", *this);
157 
158  IStringStream is(os.str());
159  dictionary dict(is);
160 
161  redirectPatchFieldPtr_.set
162  (
164  (
165  this->patch(),
166  this->internalField(),
167  dict
168  ).ptr()
169  );
170  }
171 
172  return redirectPatchFieldPtr_();
173 }
174 
175 
176 template<class Type>
178 {
179  if (this->updated())
180  {
181  return;
182  }
183 
184  const fvPatchField<Type>& fvp = redirectPatchField();
185 
186  const_cast<fvPatchField<Type>&>(fvp).updateCoeffs();
187 
188  // Copy through value
189  this->operator==(fvp);
190 
192 }
193 
194 
195 template<class Type>
197 (
198  const Pstream::commsTypes commsType
199 )
200 {
201  const fvPatchField<Type>& fvp = redirectPatchField();
202 
203  const_cast<fvPatchField<Type>&>(fvp).evaluate(commsType);
204 
206 }
207 
208 
209 template<class Type>
211 {
213  codedBase::write(os);
214 }
215 
216 
217 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Input from memory buffer stream.
Definition: IStringStream.H:52
Output to memory buffer stream.
Definition: OStringStream.H:52
string str() const
Return the string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
commsTypes
Types of communications.
Definition: UPstream.H:65
Base class for function objects and boundary conditions using dynamic code.
Definition: codedBase.H:52
void write(Ostream &os) const
Write the code for restart.
Definition: codedBase.C:461
bool updateLibrary(const dictionary &dict) const
Update library from given updated dictionary as required.
Definition: codedBase.C:398
Constructs on-the-fly a new boundary condition (derived from fixedValueFvPatchField) which is then us...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
virtual void write(Ostream &) const
Write.
codedFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const fvPatchField< Type > & redirectPatchField() const
Get reference to the underlying patch.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual void write(Ostream &) const
Write.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:217
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:210
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void evaluate(GeometricField< Type, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, GeoMesh > &x)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p