nearWallFields.H
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 Class
25  Foam::functionObjects::nearWallFields
26 
27 Description
28  Samples near-patch volume fields.
29 
30  Fields are stored
31  - every time step the field is updated with new values
32  - at output it writes the fields
33 
34  This functionObject can either be used
35  - to calculate a new field as a post-processing step or
36  - since the fields are registered, used in another functionObject
37 
38  Example of function object specification:
39  \verbatim
40  nearWallFields1
41  {
42  type nearWallFields;
43  libs ("libfieldFunctionObjects.so");
44 
45  writeControl writeTime;
46 
47  fields
48  (
49  (p pNear)
50  (U UNear)
51  );
52 
53  patch movingWall;
54 
55  distance 0.13;
56  }
57  \endverbatim
58 
59 Usage
60  \table
61  Property | Description | Required | Default value
62  type | type name: nearWallFields | yes |
63  fields | list of fields with corresponding output field names | yes |
64  patch | patch to sample | no |
65  patches | list of patches to sample | no |
66  distance | distance from patch to sample | yes |
67  \endtable
68 
69 See also
70  Foam::functionObjects::fvMeshFunctionObject
71 
72 SourceFiles
73  nearWallFields.C
74 
75 \*---------------------------------------------------------------------------*/
76 
77 #ifndef functionObjects_nearWallFields_H
78 #define functionObjects_nearWallFields_H
79 
80 #include "fvMeshFunctionObject.H"
81 #include "volFields.H"
82 #include "Tuple2.H"
83 #include "interpolationCellPoint.H"
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 namespace Foam
88 {
89 namespace functionObjects
90 {
91 
92 /*---------------------------------------------------------------------------*\
93  Class nearWallFields Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 class nearWallFields
97 :
98  public fvMeshFunctionObject
99 {
100 protected:
101 
102  // Protected member data
103 
104  // Read from dictionary
105 
106  //- Fields to process
107  List<Tuple2<word, word>> fieldSet_;
108 
109  //- Patches to sample
111 
112  //- Distance away from wall
113  scalar distance_;
114 
115  //- From original field to sampled result
116  HashTable<word> fieldMap_;
117 
118  //- From resulting back to original field
119  HashTable<word> reverseFieldMap_;
120 
121 
122  // Calculated addressing
123 
124  //- From cell to seed patch faces
126 
127  //- From cell to tracked end point
129 
130  //- Map from cell based data back to patch based data
132 
133 
134  // Locally constructed fields
135 
141 
142 
143  // Protected Member Functions
144 
145  //- Calculate addressing from cells back to patch faces
146  void calcAddressing();
147 
148  template<class Type>
149  void createFields
150  (
152  ) const;
153 
154  //- Override boundary fields with sampled values
155  template<class Type>
157  (
158  const interpolationCellPoint<Type>& interpolator,
160  ) const;
161 
162  template<class Type>
163  void sampleFields
164  (
166  ) const;
167 
168 
169 public:
170 
171  //- Runtime type information
172  TypeName("nearWallFields");
173 
174 
175  // Constructors
176 
177  //- Construct for given objectRegistry and dictionary.
178  // Allow the possibility to load fields from files
180  (
181  const word& name,
182  const Time& runTime,
183  const dictionary& dict
184  );
185 
186  //- Disallow default bitwise copy construction
187  nearWallFields(const nearWallFields&) = delete;
188 
189 
190  //- Destructor
191  virtual ~nearWallFields();
192 
193 
194  // Member Functions
195 
196  //- Read the controls
197  virtual bool read(const dictionary&);
198 
199  //- Return the list of fields required
200  virtual wordList fields() const;
201 
202  //- Calculate the near-wall fields
203  virtual bool execute();
204 
205  //- Write the near-wall fields
206  virtual bool write();
207 
208 
209  // Member Operators
210 
211  //- Disallow default bitwise assignment
212  void operator=(const nearWallFields&) = delete;
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace functionObjects
219 } // End namespace Foam
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #ifdef NoRepository
224  #include "nearWallFieldsTemplates.C"
225 #endif
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
Generic GeometricField class.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return the name of this functionObject.
Samples near-patch volume fields.
void sampleFields(PtrList< VolField< Type >> &) const
scalar distance_
Distance away from wall.
HashTable< word > reverseFieldMap_
From resulting back to original field.
void createFields(PtrList< VolField< Type >> &) const
virtual wordList fields() const
Return the list of fields required.
PtrList< volScalarField > vsf_
List< Tuple2< word, word > > fieldSet_
Fields to process.
List< List< point > > cellToSamples_
From cell to tracked end point.
labelListList cellToWalls_
From cell to seed patch faces.
void calcAddressing()
Calculate addressing from cells back to patch faces.
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, VolField< Type > &fld) const
Override boundary fields with sampled values.
PtrList< volTensorField > vtf_
TypeName("nearWallFields")
Runtime type information.
labelHashSet patchSet_
Patches to sample.
PtrList< volSymmTensorField > vSymmtf_
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
PtrList< volVectorField > vvf_
autoPtr< distributionMap > getPatchDataMapPtr_
Map from cell based data back to patch based data.
virtual bool execute()
Calculate the near-wall fields.
PtrList< volSphericalTensorField > vSpheretf_
void operator=(const nearWallFields &)=delete
Disallow default bitwise assignment.
virtual bool write()
Write the near-wall fields.
virtual bool read(const dictionary &)
Read the controls.
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 string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
dictionary dict