sampledSurfaces.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-2016 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::sampledSurfaces
29
30Description
31 Set of surfaces 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 surfaces;
41 libs (sampling);
42
43 // Write at same frequency as fields
44 writeControl outputTime;
45 writeInterval 1;
46
47 // Fields to be sampled
48 fields (p U);
49
50 // Scheme to obtain face centre value
51 sampleScheme cell;
52
53 // Scheme to obtain node values
54 // (only used if interpolate=true for the surfaces below)
55 interpolationScheme cell;
56
57 // Optional: registry storage
58 store true
59
60 // Output surface format
61 surfaceFormat vtk;
62
63 formatOptions
64 {
65 vtk
66 {
67 precision 10;
68 }
69 }
70
71 surfaces
72 {
73 f0surf
74 {
75 type meshedSurface;
76 surface f0surf.obj;
77 source cells;
78
79 // Optional: keep original regions
80 keepIds true;
81
82 // Optional: generate values on points instead of faces
83 interpolate true;
84
85 // Optional: alternative output type
86 surfaceFormat ensight;
87
88 // Optional: registry storage
89 store true
90 }
91 }
92 }
93 \endverbatim
94
95 Entries:
96 \table
97 Property | Description | Required | Default
98 type | Type-name: surfaces | yes |
99 surfaces | Dictionary or list of sample surfaces | expected |
100 fields | word/regex list of fields to sample | yes |
101 sampleScheme | scheme to obtain face centre value | no | cell
102 interpolationScheme | scheme to obtain node values | no | cellPoint
103 surfaceFormat | output surface format | yes |
104 formatOptions | dictionary of format options | no |
105 sampleOnExecute | Sample (store) on execution as well | no | false
106 store | Store surface/fields on registry | no | false
107 \endtable
108
109 Additional per-surface entries:
110 \table
111 Property | Description | Required | Default
112 store | Store surface/fields on registry | no |
113 surfaceFormat | output surface format | no |
114 formatOptions | dictionary of format options | no |
115 \endtable
116
117Note
118 The interpolationScheme is only used if interpolate=true is used by any
119 of the surfaces.
120
121SourceFiles
122 sampledSurfaces.C
123 sampledSurfacesTemplates.C
124
125\*---------------------------------------------------------------------------*/
126
127#ifndef Foam_sampledSurfaces_H
128#define Foam_sampledSurfaces_H
129
130#include "fvMeshFunctionObject.H"
131#include "sampledSurface.H"
132#include "surfaceWriter.H"
133#include "volFieldsFwd.H"
134#include "surfaceFieldsFwd.H"
135#include "IOobjectList.H"
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139namespace Foam
140{
141
142// Forward Declarations
143class polySurface;
144
145/*---------------------------------------------------------------------------*\
146 Class sampledSurfaces Declaration
147\*---------------------------------------------------------------------------*/
148
149class sampledSurfaces
150:
151 public functionObjects::fvMeshFunctionObject,
152 public PtrList<sampledSurface>
153{
154 // Static Data Members
155
156 //- Tolerance for merging points (fraction of mesh bounding box)
157 static scalar mergeTol_;
158
159 //- Local control for sampling actions
160 enum sampleActionType : unsigned
161 {
162 ACTION_NONE = 0,
163 ACTION_WRITE = 0x1,
164 ACTION_STORE = 0x2,
165 ACTION_ALL = 0xF
166 };
167
168
169 // Private Data
170
171 //- Load fields from files (not from objectRegistry)
172 const bool loadFromFiles_;
173
174 //- Output verbosity
175 bool verbose_;
176
177 //- Perform sample/store actions on execute as well
178 bool onExecute_;
179
180 //- Output path
181 fileName outputPath_;
182
183
184 // Read from dictionary
185
186 //- Names of fields to sample
187 wordRes fieldSelection_;
188
189 //- Interpolation/sample scheme to obtain face values
190 word sampleFaceScheme_;
191
192 //- Interpolation/sample scheme to obtain node values
193 word sampleNodeScheme_;
194
195
196 // Output control
197
198 //- Surface writers (one per surface)
199 PtrList<surfaceWriter> writers_;
200
201 //- Per-surface selection of store/write actions
202 List<unsigned> actions_;
203
204 //- Cached values of the global number of faces per-surface
205 labelList nFaces_;
206
207
208 // Private Member Functions
209
210 //- Return the surfaces
211 const PtrList<sampledSurface>& surfaces() const noexcept
212 {
213 return *this;
214 }
215
216 //- Return the surfaces
217 PtrList<sampledSurface>& surfaces() noexcept
219 return *this;
220 }
221
222 //- A new surfaceWriter, with per-surface formatOptions
223 static autoPtr<surfaceWriter> newWriter
224 (
225 word writeType,
227 const dictionary& surfDict
228 );
229
230
231 //- Perform sampling action with store/write
232 bool performAction(unsigned request);
233
234 //- Count selected/sampled fields per surface
235 // Returns empty IOobjectList if loadFromFiles_ is not active
236 IOobjectList preCheckFields();
237
238 //- Write sampled fieldName on surface and on outputDir path
239 template<class Type>
240 void writeSurface
241 (
243 const Field<Type>& values,
244 const word& fieldName
245 );
246
247 //- Sample and store/write a specific volume field
248 template<class Type>
249 void performAction(const VolumeField<Type>& fld, unsigned request);
250
251 //- Sample and store/write a specific surface field
252 template<class Type>
253 void performAction(const SurfaceField<Type>& fld, unsigned request);
254
255 //- Sample and write all applicable sampled fields
256 // Only uses IOobjectList when loadFromFiles_ is active
257 template<class GeoField>
258 void performAction
259 (
260 const IOobjectList& objects,
261 unsigned request
262 );
263
264
265 //- Get surface from registry if available.
266 // \return surface or nullptr
267 polySurface* getRegistrySurface(const sampledSurface& s) const;
268
269 //- Put surface onto registry, when enabled.
270 // \return surface or nullptr it surface should not be stored
271 polySurface* storeRegistrySurface(const sampledSurface& s);
272
273 //- Remove surface from registry.
274 // \return True if surface existed and was removed
275 bool removeRegistrySurface(const sampledSurface& s);
276
277 //- Store sampled field onto surface registry if it exists
278 template<class Type, class GeoMeshType>
279 bool storeRegistryField
280 (
281 const sampledSurface& s,
282 const word& fieldName,
283 const dimensionSet& dims,
284 Field<Type>&& values
285 );
286
287 //- Test surfaces for condition
288 template<class Container, class Predicate>
289 static bool testAny(const Container& items, const Predicate& pred);
290
291 //- Do any surfaces need an update?
292 virtual bool needsUpdate() const;
293
294 //- Mark the surfaces as needing an update.
295 // Return false if all surfaces were already marked as expired.
296 // Optionally force expire, even if a surface has been marked as
297 // invariant.
298 virtual bool expire(const bool force=false);
299
300 //- Update the surfaces as required.
301 // Return false if no surfaces required an update.
302 virtual bool update();
303
304 //- No copy construct
305 sampledSurfaces(const sampledSurfaces&) = delete;
306
307 //- No copy assignment
308 void operator=(const sampledSurfaces&) = delete;
309
310
311public:
312
313 //- Runtime type information
314 TypeName("surfaces");
315
316
317 // Constructors
318
319 //- Construct from Time and dictionary
321 (
322 const word& name,
323 const Time& runTime,
324 const dictionary& dict
325 );
326
327 //- Construct for given objectRegistry and dictionary
328 // allow the possibility to load fields from files
330 (
331 const word& name,
332 const objectRegistry& obr,
333 const dictionary& dict,
334 const bool loadFromFiles = false
335 );
336
337
338 //- Destructor
339 virtual ~sampledSurfaces() = default;
340
341
342 // Member Functions
343
344 //- Enable/disable verbose output
345 // \return old value
346 bool verbose(const bool on) noexcept;
347
348 //- Read the sampledSurfaces dictionary
349 virtual bool read(const dictionary& dict);
350
351 //- Sample and store if the sampleOnExecute is enabled.
352 virtual bool execute();
353
354 //- Sample and write
355 virtual bool write();
356
357 //- Update for changes of mesh - expires the surfaces
358 virtual void updateMesh(const mapPolyMesh& mpm);
359
360 //- Update for mesh point-motion - expires the surfaces
361 virtual void movePoints(const polyMesh& mesh);
362
363 //- Update for changes of mesh due to readUpdate - expires the surfaces
364 virtual void readUpdate(const polyMesh::readUpdateState state);
365
366 //- Get merge tolerance
367 static scalar mergeTol() noexcept;
368
369 //- Set merge tolerance and return old value
370 static scalar mergeTol(const scalar tol) noexcept;
371};
372
373
374// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375
376} // End namespace Foam
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380#ifdef NoRepository
382#endif
384// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
385
386#endif
387
388// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
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
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
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
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
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:71
An abstract class for surfaces with sampling.
Set of surfaces to sample.
bool verbose(const bool on) noexcept
Enable/disable verbose output.
virtual void movePoints(const polyMesh &mesh)
Update for mesh point-motion - expires the surfaces.
virtual bool read(const dictionary &dict)
Read the sampledSurfaces dictionary.
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh - expires the surfaces.
static scalar mergeTol() noexcept
Get merge tolerance.
virtual ~sampledSurfaces()=default
Destructor.
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
virtual bool execute()
Sample and store if the sampleOnExecute is enabled.
TypeName("surfaces")
Runtime type information.
virtual bool write()
Sample and write.
Base class for surface writers.
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
mesh update()
dynamicFvMesh & mesh
engineTime & runTime
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
const direction noexcept
Definition: Scalar.H:223
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))