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