Reference documentation for deal.II version 8.5.1
eigen.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2015 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__eigen_h
17 #define dealii__eigen_h
18 
19 
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>
29 
30 #include <cmath>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
37 
53 template <typename VectorType = Vector<double> >
54 class EigenPower : private Solver<VectorType>
55 {
56 public:
61 
66  {
71  double shift;
75  AdditionalData (const double shift = 0.)
76  :
77  shift(shift)
78  {}
79 
80  };
81 
87  const AdditionalData &data=AdditionalData());
88 
92  virtual ~EigenPower ();
93 
100  template <typename MatrixType>
101  void
102  solve (double &value,
103  const MatrixType &A,
104  VectorType &x);
105 
106 protected:
111 };
112 
136 template <typename VectorType = Vector<double> >
137 class EigenInverse : private Solver<VectorType>
138 {
139 public:
144 
149  {
153  double relaxation;
154 
158  unsigned int start_adaption;
167  unsigned int start_adaption = 6,
168  bool use_residual = true)
169  :
173  {}
174 
175  };
176 
182  const AdditionalData &data=AdditionalData());
183 
184 
188  virtual ~EigenInverse ();
189 
197  template <typename MatrixType>
198  void
199  solve (double &value,
200  const MatrixType &A,
201  VectorType &x);
202 
203 protected:
208 };
209 
211 //---------------------------------------------------------------------------
212 
213 
214 template <class VectorType>
217  const AdditionalData &data)
218  :
219  Solver<VectorType>(cn, mem),
220  additional_data(data)
221 {}
222 
223 
224 
225 template <class VectorType>
227 {}
228 
229 
230 
231 template <class VectorType>
232 template <typename MatrixType>
233 void
235  const MatrixType &A,
236  VectorType &x)
237 {
239 
240  deallog.push("Power method");
241 
242  VectorType *Vy = this->memory.alloc ();
243  VectorType &y = *Vy;
244  y.reinit (x);
245  VectorType *Vr = this->memory.alloc ();
246  VectorType &r = *Vr;
247  r.reinit (x);
248 
249  double length = x.l2_norm ();
250  double old_length = 0.;
251  x *= 1./length;
252 
253  A.vmult (y,x);
254 
255  // Main loop
256  int iter=0;
257  for (; conv==SolverControl::iterate; iter++)
258  {
259  y.add(additional_data.shift, x);
260 
261  // Compute absolute value of eigenvalue
262  old_length = length;
263  length = y.l2_norm ();
264 
265  // do a little trick to compute the sign
266  // with not too much effect of round-off errors.
267  double entry = 0.;
268  size_type i = 0;
269  double thresh = length/x.size();
270  do
271  {
272  Assert (i<x.size(), ExcInternalError());
273  entry = y (i++);
274  }
275  while (std::fabs(entry) < thresh);
276 
277  --i;
278 
279  // Compute unshifted eigenvalue
280  value = (entry * x (i) < 0.) ? -length : length;
281  value -= additional_data.shift;
282 
283  // Update normalized eigenvector
284  x.equ (1/length, y);
285 
286  // Compute residual
287  A.vmult (y,x);
288 
289  // Check the change of the eigenvalue
290  // Brrr, this is not really a good criterion
291  conv = this->iteration_status (iter, std::fabs(1./length-1./old_length), x);
292  }
293 
294  this->memory.free(Vy);
295  this->memory.free(Vr);
296 
297  deallog.pop();
298 
299  // in case of failure: throw exception
301  std::fabs(1./length-1./old_length)));
302 
303  // otherwise exit as normal
304 }
305 
306 //---------------------------------------------------------------------------
307 
308 template <class VectorType>
311  const AdditionalData &data)
312  :
313  Solver<VectorType>(cn, mem),
314  additional_data(data)
315 {}
316 
317 
318 
319 template <class VectorType>
321 {}
322 
323 
324 
325 template <class VectorType>
326 template <typename MatrixType>
327 void
329  const MatrixType &A,
330  VectorType &x)
331 {
332  deallog.push("Wielandt");
333 
335 
336  // Prepare matrix for solver
337  ShiftedMatrix <MatrixType> A_s(A, -value);
338 
339  // Define solver
340  ReductionControl inner_control (5000, 1.e-16, 1.e-5, false, false);
343  solver(inner_control, this->memory);
344 
345  // Next step for recomputing the shift
346  unsigned int goal = additional_data.start_adaption;
347 
348  // Auxiliary vector
349  VectorType *Vy = this->memory.alloc ();
350  VectorType &y = *Vy;
351  y.reinit (x);
352  VectorType *Vr = this->memory.alloc ();
353  VectorType &r = *Vr;
354  r.reinit (x);
355 
356  double length = x.l2_norm ();
357  double old_value = value;
358 
359  x *= 1./length;
360 
361  // Main loop
362  double res = -std::numeric_limits<double>::max();
363  size_type iter=0;
364  for (; conv==SolverControl::iterate; iter++)
365  {
366  solver.solve (A_s, y, x, prec);
367 
368  // Compute absolute value of eigenvalue
369  length = y.l2_norm ();
370 
371  // do a little trick to compute the sign
372  // with not too much effect of round-off errors.
373  double entry = 0.;
374  size_type i = 0;
375  double thresh = length/x.size();
376  do
377  {
378  Assert (i<x.size(), ExcInternalError());
379  entry = y (i++);
380  }
381  while (std::fabs(entry) < thresh);
382 
383  --i;
384 
385  // Compute unshifted eigenvalue
386  value = (entry * x (i) < 0.) ? -length : length;
387  value = 1./value;
388  value -= A_s.shift ();
389 
390  if (iter==goal)
391  {
392  const double new_shift = - additional_data.relaxation * value
393  + (1.-additional_data.relaxation) * A_s.shift();
394  A_s.shift(new_shift);
395  ++goal;
396  }
397 
398  // Update normalized eigenvector
399  x.equ (1./length, y);
400  // Compute residual
401  if (additional_data.use_residual)
402  {
403  y.equ (value, x);
404  A.vmult(r,x);
405  r.sadd(-1., value, x);
406  res = r.l2_norm();
407  // Check the residual
408  conv = this->iteration_status (iter, res, x);
409  }
410  else
411  {
412  res = std::fabs(1./value-1./old_value);
413  conv = this->iteration_status (iter, res, x);
414  }
415  old_value = value;
416  }
417 
418  this->memory.free(Vy);
419  this->memory.free(Vr);
420 
421  deallog.pop();
422 
423  // in case of failure: throw
424  // exception
427  res));
428  // otherwise exit as normal
429 }
430 
431 DEAL_II_NAMESPACE_CLOSE
432 
433 #endif
void pop()
Definition: logstream.cc:306
types::global_dof_index size_type
Definition: eigen.h:143
Continue iteration.
AdditionalData additional_data
Definition: eigen.h:110
void solve(double &value, const MatrixType &A, VectorType &x)
Definition: eigen.h:234
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void solve(double &value, const MatrixType &A, VectorType &x)
Definition: eigen.h:328
unsigned int global_dof_index
Definition: types.h:88
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:313
AdditionalData additional_data
Definition: eigen.h:207
EigenPower(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Definition: eigen.h:215
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
AdditionalData(const double shift=0.)
Definition: eigen.h:75
EigenInverse(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Definition: eigen.h:309
virtual ~EigenInverse()
Definition: eigen.h:320
virtual ~EigenPower()
Definition: eigen.h:226
unsigned int start_adaption
Definition: eigen.h:158
Definition: solver.h:325
void push(const std::string &text)
Definition: logstream.cc:294
types::global_dof_index size_type
Definition: eigen.h:60
AdditionalData(double relaxation=1., unsigned int start_adaption=6, bool use_residual=true)
Definition: eigen.h:166
void shift(const double sigma)
static ::ExceptionBase & ExcInternalError()