probes.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::probes
26 
27 Group
28  grpUtilitiesFunctionObjects
29 
30 Description
31  Set of locations to sample.
32 
33  Call write() to sample and write files.
34 
35 SourceFiles
36  probes.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef probes_H
41 #define probes_H
42 
43 #include "functionObject.H"
44 #include "HashPtrTable.H"
45 #include "OFstream.H"
46 #include "polyMesh.H"
47 #include "pointField.H"
48 #include "volFieldsFwd.H"
49 #include "surfaceFieldsFwd.H"
50 #include "surfaceMesh.H"
51 #include "wordReList.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 // Forward declaration of classes
59 class Time;
60 class objectRegistry;
61 class dictionary;
62 class fvMesh;
63 class mapPolyMesh;
64 
65 /*---------------------------------------------------------------------------*\
66  Class probes Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 class probes
70 :
71  public functionObject,
72  public pointField
73 {
74 protected:
75 
76  // Protected classes
77 
78  //- Class used for grouping field types
79  template<class Type>
80  class fieldGroup
81  :
82  public DynamicList<word>
83  {
84  public:
85  //- Construct null
86  fieldGroup()
87  :
88  DynamicList<word>(0)
89  {}
90  };
91 
92 
93  // Protected member data
94 
95  //- Const reference to fvMesh
96  const fvMesh& mesh_;
97 
98  //- Load fields from files (not from objectRegistry)
99  bool loadFromFiles_;
100 
101 
102  // Read from dictonary
103 
104  //- Names of fields to probe
106 
107  //- Fixed locations, default = yes
108  // Note: set to false for moving mesh calations where locations
109  // should move with the mesh
110  bool fixedLocations_;
111 
112  //- Interpolation scheme name
113  // Note: only possible when fixedLocations_ is true
115 
116 
117  // Calculated
118 
119  //- Categorized scalar/vector/tensor vol fields
125 
126  //- Categorized scalar/vector/tensor surf fields
132 
133  // Cells to be probed (obtained from the locations)
135 
136  // Faces to be probed
138 
139  //- Current open files
141 
142 
143  // Protected Member Functions
144 
145  //- Clear old field groups
146  void clearFieldGroups();
147 
148  //- Append fieldName to the appropriate group
149  label appendFieldGroup(const word& fieldName, const word& fieldType);
150 
151  //- Classify field types, returns the number of fields
153 
154  //- Find cells and faces containing probes
155  virtual void findElements(const fvMesh&);
156 
157  //- Classify field type and Open/close file streams,
158  // returns number of fields to sample
159  label prepare();
160 
161 
162 private:
163 
164  //- Sample and write a particular volume field
165  template<class Type>
166  void sampleAndWrite
167  (
169  );
170 
171 
172  //- Sample and write a particular surface field
173  template<class Type>
174  void sampleAndWrite
175  (
177  );
178 
179  //- Sample and write all the fields of the given type
180  template<class Type>
181  void sampleAndWrite(const fieldGroup<Type>&);
182 
183  //- Sample and write all the surface fields of the given type
184  template<class Type>
185  void sampleAndWriteSurfaceFields(const fieldGroup<Type>&);
186 
187  //- Disallow default bitwise copy construct
188  probes(const probes&);
189 
190  //- Disallow default bitwise assignment
191  void operator=(const probes&);
192 
193 
194 public:
195 
196  //- Runtime type information
197  TypeName("probes");
198 
199 
200  // Constructors
201 
202  //- Construct from Time and dictionary
203  probes
204  (
205  const word& name,
206  const Time& time,
207  const dictionary& dict
208  );
209 
210  //- Construct for given objectRegistry and dictionary.
211  // Allow the possibility to load fields from files
212  probes
213  (
214  const word& name,
215  const objectRegistry& obr,
216  const dictionary& dict,
217  const bool loadFromFiles = false
218  );
219 
220 
221  //- Destructor
222  virtual ~probes();
223 
224 
225  // Member Functions
226 
227  //- Return names of fields to probe
228  virtual const wordReList& fieldNames() const
229  {
230  return fieldSelection_;
231  }
232 
233  //- Return locations to probe
234  virtual const pointField& probeLocations() const
235  {
236  return *this;
237  }
238 
239  //- Return location for probe i
240  virtual const point& probe(const label i) const
241  {
242  return operator[](i);
243  }
244 
245  //- Cells to be probed (obtained from the locations)
246  const labelList& elements() const
247  {
248  return elementList_;
249  }
250 
251  //- Read the probes
252  virtual bool read(const dictionary&);
253 
254  //- Execute, currently does nothing
255  virtual bool execute();
256 
257  //- Sample and write
258  virtual bool write();
259 
260  //- Update for changes of mesh
261  virtual void updateMesh(const mapPolyMesh&);
262 
263  //- Update for changes of mesh
264  virtual void movePoints(const polyMesh&);
265 
266  //- Update for changes of mesh due to readUpdate
267  virtual void readUpdate(const polyMesh::readUpdateState state)
268  {}
269 
270  //- Sample a volume field at all locations
271  template<class Type>
273  (
275  ) const;
276 
277  //- Sample a single vol field on all sample locations
278  template<class Type>
279  tmp<Field<Type>> sample(const word& fieldName) const;
280 
281  //- Sample a single scalar field on all sample locations
282  template<class Type>
283  tmp<Field<Type>> sampleSurfaceFields(const word& fieldName) const;
284 
285  //- Sample a surface field at all locations
286  template<class Type>
288  (
290  ) const;
291 };
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace Foam
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 #ifdef NoRepository
301  #include "probesTemplates.C"
302 #endif
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 #endif
307 
308 // ************************************************************************* //
bool fixedLocations_
Fixed locations, default = yes.
Definition: probes.H:109
dictionary dict
fieldGroup< vector > surfaceVectorFields_
Definition: probes.H:127
void clearFieldGroups()
Clear old field groups.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual bool write()
Sample and write.
Definition: probes.C:368
fieldGroup< scalar > surfaceScalarFields_
Categorized scalar/vector/tensor surf fields.
Definition: probes.H:126
virtual const point & probe(const label i) const
Return location for probe i.
Definition: probes.H:239
fieldGroup< sphericalTensor > surfaceSphericalTensorFields_
Definition: probes.H:128
const word & name() const
Return the name of this functionObject.
label prepare()
Classify field type and Open/close file streams,.
Definition: probes.C:168
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: probes.C:389
virtual bool execute()
Execute, currently does nothing.
Definition: probes.C:362
TypeName("probes")
Runtime type information.
Set of locations to sample.
Definition: probes.H:68
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: probes.H:266
fieldGroup< symmTensor > surfaceSymmTensorFields_
Definition: probes.H:129
Generic GeometricField class.
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
Abstract base-class for Time/database function objects.
tmp< Field< Type > > sampleSurfaceFields(const word &fieldName) const
Sample a single scalar field on all sample locations.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
fieldGroup< tensor > surfaceTensorFields_
Definition: probes.H:130
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
labelList faceList_
Definition: probes.H:136
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
virtual void findElements(const fvMesh &)
Find cells and faces containing probes.
Definition: probes.C:51
fieldGroup()
Construct null.
Definition: probes.H:85
fieldGroup< scalar > scalarFields_
Categorized scalar/vector/tensor vol fields.
Definition: probes.H:119
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:95
const labelList & elements() const
Cells to be probed (obtained from the locations)
Definition: probes.H:245
label appendFieldGroup(const word &fieldName, const word &fieldType)
Append fieldName to the appropriate group.
fieldGroup< tensor > tensorFields_
Definition: probes.H:123
fieldGroup< symmTensor > symmTensorFields_
Definition: probes.H:122
wordReList fieldSelection_
Names of fields to probe.
Definition: probes.H:104
HashPtrTable< OFstream > probeFilePtrs_
Current open files.
Definition: probes.H:139
virtual void movePoints(const polyMesh &)
Update for changes of mesh.
Definition: probes.C:468
fieldGroup< vector > vectorFields_
Definition: probes.H:120
bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: probes.H:98
label classifyFields()
Classify field types, returns the number of fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
labelList elementList_
Definition: probes.H:133
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:336
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
word interpolationScheme_
Interpolation scheme name.
Definition: probes.H:113
virtual ~probes()
Destructor.
Definition: probes.C:330
fieldGroup< sphericalTensor > sphericalTensorFields_
Definition: probes.H:121
virtual const pointField & probeLocations() const
Return locations to probe.
Definition: probes.H:233
Class used for grouping field types.
Definition: probes.H:79
Namespace for OpenFOAM.
virtual const wordReList & fieldNames() const
Return names of fields to probe.
Definition: probes.H:227
void operator=(const word &)
Assignment of all addressed entries to the given value.
Definition: DynamicListI.H:384