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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::sampledSets
29
30Description
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
96Note
97 Special setFormat \c probes can be used to output ensemble results
98 in a format similar to the probes function object.
99
100SourceFiles
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
121namespace Foam
122{
123
124// Forward Declarations
125class Time;
126class dictionary;
127class globalIndex;
128
129/*---------------------------------------------------------------------------*\
130 Class sampledSets Declaration
131\*---------------------------------------------------------------------------*/
132
133class 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 //- Dictionary containing writer options
188 dictionary writeFormatOptions_;
189
190
191 // Output control
193 //- Current list of field names selected for sampling
194 DynamicList<word> selectedFieldNames_;
195
196 //- The coordSet writers (one per sampled set)
198
199 //- Current open files (non-empty on master only) when writing as probes
200 HashPtrTable<OFstream> probeFilePtrs_;
201
202
203 // Merging
204
205 //- Gathered coordSet. (Content only meaningful on master)
206 PtrList<coordSet> gatheredSets_;
207
208 //- Indexed sort order for gathered coordSet.
209 // (Content only meaningful on master)
210 List<labelList> gatheredSorting_;
211
212 //- The globalIndex for gathering. (Content only meaningful on master)
213 List<globalIndex> globalIndices_;
214
215
216 // Private Member Functions
217
218 //- Get or create new stream as required (on master)
219 OFstream* createProbeFile(const word& fieldName);
220
221 //- A new coordSet writer, with per-set formatOptions
222 static autoPtr<coordSetWriter> newWriter
223 (
224 word writeType,
226 const dictionary& surfDict
227 );
228
229 //- Perform sampling action with store/write
230 bool performAction(unsigned request);
231
232 //- Count selected/sampled fields
233 // Returns empty IOobjectList if loadFromFiles_ is not active
234 //
235 // Adjusts selectedFieldNames_
236 IOobjectList preCheckFields(unsigned request);
237
238 //- Setup the sets (optional writers)
239 void initDict(const dictionary& dict, const bool initial);
240
241 //- Combine points from all processors.
242 //- Sort by curve distance and produce index list.
243 //- Valid result only on master processor.
244 void gatherAllSets();
245
246 //- Write sampled fieldName on coordSet and on outputDir path
247 template<class Type>
248 void writeCoordSet
249 (
251 const Field<Type>& values,
252 const word& fieldName
253 );
254
255 //- Get from registry or load from disk
256 template<class GeoField>
257 tmp<GeoField> getOrLoadField(const word& fieldName) const;
258
259 //- Sample and store/write a specific volume field
260 template<class Type>
261 void performAction(const VolumeField<Type>& field, unsigned request);
262
263 //- Sample and write all applicable sampled fields
264 // Only uses IOobjectList when loadFromFiles_ is active
265 template<class GeoField>
266 void performAction(const IOobjectList& objects, unsigned request);
267
268 //- No copy construct
269 sampledSets(const sampledSets&) = delete;
270
271 //- No copy assignment
272 void operator=(const sampledSets&) = delete;
273
274
275public:
276
277 //- Runtime type information
278 TypeName("sets");
279
280
281 // Constructors
282
283 //- Construct from Time and dictionary
285 (
286 const word& name,
287 const Time& runTime,
288 const dictionary& dict
289 );
290
291 //- Construct for given objectRegistry and dictionary
292 // allow the possibility to load fields from files
294 (
295 const word& name,
296 const objectRegistry& obr,
297 const dictionary& dict,
298 const bool loadFromFiles = false
299 );
300
301
302 //- Destructor
303 virtual ~sampledSets() = default;
304
305
306 // Member Functions
307
308 //- Enable/disable verbose output
309 // \return old value
310 bool verbose(const bool on) noexcept;
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// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:68
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Output to file stream, using an OSstream.
Definition: OFstream.H:57
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Base class for writing coordSet(s) and tracks with fields.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
const word & name() const noexcept
Return the name of this functionObject.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
Set of sets to sample.
Definition: sampledSets.H:196
bool verbose(const bool on) noexcept
Enable/disable verbose output.
Definition: sampledSets.C:503
virtual ~sampledSets()=default
Destructor.
void correct()
Correct for mesh changes.
Definition: sampledSets.C:741
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:756
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:765
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:747
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:724
virtual bool write()
Sample and write.
Definition: sampledSets.C:735
TypeName("sets")
Runtime type information.
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:511
A class for managing temporary objects.
Definition: tmp.H:65
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
engineTime & runTime
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
const dictionary formatOptions(propsDict.subOrEmptyDict("formatOptions", keyType::LITERAL))