MFEM  v3.4
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  std::string::size_type pos = collection_name.find_last_of('/');
75  if (pos == std::string::npos)
76  {
77  name = collection_name;
78  // leave prefix_path empty
79  }
80  else
81  {
82  prefix_path = collection_name.substr(0, pos+1);
83  name = collection_name.substr(pos+1);
84  }
85  mesh = mesh_;
86  myid = 0;
87  num_procs = 1;
88  serial = true;
89  appendRankToFileName = false;
90 
91 #ifdef MFEM_USE_MPI
92  m_comm = MPI_COMM_NULL;
93  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
94  if (par_mesh)
95  {
96  myid = par_mesh->GetMyRank();
97  num_procs = par_mesh->GetNRanks();
98  m_comm = par_mesh->GetComm();
99  serial = false;
100  appendRankToFileName = true;
101  }
102 #endif
103  own_data = false;
104  cycle = -1;
105  time = 0.0;
106  time_step = 0.0;
109  format = SERIAL_FORMAT; // use serial mesh format
110  error = NO_ERROR;
111 }
112 
114 {
115  if (own_data && new_mesh != mesh) { delete mesh; }
116  mesh = new_mesh;
117  myid = 0;
118  num_procs = 1;
119  serial = true;
120  appendRankToFileName = false;
121 
122 #ifdef MFEM_USE_MPI
123  m_comm = MPI_COMM_NULL;
124  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
125  if (par_mesh)
126  {
127  myid = par_mesh->GetMyRank();
128  num_procs = par_mesh->GetNRanks();
129  m_comm = par_mesh->GetComm();
130  serial = false;
131  appendRankToFileName = true;
132  }
133 #endif
134 }
135 
136 #ifdef MFEM_USE_MPI
137 void DataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
138 {
139  // This seems to be the cleanest way to accomplish this
140  // and avoid duplicating fine grained details:
141 
142  SetMesh(new_mesh);
143 
144  m_comm = comm;
145  MPI_Comm_rank(comm, &myid);
146  MPI_Comm_size(comm, &num_procs);
147 }
148 #endif
149 
151 {
152  switch (fmt)
153  {
154  case SERIAL_FORMAT: break;
155 #ifdef MFEM_USE_MPI
156  case PARALLEL_FORMAT: break;
157 #endif
158  default: MFEM_ABORT("unknown format: " << fmt);
159  }
160  format = fmt;
161 }
162 
163 void DataCollection::SetPrefixPath(const std::string& prefix)
164 {
165  if (!prefix.empty())
166  {
167  prefix_path = prefix;
168  if (!prefix_path.empty() && prefix_path[prefix_path.size()-1] != '/')
169  {
170  prefix_path += '/';
171  }
172  }
173  else
174  {
175  prefix_path.clear();
176  }
177 }
178 
179 void DataCollection::Load(int cycle)
180 {
181  MFEM_ABORT("this method is not implemented");
182 }
183 
185 {
186  SaveMesh();
187 
188  if (error) { return; }
189 
190  for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
191  {
192  SaveOneField(it);
193  // Even if there is an error, try saving the other fields
194  }
195 
196  for (QFieldMapIterator it = q_field_map.begin(); it != q_field_map.end();
197  ++it)
198  {
199  SaveOneQField(it);
200  }
201 }
202 
204 {
205  int err;
206 
207  std::string dir_name = prefix_path + name;
208  if (cycle != -1)
209  {
210  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
211  }
212  err = create_directory(dir_name, mesh, myid);
213  if (err)
214  {
215  error = WRITE_ERROR;
216  MFEM_WARNING("Error creating directory: " << dir_name);
217  return; // do not even try to write the mesh
218  }
219 
220  std::string mesh_name = GetMeshFileName();
221  std::ofstream mesh_file(mesh_name.c_str());
222  mesh_file.precision(precision);
223 #ifdef MFEM_USE_MPI
224  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
225  if (pmesh && format == PARALLEL_FORMAT)
226  {
227  pmesh->ParPrint(mesh_file);
228  }
229  else
230 #endif
231  {
232  mesh->Print(mesh_file);
233  }
234  if (!mesh_file)
235  {
236  error = WRITE_ERROR;
237  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
238  }
239 }
240 
242 {
243  return (serial || format == SERIAL_FORMAT) ? "mesh" : "pmesh";
244 }
245 
247 {
249 }
250 
251 std::string DataCollection::GetFieldFileName(const std::string &field_name)
252 const
253 {
254  std::string dir_name = prefix_path + name;
255  if (cycle != -1)
256  {
257  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
258  }
259  std::string file_name = dir_name + "/" + field_name;
261  {
262  file_name += "." + to_padded_string(myid, pad_digits_rank);
263  }
264  return file_name;
265 }
266 
268 {
269  std::ofstream field_file(GetFieldFileName(it->first).c_str());
270  field_file.precision(precision);
271  (it->second)->Save(field_file);
272  if (!field_file)
273  {
274  error = WRITE_ERROR;
275  MFEM_WARNING("Error writing field to file: " << it->first);
276  }
277 }
278 
280 {
281  std::ofstream q_field_file(GetFieldFileName(it->first).c_str());
282  q_field_file.precision(precision);
283  (it->second)->Save(q_field_file);
284  if (!q_field_file)
285  {
286  error = WRITE_ERROR;
287  MFEM_WARNING("Error writing q-field to file: " << it->first);
288  }
289 }
290 
291 void DataCollection::SaveField(const std::string &field_name)
292 {
293  FieldMapIterator it = field_map.find(field_name);
294  if (it != field_map.end())
295  {
296  SaveOneField(it);
297  }
298 }
299 
300 void DataCollection::SaveQField(const std::string &q_field_name)
301 {
302  QFieldMapIterator it = q_field_map.find(q_field_name);
303  if (it != q_field_map.end())
304  {
305  SaveOneQField(it);
306  }
307 }
308 
310 {
311  if (own_data) { delete mesh; }
312  mesh = NULL;
313 
316  own_data = false;
317 }
318 
320 {
321  DeleteData();
322  field_map.clear();
323  q_field_map.clear();
324 }
325 
327 {
328  DeleteData();
329 }
330 
331 
332 // class VisItDataCollection implementation
333 
334 VisItDataCollection::VisItDataCollection(const std::string& collection_name,
335  Mesh *mesh)
336  : DataCollection(collection_name, mesh)
337 {
338  appendRankToFileName = true; // always include rank in file names
339  cycle = 0; // always include cycle in directory names
340 
341  if (mesh)
342  {
343  spatial_dim = mesh->SpaceDimension();
344  topo_dim = mesh->Dimension();
345  }
346  else
347  {
348  spatial_dim = 0;
349  topo_dim = 0;
350  }
352 }
353 
354 #ifdef MFEM_USE_MPI
356  const std::string& collection_name,
357  Mesh *mesh)
358  : DataCollection(collection_name, mesh)
359 {
360  m_comm = comm;
361  MPI_Comm_rank(comm, &myid);
362  MPI_Comm_size(comm, &num_procs);
363  appendRankToFileName = true; // always include rank in file names
364  cycle = 0; // always include cycle in directory names
365  spatial_dim = 0;
366  topo_dim = 0;
368 }
369 #endif
370 
372 {
373  DataCollection::SetMesh(new_mesh);
374  appendRankToFileName = true;
376  topo_dim = mesh->Dimension();
377 }
378 
379 #ifdef MFEM_USE_MPI
380 void VisItDataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
381 {
382  // use VisItDataCollection's custom SetMesh, then set MPI info
383  SetMesh(new_mesh);
384  m_comm = comm;
385  MPI_Comm_rank(comm, &myid);
386  MPI_Comm_size(comm, &num_procs);
387 }
388 #endif
389 
390 void VisItDataCollection::RegisterField(const std::string& name,
391  GridFunction *gf)
392 {
394  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
395 }
396 
397 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
398 {
399  visit_max_levels_of_detail = max_levels_of_detail;
400 }
401 
403 {
404  field_info_map.clear();
406 }
407 
409 {
411  SaveRootFile();
412 }
413 
415 {
416  if (myid != 0) { return; }
417 
418  std::string root_name = prefix_path + name + "_" +
420  ".mfem_root";
421  std::ofstream root_file(root_name.c_str());
422  root_file << GetVisItRootString();
423  if (!root_file)
424  {
425  error = WRITE_ERROR;
426  MFEM_WARNING("Error writing VisIt root file: " << root_name);
427  }
428 }
429 
431 {
432  DeleteAll();
433  time_step = 0.0;
434  error = NO_ERROR;
435  cycle = cycle_;
436  std::string root_name = prefix_path + name + "_" +
438  ".mfem_root";
439  LoadVisItRootFile(root_name);
440  if (format != SERIAL_FORMAT || num_procs > 1)
441  {
442 #ifndef MFEM_USE_MPI
443  MFEM_WARNING("Cannot load parallel VisIt root file in serial.");
444  error = READ_ERROR;
445 #else
446  if (m_comm == MPI_COMM_NULL)
447  {
448  MFEM_WARNING("Cannot load parallel VisIt root file without MPI"
449  " communicator");
450  error = READ_ERROR;
451  }
452  else
453  {
454  // num_procs was read from the root file, check for consistency with
455  // the associated MPI_Comm, m_comm:
456  int comm_size;
457  MPI_Comm_size(m_comm, &comm_size);
458  if (comm_size != num_procs)
459  {
460  MFEM_WARNING("Processor number mismatch: VisIt root file: "
461  << num_procs << ", MPI_comm: " << comm_size);
462  error = READ_ERROR;
463  }
464  else
465  {
466  // myid was set when setting m_comm
467  }
468  }
469 #endif
470  }
471  if (!error)
472  {
473  LoadMesh(); // sets own_data to true, when there is no error
474  }
475  if (!error)
476  {
477  LoadFields();
478  }
479  if (error)
480  {
481  DeleteAll();
482  }
483 }
484 
485 void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
486 {
487  std::ifstream root_file(root_name.c_str());
488  std::stringstream buffer;
489  buffer << root_file.rdbuf();
490  if (!buffer)
491  {
492  error = READ_ERROR;
493  MFEM_WARNING("Error reading the VisIt root file: " << root_name);
494  }
495  else
496  {
497  ParseVisItRootString(buffer.str());
498  }
499 }
500 
502 {
503  std::string mesh_fname = GetMeshFileName();
504  named_ifgzstream file(mesh_fname.c_str());
505  // TODO: in parallel, check for errors on all processors
506  if (!file)
507  {
508  error = READ_ERROR;
509  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
510  return;
511  }
512  // TODO: 1) load parallel mesh on one processor
513  if (format == SERIAL_FORMAT)
514  {
515  mesh = new Mesh(file, 1, 0, false);
516  serial = true;
517  }
518  else
519  {
520 #ifdef MFEM_USE_MPI
521  mesh = new ParMesh(m_comm, file);
522  serial = false;
523 #else
524  error = READ_ERROR;
525  MFEM_WARNING("Reading parallel format in serial is not supported");
526  return;
527 #endif
528  }
530  topo_dim = mesh->Dimension();
531  own_data = true;
532 }
533 
535 {
536  std::string path_left = prefix_path + name + "_" +
538  std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
539 
540  field_map.clear();
541  for (FieldInfoMapIterator it = field_info_map.begin();
542  it != field_info_map.end(); ++it)
543  {
544  std::string fname = path_left + it->first + path_right;
545  std::ifstream file(fname.c_str());
546  // TODO: in parallel, check for errors on all processors
547  if (!file)
548  {
549  error = READ_ERROR;
550  MFEM_WARNING("Unable to open field file: " << fname);
551  return;
552  }
553  // TODO: 1) load parallel GridFunction on one processor
554  if (serial)
555  {
556  field_map.Register(it->first, new GridFunction(mesh, file), own_data);
557  }
558  else
559  {
560 #ifdef MFEM_USE_MPI
562  it->first,
563  new ParGridFunction(dynamic_cast<ParMesh*>(mesh), file), own_data);
564 #else
565  error = READ_ERROR;
566  MFEM_WARNING("Reading parallel format in serial is not supported");
567  return;
568 #endif
569  }
570  }
571 }
572 
574 {
575  // Get the path string (relative to where the root file is, i.e. no prefix).
576  std::string path_str =
578 
579  // We have to build the json tree inside out to get all the values in there
580  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
581 
582  // Build the mesh data
583  std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
584  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
585  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
586  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
587  mesh["path"] = picojson::value(path_str + GetMeshShortFileName() +
588  file_ext_format);
589  mesh["tags"] = picojson::value(mtags);
590  mesh["format"] = picojson::value(to_string(format));
591 
592  // Build the fields data entries
593  for (FieldInfoMapIterator it = field_info_map.begin();
594  it != field_info_map.end(); ++it)
595  {
596  ftags["assoc"] = picojson::value((it->second).association);
597  ftags["comps"] = picojson::value(to_string((it->second).num_components));
598  field["path"] = picojson::value(path_str + it->first + file_ext_format);
599  field["tags"] = picojson::value(ftags);
600  fields[it->first] = picojson::value(field);
601  }
602 
603  main["cycle"] = picojson::value(double(cycle));
604  main["time"] = picojson::value(time);
605  main["time_step"] = picojson::value(time_step);
606  main["domains"] = picojson::value(double(num_procs));
607  main["mesh"] = picojson::value(mesh);
608  if (!field_info_map.empty())
609  {
610  main["fields"] = picojson::value(fields);
611  }
612 
613  dsets["main"] = picojson::value(main);
614  top["dsets"] = picojson::value(dsets);
615 
616  return picojson::value(top).serialize(true);
617 }
618 
619 void VisItDataCollection::ParseVisItRootString(const std::string& json)
620 {
621  picojson::value top, dsets, main, mesh, fields;
622  std::string parse_err = picojson::parse(top, json);
623  if (!parse_err.empty())
624  {
625  error = READ_ERROR;
626  MFEM_WARNING("Unable to parse VisIt root data.");
627  return;
628  }
629 
630  // Process "main"
631  dsets = top.get("dsets");
632  main = dsets.get("main");
633  cycle = int(main.get("cycle").get<double>());
634  time = main.get("time").get<double>();
635  if (main.contains("time_step"))
636  {
637  time_step = main.get("time_step").get<double>();
638  }
639  num_procs = int(main.get("domains").get<double>());
640  mesh = main.get("mesh");
641  fields = main.get("fields");
642 
643  // ... Process "mesh"
644 
645  // Set the DataCollection::name using the mesh path
646  std::string path = mesh.get("path").get<std::string>();
647  size_t right_sep = path.find('_');
648  if (right_sep == std::string::npos)
649  {
650  error = READ_ERROR;
651  MFEM_WARNING("Unable to parse VisIt root data.");
652  return;
653  }
654  name = path.substr(0, right_sep);
655 
656  if (mesh.contains("format"))
657  {
658  format = to_int(mesh.get("format").get<std::string>());
659  }
660  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
661  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
663  to_int(mesh.get("tags").get("max_lods").get<std::string>());
664 
665  // ... Process "fields"
666  field_info_map.clear();
667  if (fields.is<picojson::object>())
668  {
669  picojson::object fields_obj = fields.get<picojson::object>();
670  for (picojson::object::iterator it = fields_obj.begin();
671  it != fields_obj.end(); ++it)
672  {
673  picojson::value tags = it->second.get("tags");
674  field_info_map[it->first] =
675  VisItFieldInfo(tags.get("assoc").get<std::string>(),
676  to_int(tags.get("comps").get<std::string>()));
677  }
678  }
679 }
680 
681 } // 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:1034
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int format
Output mesh format: see the Format enumeration.
int error
Error state.
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
void DeleteAll()
Delete data owned by the DataCollection including field information.
std::string to_string(int i)
Definition: text.hpp:50
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
int VectorDim() const
Definition: gridfunc.cpp:291
Helper class for VisIt visualization data.
virtual void Save()
Save the collection and a VisIt root file.
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.
int main(int argc, char *argv[])
Definition: ex1.cpp:45
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.
MPI_Comm m_comm
Associated MPI communicator.
double time
Physical time (for time-dependent simulations)
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
int GetNRanks() const
Definition: pmesh.hpp:124
iterator find(const std::string &fname)
Returns an iterator to the field fname.
std::string GetMeshShortFileName() const
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.
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 to_padded_string(int i, int digits)
Definition: text.hpp:62
int to_int(const std::string &str)
Definition: text.hpp:70
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.
static const int pad_digits_default
Default value for pad_digits_*.
int Dimension() const
Definition: mesh.hpp:645
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:646
std::map< std::string, VisItFieldInfo > field_info_map
int GetMyRank() const
Definition: pmesh.hpp:125
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:123
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.
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.
std::string GetFieldFileName(const std::string &field_name) const
void clear()
Clears the map of registered fields without reclaiming memory.
void SaveRootFile()
Save a VisIt root file for the collection.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has &#39;/&#39; at the end...
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
std::string GetMeshFileName() const
void LoadVisItRootFile(const std::string &root_name)
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)
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Definition: pmesh.cpp:4494
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.