MFEM  v3.0
 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.googlecode.com.
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  pad_digits = pad_digits_default;
69  error = NO_ERROR;
70 }
71 
72 DataCollection::DataCollection(const char *collection_name, Mesh *_mesh)
73 {
74  name = collection_name;
75  mesh = _mesh;
76  myid = 0;
77  num_procs = 1;
78  serial = true;
79 #ifdef MFEM_USE_MPI
80  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
81  if (par_mesh)
82  {
83  myid = par_mesh->GetMyRank();
84  num_procs = par_mesh->GetNRanks();
85  serial = false;
86  }
87 #endif
88  own_data = false;
89  cycle = -1;
90  time = 0.0;
91  pad_digits = pad_digits_default;
92  error = NO_ERROR;
93 }
94 
95 void DataCollection::RegisterField(const char* name, GridFunction *gf)
96 {
97  field_map[name] = gf;
98 }
99 
100 GridFunction *DataCollection::GetField(const char *field_name)
101 {
102  if (HasField(field_name))
103  return field_map[field_name];
104  else
105  return NULL;
106 }
107 
109 {
110  string dir_name;
111  if (cycle == -1)
112  dir_name = name;
113  else
114  dir_name = name + "_" + to_padded_string(cycle, pad_digits);
115  int err;
116 #ifndef MFEM_USE_MPI
117  err = mkdir(dir_name.c_str(), 0777);
118  err = (err && (errno != EEXIST)) ? 1 : 0;
119 #else
120  ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
121  if (myid == 0 || pmesh == NULL)
122  {
123  err = mkdir(dir_name.c_str(), 0777);
124  err = (err && (errno != EEXIST)) ? 1 : 0;
125  if (pmesh)
126  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
127  }
128  else
129  {
130  // Wait for rank 0 to create the directory
131  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
132  }
133 #endif
134  if (err)
135  {
136  error = WRITE_ERROR;
137  MFEM_WARNING("Error creating directory: " << dir_name);
138  return; // do not even try to write mesh or fields
139  }
140 
141  string mesh_name;
142  if (serial)
143  mesh_name = dir_name + "/mesh";
144  else
145  mesh_name = dir_name + "/mesh." + to_padded_string(myid, pad_digits);
146  ofstream mesh_file(mesh_name.c_str());
147  mesh->Print(mesh_file);
148  if (!mesh_file)
149  {
150  error = WRITE_ERROR;
151  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
152  }
153  mesh_file.close();
154 
155  string field_name;
156  for (map<string,GridFunction*>::iterator it = field_map.begin();
157  it != field_map.end(); ++it)
158  {
159  if (serial)
160  field_name = dir_name + "/" + it->first;
161  else
162  field_name = dir_name + "/" + it->first + "." +
163  to_padded_string(myid, pad_digits);
164  ofstream field_file(field_name.c_str());
165  (it->second)->Save(field_file);
166  if (!field_file)
167  {
168  error = WRITE_ERROR;
169  MFEM_WARNING("Error writting field to file: " << field_name);
170  }
171  }
172 }
173 
175 {
176  if (own_data)
177  delete mesh;
178  mesh = NULL;
179  for (map<string,GridFunction*>::iterator it = field_map.begin();
180  it != field_map.end(); ++it)
181  {
182  if (own_data)
183  delete it->second;
184  it->second = NULL;
185  }
186  own_data = false;
187 }
188 
190 {
191  DeleteData();
192  field_map.clear();
193 }
194 
196 {
197  if (own_data)
198  {
199  delete mesh;
200  for (map<string,GridFunction*>::iterator it = field_map.begin();
201  it != field_map.end(); ++it)
202  delete it->second;
203  }
204 }
205 
206 
207 // class VisItDataCollection implementation
208 
209 VisItDataCollection::VisItDataCollection(const char *collection_name)
210  : DataCollection(collection_name)
211 {
212  serial = false; // always include rank in file names
213  cycle = 0; // always include cycle in directory names
214 
215  spatial_dim = 0;
216  topo_dim = 0;
218 }
219 
220 VisItDataCollection::VisItDataCollection(const char *collection_name,
221  Mesh *mesh)
222  : DataCollection(collection_name, mesh)
223 {
224  serial = false; // always include rank in file names
225  cycle = 0; // always include cycle in directory names
226 
227  spatial_dim = mesh->SpaceDimension();
228  topo_dim = mesh->Dimension();
230 }
231 
233 {
235  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
236 }
237 
238 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
239 {
240  visit_max_levels_of_detail = max_levels_of_detail;
241 }
242 
244 {
245  field_info_map.clear();
247 }
248 
250 {
251  if (myid == 0)
252  {
253  string root_name = name + "_" + to_padded_string(cycle, pad_digits) +
254  ".mfem_root";
255  ofstream root_file(root_name.c_str());
256  root_file << GetVisItRootString();
257  if (!root_file)
258  {
259  error = WRITE_ERROR;
260  MFEM_WARNING("Error writting VisIt Root file: " << root_name);
261  }
262  }
264 }
265 
267 {
268  DeleteAll();
269  cycle = _cycle;
270  string root_name = name + "_" + to_padded_string(cycle, pad_digits) +
271  ".mfem_root";
272  LoadVisItRootFile(root_name);
273  if (!error)
274  LoadMesh();
275  if (!error)
276  LoadFields();
277  if (!error)
278  own_data = true;
279  else
280  DeleteAll();
281 }
282 
284 {
285  ifstream root_file(root_name.c_str());
286  stringstream buffer;
287  buffer << root_file.rdbuf();
288  if (!buffer)
289  {
290  error = READ_ERROR;
291  MFEM_WARNING("Error reading the VisIt Root file: " << root_name);
292  }
293  else
294  {
295  ParseVisItRootString(buffer.str());
296  }
297 }
298 
300 {
301  string mesh_fname = name + "_" + to_padded_string(cycle, pad_digits) +
302  "/mesh." + to_padded_string(myid, pad_digits);
303  ifstream file(mesh_fname.c_str());
304  if (!file)
305  {
306  error = READ_ERROR;
307  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
308  return;
309  }
310  // TODO: 1) load parallel mesh on one processor
311  // 2) load parallel mesh on the same number of processors
312  mesh = new Mesh(file, 1, 1);
314  topo_dim = mesh->Dimension();
315 }
316 
318 {
319  string path_left = name + "_" + to_padded_string(cycle, pad_digits) + "/";
320  string path_right = "." + to_padded_string(myid, pad_digits);
321 
322  field_map.clear();
323  for (map<string,VisItFieldInfo>::iterator it = field_info_map.begin();
324  it != field_info_map.end(); ++it)
325  {
326  string fname = path_left + it->first + path_right;
327  ifstream file(fname.c_str());
328  if (!file)
329  {
330  error = READ_ERROR;
331  MFEM_WARNING("Unable to open field file: " << fname);
332  return;
333  }
334  // TODO: 1) load parallel GridFunction on one processor
335  // 2) load parallel GridFunction on the same number of processors
336  field_map[it->first] = new GridFunction(mesh, file);
337  }
338 }
339 
341 {
342  // Get the path string
343  string path_str = name + "_" + to_padded_string(cycle, pad_digits) + "/";
344 
345  // We have to build the json tree inside out to get all the values in there
346  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
347 
348  // Build the mesh data
349  string file_ext_format = ".%0" + to_string(pad_digits) + "d";
350  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
351  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
352  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
353  mesh["path"] = picojson::value(path_str + "mesh" + file_ext_format);
354  mesh["tags"] = picojson::value(mtags);
355 
356  // Build the fields data entries
357  for (map<string,VisItFieldInfo>::iterator it = field_info_map.begin();
358  it != field_info_map.end(); ++it)
359  {
360  ftags["assoc"] = picojson::value((it->second).association);
361  ftags["comps"] = picojson::value(to_string((it->second).num_components));
362  field["path"] = picojson::value(path_str + it->first + file_ext_format);
363  field["tags"] = picojson::value(ftags);
364  fields[it->first] = picojson::value(field);
365  }
366 
367  main["cycle"] = picojson::value(double(cycle));
368  main["time"] = picojson::value(time);
369  main["domains"] = picojson::value(double(num_procs));
370  main["mesh"] = picojson::value(mesh);
371  if (!field_info_map.empty())
372  main["fields"] = picojson::value(fields);
373 
374  dsets["main"] = picojson::value(main);
375  top["dsets"] = picojson::value(dsets);
376 
377  return picojson::value(top).serialize(true);
378 }
379 
381 {
382  picojson::value top, dsets, main, mesh, fields;
383  string parse_err = picojson::parse(top, json);
384  if (!parse_err.empty())
385  {
386  error = READ_ERROR;
387  MFEM_WARNING("Unable to parse visit root data.");
388  return;
389  }
390 
391  // Process "main"
392  dsets = top.get("dsets");
393  main = dsets.get("main");
394  cycle = int(main.get("cycle").get<double>());
395  time = main.get("time").get<double>();
396  num_procs = int(main.get("domains").get<double>());
397  mesh = main.get("mesh");
398  fields = main.get("fields");
399 
400  // ... Process "mesh"
401 
402  // Set the DataCollection::name using the mesh path
403  string path = mesh.get("path").get<string>();
404  size_t right_sep = path.find('_');
405  if (right_sep == string::npos)
406  {
407  error = READ_ERROR;
408  MFEM_WARNING("Unable to parse visit root data.");
409  return;
410  }
411  name = path.substr(0, right_sep);
412 
413  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<string>());
414  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<string>());
416  to_int(mesh.get("tags").get("max_lods").get<string>());
417 
418  // ... Process "fields"
419  field_info_map.clear();
420  if (fields.is<picojson::object>())
421  {
422  picojson::object fields_obj = fields.get<picojson::object>();
423  for (picojson::object::iterator it = fields_obj.begin();
424  it != fields_obj.end(); ++it)
425  {
426  picojson::value tags = it->second.get("tags");
427  field_info_map[it->first] =
428  VisItFieldInfo(tags.get("assoc").get<string>(),
429  to_int(tags.get("comps").get<string>()));
430  }
431  }
432 }
433 
434 }
int to_int(string str)
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:26
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.
MPI_Comm GetComm()
Definition: pmesh.hpp:68
int VectorDim() const
Definition: gridfunc.cpp:160
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)
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)
Get a pointer to a grid function in the collection.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
void Load(int _cycle=0)
Load the collection based on its VisIt data (described in its root file)
int Dimension() const
Definition: mesh.hpp:417
int SpaceDimension() const
Definition: mesh.hpp:418
std::map< std::string, VisItFieldInfo > field_info_map
void LoadVisItRootFile(std::string root_name)
int main(int argc, char *argv[])
Definition: ex1.cpp:39
int GetNRanks()
Definition: pmesh.hpp:69
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.
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.
std::string name
Name of the collection, used as a directory name when saving.
string to_string(int i)
int num_procs
Number of MPI ranks (in parallel)
int GetMyRank()
Definition: pmesh.hpp:70
Class for parallel meshes.
Definition: pmesh.hpp:27
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.