MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
datacollection.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_DATACOLLECTION
13#define MFEM_DATACOLLECTION
14
15#include "../config/config.hpp"
16#include "gridfunc.hpp"
17#include "qfunction.hpp"
18#ifdef MFEM_USE_MPI
19#include "pgridfunc.hpp"
20#endif
21#include <string>
22#include <map>
23#include <fstream>
24
25namespace mfem
26{
27
28/// Lightweight adaptor over an std::map from strings to pointer to T
29template<typename T>
31{
32public:
33 typedef std::map<std::string, T*> MapType;
34 typedef typename MapType::iterator iterator;
35 typedef typename MapType::const_iterator const_iterator;
36
37 /// Register field @a field with name @a fname
38 /** Replace existing field associated with @a fname (and optionally
39 delete associated pointer if @a own_data is true) */
40 void Register(const std::string& fname, T* field, bool own_data)
41 {
42 T*& ref = field_map[fname];
43 if (own_data)
44 {
45 delete ref; // if newly allocated -> ref is null -> OK
46 }
47 ref = field;
48 }
49
50 /// Unregister association between field @a field and name @a fname
51 /** Optionally delete associated pointer if @a own_data is true */
52 void Deregister(const std::string& fname, bool own_data)
53 {
54 iterator it = field_map.find(fname);
55 if ( it != field_map.end() )
56 {
57 if (own_data)
58 {
59 delete it->second;
60 }
61 field_map.erase(it);
62 }
63 }
64
65 /// Clear all associations between names and fields
66 /** Delete associated pointers when @a own_data is true */
67 void DeleteData(bool own_data)
68 {
69 for (iterator it = field_map.begin(); it != field_map.end(); ++it)
70 {
71 if (own_data)
72 {
73 delete it->second;
74 }
75 it->second = NULL;
76 }
77 }
78
79 /// Predicate to check if a field is associated with name @a fname
80 bool Has(const std::string& fname) const
81 {
82 return field_map.find(fname) != field_map.end();
83 }
84
85 /// Get a pointer to the field associated with name @a fname
86 /** @return Pointer to field associated with @a fname or NULL */
87 T* Get(const std::string& fname) const
88 {
89 const_iterator it = field_map.find(fname);
90 return it != field_map.end() ? it->second : NULL;
91 }
92
93 /// Returns a const reference to the underlying map
94 const MapType& GetMap() const { return field_map; }
95
96 /// Returns the number of registered fields
97 int NumFields() const { return field_map.size(); }
98
99 /// Returns a begin iterator to the registered fields
100 iterator begin() { return field_map.begin(); }
101 /// Returns a begin const iterator to the registered fields
102 const_iterator begin() const { return field_map.begin(); }
103
104 /// Returns an end iterator to the registered fields
105 iterator end() { return field_map.end(); }
106 /// Returns an end const iterator to the registered fields
107 const_iterator end() const { return field_map.end(); }
108
109 /// Returns an iterator to the field @a fname
110 iterator find(const std::string& fname)
111 { return field_map.find(fname); }
112
113 /// Returns a const iterator to the field @a fname
114 const_iterator find(const std::string& fname) const
115 { return field_map.find(fname); }
116
117 /// Clears the map of registered fields without reclaiming memory
118 void clear() { field_map.clear(); }
119
120protected:
122};
123
124
125/** A class for collecting finite element data that is part of the same
126 simulation. Currently, this class groups together grid functions (fields),
127 quadrature functions (q-fields), and the mesh that they are defined on. */
129{
130private:
131 /// A collection of named GridFunctions
133
134 /// A collection of named QuadratureFunctions
136
137public:
141
145
146 /// Format constants to be used with SetFormat().
147 /** Derived classes can define their own format enumerations and override the
148 method SetFormat() to perform input validation. */
150 {
151 SERIAL_FORMAT = 0, /**<
152 MFEM's serial ascii format, using the methods Mesh::Print() /
153 ParMesh::Print(), and GridFunction::Save() / ParGridFunction::Save().*/
154 PARALLEL_FORMAT = 1 /**<
155 MFEM's parallel ascii format, using the methods ParMesh::ParPrint() and
156 GridFunction::Save() / ParGridFunction::Save(). */
157 };
158
159protected:
160 /// Name of the collection, used as a directory name when saving
161 std::string name;
162
163 /** @brief A path where the directory with results is saved.
164 If not empty, it has '/' at the end. */
165 std::string prefix_path;
166
167 /** A FieldMap mapping registered field names to GridFunction pointers. */
169
170 /** A FieldMap mapping registered names to QuadratureFunction pointers. */
172
173 /// The (common) mesh for the collected fields
175
176 /// Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
177 /** When cycle >= 0, it is appended to directory names. */
178 int cycle;
179 /// Physical time (for time-dependent simulations)
181
182 /// Time step i.e. delta_t (for time-dependent simulations)
184
185 /// Serial or parallel run? False iff mesh is a ParMesh
186 bool serial;
187 /// Append rank to any output file names.
189
190 /// MPI rank (in parallel)
191 int myid;
192 /// Number of MPI ranks (in parallel)
194#ifdef MFEM_USE_MPI
195 /// Associated MPI communicator
196 MPI_Comm m_comm;
197#endif
198
199 /// Precision (number of digits) used for the text output of doubles
201 /// Number of digits used for the cycle and MPI rank in filenames
203
204 /// Default value for precision
205 static const int precision_default = 6;
206 /// Default value for pad_digits_*
207 static const int pad_digits_default = 6;
208
209 /// Output mesh format: see the #Format enumeration
212
213 /// Should the collection delete its mesh and fields
215
216 /// Error state
217 int error;
218
219 /// Delete data owned by the DataCollection keeping field information
220 void DeleteData();
221 /// Delete data owned by the DataCollection including field information
222 void DeleteAll();
223
224 std::string GetMeshShortFileName() const;
225 std::string GetMeshFileName() const;
226 std::string GetFieldFileName(const std::string &field_name) const;
227
228 /// Save one field to disk, assuming the collection directory exists
229 void SaveOneField(const FieldMapIterator &it);
230
231 /// Save one q-field to disk, assuming the collection directory exists
232 void SaveOneQField(const QFieldMapIterator &it);
233
234 // Helper method
235 static int create_directory(const std::string &dir_name,
236 const Mesh *mesh, int myid);
237
238public:
239 /// Initialize the collection with its name and Mesh.
240 /** When @a mesh_ is NULL, then the real mesh can be set with SetMesh(). */
241 explicit DataCollection(const std::string& collection_name,
242 Mesh *mesh_ = NULL);
243
244 /// Add a grid function to the collection
245 virtual void RegisterField(const std::string& field_name, GridFunction *gf)
246 { field_map.Register(field_name, gf, own_data); }
247
248 /// Remove a grid function from the collection
249 virtual void DeregisterField(const std::string& field_name)
250 { field_map.Deregister(field_name, own_data); }
251
252 /// Add a QuadratureFunction to the collection.
253 virtual void RegisterQField(const std::string& field_name,
255 { q_field_map.Register(field_name, qf, own_data); }
256
257 /// Remove a QuadratureFunction from the collection
258 virtual void DeregisterQField(const std::string& field_name)
259 { q_field_map.Deregister(field_name, own_data); }
260
261 /// Check if a grid function is part of the collection
262 bool HasField(const std::string& field_name) const
263 { return field_map.Has(field_name); }
264
265 /// Get a pointer to a grid function in the collection.
266 /** Returns NULL if @a field_name is not in the collection. */
267 GridFunction *GetField(const std::string& field_name)
268 { return field_map.Get(field_name); }
269
270#ifdef MFEM_USE_MPI
271 /// Return the associated MPI communicator or MPI_COMM_NULL.
272 MPI_Comm GetComm() const { return m_comm; }
273
274 /// Get a pointer to a parallel grid function in the collection.
275 /** Returns NULL if @a field_name is not in the collection.
276 @note The GridFunction pointer stored in the collection is statically
277 cast to ParGridFunction pointer. */
278 ParGridFunction *GetParField(const std::string& field_name)
279 { return static_cast<ParGridFunction*>(GetField(field_name)); }
280#endif
281
282 /// Check if a QuadratureFunction with the given name is in the collection.
283 bool HasQField(const std::string& field_name) const
284 { return q_field_map.Has(field_name); }
285
286 /// Get a pointer to a QuadratureFunction in the collection.
287 /** Returns NULL if @a field_name is not in the collection. */
288 QuadratureFunction *GetQField(const std::string& field_name)
289 { return q_field_map.Get(field_name); }
290
291 /// Get a const reference to the internal field map.
292 /** The keys in the map are the field names and the values are pointers to
293 GridFunction%s. */
295 { return field_map.GetMap(); }
296
297 /// Get a const reference to the internal q-field map.
298 /** The keys in the map are the q-field names and the values are pointers to
299 QuadratureFunction%s. */
301 { return q_field_map.GetMap(); }
302
303 /// Get a pointer to the mesh in the collection
304 Mesh *GetMesh() { return mesh; }
305
306 /// Set/change the mesh associated with the collection
307 /** When passed a Mesh, assumes the serial case: MPI rank id is set to 0 and
308 MPI num_procs is set to 1. When passed a ParMesh, MPI info from the
309 ParMesh is used to set the DataCollection's MPI rank and num_procs. */
310 virtual void SetMesh(Mesh *new_mesh);
311
312#ifdef MFEM_USE_MPI
313 /// Set/change the mesh associated with the collection.
314 /** For this case, @a comm is used to set the DataCollection's MPI rank id
315 and MPI num_procs, which influences the how files are saved for domain
316 decomposed meshes. */
317 virtual void SetMesh(MPI_Comm comm, Mesh *new_mesh);
318#endif
319
320 /// Set time cycle (for time-dependent simulations)
321 void SetCycle(int c) { cycle = c; }
322 /// Set physical time (for time-dependent simulations)
323 void SetTime(real_t t) { time = t; }
324
325 /// Set the simulation time step (for time-dependent simulations)
326 void SetTimeStep(real_t ts) { time_step = ts; }
327
328 /// Get time cycle (for time-dependent simulations)
329 int GetCycle() const { return cycle; }
330 /// Get physical time (for time-dependent simulations)
331 real_t GetTime() const { return time; }
332 /// Get the simulation time step (for time-dependent simulations)
333 real_t GetTimeStep() const { return time_step; }
334
335 /// Get the name of the collection
336 const std::string& GetCollectionName() const { return name; }
337 /// Set the ownership of collection data
338 void SetOwnData(bool o) { own_data = o; }
339
340 /// Set the precision (number of digits) used for the text output of doubles
341 void SetPrecision(int prec) { precision = prec; }
342 /// Set the number of digits used for both the cycle and the MPI rank
343 virtual void SetPadDigits(int digits)
345 /// Set the number of digits used for the cycle
346 virtual void SetPadDigitsCycle(int digits) { pad_digits_cycle = digits; }
347 /// Set the number of digits used for the MPI rank in filenames
348 virtual void SetPadDigitsRank(int digits) { pad_digits_rank = digits; }
349 /// Set the desired output mesh and data format.
350 /** See the enumeration #Format for valid options. Derived classes can define
351 their own format enumerations and override this method to perform input
352 validation. */
353 virtual void SetFormat(int fmt);
354
355 /// Set the flag for use of gz compressed files
356 virtual void SetCompression(bool comp);
357
358 /// Set the path where the DataCollection will be saved.
359 void SetPrefixPath(const std::string &prefix);
360
361 /// Get the path where the DataCollection will be saved.
362 const std::string &GetPrefixPath() const { return prefix_path; }
363
364 /// Save the collection to disk.
365 /** By default, everything is saved in the "prefix_path" directory with
366 subdirectory name "collection_name" or "collection_name_cycle" for
367 time-dependent simulations. */
368 virtual void Save();
369 /// Save the mesh, creating the collection directory.
370 virtual void SaveMesh();
371 /// Save one field, assuming the collection directory already exists.
372 virtual void SaveField(const std::string &field_name);
373 /// Save one q-field, assuming the collection directory already exists.
374 virtual void SaveQField(const std::string &field_name);
375 /// Load the collection. Not implemented in the base class DataCollection.
376 virtual void Load(int cycle_ = 0);
377
378 /// Delete the mesh and fields if owned by the collection
379 virtual ~DataCollection();
380
381 /// Errors returned by Error()
382 enum
383 {
384 // Workaround for use with headers that define NO_ERROR as a macro,
385 // e.g. winerror.h (which is included by Windows.h):
386#ifndef NO_ERROR
388#endif
389 // Use the following identifier if NO_ERROR is defined as a macro,
390 // e.g. winerror.h (which is included by Windows.h):
393 WRITE_ERROR = 2
394 };
395
396 /// Get the current error state
397 int Error() const { return error; }
398 /// Reset the error state
399 void ResetError(int err_state = No_Error) { error = err_state; }
400
401#ifdef MFEM_USE_MPI
402 friend class ParMesh;
403#endif
404};
405
406
407/// Helper class for VisIt visualization data
409{
410public:
411 std::string association = "";
413 int lod = 1;
414 std::string basis = "";
415 int order = -1;
416 VisItFieldInfo() = default;
417 VisItFieldInfo(std::string association_, int num_components_, int lod_ = 1,
418 std::string basis_ = "", int order_ = -1)
419 {
420 association = association_; num_components = num_components_; lod =lod_;
421 basis = basis_; order = order_;
422 }
423};
424
425/// Data collection with VisIt I/O routines
427{
428protected:
429 // Additional data needed in the VisIt root file, which describes the mesh
430 // and all the fields in the collection
434 std::map<std::string, VisItFieldInfo> field_info_map;
435 typedef std::map<std::string, VisItFieldInfo>::iterator FieldInfoMapIterator;
436
437 /// Prepare the VisIt root file in JSON format for the current collection
438 std::string GetVisItRootString();
439 /// Read in a VisIt root file in JSON format
440 void ParseVisItRootString(const std::string& json);
441
442 void UpdateMeshInfo();
443
444 // Helper functions for Load()
445 void LoadVisItRootFile(const std::string& root_name);
446 void LoadMesh();
447 void LoadFields();
448
449public:
450 /// Constructor. The collection name is used when saving the data.
451 /** If @a mesh_ is NULL, then the mesh can be set later by calling either
452 SetMesh() or Load(). The latter works only in serial. */
453 VisItDataCollection(const std::string& collection_name, Mesh *mesh_ = NULL);
454
455#ifdef MFEM_USE_MPI
456 /// Construct a parallel VisItDataCollection to be loaded from files.
457 /** Before loading the collection with Load(), some parameters in the
458 collection can be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc. */
459 VisItDataCollection(MPI_Comm comm, const std::string& collection_name,
460 Mesh *mesh_ = NULL);
461#endif
462
463 /// Set/change the mesh associated with the collection
464 void SetMesh(Mesh *new_mesh) override;
465
466#ifdef MFEM_USE_MPI
467 /// Set/change the mesh associated with the collection.
468 void SetMesh(MPI_Comm comm, Mesh *new_mesh) override;
469#endif
470
471 /// Add a grid function to the collection and update the root file
472 void RegisterField(const std::string& field_name,
473 GridFunction *gf) override;
474
475 /// Add a quadrature function to the collection and update the root file.
476 /** Visualization of quadrature function is not supported in VisIt(3.12).
477 A patch has been sent to VisIt developers in June 2020. */
478 void RegisterQField(const std::string& q_field_name,
479 QuadratureFunction *qf) override;
480
481 /// Set the number of digits used for both the cycle and the MPI rank
482 /// @note VisIt seems to require 6 pad digits for the MPI rank. Therefore,
483 /// this function uses this default value. This behavior can be overridden
484 /// by calling SetPadDigitsCycle() and SetPadDigitsRank() instead.
485 void SetPadDigits(int digits) override
486 { pad_digits_cycle=digits; pad_digits_rank=6; }
487
488 /// Set VisIt parameter: default levels of detail for the MultiresControl
489 void SetLevelsOfDetail(int levels_of_detail);
490
491 /// Set VisIt parameter: maximum levels of detail for the MultiresControl
492 void SetMaxLevelsOfDetail(int max_levels_of_detail);
493
494 /** @brief Delete all data owned by VisItDataCollection including field data
495 information. */
496 void DeleteAll();
497
498 /// Save the collection and a VisIt root file
499 void Save() override;
500
501 /// Save a VisIt root file for the collection
502 void SaveRootFile();
503
504 /// Load the collection based on its VisIt data (described in its root file)
505 void Load(int cycle_ = 0) override;
506
507 /// We will delete the mesh and fields if we own them
509};
510
511
512/// Abstract base class for ParaViewDataCollection and ParaViewHDFDataCollection
514{
515protected:
518 bool high_order_output = false;
519 bool restart_mode = false;
520 bool bdr_output = false;
522
523public:
524 ParaViewDataCollectionBase(const std::string &name, Mesh *mesh);
525
526 /// @brief Set the refinement level.
527 ///
528 /// In "low-order mode", every element is uniformly split based on the levels
529 /// of detail. In "high-order mode", this sets the polynomial degree for the
530 /// element transformations.
531 ///
532 /// The initial value is 1.
533 void SetLevelsOfDetail(int levels_of_detail_);
534
535 /// @brief Set the zlib compression level.
536 ///
537 /// 0 indicates no compression, -1 indicates the default compression level.
538 /// Otherwise, specify a number between 1 and 9, 1 being the fastest, and 9
539 /// being the best compression. Compression only takes effect if the output
540 /// format is BINARY or BINARY32. MFEM must be compiled with MFEM_USE_ZLIB =
541 /// YES.
542 ///
543 /// The initial compression level is 0 if MFEM is compiled with MFEM_USE_ZLIB
544 /// turned off, and -1 otherwise.
545 ///
546 /// Any nonzero compression level will enable compression.
547 void SetCompressionLevel(int compression_level_);
548
549 /// @brief Sets whether or not to output the data as high-order elements
550 /// (false by default).
551 ///
552 /// Reading high-order data requires ParaView 5.5 or later.
553 void SetHighOrderOutput(bool high_order_output_);
554
555 /// @brief Configures collection to save only fields evaluated on boundaries of
556 /// the mesh.
557 void SetBoundaryOutput(bool bdr_output_);
558
559 /// If compression is enabled, return the compression level, else return 0.
560 int GetCompressionLevel() const;
561
562 /// @brief Set the data format for the ParaView output files.
563 ///
564 /// Possible options are VTKFormat::ASCII, VTKFormat::BINARY, and
565 /// VTKFormat::BINARY32. The ASCII and BINARY options output double precision
566 /// data, whereas the BINARY32 option outputs single precision data.
567 ///
568 /// The initial format is VTKFormat::BINARY.
569 ///
570 /// VTKFormat::ASCII is not supported by ParaViewHDFDataCollection.
571 void SetDataFormat(VTKFormat fmt);
572
573 /// Returns true if the output format is BINARY or BINARY32, false if ASCII.
574 bool IsBinaryFormat() const;
575
576 /// @brief Enable or disable restart mode.
577 ///
578 /// If restart is enabled, new writes will preserve timestep metadata for any
579 /// solutions prior to the currently defined time.
580 void UseRestartMode(bool restart_mode_);
581};
582
583/// Writer for ParaView visualization (PVD and VTU format)
585{
586private:
587 std::fstream pvd_stream;
588
589 /// A collection of named Coefficients and VectorCoefficients
592
593 /** A FieldMap mapping registered names to Coefficient and VectorCoefficient
594 pointers. */
595 CoeffFieldMap coeff_field_map;
596 VCoeffFieldMap vcoeff_field_map;
597protected:
598 void WritePVTUHeader(std::ostream &out);
599 void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix);
600 void SaveDataVTU(std::ostream &out, int ref);
601 void SaveGFieldVTU(std::ostream& out, int ref_, const FieldMapIterator& it);
602 void SaveCoeffFieldVTU(std::ostream& out, int ref_, const std::string &name,
603 Coefficient &coeff);
604 void SaveVCoeffFieldVTU(std::ostream& out, int ref_, const std::string &name,
605 VectorCoefficient& coeff);
606 const char *GetDataFormatString() const;
607 const char *GetDataTypeString() const;
608
609 std::string GenerateCollectionPath();
610 std::string GenerateVTUFileName(const std::string &prefix, int rank);
611 std::string GenerateVTUPath();
612 std::string GeneratePVDFileName();
613 std::string GeneratePVTUFileName(const std::string &prefix);
614 std::string GeneratePVTUPath();
615
616public:
617 /// Constructor. The collection name is used when saving the data.
618 /** If @a mesh_ is NULL, then the mesh can be set later by calling SetMesh().
619 Before saving the data collection, some parameters in the collection can
620 be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc. */
621 ParaViewDataCollection(const std::string& collection_name,
622 Mesh *mesh_ = nullptr);
623
624 /// Get a const reference to the internal coefficient-field map.
626 { return coeff_field_map.GetMap(); }
628 { return vcoeff_field_map.GetMap(); }
629
630 /// Add a Coefficient or VectorCoefficient to the collection.
631 void RegisterCoeffField(const std::string& field_name, Coefficient *coeff)
632 { coeff_field_map.Register(field_name, coeff, own_data); }
633 void RegisterVCoeffField(const std::string& field_name,
634 VectorCoefficient *vcoeff)
635 { vcoeff_field_map.Register(field_name, vcoeff, own_data); }
636
637 /// Remove a Coefficient or VectorCoefficient from the collection
638 void DeregisterCoeffField(const std::string& field_name)
639 { coeff_field_map.Deregister(field_name, own_data); }
640 void DeregisterVCoeffField(const std::string& field_name)
641 { vcoeff_field_map.Deregister(field_name, own_data); }
642
643 /// Save the collection - the directory name is constructed based on the
644 /// cycle value
645 void Save() override;
646};
647
648#ifdef MFEM_USE_HDF5
649
650/// Writer for ParaView visualization (%VTKHDF format)
652{
653 /// The low-level VTKHDF object for I/O (pointer to implementation idiom).
654 std::unique_ptr<class VTKHDF> vtkhdf;
655
656 /// Create the VTKHDF object if it doesn't exist already.
657 void EnsureVTKHDF();
658
659 /// Save the collection (templated on floating point type).
660 template <typename FP_T> void TSave();
661
662public:
663 /// @brief Constructor. The collection name is used when saving the data.
664 ///
665 /// If @a mesh_ is NULL, then the mesh can be set later by calling SetMesh().
666 /// Before saving the data collection, some parameters in the collection can
667 /// be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc.
668 ParaViewHDFDataCollection(const std::string& collection_name,
669 Mesh *mesh_ = nullptr);
670
671 /// @brief Enable or disable compression.
672 ///
673 /// The compression level can be set with SetCompressionLevel()). VTKHDF
674 /// compression does not require MFEM to be compiled with zlib support.
675 void SetCompression(bool compression_) override;
676
677 /// Save the collection.
678 void Save() override;
679
680 /// Destructor.
682};
683
684#endif
685
686}
687#endif
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
real_t time
Physical time (for time-dependent simulations)
virtual void RegisterQField(const std::string &field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
MPI_Comm GetComm() const
Return the associated MPI communicator or MPI_COMM_NULL.
void ResetError(int err_state=No_Error)
Reset the error state.
QFieldMap::iterator QFieldMapIterator
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
const std::string & GetCollectionName() const
Get the name of the collection.
bool own_data
Should the collection delete its mesh and fields.
DataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Initialize the collection with its name and Mesh.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
int GetCycle() const
Get time cycle (for time-dependent simulations)
GFieldMap::iterator FieldMapIterator
void SetTimeStep(real_t ts)
Set the simulation time step (for time-dependent simulations)
void SaveOneField(const FieldMapIterator &it)
Save one field to disk, assuming the collection directory exists.
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
void DeleteAll()
Delete data owned by the DataCollection including field information.
GFieldMap::const_iterator FieldMapConstIterator
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
QFieldMap::const_iterator QFieldMapConstIterator
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
static const int precision_default
Default value for precision.
bool HasField(const std::string &field_name) const
Check if a grid function is part of the collection.
virtual void SetPadDigitsRank(int digits)
Set the number of digits used for the MPI rank in filenames.
virtual void SaveQField(const std::string &field_name)
Save one q-field, assuming the collection directory already exists.
int Error() const
Get the current error state.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
bool HasQField(const std::string &field_name) const
Check if a QuadratureFunction with the given name is in the collection.
bool serial
Serial or parallel run? False iff mesh is a ParMesh.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has '/' at the end.
virtual void SetPadDigits(int digits)
Set the number of digits used for both the cycle and the MPI rank.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
std::string GetFieldFileName(const std::string &field_name) const
virtual void DeregisterQField(const std::string &field_name)
Remove a QuadratureFunction from the collection.
int myid
MPI rank (in parallel)
ParGridFunction * GetParField(const std::string &field_name)
Get a pointer to a parallel grid function in the collection.
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
static const int pad_digits_default
Default value for pad_digits_*.
real_t time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
const std::string & GetPrefixPath() const
Get the path where the DataCollection will be saved.
const QFieldMapType & GetQFieldMap() const
Get a const reference to the internal q-field map.
std::string GetMeshShortFileName() const
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
virtual ~DataCollection()
Delete the mesh and fields if owned by the collection.
real_t GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
std::string GetMeshFileName() const
Format
Format constants to be used with SetFormat().
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
bool appendRankToFileName
Append rank to any output file names.
void SetOwnData(bool o)
Set the ownership of collection data.
std::string name
Name of the collection, used as a directory name when saving.
virtual void Save()
Save the collection to disk.
virtual void SaveField(const std::string &field_name)
Save one field, assuming the collection directory already exists.
QuadratureFunction * GetQField(const std::string &field_name)
Get a pointer to a QuadratureFunction in the collection.
MPI_Comm m_comm
Associated MPI communicator.
GFieldMap::MapType FieldMapType
QFieldMap::MapType QFieldMapType
virtual void SetPadDigitsCycle(int digits)
Set the number of digits used for the cycle.
Mesh * mesh
The (common) mesh for the collected fields.
virtual void SaveMesh()
Save the mesh, creating the collection directory.
int precision
Precision (number of digits) used for the text output of doubles.
int format
Output mesh format: see the Format enumeration.
real_t GetTime() const
Get physical time (for time-dependent simulations)
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
Mesh data type.
Definition mesh.hpp:65
Lightweight adaptor over an std::map from strings to pointer to T.
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
iterator end()
Returns an end iterator to the registered fields.
MapType::iterator iterator
iterator begin()
Returns a begin iterator to the registered fields.
const_iterator end() const
Returns an end const iterator to the registered fields.
std::map< std::string, T * > MapType
void DeleteData(bool own_data)
Clear all associations between names and fields.
MapType::const_iterator const_iterator
bool Has(const std::string &fname) const
Predicate to check if a field is associated with name fname.
const_iterator find(const std::string &fname) const
Returns a const iterator to the field fname.
void clear()
Clears the map of registered fields without reclaiming memory.
void Deregister(const std::string &fname, bool own_data)
Unregister association between field field and name fname.
T * Get(const std::string &fname) const
Get a pointer to the field associated with name fname.
const_iterator begin() const
Returns a begin const iterator to the registered fields.
iterator find(const std::string &fname)
Returns an iterator to the field fname.
int NumFields() const
Returns the number of registered fields.
const MapType & GetMap() const
Returns a const reference to the underlying map.
Class for parallel grid function.
Definition pgridfunc.hpp:50
Class for parallel meshes.
Definition pmesh.hpp:34
Abstract base class for ParaViewDataCollection and ParaViewHDFDataCollection.
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
ParaViewDataCollectionBase(const std::string &name, Mesh *mesh)
void SetLevelsOfDetail(int levels_of_detail_)
Set the refinement level.
void UseRestartMode(bool restart_mode_)
Enable or disable restart mode.
void SetBoundaryOutput(bool bdr_output_)
Configures collection to save only fields evaluated on boundaries of the mesh.
void SetHighOrderOutput(bool high_order_output_)
Sets whether or not to output the data as high-order elements (false by default).
void SetDataFormat(VTKFormat fmt)
Set the data format for the ParaView output files.
void SetCompressionLevel(int compression_level_)
Set the zlib compression level.
int GetCompressionLevel() const
If compression is enabled, return the compression level, else return 0.
Writer for ParaView visualization (PVD and VTU format)
void SaveVCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, VectorCoefficient &coeff)
void WritePVTUHeader(std::ostream &out)
const VCoeffFieldMap::MapType & GetVCoeffFieldMap() const
void DeregisterVCoeffField(const std::string &field_name)
const CoeffFieldMap::MapType & GetCoeffFieldMap() const
Get a const reference to the internal coefficient-field map.
void SaveGFieldVTU(std::ostream &out, int ref_, const FieldMapIterator &it)
const char * GetDataFormatString() const
void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix)
void SaveDataVTU(std::ostream &out, int ref)
void DeregisterCoeffField(const std::string &field_name)
Remove a Coefficient or VectorCoefficient from the collection.
void RegisterCoeffField(const std::string &field_name, Coefficient *coeff)
Add a Coefficient or VectorCoefficient to the collection.
ParaViewDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
void RegisterVCoeffField(const std::string &field_name, VectorCoefficient *vcoeff)
void SaveCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, Coefficient &coeff)
std::string GenerateVTUFileName(const std::string &prefix, int rank)
const char * GetDataTypeString() const
std::string GeneratePVTUFileName(const std::string &prefix)
Writer for ParaView visualization (VTKHDF format)
void Save() override
Save the collection.
ParaViewHDFDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
void SetCompression(bool compression_) override
Enable or disable compression.
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
Base class for vector Coefficients that optionally depend on time and space.
Data collection with VisIt I/O routines.
void SaveRootFile()
Save a VisIt root file for the collection.
void SetPadDigits(int digits) override
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf) override
Add a quadrature function to the collection and update the root file.
void LoadVisItRootFile(const std::string &root_name)
virtual ~VisItDataCollection()
We will delete the mesh and fields if we own them.
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
void Load(int cycle_=0) override
Load the collection based on its VisIt data (described in its root file)
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
void Save() override
Save the collection and a VisIt root file.
void SetMaxLevelsOfDetail(int max_levels_of_detail)
Set VisIt parameter: maximum levels of detail for the MultiresControl.
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
void SetMesh(Mesh *new_mesh) override
Set/change the mesh associated with the collection.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
std::map< std::string, VisItFieldInfo > field_info_map
void DeleteAll()
Delete all data owned by VisItDataCollection including field data information.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
Helper class for VisIt visualization data.
VisItFieldInfo(std::string association_, int num_components_, int lod_=1, std::string basis_="", int order_=-1)
VisItFieldInfo()=default
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
VTKFormat
Data array format for VTK and VTU files.
Definition vtk.hpp:100
float real_t
Definition config.hpp:46