MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
datacollection.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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 "../mesh/vtk.hpp"
15#include "../mesh/vtkhdf.hpp"
17#include "../general/text.hpp"
18#include "picojson.h"
19
20#include <cerrno> // errno
21#include <sstream>
22#include <regex>
23
24#ifndef _WIN32
25#include <sys/stat.h> // mkdir
26#else
27#include <direct.h> // _mkdir
28#define mkdir(dir, mode) _mkdir(dir)
29#endif
30
31namespace mfem
32{
33
34// static method
35int DataCollection::create_directory(const std::string &dir_name,
36 const Mesh *mesh, int myid)
37{
38 // create directories recursively
39 const char path_delim = '/';
40 std::string::size_type pos = 0;
41 int err_flag;
42#ifdef MFEM_USE_MPI
43 const ParMesh *pmesh = dynamic_cast<const ParMesh*>(mesh);
44#endif
45
46 do
47 {
48 pos = dir_name.find(path_delim, pos+1);
49 std::string subdir = dir_name.substr(0, pos);
50
51#ifndef MFEM_USE_MPI
52 err_flag = mkdir(subdir.c_str(), 0777);
53 err_flag = (err_flag && (errno != EEXIST)) ? 1 : 0;
54#else
55 if (myid == 0 || pmesh == NULL)
56 {
57 err_flag = mkdir(subdir.c_str(), 0777);
58 err_flag = (err_flag && (errno != EEXIST)) ? 1 : 0;
59 }
60#endif
61 }
62 while ( pos != std::string::npos );
63
64#ifdef MFEM_USE_MPI
65 if (pmesh)
66 {
67 MPI_Bcast(&err_flag, 1, MPI_INT, 0, pmesh->GetComm());
68 }
69#endif
70
71 return err_flag;
72}
73
74// class DataCollection implementation
75
76DataCollection::DataCollection(const std::string& collection_name, Mesh *mesh_)
77{
78 std::string::size_type pos = collection_name.find_last_of('/');
79 if (pos == std::string::npos)
80 {
81 name = collection_name;
82 // leave prefix_path empty
83 }
84 else
85 {
86 prefix_path = collection_name.substr(0, pos+1);
87 name = collection_name.substr(pos+1);
88 }
89 mesh = mesh_;
90 myid = 0;
91 num_procs = 1;
92 serial = true;
94
95#ifdef MFEM_USE_MPI
96 m_comm = MPI_COMM_NULL;
97 ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
98 if (par_mesh)
99 {
100 myid = par_mesh->GetMyRank();
101 num_procs = par_mesh->GetNRanks();
102 m_comm = par_mesh->GetComm();
103 serial = false;
105 }
106#endif
107 own_data = false;
108 cycle = -1;
109 time = 0.0;
110 time_step = 0.0;
113 format = SERIAL_FORMAT; // use serial mesh format
114 compression = 0;
115 error = No_Error;
116}
117
119{
120 if (own_data && new_mesh != mesh) { delete mesh; }
121 mesh = new_mesh;
122 myid = 0;
123 num_procs = 1;
124 serial = true;
125 appendRankToFileName = false;
126
127#ifdef MFEM_USE_MPI
128 m_comm = MPI_COMM_NULL;
129 ParMesh *par_mesh = dynamic_cast<ParMesh*>(mesh);
130 if (par_mesh)
131 {
132 myid = par_mesh->GetMyRank();
133 num_procs = par_mesh->GetNRanks();
134 m_comm = par_mesh->GetComm();
135 serial = false;
137 }
138#endif
139}
140
141#ifdef MFEM_USE_MPI
142void DataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
143{
144 // This seems to be the cleanest way to accomplish this
145 // and avoid duplicating fine grained details:
146
147 SetMesh(new_mesh);
148
149 m_comm = comm;
150 MPI_Comm_rank(comm, &myid);
151 MPI_Comm_size(comm, &num_procs);
152}
153#endif
154
156{
157 switch (fmt)
158 {
159 case SERIAL_FORMAT: break;
160#ifdef MFEM_USE_MPI
161 case PARALLEL_FORMAT: break;
162#endif
163 default: MFEM_ABORT("unknown format: " << fmt);
164 }
165 format = fmt;
166}
167
169{
170 compression = comp;
171#ifndef MFEM_USE_ZLIB
172 MFEM_VERIFY(!compression, "ZLib not enabled in MFEM build.");
173#endif
174}
175
176void DataCollection::SetPrefixPath(const std::string& prefix)
177{
178 if (!prefix.empty())
179 {
180 prefix_path = prefix;
181 if (!prefix_path.empty() && prefix_path[prefix_path.size()-1] != '/')
182 {
183 prefix_path += '/';
184 }
185 }
186 else
187 {
188 prefix_path.clear();
189 }
190}
191
192void DataCollection::Load(int cycle_)
193{
194 MFEM_ABORT("this method is not implemented");
195}
196
198{
199 SaveMesh();
200
201 if (error) { return; }
202
203 for (FieldMapIterator it = field_map.begin(); it != field_map.end(); ++it)
204 {
205 SaveOneField(it);
206 // Even if there is an error, try saving the other fields
207 }
208
210 ++it)
211 {
212 SaveOneQField(it);
213 }
214}
215
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 int error_code = create_directory(dir_name, mesh, myid);
224 if (error_code)
225 {
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 {
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
262std::string DataCollection::GetFieldFileName(const std::string &field_name)
263const
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 {
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 {
300 MFEM_WARNING("Error writing q-field to file: " << it->first);
301 }
302}
303
304void 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
313void DataCollection::SaveQField(const std::string &field_name)
314{
315 QFieldMapIterator it = q_field_map.find(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
338
343
344
345// class VisItDataCollection implementation
346
348{
349 if (mesh)
350 {
353 if (mesh->NURBSext)
354 {
357 }
358 }
359 else
360 {
361 spatial_dim = 0;
362 topo_dim = 0;
363 }
364}
365
366VisItDataCollection::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
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
395}
396#endif
397
399{
400 DataCollection::SetMesh(new_mesh);
403}
404
405#ifdef MFEM_USE_MPI
406void 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
416void 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 gf->FESpace()->FEColl()->Name(),
435 gf->FESpace()->FEColl()->GetOrder());
437}
438
439void VisItDataCollection::RegisterQField(const std::string& name,
441{
442 int LOD = -1;
443 Mesh *mesh = qf->GetSpace()->GetMesh();
444 for (int e=0; e<qf->GetSpace()->GetNE(); e++)
445 {
448 qf->GetIntRule(e).GetNPoints());
449
450 LOD = std::max(LOD,locLOD);
451 }
452
454 // For quadrature functions, use basis pattern:
455 // QF_{ORDER}_{VDIM}
456 int qf_vdim = qf->GetVDim();
457 int qf_order = qf->GetSpace()->GetOrder();
458 std::ostringstream oss;
459 oss << "QF_" << qf_order << "_" << qf_vdim;
460 field_info_map[name] = VisItFieldInfo("quadrature", qf->GetVDim(), LOD,
461 oss.str(), qf_order);
463}
464
466{
467 visit_levels_of_detail = levels_of_detail;
468}
469
470void VisItDataCollection::SetMaxLevelsOfDetail(int max_levels_of_detail)
471{
472 visit_max_levels_of_detail = max_levels_of_detail;
473}
474
480
486
488{
489 if (myid != 0) { return; }
490
491 std::string root_name = prefix_path + name + "_" +
493 ".mfem_root";
494 std::ofstream root_file(root_name);
495 root_file << GetVisItRootString();
496 if (!root_file)
497 {
499 MFEM_WARNING("Error writing VisIt root file: " << root_name);
500 }
501}
502
504{
505 DeleteAll();
506 time_step = 0.0;
507 error = No_Error;
508 cycle = cycle_;
509 std::string root_name = prefix_path + name + "_" +
511 ".mfem_root";
512 LoadVisItRootFile(root_name);
513 if (format != SERIAL_FORMAT || num_procs > 1)
514 {
515#ifndef MFEM_USE_MPI
516 MFEM_WARNING("Cannot load parallel VisIt root file in serial.");
518#else
519 if (m_comm == MPI_COMM_NULL)
520 {
521 MFEM_WARNING("Cannot load parallel VisIt root file without MPI"
522 " communicator");
524 }
525 else
526 {
527 // num_procs was read from the root file, check for consistency with
528 // the associated MPI_Comm, m_comm:
529 int comm_size;
530 MPI_Comm_size(m_comm, &comm_size);
531 if (comm_size != num_procs)
532 {
533 MFEM_WARNING("Processor number mismatch: VisIt root file: "
534 << num_procs << ", MPI_comm: " << comm_size);
536 }
537 else
538 {
539 // myid was set when setting m_comm
540 }
541 }
542#endif
543 }
544 if (!error)
545 {
546 LoadMesh(); // sets own_data to true, when there is no error
547 }
548 if (!error)
549 {
550 LoadFields();
551 }
552 if (error)
553 {
554 DeleteAll();
555 }
556}
557
558void VisItDataCollection::LoadVisItRootFile(const std::string& root_name)
559{
560 std::ifstream root_file(root_name);
561 std::stringstream buffer;
562 buffer << root_file.rdbuf();
563 if (!buffer)
564 {
566 MFEM_WARNING("Error reading the VisIt root file: " << root_name);
567 }
568 else
569 {
570 ParseVisItRootString(buffer.str());
571 }
572}
573
575{
576 // GetMeshFileName() uses 'serial', so we need to set it in advance.
578 std::string mesh_fname = GetMeshFileName();
579 named_ifgzstream file(mesh_fname);
580 // TODO: in parallel, check for errors on all processors
581 if (!file)
582 {
584 MFEM_WARNING("Unable to open mesh file: " << mesh_fname);
585 return;
586 }
587 // TODO: 1) load parallel mesh on one processor
588 if (format == SERIAL_FORMAT)
589 {
590 mesh = new Mesh(file, 1, 0, false);
591 serial = true;
592 }
593 else
594 {
595#ifdef MFEM_USE_MPI
596 mesh = new ParMesh(m_comm, file);
597 serial = false;
598#else
600 MFEM_WARNING("Reading parallel format in serial is not supported");
601 return;
602#endif
603 }
606 own_data = true;
607}
608
610{
611 std::string path_left = prefix_path + name + "_" +
613 std::string path_right = "." + to_padded_string(myid, pad_digits_rank);
614
616 for (FieldInfoMapIterator it = field_info_map.begin();
617 it != field_info_map.end(); ++it)
618 {
619 std::string fname = path_left + it->first + path_right;
620 mfem::ifgzstream file(fname);
621 // TODO: in parallel, check for errors on all processors
622 if (!file)
623 {
625 MFEM_WARNING("Unable to open field file: " << fname);
626 return;
627 }
628 // TODO: 1) load parallel GridFunction on one processor
629 if (serial)
630 {
631 if ((it->second).association == "nodes")
632 {
633 field_map.Register(it->first, new GridFunction(mesh, file), own_data);
634 }
635 else if ((it->second).association == "elements" || // old style
636 (it->second).association == "quadrature") // new style
637 {
638 q_field_map.Register(it->first, new QuadratureFunction(mesh, file), own_data);
639 }
640 }
641 else
642 {
643#ifdef MFEM_USE_MPI
644 if ((it->second).association == "nodes")
645 {
647 it->first,
648 new ParGridFunction(dynamic_cast<ParMesh*>(mesh), file), own_data);
649 }
650 else if ((it->second).association == "elements" || // old style
651 (it->second).association == "quadrature") // new style
652 {
653 q_field_map.Register(it->first, new QuadratureFunction(mesh, file), own_data);
654 }
655#else
657 MFEM_WARNING("Reading parallel format in serial is not supported");
658 return;
659#endif
660 }
661 }
662}
663
665{
666 // Get the path string (relative to where the root file is, i.e. no prefix).
667 std::string path_str =
669
670 // We have to build the json tree inside out to get all the values in there
671 picojson::object top, dsets, main, mesh, fields, field, mtags, ftags;
672
673 // Build the mesh data
674 std::string file_ext_format = ".%0" + to_string(pad_digits_rank) + "d";
675 mtags["spatial_dim"] = picojson::value(to_string(spatial_dim));
676 mtags["topo_dim"] = picojson::value(to_string(topo_dim));
677 mtags["max_lods"] = picojson::value(to_string(visit_max_levels_of_detail));
678 mesh["path"] = picojson::value(path_str + GetMeshShortFileName() +
679 file_ext_format);
680 mesh["tags"] = picojson::value(mtags);
681 mesh["format"] = picojson::value(to_string(format));
682
683 // Build the fields data entries
684 for (FieldInfoMapIterator it = field_info_map.begin();
685 it != field_info_map.end(); ++it)
686 {
687 ftags["assoc"] = picojson::value((it->second).association);
688 ftags["comps"] = picojson::value(to_string((it->second).num_components));
689 ftags["lod"] = picojson::value(to_string((it->second).lod));
690 ftags["basis"] = picojson::value((it->second).basis);
691 ftags["order"] = picojson::value(to_string((it->second).order));
692 field["path"] = picojson::value(path_str + it->first + file_ext_format);
693 field["tags"] = picojson::value(ftags);
694 fields[it->first] = picojson::value(field);
695 }
696
697 main["cycle"] = picojson::value(double(cycle));
698 main["time"] = picojson::value(time);
699 main["time_step"] = picojson::value(time_step);
700 main["domains"] = picojson::value(double(num_procs));
701 main["mesh"] = picojson::value(mesh);
702 if (!field_info_map.empty())
703 {
704 main["fields"] = picojson::value(fields);
705 }
706
707 dsets["main"] = picojson::value(main);
708 top["dsets"] = picojson::value(dsets);
709
710 return picojson::value(top).serialize(true);
711}
712
713void VisItDataCollection::ParseVisItRootString(const std::string& json)
714{
715 picojson::value top, dsets, main, mesh, fields;
716 std::string parse_err = picojson::parse(top, json);
717 if (!parse_err.empty())
718 {
720 MFEM_WARNING("Unable to parse VisIt root data.");
721 return;
722 }
723
724 // Process "main"
725 dsets = top.get("dsets");
726 main = dsets.get("main");
727 cycle = int(main.get("cycle").get<double>());
728 time = main.get("time").get<double>();
729 if (main.contains("time_step"))
730 {
731 time_step = main.get("time_step").get<double>();
732 }
733 num_procs = int(main.get("domains").get<double>());
734 mesh = main.get("mesh");
735 fields = main.get("fields");
736
737 // ... Process "mesh"
738
739 // Set the DataCollection::name using the mesh path
740 std::string path = mesh.get("path").get<std::string>();
741 size_t right_sep = path.rfind('_');
742 if (right_sep == std::string::npos)
743 {
745 MFEM_WARNING("Unable to parse VisIt root data.");
746 return;
747 }
748 name = path.substr(0, right_sep);
749
750 if (mesh.contains("format"))
751 {
752 format = to_int(mesh.get("format").get<std::string>());
753 }
754 spatial_dim = to_int(mesh.get("tags").get("spatial_dim").get<std::string>());
755 topo_dim = to_int(mesh.get("tags").get("topo_dim").get<std::string>());
757 to_int(mesh.get("tags").get("max_lods").get<std::string>());
758
759 // ... Process "fields"
760 field_info_map.clear();
761 if (fields.is<picojson::object>())
762 {
763 picojson::object fields_obj = fields.get<picojson::object>();
764 for (picojson::object::iterator it = fields_obj.begin();
765 it != fields_obj.end(); ++it)
766 {
767 picojson::value tags = it->second.get("tags");
768
769 // defaults that allow us to parse older mfem_root files
770 int lod = 1;
771 std::string basis = "";
772 int order = -1;
773
774 if (tags.contains("lod"))
775 {
776 lod = to_int(tags.get("lod").get<std::string>());
777 }
778
779 if (tags.contains("basis"))
780 {
781 basis = tags.get("comps").get<std::string>();
782 }
783
784 if (tags.contains("order"))
785 {
786 order = to_int(tags.get("comps").get<std::string>());
787 }
788
789 field_info_map[it->first] =
790 VisItFieldInfo(tags.get("assoc").get<std::string>(),
791 to_int(tags.get("comps").get<std::string>()),
792 lod, basis, order);
793 }
794 }
795}
796
798 const std::string &name, Mesh *mesh) : DataCollection(name, mesh)
799{
800 cycle = 0;
801#ifdef MFEM_USE_ZLIB
802 // If we have zlib, enable compression. Otherwise, compression is disabled in
803 // the DataCollection base class constructor.
804 compression = true;
805#endif
806}
807
809{
810 levels_of_detail = levels_of_detail_;
811}
812
814{
815 high_order_output = high_order_output_;
816}
817
819{
820 bdr_output = bdr_output_;
821}
822
824{
825 MFEM_ASSERT(compression_level_ >= -1 && compression_level_ <= 9,
826 "Compression level must be between -1 and 9 (inclusive).");
827 if (compression_level_ != 0) { SetCompression(true);}
828 compression_level = compression_level_;
829}
830
835
840
845
847{
848 restart_mode = restart_mode_;
849}
850
852 const std::string& collection_name, Mesh *mesh_)
853 : ParaViewDataCollectionBase(collection_name, mesh_) { }
854
859
864
866{
867 return GeneratePVTUPath();
868}
869
871{
872 return GetCollectionName() + ".pvd";
873}
874
876 const std::string &prefix)
877{
878 return prefix + ".pvtu";
879}
880
882 const std::string &prefix, int rank)
883{
884 return prefix + to_padded_string(rank, pad_digits_rank) + ".vtu";
885}
886
888{
889 // add a new collection to the PDV file
890
891 std::string col_path = GenerateCollectionPath();
892 // check if the directories are created
893 {
894 std::string path = col_path + "/" + GenerateVTUPath();
895 int error_code = create_directory(path, mesh, myid);
896 if (error_code)
897 {
899 MFEM_WARNING("Error creating directory: " << path);
900 return; // do not even try to write the mesh
901 }
902 }
903 // the directory is created
904
905 // create pvd file if needed. If we are not in restart mode, a new pvd file
906 // is always created. In restart mode, we keep any previously defined
907 // timestep values as long as they are less than the currently defined time.
908
909 if (myid == 0 && !pvd_stream.is_open())
910 {
911 std::string pvdname = col_path + "/" + GeneratePVDFileName();
912
913 bool write_header = true;
914 std::ifstream pvd_in;
915 if (restart_mode && (pvd_in.open(pvdname,std::ios::binary),pvd_in.good()))
916 {
917 // PVD file exists and restart mode enabled: preserve existing time
918 // steps less than the current time.
919 std::fstream::pos_type pos_begin = pvd_in.tellg();
920 std::fstream::pos_type pos_end = pos_begin;
921
922 std::regex regexp("timestep=\"([^[:space:]]+)\".*file=\"Cycle(\\d+)");
923 std::smatch match;
924
925 std::string line;
926 while (getline(pvd_in,line))
927 {
928 if (regex_search(line,match,regexp))
929 {
930 MFEM_ASSERT(match.size() == 3, "Unable to parse DataSet");
931 double tvalue = std::stod(match[1]);
932 if (tvalue >= GetTime()) { break; }
933 int cvalue = std::stoi(match[2]);
934 MFEM_VERIFY(cvalue < GetCycle(), "Cycle " << GetCycle() <<
935 " is too small for restart mode: trying to overwrite"
936 " existing data.");
937 pos_end = pvd_in.tellg();
938 }
939 }
940 // Since pvd_in is opened in binary mode, count will store the number
941 // of bytes from the beginning of the file until the desired insertion
942 // point (in text mode on Windows this is not the case).
943 size_t count = pos_end - pos_begin;
944 if (count != 0)
945 {
946 write_header = false;
947 std::vector<char> buf(count);
948 // Read the contents of the PVD file, from the beginning to the
949 // insertion point.
950 pvd_in.clear();
951 pvd_in.seekg(pos_begin);
952 pvd_in.read(buf.data(), count);
953 pvd_in.close();
954 // Open the PVD file in truncate mode to delete the previous
955 // contents. Open in binary mode to write the data buffer without
956 // converting \r\n to \r\r\n on Windows.
957 pvd_stream.open(pvdname,std::ios::out|std::ios::trunc|std::ios::binary);
958 pvd_stream.write(buf.data(), count);
959 // Close and reopen the file in text mode, appending to the end.
960 pvd_stream.close();
961 pvd_stream.open(pvdname,std::ios::in|std::ios::out|std::ios::ate);
962 }
963 }
964 if (write_header)
965 {
966 // Initialize new pvd file.
967 pvd_stream.open(pvdname,std::ios::out|std::ios::trunc);
968 pvd_stream << "<?xml version=\"1.0\"?>\n";
969 pvd_stream << "<VTKFile type=\"Collection\" version=\"2.2\"";
970 pvd_stream << " byte_order=\"" << VTKByteOrder() << "\">\n";
971 pvd_stream << "<Collection>" << std::endl;
972 }
973 }
974
975 std::string vtu_prefix = col_path + "/" + GenerateVTUPath() + "/";
976
977 // Save the local part of the mesh and grid functions fields to the local
978 // VTU file. Also save coefficient fields.
979 {
980 std::ofstream os(vtu_prefix + GenerateVTUFileName("proc", myid));
981 os.precision(precision);
983 }
984
985 // Save the local part of the quadrature function fields.
986 for (const auto &qfield : q_field_map)
987 {
988 MFEM_VERIFY(!bdr_output,
989 "QuadratureFunction output is not supported for "
990 "ParaViewDataCollection on domain boundary!");
991 const std::string &field_name = qfield.first;
992 std::ofstream os(vtu_prefix + GenerateVTUFileName(field_name, myid));
993 qfield.second->SaveVTU(os, pv_data_format, GetCompressionLevel(), field_name);
994 }
995
996 // MPI rank 0 also creates a "PVTU" file that points to all of the separately
997 // written VTU files.
998 // This file path is then appended to the PVD file.
999 if (myid == 0)
1000 {
1001 // Create the main PVTU file
1002 {
1003 std::ofstream pvtu_out(vtu_prefix + GeneratePVTUFileName("data"));
1004 WritePVTUHeader(pvtu_out);
1005
1006 // Grid function fields and coefficient fields
1007 pvtu_out << "<PPointData>\n";
1008 for (auto &field_it : field_map)
1009 {
1010 int vec_dim = field_it.second->VectorDim();
1011 pvtu_out << "<PDataArray type=\"" << GetDataTypeString()
1012 << "\" Name=\"" << field_it.first
1013 << "\" NumberOfComponents=\"" << vec_dim << "\" "
1014 << VTKComponentLabels(vec_dim) << " "
1015 << "format=\"" << GetDataFormatString() << "\" />\n";
1016 }
1017 for (auto &field_it : coeff_field_map)
1018 {
1019 int vec_dim = 1;
1020 pvtu_out << "<PDataArray type=\"" << GetDataTypeString()
1021 << "\" Name=\"" << field_it.first
1022 << "\" NumberOfComponents=\"" << vec_dim << "\" "
1023 << "format=\"" << GetDataFormatString() << "\" />\n";
1024 }
1025 for (auto &field_it : vcoeff_field_map)
1026 {
1027 int vec_dim = field_it.second->GetVDim();
1028 pvtu_out << "<PDataArray type=\"" << GetDataTypeString()
1029 << "\" Name=\"" << field_it.first
1030 << "\" NumberOfComponents=\"" << vec_dim << "\" "
1031 << "format=\"" << GetDataFormatString() << "\" />\n";
1032 }
1033 pvtu_out << "</PPointData>\n";
1034
1035 // Element attributes
1036 pvtu_out << "<PCellData>\n";
1037 pvtu_out << "\t<PDataArray type=\"Int32\" Name=\"" << "attribute"
1038 << "\" NumberOfComponents=\"1\""
1039 << " format=\"" << GetDataFormatString() << "\"/>\n";
1040 pvtu_out << "</PCellData>\n";
1041
1042 WritePVTUFooter(pvtu_out, "proc");
1043 }
1044
1045 // Add the latest PVTU to the PVD
1046 pvd_stream << "<DataSet timestep=\"" << GetTime()
1047 << "\" group=\"\" part=\"" << 0 << "\" file=\""
1048 << GeneratePVTUPath() + "/" + GeneratePVTUFileName("data")
1049 << "\" name=\"mesh\"/>\n";
1050
1051 // Create PVTU files for each quadrature field and add them to the PVD
1052 // file
1053 for (auto &q_field : q_field_map)
1054 {
1055 const std::string &q_field_name = q_field.first;
1056 std::string q_fname = GeneratePVTUPath() + "/"
1057 + GeneratePVTUFileName(q_field_name);
1058
1059 std::ofstream pvtu_out(col_path + "/" + q_fname);
1060 WritePVTUHeader(pvtu_out);
1061 int vec_dim = q_field.second->GetVDim();
1062 pvtu_out << "<PPointData>\n";
1063 pvtu_out << "<PDataArray type=\"" << GetDataTypeString()
1064 << "\" Name=\"" << q_field_name
1065 << "\" NumberOfComponents=\"" << vec_dim << "\" "
1066 << VTKComponentLabels(vec_dim) << " "
1067 << "format=\"" << GetDataFormatString() << "\" />\n";
1068 pvtu_out << "</PPointData>\n";
1069 WritePVTUFooter(pvtu_out, q_field_name);
1070
1071 pvd_stream << "<DataSet timestep=\"" << GetTime()
1072 << "\" group=\"\" part=\"" << 0 << "\" file=\""
1073 << q_fname << "\" name=\"" << q_field_name << "\"/>\n";
1074 }
1075 pvd_stream.flush();
1076 // Move the insertion point before the closing collection tag, so that
1077 // the PVD file is valid even when writing incrementally.
1078 std::fstream::pos_type pos = pvd_stream.tellp();
1079 pvd_stream << "</Collection>\n";
1080 pvd_stream << "</VTKFile>" << std::endl;
1081 pvd_stream.seekp(pos);
1082 }
1083}
1084
1086{
1087 os << "<?xml version=\"1.0\"?>\n";
1088 os << "<VTKFile type=\"PUnstructuredGrid\"";
1089 os << " version =\"2.2\" byte_order=\"" << VTKByteOrder() << "\">\n";
1090 os << "<PUnstructuredGrid GhostLevel=\"0\">\n";
1091
1092 os << "<PPoints>\n";
1093 os << "\t<PDataArray type=\"" << GetDataTypeString() << "\" ";
1094 os << " Name=\"Points\" NumberOfComponents=\"3\""
1095 << " format=\"" << GetDataFormatString() << "\"/>\n";
1096 os << "</PPoints>\n";
1097
1098 os << "<PCells>\n";
1099 os << "\t<PDataArray type=\"Int32\" ";
1100 os << " Name=\"connectivity\" NumberOfComponents=\"1\""
1101 << " format=\"" << GetDataFormatString() << "\"/>\n";
1102 os << "\t<PDataArray type=\"Int32\" ";
1103 os << " Name=\"offsets\" NumberOfComponents=\"1\""
1104 << " format=\"" << GetDataFormatString() << "\"/>\n";
1105 os << "\t<PDataArray type=\"UInt8\" ";
1106 os << " Name=\"types\" NumberOfComponents=\"1\""
1107 << " format=\"" << GetDataFormatString() << "\"/>\n";
1108 os << "</PCells>\n";
1109}
1110
1112 const std::string &vtu_prefix)
1113{
1114 for (int ii=0; ii<num_procs; ii++)
1115 {
1116 std::string vtu_filename = GenerateVTUFileName(vtu_prefix, ii);
1117 os << "<Piece Source=\"" << vtu_filename << "\"/>\n";
1118 }
1119 os << "</PUnstructuredGrid>\n";
1120 os << "</VTKFile>\n";
1121}
1122
1123void ParaViewDataCollection::SaveDataVTU(std::ostream &os, int ref)
1124{
1125 os << "<VTKFile type=\"UnstructuredGrid\"";
1126 if (GetCompressionLevel() != 0)
1127 {
1128 os << " compressor=\"vtkZLibDataCompressor\"";
1129 }
1130 os << " version=\"2.2\" byte_order=\"" << VTKByteOrder() << "\">\n";
1131 os << "<UnstructuredGrid>\n";
1133 bdr_output);
1134
1135 // dump out the grid functions as point data
1136 os << "<PointData >\n";
1137 // save the grid functions
1138 // iterate over all grid functions
1139 for (FieldMapIterator it=field_map.begin(); it!=field_map.end(); ++it)
1140 {
1141 MFEM_VERIFY(!bdr_output,
1142 "GridFunction output is not supported for "
1143 "ParaViewDataCollection on domain boundary!");
1144 SaveGFieldVTU(os,ref,it);
1145 }
1146 // save the coefficient functions
1147 // iterate over all Coefficient and VectorCoefficient functions
1148 for (const auto &kv : coeff_field_map)
1149 {
1150 SaveCoeffFieldVTU(os, ref, kv.first, *kv.second);
1151 }
1152 for (const auto &kv : vcoeff_field_map)
1153 {
1154 SaveVCoeffFieldVTU(os, ref, kv.first, *kv.second);
1155 }
1156 os << "</PointData>\n";
1157 // close the mesh
1158 os << "</Piece>\n"; // close the piece open in the PrintVTU method
1159 os << "</UnstructuredGrid>\n";
1160 os << "</VTKFile>" << std::endl;
1161}
1162
1163void ParaViewDataCollection::SaveGFieldVTU(std::ostream &os, int ref_,
1164 const FieldMapIterator &it)
1165{
1166 RefinedGeometry *RefG;
1167 Vector val;
1168 DenseMatrix vval, pmat;
1169 std::vector<char> buf;
1170 int vec_dim = it->second->VectorDim();
1171 os << "<DataArray type=\"" << GetDataTypeString()
1172 << "\" Name=\"" << it->first
1173 << "\" NumberOfComponents=\"" << vec_dim << "\" "
1174 << VTKComponentLabels(vec_dim) << " "
1175 << "format=\"" << GetDataFormatString() << "\" >" << '\n';
1176 if (vec_dim == 1)
1177 {
1178 for (int i = 0; i < mesh->GetNE(); i++)
1179 {
1181 mesh->GetElementBaseGeometry(i), ref_, 1);
1182 it->second->GetValues(i, RefG->RefPts, val, pmat);
1183 for (int j = 0; j < val.Size(); j++)
1184 {
1185 WriteBinaryOrASCII(os, buf, val(j), "\n", pv_data_format);
1186 }
1187 }
1188 }
1189 else
1190 {
1191 // vector data
1192 for (int i = 0; i < mesh->GetNE(); i++)
1193 {
1195 mesh->GetElementBaseGeometry(i), ref_, 1);
1196 it->second->GetVectorValues(i, RefG->RefPts, vval, pmat);
1197 for (int jj = 0; jj < vval.Width(); jj++)
1198 {
1199 for (int ii = 0; ii < vval.Height(); ii++)
1200 {
1201 WriteBinaryOrASCII(os, buf, vval(ii,jj), " ", pv_data_format);
1202 }
1203 if (pv_data_format == VTKFormat::ASCII) { os << '\n'; }
1204 }
1205 }
1206 }
1208 {
1210 }
1211 os << "</DataArray>" << std::endl;
1212}
1213
1214void ParaViewDataCollection::SaveCoeffFieldVTU(std::ostream &os, int ref_,
1215 const std::string &name, Coefficient &coeff)
1216{
1217 RefinedGeometry *RefG;
1218 real_t val;
1219 std::vector<char> buf;
1220 int vec_dim = 1;
1221 os << "<DataArray type=\"" << GetDataTypeString()
1222 << "\" Name=\"" << name
1223 << "\" NumberOfComponents=\"" << vec_dim << "\""
1224 << " format=\"" << GetDataFormatString() << "\" >" << '\n';
1225 {
1226 // scalar data
1227 if (!bdr_output)
1228 {
1229 for (int i = 0; i < mesh->GetNE(); i++)
1230 {
1232 mesh->GetElementBaseGeometry(i), ref_, 1);
1233
1235 const IntegrationRule *ir = &RefG->RefPts;
1236 for (int j = 0; j < ir->GetNPoints(); j++)
1237 {
1238 const IntegrationPoint &ip = ir->IntPoint(j);
1239 eltrans->SetIntPoint(&ip);
1240 val = coeff.Eval(*eltrans, ip);
1241 WriteBinaryOrASCII(os, buf, val, "\n", pv_data_format);
1242 }
1243 }
1244 }
1245 else
1246 {
1247 for (int i = 0; i < mesh->GetNBE(); i++)
1248 {
1250 mesh->GetBdrElementBaseGeometry(i), ref_, 1);
1251
1253 const IntegrationRule *ir = &RefG->RefPts;
1254 for (int j = 0; j < ir->GetNPoints(); j++)
1255 {
1256 const IntegrationPoint &ip = ir->IntPoint(j);
1257 eltrans->SetIntPoint(&ip);
1258 val = coeff.Eval(*eltrans, ip);
1259 WriteBinaryOrASCII(os, buf, val, "\n", pv_data_format);
1260 }
1261 }
1262 }
1263 }
1265 {
1267 }
1268 os << "</DataArray>" << std::endl;
1269}
1270
1271void ParaViewDataCollection::SaveVCoeffFieldVTU(std::ostream &os, int ref_,
1272 const std::string &name, VectorCoefficient &coeff)
1273{
1274 RefinedGeometry *RefG;
1275 Vector val;
1276 std::vector<char> buf;
1277 int vec_dim = coeff.GetVDim();
1278 os << "<DataArray type=\"" << GetDataTypeString()
1279 << "\" Name=\"" << name
1280 << "\" NumberOfComponents=\"" << vec_dim << "\""
1281 << " format=\"" << GetDataFormatString() << "\" >" << '\n';
1282 {
1283 // vector data
1284 if (!bdr_output)
1285 {
1286 for (int i = 0; i < mesh->GetNE(); i++)
1287 {
1289 mesh->GetElementBaseGeometry(i), ref_, 1);
1290
1292 const IntegrationRule *ir = &RefG->RefPts;
1293 for (int j = 0; j < ir->GetNPoints(); j++)
1294 {
1295 const IntegrationPoint &ip = ir->IntPoint(j);
1296 eltrans->SetIntPoint(&ip);
1297 coeff.Eval(val, *eltrans, ip);
1298 for (int jj = 0; jj < val.Size(); jj++)
1299 {
1300 WriteBinaryOrASCII(os, buf, val(jj), " ", pv_data_format);
1301 }
1302 if (pv_data_format == VTKFormat::ASCII) { os << '\n'; }
1303 }
1304 }
1305 }
1306 else
1307 {
1308 for (int i = 0; i < mesh->GetNBE(); i++)
1309 {
1311 mesh->GetBdrElementBaseGeometry(i), ref_, 1);
1312
1314 const IntegrationRule *ir = &RefG->RefPts;
1315 for (int j = 0; j < ir->GetNPoints(); j++)
1316 {
1317 const IntegrationPoint &ip = ir->IntPoint(j);
1318 eltrans->SetIntPoint(&ip);
1319 coeff.Eval(val, *eltrans, ip);
1320 for (int jj = 0; jj < val.Size(); jj++)
1321 {
1322 WriteBinaryOrASCII(os, buf, val(jj), " ", pv_data_format);
1323 }
1324 if (pv_data_format == VTKFormat::ASCII) { os << '\n'; }
1325 }
1326 }
1327 }
1328 }
1330 {
1332 }
1333 os << "</DataArray>" << std::endl;
1334}
1335
1337{
1339 {
1340 return "ascii";
1341 }
1342 else
1343 {
1344 return "binary";
1345 }
1346}
1347
1349{
1351 {
1352 return "Float64";
1353 }
1354 else
1355 {
1356 return "Float32";
1357 }
1358}
1359
1360#ifdef MFEM_USE_HDF5
1361
1363 const std::string &collection_name, Mesh *mesh)
1364 : ParaViewDataCollectionBase(collection_name, mesh)
1365{
1366 compression = true;
1367}
1368
1370{
1371 compression = compression_;
1372}
1373
1374void ParaViewHDFDataCollection::EnsureVTKHDF()
1375{
1376 if (!vtkhdf)
1377 {
1378 if (!prefix_path.empty())
1379 {
1380 const int error_code = create_directory(prefix_path, mesh, myid);
1381 MFEM_VERIFY(error_code == 0, "Error creating directory " << prefix_path);
1382 }
1383
1384 std::string fname = prefix_path + name + ".vtkhdf";
1385 bool use_mpi = false;
1386#ifdef MFEM_USE_MPI
1387 if (ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh))
1388 {
1389 use_mpi = true;
1390#ifdef MFEM_PARALLEL_HDF5
1391 vtkhdf.reset(new VTKHDF(fname, pmesh->GetComm(), {restart_mode, time}));
1392#else
1393 MFEM_ABORT("Requires HDF5 library with parallel support enabled");
1394#endif
1395 }
1396#endif
1397 if (!use_mpi)
1398 {
1399 vtkhdf.reset(new VTKHDF(fname, {restart_mode, time}));
1400 }
1401 }
1402}
1403
1404template <typename FP_T>
1405void ParaViewHDFDataCollection::TSave()
1406{
1407 EnsureVTKHDF();
1408
1409 if (compression)
1410 {
1411 vtkhdf->EnableCompression(compression_level >= 0 ? compression_level : 6);
1412 }
1413 else
1414 {
1415 vtkhdf->DisableCompression();
1416 }
1417
1418 vtkhdf->SaveMesh<FP_T>(*mesh, high_order_output, levels_of_detail);
1419 for (const auto &field : field_map)
1420 {
1421 vtkhdf->SaveGridFunction<FP_T>(*field.second, field.first);
1422 }
1423 vtkhdf->UpdateSteps(time);
1424 vtkhdf->Flush();
1425}
1426
1428{
1429 switch (pv_data_format)
1430 {
1431 case VTKFormat::BINARY32: TSave<float>(); break;
1432 case VTKFormat::BINARY: TSave<double>(); break;
1433 default: MFEM_ABORT("Unsupported VTK format.");
1434 }
1435}
1436
1438
1439#endif
1440
1441} // end namespace MFEM
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
real_t time
Physical time (for time-dependent simulations)
virtual void RegisterQField(const std::string &field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
QFieldMap::iterator QFieldMapIterator
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
const std::string & GetCollectionName() const
Get the name of the collection.
bool own_data
Should the collection delete its mesh and fields.
DataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Initialize the collection with its name and Mesh.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
int GetCycle() const
Get time cycle (for time-dependent simulations)
GFieldMap::iterator FieldMapIterator
void SaveOneField(const FieldMapIterator &it)
Save one field to disk, assuming the collection directory exists.
void DeleteAll()
Delete data owned by the DataCollection including field information.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
static const int precision_default
Default value for precision.
virtual void SaveQField(const std::string &field_name)
Save one q-field, assuming the collection directory already exists.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
bool serial
Serial or parallel run? False iff mesh is a ParMesh.
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has '/' at the end.
std::string GetFieldFileName(const std::string &field_name) const
int myid
MPI rank (in parallel)
static const int pad_digits_default
Default value for pad_digits_*.
real_t time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
std::string GetMeshShortFileName() const
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
virtual ~DataCollection()
Delete the mesh and fields if owned by the collection.
std::string GetMeshFileName() const
void DeleteData()
Delete data owned by the DataCollection keeping field information.
bool appendRankToFileName
Append rank to any output file names.
std::string name
Name of the collection, used as a directory name when saving.
virtual void Save()
Save the collection to disk.
virtual void SaveField(const std::string &field_name)
Save one field, assuming the collection directory already exists.
MPI_Comm m_comm
Associated MPI communicator.
Mesh * mesh
The (common) mesh for the collected fields.
virtual void SaveMesh()
Save the mesh, creating the collection directory.
int precision
Precision (number of digits) used for the text output of doubles.
int format
Output mesh format: see the Format enumeration.
real_t GetTime() const
Get physical time (for time-dependent simulations)
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:244
virtual const char * Name() const
Definition fe_coll.hpp:79
const NURBSExtension * GetNURBSext() const
Definition fespace.hpp:646
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3824
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:856
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:341
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition geom.cpp:1136
static int GetRefinementLevelFromElems(Geometry::Type geom, int Npts)
Get the Refinement level based on number of elements.
Definition geom.cpp:1972
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
FiniteElementSpace * FESpace()
int VectorDim() const
Shortcut for calling FiniteElementSpace::GetVectorDim() on the underlying fes.
Definition gridfunc.cpp:347
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Mesh data type.
Definition mesh.hpp:65
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:312
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Print the mesh to the given stream using the default MFEM mesh format.
Definition mesh.hpp:2567
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:360
ElementTransformation * GetBdrElementTransformation(int i)
Returns a pointer to the transformation defining the i-th boundary element.
Definition mesh.cpp:533
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1309
void PrintVTU(std::ostream &os, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
Definition mesh.cpp:12392
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1380
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1556
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition mesh.hpp:1559
int GetOrder() const
If all KnotVector orders are identical, return that number. Otherwise, return NURBSFECollection::Vari...
Definition nurbs.hpp:834
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
iterator end()
Returns an end iterator to the registered fields.
iterator begin()
Returns a begin iterator to the registered fields.
void DeleteData(bool own_data)
Clear all associations between names and fields.
void clear()
Clears the map of registered fields without reclaiming memory.
iterator find(const std::string &fname)
Returns an iterator to the field fname.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Class for parallel grid function.
Definition pgridfunc.hpp:50
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:405
int GetMyRank() const
Definition pmesh.hpp:407
void ParPrint(std::ostream &out, const std::string &comments="") const
Definition pmesh.cpp:6341
int GetNRanks() const
Definition pmesh.hpp:406
Abstract base class for ParaViewDataCollection and ParaViewHDFDataCollection.
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
ParaViewDataCollectionBase(const std::string &name, Mesh *mesh)
void SetLevelsOfDetail(int levels_of_detail_)
Set the refinement level.
void UseRestartMode(bool restart_mode_)
Enable or disable restart mode.
void SetBoundaryOutput(bool bdr_output_)
Configures collection to save only fields evaluated on boundaries of the mesh.
void SetHighOrderOutput(bool high_order_output_)
Sets whether or not to output the data as high-order elements (false by default).
void SetDataFormat(VTKFormat fmt)
Set the data format for the ParaView output files.
void SetCompressionLevel(int compression_level_)
Set the zlib compression level.
int GetCompressionLevel() const
If compression is enabled, return the compression level, else return 0.
void SaveVCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, VectorCoefficient &coeff)
void WritePVTUHeader(std::ostream &out)
void SaveGFieldVTU(std::ostream &out, int ref_, const FieldMapIterator &it)
const char * GetDataFormatString() const
void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix)
void SaveDataVTU(std::ostream &out, int ref)
ParaViewDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
void SaveCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, Coefficient &coeff)
std::string GenerateVTUFileName(const std::string &prefix, int rank)
const char * GetDataTypeString() const
std::string GeneratePVTUFileName(const std::string &prefix)
void Save() override
Save the collection.
ParaViewHDFDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
void SetCompression(bool compression_) override
Enable or disable compression.
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition qfunction.hpp:94
int GetVDim() const
Get the vector dimension.
Definition qfunction.hpp:87
const IntegrationRule & GetIntRule(int idx) const
Get the IntegrationRule associated with entity (element or face) idx.
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition qspace.hpp:110
Mesh * GetMesh() const
Returns the mesh.
Definition qspace.hpp:116
IntegrationRule RefPts
Definition geom.hpp:321
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void SaveRootFile()
Save a VisIt root file for the collection.
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf) override
Add a quadrature function to the collection and update the root file.
void LoadVisItRootFile(const std::string &root_name)
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
void Load(int cycle_=0) override
Load the collection based on its VisIt data (described in its root file)
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
void Save() override
Save the collection and a VisIt root file.
void SetMaxLevelsOfDetail(int max_levels_of_detail)
Set VisIt parameter: maximum levels of detail for the MultiresControl.
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
void SetMesh(Mesh *new_mesh) override
Set/change the mesh associated with the collection.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
std::map< std::string, VisItFieldInfo > field_info_map
void DeleteAll()
Delete all data owned by VisItDataCollection including field data information.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
Helper class for VisIt visualization data.
int main()
GeometryRefiner GlobGeometryRefiner
Definition geom.cpp:2014
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
void WriteBase64WithSizeAndClear(std::ostream &os, std::vector< char > &buf, int compression_level)
Encode in base 64 (and potentially compress) the given data, write it to the output stream (with a he...
Definition vtk.cpp:654
VTKFormat
Data array format for VTK and VTU files.
Definition vtk.hpp:100
@ ASCII
Data arrays will be written in ASCII format.
int to_int(const std::string &str)
Convert a string to an int.
Definition text.hpp:62
void WriteBinaryOrASCII(std::ostream &os, std::vector< char > &buf, const T &val, const char *suffix, VTKFormat format)
Write either ASCII data to the stream or binary data to the buffer depending on the given format.
Definition vtk.hpp:148
float real_t
Definition config.hpp:46
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
Definition vtk.cpp:602
std::string VTKComponentLabels(int vdim)
Returns a string defining the component labels for vector-valued data arrays for use in XML VTU files...
Definition vtk.cpp:662