16 #ifndef dealii_mg_coarse_h 17 #define dealii_mg_coarse_h 20 #include <deal.II/lac/full_matrix.h> 21 #include <deal.II/lac/householder.h> 22 #include <deal.II/lac/linear_operator.h> 23 #include <deal.II/lac/matrix_lib.h> 24 #include <deal.II/multigrid/mg_base.h> 26 DEAL_II_NAMESPACE_OPEN
37 template <
class VectorType = Vector<
double> >
66 const VectorType &src)
const;
89 template <
typename SolverType,
class VectorType = Vector<
double> >
102 template <
typename MatrixType,
typename PreconditionerType>
105 const PreconditionerType &);
115 template <
typename MatrixType,
typename PreconditionerType>
116 void initialize (SolverType &,
118 const PreconditionerType &);
131 const VectorType &src)
const;
137 template <
typename MatrixType>
138 void set_matrix (
const MatrixType &);
166 template <
class VectorType,
169 class PreconditionerType>
184 const PreconditionerType &precondition);
192 const PreconditionerType &precondition);
203 virtual void operator() (
const unsigned int level,
205 const VectorType &src)
const;
247 template <
typename number =
double,
class VectorType = Vector<number> >
263 const VectorType &src)
const;
280 template <
typename number =
double,
class VectorType = Vector<number> >
296 const VectorType &src)
const;
315 template <
class VectorType>
317 : coarse_smooth(nullptr)
321 template <
class VectorType>
323 : coarse_smooth(nullptr)
329 template <
class VectorType>
335 (&coarse_smooth_,
typeid(*this).name());
339 template <
class VectorType>
343 coarse_smooth =
nullptr;
347 template <
class VectorType>
351 const VectorType &src)
const 353 coarse_smooth->smooth(level, dst, src);
359 template <
typename SolverType,
class VectorType>
363 solver(0, typeid(*this).name()),
369 template <
typename SolverType,
class VectorType>
370 template <
typename MatrixType,
typename PreconditionerType>
374 const PreconditionerType &p)
376 solver(&s, typeid(*this).name())
382 precondition = linear_operator<VectorType>(
matrix, p);
386 template <
typename SolverType,
class VectorType>
394 template <
typename SolverType,
class VectorType>
395 template <
typename MatrixType,
typename PreconditionerType>
400 const PreconditionerType &p)
407 precondition = linear_operator<VectorType>(
matrix, p);
411 template <
typename SolverType,
class VectorType>
423 template <
typename SolverType,
class VectorType>
428 const VectorType &src)
const 431 solver->solve(matrix, dst, src, precondition);
435 template <
typename SolverType,
class VectorType>
436 template <
typename MatrixType>
451 template <
class VectorType,
454 class PreconditionerType>
458 solver(0, typeid(*this).name()),
459 matrix(0, typeid(*this).name()),
460 preconditioner(0, typeid(*this).name())
465 template <
class VectorType,
468 class PreconditionerType>
471 const MatrixType &matrix,
472 const PreconditionerType &preconditioner)
474 solver (&solver, typeid(*this).name()),
476 preconditioner (&preconditioner, typeid(*this).name())
481 template <
class VectorType,
484 class PreconditionerType>
488 const MatrixType &matrix_,
489 const PreconditionerType &preconditioner_)
493 preconditioner = &preconditioner_;
499 template <
class VectorType,
502 class PreconditionerType>
514 template <
class VectorType,
517 class PreconditionerType>
522 const VectorType &src)
const 527 solver->solve(*matrix, dst, src, *preconditioner);
534 template <
typename number,
class VectorType>
538 if (A !=
nullptr) householder.initialize(*A);
543 template <
typename number,
class VectorType>
547 householder.initialize(A);
552 template <
typename number,
class VectorType>
556 const VectorType &src)
const 558 householder.least_squares(dst, src);
565 template <
typename number,
class VectorType>
570 matrix.reinit(A.n_rows(), A.n_cols());
572 matrix.compute_inverse_svd(threshold);
576 template <
typename number,
class VectorType>
581 const VectorType &src)
const 587 template <
typename number,
class VectorType>
591 const unsigned int n = std::min(
matrix.n_rows(),
matrix.n_cols());
593 for (
unsigned int i=0; i<n; ++i)
594 deallog <<
' ' <<
matrix.singular_value(i);
595 deallog << std::endl;
601 DEAL_II_NAMESPACE_CLOSE
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const =0
void initialize(const FullMatrix< number > &A, const double threshold=0)
Contents is actually a matrix.
LinearOperator< VectorType > precondition
static ::ExceptionBase & ExcNotInitialized()
SmartPointer< const PreconditionerType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > preconditioner
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
MGCoarseGridApplySmoother()
void initialize(SolverType &, const MatrixType &, const PreconditionerType &)
LAPACKFullMatrix< number > matrix
SmartPointer< SolverType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > solver
SmartPointer< SolverType, MGCoarseGridLACIteration< SolverType, VectorType > > solver
#define Assert(cond, exc)
MGCoarseGridSVD()=default
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
MGCoarseGridIterativeSolver()
void initialize(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
LinearOperator< VectorType > matrix
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
MGCoarseGridLACIteration()
MGCoarseGridHouseholder(const FullMatrix< number > *A=nullptr)
SmartPointer< const MatrixType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > matrix
Householder< number > householder
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
SmartPointer< const MGSmootherBase< VectorType >, MGCoarseGridApplySmoother< VectorType > > coarse_smooth
void initialize(const FullMatrix< number > &A)
~MGCoarseGridLACIteration()
void set_matrix(const MatrixType &)