16#ifndef dealii_solver_fire_h
17#define dealii_solver_fire_h
90template <
typename VectorType = Vector<
double>>
147 template <
typename PreconditionerType = DiagonalMatrix<VectorType>>
149 solve(
const std::function<
double(VectorType &,
const VectorType &)> &compute,
151 const PreconditionerType &inverse_mass_matrix);
158 template <
typename MatrixType,
typename PreconditionerType>
162 const VectorType & b,
163 const PreconditionerType &preconditioner);
176 const VectorType &g)
const;
190template <
typename VectorType>
192 const double initial_timestep,
193 const double maximum_timestep,
194 const double maximum_linfty_norm)
195 : initial_timestep(initial_timestep)
196 , maximum_timestep(maximum_timestep)
197 , maximum_linfty_norm(maximum_linfty_norm)
199 AssertThrow(initial_timestep > 0. && maximum_timestep > 0. &&
200 maximum_linfty_norm > 0.,
201 ExcMessage(
"Expected positive values for initial_timestep, "
202 "maximum_timestep and maximum_linfty_norm but one "
203 "or more of the these values are not positive."));
208template <
typename VectorType>
211 const AdditionalData & data)
212 :
SolverBase<VectorType>(solver_control, vector_memory)
213 , additional_data(data)
218template <
typename VectorType>
220 const AdditionalData &data)
222 , additional_data(data)
227template <
typename VectorType>
228template <
typename PreconditionerType>
231 const std::function<
double(VectorType &,
const VectorType &)> &compute,
233 const PreconditionerType &inverse_mass_matrix)
238 const double DELAYSTEP = 5;
239 const double TIMESTEP_GROW = 1.1;
240 const double TIMESTEP_SHRINK = 0.5;
241 const double ALPHA_0 = 0.1;
242 const double ALPHA_SHRINK = 0.99;
244 using real_type =
typename VectorType::real_type;
255 VectorType &velocities = *v;
259 compute(gradients, x);
261 unsigned int iter = 0;
264 conv = this->iteration_status(iter, gradients * gradients, x);
269 const auto &maximum_timestep = additional_data.maximum_timestep;
270 double timestep = additional_data.initial_timestep;
273 double alpha = ALPHA_0;
275 unsigned int previous_iter_with_positive_v_dot_g = 0;
281 x.add(timestep, velocities);
282 inverse_mass_matrix.vmult(gradients, gradients);
283 velocities.add(-timestep, gradients);
286 compute(gradients, x);
289 conv = this->iteration_status(iter, gradient_norm_squared, x);
294 const real_type v_dot_g = velocities *
gradients;
298 const real_type velocities_norm_squared = velocities * velocities;
304 const real_type beta =
305 -alpha *
std::sqrt(velocities_norm_squared / gradient_norm_squared);
308 velocities.sadd(1. - alpha, beta, gradients);
310 if (iter - previous_iter_with_positive_v_dot_g > DELAYSTEP)
313 timestep =
std::min(timestep * TIMESTEP_GROW, maximum_timestep);
314 alpha *= ALPHA_SHRINK;
320 previous_iter_with_positive_v_dot_g = iter;
321 timestep *= TIMESTEP_SHRINK;
326 real_type vmax = velocities.linfty_norm();
331 const double minimal_timestep =
332 additional_data.maximum_linfty_norm / vmax;
333 if (minimal_timestep < timestep)
334 timestep = minimal_timestep;
337 print_vectors(iter, x, velocities, gradients);
349template <
typename VectorType>
350template <
typename MatrixType,
typename PreconditionerType>
354 const VectorType & b,
355 const PreconditionerType &preconditioner)
357 std::function<double(VectorType &,
const VectorType &)> compute_func =
358 [&](VectorType &g,
const VectorType &x) ->
double {
367 return 0.5 * A.matrix_norm_square(x) - x *
b;
370 this->solve(compute_func, x, preconditioner);
375template <
typename VectorType>
380 const VectorType &)
const
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
SolverFIRE(SolverControl &solver_control, VectorMemory< VectorType > &vector_memory, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void print_vectors(const unsigned int, const VectorType &x, const VectorType &v, const VectorType &g) const
SolverFIRE(SolverControl &solver_control, const AdditionalData &data=AdditionalData())
void solve(const std::function< double(VectorType &, const VectorType &)> &compute, VectorType &x, const PreconditionerType &inverse_mass_matrix)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
const double maximum_timestep
const double initial_timestep
const double maximum_linfty_norm
AdditionalData(const double initial_timestep=0.1, const double maximum_timestep=1, const double maximum_linfty_norm=1)