MFEM  v3.1
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 "picojson.h"
14 
15 #include <fstream>
16 #include <cerrno> // errno
17 #ifndef _WIN32
18 #include <sys/stat.h> // mkdir
19 #else
20 #include <direct.h> // _mkdir
21 #define mkdir(dir, mode) _mkdir(dir)
22 #endif
23 
24 namespace mfem
25 {
26 
27 using namespace std;
28 
29 // Helper string functions. Will go away in C++11
30 
31 string to_string(int i)
32 {
33  stringstream ss;
34  ss << i;
35 
36  // trim leading spaces
37  string out_str = ss.str();
38  out_str = out_str.substr(out_str.find_first_not_of(" \t"));
39  return out_str;
40 }
41 
42 string to_padded_string(int i, int digits)
43 {
44  ostringstream oss;
45  oss << setw(digits) << setfill('0') << i;
46  return oss.str();
47 }
48 
49 int to_int(string str)
50 {
51  int i;
52  stringstream(str) >> i;
53  return i;
54 }
55 
56 // class DataCollection implementation
57 
58 DataCollection::DataCollection(const char *collection_name)
59 {
60  name = collection_name;
61  mesh = NULL;
62  myid = 0;
63  num_procs = 1;
64  serial = true;
65  own_data = false;
66  cycle = -1;
67  time = 0.0;
68  precision = precision_default;
69  pad_digits = pad_digits_default;
70  error = NO_ERROR;
71 }
72 
73 DataCollection::DataCollection(const char *collection_name, Mesh *_mesh)
74 {
75  name = collection_name;
76  mesh = _mesh;
77  myid = 0;
78  num_procs = 1;
79  serial = true;
80 #ifdef MFEM_USE_MPI
81  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
82  if (par_mesh)
83  {
84  myid = par_mesh->GetMyRank();
85  num_procs = par_mesh->GetNRanks();
86  serial = false;
87  }
88 #endif
89  own_data = false;
90  cycle = -1;
91  time = 0.0;
92  precision = precision_default;
93  pad_digits = pad_digits_default;
94  error = NO_ERROR;
95 }
96 
98 {
99  if (own_data) { delete mesh; }
100  mesh = new_mesh;
101  myid = 0;
102  num_procs = 1;
103  serial = true;
104 #ifdef MFEM_USE_MPI
105  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
106  if (par_mesh)
107  {
108  myid = par_mesh->GetMyRank();
109  num_procs = par_mesh->GetNRanks();
110  serial = false;
111  }
112 #endif
113 }
114 
115 void DataCollection::RegisterField(const char* name, GridFunction *gf)
116 {
117  if (own_data && HasField(name))
118  {
119  delete field_map[name];
120  }
121  field_map[name] = gf;
122 }
123 
124 GridFunction *DataCollection::GetField(const char *field_name)
125 {
126  if (HasField(field_name))
127  {
128  return field_map[field_name];
129  }
130  else
131  {
132  return NULL;
133  }
134 }
135 
137 {
138  SaveMesh();
139 
140  if (error) { return; }
141 
142  for (map<string,GridFunction*>::iterator it = field_map.begin();
143  it != field_map.end(); ++it)
144  {
145  SaveOneField(it);
146  }
147 }
148 
150 {
151  string dir_name;
152  if (cycle == -1)
153  {
154  dir_name = name;
155  }
156  else
157  {
158  dir_name = name + "_" + to_padded_string(cycle, pad_digits);
159  }
160  int err;
161 #ifndef MFEM_USE_MPI
162  err = mkdir(dir_name.c_str(), 0777);
163  err = (err && (errno != EEXIST)) ? 1 : 0;
164 #else
165  ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
166  if (myid == 0 || pmesh == NULL)
167  {
168  err = mkdir(dir_name.c_str(), 0777);
169  err = (err && (errno != EEXIST)) ? 1 : 0;
170  if (pmesh)
171  {
172  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
173  }
174  }
175  else
176  {
177  // Wait for rank 0 to create the directory
178  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
179  }
180 #endif
181  if (err)
182  {
183  error = WRITE_ERROR;
184  MFEM_WARNING("Error creating directory: " << dir_name);
185  return; // do not even try to write the mesh
186  }
187 
188  string mesh_name;
189  if (serial)
190  {
191  mesh_name = dir_name + "/mesh";
192  }
193  else
194  {
195  mesh_name = dir_name + "/mesh." + to_padded_string(myid, pad_digits);
196  }
197  ofstream mesh_file(mesh_name.c_str());
198  mesh_file.precision(precision);
199  mesh->Print(mesh_file);
200  if (!mesh_file)
201  {
202  error = WRITE_ERROR;
203  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
204  }
205 }
206 
208  const std::map<std::string,GridFunction*>::iterator &it)
209 {
210  string dir_name;
211  if (cycle == -1)
212  {
213  dir_name = name;
214  }
215  else
216  {
217  dir_name = name + "_" + to_padded_string(cycle, pad_digits);
218  }
219 
220  string file_name;
221  if (serial)
222  {
223  file_name = dir_name + "/" + it->first;
224  }
225  else
226  file_name = dir_name + "/" + it->first + "." +
227  to_padded_string(myid, pad_digits);
228  ofstream field_file(file_name.c_str());
229  field_file.precision(precision);
230  (it->second)->Save(field_file);
231  if (!field_file)
232  {
233  error = WRITE_ERROR;
234  MFEM_WARNING("Error writting field to file: " << it->first);
235  }
236 }
237 
238 void DataCollection::SaveField(const char *field_name)
239 {
240  const map<string,GridFunction*>::iterator it = field_map.find(field_name);
241  if (it != field_map.end())
242  {
243  SaveOneField(it);
244  }
245 }
246 
248 {
249  if (own_data)
250  {
251  delete mesh;
252  }
253  mesh = NULL;
254  for (map<string,GridFunction*>::iterator it = field_map.begin();
255  it != field_map.end(); ++it)
256  {
257  if (own_data)
258  {
259  delete it->second;
260  }
261  it->second = NULL;
262  }
263  own_data = false;
264 }
265 
267 {
268  DeleteData();
269  field_map.clear();
270 }
271 
273 {
274  if (own_data)
275  {
276  delete mesh;
277  for (map<string,GridFunction*>::iterator it = field_map.begin();
278  it != field_map.end(); ++it)
279  {
280  delete it->second;
281  }
282  }
283 }
284 
285 
286 // class VisItDataCollection implementation
287 
288 VisItDataCollection::VisItDataCollection(const char *collection_name)
289  : DataCollection(collection_name)
290 {
291  serial = false; // always include rank in file names
292  cycle = 0; // always include cycle in directory names
293 
294  spatial_dim = 0;
295  topo_dim = 0;
297 }
298 
299 VisItDataCollection::VisItDataCollection(const char *collection_name,
300  Mesh *mesh)
301  : DataCollection(collection_name, mesh)
302 {
303  serial = false; // always include rank in file names
304  cycle = 0; // always include cycle in directory names
305 
306  spatial_dim = mesh->SpaceDimension();
307  topo_dim = mesh->Dimension();
309 }
310 
312 {
313  DataCollection::SetMesh(new_mesh);
314  serial = false; // always include rank in file names
315 
317  topo_dim = mesh->Dimension();
318 }
319 
321 {
323  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
324 }
325 
326 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
327 {
328  visit_max_levels_of_detail = max_levels_of_detail;
329 }
330 
332 {
333  field_info_map.clear();
335 }
336 
338 {
340  SaveRootFile();
341 }
342 
344 {
345  if (myid == 0)
346  {
347  string root_name = name + "_" + to_padded_string(cycle, pad_digits) +
348  ".mfem_root";
349  ofstream root_file(root_name.c_str());
350  root_file << GetVisItRootString();
351  if (!root_file)
352  {
353  error = WRITE_ERROR;
354  MFEM_WARNING("Error writting VisIt Root file: " << root_name);
355  }
356  }
357 }
358 
360 {
361  DeleteAll();
362  cycle = _cycle;
363  string root_name = name + "_" + to_padded_string(cycle, pad_digits) +
364  ".mfem_root";
365  LoadVisItRootFile(root_name);
366  if (!error)
367  {
368  LoadMesh();
369  }
370  if (!error)
371  {
372  LoadFields();
373  }
374  if (!error)
375  {
376  own_data = true;
377  }
378  else
379  {
380  DeleteAll();
381  }
382 }
383 
385 {
386  ifstream root_file(root_name.c_str());
387  stringstream buffer;
388  buffer << root_file.rdbuf();
389  if (!buffer)
390  {
391  error = READ_ERROR;
392  MFEM_WARNING("Error reading the VisIt Root file: " << root_name);
393  }
394  else
395  {
396  ParseVisItRootString(buffer.str());
397  }
398 }
399 
401 {
402  string mesh_fname = name + "_" + to_padded_string(cycle, pad_digits) +
403  "/mesh." + to_padded_string(myid, pad_digits);
404  ifstream file(mesh_fname.c_str());
405  if (!file)
406  {
407  error = READ_ERROR;
408  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
409  return;
410  }
411  // TODO: 1) load parallel mesh on one processor
412  // 2) load parallel mesh on the same number of processors
413  mesh = new Mesh(file, 1, 1);
415  topo_dim = mesh->Dimension();
416 }
417 
419 {
420  string path_left = name + "_" + to_padded_string(cycle, pad_digits) + "/";
421  string path_right = "." + to_padded_string(myid, pad_digits);
422 
423  field_map.clear();
424  for (map<string,VisItFieldInfo>::iterator it = field_info_map.begin();
425  it != field_info_map.end(); ++it)
426  {
427  string fname = path_left + it->first + path_right;
428  ifstream file(fname.c_str());
429  if (!file)
430  {
431  error = READ_ERROR;
432  MFEM_WARNING("Unable to open field file: " << fname);
433  return;
434  }
435  // TODO: 1) load parallel GridFunction on one processor
436  // 2) load parallel GridFunction on the same number of processors
437  field_map[it->first] = new GridFunction(mesh, file);
438  }
439 }
440 
442 {
443  // Get the path string
444  string path_str = name + "_" + to_padded_string(cycle, pad_digits) + "/";
445 
446  // We have to build the json tree inside out to get all the values in there
447  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
448 
449  // Build the mesh data
450  string file_ext_format = ".%0" + to_string(pad_digits) + "d";
451  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
452  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
453  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
454  mesh["path"] = picojson::value(path_str + "mesh" + file_ext_format);
455  mesh["tags"] = picojson::value(mtags);
456 
457  // Build the fields data entries
458  for (map<string,VisItFieldInfo>::iterator it = field_info_map.begin();
459  it != field_info_map.end(); ++it)
460  {
461  ftags["assoc"] = picojson::value((it->second).association);
462  ftags["comps"] = picojson::value(to_string((it->second).num_components));
463  field["path"] = picojson::value(path_str + it->first + file_ext_format);
464  field["tags"] = picojson::value(ftags);
465  fields[it->first] = picojson::value(field);
466  }
467 
468  main["cycle"] = picojson::value(double(cycle));
469  main["time"] = picojson::value(time);
470  main["domains"] = picojson::value(double(num_procs));
471  main["mesh"] = picojson::value(mesh);
472  if (!field_info_map.empty())
473  {
474  main["fields"] = picojson::value(fields);
475  }
476 
477  dsets["main"] = picojson::value(main);
478  top["dsets"] = picojson::value(dsets);
479 
480  return picojson::value(top).serialize(true);
481 }
482 
484 {
485  picojson::value top, dsets, main, mesh, fields;
486  string parse_err = picojson::parse(top, json);
487  if (!parse_err.empty())
488  {
489  error = READ_ERROR;
490  MFEM_WARNING("Unable to parse visit root data.");
491  return;
492  }
493 
494  // Process "main"
495  dsets = top.get("dsets");
496  main = dsets.get("main");
497  cycle = int(main.get("cycle").get<double>());
498  time = main.get("time").get<double>();
499  num_procs = int(main.get("domains").get<double>());
500  mesh = main.get("mesh");
501  fields = main.get("fields");
502 
503  // ... Process "mesh"
504 
505  // Set the DataCollection::name using the mesh path
506  string path = mesh.get("path").get<string>();
507  size_t right_sep = path.find('_');
508  if (right_sep == string::npos)
509  {
510  error = READ_ERROR;
511  MFEM_WARNING("Unable to parse visit root data.");
512  return;
513  }
514  name = path.substr(0, right_sep);
515 
516  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<string>());
517  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<string>());
519  to_int(mesh.get("tags").get("max_lods").get<string>());
520 
521  // ... Process "fields"
522  field_info_map.clear();
523  if (fields.is<picojson::object>())
524  {
525  picojson::object fields_obj = fields.get<picojson::object>();
526  for (picojson::object::iterator it = fields_obj.begin();
527  it != fields_obj.end(); ++it)
528  {
529  picojson::value tags = it->second.get("tags");
530  field_info_map[it->first] =
531  VisItFieldInfo(tags.get("assoc").get<string>(),
532  to_int(tags.get("comps").get<string>()));
533  }
534  }
535 }
536 
537 }
int to_int(string str)
virtual void SaveMesh()
Save the mesh, creating the collection directory.
std::map< std::string, GridFunction * > field_map
The fields and their names (used when saving)
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int error
Error state.
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
string to_padded_string(int i, int digits)
void DeleteAll()
Delete data owned by the DataCollection including field information.
int VectorDim() const
Definition: gridfunc.cpp:224
Helper class for VisIt visualization data.
virtual void Save()
Save the collection and a VisIt root file.
bool own_data
Should the collection delete its mesh and fields.
double time
Physical time (for time-dependent simulations)
int GetNRanks() const
Definition: pmesh.hpp:87
Mesh * mesh
The (common) mesh for the collected fields.
bool serial
Serial or parallel run? If false, append rank (myid) to file names.
void ParseVisItRootString(std::string json)
Read in a VisIt root file in JSON format.
GridFunction * GetField(const char *field_name)
void DeleteData()
Delete data owned by the DataCollection keeping field information.
virtual void SaveField(const char *field_name)
Save one field, assuming the collection directory already exists.
void Load(int _cycle=0)
Load the collection based on its VisIt data (described in its root file)
int Dimension() const
Definition: mesh.hpp:475
int SpaceDimension() const
Definition: mesh.hpp:476
std::map< std::string, VisItFieldInfo > field_info_map
void LoadVisItRootFile(std::string root_name)
int main(int argc, char *argv[])
Definition: ex1.cpp:45
int GetMyRank() const
Definition: pmesh.hpp:88
MPI_Comm GetComm() const
Definition: pmesh.hpp:86
int myid
MPI rank (in parallel)
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 SaveRootFile()
Save a VisIt root file for the collection.
VisItDataCollection(const char *collection_name)
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection.
DataCollection(const char *collection_name)
Create an empty collection with the given name.
int pad_digits
Number of digits used for the cycle and MPI rank in filenames.
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 SaveOneField(const std::map< std::string, GridFunction * >::iterator &it)
Save one field to disk, assuming the collection directory exists.
string to_string(int i)
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
Definition: pmesh.hpp:28
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.