16#ifndef dealii_mg_coarse_h
17#define dealii_mg_coarse_h
37template <
class VectorType = Vector<
double>>
69 const VectorType & src)
const override;
88template <
class VectorType,
91 class PreconditionerType>
105 const MatrixType &
matrix,
106 const PreconditionerType &precondition);
114 const MatrixType &
matrix,
115 const PreconditionerType &precondition);
130 const VectorType & src)
const override;
174template <
typename number =
double,
class VectorType = Vector<number>>
192 const VectorType & src)
const override;
207template <
typename number =
double,
class VectorType = Vector<number>>
225 const VectorType & src)
const;
244template <
class VectorType>
246 : coarse_smooth(nullptr)
249template <
class VectorType>
252 : coarse_smooth(nullptr)
258template <
class VectorType>
266 typeid(*this).name());
270template <
class VectorType>
274 coarse_smooth =
nullptr;
278template <
class VectorType>
282 const VectorType & src)
const
284 coarse_smooth->apply(
level, dst, src);
289template <
class VectorType,
292 class PreconditionerType>
297 : solver(0, typeid(*this).name())
298 ,
matrix(0, typeid(*this).name())
299 , preconditioner(0, typeid(*this).name())
304template <
class VectorType,
307 class PreconditionerType>
311 PreconditionerType>::
312 MGCoarseGridIterativeSolver(
SolverType & solver,
313 const MatrixType & matrix,
314 const PreconditionerType &preconditioner)
315 : solver(&solver, typeid(*this).name())
317 , preconditioner(&preconditioner, typeid(*this).name())
322template <
class VectorType,
325 class PreconditionerType>
331 PreconditionerType>::initialize(
SolverType & solver_,
332 const MatrixType & matrix_,
333 const PreconditionerType &preconditioner_)
337 preconditioner = &preconditioner_;
342template <
class VectorType,
345 class PreconditionerType>
350 PreconditionerType>::clear()
367 class PreconditionerType,
368 typename std::enable_if<
369 std::is_same<VectorType, typename SolverType::vector_type>::value,
370 VectorType>::type * =
nullptr>
373 const MatrixType & matrix,
374 const PreconditionerType &preconditioner,
376 const VectorType & src)
378 solver.solve(matrix, dst, src, preconditioner);
385 class PreconditionerType,
386 typename std::enable_if<
387 !std::is_same<VectorType, typename SolverType::vector_type>::value,
388 VectorType>::type * =
nullptr>
391 const MatrixType & matrix,
392 const PreconditionerType &preconditioner,
394 const VectorType & src)
396 typename SolverType::vector_type src_;
397 typename SolverType::vector_type dst_;
402 solver.solve(matrix, dst_, src_, preconditioner);
411template <
class VectorType,
414 class PreconditionerType>
420 PreconditionerType>::operator()(
const unsigned int ,
422 const VectorType &src)
const
429 internal::MGCoarseGridIterativeSolver::solve(
430 *solver, *matrix, *preconditioner, dst, src);
437template <
typename number,
class VectorType>
442 householder.initialize(*A);
447template <
typename number,
class VectorType>
452 householder.initialize(A);
457template <
typename number,
class VectorType>
462 const VectorType &src)
const
464 householder.least_squares(dst, src);
471template <
typename number,
class VectorType>
476 matrix.reinit(A.n_rows(), A.n_cols());
478 matrix.compute_inverse_svd(threshold);
482template <
typename number,
class VectorType>
486 const VectorType &src)
const
492template <
typename number,
class VectorType>
498 for (
unsigned int i = 0; i < n; ++i)
MGCoarseGridApplySmoother()
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
SmartPointer< const MGSmootherBase< VectorType >, MGCoarseGridApplySmoother< VectorType > > coarse_smooth
MGCoarseGridApplySmoother(const MGSmootherBase< VectorType > &coarse_smooth)
void initialize(const FullMatrix< number > &A)
Householder< number > householder
MGCoarseGridHouseholder(const FullMatrix< number > *A=nullptr)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
SmartPointer< const MatrixType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > matrix
SmartPointer< SolverType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > solver
MGCoarseGridIterativeSolver()
SmartPointer< const PreconditionerType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > preconditioner
void initialize(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
MGCoarseGridIterativeSolver(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
MGCoarseGridSVD()=default
LAPACKFullMatrix< number > matrix
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(const FullMatrix< number > &A, const double threshold=0)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotInitialized()
@ matrix
Contents is actually a matrix.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)