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-2019 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  patches (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  patches | list of patches to sample | yes |
65  distance | distance from patch to sample | yes |
66  \endtable
67 
68 See also
69  Foam::functionObjects::fvMeshFunctionObject
70 
71 SourceFiles
72  nearWallFields.C
73 
74 \*---------------------------------------------------------------------------*/
75 
76 #ifndef functionObjects_nearWallFields_H
77 #define functionObjects_nearWallFields_H
78 
79 #include "fvMeshFunctionObject.H"
80 #include "volFields.H"
81 #include "Tuple2.H"
82 #include "interpolationCellPoint.H"
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 namespace Foam
87 {
88 namespace functionObjects
89 {
90 
91 /*---------------------------------------------------------------------------*\
92  Class nearWallFields Declaration
93 \*---------------------------------------------------------------------------*/
94 
95 class nearWallFields
96 :
97  public fvMeshFunctionObject
98 {
99 protected:
100 
101  // Protected member data
102 
103  // Read from dictionary
104 
105  //- Fields to process
106  List<Tuple2<word, word>> fieldSet_;
107 
108  //- Patches to sample
110 
111  //- Distance away from wall
112  scalar distance_;
113 
114  //- From original field to sampled result
115  HashTable<word> fieldMap_;
116 
117  //- From resulting back to original field
118  HashTable<word> reverseFieldMap_;
120 
121  // Calculated addressing
122 
123  //- From cell to seed patch faces
125 
126  //- From cell to tracked end point
128 
129  //- Map from cell based data back to patch based data
131 
132 
133  // Locally constructed fields
134 
140 
141 
142  // Protected Member Functions
143 
144  //- Calculate addressing from cells back to patch faces
145  void calcAddressing();
146 
147  template<class Type>
149  (
151  ) const;
152 
153  //- Override boundary fields with sampled values
154  template<class Type>
156  (
157  const interpolationCellPoint<Type>& interpolator,
159  ) const;
161  template<class Type>
163  (
165  ) const;
166 
167 
168 public:
169 
170  //- Runtime type information
171  TypeName("nearWallFields");
172 
173 
174  // Constructors
175 
176  //- Construct for given objectRegistry and dictionary.
177  // Allow the possibility to load fields from files
179  (
180  const word& name,
181  const Time& runTime,
182  const dictionary& dict
183  );
184 
185  //- Disallow default bitwise copy construction
186  nearWallFields(const nearWallFields&) = delete;
187 
188 
189  //- Destructor
190  virtual ~nearWallFields();
191 
192 
193  // Member Functions
194 
195  //- Read the controls
196  virtual bool read(const dictionary&);
197 
198  //- Calculate the near-wall fields
199  virtual bool execute();
200 
201  //- Write the near-wall fields
202  virtual bool write();
203 
204 
205  // Member Operators
206 
207  //- Disallow default bitwise assignment
208  void operator=(const nearWallFields&) = delete;
209 };
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace functionObjects
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #ifdef NoRepository
220  #include "nearWallFieldsTemplates.C"
221 #endif
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
dictionary dict
virtual bool write()
Write the near-wall fields.
const word & name() const
Return the name of this functionObject.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
HashTable< word > reverseFieldMap_
From resulting back to original field.
void calcAddressing()
Calculate addressing from cells back to patch faces.
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
Override boundary fields with sampled values.
engineTime & runTime
Generic GeometricField class.
List< List< point > > cellToSamples_
From cell to tracked end point.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
List< Tuple2< word, word > > fieldSet_
Fields to process.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
virtual bool execute()
Calculate the near-wall fields.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
PtrList< volScalarField > vsf_
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
A class for handling words, derived from string.
Definition: word.H:59
HashTable< word > fieldMap_
From original field to sampled result.
PtrList< volTensorField > vtf_
labelListList cellToWalls_
From cell to seed patch faces.
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
void operator=(const nearWallFields &)=delete
Disallow default bitwise assignment.
Samples near-patch volume fields.
PtrList< volVectorField > vvf_
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
labelHashSet patchSet_
Patches to sample.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
PtrList< volSymmTensorField > vSymmtf_
TypeName("nearWallFields")
Runtime type information.
virtual bool read(const dictionary &)
Read the controls.
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
scalar distance_
Distance away from wall.
PtrList< volSphericalTensorField > vSpheretf_
Namespace for OpenFOAM.