Loading...
Searching...
No Matches
patchExprDriver.H
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) 2019-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::expressions::patchExpr::parseDriver
28
29Description
30 Driver for patch expressions
31
32 In addition to the standard mathematical functions, operations and
33 logical and relational operations, the patch expressions support the
34 following driver-specific functions:
35
36 Functions
37 \table
38 Function | Description | Number of arguments |
39 pos | The face centres | 0 |
40 pts | The face points | 0 |
41 area | The face area magnitudes | 0 |
42 weightAverage| Area weighted average | 1 |
43 weightSum | Area weighted sum | 1 |
44 face | The face area normal vectors | 0 |
45 normal | The face unit normal vectors | 0 |
46 point | A point-field point value | 1 |
47 faceToPoint | Interpolate face values onto points | 1 |
48 pointToFace | Interpolate point values onto faces | 1 |
49 rand | Random field | 0/1 |
50 snGrad | Surface normal field | 0 |
51 internalField | Internal field next to patch | 0 |
52 neighbourField | patch field on opposite side of coupled patch | 0 |
53 \endtable
54
55Note
56 Use namespace debug switch \c patchExpr for scanner (2), parser (4)
57 or dictionary controls as per Foam::expressions::exprDriver.
58
59SourceFiles
60 patchExprDriver.C
61 patchExprDriverFields.C
62 patchExprDriverTemplates.C
63
64\*---------------------------------------------------------------------------*/
65
66#ifndef expressions_patchExprDriver_H
67#define expressions_patchExprDriver_H
68
69#include "patchExprFwd.H"
70#include "fvExprDriver.H"
72#include "Enum.H"
73#include "volFields.H"
74#include "surfaceFields.H"
75#include "pointFields.H"
77
78// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79
80namespace Foam
81{
82namespace expressions
83{
84namespace patchExpr
85{
86
87/*---------------------------------------------------------------------------*\
88 Class parseDriver Declaration
89\*---------------------------------------------------------------------------*/
90
91class parseDriver
92:
95{
96 // Private Member Functions
97
98 static const fvPatch& getFvPatch
99 (
100 const fvMesh& fvm,
101 const dictionary& dict
102 );
103
104
105protected:
106
107 // Protected Data
108
109 //- The referenced patch
110 const fvPatch& patch_;
111
112
113 // Protected Member Functions
114
115 //- Cell selections (as logical)
117 (
118 const word& name,
119 enum topoSetSource::sourceType setType
120 ) const;
121
122 //- Face selections (as logical)
124 (
125 const word& name,
126 enum topoSetSource::sourceType setType
127 ) const;
128
129
130public:
131
132 ClassName("patchExpr::driver");
133
134 // Generated Methods
135
136 // No copy copy construct
137 parseDriver(const parseDriver&) = delete;
138
139 // No copy assignment
140 void operator=(const parseDriver&) = delete;
141
142
143 // Constructors
144
145 //- Construct for specified patch, with dictionary information
146 explicit parseDriver
147 (
148 const fvPatch& p,
150 );
151
152 //- Construct for specified patch with copy of driver context
154 (
155 const fvPatch& p,
156 const parseDriver& driver,
157 const dictionary& dict // = dictionary::null
158 );
159
160 //- Construct with patchName for the given mesh
161 parseDriver(const word& patchName, const fvMesh& mesh);
162
163 //- Construct with "patch" (mandatory) and "region" (optional)
164 //- specified in dictionary
165 parseDriver(const dictionary& dict, const fvMesh& mesh);
166
167
168 // Not generally clonable
169
170
171 //- Destructor
172 virtual ~parseDriver() = default;
173
174
175 // Public Member Functions
176
177 //- The mesh we are attached to
178 virtual const fvMesh& mesh() const
179 {
180 return patch_.boundaryMesh().mesh();
181 }
183 //- The natural field size for the expression
184 virtual label size() const
185 {
186 return patch_.patch().size();
187 }
188
189 //- The point field size for the expression
190 virtual label pointSize() const
191 {
192 return patch_.patch().nPoints();
193 }
194
195 //- Field size associated with different geometric field types
196 inline label size(const FieldAssociation geoType) const;
197
198
199 // Evaluation
200
201 //- Perform parsing on (sub) string
202 using genericRagelLemonDriver::content;
203
204 //- Execute the parser.
205 // The return value currently has no meaning.
206 virtual unsigned parse
207 (
208 const std::string& expr,
209 size_t pos = 0,
210 size_t len = std::string::npos
211 );
212
214 // General
215
216 //- Set result
217 template<class Type>
218 void setResult(Field<Type>* ptr, bool pointVal = false)
219 {
220 result().setResult<Type>(ptr, pointVal);
221 }
222
223
224 //- Retrieve variable as field if possible.
225 // Test tmp for validity to determine success of the operation.
226 template<class Type>
227 tmp<Field<Type>> getVariableIfAvailable(const word& fldName) const;
228
229
230 // Field Retrieval
231
232 //- Return named field
233 template<class Type>
234 tmp<Field<Type>> getField(const word& fldName);
235
236 //- Retrieve field (vol field)
237 template<class Type>
238 tmp<Field<Type>> getVolField(const word& fldName);
239
240 //- Retrieve field (surface field)
241 template<class Type>
242 tmp<Field<Type>> getSurfaceField(const word& fldName);
243
244 //- Retrieve field (point field)
245 template<class Type>
246 tmp<Field<Type>> getPointField(const word& fldName);
247
248 //- Return internal field next to patch
249 template<class Type>
250 tmp<Field<Type>> patchInternalField(const word& fldName);
251
252 //- Return patchField on the opposite patch of a coupled patch
253 template<class Type>
255
256 //- Return surface normal field (snGrad)
257 template<class Type>
259
260
261 // Field "shape" conversions
262
263 //- Interpolate face to point
264 template<class Type>
267 //- Interpolate point to face values
268 template<class Type>
270
271
272 // Custom Field Functions
273
274 //- The area-weighted average of a field
275 template<class Type>
276 Type areaAverage(const Field<Type>& fld) const
277 {
278 return gWeightedAverage(patch_.magSf(), fld);
279 }
280
281 //- The area-weighted sum of a field
282 template<class Type>
283 Type areaSum(const Field<Type>& fld) const
284 {
285 return gWeightedSum(patch_.magSf(), fld);
286 }
287
288 //- The face area magnitudes [magSf] - (swak = area)
290
291 //- The face centres - (swak = pos)
293
294 //- The face areas with their vector direction [Sf] - (swak = face)
296
297 //- The face unit normal direction [nf] - (expression: normal)
299
300 //- The patch point locations - (swak = pts)
302
303
304 //- Cell selection (set)
305 inline tmp<boolField> field_cellSet(const word& name) const;
306
307 //- Cell selection (zone)
308 inline tmp<boolField> field_cellZone(const word& name) const;
309
310 //- Face selection (set)
311 inline tmp<boolField> field_faceSet(const word& name) const;
312
313 //- Face selection (zone)
314 inline tmp<boolField> field_faceZone(const word& name) const;
315
316
317 //- A uniform random field
318 tmp<scalarField> field_rand(label seed=0, bool gaussian=false) const;
320 //- A Gaussian random field
321 tmp<scalarField> field_randGaussian(label seed=0) const
322 {
323 return field_rand(seed, true);
324 }
325};
326
327
328// Template specializations
329
330//- Retrieve field (surface field: bool)
331template<>
333
334//- Retrieve field (point field: bool)
335template<>
337
338
339// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341} // End namespace patchExpr
342} // End namespace expressions
343} // End namespace Foam
344
345
346// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347
348#include "patchExprDriverI.H"
349
350#ifdef NoRepository
352#endif
353
354// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355
356#endif
357
358// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
label nPoints() const
Number of points supporting patch faces.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
const exprResult & result() const noexcept
Const access to expression result.
Definition exprDriver.H:463
const dictionary & dict() const noexcept
The dictionary with all input data/specification.
Definition exprDriver.H:455
void setResult(Field< Type > *, bool wantPointData=false)
Set result field, taking ownership of the pointer.
Base driver for parsing value expressions associated with an fvMesh.
virtual label size() const
The natural field size for the expression.
tmp< Field< Type > > getVariableIfAvailable(const word &fldName) const
Retrieve variable as field if possible.
tmp< boolField > field_faceZone(const word &name) const
Face selection (zone).
tmp< Field< Type > > patchNormalField(const word &fldName)
Return surface normal field (snGrad).
tmp< Field< Type > > faceToPoint(const Field< Type > &field) const
Interpolate face to point.
tmp< Field< Type > > patchNeighbourField(const word &fldName)
Return patchField on the opposite patch of a coupled patch.
virtual label pointSize() const
The point field size for the expression.
tmp< boolField > field_cellZone(const word &name) const
Cell selection (zone).
tmp< vectorField > field_faceCentre() const
The face centres - (swak = pos).
tmp< Field< Type > > getVolField(const word &fldName)
Retrieve field (vol field).
tmp< Field< Type > > getPointField(const word &fldName)
Retrieve field (point field).
virtual unsigned parse(const std::string &expr, size_t pos=0, size_t len=std::string::npos)
Execute the parser.
void setResult(Field< Type > *ptr, bool pointVal=false)
Set result.
void operator=(const parseDriver &)=delete
tmp< boolField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
Cell selections (as logical).
tmp< Field< Type > > pointToFace(const Field< Type > &field) const
Interpolate point to face values.
tmp< vectorField > field_unitNormal() const
The face unit normal direction [nf] - (expression: normal).
parseDriver(const parseDriver &)=delete
tmp< Field< Type > > getSurfaceField(const word &fldName)
Retrieve field (surface field).
tmp< boolField > field_faceSet(const word &name) const
Face selection (set).
tmp< scalarField > field_rand(label seed=0, bool gaussian=false) const
A uniform random field.
tmp< boolField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
Face selections (as logical).
tmp< Field< Type > > patchInternalField(const word &fldName)
Return internal field next to patch.
const fvPatch & patch_
The referenced patch.
tmp< boolField > field_cellSet(const word &name) const
Cell selection (set).
Type areaAverage(const Field< Type > &fld) const
The area-weighted average of a field.
tmp< scalarField > field_faceArea() const
The face area magnitudes [magSf] - (swak = area).
tmp< vectorField > field_areaNormal() const
The face areas with their vector direction [Sf] - (swak = face).
tmp< vectorField > field_pointField() const
The patch point locations - (swak = pts).
virtual const fvMesh & mesh() const
The mesh we are attached to.
tmp< Field< Type > > getField(const word &fldName)
Return named field.
Type areaSum(const Field< Type > &fld) const
The area-weighted sum of a field.
tmp< scalarField > field_randGaussian(label seed=0) const
A Gaussian random field.
virtual ~parseDriver()=default
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
Generic interface code for Ragel/Lemon combination Subclasses should implement one or more process() ...
A class for managing temporary objects.
Definition tmp.H:75
sourceType
Enumeration defining the types of sources.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
volScalarField & p
rDeltaTY field()
Namespace for patch expressions parsing and evaluation.
A namespace for expression-related classes/traits etc.
Namespace of functions to calculate implicit derivatives returning a matrix.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Foam::surfaceFields.