MFEM  v4.0
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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_DATACOLLECTION
13 #define MFEM_DATACOLLECTION
14 
15 #include "../config/config.hpp"
16 #include "gridfunc.hpp"
17 #ifdef MFEM_USE_MPI
18 #include "pgridfunc.hpp"
19 #endif
20 #include <string>
21 #include <map>
22 
23 namespace mfem
24 {
25 
26 /// Lightweight adaptor over an std::map from strings to pointer to T
27 template<typename T>
29 {
30 public:
31  typedef std::map<std::string, T*> MapType;
32  typedef typename MapType::iterator iterator;
33  typedef typename MapType::const_iterator const_iterator;
34 
35  /// Register field @a field with name @a fname
36  /** Replace existing field associated with @a fname (and optionally
37  delete associated pointer if @a own_data is true) */
38  void Register(const std::string& fname, T* field, bool own_data)
39  {
40  T*& ref = field_map[fname];
41  if (own_data)
42  {
43  delete ref; // if newly allocated -> ref is null -> OK
44  }
45  ref = field;
46  }
47 
48  /// Unregister association between field @a field and name @a fname
49  /** Optionally delete associated pointer if @a own_data is true */
50  void Deregister(const std::string& fname, bool own_data)
51  {
52  iterator it = field_map.find(fname);
53  if ( it != field_map.end() )
54  {
55  if (own_data)
56  {
57  delete it->second;
58  }
59  field_map.erase(it);
60  }
61  }
62 
63  /// Clear all associations between names and fields
64  /** Delete associated pointers when @a own_data is true */
65  void DeleteData(bool own_data)
66  {
67  for (iterator it = field_map.begin(); it != field_map.end(); ++it)
68  {
69  if (own_data)
70  {
71  delete it->second;
72  }
73  it->second = NULL;
74  }
75  }
76 
77  /// Predicate to check if a field is associated with name @a fname
78  bool Has(const std::string& fname) const
79  {
80  return field_map.find(fname) != field_map.end();
81  }
82 
83  /// Get a pointer to the field associated with name @a fname
84  /** @return Pointer to field associated with @a fname or NULL */
85  T* Get(const std::string& fname) const
86  {
87  const_iterator it = field_map.find(fname);
88  return it != field_map.end() ? it->second : NULL;
89  }
90 
91  /// Returns a const reference to the underlying map
92  const MapType& GetMap() const { return field_map; }
93 
94  /// Returns the number of registered fields
95  int NumFields() const { return field_map.size(); }
96 
97  /// Returns a begin iterator to the registered fields
98  iterator begin() { return field_map.begin(); }
99  /// Returns a begin const iterator to the registered fields
100  const_iterator begin() const { return field_map.begin(); }
101 
102  /// Returns an end iterator to the registered fields
103  iterator end() { return field_map.end(); }
104  /// Returns an end const iterator to the registered fields
105  const_iterator end() const { return field_map.end(); }
106 
107  /// Returns an iterator to the field @a fname
108  iterator find(const std::string& fname)
109  { return field_map.find(fname); }
110 
111  /// Returns a const iterator to the field @a fname
112  const_iterator find(const std::string& fname) const
113  { return field_map.find(fname); }
114 
115  /// Clears the map of registered fields without reclaiming memory
116  void clear() { field_map.clear(); }
117 
118 protected:
120 };
121 
122 
123 /** A class for collecting finite element data that is part of the same
124  simulation. Currently, this class groups together grid functions (fields),
125  quadrature functions (q-fields), and the mesh that they are defined on. */
127 {
128 private:
129  /// A collection of named GridFunctions
131 
132  /// A collection of named QuadratureFunctions
134 public:
138 
142 
143  /// Format constants to be used with SetFormat().
144  /** Derived classes can define their own format enumerations and override the
145  method SetFormat() to perform input validation. */
146  enum Format
147  {
148  SERIAL_FORMAT = 0, /**<
149  MFEM's serial ascii format, using the methods Mesh::Print() /
150  ParMesh::Print(), and GridFunction::Save() / ParGridFunction::Save().*/
152  MFEM's parallel ascii format, using the methods ParMesh::ParPrint() and
153  GridFunction::Save() / ParGridFunction::Save(). */
154  };
155 
156 protected:
157  /// Name of the collection, used as a directory name when saving
158  std::string name;
159 
160  /** @brief A path where the directory with results is saved.
161  If not empty, it has '/' at the end. */
162  std::string prefix_path;
163 
164  /** A FieldMap mapping registered field names to GridFunction pointers. */
166 
167  /** A FieldMap mapping registered names to QuadratureFunction pointers. */
169 
170  /// The (common) mesh for the collected fields
172 
173  /// Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
174  /** When cycle >= 0, it is appended to directory names. */
175  int cycle;
176  /// Physical time (for time-dependent simulations)
177  double time;
178 
179  /// Time step i.e. delta_t (for time-dependent simulations)
180  double time_step;
181 
182  /// Serial or parallel run? False iff mesh is a ParMesh
183  bool serial;
184  /// Append rank to any output file names.
186 
187  /// MPI rank (in parallel)
188  int myid;
189  /// Number of MPI ranks (in parallel)
191 #ifdef MFEM_USE_MPI
192  /// Associated MPI communicator
193  MPI_Comm m_comm;
194 #endif
195 
196  /// Precision (number of digits) used for the text output of doubles
198  /// Number of digits used for the cycle and MPI rank in filenames
200 
201  /// Default value for precision
202  static const int precision_default = 6;
203  /// Default value for pad_digits_*
204  static const int pad_digits_default = 6;
205 
206  /// Output mesh format: see the #Format enumeration
207  int format;
209 
210  /// Should the collection delete its mesh and fields
211  bool own_data;
212 
213  /// Error state
214  int error;
215 
216  /// Delete data owned by the DataCollection keeping field information
217  void DeleteData();
218  /// Delete data owned by the DataCollection including field information
219  void DeleteAll();
220 
221  std::string GetMeshShortFileName() const;
222  std::string GetMeshFileName() const;
223  std::string GetFieldFileName(const std::string &field_name) const;
224 
225  /// Save one field to disk, assuming the collection directory exists
226  void SaveOneField(const FieldMapIterator &it);
227 
228  /// Save one q-field to disk, assuming the collection directory exists
229  void SaveOneQField(const QFieldMapIterator &it);
230 
231  // Helper method
232  static int create_directory(const std::string &dir_name,
233  const Mesh *mesh, int myid);
234 
235 public:
236  /// Initialize the collection with its name and Mesh.
237  /** When @a mesh_ is NULL, then the real mesh can be set with SetMesh(). */
238  explicit DataCollection(const std::string& collection_name,
239  Mesh *mesh_ = NULL);
240 
241  /// Add a grid function to the collection
242  virtual void RegisterField(const std::string& field_name, GridFunction *gf)
243  { field_map.Register(field_name, gf, own_data); }
244 
245  /// Remove a grid function from the collection
246  virtual void DeregisterField(const std::string& field_name)
247  { field_map.Deregister(field_name, own_data); }
248 
249  /// Add a QuadratureFunction to the collection.
250  virtual void RegisterQField(const std::string& q_field_name,
251  QuadratureFunction *qf)
252  { q_field_map.Register(q_field_name, qf, own_data); }
253 
254 
255  /// Remove a QuadratureFunction from the collection
256  virtual void DeregisterQField(const std::string& field_name)
257  { q_field_map.Deregister(field_name, own_data); }
258 
259  /// Check if a grid function is part of the collection
260  bool HasField(const std::string& name) const
261  { return field_map.Has(name); }
262 
263  /// Get a pointer to a grid function in the collection.
264  /** Returns NULL if @a field_name is not in the collection. */
265  GridFunction *GetField(const std::string& field_name)
266  { return field_map.Get(field_name); }
267 
268 #ifdef MFEM_USE_MPI
269  /// Return the associated MPI communicator or MPI_COMM_NULL.
270  MPI_Comm GetComm() const { return m_comm; }
271 
272  /// Get a pointer to a parallel grid function in the collection.
273  /** Returns NULL if @a field_name is not in the collection.
274  @note The GridFunction pointer stored in the collection is statically
275  cast to ParGridFunction pointer. */
276  ParGridFunction *GetParField(const std::string& field_name)
277  { return static_cast<ParGridFunction*>(GetField(field_name)); }
278 #endif
279 
280  /// Check if a QuadratureFunction with the given name is in the collection.
281  bool HasQField(const std::string& q_field_name) const
282  { return q_field_map.Has(q_field_name); }
283 
284  /// Get a pointer to a QuadratureFunction in the collection.
285  /** Returns NULL if @a field_name is not in the collection. */
286  QuadratureFunction *GetQField(const std::string& q_field_name)
287  { return q_field_map.Get(q_field_name); }
288 
289  /// Get a const reference to the internal field map.
290  /** The keys in the map are the field names and the values are pointers to
291  GridFunction%s. */
292  const FieldMapType &GetFieldMap() const
293  { return field_map.GetMap(); }
294 
295  /// Get a const reference to the internal q-field map.
296  /** The keys in the map are the q-field names and the values are pointers to
297  QuadratureFunction%s. */
299  { return q_field_map.GetMap(); }
300 
301  /// Get a pointer to the mesh in the collection
302  Mesh *GetMesh() { return mesh; }
303  /// Set/change the mesh associated with the collection
304  /** When passed a Mesh, assumes the serial case: MPI rank id is set to 0 and
305  MPI num_procs is set to 1. When passed a ParMesh, MPI info from the
306  ParMesh is used to set the DataCollection's MPI rank and num_procs. */
307  virtual void SetMesh(Mesh *new_mesh);
308 #ifdef MFEM_USE_MPI
309  /// Set/change the mesh associated with the collection.
310  /** For this case, @a comm is used to set the DataCollection's MPI rank id
311  and MPI num_procs, which influences the how files are saved for domain
312  decomposed meshes. */
313  virtual void SetMesh(MPI_Comm comm, Mesh *new_mesh);
314 #endif
315 
316  /// Set time cycle (for time-dependent simulations)
317  void SetCycle(int c) { cycle = c; }
318  /// Set physical time (for time-dependent simulations)
319  void SetTime(double t) { time = t; }
320 
321  /// Set the simulation time step (for time-dependent simulations)
322  void SetTimeStep(double ts) { time_step = ts; }
323 
324  /// Get time cycle (for time-dependent simulations)
325  int GetCycle() const { return cycle; }
326  /// Get physical time (for time-dependent simulations)
327  double GetTime() const { return time; }
328  /// Get the simulation time step (for time-dependent simulations)
329  double GetTimeStep() const { return time_step; }
330 
331  /// Get the name of the collection
332  const std::string& GetCollectionName() const { return name; }
333  /// Set the ownership of collection data
334  void SetOwnData(bool o) { own_data = o; }
335 
336  /// Set the precision (number of digits) used for the text output of doubles
337  void SetPrecision(int prec) { precision = prec; }
338  /// Set the number of digits used for both the cycle and the MPI rank
339  void SetPadDigits(int digits) { pad_digits_cycle=pad_digits_rank = digits; }
340  /// Set the number of digits used for the cycle
341  void SetPadDigitsCycle(int digits) { pad_digits_cycle = digits; }
342  /// Set the number of digits used for the MPI rank in filenames
343  void SetPadDigitsRank(int digits) { pad_digits_rank = digits; }
344  /// Set the desired output mesh and data format.
345  /** See the enumeration #Format for valid options. Derived classes can define
346  their own format enumerations and override this method to perform input
347  validation. */
348  virtual void SetFormat(int fmt);
349 
350  /// Set the flag for use of gz compressed files
351  void SetCompression(bool comp);
352 
353  /// Set the path where the DataCollection will be saved.
354  void SetPrefixPath(const std::string &prefix);
355 
356  /// Get the path where the DataCollection will be saved.
357  const std::string &GetPrefixPath() const { return prefix_path; }
358 
359  /// Save the collection to disk.
360  /** By default, everything is saved in the "prefix_path" directory with
361  subdirectory name "collection_name" or "collection_name_cycle" for
362  time-dependent simulations. */
363  virtual void Save();
364  /// Save the mesh, creating the collection directory.
365  virtual void SaveMesh();
366  /// Save one field, assuming the collection directory already exists.
367  virtual void SaveField(const std::string &field_name);
368  /// Save one q-field, assuming the collection directory already exists.
369  virtual void SaveQField(const std::string &q_field_name);
370 
371  /// Load the collection. Not implemented in the base class DataCollection.
372  virtual void Load(int cycle_ = 0);
373 
374  /// Delete the mesh and fields if owned by the collection
375  virtual ~DataCollection();
376 
377  /// Errors returned by Error()
378  enum { NO_ERROR = 0, READ_ERROR = 1, WRITE_ERROR = 2 };
379 
380  /// Get the current error state
381  int Error() const { return error; }
382  /// Reset the error state
383  void ResetError(int err = NO_ERROR) { error = err; }
384 };
385 
386 
387 /// Helper class for VisIt visualization data
389 {
390 public:
391  std::string association;
394  VisItFieldInfo(std::string _association, int _num_components)
395  { association = _association; num_components = _num_components; }
396 };
397 
398 /// Data collection with VisIt I/O routines
400 {
401 protected:
402  // Additional data needed in the VisIt root file, which describes the mesh
403  // and all the fields in the collection
407  std::map<std::string, VisItFieldInfo> field_info_map;
408  typedef std::map<std::string, VisItFieldInfo>::iterator FieldInfoMapIterator;
409 
410  /// Prepare the VisIt root file in JSON format for the current collection
411  std::string GetVisItRootString();
412  /// Read in a VisIt root file in JSON format
413  void ParseVisItRootString(const std::string& json);
414 
415  void UpdateMeshInfo();
416 
417  // Helper functions for Load()
418  void LoadVisItRootFile(const std::string& root_name);
419  void LoadMesh();
420  void LoadFields();
421 
422 public:
423  /// Constructor. The collection name is used when saving the data.
424  /** If @a mesh_ is NULL, then the mesh can be set later by calling either
425  SetMesh() or Load(). The latter works only in serial. */
426  VisItDataCollection(const std::string& collection_name, Mesh *mesh_ = NULL);
427 
428 #ifdef MFEM_USE_MPI
429  /// Construct a parallel VisItDataCollection to be loaded from files.
430  /** Before loading the collection with Load(), some parameters in the
431  collection can be adjusted, e.g. SetPadDigits(), SetPrefixPath(), etc. */
432  VisItDataCollection(MPI_Comm comm, const std::string& collection_name,
433  Mesh *mesh_ = NULL);
434 #endif
435 
436  /// Set/change the mesh associated with the collection
437  virtual void SetMesh(Mesh *new_mesh);
438 
439 #ifdef MFEM_USE_MPI
440  /// Set/change the mesh associated with the collection.
441  virtual void SetMesh(MPI_Comm comm, Mesh *new_mesh);
442 #endif
443 
444  /// Add a grid function to the collection and update the root file
445  virtual void RegisterField(const std::string& field_name, GridFunction *gf);
446 
447  /// Set VisIt parameter: default levels of detail for the MultiresControl
448  void SetLevelsOfDetail(int levels_of_detail);
449 
450  /// Set VisIt parameter: maximum levels of detail for the MultiresControl
451  void SetMaxLevelsOfDetail(int max_levels_of_detail);
452 
453  /** @brief Delete all data owned by VisItDataCollection including field data
454  information. */
455  void DeleteAll();
456 
457  /// Save the collection and a VisIt root file
458  virtual void Save();
459 
460  /// Save a VisIt root file for the collection
461  void SaveRootFile();
462 
463  /// Load the collection based on its VisIt data (described in its root file)
464  virtual void Load(int cycle_ = 0);
465 
466  /// We will delete the mesh and fields if we own them
467  virtual ~VisItDataCollection() {}
468 };
469 
470 }
471 
472 #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.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
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.
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.
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.
bool HasQField(const std::string &q_field_name) const
Check if a QuadratureFunction with the given name is in the collection.
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
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.
bool HasField(const std::string &name) const
Check if a grid function is part of the collection.
virtual void Save()
Save the collection and a VisIt root file.
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.
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.
void SetPadDigitsCycle(int digits)
Set the number of digits used for the cycle.
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
MPI_Comm m_comm
Associated MPI communicator.
double time
Physical time (for time-dependent simulations)
MapType::iterator iterator
void ResetError(int err=NO_ERROR)
Reset the error state.
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
iterator find(const std::string &fname)
Returns an iterator to the field fname.
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.
void SetPadDigits(int digits)
Set the number of digits used for both the cycle and the MPI rank.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
VisItFieldInfo(std::string _association, int _num_components)
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.
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.
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 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.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
static const int pad_digits_default
Default value for pad_digits_*.
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
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
void Deregister(const std::string &fname, bool own_data)
Unregister association between field field and name fname.
void SetOwnData(bool o)
Set the ownership of collection data.
virtual ~VisItDataCollection()
We will delete the mesh and fields if we own them.
void SetPadDigitsRank(int digits)
Set the number of digits used for the MPI rank in filenames.
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.
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.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:69
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 void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
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.
const_iterator end() const
Returns an end const iterator to the registered fields.
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.
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.
std::string name
Name of the collection, used as a directory name when saving.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
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)
int num_procs
Number of MPI ranks (in parallel)
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:494
std::map< std::string, T * > MapType
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.