MultigridLevel.h

class MultigridLevel

Parent classes:

  • dolfin::Variable

This class is utilized by the solver and contains all information that is needed for a multigrid level. (Although not all members are needed on each level.)

Note

Not all members are listed in this documentation. Please view the source for more information.

Methods:

void pre_smooth(uint reps)

Apply pre-smoother reps times to x.

void calc_defect()

Computes defect d = A x - b.

void restrict_to_coarser_level()

Restricts defect to coarser level: coarser\_level{\rightarrow}b = P^T d

void coarse_solve()

Compute x = A^{-1} b

void prolongate_to_finer_level()

Prolongate solution to finer level: finer\_level{\rightarrow}e = P x

void correct_error()

Correct error: x = x - e

void post_smooth(uint reps)

Apply post-smoother reps times to x

void init_system_lhs(bool galerkin_assembly=true)

Assemble matrix A0.

void init_system_rhs(bool galerkin_assembly=true)

Assemble matrix b0.

void init_boundary_conditions(const std::string& method="topological")

Init dbc_values, dbc_values_vec, x_dbc, dbc_correct and dbc_dofs. Initialize global_to_reduced_dof and reduced_to_global_dof once.

void init_reduced_system(bool update_lhs=true, bool update_rhs=true, bool update_P=true, bool keep_A0=false, bool keep_P0=false)

Create/update (if update_lhs = true) matrix A and B. Create/update (if update_rhs = true) vector br.

void init_rhs()

Create/update vector bp.

void init_prolongation(const std::string& prolongation_assembly_method="cellwise", ProlongationAssembler::BoostFlag boost_flags=ProlongationAssembler::BoostFlags::space_grouping, double eps=1e-13, bool keep_P0=false)

Initialize matrix P0/P once.

void init_vectors()

Initialize vectors b, x, d and e once.

void init_smoothers(const dolfin::Parameters& parameters)

Create/update smoothers.

void init_coarse_solver(const dolfin::Parameters& parameters)

Create/update coarse solver.

dolfin::GenericLinearAlgebraFactory& linear_algebra_factory() const

Returns a reference to the linear algebra factory used to create matrices/vectors.

bool reduce_required() const

Returns true if solve() requires a reduced vector.

boost::shared_ptr<dolfin::GenericVector> create_reduced_vector() const

Create a vector suitable for solve() (vector without boundary dofs)

void reduce_vector(const dolfin::GenericVector& x, dolfin::GenericVector& xr) const

Copy vector x to xr with boundary dofs removed.

void expand_vector(const dolfin::GenericVector& xr, dolfin::GenericVector& x) const

Copy vector xr to x with boundary values inserted.

void update_solution()

Update the solution given by u to all coarser levels (not implemented yet).

void update_solution(const dolfin::Function& u)

Update the solution of the problem to function u to all coarser levels (not implemented yet).

void update_solution(const dolfin::GenericVector& xr)

Update the solution of the problem to reduced vector xr to all coarser levels (not implemented yet).

void set_dbc_values(dolfin::GenericVector& x) const

Set x[dbc_dofs[i]] = dbc_values[i].

void zero_dbc_values(dolfin::GenericVector& x) const

Set x[dbc_dofs[i]] = 0.

Previous topic

FMG solve module

Next topic

MultigridPreconditionedKrylovSolver.h

This Page