MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
datacollection.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 
25 namespace mfem
26 {
27 
28 /// Lightweight adaptor over an std::map from strings to pointer to T
29 template<typename T>
31 {
32 public:
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 
120 protected:
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 {
130 private:
131  /// A collection of named GridFunctions
133 
134  /// A collection of named QuadratureFunctions
136 public:
140 
144 
145  /// Format constants to be used with SetFormat().
146  /** Derived classes can define their own format enumerations and override the
147  method SetFormat() to perform input validation. */
148  enum Format
149  {
150  SERIAL_FORMAT = 0, /**<
151  MFEM's serial ascii format, using the methods Mesh::Print() /
152  ParMesh::Print(), and GridFunction::Save() / ParGridFunction::Save().*/
154  MFEM's parallel ascii format, using the methods ParMesh::ParPrint() and
155  GridFunction::Save() / ParGridFunction::Save(). */
156  };
157 
158 protected:
159  /// Name of the collection, used as a directory name when saving
160  std::string name;
161 
162  /** @brief A path where the directory with results is saved.
163  If not empty, it has '/' at the end. */
164  std::string prefix_path;
165 
166  /** A FieldMap mapping registered field names to GridFunction pointers. */
168 
169  /** A FieldMap mapping registered names to QuadratureFunction pointers. */
171 
172  /// The (common) mesh for the collected fields
174 
175  /// Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
176  /** When cycle >= 0, it is appended to directory names. */
177  int cycle;
178  /// Physical time (for time-dependent simulations)
179  double time;
180 
181  /// Time step i.e. delta_t (for time-dependent simulations)
182  double time_step;
183 
184  /// Serial or parallel run? False iff mesh is a ParMesh
185  bool serial;
186  /// Append rank to any output file names.
188 
189  /// MPI rank (in parallel)
190  int myid;
191  /// Number of MPI ranks (in parallel)
193 #ifdef MFEM_USE_MPI
194  /// Associated MPI communicator
195  MPI_Comm m_comm;
196 #endif
197 
198  /// Precision (number of digits) used for the text output of doubles
200  /// Number of digits used for the cycle and MPI rank in filenames
202 
203  /// Default value for precision
204  static const int precision_default = 6;
205  /// Default value for pad_digits_*
206  static const int pad_digits_default = 6;
207 
208  /// Output mesh format: see the #Format enumeration
209  int format;
211 
212  /// Should the collection delete its mesh and fields
213  bool own_data;
214 
215  /// Error state
216  int error;
217 
218  /// Delete data owned by the DataCollection keeping field information
219  void DeleteData();
220  /// Delete data owned by the DataCollection including field information
221  void DeleteAll();
222 
223  std::string GetMeshShortFileName() const;
224  std::string GetMeshFileName() const;
225  std::string GetFieldFileName(const std::string &field_name) const;
226 
227  /// Save one field to disk, assuming the collection directory exists
228  void SaveOneField(const FieldMapIterator &it);
229 
230  /// Save one q-field to disk, assuming the collection directory exists
231  void SaveOneQField(const QFieldMapIterator &it);
232 
233  // Helper method
234  static int create_directory(const std::string &dir_name,
235  const Mesh *mesh, int myid);
236 
237 public:
238  /// Initialize the collection with its name and Mesh.
239  /** When @a mesh_ is NULL, then the real mesh can be set with SetMesh(). */
240  explicit DataCollection(const std::string& collection_name,
241  Mesh *mesh_ = NULL);
242 
243  /// Add a grid function to the collection
244  virtual void RegisterField(const std::string& field_name, GridFunction *gf)
245  { field_map.Register(field_name, gf, own_data); }
246 
247  /// Remove a grid function from the collection
248  virtual void DeregisterField(const std::string& field_name)
249  { field_map.Deregister(field_name, own_data); }
250 
251  /// Add a QuadratureFunction to the collection.
252  virtual void RegisterQField(const std::string& q_field_name,
253  QuadratureFunction *qf)
254  { q_field_map.Register(q_field_name, qf, own_data); }
255 
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& q_field_name) const
284  { return q_field_map.Has(q_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& q_field_name)
289  { return q_field_map.Get(q_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. */
294  const FieldMapType &GetFieldMap() const
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  /// Set/change the mesh associated with the collection
306  /** When passed a Mesh, assumes the serial case: MPI rank id is set to 0 and
307  MPI num_procs is set to 1. When passed a ParMesh, MPI info from the
308  ParMesh is used to set the DataCollection's MPI rank and num_procs. */
309  virtual void SetMesh(Mesh *new_mesh);
310 #ifdef MFEM_USE_MPI
311  /// Set/change the mesh associated with the collection.
312  /** For this case, @a comm is used to set the DataCollection's MPI rank id
313  and MPI num_procs, which influences the how files are saved for domain
314  decomposed meshes. */
315  virtual void SetMesh(MPI_Comm comm, Mesh *new_mesh);
316 #endif
317 
318  /// Set time cycle (for time-dependent simulations)
319  void SetCycle(int c) { cycle = c; }
320  /// Set physical time (for time-dependent simulations)
321  void SetTime(double t) { time = t; }
322 
323  /// Set the simulation time step (for time-dependent simulations)
324  void SetTimeStep(double ts) { time_step = ts; }
325 
326  /// Get time cycle (for time-dependent simulations)
327  int GetCycle() const { return cycle; }
328  /// Get physical time (for time-dependent simulations)
329  double GetTime() const { return time; }
330  /// Get the simulation time step (for time-dependent simulations)
331  double GetTimeStep() const { return time_step; }
332 
333  /// Get the name of the collection
334  const std::string& GetCollectionName() const { return name; }
335  /// Set the ownership of collection data
336  void SetOwnData(bool o) { own_data = o; }
337 
338  /// Set the precision (number of digits) used for the text output of doubles
339  void SetPrecision(int prec) { precision = prec; }
340  /// Set the number of digits used for both the cycle and the MPI rank
341  virtual void SetPadDigits(int digits)
342  { pad_digits_cycle=pad_digits_rank = digits; }
343  /// Set the number of digits used for the cycle
344  virtual void SetPadDigitsCycle(int digits) { pad_digits_cycle = digits; }
345  /// Set the number of digits used for the MPI rank in filenames
346  virtual void SetPadDigitsRank(int digits) { pad_digits_rank = digits; }
347  /// Set the desired output mesh and data format.
348  /** See the enumeration #Format for valid options. Derived classes can define
349  their own format enumerations and override this method to perform input
350  validation. */
351  virtual void SetFormat(int fmt);
352 
353  /// Set the flag for use of gz compressed files
354  virtual void SetCompression(bool comp);
355 
356  /// Set the path where the DataCollection will be saved.
357  void SetPrefixPath(const std::string &prefix);
358 
359  /// Get the path where the DataCollection will be saved.
360  const std::string &GetPrefixPath() const { return prefix_path; }
361 
362  /// Save the collection to disk.
363  /** By default, everything is saved in the "prefix_path" directory with
364  subdirectory name "collection_name" or "collection_name_cycle" for
365  time-dependent simulations. */
366  virtual void Save();
367  /// Save the mesh, creating the collection directory.
368  virtual void SaveMesh();
369  /// Save one field, assuming the collection directory already exists.
370  virtual void SaveField(const std::string &field_name);
371  /// Save one q-field, assuming the collection directory already exists.
372  virtual void SaveQField(const std::string &q_field_name);
373 
374  /// Load the collection. Not implemented in the base class DataCollection.
375  virtual void Load(int cycle_ = 0);
376 
377  /// Delete the mesh and fields if owned by the collection
378  virtual ~DataCollection();
379 
380  /// Errors returned by Error()
381  enum { NO_ERROR = 0, READ_ERROR = 1, WRITE_ERROR = 2 };
382 
383  /// Get the current error state
384  int Error() const { return error; }
385  /// Reset the error state
386  void ResetError(int err_state = NO_ERROR) { error = err_state; }
387 
388 #ifdef MFEM_USE_MPI
389  friend class ParMesh;
390 #endif
391 };
392 
393 
394 /// Helper class for VisIt visualization data
396 {
397 public:
398  std::string association;
400  int lod;
402  VisItFieldInfo(std::string association_, int num_components_, int lod_ = 1)
403  { association = association_; num_components = num_components_; lod =lod_;}
404 };
405 
406 /// Data collection with VisIt I/O routines
408 {
409 protected:
410  // Additional data needed in the VisIt root file, which describes the mesh
411  // and all the fields in the collection
415  std::map<std::string, VisItFieldInfo> field_info_map;
416  typedef std::map<std::string, VisItFieldInfo>::iterator FieldInfoMapIterator;
417 
418  /// Prepare the VisIt root file in JSON format for the current collection
419  std::string GetVisItRootString();
420  /// Read in a VisIt root file in JSON format
421  void ParseVisItRootString(const std::string& json);
422 
423  void UpdateMeshInfo();
424 
425  // Helper functions for Load()
426  void LoadVisItRootFile(const std::string& root_name);
427  void LoadMesh();
428  void LoadFields();
429 
430 public:
431  /// Constructor. The collection name is used when saving the data.
432  /** If @a mesh_ is NULL, then the mesh can be set later by calling either
433  SetMesh() or Load(). The latter works only in serial. */
434  VisItDataCollection(const std::string& collection_name, Mesh *mesh_ = NULL);
435 
436 #ifdef MFEM_USE_MPI
437  /// Construct a parallel VisItDataCollection to be loaded from files.
438  /** Before loading the collection with Load(), some parameters in the
439  collection can be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc. */
440  VisItDataCollection(MPI_Comm comm, const std::string& collection_name,
441  Mesh *mesh_ = NULL);
442 #endif
443 
444  /// Set/change the mesh associated with the collection
445  virtual void SetMesh(Mesh *new_mesh) override;
446 
447 #ifdef MFEM_USE_MPI
448  /// Set/change the mesh associated with the collection.
449  virtual void SetMesh(MPI_Comm comm, Mesh *new_mesh) override;
450 #endif
451 
452  /// Add a grid function to the collection and update the root file
453  virtual void RegisterField(const std::string& field_name,
454  GridFunction *gf) override;
455 
456  /// Add a quadrature function to the collection and update the root file.
457  /** Visualization of quadrature function is not supported in VisIt(3.12).
458  A patch has been sent to VisIt developers in June 2020. */
459  virtual void RegisterQField(const std::string& q_field_name,
460  QuadratureFunction *qf) override;
461 
462  /// Set the number of digits used for both the cycle and the MPI rank
463  /// @note VisIt seems to require 6 pad digits for the MPI rank. Therefore,
464  /// this function uses this default value. This behavior can be overridden
465  /// by calling SetPadDigitsCycle() and SetPadDigitsRank() instead.
466  virtual void SetPadDigits(int digits) override
467  { pad_digits_cycle=digits; pad_digits_rank=6; }
468 
469  /// Set VisIt parameter: default levels of detail for the MultiresControl
470  void SetLevelsOfDetail(int levels_of_detail);
471 
472  /// Set VisIt parameter: maximum levels of detail for the MultiresControl
473  void SetMaxLevelsOfDetail(int max_levels_of_detail);
474 
475  /** @brief Delete all data owned by VisItDataCollection including field data
476  information. */
477  void DeleteAll();
478 
479  /// Save the collection and a VisIt root file
480  virtual void Save() override;
481 
482  /// Save a VisIt root file for the collection
483  void SaveRootFile();
484 
485  /// Load the collection based on its VisIt data (described in its root file)
486  virtual void Load(int cycle_ = 0) override;
487 
488  /// We will delete the mesh and fields if we own them
489  virtual ~VisItDataCollection() {}
490 };
491 
492 
493 /// Helper class for ParaView visualization data
495 {
496 private:
497  int levels_of_detail;
498  std::fstream pvd_stream;
499  VTKFormat pv_data_format;
500  bool high_order_output;
501  bool restart_mode;
502 
503 protected:
504  void WritePVTUHeader(std::ostream &out);
505  void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix);
506  void SaveDataVTU(std::ostream &out, int ref);
507  void SaveGFieldVTU(std::ostream& out, int ref_, const FieldMapIterator& it);
508  const char *GetDataFormatString() const;
509  const char *GetDataTypeString() const;
510 
511  std::string GenerateCollectionPath();
512  std::string GenerateVTUFileName(const std::string &prefix, int rank);
513  std::string GenerateVTUPath();
514  std::string GeneratePVDFileName();
515  std::string GeneratePVTUFileName(const std::string &prefix);
516  std::string GeneratePVTUPath();
517 
518 
519 public:
520  /// Constructor. The collection name is used when saving the data.
521  /** If @a mesh_ is NULL, then the mesh can be set later by calling SetMesh().
522  Before saving the data collection, some parameters in the collection can
523  be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc. */
524  ParaViewDataCollection(const std::string& collection_name,
525  mfem::Mesh *mesh_ = NULL);
526 
527  /// Set refinement levels - every element is uniformly split based on
528  /// levels_of_detail_
529  void SetLevelsOfDetail(int levels_of_detail_);
530 
531  /// Save the collection - the directory name is constructed based on the
532  /// cycle value
533  virtual void Save() override;
534 
535  /// Set the data format for the ParaView output files. Possible options are
536  /// VTKFormat::ASCII, VTKFormat::BINARY, and VTKFormat::BINARY32.
537  /// The ASCII and BINARY options output double precision data, whereas the
538  /// BINARY32 option outputs single precision data.
539  void SetDataFormat(VTKFormat fmt);
540 
541  /// Set the zlib compression level. 0 indicates no compression, -1 indicates
542  /// the default compression level. Otherwise, specify a number between 1 and
543  /// 9, 1 being the fastest, and 9 being the best compression. Compression
544  /// only takes effect if the output format is BINARY or BINARY32. MFEM must
545  /// be compiled with MFEM_USE_ZLIB = YES.
546  void SetCompressionLevel(int compression_level_);
547 
548  /// Enable or disable zlib compression. If the input is true, use the default
549  /// zlib compression level (unless the compression level has previously been
550  /// set by calling SetCompressionLevel).
551  void SetCompression(bool compression_) override;
552 
553  /// Returns true if the output format is BINARY or BINARY32, false if ASCII.
554  bool IsBinaryFormat() const;
555 
556  /// Sets whether or not to output the data as high-order elements (false
557  /// by default). Reading high-order data requires ParaView 5.5 or later.
558  void SetHighOrderOutput(bool high_order_output_);
559 
560  /// Enable or disable restart mode. If restart is enabled, new writes will
561  /// preserve timestep metadata for any solutions prior to the currently
562  /// defined time.
563  void UseRestartMode(bool restart_mode_);
564 
565  /// Load the collection - not implemented in the ParaView writer
566  virtual void Load(int cycle_ = 0) override;
567 };
568 
569 }
570 #endif
virtual void SaveMesh()
Save the mesh, creating the collection directory.
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.
virtual void SetPadDigitsRank(int digits)
Set the number of digits used for the MPI rank in filenames.
void WritePVTUHeader(std::ostream &out)
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const MapType & GetMap() const
Returns a const reference to the underlying map.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
int format
Output mesh format: see the Format enumeration.
int error
Error state.
void SetDataFormat(VTKFormat fmt)
double GetTime() const
Get physical time (for time-dependent simulations)
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
void ResetError(int err_state=NO_ERROR)
Reset the error state.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
void DeleteAll()
Delete data owned by the DataCollection including field information.
const char * GetDataTypeString() const
bool HasQField(const std::string &q_field_name) const
Check if a QuadratureFunction with the given name is in the collection.
Helper class for ParaView visualization data.
const_iterator find(const std::string &fname) const
Returns a const iterator to the field fname.
iterator end()
Returns an end iterator to the registered fields.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
GFieldMap::iterator FieldMapIterator
bool HasField(const std::string &field_name) const
Check if a grid function is part of the collection.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
Format
Format constants to be used with SetFormat().
const std::string & GetCollectionName() const
Get the name of the collection.
Helper class for VisIt visualization data.
MPI_Comm GetComm() const
Return the associated MPI communicator or MPI_COMM_NULL.
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
QuadratureFunction * GetQField(const std::string &q_field_name)
Get a pointer to a QuadratureFunction in the collection.
QFieldMap::MapType QFieldMapType
bool Has(const std::string &fname) const
Predicate to check if a field is associated with name fname.
int Error() const
Get the current error state.
bool own_data
Should the collection delete its mesh and fields.
VisItFieldInfo(std::string association_, int num_components_, int lod_=1)
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
virtual void SaveQField(const std::string &q_field_name)
Save one q-field, assuming the collection directory already exists.
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
MPI_Comm m_comm
Associated MPI communicator.
virtual void Load(int cycle_=0) override
Load the collection - not implemented in the ParaView writer.
double time
Physical time (for time-dependent simulations)
MapType::iterator iterator
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
iterator find(const std::string &fname)
Returns an iterator to the field fname.
void SetCompressionLevel(int compression_level_)
std::string GetMeshShortFileName() const
MapType::const_iterator const_iterator
Mesh * mesh
The (common) mesh for the collected fields.
const QFieldMapType & GetQFieldMap() const
Get a const reference to the internal q-field map.
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
bool serial
Serial or parallel run? False iff mesh is a ParMesh.
int cycle
Time cycle; for time-dependent simulations cycle &gt;= 0, otherwise = -1.
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
virtual void Save()
Save the collection to disk.
std::string GenerateVTUFileName(const std::string &prefix, int rank)
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
int GetCycle() const
Get time cycle (for time-dependent simulations)
Data collection with VisIt I/O routines.
void SaveDataVTU(std::ostream &out, int ref)
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:96
void SaveGFieldVTU(std::ostream &out, int ref_, const FieldMapIterator &it)
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix)
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
DataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Initialize the collection with its name and Mesh.
virtual void SetMesh(Mesh *new_mesh) override
Set/change the mesh associated with the collection.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
virtual 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 SetHighOrderOutput(bool high_order_output_)
static const int pad_digits_default
Default value for pad_digits_*.
void SetCompression(bool compression_) override
const char * GetDataFormatString() const
void SetTime(double t)
Set physical time (for time-dependent simulations)
const std::string & GetPrefixPath() const
Get the path where the DataCollection will be saved.
QFieldMap::const_iterator QFieldMapConstIterator
std::string GeneratePVTUFileName(const std::string &prefix)
void Deregister(const std::string &fname, bool own_data)
Unregister association between field field and name fname.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf) override
Add a quadrature function to the collection and update the root file.
void SetOwnData(bool o)
Set the ownership of collection data.
virtual void SetPadDigits(int digits)
Set the number of digits used for both the cycle and the MPI rank.
virtual ~VisItDataCollection()
We will delete the mesh and fields if we own them.
std::map< std::string, VisItFieldInfo > field_info_map
ParGridFunction * GetParField(const std::string &field_name)
Get a pointer to a parallel grid function in the collection.
void UseRestartMode(bool restart_mode_)
int NumFields() const
Returns the number of registered fields.
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
virtual void Save() override
Save the collection and a VisIt root file.
int myid
MPI rank (in parallel)
virtual void SaveField(const std::string &field_name)
Save one field, assuming the collection directory already exists.
static const int precision_default
Default value for precision.
const_iterator begin() const
Returns a begin const iterator to the registered fields.
GFieldMap::MapType FieldMapType
void SetMaxLevelsOfDetail(int max_levels_of_detail)
Set VisIt parameter: maximum levels of detail for the MultiresControl.
virtual ~DataCollection()
Delete the mesh and fields if owned by the collection.
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
std::string GetFieldFileName(const std::string &field_name) const
GFieldMap::const_iterator FieldMapConstIterator
void clear()
Clears the map of registered fields without reclaiming memory.
void SaveRootFile()
Save a VisIt root file for the collection.
virtual void DeregisterQField(const std::string &field_name)
Remove a QuadratureFunction from the collection.
virtual void SetPadDigitsCycle(int digits)
Set the number of digits used for the cycle.
virtual void SetPadDigits(int digits) override
const_iterator end() const
Returns an end const iterator to the registered fields.
ParaViewDataCollection(const std::string &collection_name, mfem::Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has &#39;/&#39; at the end...
Lightweight adaptor over an std::map from strings to pointer to T.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
RefCoord t[3]
void SetLevelsOfDetail(int levels_of_detail_)
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
std::string GetMeshFileName() const
void LoadVisItRootFile(const std::string &root_name)
T * Get(const std::string &fname) const
Get a pointer to the field associated with name fname.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
virtual void Save() override
std::string name
Name of the collection, used as a directory name when saving.
void DeleteAll()
Delete all data owned by VisItDataCollection including field data information.
void DeleteData(bool own_data)
Clear all associations between names and fields.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
double time_step
Time step i.e. delta_t (for time-dependent simulations)
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
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
Definition: pmesh.hpp:32
Represents values or vectors of values at quadrature points on a mesh.
Definition: qfunction.hpp:23
std::map< std::string, T * > MapType
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.