Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
reduced_poisson.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2024 by Luca Heltai
4//
5// This file is part of the reduced_lagrange_multipliers application, based on
6// the deal.II library.
7//
8// The reduced_lagrange_multipliers application is free software; you can use
9// it, redistribute it, and/or modify it under the terms of the Apache-2.0
10// License WITH LLVM-exception as published by the Free Software Foundation;
11// either version 3.0 of the License, or (at your option) any later version. The
12// full text of the license can be found in the file LICENSE.md at the top level
13// of the reduced_lagrange_multipliers distribution.
14//
15// ---------------------------------------------------------------------
16
17/* ---------------------------------------------------------------------
18 *
19 * Copyright (C) 2000 - 2020 by the deal.II authors
20 *
21 * This file is part of the deal.II library.
22 *
23 * The deal.II library is free software; you can use it, redistribute
24 * it, and/or modify it under the terms of the GNU Lesser General
25 * Public License as published by the Free Software Foundation; either
26 * version 2.1 of the License, or (at your option) any later version.
27 * The full text of the license can be found in the file LICENSE.md at
28 * the top level directory of deal.II.
29 *
30 * ---------------------------------------------------------------------
31 *
32 */
33
34#include "reduced_poisson.h"
35
36#include <boost/algorithm/string.hpp>
37
38#include <type_traits>
39
41
42
43
44template <int spacedim>
46 : ParameterAcceptor("/Reduced Poisson/")
47 , rhs("/Reduced Poisson/Right hand side")
48 , bc("/Reduced Poisson/Dirichlet boundary conditions")
49 , inner_control("/Reduced Poisson/Solver/Inner control")
50 , outer_control("/Reduced Poisson/Solver/Outer control")
51{
52 add_parameter("FE degree", fe_degree, "", this->prm, Patterns::Integer(1));
53 add_parameter("Output directory", output_directory);
54 add_parameter("Output name", output_name);
55 add_parameter("Output results also before solving",
57 add_parameter("Solver type", solver_name);
58 add_parameter("Initial refinement", initial_refinement);
59 add_parameter("Dirichlet boundary ids", dirichlet_ids);
60 enter_subsection("Grid generation");
61 {
62 add_parameter("Grid generator", name_of_grid);
63 add_parameter("Grid generator arguments", arguments_for_grid);
64 }
66 enter_subsection("Refinement and remeshing");
67 {
68 add_parameter("Strategy",
70 "",
71 this->prm,
72 Patterns::Selection("fixed_fraction|fixed_number|global"));
73 add_parameter("Coarsening fraction", coarsening_fraction);
74 add_parameter("Refinement fraction", refinement_fraction);
75 add_parameter("Maximum number of cells", max_cells);
76 add_parameter("Number of refinement cycles", n_refinement_cycles);
77 }
79
80 this->prm.enter_subsection("Error");
81 convergence_table.add_parameters(this->prm);
82 this->prm.leave_subsection();
83}
84
85
86template <int dim, int spacedim>
89 : par(par)
90 , mpi_communicator(MPI_COMM_WORLD)
91 , pcout(std::cout, (Utilities::MPI::this_mpi_process(mpi_communicator) == 0))
93 pcout,
94 TimerOutput::summary,
95 TimerOutput::wall_times)
97 typename Triangulation<spacedim>::MeshSmoothing(
98 Triangulation<spacedim>::smoothing_on_refinement |
99 Triangulation<spacedim>::smoothing_on_coarsening),
100 parallel::distributed::Triangulation<
101 spacedim>::construct_multigrid_hierarchy)
102 , dh(tria)
103 , reduced_coupling(tria, par.reduced_coupling_parameters)
104 , mapping(1)
105{}
106
107
108template <int dim, int spacedim>
109void
110read_grid_and_cad_files(const std::string &grid_file_name,
111 const std::string &ids_and_cad_file_names,
113{
114 GridIn<dim, spacedim> grid_in;
115 grid_in.attach_triangulation(tria);
116 grid_in.read(grid_file_name);
117#ifdef DEAL_II_WITH_OPENCASCADE
118 using map_type = std::map<types::manifold_id, std::string>;
119 using Converter = Patterns::Tools::Convert<map_type>;
120 for (const auto &pair : Converter::to_value(ids_and_cad_file_names))
121 {
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));
126 TopoDS_Shape shape;
127 if (extension == "iges" || extension == "igs")
128 shape = OpenCASCADE::read_IGES(cad_file_name);
129 else if (extension == "step" || extension == "stp")
130 shape = OpenCASCADE::read_STEP(cad_file_name);
131 else
132 AssertThrow(false,
133 ExcNotImplemented("We found an extension that we "
134 "do not recognize as a CAD file "
135 "extension. Bailing out."));
136 const auto n_elements = OpenCASCADE::count_elements(shape);
137 if ((std::get<0>(n_elements) == 0))
138 tria.set_manifold(
139 manifold_id,
141 else if (spacedim == 3)
142 {
143 const auto t = reinterpret_cast<Triangulation<dim, 3> *>(&tria);
144 t->set_manifold(manifold_id,
146 shape));
147 }
148 else
149 tria.set_manifold(manifold_id,
151 TopoDS::Face(shape)));
152 }
153#else
154 (void)ids_and_cad_file_names;
155 AssertThrow(false, ExcNotImplemented("Generation of the grid failed."));
156#endif
157}
158
159
160
161template <int dim, int spacedim>
162void
164{
165 try
166 {
168 par.name_of_grid,
169 par.arguments_for_grid);
170 }
171 catch (...)
172 {
173 pcout << "Generating from name and argument failed." << std::endl
174 << "Trying to read from file name." << std::endl;
175 read_grid_and_cad_files(par.name_of_grid, par.arguments_for_grid, tria);
176 }
177 tria.refine_global(par.initial_refinement);
178}
179
180
181
182template <int dim, int spacedim>
183void
185{
186 TimerOutput::Scope t(computing_timer, "Initial setup");
187 fe = std::make_unique<FESystem<spacedim>>(FE_Q<spacedim>(par.fe_degree), 1);
188 quadrature = std::make_unique<QGauss<spacedim>>(par.fe_degree + 1);
189}
190
191
192template <int dim, int spacedim>
193void
195{
196 TimerOutput::Scope t(computing_timer, "Setup dofs");
197 dh.distribute_dofs(*fe);
198#ifdef MATRIX_FREE_PATH
199 dh.distribute_mg_dofs();
200#endif
201
202 owned_dofs.resize(2);
203 owned_dofs[0] = dh.locally_owned_dofs();
204 relevant_dofs.resize(2);
206 {
207 constraints.reinit(owned_dofs[0], relevant_dofs[0]);
209 for (const auto id : par.dirichlet_ids)
211 constraints.close();
212 }
213 {
214#ifdef MATRIX_FREE_PATH
215 typename MatrixFree<spacedim, double>::AdditionalData additional_data;
216 additional_data.tasks_parallel_scheme =
218 additional_data.mapping_update_flags =
220 std::shared_ptr<MatrixFree<spacedim, double>> system_mf_storage(
222 system_mf_storage->reinit(
223 mapping, dh, constraints, QGauss<1>(fe->degree + 1), additional_data);
224 stiffness_matrix.initialize(system_mf_storage);
225
226 // Perform setup for matrix-free multigrid
227 {
228 const unsigned int nlevels = tria.n_global_levels();
229 mg_matrices.resize(0, nlevels - 1);
230
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);
235
236 for (unsigned int level = 0; level < nlevels; ++level)
237 {
238 AffineConstraints<double> level_constraints(
239 dh.locally_owned_mg_dofs(level),
241 for (const types::global_dof_index dof_index :
242 mg_constrained_dofs.get_boundary_indices(level))
243 level_constraints.constrain_dof_to_zero(dof_index);
244 level_constraints.close();
245
247 additional_data_level;
248 additional_data_level.tasks_parallel_scheme =
250 additional_data_level.mapping_update_flags =
252 additional_data_level.mg_level = 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,
256 dh,
257 level_constraints,
258 QGauss<1>(fe->degree + 1),
259 additional_data_level);
260
261 mg_matrices[level].initialize(mg_mf_storage_level,
262 mg_constrained_dofs,
263 level);
264 }
265 }
266
267
268#else
269 stiffness_matrix.clear();
273 owned_dofs[0],
275 relevant_dofs[0]);
277 owned_dofs[0],
278 dsp,
280
281#endif
282 }
283 // Initialize the coupling object
284 reduced_coupling.initialize(mapping);
285
286 const auto &reduced_dh = reduced_coupling.get_dof_handler();
287 owned_dofs[1] = reduced_dh.locally_owned_dofs();
289
290 coupling_matrix.clear();
291
292 DynamicSparsityPattern dsp(dh.n_dofs(),
293 reduced_dh.n_dofs(),
294 relevant_dofs[0]);
295
296 reduced_coupling.assemble_coupling_sparsity(dsp, dh, constraints);
297
299
300 // DynamicSparsityPattern idsp(inclusions.n_dofs(),
301 // inclusions.n_dofs(),
302 // relevant_dofs[1]);
303 // for (const auto i : relevant_dofs[1])
304 // idsp.add(i, i);
305 // SparsityTools::distribute_sparsity_pattern(idsp,
306 // owned_dofs[1],
307 // mpi_communicator,
308 // relevant_dofs[1]);
309 // inclusion_matrix.reinit(owned_dofs[1],
310 // owned_dofs[1],
311 // idsp,
312 // mpi_communicator);
313 // }
314
315 // Commented out inclusions-dependent reinit
317
318#ifdef MATRIX_FREE_PATH
321#else
324#endif
325
326 pcout << " Number of degrees of freedom: " << owned_dofs[0].size() << " + "
327 << owned_dofs[1].size()
328 << " (locally owned: " << owned_dofs[0].n_elements() << " + "
329 << owned_dofs[1].n_elements() << ")" << std::endl;
330}
331
332
333#ifndef MATRIX_FREE_PATH
334template <int dim, int spacedim>
335void
337{
339 coupling_matrix = 0;
340 system_rhs = 0;
341 FEValues<spacedim> fe_values(*fe,
342 *quadrature,
345 const unsigned int dofs_per_cell = fe->n_dofs_per_cell();
346 const unsigned int n_q_points = quadrature->size();
347 FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
348 Vector<double> cell_rhs(dofs_per_cell);
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);
352 const FEValuesExtractors::Scalar scalar(0);
353 for (const auto &cell : dh.active_cell_iterators())
354 if (cell->is_locally_owned())
355 {
356 cell_matrix = 0;
357 cell_rhs = 0;
358 fe_values.reinit(cell);
359 par.rhs.value_list(fe_values.get_quadrature_points(), rhs_values);
360 for (unsigned int q = 0; q < n_q_points; ++q)
361 {
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)
365 {
366 for (unsigned int j = 0; j < dofs_per_cell; ++j)
367 {
368 cell_matrix(i, j) +=
369 grad_phi_u[i] * grad_phi_u[j] * fe_values.JxW(q);
370 }
371 cell_rhs(i) += fe_values.shape_value(i, q) * rhs_values[q] *
372 fe_values.JxW(q);
373 }
374 }
375 cell->get_dof_indices(local_dof_indices);
376 constraints.distribute_local_to_global(cell_matrix,
377 cell_rhs,
378 local_dof_indices,
380 system_rhs.block(0));
381 }
384}
385#endif
386
387// Commented out inclusions-dependent function
388/*
389template <int dim, int spacedim>
390void
391ReducedPoisson<dim, spacedim>::assemble_coupling()
392{
393 TimerOutput::Scope t(computing_timer, "Assemble Coupling matrix");
394 pcout << "Assemble coupling matrix. " << std::endl;
395
396 std::vector<types::global_dof_index>
397fe_dof_indices(fe->n_dofs_per_cell()); std::vector<types::global_dof_index>
398inclusion_dof_indices( inclusions.get_n_coefficients());
399
400 FullMatrix<double> local_coupling_matrix(fe->n_dofs_per_cell(),
401 inclusions.get_n_coefficients());
402
403 [[maybe_unused]] FullMatrix<double>
404local_bulk_matrix(fe->n_dofs_per_cell(), fe->n_dofs_per_cell());
405
406 FullMatrix<double> local_inclusion_matrix(inclusions.get_n_coefficients(),
407 inclusions.get_n_coefficients());
408
409 Vector<double> local_rhs(inclusions.get_n_coefficients());
410
411 auto particle = inclusions.inclusions_as_particles.begin();
412 while (particle != inclusions.inclusions_as_particles.end())
413 {
414 const auto &cell = particle->get_surrounding_cell();
415 const auto dh_cell =
416 typename DoFHandler<spacedim>::cell_iterator(*cell, &dh);
417 dh_cell->get_dof_indices(fe_dof_indices);
418 const auto pic =
419 inclusions.inclusions_as_particles.particles_in_cell(cell);
420 Assert(pic.begin() == particle, ExcInternalError());
421
422 auto p = pic.begin();
423 auto next_p = pic.begin();
424 while (p != pic.end())
425 {
426 const auto inclusion_id =
427inclusions.get_inclusion_id(p->get_id()); inclusion_dof_indices =
428inclusions.get_dof_indices(p->get_id()); local_coupling_matrix = 0;
429 local_inclusion_matrix = 0;
430 local_rhs = 0;
431
432 // Extract all points that refer to the same inclusion
433 std::vector<Point<spacedim>> ref_q_points;
434 for (; next_p != pic.end() &&
435 inclusions.get_inclusion_id(next_p->get_id()) ==
436inclusion_id;
437 ++next_p)
438 ref_q_points.push_back(next_p->get_reference_location());
439 FEValues<spacedim, spacedim> fev(*fe,
440 ref_q_points,
441 update_values | update_gradients);
442 fev.reinit(dh_cell);
443 for (unsigned int q = 0; q < ref_q_points.size(); ++q)
444 {
445 const auto id = p->get_id();
446 const auto &inclusion_fe_values =
447inclusions.get_fe_values(id); const auto &real_q =
448p->get_location();
449
450 // Coupling and inclusions matrix
451 for (unsigned int j = 0; j < inclusions.get_n_coefficients();
452++j)
453 {
454 for (unsigned int i = 0; i < fe->n_dofs_per_cell(); ++i)
455 local_coupling_matrix(i, j) +=
456 (fev.shape_value(i, q)) * inclusion_fe_values[j];
457 local_rhs(j) +=
458 inclusion_fe_values[j] *
459 inclusions.get_inclusion_data(inclusion_id, id, real_q);
460
461 local_inclusion_matrix(j, j) +=
462 (inclusion_fe_values[j] * inclusion_fe_values[j] /
463 inclusion_fe_values[0]);
464 }
465 ++p;
466 }
467 // I expect p and next_p to be the same now.
468 Assert(p == next_p, ExcInternalError());
469
470 // Add local matrices to global ones
471 constraints.distribute_local_to_global(local_coupling_matrix,
472 fe_dof_indices,
473 inclusion_constraints,
474 inclusion_dof_indices,
475 coupling_matrix);
476 inclusion_constraints.distribute_local_to_global(
477 local_rhs, inclusion_dof_indices, system_rhs.block(1));
478
479 inclusion_constraints.distribute_local_to_global(
480 local_inclusion_matrix, inclusion_dof_indices,
481inclusion_matrix);
482 }
483 particle = pic.end();
484 }
485 coupling_matrix.compress(VectorOperation::add);
486#ifndef MATRIX_FREE_PATH
487 stiffness_matrix.compress(VectorOperation::add);
488#endif
489 inclusion_matrix.compress(VectorOperation::add);
490 system_rhs.compress(VectorOperation::add);
491 pcout << "System rhs: " << system_rhs.l2_norm() << std::endl;
492}
493*/
494
495#ifdef MATRIX_FREE_PATH
496template <int dim, int spacedim>
497void
499{
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));
503 system_rhs = 0;
504 solution.block(0) = 0;
505 constraints.distribute(solution.block(0));
506 solution.block(0).update_ghost_values();
507
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();
511 ++cell)
512 {
513 phi.reinit(cell);
514 phi.read_dof_values_plain(solution.block(0));
515 phi.evaluate(EvaluationFlags::gradients);
516 for (unsigned int q = 0; q < phi.n_q_points; ++q)
517 {
519 phi.quadrature_point(q);
520
521 VectorizedArray<double> f_value = 0.0;
522 for (unsigned int v = 0; v < VectorizedArray<double>::size(); ++v)
523 {
525 for (unsigned int d = 0; d < spacedim; ++d)
526 p[d] = p_vect[d][v];
527 f_value[v] = par.rhs.value(p);
528 }
529 phi.submit_gradient(-phi.get_gradient(q), q);
530 phi.submit_value(f_value, q);
531 }
533 phi.distribute_local_to_global(system_rhs.block(0));
534 }
535
536 system_rhs.compress(VectorOperation::add);
537}
538#endif
539
540
541template <int dim, int spacedim>
542void
544{
546 pcout << "Preparing solve." << std::endl;
547 SolverCG<VectorType> cg_stiffness(par.inner_control);
548#ifdef MATRIX_FREE_PATH
549
550 using Payload = dealii::internal::LinearOperatorImplementation::EmptyPayload;
553
554 MGTransferMatrixFree<spacedim, float> mg_transfer(mg_constrained_dofs);
555 mg_transfer.build(dh);
556
557 using SmootherType =
558 PreconditionChebyshev<LevelMatrixType,
560 mg::SmootherRelaxation<SmootherType,
562 mg_smoother;
564 smoother_data.resize(0, tria.n_global_levels() - 1);
565 for (unsigned int level = 0; level < tria.n_global_levels(); ++level)
566 {
567 if (level > 0)
568 {
569 smoother_data[level].smoothing_range = 15.;
570 smoother_data[level].degree = 5;
571 smoother_data[level].eig_cg_n_iterations = 10;
572 }
573 else
574 {
575 smoother_data[0].smoothing_range = 1e-3;
576 smoother_data[0].degree = numbers::invalid_unsigned_int;
577 smoother_data[0].eig_cg_n_iterations = mg_matrices[0].m();
578 }
579 mg_matrices[level].compute_diagonal();
580 smoother_data[level].preconditioner =
581 mg_matrices[level].get_matrix_diagonal_inverse();
582 }
583 mg_smoother.initialize(mg_matrices, smoother_data);
584
586 mg_coarse;
587 mg_coarse.initialize(mg_smoother);
588
590
592 mg_interface_matrices;
593 mg_interface_matrices.resize(0, tria.n_global_levels() - 1);
594 for (unsigned int level = 0; level < tria.n_global_levels(); ++level)
595 mg_interface_matrices[level].initialize(mg_matrices[level]);
597 mg_interface_matrices);
598
600 mg_matrix, mg_coarse, mg_transfer, mg_smoother, mg_smoother);
601 mg.set_edge_matrices(mg_interface, mg_interface);
602
603 PreconditionMG<spacedim,
606 preconditioner(dh, mg, mg_transfer);
607
608 auto invA = A;
609 invA = inverse_operator(A, cg_stiffness, preconditioner);
610#else
611 using Payload =
615
616 LA::MPI::PreconditionAMG prec_A;
617 {
618 LA::MPI::PreconditionAMG::AdditionalData data;
619# ifdef USE_PETSC_LA
620 data.symmetric_operator = true;
621# endif
622 pcout << "Initialize AMG...";
623 prec_A.initialize(stiffness_matrix, data);
624 pcout << "done." << std::endl;
625 }
626 const auto amgA = linear_operator<VectorType, VectorType, Payload>(A, prec_A);
627 auto invA = A;
628 invA = inverse_operator(A, cg_stiffness, amgA);
629#endif
630
631
632 // Some aliases
633 auto &u = solution.block(0);
634 auto &lambda = solution.block(1);
635
636 const auto &f = system_rhs.block(0);
637 const auto &g = system_rhs.block(1);
638
639 if (reduced_coupling.get_dof_handler().n_dofs() == 0)
640 {
641 u = invA * f;
642 }
643 else
644 {
645#ifdef MATRIX_FREE_PATH
646 auto Bt =
648 Bt.reinit_range_vector = [this](VectorType &vec, const bool) {
649 vec.reinit(owned_dofs[0], relevant_dofs[0], mpi_communicator);
650 };
651 Bt.reinit_domain_vector = [this](VectorType &vec, const bool) {
652 vec.reinit(owned_dofs[1], relevant_dofs[1], mpi_communicator);
653 };
654
656#else
657 const auto Bt =
660 // const auto B = linear_operator<VectorType, VectorType, Payload>(
661 // coupling_matrix_transpose);
662#endif
663
664 if (par.solver_name == "Schur")
665 {
666 // Schur complement
667 pcout << " Prepare schur... ";
668 const auto S = B * invA * Bt;
669 pcout << "S was built." << std::endl;
670
671 // Schur complement preconditioner
672 auto invS = S;
673 SolverFGMRES<VectorType> solver_schur(par.outer_control);
674 invS =
676 solver_schur);
677
678 pcout << " f norm: " << f.l2_norm() << ", g norm: " << g.l2_norm()
679 << std::endl;
680
681 // Compute Lambda first
682 lambda = invS * (B * invA * f - g);
683 pcout << " Solved for lambda in " << par.outer_control.last_step()
684 << " iterations." << std::endl;
685
686 // Then compute u
687 u = invA * (f - Bt * lambda);
688 pcout << " u norm: " << u.l2_norm()
689 << ", lambda norm: " << lambda.l2_norm() << std::endl;
690
691 pcout << " Solved for u in " << par.inner_control.last_step()
692 << " iterations." << std::endl;
693
694 constraints.distribute(u);
695 reduced_coupling.get_coupling_constraints().distribute(lambda);
696 // solution.update_ghost_values();
698 }
699 else if (par.solver_name == "AL")
700 {
701 pcout << "Prepare AL preconditioner... " << std::endl;
702
704 (std::is_same_v<LA::MPI::Vector, TrilinosWrappers::MPI::Vector>),
705 ExcNotImplemented());
706
707 LA::MPI::SparseMatrix reduced_mass_matrix;
708 const auto &reduced_dh = reduced_coupling.get_dof_handler();
709 DynamicSparsityPattern dsp_reduced_mass(relevant_dofs[1]);
710 DoFTools::make_sparsity_pattern(reduced_dh, dsp_reduced_mass);
712 owned_dofs[1],
714 relevant_dofs[1]);
715 reduced_mass_matrix.reinit(owned_dofs[1],
716 owned_dofs[1],
717 dsp_reduced_mass,
719
720 // Create mass matrix associated with the reduced dof handler
721 reduced_coupling.assemble_coupling_mass_matrix(reduced_mass_matrix);
722
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;
726
727
728 const auto M = linear_operator<LA::MPI::Vector>(reduced_mass_matrix);
729
730 // Augmented Lagrangian solver
732 M_inv_ilu.initialize(reduced_mass_matrix);
733
734 TrilinosWrappers::MPI::Vector inverse_squares_reduced; // diag(M)^{-2}
735 inverse_squares_reduced.reinit(owned_dofs[1], mpi_communicator);
736 for (const types::global_dof_index local_idx : owned_dofs[1])
737 {
738 const double el = reduced_mass_matrix.diag_element(local_idx);
739 Assert(std::abs(el) > 1e-10,
740 ExcMessage(
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);
745 }
746
747 inverse_squares_reduced.compress(VectorOperation::insert);
748
749
750 SolverControl solver_control(100, 1e-15, false, false);
752 solver_control);
753 auto invM = inverse_operator(M, solver_mass_matrix, M_inv_ilu);
754 auto invW = invM * invM;
755
756 const double gamma = 10; // TODO: add to parameters file
757 auto Aug = A + gamma * Bt * invW * B;
758
759 TrilinosWrappers::SparseMatrix augmented_matrix;
760 pcout << "Building augmented matrix..." << std::endl;
763 inverse_squares_reduced,
764 gamma,
765 augmented_matrix);
766 pcout << "done." << std::endl;
767
768 TrilinosWrappers::PreconditionAMG prec_amg_augmented_block;
770 prec_amg_augmented_block.initialize(augmented_matrix, data);
771
772 auto Zero = M * 0.0;
774 {{{{Aug, Bt}}, {{B, Zero}}}});
775
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);
780
781 // lagrangian term
782 LA::MPI::Vector tmp;
783 tmp.reinit(system_rhs.block(0));
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); // ! augmented
787 system_rhs_block.block(1) = system_rhs.block(1);
788
789 SolverCG<LA::MPI::Vector> solver_lagrangian(par.inner_control);
790
791
792 auto Aug_inv =
793 inverse_operator(Aug, solver_lagrangian, prec_amg_augmented_block);
794 SolverFGMRES<LA::MPI::BlockVector> solver_fgmres(par.outer_control);
795
797 augmented_lagrangian_preconditioner{Aug_inv, B, Bt, invW, gamma};
798
799 solver_fgmres.solve(AA,
800 solution_block,
801 system_rhs_block,
802 augmented_lagrangian_preconditioner);
803
804 pcout << " Solved with AL preconditioner in "
805 << par.outer_control.last_step() << " iterations." << std::endl;
806
807 constraints.distribute(solution_block.block(0));
808 reduced_coupling.get_coupling_constraints().distribute(
809 solution_block.block(1));
810 // solution.update_ghost_values();
811 locally_relevant_solution = solution_block;
812
813
814#ifdef DEBUG
815 // Estimate condition number of BBt using CG
816 {
817 auto output_double_number = [this](double input,
818 const std::string &text) {
820 std::cout << text << input << std::endl;
821 };
822
823 // Estimate condition number:
824 pcout << "- - - - - - - - - - - - - - - - - - - - - - - -"
825 << std::endl;
826 pcout << "Estimate condition number of BBt using CG" << std::endl;
827 SolverControl solver_control(100000, 1e-12);
828 SolverCG<TrilinosWrappers::MPI::Vector> solver_cg(solver_control);
829
831 std::bind(output_double_number,
832 std::placeholders::_1,
833 "Condition number estimate: "));
834 using PayloadType = dealii::TrilinosWrappers::internal::
835 LinearOperatorImplementation::TrilinosPayload;
836
837 auto BBt = B * Bt;
838
840 u = 0.;
842 f = 1.;
844 try
845 {
846 solver_cg.solve(BBt, u, f, prec_no);
847 }
848 catch (...)
849 {
850 pcout
851 << "***BBt solve not successfull (see condition number above)***"
852 << std::endl;
853 }
854 }
855#endif
856 }
857 else
858 {
860 }
861 }
862}
863
864
865
866template <int dim, int spacedim>
867void
869{
871 Vector<float> error_per_cell(tria.n_active_cells());
873 QGauss<spacedim - 1>(par.fe_degree +
874 1),
875 {},
877 error_per_cell);
878 if (par.refinement_strategy == "fixed_fraction")
879 {
881 tria, error_per_cell, par.refinement_fraction, par.coarsening_fraction);
882 }
883 else if (par.refinement_strategy == "fixed_number")
884 {
886 tria,
887 error_per_cell,
888 par.refinement_fraction,
889 par.coarsening_fraction,
890 par.max_cells);
891 }
892 else if (par.refinement_strategy == "global")
893 for (const auto &cell : tria.active_cell_iterators())
894 cell->set_refine_flag();
895
897 tria.prepare_coarsening_and_refinement();
898 // inclusions.inclusions_as_particles.prepare_for_coarsening_and_refinement();
901 tria.execute_coarsening_and_refinement();
902 // inclusions.inclusions_as_particles.unpack_after_coarsening_and_refinement();
903 setup_dofs();
904 transfer.interpolate(solution.block(0));
905 constraints.distribute(solution.block(0));
906 locally_relevant_solution.block(0) = solution.block(0);
907}
908
909
910
911template <int dim, int spacedim>
912std::string
914{
915 TimerOutput::Scope t(computing_timer, "Output results");
916 std::string solution_name = "solution";
917 DataOut<spacedim> data_out;
918 data_out.attach_dof_handler(dh);
919 data_out.add_data_vector(locally_relevant_solution.block(0), solution_name);
920 Vector<float> subdomain(tria.n_active_cells());
921 for (unsigned int i = 0; i < subdomain.size(); ++i)
922 subdomain(i) = tria.locally_owned_subdomain();
923 data_out.add_data_vector(subdomain, "subdomain");
924 data_out.build_patches();
925 const std::string filename =
926 par.output_name + "_" + std::to_string(cycle) + ".vtu";
927 data_out.write_vtu_in_parallel(par.output_directory + "/" + filename,
929 return filename;
930}
931
932
933template <int dim, int spacedim>
934void
936{
937 static std::vector<std::pair<double, std::string>> cycles_and_solutions;
938 static std::vector<std::pair<double, std::string>> cycles_and_particles;
939
940 if (cycles_and_solutions.size() == cycle)
941 {
942 cycles_and_solutions.push_back({(double)cycle, output_solution()});
943
944 const std::string particles_filename =
945 par.output_name + "_particles_" + std::to_string(cycle) + ".vtu";
946 reduced_coupling.output_particles(par.output_directory + "/" +
947 particles_filename);
948
949 cycles_and_particles.push_back({(double)cycle, particles_filename});
950
951 std::ofstream pvd_solutions(par.output_directory + "/" + par.output_name +
952 ".pvd");
953 std::ofstream pvd_particles(par.output_directory + "/" + par.output_name +
954 "_particles.pvd");
955 DataOutBase::write_pvd_record(pvd_solutions, cycles_and_solutions);
956 DataOutBase::write_pvd_record(pvd_particles, cycles_and_particles);
957 }
958}
959
960template <int dim, int spacedim>
961void
963{
964#ifdef USE_PETSC_LA
965 pcout << "Running ReducedPoisson<" << Utilities::dim_string(dim, spacedim)
966 << "> using PETSc." << std::endl;
967#else
968 pcout << "Running ReducedPoisson<" << Utilities::dim_string(dim, spacedim)
969 << "> using Trilinos with "
970 << Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) << " MPI ranks."
971 << std::endl;
972#endif
974 {
975 par.prm.print_parameters(par.output_directory + "/" + "used_parameters_" +
976 std::to_string(dim) +
977 std::to_string(spacedim) + ".prm",
979 }
980}
981
982template <int dim, int spacedim>
983void
985{
987 make_grid();
988 setup_fe();
989 for (cycle = 0; cycle < par.n_refinement_cycles; ++cycle)
990 {
991 setup_dofs();
992 if (par.output_results_before_solving)
994#ifdef MATRIX_FREE_PATH
995 assemble_rhs();
996#else
998#endif
999 reduced_coupling.assemble_coupling_matrix(coupling_matrix,
1000 dh,
1001 constraints);
1002 reduced_coupling.assemble_reduced_rhs(system_rhs.block(1));
1003
1004#ifdef MATRIX_FREE_PATH
1005 // MappingQ1<spacedim> mapping;
1006 // coupling_operator =
1007 // std::make_unique<CouplingOperator<spacedim,
1008 // double>>(
1009 // inclusions, dh, constraints,
1010 // mapping, *fe);
1011#endif
1012 // return;
1013 solve();
1015 par.convergence_table.error_from_exact(dh,
1017 par.bc);
1018 if (cycle != par.n_refinement_cycles - 1)
1020 if (pcout.is_active())
1021 par.convergence_table.output_table(pcout.get_stream());
1022 }
1023}
1024
1025
1026// Template instantiations
1027template class ReducedPoissonParameters<2>;
1028template class ReducedPoissonParameters<3>;
1029
1030template class ReducedPoisson<2>;
1031template class ReducedPoisson<3>;
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 &section_name="")
static ParameterHandler prm
void enter_subsection(const std::string &subsection)
void add_parameter(const std::string &entry, ParameterType &parameter, 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
BlockVectorType solution
const ReducedPoissonParameters< spacedim > & par
std::vector< IndexSet > owned_dofs
parallel::distributed::Triangulation< spacedim > tria
AffineConstraints< double > constraints
ConditionalOStream pcout
unsigned int cycle
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
ParsedConvergenceTable convergence_table
std::list< types::boundary_id > dirichlet_ids
ParameterAcceptorProxy< ReductionControl > inner_control
ParameterAcceptorProxy< ReductionControl > outer_control
ParameterAcceptorProxy< Functions::ParsedFunction< spacedim > > rhs
#define DEAL_II_NOT_IMPLEMENTED()
unsigned int level
#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)
update_values
update_JxW_values
update_gradients
update_quadrature_points
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 > > &times_and_names)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
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)
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm mpi_comm, const IndexSet &locally_relevant_rows)
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 interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
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