codedFixedValuePointPatchField.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) 2012-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 
28 #include "fieldMapper.H"
29 #include "pointFields.H"
30 #include "dynamicCode.H"
31 #include "dynamicCodeContext.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  {"code", "codeInclude", "localCode"}
39 );
40 
41 template<class Type>
43 (
44  {word::null, word::null, word::null}
45 );
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 template<class Type>
52 (
53  dynamicCode& dynCode,
54  const dynamicCodeContext& context
55 ) const
56 {
57  // Take no chances - typeName must be identical to codeName()
58  dynCode.setFilterVariable("typeName", codeName());
59 
60  // Set TemplateType and FieldType filter variables
61  // (for pointPatchField)
62  word fieldType(pTraits<Type>::typeName);
63 
64  // Template type for pointPatchField
65  dynCode.setFilterVariable("TemplateType", fieldType);
66 
67  // Name for pointPatchField - eg, ScalarField, VectorField, ...
68  fieldType[0] = toupper(fieldType[0]);
69  dynCode.setFilterVariable("FieldType", fieldType + "Field");
70 
71  // Compile filtered C template
72  dynCode.addCompileFile(codeTemplateC("codedFixedValuePointPatchField"));
73 
74  // Copy filtered H template
75  dynCode.addCopyFile(codeTemplateH("codedFixedValuePointPatchField"));
76 
77  // Make verbose if debugging
78  dynCode.setFilterVariable("verbose", Foam::name(bool(debug)));
79 
80  if (debug)
81  {
82  Info<<"compile " << codeName() << " sha1: " << context.sha1() << endl;
83  }
84 
85  // Define Make/options
86  dynCode.setMakeOptions
87  (
88  "EXE_INC = -g \\\n"
89  "-I$(LIB_SRC)/finiteVolume/lnInclude \\\n"
90  + context.options()
91  + "\n\nLIB_LIBS = \\\n"
92  + " -lOpenFOAM \\\n"
93  + " -lfiniteVolume \\\n"
94  + context.libs()
95  );
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
101 template<class Type>
103 (
104  const pointPatch& p,
106  const dictionary& dict
107 )
108 :
109  fixedValuePointPatchField<Type>(p, iF, dict),
110  codedBase(dict, codeKeys, codeDictVars)
111 {
112  // Compile the library containing user-defined pointPatchField
114 }
115 
116 
117 template<class Type>
119 (
121  const pointPatch& p,
123  const fieldMapper& mapper
124 )
125 :
126  fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
127  codedBase(ptf)
128 {}
129 
130 
131 template<class Type>
133 (
136 )
137 :
138  fixedValuePointPatchField<Type>(ptf, iF),
139  codedBase(ptf)
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
145 template<class Type>
148 {
149  if (!redirectPatchFieldPtr_.valid())
150  {
151  OStringStream os;
152  writeEntry(os, "type", codeName());
153  writeEntry(os, "value", static_cast<const Field<Type>&>(*this));
154  IStringStream is(os.str());
155  dictionary dict(is);
156 
157  redirectPatchFieldPtr_.set
158  (
160  (
161  this->patch(),
162  this->internalField(),
163  dict
164  ).ptr()
165  );
166  }
167 
168  return redirectPatchFieldPtr_();
169 }
170 
171 
172 template<class Type>
174 {
175  if (this->updated())
176  {
177  return;
178  }
179 
180  const pointPatchField<Type>& fvp = redirectPatchField();
181 
182  const_cast<pointPatchField<Type>&>(fvp).updateCoeffs();
183 
184  // Copy through value
185  this->operator==(fvp);
186 
188 }
189 
190 
191 template<class Type>
193 (
194  const Pstream::commsTypes commsType
195 )
196 {
197  const pointPatchField<Type>& fvp = redirectPatchField();
198 
199  const_cast<pointPatchField<Type>&>(fvp).evaluate(commsType);
200 
202 }
203 
204 
205 template<class Type>
207 {
209  codedBase::write(os);
210 }
211 
212 
213 // ************************************************************************* //
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...
Pre-declare SubField and related Field type.
Definition: Field.H:83
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 fixedValuePointPatchField) which is then...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
codedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
const pointPatchField< Type > & redirectPatchField() const
Get reference to the underlying patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
A FixedValue boundary condition for pointField.
Abstract base class for point-mesh patch fields.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
dictionary dict
volScalarField & p