PRISMS-PF  v2.1
checkpoint.cc
Go to the documentation of this file.
1 #include "../../include/matrixFreePDE.h"
2 #include <boost/archive/binary_iarchive.hpp>
3 #include <boost/archive/binary_oarchive.hpp>
4 
5 #ifdef DEAL_II_WITH_ZLIB
6 # include <zlib.h>
7 #endif
8 
9 // Save a checkpoint
10 template <int dim, int degree>
12  computing_timer.enter_section("matrixFreePDE: save_checkpoint");
13  unsigned int my_id = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
14 
15  if (my_id == 0)
16  {
17  // if we have previously written a snapshot, then keep the last
18  // snapshot in case this one fails to save. Note: static variables
19  // will only be initialized once per model run.
20  static bool previous_snapshot_exists = (userInputs.resume_from_checkpoint == true);
21 
22  if (previous_snapshot_exists == true)
23  {
24  move_file ("restart.mesh","restart.mesh.old");
25  move_file ("restart.mesh.info","restart.mesh.info.old");
26  move_file ("restart.time.info","restart.time.info.old");
27  }
28  // from now on, we know that if we get into this
29  // function again that a snapshot has previously
30  // been written
31  previous_snapshot_exists = true;
32  }
33 
34  // save Triangulation and Solution vectors:
35  {
36  // Serializing all of the scalars together and all of the vectors together
37 
38  // First, get lists of scalar and vector fields
39  std::vector<unsigned int> scalar_var_indices, vector_var_indices;
40  for (unsigned int var=0; var < userInputs.number_of_variables; var++){
41  if (userInputs.var_type[var] == SCALAR){
42  scalar_var_indices.push_back(var);
43  }
44  else {
45  vector_var_indices.push_back(var);
46  }
47  }
48 
49  // Second, build one solution set list for scalars and one for vectors
50  std::vector<const vectorType *> solSet_transfer_scalars;
51  std::vector<const vectorType *> solSet_transfer_vectors;
52  for(unsigned int var = 0; var < userInputs.number_of_variables; ++var){
53  if (userInputs.var_type[var] == SCALAR){
54  solSet_transfer_scalars.push_back(solutionSet[var]);
55  }
56  else {
57  solSet_transfer_vectors.push_back(solutionSet[var]);
58  }
59  }
60 
61  // Finally, save the triangulation and the solutionTransfer objects
62  if (scalar_var_indices.size() > 0 && vector_var_indices.size() == 0){
63  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_scalars (*dofHandlersSet[scalar_var_indices[0]]);
64  system_trans_scalars.prepare_serialization(solSet_transfer_scalars);
65 
66  triangulation.save ("restart.mesh");
67  }
68  else if (scalar_var_indices.size() == 0 && vector_var_indices.size() > 0){
69  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_vectors (*dofHandlersSet[vector_var_indices[0]]);
70  system_trans_vectors.prepare_serialization(solSet_transfer_vectors);
71 
72  triangulation.save ("restart.mesh");
73  }
74  else {
75  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_scalars (*dofHandlersSet[scalar_var_indices[0]]);
76  system_trans_scalars.prepare_serialization(solSet_transfer_scalars);
77 
78  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_vectors (*dofHandlersSet[vector_var_indices[0]]);
79  system_trans_vectors.prepare_serialization(solSet_transfer_vectors);
80 
81  triangulation.save ("restart.mesh");
82  }
83 
84  }
85 
86  // Save information about the current increment and current time
87  if (my_id == 0){
88  std::ofstream time_info_file;
89  time_info_file.open("restart.time.info");
90  time_info_file << currentIncrement << " (currentIncrement)\n";
91  time_info_file << currentTime << " (currentTime)\n";
92  time_info_file.close();
93  }
94 
95  pcout << "*** Checkpoint created! ***" << std::endl << std::endl;
96  computing_timer.exit_section("matrixFreePDE: save_checkpoint");
97 
98 }
99 
100 
101 // Load from a previously created checkpoint
102 template <int dim, int degree>
104 
105  // First check existence of the two restart files for the mesh and field variables
106  verify_checkpoint_file_exists("restart.mesh");
107  verify_checkpoint_file_exists("restart.mesh.info");
108 
109  pcout << std::endl << "*** Resuming from a checkpoint! ***" << std::endl << std::endl;
110 
111  try
112  {
113  triangulation.load ("restart.mesh");
114  }
115  catch (...)
116  {
117  AssertThrow(false, ExcMessage("PRISMS-PF Error: Cannot open snapshot mesh file or read the triangulation stored there."));
118  }
119 
120 }
121 
122 // Load from a previously saved checkpoint
123 template <int dim, int degree>
125 
126  // Serializing all of the scalars together and all of the vectors together
127 
128  // First, get lists of scalar and vector fields
129  std::vector<unsigned int> scalar_var_indices, vector_var_indices;
130  for (unsigned int var=0; var < userInputs.number_of_variables; var++){
131  if (userInputs.var_type[var] == SCALAR){
132  scalar_var_indices.push_back(var);
133  }
134  else {
135  vector_var_indices.push_back(var);
136  }
137  }
138 
139  // Second, build one solution set list for scalars and one for vectors
140  std::vector<vectorType *> solSet_transfer_scalars;
141  std::vector<vectorType *> solSet_transfer_vectors;
142  for(unsigned int var = 0; var < userInputs.number_of_variables; ++var){
143  if (userInputs.var_type[var] == SCALAR){
144  solSet_transfer_scalars.push_back(solutionSet[var]);
145  }
146  else {
147  solSet_transfer_vectors.push_back(solutionSet[var]);
148  }
149  }
150 
151  // Finally, deserialize the fields to the solSet_transfer objects, which contain pointers to solutionSet
152  if (scalar_var_indices.size() > 0){
153  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_scalars (*dofHandlersSet[scalar_var_indices[0]]);
154  system_trans_scalars.deserialize(solSet_transfer_scalars);
155  }
156  if (vector_var_indices.size() > 0){
157  parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_vectors (*dofHandlersSet[vector_var_indices[0]]);
158  system_trans_vectors.deserialize(solSet_transfer_vectors);
159  }
160 
161 }
162 
163 // Load from a previously saved checkpoint
164 template <int dim, int degree>
166 
167  // Make sure that restart.time.info exists
168  verify_checkpoint_file_exists("restart.time.info");
169 
170  std::ifstream time_info_file;
171  time_info_file.open("restart.time.info");
172  std::string line;
173  std::getline(time_info_file, line);
174  line.erase(line.end()-19,line.end());
175  currentIncrement = dealii::Utilities::string_to_int(line);
176 
177  std::getline(time_info_file, line);
178  line.erase(line.end()-14,line.end());
179  currentTime = dealii::Utilities::string_to_double(line);
180  time_info_file.close();
181 
182 }
183 
184 
185 // Move/rename a checkpoint file
186 template <int dim, int degree>
187 void MatrixFreePDE<dim,degree>::move_file (const std::string &old_name, const std::string &new_name){
188 
189  int error = system (("mv " + old_name + " " + new_name).c_str());
190 
191  // If the above call failed, e.g. because there is no command-line
192  // available, try with internal functions.
193  if (error != 0)
194  {
195  std::ifstream ifile(new_name);
196  if (static_cast<bool>(ifile))
197  {
198  error = remove(new_name.c_str());
199  AssertThrow (error == 0, ExcMessage(std::string ("Unable to remove file: "
200  + new_name
201  + ", although it seems to exist. "
202  + "The error code is "
203  + dealii::Utilities::to_string(error) + ".")));
204  }
205 
206  error = rename(old_name.c_str(),new_name.c_str());
207  AssertThrow (error == 0, ExcMessage(std::string ("Unable to rename files: ")
208  +
209  old_name + " -> " + new_name
210  + ". The error code is "
211  + dealii::Utilities::to_string(error) + "."));
212  }
213 }
214 
215 template <int dim, int degree>
217  std::ifstream in (filename);
218  if (!in){
219  AssertThrow (false,
220  ExcMessage (std::string("PRISMS-PF Error: You are trying to restart a previous computation, "
221  "but the restart file <")
222  +
223  filename
224  +
225  "> does not appear to exist!"));
226  }
227 }
228 
229 #include "../../include/matrixFreePDE_template_instantiations.h"
void save_checkpoint()
Definition: checkpoint.cc:11
void verify_checkpoint_file_exists(const std::string filename)
Definition: checkpoint.cc:216
void load_checkpoint_triangulation()
Definition: checkpoint.cc:103
void load_checkpoint_fields()
Definition: checkpoint.cc:124
void load_checkpoint_time_info()
Definition: checkpoint.cc:165
void move_file(const std::string &, const std::string &)
Definition: checkpoint.cc:187