nearWallFields.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::functionObjects::nearWallFields
29
30Group
31 grpFieldFunctionObjects
32
33Description
34 Samples near-patch volume fields within an input distance range.
35
36 Operands:
37 \table
38 Operand | Type | Location
39 input | vol<Type>Field | $FOAM_CASE/<time>/<inpField>
40 output file | - | -
41 output field | vol<Type>Field | $FOAM_CASE/<time>/<outField>
42 \endtable
43
44 where \c <Type>=Scalar/Vector/SphericalTensor/SymmTensor/Tensor.
45
46Usage
47 Minimal example by using \c system/controlDict.functions:
48 \verbatim
49 nearWallFields1
50 {
51 // Mandatory entries (unmodifiable)
52 type nearWallFields;
53 libs (fieldFunctionObjects);
54 fields
55 (
56 (<field1> <outField1>)
57 (<field2> <outField2>)
58 );
59 patches (<patch1> <patch2> ... <patchN>);
60 distance 0.01;
61
62 // Optional (inherited) entries
63 ...
64 }
65 \endverbatim
66
67 where the entries mean:
68 \table
69 Property | Description | Type | Req'd | Dflt
70 type | Type name: nearWallFields | word | yes | -
71 libs | Library name: fieldFunctionObjects | word | yes | -
72 fields | Names of input-output fields | wordHashTable | yes | -
73 patches | Names of patches to sample | wordList | yes | -
74 distance | Wall-normal distance from patch to sample | scalar | yes | -
75 \endtable
76
77 The inherited entries are elaborated in:
78 - \link functionObject.H \endlink
79
80 Usage by the \c postProcess utility is not available.
81
82See also
83 - Foam::functionObject
84 - Foam::functionObjects::fvMeshFunctionObject
85 - ExtendedCodeGuide::functionObjects::field::nearWallFields
86
87SourceFiles
88 nearWallFields.C
89 nearWallFieldsTemplates.C
90
91\*---------------------------------------------------------------------------*/
92
93#ifndef functionObjects_nearWallFields_H
94#define functionObjects_nearWallFields_H
95
97#include "volFields.H"
98#include "Tuple2.H"
100
101// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
102
103namespace Foam
104{
105namespace functionObjects
106{
107
108/*---------------------------------------------------------------------------*\
109 Class nearWallFields Declaration
110\*---------------------------------------------------------------------------*/
111
112class nearWallFields
113:
114 public fvMeshFunctionObject
115{
116protected:
117
118 // Protected Data
119
120 // Read from dictionary
121
122 //- Fields to process (input-name output-name)
123 List<Tuple2<word, word>> fieldSet_;
124
125 //- Patches to sample
127
128 //- Distance away from wall
129 scalar distance_;
130
131 //- From original field to sampled result
132 HashTable<word> fieldMap_;
133
134 //- From resulting back to original field
135 HashTable<word> reverseFieldMap_;
136
137
138 // Calculated addressing
139
140 //- From cell to seed patch faces
142
143 //- From cell to tracked end point
144 List<List<point>> cellToSamples_;
145
146 //- Map from cell based data back to patch based data
147 autoPtr<mapDistribute> getPatchDataMapPtr_;
148
149
150 // Locally constructed fields
151
152 PtrList<volScalarField> vsf_;
153 PtrList<volVectorField> vvf_;
154 PtrList<volSphericalTensorField> vSpheretf_;
155 PtrList<volSymmTensorField> vSymmtf_;
156 PtrList<volTensorField> vtf_;
157
158
159 // Protected Member Functions
160
161 //- Calculate addressing from cells back to patch faces
162 void calcAddressing();
164 template<class Type>
165 void createFields
166 (
168 ) const;
169
170 //- Override boundary fields with sampled values
171 template<class Type>
173 (
174 const interpolationCellPoint<Type>& interpolator,
176 ) const;
178 template<class Type>
179 void sampleFields
182 ) const;
184
185public:
187 //- Runtime type information
188 TypeName("nearWallFields");
189
190
191 // Constructors
193 //- Construct for given objectRegistry and dictionary.
194 // Allow the possibility to load fields from files
196 (
197 const word& name,
198 const Time& runTime,
199 const dictionary& dict
200 );
201
202 //- No copy construct
205 //- No copy assignment
206 void operator=(const nearWallFields&) = delete;
208
209 //- Destructor
210 virtual ~nearWallFields() = default;
211
212
213 // Member Functions
214
215 //- Read the controls
216 virtual bool read(const dictionary& dict);
217
218 //- Calculate the near-wall fields
219 virtual bool execute();
220
221 //- Write the near-wall fields
222 virtual bool write();
223};
224
225
226// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227
228} // End namespace functionObjects
229} // End namespace Foam
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233#ifdef NoRepository
235#endif
236
237// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238
239#endif
240
241// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic GeometricField class.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const word & name() const noexcept
Return the name of this functionObject.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Samples near-patch volume fields within an input distance range.
scalar distance_
Distance away from wall.
HashTable< word > reverseFieldMap_
From resulting back to original field.
virtual ~nearWallFields()=default
Destructor.
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
PtrList< volScalarField > vsf_
List< Tuple2< word, word > > fieldSet_
Fields to process (input-name output-name)
List< List< point > > cellToSamples_
From cell to tracked end point.
virtual bool read(const dictionary &dict)
Read the controls.
labelListList cellToWalls_
From cell to seed patch faces.
void calcAddressing()
Calculate addressing from cells back to patch faces.
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
Override boundary fields with sampled values.
PtrList< volTensorField > vtf_
TypeName("nearWallFields")
Runtime type information.
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
labelHashSet patchSet_
Patches to sample.
PtrList< volSymmTensorField > vSymmtf_
nearWallFields(const nearWallFields &)=delete
No copy construct.
PtrList< volVectorField > vvf_
virtual bool execute()
Calculate the near-wall fields.
PtrList< volSphericalTensorField > vSpheretf_
void operator=(const nearWallFields &)=delete
No copy assignment.
virtual bool write()
Write the near-wall fields.
HashTable< word > fieldMap_
From original field to sampled result.
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
A class for handling words, derived from Foam::string.
Definition: word.H:68
engineTime & runTime
Namespace for OpenFOAM.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73