probes.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::probes
26 
27 Description
28  Set of locations to sample.
29 
30  Call write() to sample and write files.
31 
32 SourceFiles
33  probes.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef probes_H
38 #define probes_H
39 
40 #include "functionObject.H"
41 #include "HashPtrTable.H"
42 #include "OFstream.H"
43 #include "polyMesh.H"
44 #include "pointField.H"
45 #include "volFieldsFwd.H"
46 #include "surfaceFieldsFwd.H"
47 #include "surfaceMesh.H"
48 #include "wordReList.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward declaration of classes
56 class Time;
57 class objectRegistry;
58 class dictionary;
59 class fvMesh;
60 class mapPolyMesh;
61 
62 /*---------------------------------------------------------------------------*\
63  Class probes Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class probes
67 :
68  public functionObject,
69  public pointField
70 {
71 protected:
72 
73  // Protected classes
74 
75  //- Class used for grouping field types
76  template<class Type>
77  class fieldGroup
78  :
79  public DynamicList<word>
80  {
81  public:
82  //- Construct null
83  fieldGroup()
84  :
85  DynamicList<word>(0)
86  {}
87  };
88 
89 
90  // Protected member data
91 
92  //- Const reference to fvMesh
93  const fvMesh& mesh_;
94 
95  //- Load fields from files (not from objectRegistry)
96  bool loadFromFiles_;
97 
98 
99  // Read from dictonary
100 
101  //- Names of fields to probe
103 
104  //- Fixed locations, default = yes
105  // Note: set to false for moving mesh calations where locations
106  // should move with the mesh
107  bool fixedLocations_;
108 
109  //- Interpolation scheme name
110  // Note: only possible when fixedLocations_ is true
112 
113 
114  // Calculated
115 
116  //- Categorized scalar/vector/tensor vol fields
122 
123  //- Categorized scalar/vector/tensor surf fields
129 
130  // Cells to be probed (obtained from the locations)
132 
133  // Faces to be probed
135 
136  //- Current open files
138 
139 
140  // Protected Member Functions
141 
142  //- Clear old field groups
143  void clearFieldGroups();
144 
145  //- Append fieldName to the appropriate group
146  label appendFieldGroup(const word& fieldName, const word& fieldType);
147 
148  //- Classify field types, returns the number of fields
150 
151  //- Find cells and faces containing probes
152  virtual void findElements(const fvMesh&);
153 
154  //- Classify field type and Open/close file streams,
155  // returns number of fields to sample
156  label prepare();
157 
158 
159 private:
160 
161  //- Sample and write a particular volume field
162  template<class Type>
163  void sampleAndWrite
164  (
166  );
167 
168 
169  //- Sample and write a particular surface field
170  template<class Type>
171  void sampleAndWrite
172  (
174  );
175 
176  //- Sample and write all the fields of the given type
177  template<class Type>
178  void sampleAndWrite(const fieldGroup<Type>&);
179 
180  //- Sample and write all the surface fields of the given type
181  template<class Type>
182  void sampleAndWriteSurfaceFields(const fieldGroup<Type>&);
183 
184 
185 public:
186 
187  //- Runtime type information
188  TypeName("probes");
189 
190 
191  // Constructors
192 
193  //- Construct from Time and dictionary
194  probes
195  (
196  const word& name,
197  const Time& time,
198  const dictionary& dict
199  );
200 
201  //- Construct for given objectRegistry and dictionary.
202  // Allow the possibility to load fields from files
203  probes
204  (
205  const word& name,
206  const objectRegistry& obr,
207  const dictionary& dict,
208  const bool loadFromFiles = false
209  );
210 
211  //- Disallow default bitwise copy construction
212  probes(const probes&) = delete;
213 
214 
215  //- Destructor
216  virtual ~probes();
217 
218 
219  // Member Functions
220 
221  //- Return names of fields to probe
222  virtual const wordReList& fieldNames() const
223  {
224  return fieldSelection_;
225  }
226 
227  //- Return locations to probe
228  virtual const pointField& probeLocations() const
229  {
230  return *this;
231  }
232 
233  //- Return location for probe i
234  virtual const point& probe(const label i) const
235  {
236  return operator[](i);
237  }
238 
239  //- Cells to be probed (obtained from the locations)
240  const labelList& elements() const
241  {
242  return elementList_;
243  }
244 
245  //- Read the probes
246  virtual bool read(const dictionary&);
247 
248  //- Execute, currently does nothing
249  virtual bool execute();
250 
251  //- Sample and write
252  virtual bool write();
253 
254  //- Update for changes of mesh
255  virtual void updateMesh(const mapPolyMesh&);
256 
257  //- Update for changes of mesh
258  virtual void movePoints(const polyMesh&);
259 
260  //- Update for changes of mesh due to readUpdate
261  virtual void readUpdate(const polyMesh::readUpdateState state)
262  {}
263 
264  //- Sample a volume field at all locations
265  template<class Type>
267  (
269  ) const;
270 
271  //- Sample a single vol field on all sample locations
272  template<class Type>
273  tmp<Field<Type>> sample(const word& fieldName) const;
274 
275  //- Sample a single scalar field on all sample locations
276  template<class Type>
277  tmp<Field<Type>> sampleSurfaceFields(const word& fieldName) const;
278 
279  //- Sample a surface field at all locations
280  template<class Type>
282  (
284  ) const;
285 
286 
287  // Member Operators
288 
289  //- Disallow default bitwise assignment
290  void operator=(const probes&) = delete;
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:106
dictionary dict
fieldGroup< vector > surfaceVectorFields_
Definition: probes.H:124
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:123
virtual const point & probe(const label i) const
Return location for probe i.
Definition: probes.H:233
fieldGroup< sphericalTensor > surfaceSphericalTensorFields_
Definition: probes.H:125
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:158
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:65
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: probes.H:260
fieldGroup< symmTensor > surfaceSymmTensorFields_
Definition: probes.H:126
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:127
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
labelList faceList_
Definition: probes.H:133
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:82
fieldGroup< scalar > scalarFields_
Categorized scalar/vector/tensor vol fields.
Definition: probes.H:116
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:92
const labelList & elements() const
Cells to be probed (obtained from the locations)
Definition: probes.H:239
label appendFieldGroup(const word &fieldName, const word &fieldType)
Append fieldName to the appropriate group.
fieldGroup< tensor > tensorFields_
Definition: probes.H:120
fieldGroup< symmTensor > symmTensorFields_
Definition: probes.H:119
wordReList fieldSelection_
Names of fields to probe.
Definition: probes.H:101
HashPtrTable< OFstream > probeFilePtrs_
Current open files.
Definition: probes.H:136
virtual void movePoints(const polyMesh &)
Update for changes of mesh.
Definition: probes.C:468
fieldGroup< vector > vectorFields_
Definition: probes.H:117
bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: probes.H:95
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:130
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:336
probes(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
Definition: probes.C:281
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:110
virtual ~probes()
Destructor.
Definition: probes.C:330
fieldGroup< sphericalTensor > sphericalTensorFields_
Definition: probes.H:118
virtual const pointField & probeLocations() const
Return locations to probe.
Definition: probes.H:227
Class used for grouping field types.
Definition: probes.H:76
Namespace for OpenFOAM.
virtual const wordReList & fieldNames() const
Return names of fields to probe.
Definition: probes.H:221
void operator=(const word &)
Assignment of all addressed entries to the given value.
Definition: DynamicListI.H:388