1 #include "../../include/matrixFreePDE.h" 2 #include <boost/archive/binary_iarchive.hpp> 3 #include <boost/archive/binary_oarchive.hpp> 5 #ifdef DEAL_II_WITH_ZLIB 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);
20 static bool previous_snapshot_exists = (userInputs.resume_from_checkpoint ==
true);
22 if (previous_snapshot_exists ==
true)
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");
31 previous_snapshot_exists =
true;
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);
45 vector_var_indices.push_back(var);
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]);
57 solSet_transfer_vectors.push_back(solutionSet[var]);
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);
66 triangulation.save (
"restart.mesh");
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);
72 triangulation.save (
"restart.mesh");
75 parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_scalars (*dofHandlersSet[scalar_var_indices[0]]);
76 system_trans_scalars.prepare_serialization(solSet_transfer_scalars);
78 parallel::distributed::SolutionTransfer<dim, vectorType> system_trans_vectors (*dofHandlersSet[vector_var_indices[0]]);
79 system_trans_vectors.prepare_serialization(solSet_transfer_vectors);
81 triangulation.save (
"restart.mesh");
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();
95 pcout <<
"*** Checkpoint created! ***" << std::endl << std::endl;
96 computing_timer.exit_section(
"matrixFreePDE: save_checkpoint");
102 template <
int dim,
int degree>
106 verify_checkpoint_file_exists(
"restart.mesh");
107 verify_checkpoint_file_exists(
"restart.mesh.info");
109 pcout << std::endl <<
"*** Resuming from a checkpoint! ***" << std::endl << std::endl;
113 triangulation.load (
"restart.mesh");
117 AssertThrow(
false, ExcMessage(
"PRISMS-PF Error: Cannot open snapshot mesh file or read the triangulation stored there."));
123 template <
int dim,
int degree>
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);
135 vector_var_indices.push_back(var);
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]);
147 solSet_transfer_vectors.push_back(solutionSet[var]);
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);
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);
164 template <
int dim,
int degree>
168 verify_checkpoint_file_exists(
"restart.time.info");
170 std::ifstream time_info_file;
171 time_info_file.open(
"restart.time.info");
173 std::getline(time_info_file, line);
174 line.erase(line.end()-19,line.end());
175 currentIncrement = dealii::Utilities::string_to_int(line);
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();
186 template <
int dim,
int degree>
189 int error = system ((
"mv " + old_name +
" " + new_name).c_str());
195 std::ifstream ifile(new_name);
196 if (static_cast<bool>(ifile))
198 error =
remove(new_name.c_str());
199 AssertThrow (error == 0, ExcMessage(std::string (
"Unable to remove file: " 201 +
", although it seems to exist. " 202 +
"The error code is " 203 + dealii::Utilities::to_string(error) +
".")));
206 error = rename(old_name.c_str(),new_name.c_str());
207 AssertThrow (error == 0, ExcMessage(std::string (
"Unable to rename files: ")
209 old_name +
" -> " + new_name
210 +
". The error code is " 211 + dealii::Utilities::to_string(error) +
"."));
215 template <
int dim,
int degree>
217 std::ifstream in (filename);
220 ExcMessage (std::string(
"PRISMS-PF Error: You are trying to restart a previous computation, " 221 "but the restart file <")
225 "> does not appear to exist!"));
229 #include "../../include/matrixFreePDE_template_instantiations.h"
void verify_checkpoint_file_exists(const std::string filename)
void load_checkpoint_triangulation()
void load_checkpoint_fields()
void load_checkpoint_time_info()
void move_file(const std::string &, const std::string &)