BackupRestore.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014 IRIS AS
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef OPM_BACKUPRESTORE_HEADER_INCLUDED
21 #define OPM_BACKUPRESTORE_HEADER_INCLUDED
22 
23 #include <iostream>
24 
25 #include <opm/core/simulator/SimulatorState.hpp>
26 #include <opm/core/simulator/BlackoilState.hpp>
28 
29 namespace Opm {
30 
31  namespace {
32  template <class T>
33  void writeValue( std::ostream& stream, const T& value )
34  {
35  stream.write( (const char *) &value, sizeof( T ) );
36  }
37 
38  template <class T>
39  void readValue( std::istream& stream, T& value )
40  {
41  stream.read( (char *) &value, sizeof( T ) );
42  }
43 
44  // write any container that provides a size method and a data method
45  template <class Vector>
46  void writeContainer( std::ostream& stream, const Vector& vector)
47  {
48  typedef typename Vector :: value_type T;
49  unsigned int size = vector.size() * sizeof( T );
50  writeValue( stream, size );
51  if( size > 0 ) {
52  stream.write( (const char *) vector.data(), size );
53  }
54  }
55 
56  // write map where the key is a std::string
57  // and the value is a std::vector or std::array
58  template <class Map>
59  void writeMap( std::ostream& stream, const Map& m)
60  {
61  const unsigned int mapEntries = m.size();
62  writeValue( stream, mapEntries );
63  if( mapEntries > 0 )
64  {
65  const auto& end = m.end();
66  for(auto it = m.begin(); it != end; ++it )
67  {
68  // write key (assume that key is a container)
69  writeContainer( stream, (*it).first );
70  // write value (assume that value is a container)
71  writeContainer( stream, (*it).second );
72  }
73  }
74  }
75 
76  template <class Container>
77  void resizeContainer( Container& container, size_t size )
78  {
79  container.resize( size );
80  }
81 
82  template <class T, size_t n>
83  void resizeContainer( std::array<T, n>& /* a */, size_t size )
84  {
85  static_cast<void>(size);
86  assert( int(size) == int(n) );
87  }
88 
89  template <class Container>
90  void readData(std::istream& stream, Container& container, size_t datasize)
91  {
92  stream.read( reinterpret_cast<char*>( container.data() ), datasize );
93  }
94 
95  //We need to be careful with string, because string.data() returns something that
96  //"A program shall not alter any of the characters in this sequence."
97  template <>
98  void readData<std::string>(std::istream& stream, std::string& string, size_t datasize)
99  {
100  std::vector<char> raw_data(datasize);
101  readData(stream, raw_data, datasize);
102  string = std::string(raw_data.data(), datasize);
103  }
104 
105  template <class Container>
106  void readContainerImpl( std::istream& stream, Container& container, const bool adjustSize )
107  {
108  typedef typename Container :: value_type T;
109  size_t dataSize = 0;
110  readValue( stream, dataSize );
111  if( adjustSize && dataSize > 0 ) {
112  resizeContainer( container, dataSize/sizeof(T) );
113  }
114 
115  if( dataSize != container.size() * sizeof( T ) )
116  {
117  OPM_THROW(std::logic_error,
118  "Size of stored data and simulation data does not match"
119  << dataSize << " " << (container.size() * sizeof( T )) );
120  }
121  if( dataSize > 0 ) {
122  readData(stream, container, dataSize);
123  }
124  }
125 
126  template <class Container>
127  void readAndResizeContainer( std::istream& stream, Container& container )
128  {
129  readContainerImpl( stream, container, true );
130  }
131 
132  template <class Container>
133  void readContainer( std::istream& stream, Container& container )
134  {
135  readContainerImpl( stream, container, false );
136  }
137 
138  template <class Map>
139  void readMap( std::istream& stream, Map& m)
140  {
141  m.clear();
142  unsigned int mapEntries = 0;
143  readValue( stream, mapEntries );
144  for( unsigned int entry = 0; entry<mapEntries; ++entry )
145  {
146  std::pair< typename Map :: key_type, typename Map :: mapped_type > mapEntry;
147  // read key
148  readAndResizeContainer( stream, mapEntry.first );
149  // read values
150  readContainer( stream, mapEntry.second );
151 
152  // insert entry into map
153  m.insert( mapEntry );
154  }
155  }
156 
157  enum { SimulatorStateId = 0,
158  WellStateId = 1,
159  BlackoilStateId = 2,
160  WellStateFullyImplicitBackoilId = 3 };
161 
162  inline int objectId( const SimulatorState& /* state */) {
163  return SimulatorStateId;
164  }
165  inline int objectId( const WellState& /* state */) {
166  return WellStateId;
167  }
168  inline int objectId( const BlackoilState& /* state */) {
169  return BlackoilStateId;
170  }
171  inline int objectId( const WellStateFullyImplicitBlackoil& /* state */) {
172  return WellStateFullyImplicitBackoilId;
173  }
174 
175  template <class State>
176  void checkObjectId( std::istream& in, const State& state )
177  {
178  // read id and compare with object
179  int id = -1;
180  readValue( in, id );
181  if( id != objectId( state ) ) {
182  OPM_THROW(std::logic_error,"backup-restore object type mismatch");
183  }
184  }
185  }
186 
187  // SimulatorState
188  inline
189  std::ostream& operator << (std::ostream& out, const SimulatorState& state )
190  {
191  // write id of object to stream
192  writeValue( out, objectId( state ) );
193 
194  const int numPhases = state.numPhases();
195  writeValue( out, numPhases );
196 
197  // write variables
198  writeContainer( out, state.pressure() );
199  writeContainer( out, state.temperature() );
200  writeContainer( out, state.facepressure() );
201  writeContainer( out, state.faceflux() );
202  writeContainer( out, state.saturation() );
203 
204  return out;
205  }
206 
207  inline
208  std::istream& operator >> (std::istream& in, SimulatorState& state )
209  {
210  // check id of stored object
211  checkObjectId( in, state );
212 
213  int numPhases = 0;
214  readValue( in, numPhases );
215 
216  if( numPhases != state.numPhases() )
217  OPM_THROW(std::logic_error,"num phases wrong");
218 
219  // read variables
220  readContainer( in, state.pressure() );
221  readContainer( in, state.temperature() );
222  readContainer( in, state.facepressure() );
223  readContainer( in, state.faceflux() );
224  readContainer( in, state.saturation() );
225 
226  return in;
227  }
228 
229  // BlackoilState
230  inline
231  std::ostream& operator << (std::ostream& out, const BlackoilState& state )
232  {
233  // write id of object to stream
234  writeValue( out, objectId( state ) );
235 
236  // backup simulator state
237  const SimulatorState& simstate = static_cast< const SimulatorState& > (state);
238  out << simstate;
239 
240  // backup additional blackoil state variables
241  writeContainer( out, state.surfacevol() );
242  writeContainer( out, state.gasoilratio() );
243  writeContainer( out, state.rv() );
244 
245  return out;
246  }
247 
248  inline
249  std::istream& operator >> (std::istream& in, BlackoilState& state )
250  {
251  // check id of stored object
252  checkObjectId( in, state );
253 
254  // restore simulator state
255  SimulatorState& simstate = static_cast< SimulatorState& > (state);
256  in >> simstate;
257 
258  // restore additional blackoil state variables
259  readContainer( in, state.surfacevol() );
260  readContainer( in, state.gasoilratio() );
261  readContainer( in, state.rv() );
262 
263  return in;
264  }
265 
266  // WellState
267  inline
268  std::ostream& operator << (std::ostream& out, const WellState& state )
269  {
270  // write id of object to stream
271  writeValue( out, objectId( state ) );
272 
273  // backup well state
274  writeContainer( out, state.bhp() );
275  writeContainer( out, state.temperature() );
276  writeContainer( out, state.wellRates() );
277  writeContainer( out, state.perfRates() );
278  writeContainer( out, state.perfPress() );
279 
280  return out;
281  }
282 
283  inline
284  std::istream& operator >> (std::istream& in, WellState& state )
285  {
286  // check id of stored object
287  checkObjectId( in, state );
288 
289  // restore well state
290  readAndResizeContainer( in, state.bhp() );
291  readAndResizeContainer( in, state.temperature() );
292  readAndResizeContainer( in, state.wellRates() );
293  readAndResizeContainer( in, state.perfRates() );
294  readAndResizeContainer( in, state.perfPress() );
295 
296  return in;
297  }
298 
299  // WellStateFullyImplicitBlackoil
300  inline
301  std::ostream& operator << (std::ostream& out, const WellStateFullyImplicitBlackoil& state )
302  {
303  // write id of object to stream
304  writeValue( out, objectId( state ) );
305 
306  // backup well state
307  const WellState& wellState = static_cast< const WellState& > (state);
308  out << wellState;
309 
310  const int numWells = state.numWells();
311  writeValue( out, numWells );
312  if( numWells > 0 )
313  {
314  const int numPhases = state.numPhases();
315  writeValue( out, numPhases );
316 
317  // backup additional variables
318  writeContainer( out, state.perfPhaseRates() );
319  writeContainer( out, state.currentControls() );
320  writeMap( out, state.wellMap() );
321  }
322 
323  return out;
324  }
325 
326  inline
327  std::istream& operator >> (std::istream& in, WellStateFullyImplicitBlackoil& state )
328  {
329  // check id of stored object
330  checkObjectId( in, state );
331 
332  // restore well state
333  WellState& wellState = static_cast< WellState& > (state);
334  in >> wellState;
335 
336  int numWells = 0;
337  readValue( in, numWells );
338  if( numWells != state.numWells() )
339  OPM_THROW(std::logic_error,"wrong numWells");
340 
341  if( numWells > 0 )
342  {
343  int numPhases = 0;
344  readValue( in, numPhases );
345  if( numPhases != state.numPhases() )
346  OPM_THROW(std::logic_error,"wrong numPhases");
347 
348  // restore additional variables
349  readAndResizeContainer( in, state.perfPhaseRates() );
350  readAndResizeContainer( in, state.currentControls() );
351  readMap( in, state.wellMap() );
352  }
353 
354  return in;
355  }
356 
357 } // namespace Opm
358 
359 #endif // OPM_BACKUPRESTORE_HEADER_INCLUDED
Definition: AdditionalObjectDeleter.hpp:22
Definition: WellStateFullyImplicitBlackoil.hpp:41
std::vector< double > & perfPhaseRates()
One rate per phase and well connection.
Definition: WellStateFullyImplicitBlackoil.hpp:173
Ostream & operator<<(Ostream &os, const AutoDiff< Scalar > &fw)
Definition: AutoDiff.hpp:141
std::vector< int > & currentControls()
One current control per well.
Definition: WellStateFullyImplicitBlackoil.hpp:177
std::istream & operator>>(std::istream &in, SimulatorState &state)
Definition: BackupRestore.hpp:208