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 | Copyright (C) 2011-2016 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 Group
28  grpFieldFunctionObjects
29 
30 Description
31  This function object samples near-patch volume fields
32 
33  Fields are stored
34  - every time step the field is updated with new values
35  - at output it writes the fields
36 
37  This functionObject can either be used
38  - to calculate a new field as a post-processing step or
39  - since the fields are registered, used in another functionObject
40 
41  Example of function object specification:
42  \verbatim
43  nearWallFields1
44  {
45  type nearWallFields;
46  libs ("libfieldFunctionObjects.so");
47  ...
48  fields ((p pNear)(U UNear));
49  patches (movingWall);
50  distance 0.13;
51  }
52  \endverbatim
53 
54 Usage
55  \table
56  Property | Description | Required | Default value
57  type | type name: nearWallFields | yes |
58  fields | list of fields with correspoding output field names | yes |
59  patches | list of patches to sample | yes |
60  distance | distance from patch to sample | yes |
61  \endtable
62 
63 See also
64  Foam::functionObjects::fvMeshFunctionObject
65 
66 SourceFiles
67  nearWallFields.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef functionObjects_nearWallFields_H
72 #define functionObjects_nearWallFields_H
73 
74 #include "fvMeshFunctionObject.H"
75 #include "volFields.H"
76 #include "Tuple2.H"
77 #include "interpolationCellPoint.H"
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 namespace Foam
82 {
83 namespace functionObjects
84 {
85 
86 /*---------------------------------------------------------------------------*\
87  Class nearWallFields Declaration
88 \*---------------------------------------------------------------------------*/
89 
90 class nearWallFields
91 :
92  public fvMeshFunctionObject
93 {
94 protected:
95 
96  // Protected member data
97 
98  // Read from dictionary
99 
100  //- Fields to process
101  List<Tuple2<word, word>> fieldSet_;
102 
103  //- Patches to sample
105 
106  //- Distance away from wall
107  scalar distance_;
108 
109  //- From original field to sampled result
110  HashTable<word> fieldMap_;
111 
112  //- From resulting back to original field
113  HashTable<word> reverseFieldMap_;
115 
116  // Calculated addressing
117 
118  //- From cell to seed patch faces
120 
121  //- From cell to tracked end point
123 
124  //- Map from cell based data back to patch based data
126 
127 
128  // Locally constructed fields
129 
135 
136 
137  // Protected Member Functions
138 
139  //- Calculate addressing from cells back to patch faces
140  void calcAddressing();
141 
142  template<class Type>
144  (
146  ) const;
147 
148  //- Override boundary fields with sampled values
149  template<class Type>
151  (
152  const interpolationCellPoint<Type>& interpolator,
154  ) const;
156  template<class Type>
158  (
160  ) const;
161 
162 
163 private:
164 
165  //- Disallow default bitwise copy construct
167 
168  //- Disallow default bitwise assignment
169  void operator=(const nearWallFields&);
170 
171 public:
172 
173  //- Runtime type information
174  TypeName("nearWallFields");
175 
176 
177  // Constructors
178 
179  //- Construct for given objectRegistry and dictionary.
180  // Allow the possibility to load fields from files
182  (
183  const word& name,
184  const Time& runTime,
185  const dictionary& dict
186  );
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 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace functionObjects
209 } // End namespace Foam
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 #ifdef NoRepository
214  #include "nearWallFieldsTemplates.C"
215 #endif
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #endif
220 
221 // ************************************************************************* //
dictionary dict
virtual bool write()
Write the near-wall fields.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
HashTable< word > reverseFieldMap_
From resulting back to original field.
void calcAddressing()
Calculate addressing from cells back to patch faces.
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.
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
virtual bool execute()
Calculate the near-wall fields.
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
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))
const word & name() const
Return the name of this functionObject.
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.
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
Override boundary fields with sampled values.
This function object 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:62
PtrList< volSymmTensorField > vSymmtf_
TypeName("nearWallFields")
Runtime type information.
virtual bool read(const dictionary &)
Read the controls.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
scalar distance_
Distance away from wall.
PtrList< volSphericalTensorField > vSpheretf_
Namespace for OpenFOAM.