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  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  writeControl writeTime;
49 
50  fields
51  (
52  (p pNear)
53  (U UNear)
54  );
55 
56  patches (movingWall);
57 
58  distance 0.13;
59  }
60  \endverbatim
61 
62 Usage
63  \table
64  Property | Description | Required | Default value
65  type | type name: nearWallFields | yes |
66  fields | list of fields with correspoding output field names | yes |
67  patches | list of patches to sample | yes |
68  distance | distance from patch to sample | yes |
69  \endtable
70 
71 See also
72  Foam::functionObjects::fvMeshFunctionObject
73 
74 SourceFiles
75  nearWallFields.C
76 
77 \*---------------------------------------------------------------------------*/
78 
79 #ifndef functionObjects_nearWallFields_H
80 #define functionObjects_nearWallFields_H
81 
82 #include "fvMeshFunctionObject.H"
83 #include "volFields.H"
84 #include "Tuple2.H"
85 #include "interpolationCellPoint.H"
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 namespace Foam
90 {
91 namespace functionObjects
92 {
93 
94 /*---------------------------------------------------------------------------*\
95  Class nearWallFields Declaration
96 \*---------------------------------------------------------------------------*/
97 
98 class nearWallFields
99 :
100  public fvMeshFunctionObject
101 {
102 protected:
103 
104  // Protected member data
105 
106  // Read from dictionary
107 
108  //- Fields to process
109  List<Tuple2<word, word>> fieldSet_;
110 
111  //- Patches to sample
113 
114  //- Distance away from wall
115  scalar distance_;
116 
117  //- From original field to sampled result
118  HashTable<word> fieldMap_;
119 
120  //- From resulting back to original field
121  HashTable<word> reverseFieldMap_;
123 
124  // Calculated addressing
125 
126  //- From cell to seed patch faces
128 
129  //- From cell to tracked end point
131 
132  //- Map from cell based data back to patch based data
134 
135 
136  // Locally constructed fields
137 
143 
144 
145  // Protected Member Functions
146 
147  //- Calculate addressing from cells back to patch faces
148  void calcAddressing();
149 
150  template<class Type>
152  (
154  ) const;
155 
156  //- Override boundary fields with sampled values
157  template<class Type>
159  (
160  const interpolationCellPoint<Type>& interpolator,
162  ) const;
164  template<class Type>
166  (
168  ) const;
169 
170 
171 private:
172 
173  //- Disallow default bitwise copy construct
175 
176  //- Disallow default bitwise assignment
177  void operator=(const nearWallFields&);
178 
179 public:
180 
181  //- Runtime type information
182  TypeName("nearWallFields");
183 
184 
185  // Constructors
186 
187  //- Construct for given objectRegistry and dictionary.
188  // Allow the possibility to load fields from files
190  (
191  const word& name,
192  const Time& runTime,
193  const dictionary& dict
194  );
195 
196 
197  //- Destructor
198  virtual ~nearWallFields();
199 
200 
201  // Member Functions
202 
203  //- Read the controls
204  virtual bool read(const dictionary&);
205 
206  //- Calculate the near-wall fields
207  virtual bool execute();
208 
209  //- Write the near-wall fields
210  virtual bool write();
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace functionObjects
217 } // End namespace Foam
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #ifdef NoRepository
222  #include "nearWallFieldsTemplates.C"
223 #endif
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 #endif
228 
229 // ************************************************************************* //
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:137
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.
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:210
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.
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:63
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.