16 #ifndef dealii__eigen_h 17 #define dealii__eigen_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/lac/shifted_matrix.h> 22 #include <deal.II/lac/solver.h> 23 #include <deal.II/lac/solver_control.h> 24 #include <deal.II/lac/solver_cg.h> 25 #include <deal.II/lac/solver_gmres.h> 26 #include <deal.II/lac/solver_minres.h> 27 #include <deal.II/lac/vector_memory.h> 28 #include <deal.II/lac/precondition.h> 32 DEAL_II_NAMESPACE_OPEN
53 template <
typename VectorType = Vector<
double> >
87 const AdditionalData &data=AdditionalData());
100 template <
typename MatrixType>
102 solve (
double &value,
136 template <
typename VectorType = Vector<
double> >
182 const AdditionalData &data=AdditionalData());
197 template <
typename MatrixType>
199 solve (
double &value,
214 template <
class VectorType>
219 Solver<VectorType>(cn, mem),
220 additional_data(data)
225 template <
class VectorType>
231 template <
class VectorType>
232 template <
typename MatrixType>
240 deallog.
push(
"Power method");
242 VectorType *Vy = this->memory.alloc ();
245 VectorType *Vr = this->memory.alloc ();
249 double length = x.l2_norm ();
250 double old_length = 0.;
259 y.add(additional_data.shift, x);
263 length = y.l2_norm ();
269 double thresh = length/x.size();
275 while (std::fabs(entry) < thresh);
280 value = (entry * x (i) < 0.) ? -length : length;
281 value -= additional_data.shift;
291 conv = this->iteration_status (iter, std::fabs(1./length-1./old_length), x);
294 this->memory.free(Vy);
295 this->memory.free(Vr);
301 std::fabs(1./length-1./old_length)));
308 template <
class VectorType>
313 Solver<VectorType>(cn, mem),
314 additional_data(data)
319 template <
class VectorType>
325 template <
class VectorType>
326 template <
typename MatrixType>
332 deallog.
push(
"Wielandt");
343 solver(inner_control, this->memory);
346 unsigned int goal = additional_data.start_adaption;
349 VectorType *Vy = this->memory.alloc ();
352 VectorType *Vr = this->memory.alloc ();
356 double length = x.l2_norm ();
357 double old_value = value;
362 double res = -std::numeric_limits<double>::max();
366 solver.
solve (A_s, y, x, prec);
369 length = y.l2_norm ();
375 double thresh = length/x.size();
381 while (std::fabs(entry) < thresh);
386 value = (entry * x (i) < 0.) ? -length : length;
388 value -= A_s.
shift ();
392 const double new_shift = - additional_data.relaxation * value
393 + (1.-additional_data.relaxation) * A_s.
shift();
394 A_s.
shift(new_shift);
399 x.equ (1./length, y);
401 if (additional_data.use_residual)
405 r.sadd(-1., value, x);
408 conv = this->iteration_status (iter, res, x);
412 res = std::fabs(1./value-1./old_value);
413 conv = this->iteration_status (iter, res, x);
418 this->memory.free(Vy);
419 this->memory.free(Vr);
431 DEAL_II_NAMESPACE_CLOSE
types::global_dof_index size_type
AdditionalData additional_data
void solve(double &value, const MatrixType &A, VectorType &x)
#define AssertThrow(cond, exc)
void solve(double &value, const MatrixType &A, VectorType &x)
unsigned int global_dof_index
Stop iteration, goal reached.
#define Assert(cond, exc)
AdditionalData additional_data
EigenPower(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
AdditionalData(const double shift=0.)
EigenInverse(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
unsigned int start_adaption
void push(const std::string &text)
types::global_dof_index size_type
AdditionalData(double relaxation=1., unsigned int start_adaption=6, bool use_residual=true)
void shift(const double sigma)
static ::ExceptionBase & ExcInternalError()