MFEM  v4.3.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-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "fem.hpp"
13 #include "../mesh/nurbs.hpp"
14 #include "../general/binaryio.hpp"
15 #include "../general/text.hpp"
16 #include "picojson.h"
17 
18 #include <cerrno> // errno
19 #include <sstream>
20 #include <regex>
21 
22 #ifndef _WIN32
23 #include <sys/stat.h> // mkdir
24 #else
25 #include <direct.h> // _mkdir
26 #define mkdir(dir, mode) _mkdir(dir)
27 #endif
28 
29 namespace mfem
30 {
31 
32 // static method
33 int DataCollection::create_directory(const std::string &dir_name,
34  const Mesh *mesh, int myid)
35 {
36  // create directories recursively
37  const char path_delim = '/';
38  std::string::size_type pos = 0;
39  int err;
40 #ifdef MFEM_USE_MPI
41  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
42 #endif
43 
44  do
45  {
46  pos = dir_name.find(path_delim, pos+1);
47  std::string subdir = dir_name.substr(0, pos);
48 
49 #ifndef MFEM_USE_MPI
50  err = mkdir(subdir.c_str(), 0777);
51  err = (err && (errno != EEXIST)) ? 1 : 0;
52 #else
53  if (myid == 0 || pmesh == NULL)
54  {
55  err = mkdir(subdir.c_str(), 0777);
56  err = (err && (errno != EEXIST)) ? 1 : 0;
57  }
58 #endif
59  }
60  while ( pos != std::string::npos );
61 
62 #ifdef MFEM_USE_MPI
63  if (pmesh)
64  {
65  MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->GetComm());
66  }
67 #endif
68 
69  return err;
70 }
71 
72 // class DataCollection implementation
73 
74 DataCollection::DataCollection(const std::string& collection_name, Mesh *mesh_)
75 {
76  std::string::size_type pos = collection_name.find_last_of('/');
77  if (pos == std::string::npos)
78  {
79  name = collection_name;
80  // leave prefix_path empty
81  }
82  else
83  {
84  prefix_path = collection_name.substr(0, pos+1);
85  name = collection_name.substr(pos+1);
86  }
87  mesh = mesh_;
88  myid = 0;
89  num_procs = 1;
90  serial = true;
91  appendRankToFileName = false;
92 
93 #ifdef MFEM_USE_MPI
94  m_comm = MPI_COMM_NULL;
95  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
96  if (par_mesh)
97  {
98  myid = par_mesh->GetMyRank();
99  num_procs = par_mesh->GetNRanks();
100  m_comm = par_mesh->GetComm();
101  serial = false;
102  appendRankToFileName = true;
103  }
104 #endif
105  own_data = false;
106  cycle = -1;
107  time = 0.0;
108  time_step = 0.0;
111  format = SERIAL_FORMAT; // use serial mesh format
112  compression = false;
113  error = NO_ERROR;
114 }
115 
117 {
118  if (own_data && new_mesh != mesh) { delete mesh; }
119  mesh = new_mesh;
120  myid = 0;
121  num_procs = 1;
122  serial = true;
123  appendRankToFileName = false;
124 
125 #ifdef MFEM_USE_MPI
126  m_comm = MPI_COMM_NULL;
127  ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
128  if (par_mesh)
129  {
130  myid = par_mesh->GetMyRank();
131  num_procs = par_mesh->GetNRanks();
132  m_comm = par_mesh->GetComm();
133  serial = false;
134  appendRankToFileName = true;
135  }
136 #endif
137 }
138 
139 #ifdef MFEM_USE_MPI
140 void DataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
141 {
142  // This seems to be the cleanest way to accomplish this
143  // and avoid duplicating fine grained details:
144 
145  SetMesh(new_mesh);
146 
147  m_comm = comm;
148  MPI_Comm_rank(comm, &myid);
149  MPI_Comm_size(comm, &num_procs);
150 }
151 #endif
152 
154 {
155  switch (fmt)
156  {
157  case SERIAL_FORMAT: break;
158 #ifdef MFEM_USE_MPI
159  case PARALLEL_FORMAT: break;
160 #endif
161  default: MFEM_ABORT("unknown format: " << fmt);
162  }
163  format = fmt;
164 }
165 
167 {
168  compression = comp;
169 #ifndef MFEM_USE_ZLIB
170  MFEM_VERIFY(!compression, "ZLib not enabled in MFEM build.");
171 #endif
172 }
173 
174 void DataCollection::SetPrefixPath(const std::string& prefix)
175 {
176  if (!prefix.empty())
177  {
178  prefix_path = prefix;
179  if (!prefix_path.empty() && prefix_path[prefix_path.size()-1] != '/')
180  {
181  prefix_path += '/';
182  }
183  }
184  else
185  {
186  prefix_path.clear();
187  }
188 }
189 
190 void DataCollection::Load(int cycle)
191 {
192  MFEM_ABORT("this method is not implemented");
193 }
194 
196 {
197  SaveMesh();
198 
199  if (error) { return; }
200 
201  for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
202  {
203  SaveOneField(it);
204  // Even if there is an error, try saving the other fields
205  }
206 
207  for (QFieldMapIterator it = q_field_map.begin(); it != q_field_map.end();
208  ++it)
209  {
210  SaveOneQField(it);
211  }
212 }
213 
215 {
216  int err;
217 
218  std::string dir_name = prefix_path + name;
219  if (cycle != -1)
220  {
221  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
222  }
223  err = create_directory(dir_name, mesh, myid);
224  if (err)
225  {
226  error = WRITE_ERROR;
227  MFEM_WARNING("Error creating directory: " << dir_name);
228  return; // do not even try to write the mesh
229  }
230 
231  std::string mesh_name = GetMeshFileName();
232  mfem::ofgzstream mesh_file(mesh_name, compression);
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  mfem::ofgzstream field_file(GetFieldFileName(it->first), compression);
281 
282  field_file.precision(precision);
283  (it->second)->Save(field_file);
284  if (!field_file)
285  {
286  error = WRITE_ERROR;
287  MFEM_WARNING("Error writing field to file: " << it->first);
288  }
289 }
290 
292 {
293  mfem::ofgzstream q_field_file(GetFieldFileName(it->first), compression);
294 
295  q_field_file.precision(precision);
296  (it->second)->Save(q_field_file);
297  if (!q_field_file)
298  {
299  error = WRITE_ERROR;
300  MFEM_WARNING("Error writing q-field to file: " << it->first);
301  }
302 }
303 
304 void DataCollection::SaveField(const std::string &field_name)
305 {
306  FieldMapIterator it = field_map.find(field_name);
307  if (it != field_map.end())
308  {
309  SaveOneField(it);
310  }
311 }
312 
313 void DataCollection::SaveQField(const std::string &q_field_name)
314 {
315  QFieldMapIterator it = q_field_map.find(q_field_name);
316  if (it != q_field_map.end())
317  {
318  SaveOneQField(it);
319  }
320 }
321 
323 {
324  if (own_data) { delete mesh; }
325  mesh = NULL;
326 
329  own_data = false;
330 }
331 
333 {
334  DeleteData();
335  field_map.clear();
336  q_field_map.clear();
337 }
338 
340 {
341  DeleteData();
342 }
343 
344 
345 // class VisItDataCollection implementation
346 
348 {
349  if (mesh)
350  {
352  topo_dim = mesh->Dimension();
353  if (mesh->NURBSext)
354  {
357  }
358  }
359  else
360  {
361  spatial_dim = 0;
362  topo_dim = 0;
363  }
364 }
365 
366 VisItDataCollection::VisItDataCollection(const std::string& collection_name,
367  Mesh *mesh)
368  : DataCollection(collection_name, mesh)
369 {
370  appendRankToFileName = true; // always include rank in file names
371  cycle = 0; // always include cycle in directory names
372 
375 
376  UpdateMeshInfo();
377 }
378 
379 #ifdef MFEM_USE_MPI
381  const std::string& collection_name,
382  Mesh *mesh)
383  : DataCollection(collection_name, mesh)
384 {
385  m_comm = comm;
386  MPI_Comm_rank(comm, &myid);
387  MPI_Comm_size(comm, &num_procs);
388  appendRankToFileName = true; // always include rank in file names
389  cycle = 0; // always include cycle in directory names
390 
393 
394  UpdateMeshInfo();
395 }
396 #endif
397 
399 {
400  DataCollection::SetMesh(new_mesh);
401  appendRankToFileName = true;
402  UpdateMeshInfo();
403 }
404 
405 #ifdef MFEM_USE_MPI
406 void VisItDataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
407 {
408  // use VisItDataCollection's custom SetMesh, then set MPI info
409  SetMesh(new_mesh);
410  m_comm = comm;
411  MPI_Comm_rank(comm, &myid);
412  MPI_Comm_size(comm, &num_procs);
413 }
414 #endif
415 
416 void VisItDataCollection::RegisterField(const std::string& name,
417  GridFunction *gf)
418 {
419  int LOD = 1;
420  if (gf->FESpace()->GetNURBSext())
421  {
422  LOD = gf->FESpace()->GetNURBSext()->GetOrder();
423  }
424  else
425  {
426  for (int e=0; e<gf->FESpace()->GetNE(); e++)
427  {
428  LOD = std::max(LOD,gf->FESpace()->GetFE(e)->GetOrder());
429  }
430  }
431 
433  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim(), LOD);
434  visit_levels_of_detail = std::max(visit_levels_of_detail, LOD);
435 }
436 
437 void VisItDataCollection::RegisterQField(const std::string& name,
438  QuadratureFunction *qf)
439 {
440  int LOD = -1;
441  Mesh *mesh = qf->GetSpace()->GetMesh();
442  for (int e=0; e<qf->GetSpace()->GetNE(); e++)
443  {
445  mesh->GetElementBaseGeometry(e),
446  qf->GetElementIntRule(e).GetNPoints());
447 
448  LOD = std::max(LOD,locLOD);
449  }
450 
452  field_info_map[name] = VisItFieldInfo("elements", 1, LOD);
454 }
455 
456 void VisItDataCollection::SetLevelsOfDetail(int levels_of_detail)
457 {
458  visit_levels_of_detail = levels_of_detail;
459 }
460 
461 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
462 {
463  visit_max_levels_of_detail = max_levels_of_detail;
464 }
465 
467 {
468  field_info_map.clear();
470 }
471 
473 {
475  SaveRootFile();
476 }
477 
479 {
480  if (myid != 0) { return; }
481 
482  std::string root_name = prefix_path + name + "_" +
484  ".mfem_root";
485  std::ofstream root_file(root_name.c_str());
486  root_file << GetVisItRootString();
487  if (!root_file)
488  {
489  error = WRITE_ERROR;
490  MFEM_WARNING("Error writing VisIt root file: " << root_name);
491  }
492 }
493 
495 {
496  DeleteAll();
497  time_step = 0.0;
498  error = NO_ERROR;
499  cycle = cycle_;
500  std::string root_name = prefix_path + name + "_" +
502  ".mfem_root";
503  LoadVisItRootFile(root_name);
504  if (format != SERIAL_FORMAT || num_procs > 1)
505  {
506 #ifndef MFEM_USE_MPI
507  MFEM_WARNING("Cannot load parallel VisIt root file in serial.");
508  error = READ_ERROR;
509 #else
510  if (m_comm == MPI_COMM_NULL)
511  {
512  MFEM_WARNING("Cannot load parallel VisIt root file without MPI"
513  " communicator");
514  error = READ_ERROR;
515  }
516  else
517  {
518  // num_procs was read from the root file, check for consistency with
519  // the associated MPI_Comm, m_comm:
520  int comm_size;
521  MPI_Comm_size(m_comm, &comm_size);
522  if (comm_size != num_procs)
523  {
524  MFEM_WARNING("Processor number mismatch: VisIt root file: "
525  << num_procs << ", MPI_comm: " << comm_size);
526  error = READ_ERROR;
527  }
528  else
529  {
530  // myid was set when setting m_comm
531  }
532  }
533 #endif
534  }
535  if (!error)
536  {
537  LoadMesh(); // sets own_data to true, when there is no error
538  }
539  if (!error)
540  {
541  LoadFields();
542  }
543  if (error)
544  {
545  DeleteAll();
546  }
547 }
548 
549 void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
550 {
551  std::ifstream root_file(root_name.c_str());
552  std::stringstream buffer;
553  buffer << root_file.rdbuf();
554  if (!buffer)
555  {
556  error = READ_ERROR;
557  MFEM_WARNING("Error reading the VisIt root file: " << root_name);
558  }
559  else
560  {
561  ParseVisItRootString(buffer.str());
562  }
563 }
564 
566 {
567  // GetMeshFileName() uses 'serial', so we need to set it in advance.
568  serial = (format == SERIAL_FORMAT);
569  std::string mesh_fname = GetMeshFileName();
570  named_ifgzstream file(mesh_fname);
571  // TODO: in parallel, check for errors on all processors
572  if (!file)
573  {
574  error = READ_ERROR;
575  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
576  return;
577  }
578  // TODO: 1) load parallel mesh on one processor
579  if (format == SERIAL_FORMAT)
580  {
581  mesh = new Mesh(file, 1, 0, false);
582  serial = true;
583  }
584  else
585  {
586 #ifdef MFEM_USE_MPI
587  mesh = new ParMesh(m_comm, file);
588  serial = false;
589 #else
590  error = READ_ERROR;
591  MFEM_WARNING("Reading parallel format in serial is not supported");
592  return;
593 #endif
594  }
596  topo_dim = mesh->Dimension();
597  own_data = true;
598 }
599 
601 {
602  std::string path_left = prefix_path + name + "_" +
604  std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
605 
606  field_map.clear();
607  for (FieldInfoMapIterator it = field_info_map.begin();
608  it != field_info_map.end(); ++it)
609  {
610  std::string fname = path_left + it->first + path_right;
611  mfem::ifgzstream file(fname);
612  // TODO: in parallel, check for errors on all processors
613  if (!file)
614  {
615  error = READ_ERROR;
616  MFEM_WARNING("Unable to open field file: " << fname);
617  return;
618  }
619  // TODO: 1) load parallel GridFunction on one processor
620  if (serial)
621  {
622  if ((it->second).association == "nodes")
623  {
624  field_map.Register(it->first, new GridFunction(mesh, file), own_data);
625  }
626  else if ((it->second).association == "elements")
627  {
628  q_field_map.Register(it->first, new QuadratureFunction(mesh, file), own_data);
629  }
630  }
631  else
632  {
633 #ifdef MFEM_USE_MPI
634  if ((it->second).association == "nodes")
635  {
637  it->first,
638  new ParGridFunction(dynamic_cast<ParMesh*>(mesh), file), own_data);
639  }
640  else if ((it->second).association == "elements")
641  {
642  q_field_map.Register(it->first, new QuadratureFunction(mesh, file), own_data);
643  }
644 #else
645  error = READ_ERROR;
646  MFEM_WARNING("Reading parallel format in serial is not supported");
647  return;
648 #endif
649  }
650  }
651 }
652 
654 {
655  // Get the path string (relative to where the root file is, i.e. no prefix).
656  std::string path_str =
658 
659  // We have to build the json tree inside out to get all the values in there
660  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
661 
662  // Build the mesh data
663  std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
664  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
665  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
666  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
667  mesh["path"] = picojson::value(path_str + GetMeshShortFileName() +
668  file_ext_format);
669  mesh["tags"] = picojson::value(mtags);
670  mesh["format"] = picojson::value(to_string(format));
671 
672  // Build the fields data entries
673  for (FieldInfoMapIterator it = field_info_map.begin();
674  it != field_info_map.end(); ++it)
675  {
676  ftags["assoc"] = picojson::value((it->second).association);
677  ftags["comps"] = picojson::value(to_string((it->second).num_components));
678  ftags["lod"] = picojson::value(to_string((it->second).lod));
679  field["path"] = picojson::value(path_str + it->first + file_ext_format);
680  field["tags"] = picojson::value(ftags);
681  fields[it->first] = picojson::value(field);
682  }
683 
684  main["cycle"] = picojson::value(double(cycle));
685  main["time"] = picojson::value(time);
686  main["time_step"] = picojson::value(time_step);
687  main["domains"] = picojson::value(double(num_procs));
688  main["mesh"] = picojson::value(mesh);
689  if (!field_info_map.empty())
690  {
691  main["fields"] = picojson::value(fields);
692  }
693 
694  dsets["main"] = picojson::value(main);
695  top["dsets"] = picojson::value(dsets);
696 
697  return picojson::value(top).serialize(true);
698 }
699 
700 void VisItDataCollection::ParseVisItRootString(const std::string& json)
701 {
702  picojson::value top, dsets, main, mesh, fields;
703  std::string parse_err = picojson::parse(top, json);
704  if (!parse_err.empty())
705  {
706  error = READ_ERROR;
707  MFEM_WARNING("Unable to parse VisIt root data.");
708  return;
709  }
710 
711  // Process "main"
712  dsets = top.get("dsets");
713  main = dsets.get("main");
714  cycle = int(main.get("cycle").get<double>());
715  time = main.get("time").get<double>();
716  if (main.contains("time_step"))
717  {
718  time_step = main.get("time_step").get<double>();
719  }
720  num_procs = int(main.get("domains").get<double>());
721  mesh = main.get("mesh");
722  fields = main.get("fields");
723 
724  // ... Process "mesh"
725 
726  // Set the DataCollection::name using the mesh path
727  std::string path = mesh.get("path").get<std::string>();
728  size_t right_sep = path.find('_');
729  if (right_sep == std::string::npos)
730  {
731  error = READ_ERROR;
732  MFEM_WARNING("Unable to parse VisIt root data.");
733  return;
734  }
735  name = path.substr(0, right_sep);
736 
737  if (mesh.contains("format"))
738  {
739  format = to_int(mesh.get("format").get<std::string>());
740  }
741  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
742  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
744  to_int(mesh.get("tags").get("max_lods").get<std::string>());
745 
746  // ... Process "fields"
747  field_info_map.clear();
748  if (fields.is<picojson::object>())
749  {
750  picojson::object fields_obj = fields.get<picojson::object>();
751  for (picojson::object::iterator it = fields_obj.begin();
752  it != fields_obj.end(); ++it)
753  {
754  picojson::value tags = it->second.get("tags");
755  field_info_map[it->first] =
756  VisItFieldInfo(tags.get("assoc").get<std::string>(),
757  to_int(tags.get("comps").get<std::string>()));
758  }
759  }
760 }
761 
763  collection_name,
764  Mesh *mesh_)
765  : DataCollection(collection_name, mesh_),
766  levels_of_detail(1),
767  pv_data_format(VTKFormat::BINARY),
768  high_order_output(false),
769  restart_mode(false)
770 {
771 #ifdef MFEM_USE_ZLIB
772  compression = -1; // default zlib compression level, equivalent to 6
773 #else
774  compression = 0;
775 #endif
776 }
777 
778 void ParaViewDataCollection::SetLevelsOfDetail(int levels_of_detail_)
779 {
780  levels_of_detail = levels_of_detail_;
781 }
782 
784 {
785  MFEM_WARNING("ParaViewDataCollection::Load() is not implemented!");
786 }
787 
789 {
790  std::string out = "";
792  return out;
793 }
794 
796 {
797  std::string out = "Cycle" + to_padded_string(cycle,pad_digits_cycle);
798  return out;
799 }
800 
802 {
803  std::string out = GeneratePVTUPath();
804  return out;
805 }
806 
808 {
809  std::string out = GetCollectionName()+".pvd";
810  return out;
811 }
812 
814 {
815  std::string out = "data.pvtu";
816  return out;
817 }
818 
820 {
821  std::string out = "proc" + to_padded_string(myid,pad_digits_rank)+".vtu";
822  return out;
823 }
825 {
826  std::string out = "proc" + to_padded_string(crank,pad_digits_rank)+".vtu";
827  return out;
828 }
829 
831 {
832  // add a new collection to the PDV file
833 
834  // check if the directories are created
835  {
836  std::string path = GenerateCollectionPath()+"/"+GenerateVTUPath();
837  int err = create_directory(path, mesh, myid);
838  if (err)
839  {
840  error = WRITE_ERROR;
841  MFEM_WARNING("Error creating directory: " << path);
842  return; // do not even try to write the mesh
843  }
844  }
845  // the directory is created
846 
847  // create pvd file if needed. If we are not in restart mode, a new pvd file
848  // is always created. In restart mode, we keep any previously defined
849  // timestep values as long as they are less than the currently defined time.
850 
851  if (myid == 0 && !pvd_stream.is_open())
852  {
853  std::string dpath=GenerateCollectionPath();
854  std::string pvdname=dpath+"/"+GeneratePVDFileName();
855 
856  std::ifstream pvd_in;
857  if (restart_mode && (pvd_in.open(pvdname,std::ios::binary),pvd_in.good()))
858  {
859  // PVD file exists and restart mode enabled: preserve existing time
860  // steps less than the current time.
861  std::fstream::pos_type pos_begin = pvd_in.tellg();
862  std::fstream::pos_type pos_end = pos_begin;
863 
864  std::regex regexp("timestep=\"([^[:space:]]+)\".*file=\"Cycle(\\d+)");
865  std::smatch match;
866 
867  std::string line;
868  while (getline(pvd_in,line))
869  {
870  if (regex_search(line,match,regexp))
871  {
872  MFEM_ASSERT(match.size() == 3, "Unable to parse DataSet");
873  double tvalue = std::stod(match[1]);
874  if (tvalue >= GetTime()) { break; }
875  int cvalue = std::stoi(match[2]);
876  MFEM_VERIFY(cvalue < GetCycle(), "Cycle " << GetCycle() <<
877  " is too small for restart mode: trying to overwrite"
878  " existing data.");
879  pos_end = pvd_in.tellg();
880  }
881  }
882  size_t count = pos_end - pos_begin;
883  std::vector<char> buf(count);
884  pvd_in.clear();
885  pvd_in.seekg(pos_begin);
886  pvd_in.read(buf.data(), count);
887  pvd_in.close();
888  pvd_stream.open(pvdname.c_str(),std::ios::out);
889  pvd_stream.write(buf.data(), count);
890  }
891  else
892  {
893  // initialize new pvd file
894  pvd_stream.open(pvdname.c_str(),std::ios::out);
895  // initialize the file
896  pvd_stream << "<?xml version=\"1.0\"?>\n";
897  pvd_stream << "<VTKFile type=\"Collection\" version=\"0.1\"";
898  pvd_stream << " byte_order=\"" << VTKByteOrder() << "\">\n";
899  pvd_stream << "<Collection>" << std::endl;
900  }
901  }
902 
903  // define the vtu file
904  {
905  std::string fname = GenerateCollectionPath()+"/"+GenerateVTUPath()+"/"
907  std::fstream out(fname.c_str(), std::ios::out);
908  out.precision(precision);
909  SaveDataVTU(out,levels_of_detail);
910  out.close();
911  }
912 
913  // define the pvtu file only on process 0
914  if (myid==0)
915  {
916  std::string fname = GenerateCollectionPath()+"/"+GeneratePVTUPath()+"/"
918  std::fstream out(fname.c_str(), std::ios::out);
919 
920  out << "<?xml version=\"1.0\"?>\n";
921  out << "<VTKFile type=\"PUnstructuredGrid\"";
922  out << " version =\"0.1\" byte_order=\"" << VTKByteOrder() << "\">\n";
923  out << "<PUnstructuredGrid GhostLevel=\"0\">\n";
924 
925  out << "<PPoints>\n";
926  out << "\t<PDataArray type=\"" << GetDataTypeString() << "\" ";
927  out << " Name=\"Points\" NumberOfComponents=\"3\""
928  << " format=\"" << GetDataFormatString() << "\"/>\n";
929  out << "</PPoints>\n";
930 
931  out << "<PCells>\n";
932  out << "\t<PDataArray type=\"Int32\" ";
933  out << " Name=\"connectivity\" NumberOfComponents=\"1\""
934  << " format=\"" << GetDataFormatString() << "\"/>\n";
935  out << "\t<PDataArray type=\"Int32\" ";
936  out << " Name=\"offsets\" NumberOfComponents=\"1\""
937  << " format=\"" << GetDataFormatString() << "\"/>\n";
938  out << "\t<PDataArray type=\"UInt8\" ";
939  out << " Name=\"types\" NumberOfComponents=\"1\""
940  << " format=\"" << GetDataFormatString() << "\"/>\n";
941  out << "</PCells>\n";
942 
943  out << "<PPointData>\n";
944  for (FieldMapIterator it=field_map.begin(); it!=field_map.end(); ++it)
945  {
946  int vec_dim=it->second->VectorDim();
947  out << "<PDataArray type=\"" << GetDataTypeString()
948  << "\" Name=\"" << it->first
949  << "\" NumberOfComponents=\"" << vec_dim << "\" "
950  << "format=\"" << GetDataFormatString() << "\" />\n";
951  }
952  out << "</PPointData>\n";
953 
954  // CELL DATA
955  out << "<PCellData>\n";
956  out << "\t<PDataArray type=\"Int32\" Name=\"" << "attribute"
957  << "\" NumberOfComponents=\"1\""
958  << " format=\"" << GetDataFormatString() << "\"/>\n";
959  out << "</PCellData>\n";
960 
961  for (int ii=0; ii<num_procs; ii++)
962  {
963  // this one is generated without the path
964  std::string nfname=GenerateVTUFileName(ii);
965  out << "<Piece Source=\"" << nfname << "\"/>\n";
966  }
967  out << "</PUnstructuredGrid>\n";
968  out << "</VTKFile>\n";
969  out.close();
970 
971  fname = GeneratePVTUPath()+"/"+GeneratePVTUFileName();
972  // add the pvtu file to the pvd_stream
973  pvd_stream << "<DataSet timestep=\"" << GetTime(); // GetCycle();
974  pvd_stream << "\" group=\"\" part=\"" << 0 << "\" file=\"";
975  pvd_stream << fname << "\"/>\n";
976  std::fstream::pos_type pos = pvd_stream.tellp();
977  pvd_stream << "</Collection>\n";
978  pvd_stream << "</VTKFile>" << std::endl;
979  pvd_stream.seekp(pos);
980  }
981 }
982 
983 void ParaViewDataCollection::SaveDataVTU(std::ostream &out, int ref)
984 {
985  out << "<VTKFile type=\"UnstructuredGrid\"";
986  if (compression != 0)
987  {
988  out << " compressor=\"vtkZLibDataCompressor\"";
989  }
990  out << " version=\"0.1\" byte_order=\"" << VTKByteOrder() << "\">\n";
991  out << "<UnstructuredGrid>\n";
992  mesh->PrintVTU(out,ref,pv_data_format,high_order_output,compression);
993 
994  // dump out the grid functions as point data
995  out << "<PointData >\n";
996  // save the grid functions
997  // iterate over all grid functions
998  for (FieldMapIterator it=field_map.begin(); it!=field_map.end(); ++it)
999  {
1000  SaveGFieldVTU(out,ref,it);
1001  }
1002  // iterate over all quadrature functions
1003  // if the Quadrature functions are dumped as cell data
1004  // the cycle should be moved before the grid functions
1005  // and the PrintVTU CellData section should be open in the mesh dump
1006  for (QFieldMapIterator it=q_field_map.begin(); it!=q_field_map.end(); ++it)
1007  {
1008  // save the quadrature functions
1009  // this one is not implemented yet
1010  SaveQFieldVTU(out,ref,it);
1011  }
1012  out << "</PointData>\n";
1013  // close the mesh
1014  out << "</Piece>\n"; // close the piece open in the PrintVTU method
1015  out << "</UnstructuredGrid>\n";
1016  out << "</VTKFile>" << std::endl;
1017 }
1018 
1019 void ParaViewDataCollection::SaveQFieldVTU(std::ostream &out, int ref,
1020  const QFieldMapIterator& it )
1021 {
1022  MFEM_WARNING("SaveQFieldVTU is not currently implemented - field name:"<<it->second);
1023 }
1024 
1025 void ParaViewDataCollection::SaveGFieldVTU(std::ostream &out, int ref_,
1026  const FieldMapIterator& it)
1027 {
1028  RefinedGeometry *RefG;
1029  Vector val;
1030  DenseMatrix vval, pmat;
1031  std::vector<char> buf;
1032  int vec_dim = it->second->VectorDim();
1033  if (vec_dim == 1)
1034  {
1035  // scalar data
1036  out << "<DataArray type=\"" << GetDataTypeString()
1037  << "\" Name=\"" << it->first;
1038  out << "\" NumberOfComponents=\"1\" format=\""
1039  << GetDataFormatString() << "\" >\n";
1040  for (int i = 0; i < mesh->GetNE(); i++)
1041  {
1042  RefG = GlobGeometryRefiner.Refine(
1043  mesh->GetElementBaseGeometry(i), ref_, 1);
1044  it->second->GetValues(i, RefG->RefPts, val, pmat);
1045  for (int j = 0; j < val.Size(); j++)
1046  {
1047  if (pv_data_format == VTKFormat::ASCII)
1048  {
1049  out << ZeroSubnormal(val(j)) << '\n';
1050  }
1051  else if (pv_data_format == VTKFormat::BINARY)
1052  {
1053  bin_io::AppendBytes(buf, val(j));
1054  }
1055  else
1056  {
1057  bin_io::AppendBytes<float>(buf, float(val(j)));
1058  }
1059  }
1060  }
1061  }
1062  else
1063  {
1064  // vector data
1065  out << "<DataArray type=\"" << GetDataTypeString()
1066  << "\" Name=\"" << it->first;
1067  out << "\" NumberOfComponents=\"" << vec_dim << "\""
1068  << " format=\"" << GetDataFormatString() << "\" >" << '\n';
1069  for (int i = 0; i < mesh->GetNE(); i++)
1070  {
1071  RefG = GlobGeometryRefiner.Refine(
1072  mesh->GetElementBaseGeometry(i), ref_, 1);
1073 
1074  it->second->GetVectorValues(i, RefG->RefPts, vval, pmat);
1075 
1076  for (int jj = 0; jj < vval.Width(); jj++)
1077  {
1078  for (int ii = 0; ii < vval.Height(); ii++)
1079  {
1080  if (pv_data_format == VTKFormat::ASCII)
1081  {
1082  out << ZeroSubnormal(vval(ii,jj)) << ' ';
1083  }
1084  else if (pv_data_format == VTKFormat::BINARY)
1085  {
1086  bin_io::AppendBytes(buf, vval(ii,jj));
1087  }
1088  else
1089  {
1090  bin_io::AppendBytes<float>(buf, float(vval(ii,jj)));
1091  }
1092  }
1093  if (pv_data_format == VTKFormat::ASCII) { out << '\n'; }
1094  }
1095  }
1096  }
1097 
1098  if (IsBinaryFormat())
1099  {
1100  WriteVTKEncodedCompressed(out,buf.data(),buf.size(),compression);
1101  out << '\n';
1102  }
1103  out << "</DataArray>" << std::endl;
1104 }
1105 
1107 {
1108  pv_data_format = fmt;
1109 }
1110 
1112 {
1113  return pv_data_format != VTKFormat::ASCII;
1114 }
1115 
1116 void ParaViewDataCollection::SetHighOrderOutput(bool high_order_output_)
1117 {
1118  high_order_output = high_order_output_;
1119 }
1120 
1122 {
1123  MFEM_ASSERT(compression_level_ >= -1 && compression_level_ <= 9,
1124  "Compression level must be between -1 and 9 (inclusive).");
1125  compression = compression_level_;
1126 }
1127 
1129 {
1130  // If we are enabling compression, and it was disabled previously, use the
1131  // default compression level. Otherwise, leave the compression level
1132  // unchanged.
1133  if (compression_ && compression == 0)
1134  {
1135  SetCompressionLevel(-1);
1136  }
1137 }
1138 
1140 {
1141  restart_mode = restart_mode_;
1142 }
1143 
1145 {
1146  if (pv_data_format == VTKFormat::ASCII)
1147  {
1148  return "ascii";
1149  }
1150  else
1151  {
1152  return "binary";
1153  }
1154 }
1155 
1157 {
1158  if (pv_data_format==VTKFormat::ASCII || pv_data_format==VTKFormat::BINARY)
1159  {
1160  return "Float64";
1161  }
1162  else
1163  {
1164  return "Float32";
1165  }
1166 }
1167 
1168 } // end namespace MFEM
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
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:1387
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: gridfunc.hpp:832
int format
Output mesh format: see the Format enumeration.
int error
Error state.
void SetDataFormat(VTKFormat fmt)
double GetTime() const
Get physical time (for time-dependent simulations)
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
void DeleteAll()
Delete data owned by the DataCollection including field information.
const char * GetDataTypeString() const
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 Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
const std::string & GetCollectionName() const
Get the name of the collection.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:327
int VectorDim() const
Definition: gridfunc.cpp:311
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
Helper class for VisIt visualization data.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:944
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
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.
void PrintVTU(std::ostream &out, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
Definition: mesh.cpp:9884
bool own_data
Should the collection delete its mesh and fields.
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:971
virtual void SaveQField(const std::string &q_field_name)
Save one q-field, assuming the collection directory already exists.
T ZeroSubnormal(T val)
Definition: vector.hpp:473
MPI_Comm m_comm
Associated MPI communicator.
virtual void Load(int cycle_=0) override
Load the collection - not implemented in the ParaView writer.
double time
Physical time (for time-dependent simulations)
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
int GetNRanks() const
Definition: pmesh.hpp:277
iterator find(const std::string &fname)
Returns an iterator to the field fname.
void SetCompressionLevel(int compression_level_)
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:567
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:947
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
void SaveQFieldVTU(std::ostream &out, int ref, const QFieldMapIterator &it)
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.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
virtual void Save()
Save the collection to disk.
std::string to_padded_string(int i, int digits)
Convert an integer to a 0-padded string with the given number of digits.
Definition: text.hpp:54
int to_int(const std::string &str)
Convert a string to an int.
Definition: text.hpp:62
int GetCycle() const
Get time cycle (for time-dependent simulations)
void SaveDataVTU(std::ostream &out, int ref)
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:775
VTKFormat
Definition: vtk.hpp:61
void SaveGFieldVTU(std::ostream &out, int ref_, const FieldMapIterator &it)
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)
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1518
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.
IntegrationRule RefPts
Definition: geom.hpp:262
void SetHighOrderOutput(bool high_order_output_)
static const int pad_digits_default
Default value for pad_digits_*.
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition: geom.cpp:942
void SetCompression(bool compression_) override
virtual int GetRefinementLevelFromElems(Geometry::Type geom, int Npts)
Get the Refinement level based on number of elements.
Definition: geom.cpp:1477
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:629
int Dimension() const
Definition: mesh.hpp:911
const char * GetDataFormatString() const
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:348
int SpaceDimension() const
Definition: mesh.hpp:912
std::map< std::string, VisItFieldInfo > field_info_map
void UseRestartMode(bool restart_mode_)
int GetMyRank() const
Definition: pmesh.hpp:278
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:412
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:276
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
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:206
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.
const char * VTKByteOrder()
Definition: vtk.cpp:590
ParaViewDataCollection(const std::string &collection_name, mfem::Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has &#39;/&#39; at the end...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2388
void WriteVTKEncodedCompressed(std::ostream &out, const void *bytes, uint32_t nbytes, int compression_level)
Definition: vtk.cpp:547
void AppendBytes(std::vector< char > &vec, const T &val)
Append the binary representation of val to the byte buffer vec.
Definition: binaryio.hpp:55
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a quadrature function to the collection and update the root file.
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
void SetLevelsOfDetail(int levels_of_detail_)
std::string GetMeshFileName() const
Vector data type.
Definition: vector.hpp:60
void LoadVisItRootFile(const std::string &root_name)
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
virtual void Save() override
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)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
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:5525
Class for parallel meshes.
Definition: pmesh.hpp:32
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:734
int main()
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.