Reference documentation for deal.II version 9.1.0-pre
solver_fire.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_solver_fire_h
17 #define dealii_solver_fire_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/logstream.h>
23 
24 #include <deal.II/lac/diagonal_matrix.h>
25 #include <deal.II/lac/solver.h>
26 
27 #include <functional>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
35 
90 template <typename VectorType = Vector<double>>
91 class SolverFIRE : public Solver<VectorType>
92 {
93 public:
98  {
104  explicit AdditionalData(const double initial_timestep = 0.1,
105  const double maximum_timestep = 1,
106  const double maximum_linfty_norm = 1);
107 
111  const double initial_timestep;
112 
116  const double maximum_timestep;
117 
121  const double maximum_linfty_norm;
122  };
123 
127  SolverFIRE(SolverControl & solver_control,
128  VectorMemory<VectorType> &vector_memory,
129  const AdditionalData & data = AdditionalData());
130 
135  SolverFIRE(SolverControl &solver_control, const AdditionalData &data);
136 
140  virtual ~SolverFIRE();
141 
151  template <typename PreconditionerType = DiagonalMatrix<VectorType>>
152  void
153  solve(const std::function<double(VectorType &, const VectorType &)> &compute,
154  VectorType & x,
155  const PreconditionerType &inverse_mass_matrix);
156 
162  template <typename MatrixType, typename PreconditionerType>
163  void
164  solve(const MatrixType & A,
165  VectorType & x,
166  const VectorType & b,
167  const PreconditionerType &preconditioner);
168 
169 protected:
176  virtual void
177  print_vectors(const unsigned int,
178  const VectorType &x,
179  const VectorType &v,
180  const VectorType &g) const;
181 
186 };
187 
190 /*------------------------- Implementation ----------------------------*/
191 
192 #ifndef DOXYGEN
193 
194 template <typename VectorType>
196  const double initial_timestep,
197  const double maximum_timestep,
198  const double maximum_linfty_norm)
199  : initial_timestep(initial_timestep)
200  , maximum_timestep(maximum_timestep)
201  , maximum_linfty_norm(maximum_linfty_norm)
202 {
203  AssertThrow(initial_timestep > 0. && maximum_timestep > 0. &&
204  maximum_linfty_norm > 0.,
205  ExcMessage("Expected positive values for initial_timestep, "
206  "maximum_timestep and maximum_linfty_norm but one "
207  "or more of the these values are not positive."));
208 }
209 
210 
211 
212 template <typename VectorType>
214  VectorMemory<VectorType> &vector_memory,
215  const AdditionalData & data)
216  : Solver<VectorType>(solver_control, vector_memory)
217  , additional_data(data)
218 {}
219 
220 
221 
222 template <typename VectorType>
224  const AdditionalData &data)
225  : Solver<VectorType>(solver_control)
226  , additional_data(data)
227 {}
228 
229 
230 
231 template <typename VectorType>
233 {}
234 
235 
236 
237 template <typename VectorType>
238 template <typename PreconditionerType>
239 void
241  const std::function<double(VectorType &, const VectorType &)> &compute,
242  VectorType & x,
243  const PreconditionerType &inverse_mass_matrix)
244 {
245  LogStream::Prefix prefix("FIRE");
246 
247  // FIRE algorithm constants
248  const double DELAYSTEP = 5;
249  const double TIMESTEP_GROW = 1.1;
250  const double TIMESTEP_SHRINK = 0.5;
251  const double ALPHA_0 = 0.1;
252  const double ALPHA_SHRINK = 0.99;
253 
254  using real_type = typename VectorType::real_type;
255 
256  typename VectorMemory<VectorType>::Pointer v(this->memory);
257  typename VectorMemory<VectorType>::Pointer g(this->memory);
258 
259  // Set velocities to zero but not gradients
260  // as we are going to compute them soon.
261  v->reinit(x, false);
262  g->reinit(x, true);
263 
264  // Refer to v and g with some readable names.
265  VectorType &velocities = *v;
266  VectorType &gradients = *g;
267 
268  // Update gradients for the new x.
269  compute(gradients, x);
270 
271  unsigned int iter = 0;
272 
274  conv = this->iteration_status(iter, gradients * gradients, x);
275  if (conv != SolverControl::iterate)
276  return;
277 
278  // Refer to additional data members with some readable names.
279  const auto &maximum_timestep = additional_data.maximum_timestep;
280  double timestep = additional_data.initial_timestep;
281 
282  // First scaling factor.
283  double alpha = ALPHA_0;
284 
285  unsigned int previous_iter_with_positive_v_dot_g = 0;
286 
287  while (conv == SolverControl::iterate)
288  {
289  ++iter;
290  // Euler integration step.
291  x.add(timestep, velocities); // x += dt * v
292  inverse_mass_matrix.vmult(gradients, gradients); // g = M^{-1} * g
293  velocities.add(-timestep, gradients); // v -= dt * h
294 
295  // Compute gradients for the new x.
296  compute(gradients, x);
297 
298  const real_type gradient_norm_squared = gradients * gradients;
299  conv = this->iteration_status(iter, gradient_norm_squared, x);
300  if (conv != SolverControl::iterate)
301  break;
302 
303  // v_dot_g = V * G
304  const real_type v_dot_g = velocities * gradients;
305 
306  if (v_dot_g < 0.)
307  {
308  const real_type velocities_norm_squared = velocities * velocities;
309 
310  // Check if we divide by zero in DEBUG mode.
311  Assert(gradient_norm_squared > 0., ExcInternalError());
312 
313  // beta = - alpha |V|/|G|
314  const real_type beta =
315  -alpha * std::sqrt(velocities_norm_squared / gradient_norm_squared);
316 
317  // V = (1-alpha) V + beta G.
318  velocities.sadd(1. - alpha, beta, gradients);
319 
320  if (iter - previous_iter_with_positive_v_dot_g > DELAYSTEP)
321  {
322  // Increase timestep and decrease alpha.
323  timestep = std::min(timestep * TIMESTEP_GROW, maximum_timestep);
324  alpha *= ALPHA_SHRINK;
325  }
326  }
327  else
328  {
329  // Decrease timestep, reset alpha and set V = 0.
330  previous_iter_with_positive_v_dot_g = iter;
331  timestep *= TIMESTEP_SHRINK;
332  alpha = ALPHA_0;
333  velocities = 0.;
334  }
335 
336  real_type vmax = velocities.linfty_norm();
337 
338  // Change timestep if any dof would move more than maximum_linfty_norm.
339  if (vmax > 0.)
340  {
341  const double minimal_timestep =
343  if (minimal_timestep < timestep)
344  timestep = minimal_timestep;
345  }
346 
347  print_vectors(iter, x, velocities, gradients);
348 
349  } // While we need to iterate.
350 
351  // In the case of failure: throw exception.
352  if (conv != SolverControl::success)
353  AssertThrow(false,
354  SolverControl::NoConvergence(iter, gradients * gradients));
355 }
356 
357 
358 
359 template <typename VectorType>
360 template <typename MatrixType, typename PreconditionerType>
361 void
362 SolverFIRE<VectorType>::solve(const MatrixType & A,
363  VectorType & x,
364  const VectorType & b,
365  const PreconditionerType &preconditioner)
366 {
367  std::function<double(VectorType &, const VectorType &)> compute_func =
368  [&](VectorType &g, const VectorType &x) -> double {
369  // Residual of the quadratic form @f$ \frac{1}{2} xAx - xb @f$.
370  // G = b - Ax
371  A.residual(g, x, b);
372 
373  // Gradient G = Ax -b.
374  g *= -1.;
375 
376  // The quadratic form @f$\frac{1}{2} xAx - xb @f$.
377  return 0.5 * A.matrix_norm_square(x) - x * b;
378  };
379 
380  this->solve(compute_func, x, preconditioner);
381 }
382 
383 
384 
385 template <typename VectorType>
386 void
387 SolverFIRE<VectorType>::print_vectors(const unsigned int,
388  const VectorType &,
389  const VectorType &,
390  const VectorType &) const
391 {}
392 
393 
394 
395 #endif // DOXYGEN
396 
397 DEAL_II_NAMESPACE_CLOSE
398 
399 #endif
Continue iteration.
AdditionalData(const double initial_timestep=0.1, const double maximum_timestep=1, const double maximum_linfty_norm=1)
const AdditionalData additional_data
Definition: solver_fire.h:185
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
static::ExceptionBase & ExcMessage(std::string arg1)
Stop iteration, goal reached.
const double maximum_timestep
Definition: solver_fire.h:116
#define Assert(cond, exc)
Definition: exceptions.h:1227
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate), StateCombiner > iteration_status
Definition: solver.h:458
void solve(const std::function< double(VectorType &, const VectorType &)> &compute, VectorType &x, const PreconditionerType &inverse_mass_matrix)
virtual void print_vectors(const unsigned int, const VectorType &x, const VectorType &v, const VectorType &g) const
VectorMemory< VectorType > & memory
Definition: solver.h:407
const double initial_timestep
Definition: solver_fire.h:111
Definition: solver.h:328
const double maximum_linfty_norm
Definition: solver_fire.h:121
SolverFIRE(SolverControl &solver_control, VectorMemory< VectorType > &vector_memory, const AdditionalData &data=AdditionalData())
virtual ~SolverFIRE()
static::ExceptionBase & ExcInternalError()