Reference documentation for deal.II version 8.5.1
mg_coarse.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2016 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__mg_coarse_h
17 #define dealii__mg_coarse_h
18 
19 
20 #include <deal.II/lac/full_matrix.h>
21 #include <deal.II/lac/matrix_lib.h>
22 #include <deal.II/lac/householder.h>
23 #include <deal.II/multigrid/mg_base.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
29 
36 template<class VectorType = Vector<double> >
37 class MGCoarseGridApplySmoother : public MGCoarseGridBase<VectorType>
38 {
39 public:
44 
49 
54 
58  void operator() (const unsigned int level,
59  VectorType &dst,
60  const VectorType &src) const;
61 
62 private:
67 };
68 
69 
83 template<typename SolverType, class VectorType = Vector<double> >
84 class MGCoarseGridLACIteration : public MGCoarseGridBase<VectorType>
85 {
86 public:
91 
96  template<typename MatrixType, typename PreconditionerType>
97  MGCoarseGridLACIteration (SolverType &,
98  const MatrixType &,
99  const PreconditionerType &);
100 
105 
109  template<typename MatrixType, typename PreconditionerType>
110  void initialize (SolverType &,
111  const MatrixType &,
112  const PreconditionerType &);
113 
117  void clear ();
118 
123  void operator() (const unsigned int level,
124  VectorType &dst,
125  const VectorType &src) const;
126 
131  template <typename MatrixType>
132  void set_matrix (const MatrixType &);
133 
134 private:
139 
144 
149 
150 } DEAL_II_DEPRECATED;
151 
152 
153 
160 template<class VectorType,
161  class SolverType,
162  class MatrixType,
163  class PreconditionerType>
165 {
166 public:
171 
176  MGCoarseGridIterativeSolver (SolverType &solver,
177  const MatrixType &matrix,
178  const PreconditionerType &precondition);
179 
184  void initialize (SolverType &solver,
185  const MatrixType &matrix,
186  const PreconditionerType &precondition);
187 
191  void clear ();
192 
197  virtual void operator() (const unsigned int level,
198  VectorType &dst,
199  const VectorType &src) const;
200 
201 private:
205  SmartPointer<SolverType,MGCoarseGridIterativeSolver<VectorType,
206  SolverType,
207  MatrixType,
208  PreconditionerType> > solver;
209 
213  SmartPointer<const MatrixType,MGCoarseGridIterativeSolver<VectorType,
214  SolverType,
215  MatrixType,
216  PreconditionerType> > matrix;
217 
221  SmartPointer<const PreconditionerType,MGCoarseGridIterativeSolver<VectorType,
222  SolverType,
223  MatrixType,
224  PreconditionerType> > preconditioner;
225 
226 };
227 
228 
229 
230 
241 template<typename number = double, class VectorType = Vector<number> >
242 class MGCoarseGridHouseholder : public MGCoarseGridBase<VectorType>
243 {
244 public:
249 
253  void initialize (const FullMatrix<number> &A);
254 
255  void operator() (const unsigned int level,
256  VectorType &dst,
257  const VectorType &src) const;
258 
259 private:
264 };
265 
274 template<typename number = double, class VectorType = Vector<number> >
275 class MGCoarseGridSVD : public MGCoarseGridBase<VectorType>
276 {
277 public:
281  MGCoarseGridSVD ();
282 
286  void initialize (const FullMatrix<number> &A, const double threshold = 0);
287 
288  void operator() (const unsigned int level,
289  VectorType &dst,
290  const VectorType &src) const;
291 
295  void log () const;
296 
297 private:
298 
303 };
304 
307 #ifndef DOXYGEN
308 /* ------------------ Functions for MGCoarseGridApplySmoother -----------*/
309 template<class VectorType>
311  : coarse_smooth(NULL)
312 {
313 }
314 
315 template<class VectorType>
317  : coarse_smooth(NULL)
318 {
319  initialize(coarse_smooth);
320 }
321 
322 
323 template<class VectorType>
324 void
326 {
327  coarse_smooth =
329  (&coarse_smooth_,typeid(*this).name());
330 }
331 
332 template<class VectorType>
333 void
335  VectorType &dst,
336  const VectorType &src) const
337 {
338  coarse_smooth->smooth(level, dst, src);
339 }
340 
341 /* ------------------ Functions for MGCoarseGridLACIteration ------------ */
342 
343 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
344 
345 template<typename SolverType, class VectorType>
348  :
349  solver(0, typeid(*this).name()),
350  matrix(0),
351  precondition(0)
352 {}
353 
354 
355 template<typename SolverType, class VectorType>
356 template<typename MatrixType, typename PreconditionerType>
358 ::MGCoarseGridLACIteration (SolverType &s,
359  const MatrixType &m,
360  const PreconditionerType &p)
361  :
362  solver(&s, typeid(*this).name())
363 {
364  matrix = new PointerMatrix<MatrixType, VectorType>(&m);
365  precondition = new PointerMatrix<PreconditionerType, VectorType>(&p);
366 }
367 
368 
369 template<typename SolverType, class VectorType>
372 {
373  clear();
374 }
375 
376 
377 template<typename SolverType, class VectorType>
378 template<typename MatrixType, typename PreconditionerType>
379 void
381 ::initialize (SolverType &s,
382  const MatrixType &m,
383  const PreconditionerType &p)
384 {
385  solver = &s;
386  if (matrix)
387  delete matrix;
388  matrix = new PointerMatrix<MatrixType, VectorType>(&m);
389  if (precondition)
390  delete precondition;
391  precondition = new PointerMatrix<PreconditionerType, VectorType>(&p);
392 }
393 
394 
395 template<typename SolverType, class VectorType>
396 void
398 ::clear()
399 {
400  solver = 0;
401  if (matrix)
402  delete matrix;
403  matrix = 0;
404  if (precondition)
405  delete precondition;
406  precondition = 0;
407 }
408 
409 
410 template<typename SolverType, class VectorType>
411 void
413 ::operator() (const unsigned int /* level */,
414  VectorType &dst,
415  const VectorType &src) const
416 {
417  Assert(solver!=0, ExcNotInitialized());
418  Assert(matrix!=0, ExcNotInitialized());
419  Assert(precondition!=0, ExcNotInitialized());
420  solver->solve(*matrix, dst, src, *precondition);
421 }
422 
423 
424 template<typename SolverType, class VectorType>
425 template<typename MatrixType>
426 void
428 ::set_matrix(const MatrixType &m)
429 {
430  if (matrix)
431  delete matrix;
432  matrix = new PointerMatrix<MatrixType, VectorType>(&m);
433 }
434 
435 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
436 
437 
438 /* ------------------ Functions for MGCoarseGridIterativeSolver ------------ */
439 
440 template <class VectorType,
441  class SolverType,
442  class MatrixType,
443  class PreconditionerType>
446  :
447  solver(0, typeid(*this).name()),
448  matrix(0, typeid(*this).name()),
449  preconditioner(0, typeid(*this).name())
450 {}
451 
452 
453 
454 template <class VectorType,
455  class SolverType,
456  class MatrixType,
457  class PreconditionerType>
459 ::MGCoarseGridIterativeSolver (SolverType &solver,
460  const MatrixType &matrix,
461  const PreconditionerType &preconditioner)
462  :
463  solver (&solver, typeid(*this).name()),
464  matrix (&matrix, typeid(*this).name()),
465  preconditioner (&preconditioner, typeid(*this).name())
466 {}
467 
468 
469 
470 template <class VectorType,
471  class SolverType,
472  class MatrixType,
473  class PreconditionerType>
474 void
476 ::initialize (SolverType &solver_,
477  const MatrixType &matrix_,
478  const PreconditionerType &preconditioner_)
479 {
480  solver = &solver_;
481  matrix = &matrix_;
482  preconditioner = &preconditioner_;
483 
484 }
485 
486 
487 
488 template <class VectorType,
489  class SolverType,
490  class MatrixType,
491  class PreconditionerType>
492 void
494 ::clear ()
495 {
496  solver = 0;
497  matrix = 0;
498  preconditioner = 0;
499 }
500 
501 
502 
503 template <class VectorType,
504  class SolverType,
505  class MatrixType,
506  class PreconditionerType>
507 void
509 ::operator() (const unsigned int /*level*/,
510  VectorType &dst,
511  const VectorType &src) const
512 {
513  Assert(solver!=0, ExcNotInitialized());
514  Assert(matrix!=0, ExcNotInitialized());
515  Assert(preconditioner!=0, ExcNotInitialized());
516  solver->solve(*matrix, dst, src, *preconditioner);
517 }
518 
519 
520 
521 /* ------------------ Functions for MGCoarseGridHouseholder ------------ */
522 
523 template<typename number, class VectorType>
525 (const FullMatrix<number> *A)
526 {
527  if (A != 0) householder.initialize(*A);
528 }
529 
530 
531 
532 template<typename number, class VectorType>
533 void
535 {
536  householder.initialize(A);
537 }
538 
539 
540 
541 template<typename number, class VectorType>
542 void
544  VectorType &dst,
545  const VectorType &src) const
546 {
547  householder.least_squares(dst, src);
548 }
549 
550 //---------------------------------------------------------------------------
551 
552 template<typename number, class VectorType>
553 inline
555 {}
556 
557 
558 
559 template<typename number, class VectorType>
560 void
562  double threshold)
563 {
564  matrix.reinit(A.n_rows(), A.n_cols());
565  matrix = A;
566  matrix.compute_inverse_svd(threshold);
567 }
568 
569 
570 template<typename number, class VectorType>
571 void
573  const unsigned int /*level*/,
574  VectorType &dst,
575  const VectorType &src) const
576 {
577  matrix.vmult(dst, src);
578 }
579 
580 
581 template<typename number, class VectorType>
582 void
584 {
585  const unsigned int n = std::min(matrix.n_rows(), matrix.n_cols());
586 
587  for (unsigned int i=0; i<n; ++i)
588  deallog << ' ' << matrix.singular_value(i);
589  deallog << std::endl;
590 }
591 
592 
593 #endif // DOXYGEN
594 
595 DEAL_II_NAMESPACE_CLOSE
596 
597 #endif
void initialize(const FullMatrix< number > &A, const double threshold=0)
MGCoarseGridHouseholder(const FullMatrix< number > *A=0)
static ::ExceptionBase & ExcNotInitialized()
SmartPointer< const PreconditionerType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > preconditioner
Definition: mg_coarse.h:224
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(SolverType &, const MatrixType &, const PreconditionerType &)
LAPACKFullMatrix< number > matrix
Definition: mg_coarse.h:302
PointerMatrixBase< VectorType > * matrix
Definition: mg_coarse.h:143
SmartPointer< SolverType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > solver
Definition: mg_coarse.h:208
SmartPointer< SolverType, MGCoarseGridLACIteration< SolverType, VectorType > > solver
Definition: mg_coarse.h:138
#define Assert(cond, exc)
Definition: exceptions.h:313
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
PointerMatrixBase< VectorType > * precondition
Definition: mg_coarse.h:148
SmartPointer< const MatrixType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > matrix
Definition: mg_coarse.h:216
Householder< number > householder
Definition: mg_coarse.h:263
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
SmartPointer< const MGSmootherBase< VectorType >, MGCoarseGridApplySmoother< VectorType > > coarse_smooth
Definition: mg_coarse.h:66
void initialize(const FullMatrix< number > &A)
void set_matrix(const MatrixType &)
void log() const