MFEM  v4.0
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 "../mesh/nurbs.hpp"
14 #include "../general/text.hpp"
15 #include "picojson.h"
16 
17 #include <fstream>
18 #include <cerrno> // errno
19 #include <sstream>
20 
21 #ifndef _WIN32
22 #include <sys/stat.h> // mkdir
23 #else
24 #include <direct.h> // _mkdir
25 #define mkdir(dir, mode) _mkdir(dir)
26 #endif
27 
28 namespace mfem
29 {
30 
31 // static method
32 int DataCollection::create_directory(const std::string &dir_name,
33  const Mesh *mesh, int myid)
34 {
35  // create directories recursively
36  const char path_delim = '/';
37  std::string::size_type pos = 0;
38  int err;
39 #ifdef MFEM_USE_MPI
40  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
41 #endif
42 
43  do
44  {
45  pos = dir_name.find(path_delim, pos+1);
46  std::string subdir = dir_name.substr(0, pos);
47 
48 #ifndef MFEM_USE_MPI
49  err = mkdir(subdir.c_str(), 0777);
50  err = (err && (errno != EEXIST)) ? 1 : 0;
51 #else
52  if (myid == 0 || pmesh == NULL)
53  {
54  err = mkdir(subdir.c_str(), 0777);
55  err = (err && (errno != EEXIST)) ? 1 : 0;
56  }
57 #endif
58  }
59  while ( pos != std::string::npos );
60 
61 #ifdef MFEM_USE_MPI
62  if (pmesh)
63  {
64  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
65  }
66 #endif
67 
68  return err;
69 }
70 
71 // class DataCollection implementation
72 
73 DataCollection::DataCollection(const std::string& collection_name, Mesh *mesh_)
74 {
75  std::string::size_type pos = collection_name.find_last_of('/');
76  if (pos == std::string::npos)
77  {
78  name = collection_name;
79  // leave prefix_path empty
80  }
81  else
82  {
83  prefix_path = collection_name.substr(0, pos+1);
84  name = collection_name.substr(pos+1);
85  }
86  mesh = mesh_;
87  myid = 0;
88  num_procs = 1;
89  serial = true;
90  appendRankToFileName = false;
91 
92 #ifdef MFEM_USE_MPI
93  m_comm = MPI_COMM_NULL;
94  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
95  if (par_mesh)
96  {
97  myid = par_mesh->GetMyRank();
98  num_procs = par_mesh->GetNRanks();
99  m_comm = par_mesh->GetComm();
100  serial = false;
101  appendRankToFileName = true;
102  }
103 #endif
104  own_data = false;
105  cycle = -1;
106  time = 0.0;
107  time_step = 0.0;
110  format = SERIAL_FORMAT; // use serial mesh format
111  compression = false;
112  error = NO_ERROR;
113 }
114 
116 {
117  if (own_data && new_mesh != mesh) { delete mesh; }
118  mesh = new_mesh;
119  myid = 0;
120  num_procs = 1;
121  serial = true;
122  appendRankToFileName = false;
123 
124 #ifdef MFEM_USE_MPI
125  m_comm = MPI_COMM_NULL;
126  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
127  if (par_mesh)
128  {
129  myid = par_mesh->GetMyRank();
130  num_procs = par_mesh->GetNRanks();
131  m_comm = par_mesh->GetComm();
132  serial = false;
133  appendRankToFileName = true;
134  }
135 #endif
136 }
137 
138 #ifdef MFEM_USE_MPI
139 void DataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
140 {
141  // This seems to be the cleanest way to accomplish this
142  // and avoid duplicating fine grained details:
143 
144  SetMesh(new_mesh);
145 
146  m_comm = comm;
147  MPI_Comm_rank(comm, &myid);
148  MPI_Comm_size(comm, &num_procs);
149 }
150 #endif
151 
153 {
154  switch (fmt)
155  {
156  case SERIAL_FORMAT: break;
157 #ifdef MFEM_USE_MPI
158  case PARALLEL_FORMAT: break;
159 #endif
160  default: MFEM_ABORT("unknown format: " << fmt);
161  }
162  format = fmt;
163 }
164 
166 {
167  compression = comp;
168 #ifdef MFEM_USE_GZSTREAM
169  MFEM_ASSERT(!compression, "GZStream not enabled in MFEM build.");
170 #endif
171 }
172 
173 void DataCollection::SetPrefixPath(const std::string& prefix)
174 {
175  if (!prefix.empty())
176  {
177  prefix_path = prefix;
178  if (!prefix_path.empty() && prefix_path[prefix_path.size()-1] != '/')
179  {
180  prefix_path += '/';
181  }
182  }
183  else
184  {
185  prefix_path.clear();
186  }
187 }
188 
189 void DataCollection::Load(int cycle)
190 {
191  MFEM_ABORT("this method is not implemented");
192 }
193 
195 {
196  SaveMesh();
197 
198  if (error) { return; }
199 
200  for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
201  {
202  SaveOneField(it);
203  // Even if there is an error, try saving the other fields
204  }
205 
206  for (QFieldMapIterator it = q_field_map.begin(); it != q_field_map.end();
207  ++it)
208  {
209  SaveOneQField(it);
210  }
211 }
212 
214 {
215  int err;
216 
217  std::string dir_name = prefix_path + name;
218  if (cycle != -1)
219  {
220  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
221  }
222  err = create_directory(dir_name, mesh, myid);
223  if (err)
224  {
225  error = WRITE_ERROR;
226  MFEM_WARNING("Error creating directory: " << dir_name);
227  return; // do not even try to write the mesh
228  }
229 
230  std::string mesh_name = GetMeshFileName();
231  const char *mode = (compression) ? "zwb6" : "w";
232  ofgzstream mesh_file(mesh_name.c_str(), mode);
233  mesh_file.precision(precision);
234 #ifdef MFEM_USE_MPI
235  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
236  if (pmesh && format == PARALLEL_FORMAT)
237  {
238  pmesh->ParPrint(mesh_file);
239  }
240  else
241 #endif
242  {
243  mesh->Print(mesh_file);
244  }
245  if (!mesh_file)
246  {
247  error = WRITE_ERROR;
248  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
249  }
250 }
251 
253 {
254  return (serial || format == SERIAL_FORMAT) ? "mesh" : "pmesh";
255 }
256 
258 {
260 }
261 
262 std::string DataCollection::GetFieldFileName(const std::string &field_name)
263 const
264 {
265  std::string dir_name = prefix_path + name;
266  if (cycle != -1)
267  {
268  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
269  }
270  std::string file_name = dir_name + "/" + field_name;
272  {
273  file_name += "." + to_padded_string(myid, pad_digits_rank);
274  }
275  return file_name;
276 }
277 
279 {
280  const char *mode = (compression) ? "zwb6" : "w";
281  ofgzstream field_file(GetFieldFileName(it->first).c_str(), mode);
282 
283  field_file.precision(precision);
284  (it->second)->Save(field_file);
285  if (!field_file)
286  {
287  error = WRITE_ERROR;
288  MFEM_WARNING("Error writing field to file: " << it->first);
289  }
290 }
291 
293 {
294  const char *mode = (compression) ? "zwb6" : "w";
295  ofgzstream q_field_file(GetFieldFileName(it->first).c_str(), mode);
296  q_field_file.precision(precision);
297  (it->second)->Save(q_field_file);
298  if (!q_field_file)
299  {
300  error = WRITE_ERROR;
301  MFEM_WARNING("Error writing q-field to file: " << it->first);
302  }
303 }
304 
305 void DataCollection::SaveField(const std::string &field_name)
306 {
307  FieldMapIterator it = field_map.find(field_name);
308  if (it != field_map.end())
309  {
310  SaveOneField(it);
311  }
312 }
313 
314 void DataCollection::SaveQField(const std::string &q_field_name)
315 {
316  QFieldMapIterator it = q_field_map.find(q_field_name);
317  if (it != q_field_map.end())
318  {
319  SaveOneQField(it);
320  }
321 }
322 
324 {
325  if (own_data) { delete mesh; }
326  mesh = NULL;
327 
330  own_data = false;
331 }
332 
334 {
335  DeleteData();
336  field_map.clear();
337  q_field_map.clear();
338 }
339 
341 {
342  DeleteData();
343 }
344 
345 
346 // class VisItDataCollection implementation
347 
349 {
350  if (mesh)
351  {
353  topo_dim = mesh->Dimension();
354  if (mesh->NURBSext)
355  {
358  }
359  }
360  else
361  {
362  spatial_dim = 0;
363  topo_dim = 0;
364  }
365 }
366 
367 VisItDataCollection::VisItDataCollection(const std::string& collection_name,
368  Mesh *mesh)
369  : DataCollection(collection_name, mesh)
370 {
371  appendRankToFileName = true; // always include rank in file names
372  cycle = 0; // always include cycle in directory names
373 
376 
377  UpdateMeshInfo();
378 }
379 
380 #ifdef MFEM_USE_MPI
382  const std::string& collection_name,
383  Mesh *mesh)
384  : DataCollection(collection_name, mesh)
385 {
386  m_comm = comm;
387  MPI_Comm_rank(comm, &myid);
388  MPI_Comm_size(comm, &num_procs);
389  appendRankToFileName = true; // always include rank in file names
390  cycle = 0; // always include cycle in directory names
391 
394 
395  UpdateMeshInfo();
396 }
397 #endif
398 
400 {
401  DataCollection::SetMesh(new_mesh);
402  appendRankToFileName = true;
403  UpdateMeshInfo();
404 }
405 
406 #ifdef MFEM_USE_MPI
407 void VisItDataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
408 {
409  // use VisItDataCollection's custom SetMesh, then set MPI info
410  SetMesh(new_mesh);
411  m_comm = comm;
412  MPI_Comm_rank(comm, &myid);
413  MPI_Comm_size(comm, &num_procs);
414 }
415 #endif
416 
417 void VisItDataCollection::RegisterField(const std::string& name,
418  GridFunction *gf)
419 {
421  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
422 
423  int LOD = 1;
424  if (gf->FESpace()->GetNURBSext())
425  {
426  LOD = gf->FESpace()->GetNURBSext()->GetOrder();
427  }
428  else
429  {
430  for (int e=0; e<gf->FESpace()->GetNE() ; e++)
431  {
432  LOD = std::max(LOD,gf->FESpace()->GetFE(e)->GetOrder());
433  }
434  }
435 
437 }
438 
439 void VisItDataCollection::SetLevelsOfDetail(int levels_of_detail)
440 {
441  visit_levels_of_detail = levels_of_detail;
442 }
443 
444 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
445 {
446  visit_max_levels_of_detail = max_levels_of_detail;
447 }
448 
450 {
451  field_info_map.clear();
453 }
454 
456 {
458  SaveRootFile();
459 }
460 
462 {
463  if (myid != 0) { return; }
464 
465  std::string root_name = prefix_path + name + "_" +
467  ".mfem_root";
468  std::ofstream root_file(root_name.c_str());
469  root_file << GetVisItRootString();
470  if (!root_file)
471  {
472  error = WRITE_ERROR;
473  MFEM_WARNING("Error writing VisIt root file: " << root_name);
474  }
475 }
476 
478 {
479  DeleteAll();
480  time_step = 0.0;
481  error = NO_ERROR;
482  cycle = cycle_;
483  std::string root_name = prefix_path + name + "_" +
485  ".mfem_root";
486  LoadVisItRootFile(root_name);
487  if (format != SERIAL_FORMAT || num_procs > 1)
488  {
489 #ifndef MFEM_USE_MPI
490  MFEM_WARNING("Cannot load parallel VisIt root file in serial.");
491  error = READ_ERROR;
492 #else
493  if (m_comm == MPI_COMM_NULL)
494  {
495  MFEM_WARNING("Cannot load parallel VisIt root file without MPI"
496  " communicator");
497  error = READ_ERROR;
498  }
499  else
500  {
501  // num_procs was read from the root file, check for consistency with
502  // the associated MPI_Comm, m_comm:
503  int comm_size;
504  MPI_Comm_size(m_comm, &comm_size);
505  if (comm_size != num_procs)
506  {
507  MFEM_WARNING("Processor number mismatch: VisIt root file: "
508  << num_procs << ", MPI_comm: " << comm_size);
509  error = READ_ERROR;
510  }
511  else
512  {
513  // myid was set when setting m_comm
514  }
515  }
516 #endif
517  }
518  if (!error)
519  {
520  LoadMesh(); // sets own_data to true, when there is no error
521  }
522  if (!error)
523  {
524  LoadFields();
525  }
526  if (error)
527  {
528  DeleteAll();
529  }
530 }
531 
532 void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
533 {
534  std::ifstream root_file(root_name.c_str());
535  std::stringstream buffer;
536  buffer << root_file.rdbuf();
537  if (!buffer)
538  {
539  error = READ_ERROR;
540  MFEM_WARNING("Error reading the VisIt root file: " << root_name);
541  }
542  else
543  {
544  ParseVisItRootString(buffer.str());
545  }
546 }
547 
549 {
550  std::string mesh_fname = GetMeshFileName();
551  named_ifgzstream file(mesh_fname.c_str());
552  // TODO: in parallel, check for errors on all processors
553  if (!file)
554  {
555  error = READ_ERROR;
556  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
557  return;
558  }
559  // TODO: 1) load parallel mesh on one processor
560  if (format == SERIAL_FORMAT)
561  {
562  mesh = new Mesh(file, 1, 0, false);
563  serial = true;
564  }
565  else
566  {
567 #ifdef MFEM_USE_MPI
568  mesh = new ParMesh(m_comm, file);
569  serial = false;
570 #else
571  error = READ_ERROR;
572  MFEM_WARNING("Reading parallel format in serial is not supported");
573  return;
574 #endif
575  }
577  topo_dim = mesh->Dimension();
578  own_data = true;
579 }
580 
582 {
583  std::string path_left = prefix_path + name + "_" +
585  std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
586 
587  field_map.clear();
588  for (FieldInfoMapIterator it = field_info_map.begin();
589  it != field_info_map.end(); ++it)
590  {
591  std::string fname = path_left + it->first + path_right;
592  ifgzstream file(fname.c_str());
593  // TODO: in parallel, check for errors on all processors
594  if (!file)
595  {
596  error = READ_ERROR;
597  MFEM_WARNING("Unable to open field file: " << fname);
598  return;
599  }
600  // TODO: 1) load parallel GridFunction on one processor
601  if (serial)
602  {
603  field_map.Register(it->first, new GridFunction(mesh, file), own_data);
604  }
605  else
606  {
607 #ifdef MFEM_USE_MPI
609  it->first,
610  new ParGridFunction(dynamic_cast<ParMesh*>(mesh), file), own_data);
611 #else
612  error = READ_ERROR;
613  MFEM_WARNING("Reading parallel format in serial is not supported");
614  return;
615 #endif
616  }
617  }
618 }
619 
621 {
622  // Get the path string (relative to where the root file is, i.e. no prefix).
623  std::string path_str =
625 
626  // We have to build the json tree inside out to get all the values in there
627  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
628 
629  // Build the mesh data
630  std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
631  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
632  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
633  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
634  mesh["path"] = picojson::value(path_str + GetMeshShortFileName() +
635  file_ext_format);
636  mesh["tags"] = picojson::value(mtags);
637  mesh["format"] = picojson::value(to_string(format));
638 
639  // Build the fields data entries
640  for (FieldInfoMapIterator it = field_info_map.begin();
641  it != field_info_map.end(); ++it)
642  {
643  ftags["assoc"] = picojson::value((it->second).association);
644  ftags["comps"] = picojson::value(to_string((it->second).num_components));
645  ftags["lod"] = picojson::value(to_string(visit_levels_of_detail));
646  field["path"] = picojson::value(path_str + it->first + file_ext_format);
647  field["tags"] = picojson::value(ftags);
648  fields[it->first] = picojson::value(field);
649  }
650 
651  main["cycle"] = picojson::value(double(cycle));
652  main["time"] = picojson::value(time);
653  main["time_step"] = picojson::value(time_step);
654  main["domains"] = picojson::value(double(num_procs));
655  main["mesh"] = picojson::value(mesh);
656  if (!field_info_map.empty())
657  {
658  main["fields"] = picojson::value(fields);
659  }
660 
661  dsets["main"] = picojson::value(main);
662  top["dsets"] = picojson::value(dsets);
663 
664  return picojson::value(top).serialize(true);
665 }
666 
667 void VisItDataCollection::ParseVisItRootString(const std::string& json)
668 {
669  picojson::value top, dsets, main, mesh, fields;
670  std::string parse_err = picojson::parse(top, json);
671  if (!parse_err.empty())
672  {
673  error = READ_ERROR;
674  MFEM_WARNING("Unable to parse VisIt root data.");
675  return;
676  }
677 
678  // Process "main"
679  dsets = top.get("dsets");
680  main = dsets.get("main");
681  cycle = int(main.get("cycle").get<double>());
682  time = main.get("time").get<double>();
683  if (main.contains("time_step"))
684  {
685  time_step = main.get("time_step").get<double>();
686  }
687  num_procs = int(main.get("domains").get<double>());
688  mesh = main.get("mesh");
689  fields = main.get("fields");
690 
691  // ... Process "mesh"
692 
693  // Set the DataCollection::name using the mesh path
694  std::string path = mesh.get("path").get<std::string>();
695  size_t right_sep = path.find('_');
696  if (right_sep == std::string::npos)
697  {
698  error = READ_ERROR;
699  MFEM_WARNING("Unable to parse VisIt root data.");
700  return;
701  }
702  name = path.substr(0, right_sep);
703 
704  if (mesh.contains("format"))
705  {
706  format = to_int(mesh.get("format").get<std::string>());
707  }
708  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
709  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
711  to_int(mesh.get("tags").get("max_lods").get<std::string>());
712 
713  // ... Process "fields"
714  field_info_map.clear();
715  if (fields.is<picojson::object>())
716  {
717  picojson::object fields_obj = fields.get<picojson::object>();
718  for (picojson::object::iterator it = fields_obj.begin();
719  it != fields_obj.end(); ++it)
720  {
721  picojson::value tags = it->second.get("tags");
722  field_info_map[it->first] =
723  VisItFieldInfo(tags.get("assoc").get<std::string>(),
724  to_int(tags.get("comps").get<std::string>()));
725  }
726  }
727 }
728 
729 } // 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:1156
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 GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:315
int VectorDim() const
Definition: gridfunc.cpp:302
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.
void SetCompression(bool comp)
Set the flag for use of gz compressed files.
int main(int argc, char *argv[])
Definition: ex1.cpp:58
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:224
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.
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:374
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.
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
static const int pad_digits_default
Default value for pad_digits_*.
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:430
int Dimension() const
Definition: mesh.hpp:713
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
int GetOrder() const
If all orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
Definition: nurbs.hpp:340
int SpaceDimension() const
Definition: mesh.hpp:714
std::map< std::string, VisItFieldInfo > field_info_map
int GetMyRank() const
Definition: pmesh.hpp:225
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:274
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:223
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.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
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...
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1541
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:4969
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.