MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
datacollection.cpp
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 #include "fem.hpp"
13 #include "../general/text.hpp"
14 #include "picojson.h"
15 
16 #include <fstream>
17 #include <cerrno> // errno
18 #include <sstream>
19 
20 #ifndef _WIN32
21 #include <sys/stat.h> // mkdir
22 #else
23 #include <direct.h> // _mkdir
24 #define mkdir(dir, mode) _mkdir(dir)
25 #endif
26 
27 namespace mfem
28 {
29 
30 // static method
31 int DataCollection::create_directory(const std::string &dir_name,
32  const Mesh *mesh, int myid)
33 {
34  // create directories recursively
35  const char path_delim = '/';
36  std::string::size_type pos = 0;
37  int err;
38 #ifdef MFEM_USE_MPI
39  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
40 #endif
41 
42  do
43  {
44  pos = dir_name.find(path_delim, pos+1);
45  std::string subdir = dir_name.substr(0, pos);
46 
47 #ifndef MFEM_USE_MPI
48  err = mkdir(subdir.c_str(), 0777);
49  err = (err && (errno != EEXIST)) ? 1 : 0;
50 #else
51  if (myid == 0 || pmesh == NULL)
52  {
53  err = mkdir(subdir.c_str(), 0777);
54  err = (err && (errno != EEXIST)) ? 1 : 0;
55  }
56 #endif
57  }
58  while ( pos != std::string::npos );
59 
60 #ifdef MFEM_USE_MPI
61  if (pmesh)
62  {
63  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
64  }
65 #endif
66 
67  return err;
68 }
69 
70 // class DataCollection implementation
71 
72 DataCollection::DataCollection(const std::string& collection_name, Mesh *mesh_)
73 {
74  name = collection_name;
75  // leave prefix_path empty
76  mesh = mesh_;
77  myid = 0;
78  num_procs = 1;
79  serial = true;
80  appendRankToFileName = false;
81 
82 #ifdef MFEM_USE_MPI
83  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
84  if (par_mesh)
85  {
86  myid = par_mesh->GetMyRank();
87  num_procs = par_mesh->GetNRanks();
88  serial = false;
89  appendRankToFileName = true;
90  }
91 #endif
92  own_data = false;
93  cycle = -1;
94  time = 0.0;
95  time_step = 0.0;
98  format = 0; // use older serial mesh format
99  error = NO_ERROR;
100 }
101 
103 {
104  if (own_data && new_mesh != mesh) { delete mesh; }
105  mesh = new_mesh;
106  myid = 0;
107  num_procs = 1;
108  serial = true;
109  appendRankToFileName = false;
110 
111 #ifdef MFEM_USE_MPI
112  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
113  if (par_mesh)
114  {
115  myid = par_mesh->GetMyRank();
116  num_procs = par_mesh->GetNRanks();
117  serial = false;
118  appendRankToFileName = true;
119  }
120 #endif
121 }
122 
123 void DataCollection::RegisterField(const std::string& name, GridFunction *gf)
124 {
125  GridFunction *&ref = field_map[name];
126  if (own_data)
127  {
128  delete ref; // if newly allocated -> ref is null -> OK
129  }
130  ref = gf;
131 }
132 
133 void DataCollection::DeregisterField(const std::string& name)
134 {
135  FieldMapIterator it = field_map.find(name);
136  if (it != field_map.end())
137  {
138  if (own_data)
139  {
140  delete it->second;
141  }
142  field_map.erase(it);
143  }
144 }
145 
146 void DataCollection::RegisterQField(const std::string& q_field_name,
147  QuadratureFunction *qf)
148 {
149  QuadratureFunction *&ref = q_field_map[q_field_name];
150  if (own_data)
151  {
152  delete ref; // if newly allocated -> ref is null -> OK
153  }
154  ref = qf;
155 }
156 
157 void DataCollection::DeregisterQField(const std::string& name)
158 {
159  QFieldMapIterator it = q_field_map.find(name);
160  if (it != q_field_map.end())
161  {
162  if (own_data)
163  {
164  delete it->second;
165  }
166  q_field_map.erase(it);
167  }
168 }
169 
170 GridFunction *DataCollection::GetField(const std::string& field_name)
171 {
172  FieldMapConstIterator it = field_map.find(field_name);
173 
174  return (it != field_map.end()) ? it->second : NULL;
175 }
176 
177 QuadratureFunction *DataCollection::GetQField(const std::string& q_field_name)
178 {
179  QFieldMapConstIterator it = q_field_map.find(q_field_name);
180 
181  return (it != q_field_map.end()) ? it->second : NULL;
182 }
183 
184 void DataCollection::SetPrefixPath(const std::string& prefix)
185 {
186  if (!prefix.empty())
187  {
188  prefix_path = prefix;
189  if (!prefix_path.empty() && prefix_path[prefix_path.size()-1] != '/')
190  {
191  prefix_path += '/';
192  }
193  }
194  else
195  {
196  prefix_path.clear();
197  }
198 }
199 
200 void DataCollection::Load(int cycle)
201 {
202  MFEM_ABORT("this method is not implemented");
203 }
204 
206 {
207  SaveMesh();
208 
209  if (error) { return; }
210 
211  for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
212  {
213  SaveOneField(it);
214  // Even if there is an error, try saving the other fields
215  }
216 
217  for (QFieldMapIterator it = q_field_map.begin(); it != q_field_map.end();
218  ++it)
219  {
220  SaveOneQField(it);
221  }
222 }
223 
225 {
226  int err;
227 
228  std::string dir_name = prefix_path + name;
229  if (cycle != -1)
230  {
231  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
232  }
233  err = create_directory(dir_name, mesh, myid);
234  if (err)
235  {
236  error = WRITE_ERROR;
237  MFEM_WARNING("Error creating directory: " << dir_name);
238  return; // do not even try to write the mesh
239  }
240 
241  std::string mesh_name = dir_name +
242  ((serial || format == 0 )? "/mesh" : "/pmesh");
244  {
245  mesh_name += "." + to_padded_string(myid, pad_digits_rank);
246  }
247  std::ofstream mesh_file(mesh_name.c_str());
248  mesh_file.precision(precision);
249 #ifdef MFEM_USE_MPI
250  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
251  if (pmesh && format == 1 )
252  {
253  pmesh->ParPrint(mesh_file);
254  }
255  else
256 #endif
257  {
258  mesh->Print(mesh_file);
259  }
260  if (!mesh_file)
261  {
262  error = WRITE_ERROR;
263  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
264  }
265 }
266 
267 std::string DataCollection::GetFieldFileName(const std::string &field_name)
268 {
269  std::string dir_name = prefix_path + name;
270  if (cycle != -1)
271  {
272  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
273  }
274  std::string file_name = dir_name + "/" + field_name;
276  {
277  file_name += "." + to_padded_string(myid, pad_digits_rank);
278  }
279  return file_name;
280 }
281 
283 {
284  std::ofstream field_file(GetFieldFileName(it->first).c_str());
285  field_file.precision(precision);
286  (it->second)->Save(field_file);
287  if (!field_file)
288  {
289  error = WRITE_ERROR;
290  MFEM_WARNING("Error writing field to file: " << it->first);
291  }
292 }
293 
295 {
296  std::ofstream q_field_file(GetFieldFileName(it->first).c_str());
297  q_field_file.precision(precision);
298  (it->second)->Save(q_field_file);
299  if (!q_field_file)
300  {
301  error = WRITE_ERROR;
302  MFEM_WARNING("Error writing q-field to file: " << it->first);
303  }
304 }
305 
306 void DataCollection::SaveField(const std::string &field_name)
307 {
308  FieldMapIterator it = field_map.find(field_name);
309  if (it != field_map.end())
310  {
311  SaveOneField(it);
312  }
313 }
314 
315 void DataCollection::SaveQField(const std::string &q_field_name)
316 {
317  QFieldMapIterator it = q_field_map.find(q_field_name);
318  if (it != q_field_map.end())
319  {
320  SaveOneQField(it);
321  }
322 }
323 
325 {
326  if (own_data) { delete mesh; }
327  mesh = NULL;
328 
329  for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
330  {
331  if (own_data) { delete it->second; }
332  it->second = NULL;
333  }
334  for (QFieldMapIterator it = q_field_map.begin();
335  it != q_field_map.end(); ++it)
336  {
337  if (own_data) { delete it->second; }
338  it->second = NULL;
339  }
340  own_data = false;
341 }
342 
344 {
345  DeleteData();
346  field_map.clear();
347  q_field_map.clear();
348 }
349 
351 {
352  DeleteData();
353 }
354 
355 
356 // class VisItDataCollection implementation
357 
358 VisItDataCollection::VisItDataCollection(const std::string& collection_name,
359  Mesh *mesh)
360  : DataCollection(collection_name, mesh)
361 {
362  appendRankToFileName = true; // always include rank in file names
363  cycle = 0; // always include cycle in directory names
364 
365  if (mesh)
366  {
367  spatial_dim = mesh->SpaceDimension();
368  topo_dim = mesh->Dimension();
369  }
370  else
371  {
372  spatial_dim = 0;
373  topo_dim = 0;
374  }
376 }
377 
379 {
380  DataCollection::SetMesh(new_mesh);
381  appendRankToFileName = true;
383  topo_dim = mesh->Dimension();
384 }
385 
386 void VisItDataCollection::RegisterField(const std::string& name,
387  GridFunction *gf)
388 {
390  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
391 }
392 
393 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
394 {
395  visit_max_levels_of_detail = max_levels_of_detail;
396 }
397 
399 {
400  field_info_map.clear();
402 }
403 
405 {
407  SaveRootFile();
408 }
409 
411 {
412  if (myid != 0) { return; }
413 
414  std::string root_name = prefix_path + name + "_" +
416  ".mfem_root";
417  std::ofstream root_file(root_name.c_str());
418  root_file << GetVisItRootString();
419  if (!root_file)
420  {
421  error = WRITE_ERROR;
422  MFEM_WARNING("Error writting VisIt Root file: " << root_name);
423  }
424 }
425 
427 {
428  DeleteAll();
429  cycle = cycle_;
430  std::string root_name = prefix_path + name + "_" +
432  ".mfem_root";
433  LoadVisItRootFile(root_name);
434  if (!error)
435  {
436  LoadMesh();
437  }
438  if (!error)
439  {
440  LoadFields();
441  }
442  if (!error)
443  {
444  own_data = true;
445  }
446  else
447  {
448  DeleteAll();
449  }
450 }
451 
452 void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
453 {
454  std::ifstream root_file(root_name.c_str());
455  std::stringstream buffer;
456  buffer << root_file.rdbuf();
457  if (!buffer)
458  {
459  error = READ_ERROR;
460  MFEM_WARNING("Error reading the VisIt Root file: " << root_name);
461  }
462  else
463  {
464  ParseVisItRootString(buffer.str());
465  }
466 }
467 
469 {
470  std::string mesh_fname = prefix_path + name + "_" +
473  named_ifgzstream file(mesh_fname.c_str());
474  if (!file)
475  {
476  error = READ_ERROR;
477  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
478  return;
479  }
480  // TODO: 1) load parallel mesh on one processor
481  // 2) load parallel mesh on the same number of processors
482  mesh = new Mesh(file, 1, 1);
484  topo_dim = mesh->Dimension();
485 }
486 
488 {
489  std::string path_left = prefix_path + name + "_" +
491  std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
492 
493  field_map.clear();
494  for (FieldInfoMapIterator it = field_info_map.begin();
495  it != field_info_map.end(); ++it)
496  {
497  std::string fname = path_left + it->first + path_right;
498  std::ifstream file(fname.c_str());
499  if (!file)
500  {
501  error = READ_ERROR;
502  MFEM_WARNING("Unable to open field file: " << fname);
503  return;
504  }
505  // TODO: 1) load parallel GridFunction on one processor
506  // 2) load parallel GridFunction on the same number of processors
507  field_map[it->first] = new GridFunction(mesh, file);
508  }
509 }
510 
512 {
513  // Get the path string (relative to where the root file is, i.e. no prefix).
514  std::string path_str =
516 
517  // We have to build the json tree inside out to get all the values in there
518  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
519 
520  // Build the mesh data
521  std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
522  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
523  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
524  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
525  mesh["path"] = picojson::value(path_str + ((format==0)?"":"p") + "mesh" +
526  file_ext_format);
527  mesh["tags"] = picojson::value(mtags);
528 
529  // Build the fields data entries
530  for (FieldInfoMapIterator it = field_info_map.begin();
531  it != field_info_map.end(); ++it)
532  {
533  ftags["assoc"] = picojson::value((it->second).association);
534  ftags["comps"] = picojson::value(to_string((it->second).num_components));
535  field["path"] = picojson::value(path_str + it->first + file_ext_format);
536  field["tags"] = picojson::value(ftags);
537  fields[it->first] = picojson::value(field);
538  }
539 
540  main["cycle"] = picojson::value(double(cycle));
541  main["time"] = picojson::value(time);
542  main["domains"] = picojson::value(double(num_procs));
543  main["mesh"] = picojson::value(mesh);
544  if (!field_info_map.empty())
545  {
546  main["fields"] = picojson::value(fields);
547  }
548 
549  dsets["main"] = picojson::value(main);
550  top["dsets"] = picojson::value(dsets);
551 
552  return picojson::value(top).serialize(true);
553 }
554 
555 void VisItDataCollection::ParseVisItRootString(const std::string& json)
556 {
557  picojson::value top, dsets, main, mesh, fields;
558  std::string parse_err = picojson::parse(top, json);
559  if (!parse_err.empty())
560  {
561  error = READ_ERROR;
562  MFEM_WARNING("Unable to parse visit root data.");
563  return;
564  }
565 
566  // Process "main"
567  dsets = top.get("dsets");
568  main = dsets.get("main");
569  cycle = int(main.get("cycle").get<double>());
570  time = main.get("time").get<double>();
571  num_procs = int(main.get("domains").get<double>());
572  mesh = main.get("mesh");
573  fields = main.get("fields");
574 
575  // ... Process "mesh"
576 
577  // Set the DataCollection::name using the mesh path
578  std::string path = mesh.get("path").get<std::string>();
579  size_t right_sep = path.find('_');
580  if (right_sep == std::string::npos)
581  {
582  error = READ_ERROR;
583  MFEM_WARNING("Unable to parse visit root data.");
584  return;
585  }
586  name = path.substr(0, right_sep);
587 
588  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
589  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
591  to_int(mesh.get("tags").get("max_lods").get<std::string>());
592 
593  // ... Process "fields"
594  field_info_map.clear();
595  if (fields.is<picojson::object>())
596  {
597  picojson::object fields_obj = fields.get<picojson::object>();
598  for (picojson::object::iterator it = fields_obj.begin();
599  it != fields_obj.end(); ++it)
600  {
601  picojson::value tags = it->second.get("tags");
602  field_info_map[it->first] =
603  VisItFieldInfo(tags.get("assoc").get<std::string>(),
604  to_int(tags.get("comps").get<std::string>()));
605  }
606  }
607 }
608 
609 } // end namespace MFEM
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.
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1019
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int format
Output mesh format: 0 - serial format (default), 1 - parallel format.
int error
Error state.
bool appendRankToFileName
Append rank to any output file names.
void DeleteAll()
Delete data owned by the DataCollection including field information.
std::string to_string(int i)
Definition: text.hpp:48
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
int VectorDim() const
Definition: gridfunc.cpp:261
Helper class for VisIt visualization data.
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.
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.
double time
Physical time (for time-dependent simulations)
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
FieldMapType::const_iterator FieldMapConstIterator
int GetNRanks() const
Definition: pmesh.hpp:118
Mesh * mesh
The (common) mesh for the collected fields.
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
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.
std::string to_padded_string(int i, int digits)
Definition: text.hpp:60
int to_int(const std::string &str)
Definition: text.hpp:68
std::string GetFieldFileName(const std::string &field_name)
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_*.
int Dimension() const
Definition: mesh.hpp:641
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
int SpaceDimension() const
Definition: mesh.hpp:642
int main()
std::map< std::string, VisItFieldInfo > field_info_map
int GetMyRank() const
Definition: pmesh.hpp:119
FieldMapType::iterator FieldMapIterator
The fields and their names (used when saving)
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:117
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.
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.
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.
double time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Definition: pmesh.cpp:4465
Class for parallel meshes.
Definition: pmesh.hpp:29
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.