MFEM  v3.3.2
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 /** A class for collecting finite element data that is part of the same
27  simulation. Currently, this class groups together grid functions (fields),
28  quadrature functions (q-fields), and the mesh that they are defined on. */
30 {
31 public:
32  typedef std::map<std::string, GridFunction*> FieldMapType;
33  typedef std::map<std::string, QuadratureFunction*> QFieldMapType;
34 
35 protected:
36  /// Name of the collection, used as a directory name when saving
37  std::string name;
38 
39  /// A path where the directory with results is saved.
40  /// If not empty, it has '/' at the end.
41  std::string prefix_path;
42 
43  /// The fields and their names (used when saving)
44  typedef FieldMapType::iterator FieldMapIterator;
45  typedef FieldMapType::const_iterator FieldMapConstIterator;
46  /** An std::map containing the registered fields' names as keys (std::string)
47  and their GridFunction pointers as values. */
49 
50  typedef QFieldMapType::iterator QFieldMapIterator;
51  typedef QFieldMapType::const_iterator QFieldMapConstIterator;
53 
54  /// The (common) mesh for the collected fields
56 
57  /** Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
58  When cycle >= 0, it is appended to directory names. */
59  int cycle;
60  /// Physical time (for time-dependent simulations)
61  double time;
62 
63  /// Time step i.e. delta_t (for time-dependent simulations)
64  double time_step;
65 
66  /// Serial or parallel run?
67  bool serial;
68  /// Append rank to any output file names.
70 
71  /// MPI rank (in parallel)
72  int myid;
73  /// Number of MPI ranks (in parallel)
74  int num_procs;
75 
76  /// Precision (number of digits) used for the text output of doubles
77  int precision;
78  /// Number of digits used for the cycle and MPI rank in filenames
80 
81  /// Default value for precision
82  static const int precision_default = 6;
83  /// Default value for pad_digits_*
84  static const int pad_digits_default = 6;
85 
86  /// Output mesh format: 0 - serial format (default), 1 - parallel format
87  int format;
88 
89  /// Should the collection delete its mesh and fields
90  bool own_data;
91 
92  /// Error state
93  int error;
94 
95  /// Delete data owned by the DataCollection keeping field information
96  void DeleteData();
97  /// Delete data owned by the DataCollection including field information
98  void DeleteAll();
99 
100  std::string GetFieldFileName(const std::string &field_name);
101 
102  /// Save one field to disk, assuming the collection directory exists
103  void SaveOneField(const FieldMapIterator &it);
104 
105  /// Save one q-field to disk, assuming the collection directory exists
106  void SaveOneQField(const QFieldMapIterator &it);
107 
108  // Helper method
109  static int create_directory(const std::string &dir_name,
110  const Mesh *mesh, int myid);
111 
112 public:
113  /// Initialize the collection with its name and Mesh.
114  /** When @a mesh_ is NULL, then the real mesh can be set with SetMesh(). */
115  DataCollection(const std::string& collection_name, Mesh *mesh_ = NULL);
116 
117  /// Add a grid function to the collection
118  virtual void RegisterField(const std::string& field_name, GridFunction *gf);
119 
120  /// Remove a grid function from the collection
121  virtual void DeregisterField(const std::string& field_name);
122 
123  /// Add a QuadratureFunction to the collection.
124  virtual void RegisterQField(const std::string& q_field_name,
125  QuadratureFunction *qf);
126 
127  /// Remove a QuadratureFunction from the collection
128  virtual void DeregisterQField(const std::string& field_name);
129 
130  /// Check if a grid function is part of the collection
131  bool HasField(const std::string& name) const
132  { return field_map.find(name) != field_map.end(); }
133 
134  /// Get a pointer to a grid function in the collection.
135  /** Returns NULL if @a field_name is not in the collection. */
136  GridFunction *GetField(const std::string& field_name);
137 
138 #ifdef MFEM_USE_MPI
139  /// Get a pointer to a parallel grid function in the collection.
140  /** Returns NULL if @a field_name is not in the collection.
141  @note The GridFunction pointer stored in the collection is statically
142  cast to ParGridFunction pointer. */
143  ParGridFunction *GetParField(const std::string& field_name)
144  { return static_cast<ParGridFunction*>(GetField(field_name)); }
145 #endif
146 
147  /// Check if a QuadratureFunction with the given name is in the collection.
148  bool HasQField(const std::string& q_field_name) const
149  { return q_field_map.find(q_field_name) != q_field_map.end(); }
150 
151  /// Get a pointer to a QuadratureFunction in the collection.
152  /** Returns NULL if @a field_name is not in the collection. */
153  QuadratureFunction *GetQField(const std::string& q_field_name);
154 
155  /// Get a const reference to the internal field map.
156  /** The keys in the map are the field names and the values are pointers to
157  GridFunction%s. */
158  const FieldMapType &GetFieldMap() const { return field_map; }
159 
160  /// Get a const reference to the internal q-field map.
161  /** The keys in the map are the q-field names and the values are pointers to
162  QuadratureFunction%s. */
163  const QFieldMapType &GetQFieldMap() const { return q_field_map; }
164 
165  /// Get a pointer to the mesh in the collection
166  Mesh *GetMesh() { return mesh; }
167  /// Set/change the mesh associated with the collection
168  virtual void SetMesh(Mesh *new_mesh);
169 
170  /// Set time cycle (for time-dependent simulations)
171  void SetCycle(int c) { cycle = c; }
172  /// Set physical time (for time-dependent simulations)
173  void SetTime(double t) { time = t; }
174 
175  /// Set the simulation time step (for time-dependent simulations)
176  void SetTimeStep(double ts) { time_step = ts; }
177 
178  /// Get time cycle (for time-dependent simulations)
179  int GetCycle() const { return cycle; }
180  /// Get physical time (for time-dependent simulations)
181  double GetTime() const { return time; }
182  /// Get the simulation time step (for time-dependent simulations)
183  double GetTimeStep() const { return time_step; }
184 
185  /// Get the name of the collection
186  const std::string& GetCollectionName() const { return name; }
187  /// Set the ownership of collection data
188  void SetOwnData(bool o) { own_data = o; }
189 
190  /// Set the precision (number of digits) used for the text output of doubles
191  void SetPrecision(int prec) { precision = prec; }
192  /// Set the number of digits used for both the cycle and the MPI rank
193  void SetPadDigits(int digits) { pad_digits_cycle=pad_digits_rank = digits; }
194  /// Set the number of digits used for the cycle
195  void SetPadDigitsCycle(int digits) { pad_digits_cycle = digits; }
196  /// Set the number of digits used for the MPI rank in filenames
197  void SetPadDigitsRank(int digits) { pad_digits_rank = digits; }
198  /** @brief Set the desired output mesh format: 0 - serial format (default),
199  1 - parallel format. */
200  void SetFormat(int fmt) { format = fmt; }
201 
202  /// Set the path where the DataCollection will be saved.
203  void SetPrefixPath(const std::string &prefix);
204 
205  /// Get the path where the DataCollection will be saved.
206  const std::string &GetPrefixPath() const { return prefix_path; }
207 
208  /** Save the collection to disk. By default, everything is saved in a
209  directory with name "collection_name" or "collection_name_cycle" for
210  time-dependent simulations. */
211  virtual void Save();
212  /// Save the mesh, creating the collection directory.
213  virtual void SaveMesh();
214  /// Save one field, assuming the collection directory already exists.
215  virtual void SaveField(const std::string &field_name);
216  /// Save one q-field, assuming the collection directory already exists.
217  virtual void SaveQField(const std::string &q_field_name);
218 
219  /// Load the collection. Not implemented in the base class DataCollection.
220  virtual void Load(int cycle_ = 0);
221 
222  /// Delete the mesh and fields if owned by the collection
223  virtual ~DataCollection();
224 
225  /// Errors returned by Error()
226  enum { NO_ERROR = 0, READ_ERROR = 1, WRITE_ERROR = 2 };
227 
228  /// Get the current error state
229  int Error() const { return error; }
230  /// Reset the error state
231  void ResetError(int err = NO_ERROR) { error = err; }
232 };
233 
234 
235 /// Helper class for VisIt visualization data
237 {
238 public:
239  std::string association;
242  VisItFieldInfo(std::string _association, int _num_components)
243  { association = _association; num_components = _num_components; }
244 };
245 
246 /// Data collection with VisIt I/O routines
248 {
249 protected:
250  // Additional data needed in the VisIt root file, which describes the mesh
251  // and all the fields in the collection
254  std::map<std::string, VisItFieldInfo> field_info_map;
255  typedef std::map<std::string, VisItFieldInfo>::iterator FieldInfoMapIterator;
256 
257  /// Prepare the VisIt root file in JSON format for the current collection
258  std::string GetVisItRootString();
259  /// Read in a VisIt root file in JSON format
260  void ParseVisItRootString(const std::string& json);
261 
262  // Helper functions for Load()
263  void LoadVisItRootFile(const std::string& root_name);
264  void LoadMesh();
265  void LoadFields();
266 
267 public:
268  /// Constructor. The collection name is used when saving the data.
269  /** If @a mesh_ is NULL, then the mesh can be set later by calling either
270  SetMesh() or Load(). The latter works only in serial. */
271  VisItDataCollection(const std::string& collection_name, Mesh *mesh_ = NULL);
272 
273  /// Set/change the mesh associated with the collection
274  virtual void SetMesh(Mesh *new_mesh);
275 
276  /// Add a grid function to the collection and update the root file
277  virtual void RegisterField(const std::string& field_name, GridFunction *gf);
278 
279  /// Set VisIt parameter: maximum levels of detail for the MultiresControl
280  void SetMaxLevelsOfDetail(int max_levels_of_detail);
281 
282  /** Delete all data owned by VisItDataCollection including field data
283  information. */
284  void DeleteAll();
285 
286  /// Save the collection and a VisIt root file
287  virtual void Save();
288 
289  /// Save a VisIt root file for the collection
290  void SaveRootFile();
291 
292  /// Load the collection based on its VisIt data (described in its root file)
293  virtual void Load(int cycle_ = 0);
294 
295  /// We will delete the mesh and fields if we own them
296  virtual ~VisItDataCollection() {}
297 };
298 
299 }
300 
301 #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
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
int format
Output mesh format: 0 - serial format (default), 1 - parallel format.
int error
Error state.
double GetTime() const
Get physical time (for time-dependent simulations)
bool appendRankToFileName
Append rank to any output file names.
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.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
const std::string & GetCollectionName() const
Get the name of the collection.
std::map< std::string, GridFunction * > FieldMapType
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.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
int Error() const
Get the current error state.
bool own_data
Should the collection delete its mesh and fields.
QFieldMapType q_field_map
QuadratureFunction * GetQField(const std::string &q_field_name)
Get a pointer to a QuadratureFunction in the collection.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
QFieldMapType::iterator QFieldMapIterator
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)
double time
Physical time (for time-dependent simulations)
void ResetError(int err=NO_ERROR)
Reset the error state.
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
FieldMapType::const_iterator FieldMapConstIterator
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.
VisItFieldInfo(std::string _association, int _num_components)
bool serial
Serial or parallel run?
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
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.
std::string GetFieldFileName(const std::string &field_name)
void SetFormat(int fmt)
Set the desired output mesh format: 0 - serial format (default), 1 - parallel format.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
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.
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.
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
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.
FieldMapType::iterator FieldMapIterator
The fields and their names (used when saving)
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.
QFieldMapType::const_iterator QFieldMapConstIterator
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)
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
virtual void DeregisterQField(const std::string &field_name)
Remove a QuadratureFunction from the collection.
void SaveRootFile()
Save a VisIt root file for 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.
void LoadVisItRootFile(const std::string &root_name)
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to 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.
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:354
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
std::map< std::string, QuadratureFunction * > QFieldMapType