AMDiS  0.1
The Adaptive Multi-Dimensional Simulation Toolbox
BackupRestore.hpp
1 #pragma once
2 
3 #include <cstdint>
4 #include <fstream>
5 #include <memory>
6 #include <vector>
7 
8 #include <dune/grid/common/backuprestore.hh>
9 #include <dune/grid/common/gridfactory.hh>
10 
11 namespace AMDiS
12 {
13  // Fallback Implementation for Grids not supporting direct BackupRestore
14  template <class Grid>
16  {
17  template <class GridView>
18  class Writer
19  {
20  enum { dim = GridView::dimension };
21  enum { dow = GridView::dimensionworld };
22 
23  using ct = typename GridView::ctype;
24 
25  public:
26  Writer(GridView const& gv)
27  : gridView_(gv)
28  {
29  std::int64_t num = 0;
30  indexMap_.resize(gridView_.size(dim));
31  auto const& indexSet = gridView_.indexSet();
32  for (auto const& vertex : vertices(gridView_))
33  indexMap_[indexSet.index(vertex)] = std::int64_t(num++);
34  }
35 
36  void writeVertices(std::ostream& out) const
37  {
38  for (auto const& vertex : vertices(gridView_)) {
39  auto v = vertex.geometry().center();
40  out.write((char*)&v, dow*sizeof(ct));
41  }
42  }
43 
44  void writeElements(std::ostream& out) const
45  {
46  auto const& indexSet = gridView_.indexSet();
47 
48  std::vector<std::int64_t> connectivity;
49  connectivity.reserve(8);
50  for (auto const& e : elements(gridView_)) {
51  unsigned int id = e.type().id();
52  out.write((char*)&id, sizeof(unsigned int));
53 
54  connectivity.clear();
55  for (unsigned int j = 0; j < e.subEntities(dim); ++j)
56  connectivity.emplace_back(indexMap_[indexSet.subIndex(e,j,dim)]);
57 
58  out.write((char*)connectivity.data(), connectivity.size()*sizeof(std::int64_t));
59  }
60  }
61 
62  private:
63  GridView gridView_;
64  std::vector<std::int64_t> indexMap_;
65  };
66 
67  class Reader
68  {
69  enum { dim = Grid::dimension };
70  enum { dow = Grid::dimensionworld };
71 
72  using ct = typename Grid::ctype;
73  using Factory = Dune::GridFactory<Grid>;
74  using GlobalCoordinates = Dune::FieldVector<ct,dow>;
75 
76  public:
77  Reader(Factory& factory, std::istream& in)
78  : factory_(factory)
79  {
80  in.read((char*)&numElements_, sizeof(std::int64_t));
81  in.read((char*)&numVertices_, sizeof(std::int64_t));
82  }
83 
84  void readVertices(std::istream& in) const
85  {
86  GlobalCoordinates p;
87  for (std::int64_t i = 0; i < numVertices_; ++i) {
88  in.read((char*)&p[0], dow*sizeof(ct));
89  factory_.insertVertex(p);
90  }
91  }
92 
93  void readElements(std::istream& in) const
94  {
95  std::vector<std::int64_t> connectivity(8);
96  std::vector<unsigned int> vertices;
97  vertices.reserve(8);
98 
99  for (std::int64_t i = 0; i < numElements_; ++i) {
100  unsigned int id = 0;
101  in.read((char*)&id, sizeof(unsigned int));
102 
103  Dune::GeometryType type(id,dim);
104  auto refElem = Dune::referenceElement<ct,dim>(type);
105 
106  in.read((char*)connectivity.data(), refElem.size(dim)*sizeof(std::int64_t));
107  vertices.clear();
108  std::copy_n(connectivity.begin(),refElem.size(dim),std::back_inserter(vertices));
109 
110  factory_.insertElement(type, vertices);
111  }
112  }
113 
114  private:
115  Factory& factory_;
116  std::int64_t numElements_ = 0;
117  std::int64_t numVertices_ = 0;
118  };
119 
120  public:
122  template <class GridView>
123  static void backup(GridView const& gv, std::string const& filename)
124  {
125  std::ofstream out(filename, std::ios::binary);
126  backup(gv,out);
127  }
128 
130  template <class GridView>
131  static void backup(GridView const& gv, std::ostream& out)
132  {
133  std::int32_t dim = GridView::dimension;
134  std::int32_t dow = GridView::dimensionworld;
135  std::int64_t num_elements = gv.size(0);
136  std::int64_t num_vertices = gv.size(dim);
137 
138  out.write((char*)&dim, sizeof(std::int32_t));
139  out.write((char*)&dow, sizeof(std::int32_t));
140  out.write((char*)&num_elements, sizeof(std::int64_t));
141  out.write((char*)&num_vertices, sizeof(std::int64_t));
142 
143  Writer<GridView> writer(gv);
144  writer.writeVertices(out);
145  writer.writeElements(out);
146  }
147 
149  static Grid* restore(std::string const& filename)
150  {
151  std::ifstream in(filename, std::ios::binary);
152  return restore(in);
153  }
154 
156  static Grid* restore(std::istream& in)
157  {
158  Dune::GridFactory<Grid> factory;
159 
160  std::int32_t dim = 0, dow = 0;
161 
162  in.read((char*)&dim, sizeof(std::int32_t));
163  in.read((char*)&dow, sizeof(std::int32_t));
164 
165  assert(Grid::dimension == dim);
166  assert(Grid::dimensionworld == dow);
167 
168  Reader reader(factory, in);
169  reader.readVertices(in);
170  reader.readElements(in);
171 
172  std::unique_ptr<Grid> ptr(factory.createGrid());
173  return ptr.release();
174  }
175  };
176 
177 } // end namespace AMDiS
static void backup(GridView const &gv, std::string const &filename)
Write a hierarchic grid to disk.
Definition: BackupRestore.hpp:123
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
Definition: BackupRestore.hpp:15
static Grid * restore(std::istream &in)
read a hierarchic grid from a stream
Definition: BackupRestore.hpp:156
static Grid * restore(std::string const &filename)
read a hierarchic grid from disk
Definition: BackupRestore.hpp:149
static void backup(GridView const &gv, std::ostream &out)
write a hierarchic grid into a stream
Definition: BackupRestore.hpp:131