52class InverseMassOperator
55 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
57 typedef InverseMassOperator<dim, n_components, Number> This;
59 typedef CellIntegrator<dim, n_components, Number> Integrator;
62 typedef dealii::MatrixFreeOperators::CellwiseInverseMassMatrix<dim, -1, n_components, Number>
63 InverseMassAsMatrixFreeOperator;
65 typedef std::pair<unsigned int, unsigned int> Range;
68 InverseMassOperator() : matrix_free(
nullptr), dof_index(0), quad_index(0)
73 initialize(dealii::MatrixFree<dim, Number>
const & matrix_free_in,
76 this->matrix_free = &matrix_free_in;
77 dof_index = inverse_mass_operator_data.dof_index;
78 quad_index = inverse_mass_operator_data.quad_index;
80 data = inverse_mass_operator_data.parameters;
82 dealii::FiniteElement<dim>
const & fe = matrix_free->get_dof_handler(dof_index).get_fe();
85 AssertThrow(fe.conforms(dealii::FiniteElementData<dim>::L2),
86 dealii::ExcMessage(
"InverseMassOperator only implemented for DG!"));
88 if(data.implementation_type == InverseMassType::MatrixfreeOperator)
93 fe.base_element(0).dofs_per_cell == dealii::Utilities::pow(fe.degree + 1, dim),
95 "The matrix-free inverse mass operator is currently only available for tensor-product elements."));
100 this->matrix_free->get_shape_info(0, quad_index).data[0].n_q_points_1d == fe.degree + 1,
102 "The matrix-free inverse mass operator is currently only available if n_q_points_1d = n_nodes_1d."));
106 else if(data.implementation_type == InverseMassType::ElementwiseKrylovSolver or
107 data.implementation_type == InverseMassType::BlockMatrices)
110 dealii::AffineConstraints<Number> constraint;
115 mass_operator_data.dof_index = dof_index;
116 mass_operator_data.quad_index = quad_index;
117 if(data.implementation_type == InverseMassType::ElementwiseKrylovSolver)
119 mass_operator_data.implement_block_diagonal_preconditioner_matrix_free =
true;
120 mass_operator_data.solver_block_diagonal = Elementwise::Solver::CG;
121 if(inverse_mass_operator_data.parameters.preconditioner == PreconditionerMass::None)
123 mass_operator_data.preconditioner_block_diagonal = Elementwise::Preconditioner::None;
125 else if(inverse_mass_operator_data.parameters.preconditioner ==
126 PreconditionerMass::PointJacobi)
128 mass_operator_data.preconditioner_block_diagonal =
129 Elementwise::Preconditioner::PointJacobi;
131 mass_operator_data.solver_data_block_diagonal =
132 inverse_mass_operator_data.parameters.solver_data;
135 mass_operator.initialize(*matrix_free, constraint, mass_operator_data);
137 block_jacobi_preconditioner =
138 std::make_shared<BlockJacobiPreconditioner<MassOperator<dim, n_components, Number>>>(
139 mass_operator,
true );
143 AssertThrow(
false, dealii::ExcMessage(
"The specified InverseMassType is not implemented."));
154 if(data.implementation_type == InverseMassType::MatrixfreeOperator)
159 else if(data.implementation_type == InverseMassType::ElementwiseKrylovSolver or
160 data.implementation_type == InverseMassType::BlockMatrices)
166 block_jacobi_preconditioner->update();
170 AssertThrow(
false, dealii::ExcMessage(
"The specified InverseMassType is not implemented."));
178 dst.zero_out_ghost_values();
180 if(data.implementation_type == InverseMassType::MatrixfreeOperator)
182 matrix_free->cell_loop(&This::cell_loop_matrix_free_operator,
this, dst, src);
186 block_jacobi_preconditioner->vmult(dst, src);
192 apply_scale(VectorType & dst,
double const scaling_factor, VectorType
const & src)
const
194 if(data.implementation_type == InverseMassType::MatrixfreeOperator)
200 dst.zero_out_ghost_values();
202 matrix_free->cell_loop(
203 &This::cell_loop_matrix_free_operator,
208 [&](
const unsigned int start_range,
const unsigned int end_range) {
209 for(
unsigned int i = start_range; i < end_range; ++i)
210 dst.local_element(i) *= scaling_factor;
217 dst *= scaling_factor;
224 cell_loop_matrix_free_operator(dealii::MatrixFree<dim, Number>
const &,
226 VectorType
const & src,
227 Range
const & cell_range)
const
229 Integrator integrator(*matrix_free, dof_index, quad_index);
230 InverseMassAsMatrixFreeOperator inverse_mass(integrator);
232 for(
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
234 integrator.reinit(cell);
235 integrator.read_dof_values(src, 0);
237 inverse_mass.apply(integrator.begin_dof_values(), integrator.begin_dof_values());
239 integrator.set_dof_values(dst, 0);
243 dealii::MatrixFree<dim, Number>
const * matrix_free;
245 unsigned int dof_index, quad_index;
247 InverseMassParameters data;
254 std::shared_ptr<BlockJacobiPreconditioner<MassOperator<dim, n_components, Number>>>
255 block_jacobi_preconditioner;
259 MassOperator<dim, n_components, Number> mass_operator;
286 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
290 initialize(dealii::MatrixFree<dim, Number>
const & matrix_free,
291 dealii::AffineConstraints<Number>
const & constraints,
296 mass_operator_data.dof_index = inverse_mass_operator_data.dof_index;
297 mass_operator_data.quad_index = inverse_mass_operator_data.quad_index;
298 mass_operator.initialize(matrix_free, constraints, mass_operator_data);
301 solver_data.max_iter = inverse_mass_operator_data.parameters.solver_data.max_iter;
302 solver_data.solver_tolerance_abs = inverse_mass_operator_data.parameters.solver_data.abs_tol;
303 solver_data.solver_tolerance_rel = inverse_mass_operator_data.parameters.solver_data.rel_tol;
305 if(inverse_mass_operator_data.parameters.preconditioner == PreconditionerMass::None)
307 solver_data.use_preconditioner =
false;
309 else if(inverse_mass_operator_data.parameters.preconditioner == PreconditionerMass::PointJacobi)
312 std::make_shared<JacobiPreconditioner<MassOperator<dim, n_components, Number>>>(
313 mass_operator,
true );
315 solver_data.use_preconditioner =
true;
319 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
323 std::make_shared<Krylov::SolverCG<MassOperator<dim, n_components, Number>,
325 VectorType>>(mass_operator, *preconditioner, solver_data);
333 apply(VectorType & dst, VectorType
const & src)
const
335 Assert(solver.get() != 0, dealii::ExcMessage(
"Mass solver has not been initialized."));
345 return solver->solve(dst, temp);
349 return solver->solve(dst, src);
356 std::shared_ptr<PreconditionerBase<Number>> preconditioner;
357 std::shared_ptr<Krylov::SolverBase<VectorType>> solver;