MFEM  v4.1.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-2020, 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 
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 #ifndef MFEM_USE_ZLIB
169  MFEM_VERIFY(!compression, "ZLib 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  mfem::ofgzstream mesh_file(mesh_name, compression);
232  mesh_file.precision(precision);
233 #ifdef MFEM_USE_MPI
234  const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
235  if (pmesh && format == PARALLEL_FORMAT)
236  {
237  pmesh->ParPrint(mesh_file);
238  }
239  else
240 #endif
241  {
242  mesh->Print(mesh_file);
243  }
244  if (!mesh_file)
245  {
246  error = WRITE_ERROR;
247  MFEM_WARNING("Error writing mesh to file: " << mesh_name);
248  }
249 }
250 
252 {
253  return (serial || format == SERIAL_FORMAT) ? "mesh" : "pmesh";
254 }
255 
257 {
259 }
260 
261 std::string DataCollection::GetFieldFileName(const std::string &field_name)
262 const
263 {
264  std::string dir_name = prefix_path + name;
265  if (cycle != -1)
266  {
267  dir_name += "_" + to_padded_string(cycle, pad_digits_cycle);
268  }
269  std::string file_name = dir_name + "/" + field_name;
271  {
272  file_name += "." + to_padded_string(myid, pad_digits_rank);
273  }
274  return file_name;
275 }
276 
278 {
279  mfem::ofgzstream field_file(GetFieldFileName(it->first), compression);
280 
281  field_file.precision(precision);
282  (it->second)->Save(field_file);
283  if (!field_file)
284  {
285  error = WRITE_ERROR;
286  MFEM_WARNING("Error writing field to file: " << it->first);
287  }
288 }
289 
291 {
292  mfem::ofgzstream q_field_file(GetFieldFileName(it->first), compression);
293 
294  q_field_file.precision(precision);
295  (it->second)->Save(q_field_file);
296  if (!q_field_file)
297  {
298  error = WRITE_ERROR;
299  MFEM_WARNING("Error writing q-field to file: " << it->first);
300  }
301 }
302 
303 void DataCollection::SaveField(const std::string &field_name)
304 {
305  FieldMapIterator it = field_map.find(field_name);
306  if (it != field_map.end())
307  {
308  SaveOneField(it);
309  }
310 }
311 
312 void DataCollection::SaveQField(const std::string &q_field_name)
313 {
314  QFieldMapIterator it = q_field_map.find(q_field_name);
315  if (it != q_field_map.end())
316  {
317  SaveOneQField(it);
318  }
319 }
320 
322 {
323  if (own_data) { delete mesh; }
324  mesh = NULL;
325 
328  own_data = false;
329 }
330 
332 {
333  DeleteData();
334  field_map.clear();
335  q_field_map.clear();
336 }
337 
339 {
340  DeleteData();
341 }
342 
343 
344 // class VisItDataCollection implementation
345 
347 {
348  if (mesh)
349  {
351  topo_dim = mesh->Dimension();
352  if (mesh->NURBSext)
353  {
356  }
357  }
358  else
359  {
360  spatial_dim = 0;
361  topo_dim = 0;
362  }
363 }
364 
365 VisItDataCollection::VisItDataCollection(const std::string& collection_name,
366  Mesh *mesh)
367  : DataCollection(collection_name, mesh)
368 {
369  appendRankToFileName = true; // always include rank in file names
370  cycle = 0; // always include cycle in directory names
371 
374 
375  UpdateMeshInfo();
376 }
377 
378 #ifdef MFEM_USE_MPI
380  const std::string& collection_name,
381  Mesh *mesh)
382  : DataCollection(collection_name, mesh)
383 {
384  m_comm = comm;
385  MPI_Comm_rank(comm, &myid);
386  MPI_Comm_size(comm, &num_procs);
387  appendRankToFileName = true; // always include rank in file names
388  cycle = 0; // always include cycle in directory names
389 
392 
393  UpdateMeshInfo();
394 }
395 #endif
396 
398 {
399  DataCollection::SetMesh(new_mesh);
400  appendRankToFileName = true;
401  UpdateMeshInfo();
402 }
403 
404 #ifdef MFEM_USE_MPI
405 void VisItDataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
406 {
407  // use VisItDataCollection's custom SetMesh, then set MPI info
408  SetMesh(new_mesh);
409  m_comm = comm;
410  MPI_Comm_rank(comm, &myid);
411  MPI_Comm_size(comm, &num_procs);
412 }
413 #endif
414 
415 void VisItDataCollection::RegisterField(const std::string& name,
416  GridFunction *gf)
417 {
419  field_info_map[name] = VisItFieldInfo("nodes", gf->VectorDim());
420 
421  int LOD = 1;
422  if (gf->FESpace()->GetNURBSext())
423  {
424  LOD = gf->FESpace()->GetNURBSext()->GetOrder();
425  }
426  else
427  {
428  for (int e=0; e<gf->FESpace()->GetNE(); e++)
429  {
430  LOD = std::max(LOD,gf->FESpace()->GetFE(e)->GetOrder());
431  }
432  }
433 
435 }
436 
437 void VisItDataCollection::SetLevelsOfDetail(int levels_of_detail)
438 {
439  visit_levels_of_detail = levels_of_detail;
440 }
441 
442 void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
443 {
444  visit_max_levels_of_detail = max_levels_of_detail;
445 }
446 
448 {
449  field_info_map.clear();
451 }
452 
454 {
456  SaveRootFile();
457 }
458 
460 {
461  if (myid != 0) { return; }
462 
463  std::string root_name = prefix_path + name + "_" +
465  ".mfem_root";
466  std::ofstream root_file(root_name.c_str());
467  root_file << GetVisItRootString();
468  if (!root_file)
469  {
470  error = WRITE_ERROR;
471  MFEM_WARNING("Error writing VisIt root file: " << root_name);
472  }
473 }
474 
476 {
477  DeleteAll();
478  time_step = 0.0;
479  error = NO_ERROR;
480  cycle = cycle_;
481  std::string root_name = prefix_path + name + "_" +
483  ".mfem_root";
484  LoadVisItRootFile(root_name);
485  if (format != SERIAL_FORMAT || num_procs > 1)
486  {
487 #ifndef MFEM_USE_MPI
488  MFEM_WARNING("Cannot load parallel VisIt root file in serial.");
489  error = READ_ERROR;
490 #else
491  if (m_comm == MPI_COMM_NULL)
492  {
493  MFEM_WARNING("Cannot load parallel VisIt root file without MPI"
494  " communicator");
495  error = READ_ERROR;
496  }
497  else
498  {
499  // num_procs was read from the root file, check for consistency with
500  // the associated MPI_Comm, m_comm:
501  int comm_size;
502  MPI_Comm_size(m_comm, &comm_size);
503  if (comm_size != num_procs)
504  {
505  MFEM_WARNING("Processor number mismatch: VisIt root file: "
506  << num_procs << ", MPI_comm: " << comm_size);
507  error = READ_ERROR;
508  }
509  else
510  {
511  // myid was set when setting m_comm
512  }
513  }
514 #endif
515  }
516  if (!error)
517  {
518  LoadMesh(); // sets own_data to true, when there is no error
519  }
520  if (!error)
521  {
522  LoadFields();
523  }
524  if (error)
525  {
526  DeleteAll();
527  }
528 }
529 
530 void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
531 {
532  std::ifstream root_file(root_name.c_str());
533  std::stringstream buffer;
534  buffer << root_file.rdbuf();
535  if (!buffer)
536  {
537  error = READ_ERROR;
538  MFEM_WARNING("Error reading the VisIt root file: " << root_name);
539  }
540  else
541  {
542  ParseVisItRootString(buffer.str());
543  }
544 }
545 
547 {
548  std::string mesh_fname = GetMeshFileName();
549  named_ifgzstream file(mesh_fname);
550  // TODO: in parallel, check for errors on all processors
551  if (!file)
552  {
553  error = READ_ERROR;
554  MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
555  return;
556  }
557  // TODO: 1) load parallel mesh on one processor
558  if (format == SERIAL_FORMAT)
559  {
560  mesh = new Mesh(file, 1, 0, false);
561  serial = true;
562  }
563  else
564  {
565 #ifdef MFEM_USE_MPI
566  mesh = new ParMesh(m_comm, file);
567  serial = false;
568 #else
569  error = READ_ERROR;
570  MFEM_WARNING("Reading parallel format in serial is not supported");
571  return;
572 #endif
573  }
575  topo_dim = mesh->Dimension();
576  own_data = true;
577 }
578 
580 {
581  std::string path_left = prefix_path + name + "_" +
583  std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
584 
585  field_map.clear();
586  for (FieldInfoMapIterator it = field_info_map.begin();
587  it != field_info_map.end(); ++it)
588  {
589  std::string fname = path_left + it->first + path_right;
590  mfem::ifgzstream file(fname);
591  // TODO: in parallel, check for errors on all processors
592  if (!file)
593  {
594  error = READ_ERROR;
595  MFEM_WARNING("Unable to open field file: " << fname);
596  return;
597  }
598  // TODO: 1) load parallel GridFunction on one processor
599  if (serial)
600  {
601  field_map.Register(it->first, new GridFunction(mesh, file), own_data);
602  }
603  else
604  {
605 #ifdef MFEM_USE_MPI
607  it->first,
608  new ParGridFunction(dynamic_cast<ParMesh*>(mesh), file), own_data);
609 #else
610  error = READ_ERROR;
611  MFEM_WARNING("Reading parallel format in serial is not supported");
612  return;
613 #endif
614  }
615  }
616 }
617 
619 {
620  // Get the path string (relative to where the root file is, i.e. no prefix).
621  std::string path_str =
623 
624  // We have to build the json tree inside out to get all the values in there
625  picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
626 
627  // Build the mesh data
628  std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
629  mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
630  mtags["topo_dim"] = picojson::value(to_string(topo_dim));
631  mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
632  mesh["path"] = picojson::value(path_str + GetMeshShortFileName() +
633  file_ext_format);
634  mesh["tags"] = picojson::value(mtags);
635  mesh["format"] = picojson::value(to_string(format));
636 
637  // Build the fields data entries
638  for (FieldInfoMapIterator it = field_info_map.begin();
639  it != field_info_map.end(); ++it)
640  {
641  ftags["assoc"] = picojson::value((it->second).association);
642  ftags["comps"] = picojson::value(to_string((it->second).num_components));
643  ftags["lod"] = picojson::value(to_string(visit_levels_of_detail));
644  field["path"] = picojson::value(path_str + it->first + file_ext_format);
645  field["tags"] = picojson::value(ftags);
646  fields[it->first] = picojson::value(field);
647  }
648 
649  main["cycle"] = picojson::value(double(cycle));
650  main["time"] = picojson::value(time);
651  main["time_step"] = picojson::value(time_step);
652  main["domains"] = picojson::value(double(num_procs));
653  main["mesh"] = picojson::value(mesh);
654  if (!field_info_map.empty())
655  {
656  main["fields"] = picojson::value(fields);
657  }
658 
659  dsets["main"] = picojson::value(main);
660  top["dsets"] = picojson::value(dsets);
661 
662  return picojson::value(top).serialize(true);
663 }
664 
665 void VisItDataCollection::ParseVisItRootString(const std::string& json)
666 {
667  picojson::value top, dsets, main, mesh, fields;
668  std::string parse_err = picojson::parse(top, json);
669  if (!parse_err.empty())
670  {
671  error = READ_ERROR;
672  MFEM_WARNING("Unable to parse VisIt root data.");
673  return;
674  }
675 
676  // Process "main"
677  dsets = top.get("dsets");
678  main = dsets.get("main");
679  cycle = int(main.get("cycle").get<double>());
680  time = main.get("time").get<double>();
681  if (main.contains("time_step"))
682  {
683  time_step = main.get("time_step").get<double>();
684  }
685  num_procs = int(main.get("domains").get<double>());
686  mesh = main.get("mesh");
687  fields = main.get("fields");
688 
689  // ... Process "mesh"
690 
691  // Set the DataCollection::name using the mesh path
692  std::string path = mesh.get("path").get<std::string>();
693  size_t right_sep = path.find('_');
694  if (right_sep == std::string::npos)
695  {
696  error = READ_ERROR;
697  MFEM_WARNING("Unable to parse VisIt root data.");
698  return;
699  }
700  name = path.substr(0, right_sep);
701 
702  if (mesh.contains("format"))
703  {
704  format = to_int(mesh.get("format").get<std::string>());
705  }
706  spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
707  topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
709  to_int(mesh.get("tags").get("max_lods").get<std::string>());
710 
711  // ... Process "fields"
712  field_info_map.clear();
713  if (fields.is<picojson::object>())
714  {
715  picojson::object fields_obj = fields.get<picojson::object>();
716  for (picojson::object::iterator it = fields_obj.begin();
717  it != fields_obj.end(); ++it)
718  {
719  picojson::value tags = it->second.get("tags");
720  field_info_map[it->first] =
721  VisItFieldInfo(tags.get("assoc").get<std::string>(),
722  to_int(tags.get("comps").get<std::string>()));
723  }
724  }
725 }
726 
728  collection_name,
729  Mesh *mesh_)
730  : DataCollection(collection_name, mesh_),
731  levels_of_detail(1),
732  pv_data_format(VTKFormat::BINARY),
733  high_order_output(false)
734 {
735 #ifdef MFEM_USE_ZLIB
736  compression = -1; // default zlib compression level, equivalent to 6
737 #else
738  compression = 0;
739 #endif
740 }
741 
742 void ParaViewDataCollection::RegisterField(const std::string& field_name,
743  mfem::GridFunction *gf)
744 {
745  DataCollection::RegisterField(field_name,gf);
746 }
747 
748 void ParaViewDataCollection::SetLevelsOfDetail(int levels_of_detail_)
749 {
750  levels_of_detail = levels_of_detail_;
751 }
752 
754 {
755  MFEM_WARNING("ParaViewDataCollection::Load() is not implemented!");
756 }
757 
759 {
760  std::string out = "";
762  return out;
763 }
764 
766 {
767  std::string out = "Cycle" + to_padded_string(cycle,pad_digits_cycle);
768  return out;
769 }
770 
772 {
773  std::string out = GeneratePVTUPath();
774  return out;
775 }
776 
778 {
779  std::string out = GetCollectionName()+".pvd";
780  return out;
781 }
782 
784 {
785  std::string out = "data.pvtu";
786  return out;
787 }
788 
790 {
791  std::string out = "proc" + to_padded_string(myid,pad_digits_rank)+".vtu";
792  return out;
793 }
795 {
796  std::string out = "proc" + to_padded_string(crank,pad_digits_rank)+".vtu";
797  return out;
798 }
799 
801 {
802  // add a new collection to the PDV file
803 
804  // check if the directories are created
805  {
806  std::string path = GenerateCollectionPath()+"/"+GenerateVTUPath();
807  int err = create_directory(path, mesh, myid);
808  if (err)
809  {
810  error = WRITE_ERROR;
811  MFEM_WARNING("Error creating directory: " << path);
812  return; // do not even try to write the mesh
813  }
814  }
815  // the directory is created
816 
817  // create pvd file if needed
818  if (!pvd_stream.is_open())
819  {
820  std::string dpath=GenerateCollectionPath();
821  std::string pvdname=dpath+"/"+GeneratePVDFileName();
822  pvd_stream.open(pvdname.c_str(),std::ios::out);
823  // initialize the file
824  pvd_stream << "<?xml version=\"1.0\"?>\n";
825  pvd_stream << "<VTKFile type=\"Collection\" version=\"0.1\"";
826  pvd_stream << " byte_order=\"" << VTKByteOrder() << "\">\n";
827  pvd_stream << "<Collection>" << std::endl;
828  }
829 
830  // define the vtu file
831  {
832  std::string fname = GenerateCollectionPath()+"/"+GenerateVTUPath()+"/"
834  std::fstream out(fname.c_str(), std::ios::out);
835  out.precision(precision);
836  SaveDataVTU(out,levels_of_detail);
837  out.close();
838  }
839 
840  // define the pvtu file only on process 0
841  if (myid==0)
842  {
843  std::string fname = GenerateCollectionPath()+"/"+GeneratePVTUPath()+"/"
845  std::fstream out(fname.c_str(), std::ios::out);
846 
847  out << "<?xml version=\"1.0\"?>\n";
848  out << "<VTKFile type=\"PUnstructuredGrid\"";
849  out << " version =\"0.1\" byte_order=\"" << VTKByteOrder() << "\">\n";
850  out << "<PUnstructuredGrid GhostLevel=\"0\">\n";
851 
852  out << "<PPoints>\n";
853  out << "\t<PDataArray type=\"" << GetDataTypeString() << "\" ";
854  out << " Name=\"Points\" NumberOfComponents=\"3\""
855  << " format=\"" << GetDataFormatString() << "\"/>\n";
856  out << "</PPoints>\n";
857 
858  out << "<PCells>\n";
859  out << "\t<PDataArray type=\"Int32\" ";
860  out << " Name=\"connectivity\" NumberOfComponents=\"1\""
861  << " format=\"" << GetDataFormatString() << "\"/>\n";
862  out << "\t<PDataArray type=\"Int32\" ";
863  out << " Name=\"offsets\" NumberOfComponents=\"1\""
864  << " format=\"" << GetDataFormatString() << "\"/>\n";
865  out << "\t<PDataArray type=\"UInt8\" ";
866  out << " Name=\"types\" NumberOfComponents=\"1\""
867  << " format=\"" << GetDataFormatString() << "\"/>\n";
868  out << "</PCells>\n";
869 
870  out << "<PPointData>\n";
871  for (FieldMapIterator it=field_map.begin(); it!=field_map.end(); ++it)
872  {
873  int vec_dim=it->second->VectorDim();
874  out << "<PDataArray type=\"" << GetDataTypeString()
875  << "\" Name=\"" << it->first
876  << "\" NumberOfComponents=\"" << vec_dim << "\" "
877  << "format=\"" << GetDataFormatString() << "\" />\n";
878  }
879  out << "</PPointData>\n";
880 
881  // CELL DATA
882  out << "<PCellData>\n";
883  out << "\t<PDataArray type=\"Int32\" Name=\"" << "material"
884  << "\" NumberOfComponents=\"1\""
885  << " format=\"" << GetDataFormatString() << "\"/>\n";
886  out << "</PCellData>\n";
887 
888  for (int ii=0; ii<num_procs; ii++)
889  {
890  // this one is generated without the path
891  std::string nfname=GenerateVTUFileName(ii);
892  out << "<Piece Source=\"" << nfname << "\"/>\n";
893  }
894  out << "</PUnstructuredGrid>\n";
895  out << "</VTKFile>\n";
896  out.close();
897 
898  fname = GeneratePVTUPath()+"/"+GeneratePVTUFileName();
899  // add the pvtu file to the pvd_stream
900  pvd_stream << "<DataSet timestep=\"" << GetTime(); // GetCycle();
901  pvd_stream << "\" group=\"\" part=\"" << 0 << "\" file=\"";
902  pvd_stream << fname << "\"/>\n";
903  std::fstream::pos_type pos = pvd_stream.tellp();
904  pvd_stream << "</Collection>\n";
905  pvd_stream << "</VTKFile>" << std::endl;
906  pvd_stream.seekp(pos);
907  }
908 }
909 
910 void ParaViewDataCollection::SaveDataVTU(std::ostream &out, int ref)
911 {
912  out << "<VTKFile type=\"UnstructuredGrid\"";
913  if (compression != 0)
914  {
915  out << " compressor=\"vtkZLibDataCompressor\"";
916  }
917  out << " version=\"0.1\" byte_order=\"" << VTKByteOrder() << "\">\n";
918  out << "<UnstructuredGrid>\n";
919  mesh->PrintVTU(out,ref,pv_data_format,high_order_output,compression);
920 
921  // dump out the grid functions as point data
922  out << "<PointData >\n";
923  // save the grid functions
924  // iterate over all grid functions
925  for (FieldMapIterator it=field_map.begin(); it!=field_map.end(); ++it)
926  {
927  SaveGFieldVTU(out,ref,it);
928  }
929  // iterate over all quadrature functions
930  // if the Quadrature functions are dumped as cell data
931  // the cycle should be moved before the grid functions
932  // and the PrintVTU CellData section should be open in the mesh dump
933  for (QFieldMapIterator it=q_field_map.begin(); it!=q_field_map.end(); ++it)
934  {
935  // save the quadrature functions
936  // this one is not implemented yet
937  SaveQFieldVTU(out,ref,it);
938  }
939  out << "</PointData>\n";
940  // close the mesh
941  out << "</Piece>\n"; // close the piece open in the PrintVTU method
942  out << "</UnstructuredGrid>\n";
943  out << "</VTKFile>" << std::endl;
944 }
945 
946 void ParaViewDataCollection::SaveQFieldVTU(std::ostream &out, int ref,
947  const QFieldMapIterator& it )
948 {
949  MFEM_WARNING("SaveQFieldVTU is not currently implemented - field name:"<<it->second);
950 }
951 
952 void ParaViewDataCollection::SaveGFieldVTU(std::ostream &out, int ref_,
953  const FieldMapIterator& it)
954 {
955  RefinedGeometry *RefG;
956  Vector val;
957  DenseMatrix vval, pmat;
958  std::vector<char> buf;
959  int vec_dim = it->second->VectorDim();
960  if (vec_dim == 1)
961  {
962  // scalar data
963  out << "<DataArray type=\"" << GetDataTypeString()
964  << "\" Name=\"" << it->first;
965  out << "\" NumberOfComponents=\"1\" format=\""
966  << GetDataFormatString() << "\" >\n";
967  for (int i = 0; i < mesh->GetNE(); i++)
968  {
970  mesh->GetElementBaseGeometry(i), ref_, 1);
971  it->second->GetValues(i, RefG->RefPts, val, pmat);
972  for (int j = 0; j < val.Size(); j++)
973  {
974  if (pv_data_format == VTKFormat::ASCII)
975  {
976  out << val(j) << '\n';
977  }
978  else if (pv_data_format == VTKFormat::BINARY)
979  {
980  bin_io::AppendBytes(buf, val(j));
981  }
982  else
983  {
984  bin_io::AppendBytes<float>(buf, float(val(j)));
985  }
986  }
987  }
988  }
989  else
990  {
991  // vector data
992  out << "<DataArray type=\"" << GetDataTypeString()
993  << "\" Name=\"" << it->first;
994  out << "\" NumberOfComponents=\"" << vec_dim << "\""
995  << " format=\"" << GetDataFormatString() << "\" >" << '\n';
996  for (int i = 0; i < mesh->GetNE(); i++)
997  {
999  mesh->GetElementBaseGeometry(i), ref_, 1);
1000 
1001  it->second->GetVectorValues(i, RefG->RefPts, vval, pmat);
1002 
1003  for (int jj = 0; jj < vval.Width(); jj++)
1004  {
1005  for (int ii = 0; ii < vval.Height(); ii++)
1006  {
1007  if (pv_data_format == VTKFormat::ASCII)
1008  {
1009  out << vval(ii,jj) << ' ';
1010  }
1011  else if (pv_data_format == VTKFormat::BINARY)
1012  {
1013  bin_io::AppendBytes(buf, vval(ii,jj));
1014  }
1015  else
1016  {
1017  bin_io::AppendBytes<float>(buf, float(vval(ii,jj)));
1018  }
1019  }
1020  if (pv_data_format == VTKFormat::ASCII) { out << '\n'; }
1021  }
1022  }
1023  }
1024 
1025  if (IsBinaryFormat())
1026  {
1027  WriteVTKEncodedCompressed(out,buf.data(),buf.size(),compression);
1028  out << '\n';
1029  }
1030  out << "</DataArray>" << std::endl;
1031 }
1032 
1034 {
1035  pv_data_format = fmt;
1036 }
1037 
1039 {
1040  return pv_data_format != VTKFormat::ASCII;
1041 }
1042 
1043 void ParaViewDataCollection::SetHighOrderOutput(bool high_order_output_)
1044 {
1045  high_order_output = high_order_output_;
1046 }
1047 
1049 {
1050  MFEM_ASSERT(compression_level_ >= -1 && compression_level_ <= 9,
1051  "Compression level must be between -1 and 9 (inclusive).");
1052  compression = compression_level_;
1053 }
1054 
1056 {
1057  // If we are enabling compression, and it was disabled previously, use the
1058  // default compression level. Otherwise, leave the compression level
1059  // unchanged.
1060  if (compression_ && compression == 0)
1061  {
1062  SetCompressionLevel(-1);
1063  }
1064 }
1065 
1067 {
1068  if (pv_data_format == VTKFormat::ASCII)
1069  {
1070  return "ascii";
1071  }
1072  else
1073  {
1074  return "binary";
1075  }
1076 }
1077 
1079 {
1080  if (pv_data_format==VTKFormat::ASCII || pv_data_format==VTKFormat::BINARY)
1081  {
1082  return "Float64";
1083  }
1084  else
1085  {
1086  return "Float32";
1087  }
1088 }
1089 
1090 } // 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:1188
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.
void SetDataFormat(VTKFormat fmt)
double GetTime() const
Get physical time (for time-dependent simulations)
void PrintVTU(std::ostream &out, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0)
Definition: mesh.cpp:8714
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
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:63
const std::string & GetCollectionName() const
Get the name of the collection.
virtual void RegisterField(const std::string &field_name, mfem::GridFunction *gf) override
Add a grid function to 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:318
int VectorDim() const
Definition: gridfunc.cpp:302
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
Helper class for VisIt visualization data.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
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.
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
int main(int argc, char *argv[])
Definition: ex1.cpp:62
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:790
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.
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:231
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:408
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:57
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
void SaveDataVTU(std::ostream &out, int ref)
VTKFormat
Definition: vtk.hpp:22
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:1395
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:239
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:925
void SetCompression(bool compression_) override
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:442
int Dimension() const
Definition: mesh.hpp:744
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:745
std::map< std::string, VisItFieldInfo > field_info_map
int GetMyRank() const
Definition: pmesh.hpp:232
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:297
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:230
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:191
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:461
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...
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
void WriteVTKEncodedCompressed(std::ostream &out, const void *bytes, uint32_t nbytes, int compression_level)
Definition: vtk.cpp:418
void AppendBytes(std::vector< char > &vec, const T &val)
Definition: binaryio.hpp:43
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:48
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:5218
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.