libMesh
twostep_time_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #include "libmesh/twostep_time_solver.h"
20 #include "libmesh/diff_system.h"
21 #include "libmesh/euler_solver.h"
22 #include "libmesh/numeric_vector.h"
23 #include "libmesh/auto_ptr.h" // libmesh_make_unique
24 
25 namespace libMesh
26 {
27 
28 
29 
32 
33 {
34  // We start with a reasonable time solver: implicit Euler
35  core_time_solver = libmesh_make_unique<EulerSolver>(s);
36 }
37 
38 
39 
41 {
42 }
43 
44 
45 
47 {
49 
50  // The core_time_solver will handle any first_solve actions
51  first_solve = false;
52 
53  // We may have to repeat timesteps entirely if our error is bad
54  // enough
55  bool max_tolerance_met = false;
56 
57  // Calculating error values each time
58  Real single_norm(0.), double_norm(0.), error_norm(0.),
59  relative_error(0.);
60 
61  while (!max_tolerance_met)
62  {
63  // If we've been asked to reduce deltat if necessary, make sure
64  // the core timesolver does so
65  core_time_solver->reduce_deltat_on_diffsolver_failure =
67 
68  if (!quiet)
69  {
70  libMesh::out << "\n === Computing adaptive timestep === "
71  << std::endl;
72  }
73 
74  // Use the double-length timestep first (so the
75  // old_nonlinear_solution won't have to change)
76  core_time_solver->solve();
77 
78  // Save a copy of the double-length nonlinear solution
79  // and the old nonlinear solution
80  std::unique_ptr<NumericVector<Number>> double_solution =
81  _system.solution->clone();
82  std::unique_ptr<NumericVector<Number>> old_solution =
83  _system.get_vector("_old_nonlinear_solution").clone();
84 
85  double_norm = calculate_norm(_system, *double_solution);
86  if (!quiet)
87  {
88  libMesh::out << "Double norm = " << double_norm << std::endl;
89  }
90 
91  // Then reset the initial guess for our single-length calcs
92  *(_system.solution) = _system.get_vector("_old_nonlinear_solution");
93 
94  // Call two single-length timesteps
95  // Be sure that the core_time_solver does not change the
96  // timestep here. (This is unlikely because it just succeeded
97  // with a timestep twice as large!)
98  // FIXME: even if diffsolver failure is unlikely, we ought to
99  // do *something* if it happens
100  core_time_solver->reduce_deltat_on_diffsolver_failure = 0;
101 
102  Real old_time = _system.time;
103  Real old_deltat = _system.deltat;
104  _system.deltat *= 0.5;
105  core_time_solver->solve();
106  core_time_solver->advance_timestep();
107  core_time_solver->solve();
108 
109  single_norm = calculate_norm(_system, *_system.solution);
110  if (!quiet)
111  {
112  libMesh::out << "Single norm = " << single_norm << std::endl;
113  }
114 
115  // Reset the core_time_solver's reduce_deltat... value.
116  core_time_solver->reduce_deltat_on_diffsolver_failure =
118 
119  // But then back off just in case our advance_timestep() isn't
120  // called.
121  // FIXME: this probably doesn't work with multistep methods
122  _system.get_vector("_old_nonlinear_solution") = *old_solution;
123  _system.time = old_time;
124  _system.deltat = old_deltat;
125 
126  // Find the relative error
127  *double_solution -= *(_system.solution);
128  error_norm = calculate_norm(_system, *double_solution);
129  relative_error = error_norm / _system.deltat /
130  std::max(double_norm, single_norm);
131 
132  // If the relative error makes no sense, we're done
133  if (!double_norm && !single_norm)
134  return;
135 
136  if (!quiet)
137  {
138  libMesh::out << "Error norm = " << error_norm << std::endl;
139  libMesh::out << "Local relative error = "
140  << (error_norm /
141  std::max(double_norm, single_norm))
142  << std::endl;
143  libMesh::out << "Global relative error = "
144  << (error_norm / _system.deltat /
145  std::max(double_norm, single_norm))
146  << std::endl;
147  libMesh::out << "old delta t = " << _system.deltat << std::endl;
148  }
149 
150  // If our upper tolerance is negative, that means we want to set
151  // it based on the first successful time step
152  if (this->upper_tolerance < 0)
153  this->upper_tolerance = -this->upper_tolerance * relative_error;
154 
155  // If we haven't met our upper error tolerance, we'll have to
156  // repeat this timestep entirely
157  if (this->upper_tolerance && relative_error > this->upper_tolerance)
158  {
159  // Reset the initial guess for our next try
160  *(_system.solution) =
161  _system.get_vector("_old_nonlinear_solution");
162 
163  // Chop delta t in half
164  _system.deltat /= 2.;
165 
166  if (!quiet)
167  {
168  libMesh::out << "Failed to meet upper error tolerance"
169  << std::endl;
170  libMesh::out << "Retrying with delta t = "
171  << _system.deltat << std::endl;
172  }
173  }
174  else
175  max_tolerance_met = true;
176  }
177 
178 
179  // Otherwise, compare the relative error to the tolerance
180  // and adjust deltat
182 
183  // If our target tolerance is negative, that means we want to set
184  // it based on the first successful time step
185  if (this->target_tolerance < 0)
186  this->target_tolerance = -this->target_tolerance * relative_error;
187 
188  const Real global_shrink_or_growth_factor =
189  std::pow(this->target_tolerance / relative_error,
190  static_cast<Real>(1. / core_time_solver->error_order()));
191 
192  const Real local_shrink_or_growth_factor =
193  std::pow(this->target_tolerance /
194  (error_norm/std::max(double_norm, single_norm)),
195  static_cast<Real>(1. / (core_time_solver->error_order()+1.)));
196 
197  if (!quiet)
198  {
199  libMesh::out << "The global growth/shrink factor is: "
200  << global_shrink_or_growth_factor << std::endl;
201  libMesh::out << "The local growth/shrink factor is: "
202  << local_shrink_or_growth_factor << std::endl;
203  }
204 
205  // The local s.o.g. factor is based on the expected **local**
206  // truncation error for the timestepping method, the global
207  // s.o.g. factor is based on the method's **global** truncation
208  // error. You can shrink/grow the timestep to attempt to satisfy
209  // either a global or local time-discretization error tolerance.
210 
211  Real shrink_or_growth_factor =
212  this->global_tolerance ? global_shrink_or_growth_factor :
213  local_shrink_or_growth_factor;
214 
215  if (this->max_growth && this->max_growth < shrink_or_growth_factor)
216  {
217  if (!quiet && this->global_tolerance)
218  {
219  libMesh::out << "delta t is constrained by max_growth" << std::endl;
220  }
221  shrink_or_growth_factor = this->max_growth;
222  }
223 
224  _system.deltat *= shrink_or_growth_factor;
225 
226  // Restrict deltat to max-allowable value if necessary
227  if ((this->max_deltat != 0.0) && (_system.deltat > this->max_deltat))
228  {
229  if (!quiet)
230  {
231  libMesh::out << "delta t is constrained by maximum-allowable delta t."
232  << std::endl;
233  }
234  _system.deltat = this->max_deltat;
235  }
236 
237  // Restrict deltat to min-allowable value if necessary
238  if ((this->min_deltat != 0.0) && (_system.deltat < this->min_deltat))
239  {
240  if (!quiet)
241  {
242  libMesh::out << "delta t is constrained by minimum-allowable delta t."
243  << std::endl;
244  }
245  _system.deltat = this->min_deltat;
246  }
247 
248  if (!quiet)
249  {
250  libMesh::out << "new delta t = " << _system.deltat << std::endl;
251  }
252 }
253 
254 } // namespace libMesh
libMesh::TwostepTimeSolver::TwostepTimeSolver
TwostepTimeSolver(sys_type &s)
Constructor.
Definition: twostep_time_solver.C:30
libMesh::DifferentiableSystem::deltat
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
Definition: diff_system.h:260
libMesh::AdaptiveTimeSolver::calculate_norm
virtual Real calculate_norm(System &, NumericVector< Number > &)
A helper function to calculate error norms.
Definition: adaptive_time_solver.C:155
libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure
unsigned int reduce_deltat_on_diffsolver_failure
This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat ...
Definition: time_solver.h:220
libMesh::TimeSolver::quiet
bool quiet
Print extra debugging information if quiet == false.
Definition: time_solver.h:191
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::NumericVector::clone
virtual std::unique_ptr< NumericVector< T > > clone() const =0
libMesh::AdaptiveTimeSolver::upper_tolerance
Real upper_tolerance
This tolerance is the maximum relative error between an exact time integration and a single time step...
Definition: adaptive_time_solver.h:161
libMesh::DifferentiableSystem
This class provides a specific system class.
Definition: diff_system.h:53
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::TwostepTimeSolver::~TwostepTimeSolver
~TwostepTimeSolver()
Destructor.
Definition: twostep_time_solver.C:40
libMesh::AdaptiveTimeSolver::core_time_solver
std::unique_ptr< UnsteadySolver > core_time_solver
This object is used to take timesteps.
Definition: adaptive_time_solver.h:114
libMesh::AdaptiveTimeSolver::last_deltat
Real last_deltat
We need to store the value of the last deltat used, so that advance_timestep() will increment the sys...
Definition: adaptive_time_solver.h:207
libMesh::AdaptiveTimeSolver::global_tolerance
bool global_tolerance
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accura...
Definition: adaptive_time_solver.h:198
libMesh::AdaptiveTimeSolver::max_growth
Real max_growth
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old de...
Definition: adaptive_time_solver.h:181
std::pow
double pow(double a, int b)
Definition: libmesh_augment_std_namespace.h:58
libMesh::TimeSolver::_system
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:257
libMesh::System::time
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1561
libMesh::UnsteadySolver::first_solve
bool first_solve
A bool that will be true the first time solve() is called, and false thereafter.
Definition: unsteady_solver.h:161
libMesh::AdaptiveTimeSolver
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with ...
Definition: adaptive_time_solver.h:49
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::TwostepTimeSolver::solve
virtual void solve() override
This method solves for the solution at the next timestep.
Definition: twostep_time_solver.C:46
libMesh::AdaptiveTimeSolver::max_deltat
Real max_deltat
Do not allow the adaptive time solver to select deltat > max_deltat.
Definition: adaptive_time_solver.h:167
libMesh::AdaptiveTimeSolver::min_deltat
Real min_deltat
Do not allow the adaptive time solver to select deltat < min_deltat.
Definition: adaptive_time_solver.h:173
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::out
OStreamProxy out
libMesh::System::get_vector
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
libMesh::AdaptiveTimeSolver::target_tolerance
Real target_tolerance
This tolerance is the target relative error between an exact time integration and a single time step ...
Definition: adaptive_time_solver.h:144