55template <
typename VectorType = Vector<
double>>
87 const AdditionalData & data = AdditionalData());
96 template <
typename MatrixType>
98 solve(
double &value,
const MatrixType &A, VectorType &x);
128template <
typename VectorType = Vector<
double>>
172 const AdditionalData & data = AdditionalData());
181 template <
typename MatrixType>
183 solve(
double &value,
const MatrixType &A, VectorType &x);
196template <
class VectorType>
201 , additional_data(data)
206template <
class VectorType>
207template <
typename MatrixType>
216 VectorType & y = *Vy;
219 VectorType & r = *Vr;
222 double length = x.l2_norm();
223 double old_length = 0.;
232 y.add(additional_data.shift, x);
236 length = y.l2_norm();
242 double thresh = length / x.size();
248 while (std::fabs(entry) < thresh);
253 value = (entry * x(i) < 0.) ? -length : length;
254 value -= additional_data.shift;
257 x.equ(1 / length, y);
264 conv = this->iteration_status(iter,
265 std::fabs(1. / length - 1. / old_length),
272 iter, std::fabs(1. / length - 1. / old_length)));
279template <
class VectorType>
284 , additional_data(data)
289template <
class VectorType>
290template <
typename MatrixType>
302 double current_shift = -value;
311 unsigned int goal = additional_data.start_adaption;
315 VectorType & y = *Vy;
318 VectorType & r = *Vr;
321 double length = x.l2_norm();
322 double old_value = value;
327 double res = std::numeric_limits<double>::lowest();
331 solver.
solve(A_s, y, x, prec);
334 length = y.l2_norm();
340 double thresh = length / x.size();
346 while (std::fabs(entry) < thresh);
351 value = (entry * x(i) < 0. ? -1. : 1.) / length - current_shift;
355 const auto & relaxation = additional_data.relaxation;
356 const double new_shift =
357 relaxation * (-value) + (1. - relaxation) * current_shift;
360 current_shift = new_shift;
366 x.equ(1. / length, y);
368 if (additional_data.use_residual)
372 r.sadd(-1., value, x);
375 conv = this->iteration_status(iter, res, x);
379 res = std::fabs(1. / value - 1. / old_value);
380 conv = this->iteration_status(iter, res, x);
EigenInverse(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
AdditionalData additional_data
void solve(double &value, const MatrixType &A, VectorType &x)
EigenPower(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void solve(double &value, const MatrixType &A, VectorType &x)
AdditionalData additional_data
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
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()
#define AssertThrow(cond, exc)
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
LinearOperator< Range, Domain, Payload > identity_operator(const LinearOperator< Range, Domain, Payload > &)
unsigned int global_dof_index
AdditionalData(double relaxation=1., unsigned int start_adaption=6, bool use_residual=true)
unsigned int start_adaption
AdditionalData(const double shift=0.)