sampledSets.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::sampledSets
29 
30 Description
31  Set of sets to sample.
32 
33  The write() method is used to sample and write files.
34 
35  Example of function object specification:
36 
37  \verbatim
38  surfaces
39  {
40  type sets;
41  libs (sampling);
42 
43  // Write at same frequency as fields
44  writeControl writeTime;
45  writeInterval 1;
46 
47  // Fields to be sampled
48  fields (p U);
49 
50  // Scheme to obtain values
51  interpolationScheme cellPoint;
52 
53  // Output format
54  setFormat vtk;
55 
56  formatOptions
57  {
58  vtk
59  {
60  precision 10;
61  }
62  }
63 
64  sets
65  {
66  Uref
67  {
68  type cloud;
69  axis xyz;
70  points ((-0.0508 0.0508 0.01));
71  }
72  }
73  }
74  \endverbatim
75 
76  Entries:
77  \table
78  Property | Description | Required | Default
79  type | Type-name: sets | yes |
80  sets | Dictionary or list of sample sets | expected |
81  fields | word/regex list of fields to sample | yes |
82  interpolationScheme | scheme to obtain values | no | cellPoint
83  setFormat | output format | yes |
84  sampleOnExecute | Sample (store) on execution as well | no | false
85  formatOptions | dictionary of format options | no |
86  \endtable
87 
88  Additional per-set entries:
89  \table
90  Property | Description | Required | Default
91  store | Store surface/fields on registry | no |
92  setFormat | output format | no |
93  formatOptions | dictionary of format options | no |
94  \endtable
95 
96 Note
97  Special setFormat \c probes can be used to output ensemble results
98  in a format similar to the probes function object.
99 
100 SourceFiles
101  sampledSets.C
102  sampledSetsImpl.C
103 
104 \*---------------------------------------------------------------------------*/
105 
106 #ifndef Foam_sampledSets_H
107 #define Foam_sampledSets_H
108 
109 #include "fvMeshFunctionObject.H"
110 #include "HashPtrTable.H"
111 #include "OFstream.H"
112 #include "sampledSet.H"
113 #include "meshSearch.H"
114 #include "coordSet.H"
115 #include "coordSetWriter.H"
116 #include "volFieldsFwd.H"
117 #include "IOobjectList.H"
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 namespace Foam
122 {
123 
124 // Forward Declarations
125 class Time;
126 class dictionary;
127 class globalIndex;
128 
129 /*---------------------------------------------------------------------------*\
130  Class sampledSets Declaration
131 \*---------------------------------------------------------------------------*/
132 
133 class sampledSets
134 :
135  public functionObjects::fvMeshFunctionObject,
136  public PtrList<sampledSet>
137 {
138  // Static Data Members
139 
140  //- Local control for sampling actions
141  enum sampleActionType : unsigned
142  {
143  ACTION_NONE = 0,
144  ACTION_WRITE = 0x1,
145  ACTION_STORE = 0x2,
146  ACTION_ALL = 0xF
147  };
148 
149  // Private Data
150 
151  //- Keep the dictionary to recreate sets for moving mesh cases
152  dictionary dict_;
153 
154  //- Load fields from files (not from objectRegistry)
155  bool loadFromFiles_;
156 
157  //- Output verbosity
158  bool verbose_;
159 
160  //- Perform sample/store actions on execute as well
161  bool onExecute_;
162 
163  //- Correct meshSearch and update sets
164  bool needsCorrect_;
165 
166  //- Output in raw format like 'probes' does
167  bool writeAsProbes_;
168 
169  //- Output path
170  fileName outputPath_;
171 
172  //- Mesh search engine
173  meshSearch searchEngine_;
174 
175 
176  // Read from dictionary
177 
178  //- Names of fields to sample
179  wordRes fieldSelection_;
180 
181  //- Interpolation/sample scheme to obtain values at the points
182  word samplePointScheme_;
183 
184  //- Output format to use
185  word writeFormat_;
186 
187 
188  // Output control
189 
190  //- Current list of field names selected for sampling
191  DynamicList<word> selectedFieldNames_;
192 
193  //- The coordSet writers (one per sampled set)
194  PtrList<coordSetWriter> writers_;
195 
196  //- Current open files (non-empty on master only) when writing as probes
197  HashPtrTable<OFstream> probeFilePtrs_;
198 
199 
200  // Merging
201 
202  //- Gathered coordSet. (Content only meaningful on master)
203  PtrList<coordSet> gatheredSets_;
204 
205  //- Indexed sort order for gathered coordSet.
206  // (Content only meaningful on master)
207  List<labelList> gatheredSorting_;
208 
209  //- The globalIndex for gathering. (Content only meaningful on master)
210  List<globalIndex> globalIndices_;
211 
212 
213  // Private Member Functions
214 
215  //- Get or create new stream as required (on master)
216  OFstream* createProbeFile(const word& fieldName);
217 
218  //- A new coordSet writer, with per-set formatOptions
219  static autoPtr<coordSetWriter> newWriter
220  (
221  word writerType,
222  const dictionary& topDict,
223  const dictionary& setDict
224  );
225 
226  //- Perform sampling action with store/write
227  bool performAction(unsigned request);
228 
229  //- Count selected/sampled fields
230  // Returns empty IOobjectList if loadFromFiles_ is not active
231  //
232  // Adjusts selectedFieldNames_
233  IOobjectList preCheckFields(unsigned request);
234 
235  //- Setup the sets (optional writers)
236  void initDict(const dictionary& dict, const bool initial);
237 
238  //- Combine points from all processors.
239  //- Sort by curve distance and produce index list.
240  //- Valid result only on master processor.
241  void gatherAllSets();
242 
243  //- Write sampled fieldName on coordSet and on outputDir path
244  template<class Type>
245  void writeCoordSet
246  (
248  const Field<Type>& values,
249  const word& fieldName
250  );
251 
252  //- Get from registry or load from disk
253  template<class GeoField>
254  tmp<GeoField> getOrLoadField(const word& fieldName) const;
255 
256  //- Sample and store/write a specific volume field
257  template<class Type>
258  void performAction(const VolumeField<Type>& field, unsigned request);
259 
260  //- Sample and write all applicable sampled fields
261  // Only uses IOobjectList when loadFromFiles_ is active
262  template<class GeoField>
263  void performAction(const IOobjectList& objects, unsigned request);
264 
265  //- No copy construct
266  sampledSets(const sampledSets&) = delete;
267 
268  //- No copy assignment
269  void operator=(const sampledSets&) = delete;
270 
271 
272 public:
273 
274  //- Runtime type information
275  TypeName("sets");
276 
277 
278  // Constructors
279 
280  //- Construct from Time and dictionary
282  (
283  const word& name,
284  const Time& runTime,
285  const dictionary& dict
286  );
287 
288  //- Construct for given objectRegistry and dictionary
289  // allow the possibility to load fields from files
291  (
292  const word& name,
293  const objectRegistry& obr,
294  const dictionary& dict,
295  const bool loadFromFiles = false
296  );
297 
298 
299  //- Destructor
300  virtual ~sampledSets() = default;
301 
302 
303  // Member Functions
304 
305  //- Enable/disable verbose output
306  // \return old value
307  bool verbose(const bool on) noexcept;
308 
309  //- Return names of fields to sample
310  const wordRes& fieldNames() const noexcept { return fieldSelection_; }
311 
312  //- Read the sampledSets
313  virtual bool read(const dictionary&);
314 
315  //- Execute, currently does nothing
316  virtual bool execute();
317 
318  //- Sample and write
319  virtual bool write();
320 
321  //- Correct for mesh changes
322  void correct();
323 
324  //- Update for changes of mesh
325  virtual void updateMesh(const mapPolyMesh&);
326 
327  //- Update for mesh point-motion
328  virtual void movePoints(const polyMesh&);
329 
330  //- Update for changes of mesh due to readUpdate
331  virtual void readUpdate(const polyMesh::readUpdateState state);
332 };
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 } // End namespace Foam
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 #endif
342 
343 // ************************************************************************* //
dictionary dict
const wordRes & fieldNames() const noexcept
Return names of fields to sample.
Definition: sampledSets.H:442
rDeltaTY field()
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:738
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Output to file stream, using an OSstream.
Definition: OFstream.H:49
engineTime & runTime
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:499
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers...
Definition: HashPtrTable.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
const word & name() const noexcept
Return the name of this functionObject.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
Base class for writing coordSet(s) and tracks with fields.
void correct()
Correct for mesh changes.
Definition: sampledSets.C:723
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:729
bool verbose(const bool on) noexcept
Enable/disable verbose output.
Definition: sampledSets.C:491
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:706
const direction noexcept
Definition: Scalar.H:258
virtual ~sampledSets()=default
Destructor.
TypeName("sets")
Runtime type information.
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:747
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual bool write()
Sample and write.
Definition: sampledSets.C:717
Set of sets to sample.
Definition: sampledSets.H:188
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
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:90
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
Namespace for OpenFOAM.