36#include <boost/algorithm/string.hpp>
44template <
int spacedim>
47 ,
rhs(
"/Reduced Poisson/Right hand side")
48 ,
bc(
"/Reduced Poisson/Dirichlet boundary conditions")
80 this->
prm.enter_subsection(
"Error");
82 this->
prm.leave_subsection();
86template <
int dim,
int spacedim>
101 spacedim>::construct_multigrid_hierarchy)
108template <
int dim,
int spacedim>
111 const std::string &ids_and_cad_file_names,
116 grid_in.
read(grid_file_name);
117#ifdef DEAL_II_WITH_OPENCASCADE
118 using map_type = std::map<types::manifold_id, std::string>;
120 for (
const auto &pair : Converter::to_value(ids_and_cad_file_names))
122 const auto &manifold_id = pair.first;
123 const auto &cad_file_name = pair.second;
124 const auto extension = boost::to_lower_copy(
125 cad_file_name.substr(cad_file_name.find_last_of(
'.') + 1));
127 if (extension ==
"iges" || extension ==
"igs")
129 else if (extension ==
"step" || extension ==
"stp")
133 ExcNotImplemented(
"We found an extension that we "
134 "do not recognize as a CAD file "
135 "extension. Bailing out."));
137 if ((std::get<0>(n_elements) == 0))
141 else if (spacedim == 3)
151 TopoDS::Face(shape)));
154 (void)ids_and_cad_file_names;
155 AssertThrow(
false, ExcNotImplemented(
"Generation of the grid failed."));
161template <
int dim,
int spacedim>
169 par.arguments_for_grid);
173 pcout <<
"Generating from name and argument failed." << std::endl
174 <<
"Trying to read from file name." << std::endl;
177 tria.refine_global(
par.initial_refinement);
182template <
int dim,
int spacedim>
188 quadrature = std::make_unique<QGauss<spacedim>>(
par.fe_degree + 1);
192template <
int dim,
int spacedim>
197 dh.distribute_dofs(*
fe);
198#ifdef MATRIX_FREE_PATH
199 dh.distribute_mg_dofs();
209 for (
const auto id :
par.dirichlet_ids)
214#ifdef MATRIX_FREE_PATH
220 std::shared_ptr<MatrixFree<spacedim, double>> system_mf_storage(
222 system_mf_storage->reinit(
228 const unsigned int nlevels =
tria.n_global_levels();
229 mg_matrices.resize(0, nlevels - 1);
231 const std::set<types::boundary_id> dirichlet_boundary_ids = {0};
232 mg_constrained_dofs.initialize(
dh);
233 mg_constrained_dofs.make_zero_boundary_constraints(
234 dh, dirichlet_boundary_ids);
239 dh.locally_owned_mg_dofs(
level),
242 mg_constrained_dofs.get_boundary_indices(
level))
244 level_constraints.
close();
247 additional_data_level;
253 std::shared_ptr<MatrixFree<spacedim, float>> mg_mf_storage_level =
254 std::make_shared<MatrixFree<spacedim, float>>();
255 mg_mf_storage_level->reinit(
mapping,
259 additional_data_level);
261 mg_matrices[
level].initialize(mg_mf_storage_level,
287 owned_dofs[1] = reduced_dh.locally_owned_dofs();
318#ifdef MATRIX_FREE_PATH
326 pcout <<
" Number of degrees of freedom: " <<
owned_dofs[0].size() <<
" + "
328 <<
" (locally owned: " <<
owned_dofs[0].n_elements() <<
" + "
329 <<
owned_dofs[1].n_elements() <<
")" << std::endl;
333#ifndef MATRIX_FREE_PATH
334template <
int dim,
int spacedim>
345 const unsigned int dofs_per_cell =
fe->n_dofs_per_cell();
346 const unsigned int n_q_points =
quadrature->size();
349 std::vector<double> rhs_values(n_q_points);
350 std::vector<Tensor<1, spacedim>> grad_phi_u(dofs_per_cell);
351 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
353 for (
const auto &cell :
dh.active_cell_iterators())
354 if (cell->is_locally_owned())
360 for (
unsigned int q = 0; q < n_q_points; ++q)
362 for (
unsigned int k = 0; k < dofs_per_cell; ++k)
363 grad_phi_u[k] = fe_values[scalar].gradient(k, q);
364 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
366 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
369 grad_phi_u[i] * grad_phi_u[j] * fe_values.
JxW(q);
371 cell_rhs(i) += fe_values.
shape_value(i, q) * rhs_values[q] *
375 cell->get_dof_indices(local_dof_indices);
376 constraints.distribute_local_to_global(cell_matrix,
495#ifdef MATRIX_FREE_PATH
496template <
int dim,
int spacedim>
500 stiffness_matrix.get_matrix_free()->initialize_dof_vector(solution.block(0));
501 stiffness_matrix.get_matrix_free()->initialize_dof_vector(
502 system_rhs.block(0));
504 solution.block(0) = 0;
505 constraints.distribute(solution.block(0));
506 solution.block(0).update_ghost_values();
508 FEEvaluation<spacedim, -1> phi(*stiffness_matrix.get_matrix_free());
509 for (
unsigned int cell = 0;
510 cell < stiffness_matrix.get_matrix_free()->n_cell_batches();
514 phi.read_dof_values_plain(solution.block(0));
516 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
519 phi.quadrature_point(q);
522 for (
unsigned int v = 0; v < VectorizedArray<double>::size(); ++v)
525 for (
unsigned int d = 0;
d < spacedim; ++
d)
527 f_value[v] = par.rhs.value(p);
529 phi.submit_gradient(-phi.get_gradient(q), q);
530 phi.submit_value(f_value, q);
533 phi.distribute_local_to_global(system_rhs.block(0));
541template <
int dim,
int spacedim>
546 pcout <<
"Preparing solve." << std::endl;
548#ifdef MATRIX_FREE_PATH
550 using Payload = dealii::internal::LinearOperatorImplementation::EmptyPayload;
564 smoother_data.
resize(0,
tria.n_global_levels() - 1);
569 smoother_data[
level].smoothing_range = 15.;
570 smoother_data[
level].degree = 5;
571 smoother_data[
level].eig_cg_n_iterations = 10;
575 smoother_data[0].smoothing_range = 1
e-3;
577 smoother_data[0].eig_cg_n_iterations = mg_matrices[0].m();
579 mg_matrices[
level].compute_diagonal();
580 smoother_data[
level].preconditioner =
581 mg_matrices[
level].get_matrix_diagonal_inverse();
583 mg_smoother.initialize(mg_matrices, smoother_data);
592 mg_interface_matrices;
593 mg_interface_matrices.
resize(0,
tria.n_global_levels() - 1);
595 mg_interface_matrices[
level].initialize(mg_matrices[
level]);
597 mg_interface_matrices);
600 mg_matrix, mg_coarse, mg_transfer, mg_smoother, mg_smoother);
601 mg.set_edge_matrices(mg_interface, mg_interface);
606 preconditioner(
dh,
mg, mg_transfer);
616 LA::MPI::PreconditionAMG prec_A;
618 LA::MPI::PreconditionAMG::AdditionalData
data;
620 data.symmetric_operator =
true;
622 pcout <<
"Initialize AMG...";
624 pcout <<
"done." << std::endl;
645#ifdef MATRIX_FREE_PATH
648 Bt.reinit_range_vector = [
this](
VectorType &vec,
const bool) {
651 Bt.reinit_domain_vector = [
this](
VectorType &vec,
const bool) {
664 if (
par.solver_name ==
"Schur")
667 pcout <<
" Prepare schur... ";
668 const auto S = B * invA * Bt;
669 pcout <<
"S was built." << std::endl;
678 pcout <<
" f norm: " << f.l2_norm() <<
", g norm: " << g.l2_norm()
682 lambda = invS * (B * invA * f - g);
683 pcout <<
" Solved for lambda in " <<
par.outer_control.last_step()
684 <<
" iterations." << std::endl;
687 u = invA * (f - Bt * lambda);
688 pcout <<
" u norm: " << u.l2_norm()
689 <<
", lambda norm: " << lambda.l2_norm() << std::endl;
691 pcout <<
" Solved for u in " <<
par.inner_control.last_step()
692 <<
" iterations." << std::endl;
699 else if (
par.solver_name ==
"AL")
701 pcout <<
"Prepare AL preconditioner... " << std::endl;
704 (std::is_same_v<LA::MPI::Vector, TrilinosWrappers::MPI::Vector>),
705 ExcNotImplemented());
707 LA::MPI::SparseMatrix reduced_mass_matrix;
723 pcout <<
"Reduced mass matrix size: " << reduced_mass_matrix.m()
724 <<
" x " << reduced_mass_matrix.n()
725 <<
", norm: " << reduced_mass_matrix.linfty_norm() << std::endl;
738 const double el = reduced_mass_matrix.diag_element(local_idx);
741 "Diagonal element " + std::to_string(local_idx) +
742 " of reduced mass matrix (" + std::to_string(el) +
743 ") is close to zero. Cannot compute inverse square."));
744 inverse_squares_reduced(local_idx) = 1. / (el * el);
754 auto invW = invM * invM;
756 const double gamma = 10;
757 auto Aug = A + gamma * Bt * invW * B;
760 pcout <<
"Building augmented matrix..." << std::endl;
763 inverse_squares_reduced,
766 pcout <<
"done." << std::endl;
774 {{{{Aug, Bt}}, {{B, Zero}}}});
776 LA::MPI::BlockVector solution_block;
777 LA::MPI::BlockVector system_rhs_block;
778 AA.reinit_domain_vector(solution_block,
false);
779 AA.reinit_range_vector(system_rhs_block,
false);
784 tmp = gamma * Bt * invW *
system_rhs.block(1);
785 system_rhs_block.block(0) =
system_rhs.block(0);
786 system_rhs_block.block(0).add(1., tmp);
787 system_rhs_block.block(1) =
system_rhs.block(1);
797 augmented_lagrangian_preconditioner{Aug_inv, B, Bt, invW, gamma};
799 solver_fgmres.
solve(AA,
802 augmented_lagrangian_preconditioner);
804 pcout <<
" Solved with AL preconditioner in "
805 <<
par.outer_control.last_step() <<
" iterations." << std::endl;
809 solution_block.block(1));
817 auto output_double_number = [
this](
double input,
818 const std::string &text) {
820 std::cout << text << input << std::endl;
824 pcout <<
"- - - - - - - - - - - - - - - - - - - - - - - -"
826 pcout <<
"Estimate condition number of BBt using CG" << std::endl;
831 std::bind(output_double_number,
832 std::placeholders::_1,
833 "Condition number estimate: "));
834 using PayloadType = dealii::TrilinosWrappers::internal::
835 LinearOperatorImplementation::TrilinosPayload;
846 solver_cg.
solve(BBt, u, f, prec_no);
851 <<
"***BBt solve not successfull (see condition number above)***"
866template <
int dim,
int spacedim>
878 if (
par.refinement_strategy ==
"fixed_fraction")
881 tria, error_per_cell,
par.refinement_fraction,
par.coarsening_fraction);
883 else if (
par.refinement_strategy ==
"fixed_number")
888 par.refinement_fraction,
889 par.coarsening_fraction,
892 else if (
par.refinement_strategy ==
"global")
893 for (
const auto &cell :
tria.active_cell_iterators())
894 cell->set_refine_flag();
897 tria.prepare_coarsening_and_refinement();
901 tria.execute_coarsening_and_refinement();
911template <
int dim,
int spacedim>
916 std::string solution_name =
"solution";
921 for (
unsigned int i = 0; i < subdomain.
size(); ++i)
922 subdomain(i) =
tria.locally_owned_subdomain();
925 const std::string filename =
926 par.output_name +
"_" + std::to_string(
cycle) +
".vtu";
933template <
int dim,
int spacedim>
937 static std::vector<std::pair<double, std::string>> cycles_and_solutions;
938 static std::vector<std::pair<double, std::string>> cycles_and_particles;
940 if (cycles_and_solutions.size() ==
cycle)
944 const std::string particles_filename =
945 par.output_name +
"_particles_" + std::to_string(
cycle) +
".vtu";
949 cycles_and_particles.push_back({(double)
cycle, particles_filename});
951 std::ofstream pvd_solutions(
par.output_directory +
"/" +
par.output_name +
953 std::ofstream pvd_particles(
par.output_directory +
"/" +
par.output_name +
960template <
int dim,
int spacedim>
966 <<
"> using PETSc." << std::endl;
969 <<
"> using Trilinos with "
975 par.prm.print_parameters(
par.output_directory +
"/" +
"used_parameters_" +
976 std::to_string(dim) +
977 std::to_string(spacedim) +
".prm",
982template <
int dim,
int spacedim>
992 if (
par.output_results_before_solving)
994#ifdef MATRIX_FREE_PATH
1004#ifdef MATRIX_FREE_PATH
1015 par.convergence_table.error_from_exact(
dh,
1018 if (
cycle !=
par.n_refinement_cycles - 1)
1020 if (
pcout.is_active())
1021 par.convergence_table.output_table(
pcout.get_stream());
void constrain_dof_to_zero(const size_type constrained_dof)
void write_vtu_in_parallel(const std::string &filename, const MPI_Comm comm) const
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
const std::vector< Point< spacedim > > & get_quadrature_points() const
double JxW(const unsigned int q_point) const
const double & shape_value(const unsigned int i, const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
void attach_triangulation(Triangulation< dim, spacedim > &tria)
void read(std::istream &in, Format format=Default)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&...args)
void build(const DoFHandler< dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
ParameterAcceptor(const std::string §ion_name="")
static ParameterHandler prm
void enter_subsection(const std::string &subsection)
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType * > &all_in)
void interpolate(std::vector< VectorType * > &all_out)
boost::signals2::connection connect_condition_number_slot(const std::function< void(double)> &slot, const bool every_iteration=false)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerTypes &...preconditioners)
void compress(VectorOperation::values operation)
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
virtual size_type size() const override
BlockVectorType system_rhs
std::unique_ptr< FiniteElement< spacedim > > fe
BlockVectorType locally_relevant_solution
MappingQ< spacedim > mapping
std::vector< IndexSet > relevant_dofs
std::unique_ptr< Quadrature< spacedim > > quadrature
const ReducedPoissonParameters< spacedim > & par
std::vector< IndexSet > owned_dofs
void refine_and_transfer()
parallel::distributed::Triangulation< spacedim > tria
AffineConstraints< double > constraints
DoFHandler< spacedim > dh
LA::MPI::SparseMatrix coupling_matrix
MPI_Comm mpi_communicator
ReducedPoisson(const ReducedPoissonParameters< spacedim > &par)
void print_parameters() const
LA::MPI::Vector VectorType
void output_results() const
TimerOutput computing_timer
ReducedCoupling< 1, 2, spacedim, 1 > reduced_coupling
LA::MPI::SparseMatrix stiffness_matrix
void assemble_poisson_system()
std::string output_solution() const
ParameterAcceptorProxy< Functions::ParsedFunction< spacedim > > bc
std::string refinement_strategy
ParsedConvergenceTable convergence_table
double refinement_fraction
ReducedPoissonParameters()
bool output_results_before_solving
std::list< types::boundary_id > dirichlet_ids
std::string arguments_for_grid
ParameterAcceptorProxy< ReductionControl > inner_control
std::string output_directory
unsigned int initial_refinement
ParameterAcceptorProxy< ReductionControl > outer_control
double coarsening_fraction
unsigned int n_refinement_cycles
ParameterAcceptorProxy< Functions::ParsedFunction< spacedim > > rhs
#define DEAL_II_NOT_IMPLEMENTED()
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
LinearOperator< Domain, Range, Payload > transpose_operator(const LinearOperator< Range, Domain, Payload > &op)
LinearOperator< Domain, Range, Payload > inverse_operator(const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const Preconditioner &preconditioner)
BlockLinearOperator< Range, Domain, BlockPayload > block_operator(const BlockMatrixType &matrix)
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
std::vector< index_type > data
void write_pvd_record(std::ostream &out, const std::vector< std::pair< double, std::string > > ×_and_names)
void generate_from_name_and_arguments(Triangulation< dim, spacedim > &tria, const std::string &grid_generator_function_name, const std::string &grid_generator_function_arguments)
std::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
std::string dim_string(const int dim, const int spacedim)
void create_augmented_block(const MatrixType &A, const MatrixType &Ct, const VectorType &scaling_vector, const double gamma, MatrixType &augmented_matrix)
constexpr unsigned int invalid_unsigned_int
void refine_and_coarsen_fixed_fraction(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error, const VectorTools::NormType norm_type=VectorTools::L1_norm)
void refine_and_coarsen_fixed_number(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const types::global_cell_index max_n_cells=std::numeric_limits< types::global_cell_index >::max())
::SolutionTransfer< dim, VectorType, spacedim > SolutionTransfer
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
void read_grid_and_cad_files(const std::string &grid_file_name, const std::string &ids_and_cad_file_names, Triangulation< dim, spacedim > &tria)
TasksParallelScheme tasks_parallel_scheme
UpdateFlags mapping_update_flags