Processing math: 25%
 Reference documentation for deal.II version 9.4.1
\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}} \newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=} \newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]} \newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
step-75.h
Go to the documentation of this file.
1) const override
465 * {
466 * const std::array<double, dim> p_sphere =
468 *
469 * constexpr const double alpha = 2. / 3.;
470 * return std::pow(p_sphere[0], alpha) * std::sin(alpha * p_sphere[1]);
471 * }
472 * };
473 *
474 *
475 *
476 * @endcode
477 *
478 *
479 * <a name="Parameters"></a>
480 * <h3>Parameters</h3>
481 *
482
483 *
484 * For this tutorial, we will use a simplified set of parameters. It is also
485 * possible to use a ParameterHandler class here, but to keep this tutorial
486 * short we decided on using simple structs. The actual intention of all these
487 * parameters will be described in the upcoming classes at their respective
488 * location where they are used.
489 *
490
491 *
492 * The following parameter set controls the coarse-grid solver, the smoothers,
493 * and the inter-grid transfer scheme of the multigrid mechanism.
494 * We populate it with default parameters.
495 *
496 * @code
497 * struct MultigridParameters
498 * {
499 * struct
500 * {
501 * std::string type = "cg_with_amg";
502 * unsigned int maxiter = 10000;
503 * double abstol = 1e-20;
504 * double reltol = 1e-4;
505 * unsigned int smoother_sweeps = 1;
506 * unsigned int n_cycles = 1;
507 * std::string smoother_type = "ILU";
508 * } coarse_solver;
509 *
510 * struct
511 * {
512 * std::string type = "chebyshev";
513 * double smoothing_range = 20;
514 * unsigned int degree = 5;
515 * unsigned int eig_cg_n_iterations = 20;
516 * } smoother;
517 *
518 * struct
519 * {
521 * p_sequence = MGTransferGlobalCoarseningTools::
522 * PolynomialCoarseningSequenceType::decrease_by_one;
523 * bool perform_h_transfer = true;
524 * } transfer;
525 * };
526 *
527 *
528 *
529 * @endcode
530 *
531 * This is the general parameter struct for the problem class. You will find
532 * this struct divided into several categories, including general runtime
533 * parameters, level limits, refine and coarsen fractions, as well as
534 * parameters for cell weighting. It also contains an instance of the above
535 * struct for multigrid parameters which will be passed to the multigrid
536 * algorithm.
537 *
538 * @code
539 * struct Parameters
540 * {
541 * unsigned int n_cycles = 8;
542 * double tolerance_factor = 1e-12;
543 *
544 * MultigridParameters mg_data;
545 *
546 * unsigned int min_h_level = 5;
547 * unsigned int max_h_level = 12;
548 * unsigned int min_p_degree = 2;
549 * unsigned int max_p_degree = 6;
550 * unsigned int max_p_level_difference = 1;
551 *
552 * double refine_fraction = 0.3;
553 * double coarsen_fraction = 0.03;
554 * double p_refine_fraction = 0.9;
555 * double p_coarsen_fraction = 0.9;
556 *
557 * double weighting_factor = 1.;
558 * double weighting_exponent = 1.;
559 * };
560 *
561 *
562 *
563 * @endcode
564 *
565 *
566 * <a name="MatrixfreeLaplaceoperator"></a>
567 * <h3>Matrix-free Laplace operator</h3>
568 *
569
570 *
571 * This is a matrix-free implementation of the Laplace operator that will
572 * basically take over the part of the `assemble_system()` function from other
573 * tutorials. The meaning of all member functions will be explained at their
574 * definition later.
575 *
576
577 *
578 * We will use the FEEvaluation class to evaluate the solution vector
579 * at the quadrature points and to perform the integration. In contrast to
580 * other tutorials, the template arguments `degree` is set to @f$-1@f$ and
581 * `number of quadrature in 1D` to @f$0@f$. In this case, FEEvaluation selects
582 * dynamically the correct polynomial degree and number of quadrature
583 * points. Here, we introduce an alias to FEEvaluation with the correct
584 * template parameters so that we do not have to worry about them later on.
585 *
586 * @code
587 * template <int dim, typename number>
588 * class LaplaceOperator : public Subscriptor
589 * {
590 * public:
592 *
593 * using FECellIntegrator = FEEvaluation<dim, -1, 0, 1, number>;
594 *
595 * LaplaceOperator() = default;
596 *
597 * LaplaceOperator(const hp::MappingCollection<dim> &mapping,
598 * const DoFHandler<dim> & dof_handler,
599 * const hp::QCollection<dim> & quad,
600 * const AffineConstraints<number> & constraints,
601 * VectorType & system_rhs);
602 *
603 * void reinit(const hp::MappingCollection<dim> &mapping,
604 * const DoFHandler<dim> & dof_handler,
605 * const hp::QCollection<dim> & quad,
606 * const AffineConstraints<number> & constraints,
607 * VectorType & system_rhs);
608 *
609 * types::global_dof_index m() const;
610 *
611 * number el(unsigned int, unsigned int) const;
612 *
613 * void initialize_dof_vector(VectorType &vec) const;
614 *
615 * void vmult(VectorType &dst, const VectorType &src) const;
616 *
617 * void Tvmult(VectorType &dst, const VectorType &src) const;
618 *
619 * const TrilinosWrappers::SparseMatrix &get_system_matrix() const;
620 *
621 * void compute_inverse_diagonal(VectorType &diagonal) const;
622 *
623 * private:
624 * void do_cell_integral_local(FECellIntegrator &integrator) const;
625 *
626 * void do_cell_integral_global(FECellIntegrator &integrator,
627 * VectorType & dst,
628 * const VectorType &src) const;
629 *
630 *
631 * void do_cell_integral_range(
632 * const MatrixFree<dim, number> & matrix_free,
633 * VectorType & dst,
634 * const VectorType & src,
635 * const std::pair<unsigned int, unsigned int> &range) const;
636 *
637 * MatrixFree<dim, number> matrix_free;
638 *
639 * @endcode
640 *
641 * To solve the equation system on the coarsest level with an AMG
642 * preconditioner, we need an actual system matrix on the coarsest level.
643 * For this purpose, we provide a mechanism that optionally computes a
644 * matrix from the matrix-free formulation, for which we introduce a
645 * dedicated SparseMatrix object. In the default case, this matrix stays
646 * empty. Once `get_system_matrix()` is called, this matrix is filled (lazy
647 * allocation). Since this is a `const` function, we need the "mutable"
648 * keyword here. We also need a the constraints object to build the matrix.
649 *
650 * @code
651 * AffineConstraints<number> constraints;
652 * mutable TrilinosWrappers::SparseMatrix system_matrix;
653 * };
654 *
655 *
656 *
657 * @endcode
658 *
659 * The following section contains functions to initialize and reinitialize
660 * the class. In particular, these functions initialize the internal
661 * MatrixFree instance. For sake of simplicity, we also compute the system
662 * right-hand-side vector.
663 *
664 * @code
665 * template <int dim, typename number>
666 * LaplaceOperator<dim, number>::LaplaceOperator(
667 * const hp::MappingCollection<dim> &mapping,
668 * const DoFHandler<dim> & dof_handler,
669 * const hp::QCollection<dim> & quad,
670 * const AffineConstraints<number> & constraints,
671 * VectorType & system_rhs)
672 * {
673 * this->reinit(mapping, dof_handler, quad, constraints, system_rhs);
674 * }
675 *
676 *
677 *
678 * template <int dim, typename number>
679 * void LaplaceOperator<dim, number>::reinit(
680 * const hp::MappingCollection<dim> &mapping,
681 * const DoFHandler<dim> & dof_handler,
682 * const hp::QCollection<dim> & quad,
683 * const AffineConstraints<number> & constraints,
684 * VectorType & system_rhs)
685 * {
686 * @endcode
687 *
688 * Clear internal data structures (in the case that the operator is reused).
689 *
690 * @code
691 * this->system_matrix.clear();
692 *
693 * @endcode
694 *
695 * Copy the constraints, since they might be needed for computation of the
696 * system matrix later on.
697 *
698 * @code
699 * this->constraints.copy_from(constraints);
700 *
701 * @endcode
702 *
703 * Set up MatrixFree. At the quadrature points, we only need to evaluate
704 * the gradient of the solution and test with the gradient of the shape
705 * functions so that we only need to set the flag `update_gradients`.
706 *
707 * @code
710 *
711 * matrix_free.reinit(mapping, dof_handler, constraints, quad, data);
712 *
713 * @endcode
714 *
715 * Compute the right-hand side vector. For this purpose, we set up a second
716 * MatrixFree instance that uses a modified AffineConstraints not containing
717 * the constraints due to Dirichlet-boundary conditions. This modified
718 * operator is applied to a vector with only the Dirichlet values set. The
719 * result is the negative right-hand-side vector.
720 *
721 * @code
722 * {
723 * AffineConstraints<number> constraints_without_dbc;
724 *
725 * const IndexSet locally_relevant_dofs =
727 * constraints_without_dbc.reinit(locally_relevant_dofs);
728 *
730 * constraints_without_dbc);
731 * constraints_without_dbc.close();
732 *
733 * VectorType b, x;
734 *
735 * this->initialize_dof_vector(system_rhs);
736 *
737 * MatrixFree<dim, number> matrix_free;
738 * matrix_free.reinit(
739 * mapping, dof_handler, constraints_without_dbc, quad, data);
740 *
741 * matrix_free.initialize_dof_vector(b);
742 * matrix_free.initialize_dof_vector(x);
743 *
744 * constraints.distribute(x);
745 *
746 * matrix_free.cell_loop(&LaplaceOperator::do_cell_integral_range,
747 * this,
748 * b,
749 * x);
750 *
751 * constraints.set_zero(b);
752 *
753 * system_rhs -= b;
754 * }
755 * }
756 *
757 *
758 *
759 * @endcode
760 *
761 * The following functions are implicitly needed by the multigrid algorithm,
762 * including the smoothers.
763 *
764
765 *
766 * Since we do not have a matrix, query the DoFHandler for the number of
767 * degrees of freedom.
768 *
769 * @code
770 * template <int dim, typename number>
771 * types::global_dof_index LaplaceOperator<dim, number>::m() const
772 * {
773 * return matrix_free.get_dof_handler().n_dofs();
774 * }
775 *
776 *
777 *
778 * @endcode
779 *
780 * Access a particular element in the matrix. This function is neither
781 * needed nor implemented, however, is required to compile the program.
782 *
783 * @code
784 * template <int dim, typename number>
785 * number LaplaceOperator<dim, number>::el(unsigned int, unsigned int) const
786 * {
787 * Assert(false, ExcNotImplemented());
788 * return 0;
789 * }
790 *
791 *
792 *
793 * @endcode
794 *
795 * Initialize the given vector. We simply delegate the task to the
796 * MatrixFree function with the same name.
797 *
798 * @code
799 * template <int dim, typename number>
800 * void
801 * LaplaceOperator<dim, number>::initialize_dof_vector(VectorType &vec) const
802 * {
803 * matrix_free.initialize_dof_vector(vec);
804 * }
805 *
806 *
807 *
808 * @endcode
809 *
810 * Perform an operator evaluation by looping with the help of MatrixFree
811 * over all cells and evaluating the effect of the cell integrals (see also:
812 * `do_cell_integral_local()` and `do_cell_integral_global()`).
813 *
814 * @code
815 * template <int dim, typename number>
816 * void LaplaceOperator<dim, number>::vmult(VectorType & dst,
817 * const VectorType &src) const
818 * {
819 * this->matrix_free.cell_loop(
820 * &LaplaceOperator::do_cell_integral_range, this, dst, src, true);
821 * }
822 *
823 *
824 *
825 * @endcode
826 *
827 * Perform the transposed operator evaluation. Since we are considering
828 * symmetric "matrices", this function can simply delegate it task to vmult().
829 *
830 * @code
831 * template <int dim, typename number>
832 * void LaplaceOperator<dim, number>::Tvmult(VectorType & dst,
833 * const VectorType &src) const
834 * {
835 * this->vmult(dst, src);
836 * }
837 *
838 *
839 *
840 * @endcode
841 *
842 * Since we do not have a system matrix, we cannot loop over the the
843 * diagonal entries of the matrix. Instead, we compute the diagonal by
844 * performing a sequence of operator evaluations to unit basis vectors.
845 * For this purpose, an optimized function from the MatrixFreeTools
846 * namespace is used. The inversion is performed manually afterwards.
847 *
848 * @code
849 * template <int dim, typename number>
850 * void LaplaceOperator<dim, number>::compute_inverse_diagonal(
851 * VectorType &diagonal) const
852 * {
853 * this->matrix_free.initialize_dof_vector(diagonal);
855 * diagonal,
856 * &LaplaceOperator::do_cell_integral_local,
857 * this);
858 *
859 * for (auto &i : diagonal)
860 * i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
861 * }
862 *
863 *
864 *
865 * @endcode
866 *
867 * In the matrix-free context, no system matrix is set up during
868 * initialization of this class. As a consequence, it has to be computed
869 * here if it should be requested. Since the matrix is only computed in
870 * this tutorial for linear elements (on the coarse grid), this is
871 * acceptable.
872 * The matrix entries are obtained via sequence of operator evaluations.
873 * For this purpose, the optimized function MatrixFreeTools::compute_matrix()
874 * is used. The matrix will only be computed if it has not been set up yet
875 * (lazy allocation).
876 *
877 * @code
878 * template <int dim, typename number>
880 * LaplaceOperator<dim, number>::get_system_matrix() const
881 * {
882 * if (system_matrix.m() == 0 && system_matrix.n() == 0)
883 * {
884 * const auto &dof_handler = this->matrix_free.get_dof_handler();
885 *
887 * dof_handler.locally_owned_dofs(),
888 * dof_handler.get_triangulation().get_communicator());
889 *
890 * DoFTools::make_sparsity_pattern(dof_handler, dsp, this->constraints);
891 *
892 * dsp.compress();
893 * system_matrix.reinit(dsp);
894 *
896 * matrix_free,
897 * constraints,
898 * system_matrix,
899 * &LaplaceOperator::do_cell_integral_local,
900 * this);
901 * }
902 *
903 * return this->system_matrix;
904 * }
905 *
906 *
907 *
908 * @endcode
909 *
910 * Perform cell integral on a cell batch without gathering and scattering
911 * the values. This function is needed for the MatrixFreeTools functions
912 * since these functions operate directly on the buffers of FEEvaluation.
913 *
914 * @code
915 * template <int dim, typename number>
916 * void LaplaceOperator<dim, number>::do_cell_integral_local(
917 * FECellIntegrator &integrator) const
918 * {
920 *
921 * for (unsigned int q = 0; q < integrator.n_q_points; ++q)
922 * integrator.submit_gradient(integrator.get_gradient(q), q);
923 *
924 * integrator.integrate(EvaluationFlags::gradients);
925 * }
926 *
927 *
928 *
929 * @endcode
930 *
931 * Same as above but with access to the global vectors.
932 *
933 * @code
934 * template <int dim, typename number>
935 * void LaplaceOperator<dim, number>::do_cell_integral_global(
936 * FECellIntegrator &integrator,
937 * VectorType & dst,
938 * const VectorType &src) const
939 * {
940 * integrator.gather_evaluate(src, EvaluationFlags::gradients);
941 *
942 * for (unsigned int q = 0; q < integrator.n_q_points; ++q)
943 * integrator.submit_gradient(integrator.get_gradient(q), q);
944 *
945 * integrator.integrate_scatter(EvaluationFlags::gradients, dst);
946 * }
947 *
948 *
949 *
950 * @endcode
951 *
952 * This function loops over all cell batches within a cell-batch range and
953 * calls the above function.
954 *
955 * @code
956 * template <int dim, typename number>
957 * void LaplaceOperator<dim, number>::do_cell_integral_range(
958 * const MatrixFree<dim, number> & matrix_free,
959 * VectorType & dst,
960 * const VectorType & src,
961 * const std::pair<unsigned int, unsigned int> &range) const
962 * {
963 * FECellIntegrator integrator(matrix_free, range);
964 *
965 * for (unsigned cell = range.first; cell < range.second; ++cell)
966 * {
967 * integrator.reinit(cell);
968 *
969 * do_cell_integral_global(integrator, dst, src);
970 * }
971 * }
972 *
973 *
974 *
975 * @endcode
976 *
977 *
978 * <a name="Solverandpreconditioner"></a>
979 * <h3>Solver and preconditioner</h3>
980 *
981
982 *
983 *
984 * <a name="Conjugategradientsolverwithmultigridpreconditioner"></a>
985 * <h4>Conjugate-gradient solver with multigrid preconditioner</h4>
986 *
987
988 *
989 * This function solves the equation system with a sequence of provided
990 * multigrid objects. It is meant to be treated as general as possible, hence
991 * the multitude of template parameters.
992 *
993 * @code
994 * template <typename VectorType,
995 * int dim,
996 * typename SystemMatrixType,
997 * typename LevelMatrixType,
998 * typename MGTransferType>
999 * static void
1000 * mg_solve(SolverControl & solver_control,
1001 * VectorType & dst,
1002 * const VectorType & src,
1003 * const MultigridParameters &mg_data,
1004 * const DoFHandler<dim> & dof,
1005 * const SystemMatrixType & fine_matrix,
1006 * const MGLevelObject<std::unique_ptr<LevelMatrixType>> &mg_matrices,
1007 * const MGTransferType & mg_transfer)
1008 * {
1009 * AssertThrow(mg_data.coarse_solver.type == "cg_with_amg",
1010 * ExcNotImplemented());
1011 * AssertThrow(mg_data.smoother.type == "chebyshev", ExcNotImplemented());
1012 *
1013 * const unsigned int min_level = mg_matrices.min_level();
1014 * const unsigned int max_level = mg_matrices.max_level();
1015 *
1016 * using SmootherPreconditionerType = DiagonalMatrix<VectorType>;
1017 * using SmootherType = PreconditionChebyshev<LevelMatrixType,
1018 * VectorType,
1019 * SmootherPreconditionerType>;
1020 * using PreconditionerType = PreconditionMG<dim, VectorType, MGTransferType>;
1021 *
1022 * @endcode
1023 *
1024 * We initialize level operators and Chebyshev smoothers here.
1025 *
1026 * @code
1027 * mg::Matrix<VectorType> mg_matrix(mg_matrices);
1028 *
1030 * min_level, max_level);
1031 *
1032 * for (unsigned int level = min_level; level <= max_level; ++level)
1033 * {
1034 * smoother_data[level].preconditioner =
1035 * std::make_shared<SmootherPreconditionerType>();
1036 * mg_matrices[level]->compute_inverse_diagonal(
1037 * smoother_data[level].preconditioner->get_vector());
1038 * smoother_data[level].smoothing_range = mg_data.smoother.smoothing_range;
1039 * smoother_data[level].degree = mg_data.smoother.degree;
1040 * smoother_data[level].eig_cg_n_iterations =
1041 * mg_data.smoother.eig_cg_n_iterations;
1042 * }
1043 *
1045 * mg_smoother;
1046 * mg_smoother.initialize(mg_matrices, smoother_data);
1047 *
1048 * @endcode
1049 *
1050 * Next, we initialize the coarse-grid solver. We use conjugate-gradient
1051 * method with AMG as preconditioner.
1052 *
1053 * @code
1054 * ReductionControl coarse_grid_solver_control(mg_data.coarse_solver.maxiter,
1055 * mg_data.coarse_solver.abstol,
1056 * mg_data.coarse_solver.reltol,
1057 * false,
1058 * false);
1059 * SolverCG<VectorType> coarse_grid_solver(coarse_grid_solver_control);
1060 *
1061 * std::unique_ptr<MGCoarseGridBase<VectorType>> mg_coarse;
1062 *
1063 * TrilinosWrappers::PreconditionAMG precondition_amg;
1065 * amg_data.smoother_sweeps = mg_data.coarse_solver.smoother_sweeps;
1066 * amg_data.n_cycles = mg_data.coarse_solver.n_cycles;
1067 * amg_data.smoother_type = mg_data.coarse_solver.smoother_type.c_str();
1068 *
1069 * precondition_amg.initialize(mg_matrices[min_level]->get_system_matrix(),
1070 * amg_data);
1071 *
1072 * mg_coarse =
1073 * std::make_unique<MGCoarseGridIterativeSolver<VectorType,
1075 * LevelMatrixType,
1076 * decltype(precondition_amg)>>(
1077 * coarse_grid_solver, *mg_matrices[min_level], precondition_amg);
1078 *
1079 * @endcode
1080 *
1081 * Finally, we create the Multigrid object, convert it to a preconditioner,
1082 * and use it inside of a conjugate-gradient solver to solve the linear
1083 * system of equations.
1084 *
1085 * @code
1087 * mg_matrix, *mg_coarse, mg_transfer, mg_smoother, mg_smoother);
1088 *
1089 * PreconditionerType preconditioner(dof, mg, mg_transfer);
1090 *
1091 * SolverCG<VectorType>(solver_control)
1092 * .solve(fine_matrix, dst, src, preconditioner);
1093 * }
1094 *
1095 *
1096 *
1097 * @endcode
1098 *
1099 *
1100 * <a name="Hybridpolynomialgeometricglobalcoarseningmultigridpreconditioner"></a>
1101 * <h4>Hybrid polynomial/geometric-global-coarsening multigrid preconditioner</h4>
1102 *
1103
1104 *
1105 * The above function deals with the actual solution for a given sequence of
1106 * multigrid objects. This functions creates the actual multigrid levels, in
1107 * particular the operators, and the transfer operator as a
1109 *
1110 * @code
1111 * template <typename VectorType, typename OperatorType, int dim>
1112 * void solve_with_gmg(SolverControl & solver_control,
1113 * const OperatorType & system_matrix,
1114 * VectorType & dst,
1115 * const VectorType & src,
1116 * const MultigridParameters & mg_data,
1117 * const hp::MappingCollection<dim> mapping_collection,
1118 * const DoFHandler<dim> & dof_handler,
1119 * const hp::QCollection<dim> & quadrature_collection)
1120 * {
1121 * @endcode
1122 *
1123 * Create a DoFHandler and operator for each multigrid level,
1124 * as well as, create transfer operators. To be able to
1125 * set up the operators, we need a set of DoFHandler that we create
1126 * via global coarsening of p or h. For latter, we need also a sequence
1127 * of Triangulation objects that are obtained by
1129 *
1130
1131 *
1132 * In case no h-transfer is requested, we provide an empty deleter for the
1133 * `emplace_back()` function, since the Triangulation of our DoFHandler is
1134 * an external field and its destructor is called somewhere else.
1135 *
1136 * @code
1137 * MGLevelObject<DoFHandler<dim>> dof_handlers;
1140 *
1141 * std::vector<std::shared_ptr<const Triangulation<dim>>>
1142 * coarse_grid_triangulations;
1143 * if (mg_data.transfer.perform_h_transfer)
1144 * coarse_grid_triangulations =
1146 * dof_handler.get_triangulation());
1147 * else
1148 * coarse_grid_triangulations.emplace_back(
1149 * &(dof_handler.get_triangulation()), [](auto *) {});
1150 *
1151 * @endcode
1152 *
1153 * Determine the total number of levels for the multigrid operation and
1154 * allocate sufficient memory for all levels.
1155 *
1156 * @code
1157 * const unsigned int n_h_levels = coarse_grid_triangulations.size() - 1;
1158 *
1159 * const auto get_max_active_fe_degree = [&](const auto &dof_handler) {
1160 * unsigned int max = 0;
1161 *
1162 * for (auto &cell : dof_handler.active_cell_iterators())
1163 * if (cell->is_locally_owned())
1164 * max =
1165 * std::max(max, dof_handler.get_fe(cell->active_fe_index()).degree);
1166 *
1167 * return Utilities::MPI::max(max, MPI_COMM_WORLD);
1168 * };
1169 *
1170 * const unsigned int n_p_levels =
1172 * get_max_active_fe_degree(dof_handler), mg_data.transfer.p_sequence)
1173 * .size();
1174 *
1175 * std::map<unsigned int, unsigned int> fe_index_for_degree;
1176 * for (unsigned int i = 0; i < dof_handler.get_fe_collection().size(); ++i)
1177 * {
1178 * const unsigned int degree = dof_handler.get_fe(i).degree;
1179 * Assert(fe_index_for_degree.find(degree) == fe_index_for_degree.end(),
1180 * ExcMessage("FECollection does not contain unique degrees."));
1181 * fe_index_for_degree[degree] = i;
1182 * }
1183 *
1184 * unsigned int minlevel = 0;
1185 * unsigned int maxlevel = n_h_levels + n_p_levels - 1;
1186 *
1187 * dof_handlers.resize(minlevel, maxlevel);
1188 * operators.resize(minlevel, maxlevel);
1189 * transfers.resize(minlevel, maxlevel);
1190 *
1191 * @endcode
1192 *
1193 * Loop from the minimum (coarsest) to the maximum (finest) level and set up
1194 * DoFHandler accordingly. We start with the h-levels, where we distribute
1195 * on increasingly finer meshes linear elements.
1196 *
1197 * @code
1198 * for (unsigned int l = 0; l < n_h_levels; ++l)
1199 * {
1200 * dof_handlers[l].reinit(*coarse_grid_triangulations[l]);
1201 * dof_handlers[l].distribute_dofs(dof_handler.get_fe_collection());
1202 * }
1203 *
1204 * @endcode
1205 *
1206 * After we reached the finest mesh, we will adjust the polynomial degrees
1207 * on each level. We reverse iterate over our data structure and start at
1208 * the finest mesh that contains all information about the active FE
1209 * indices. We then lower the polynomial degree of each cell level by level.
1210 *
1211 * @code
1212 * for (unsigned int i = 0, l = maxlevel; i < n_p_levels; ++i, --l)
1213 * {
1214 * dof_handlers[l].reinit(dof_handler.get_triangulation());
1215 *
1216 * if (l == maxlevel) // finest level
1217 * {
1218 * auto &dof_handler_mg = dof_handlers[l];
1219 *
1220 * auto cell_other = dof_handler.begin_active();
1221 * for (auto &cell : dof_handler_mg.active_cell_iterators())
1222 * {
1223 * if (cell->is_locally_owned())
1224 * cell->set_active_fe_index(cell_other->active_fe_index());
1225 * cell_other++;
1226 * }
1227 * }
1228 * else // coarse level
1229 * {
1230 * auto &dof_handler_fine = dof_handlers[l + 1];
1231 * auto &dof_handler_coarse = dof_handlers[l + 0];
1232 *
1233 * auto cell_other = dof_handler_fine.begin_active();
1234 * for (auto &cell : dof_handler_coarse.active_cell_iterators())
1235 * {
1236 * if (cell->is_locally_owned())
1237 * {
1238 * const unsigned int next_degree =
1241 * cell_other->get_fe().degree,
1242 * mg_data.transfer.p_sequence);
1243 * Assert(fe_index_for_degree.find(next_degree) !=
1244 * fe_index_for_degree.end(),
1245 * ExcMessage("Next polynomial degree in sequence "
1246 * "does not exist in FECollection."));
1247 *
1248 * cell->set_active_fe_index(fe_index_for_degree[next_degree]);
1249 * }
1250 * cell_other++;
1251 * }
1252 * }
1253 *
1254 * dof_handlers[l].distribute_dofs(dof_handler.get_fe_collection());
1255 * }
1256 *
1257 * @endcode
1258 *
1259 * Next, we will create all data structures additionally needed on each
1260 * multigrid level. This involves determining constraints with homogeneous
1261 * Dirichlet boundary conditions, and building the operator just like on the
1262 * active level.
1263 *
1264 * @code
1266 * constraints(minlevel, maxlevel);
1267 *
1268 * for (unsigned int level = minlevel; level <= maxlevel; ++level)
1269 * {
1270 * const auto &dof_handler = dof_handlers[level];
1271 * auto & constraint = constraints[level];
1272 *
1273 * const IndexSet locally_relevant_dofs =
1275 * constraint.reinit(locally_relevant_dofs);
1276 *
1277 * DoFTools::make_hanging_node_constraints(dof_handler, constraint);
1278 * VectorTools::interpolate_boundary_values(mapping_collection,
1279 * dof_handler,
1280 * 0,
1282 * constraint);
1283 * constraint.close();
1284 *
1285 * VectorType dummy;
1286 *
1287 * operators[level] = std::make_unique<OperatorType>(mapping_collection,
1288 * dof_handler,
1289 * quadrature_collection,
1290 * constraint,
1291 * dummy);
1292 * }
1293 *
1294 * @endcode
1295 *
1296 * Set up intergrid operators and collect transfer operators within a single
1297 * operator as needed by the Multigrid solver class.
1298 *
1299 * @code
1300 * for (unsigned int level = minlevel; level < maxlevel; ++level)
1301 * transfers[level + 1].reinit(dof_handlers[level + 1],
1302 * dof_handlers[level],
1303 * constraints[level + 1],
1304 * constraints[level]);
1305 *
1307 * transfers, [&](const auto l, auto &vec) {
1308 * operators[l]->initialize_dof_vector(vec);
1309 * });
1310 *
1311 * @endcode
1312 *
1313 * Finally, proceed to solve the problem with multigrid.
1314 *
1315 * @code
1316 * mg_solve(solver_control,
1317 * dst,
1318 * src,
1319 * mg_data,
1320 * dof_handler,
1321 * system_matrix,
1322 * operators,
1323 * transfer);
1324 * }
1325 *
1326 *
1327 *
1328 * @endcode
1329 *
1330 *
1331 * <a name="ThecodeLaplaceProblemcodeclasstemplate"></a>
1332 * <h3>The <code>LaplaceProblem</code> class template</h3>
1333 *
1334
1335 *
1336 * Now we will finally declare the main class of this program, which solves
1337 * the Laplace equation on subsequently refined function spaces. Its structure
1338 * will look familiar as it is similar to the main classes of @ref step_27 "step-27" and
1339 * @ref step_40 "step-40". There are basically just two additions:
1340 * - The SparseMatrix object that would hold the system matrix has been
1341 * replaced by an object of the LaplaceOperator class for the MatrixFree
1342 * formulation.
1343 * - An object of parallel::CellWeights, which will help us with load
1344 * balancing, has been added.
1345 *
1346 * @code
1347 * template <int dim>
1348 * class LaplaceProblem
1349 * {
1350 * public:
1351 * LaplaceProblem(const Parameters &parameters);
1352 *
1353 * void run();
1354 *
1355 * private:
1356 * void initialize_grid();
1357 * void setup_system();
1358 * void print_diagnostics();
1359 * void solve_system();
1360 * void compute_indicators();
1361 * void adapt_resolution();
1362 * void output_results(const unsigned int cycle);
1363 *
1364 * MPI_Comm mpi_communicator;
1365 *
1366 * const Parameters prm;
1367 *
1369 * DoFHandler<dim> dof_handler;
1370 *
1371 * hp::MappingCollection<dim> mapping_collection;
1372 * hp::FECollection<dim> fe_collection;
1373 * hp::QCollection<dim> quadrature_collection;
1374 * hp::QCollection<dim - 1> face_quadrature_collection;
1375 *
1376 * IndexSet locally_owned_dofs;
1377 * IndexSet locally_relevant_dofs;
1378 *
1379 * AffineConstraints<double> constraints;
1380 *
1381 * LaplaceOperator<dim, double> laplace_operator;
1382 * LinearAlgebra::distributed::Vector<double> locally_relevant_solution;
1384 *
1385 * std::unique_ptr<FESeries::Legendre<dim>> legendre;
1386 * std::unique_ptr<parallel::CellWeights<dim>> cell_weights;
1387 *
1388 * Vector<float> estimated_error_per_cell;
1389 * Vector<float> hp_decision_indicators;
1390 *
1391 * ConditionalOStream pcout;
1392 * TimerOutput computing_timer;
1393 * };
1394 *
1395 *
1396 *
1397 * @endcode
1398 *
1399 *
1400 * <a name="ThecodeLaplaceProblemcodeclassimplementation"></a>
1401 * <h3>The <code>LaplaceProblem</code> class implementation</h3>
1402 *
1403
1404 *
1405 *
1406 * <a name="Constructor"></a>
1407 * <h4>Constructor</h4>
1408 *
1409
1410 *
1411 * The constructor starts with an initializer list that looks similar to the
1412 * one of @ref step_40 "step-40". We again prepare the ConditionalOStream object to allow
1413 * only the first process to output anything over the console, and initialize
1414 * the computing timer properly.
1415 *
1416 * @code
1417 * template <int dim>
1418 * LaplaceProblem<dim>::LaplaceProblem(const Parameters &parameters)
1419 * : mpi_communicator(MPI_COMM_WORLD)
1420 * , prm(parameters)
1421 * , triangulation(mpi_communicator)
1422 * , dof_handler(triangulation)
1423 * , pcout(std::cout,
1424 * (Utilities::MPI::this_mpi_process(mpi_communicator) == 0))
1425 * , computing_timer(mpi_communicator,
1426 * pcout,
1429 * {
1430 * Assert(prm.min_h_level <= prm.max_h_level,
1431 * ExcMessage(
1432 * "Triangulation level limits have been incorrectly set up."));
1433 * Assert(prm.min_p_degree <= prm.max_p_degree,
1434 * ExcMessage("FECollection degrees have been incorrectly set up."));
1435 *
1436 * @endcode
1437 *
1438 * We need to prepare the data structures for the hp-functionality in the
1439 * actual body of the constructor, and create corresponding objects for
1440 * every degree in the specified range from the parameter struct. As we are
1441 * only dealing with non-distorted rectangular cells, a linear mapping
1442 * object is sufficient in this context.
1443 *
1444
1445 *
1446 * In the Parameters struct, we provide ranges for levels on which the
1447 * function space is operating with a reasonable resolution. The multigrid
1448 * algorithm requires linear elements on the coarsest possible level. So we
1449 * start with the lowest polynomial degree and fill the collection with
1450 * consecutively higher degrees until the user-specified maximum is
1451 * reached.
1452 *
1453 * @code
1454 * mapping_collection.push_back(MappingQ1<dim>());
1455 *
1456 * for (unsigned int degree = 1; degree <= prm.max_p_degree; ++degree)
1457 * {
1458 * fe_collection.push_back(FE_Q<dim>(degree));
1459 * quadrature_collection.push_back(QGauss<dim>(degree + 1));
1460 * face_quadrature_collection.push_back(QGauss<dim - 1>(degree + 1));
1461 * }
1462 *
1463 * @endcode
1464 *
1465 * As our FECollection contains more finite elements than we want to use for
1466 * the finite element approximation of our solution, we would like to limit
1467 * the range on which active FE indices can operate on. For this, the
1468 * FECollection class allows to register a hierarchy that determines the
1469 * succeeding and preceding finite element in case of of p-refinement and
1470 * p-coarsening, respectively. All functions in the hp::Refinement namespace
1471 * consult this hierarchy to determine future FE indices. We will register
1472 * such a hierarchy that only works on finite elements with polynomial
1473 * degrees in the proposed range <code>[min_p_degree, max_p_degree]</code>.
1474 *
1475 * @code
1476 * const unsigned int min_fe_index = prm.min_p_degree - 1;
1477 * fe_collection.set_hierarchy(
1478 * /*next_index=*/
1479 * [](const typename hp::FECollection<dim> &fe_collection,
1480 * const unsigned int fe_index) -> unsigned int {
1481 * return ((fe_index + 1) < fe_collection.size()) ? fe_index + 1 :
1482 * fe_index;
1483 * },
1484 * /*previous_index=*/
1485 * [min_fe_index](const typename hp::FECollection<dim> &,
1486 * const unsigned int fe_index) -> unsigned int {
1487 * Assert(fe_index >= min_fe_index,
1488 * ExcMessage("Finite element is not part of hierarchy!"));
1489 * return (fe_index > min_fe_index) ? fe_index - 1 : fe_index;
1490 * });
1491 *
1492 * @endcode
1493 *
1494 * We initialize the FESeries::Legendre object in the default configuration
1495 * for smoothness estimation.
1496 *
1497 * @code
1498 * legendre = std::make_unique<FESeries::Legendre<dim>>(
1500 *
1501 * @endcode
1502 *
1503 * The next part is going to be tricky. During execution of refinement, a
1504 * few hp-algorithms need to interfere with the actual refinement process on
1505 * the Triangulation object. We do this by connecting several functions to
1506 * Triangulation::Signals: signals will be called at different stages during
1507 * the actual refinement process and trigger all connected functions. We
1508 * require this functionality for load balancing and to limit the polynomial
1509 * degrees of neighboring cells.
1510 *
1511
1512 *
1513 * For the former, we would like to assign a weight to every cell that is
1514 * proportional to the number of degrees of freedom of its future finite
1515 * element. The library offers a class parallel::CellWeights that allows to
1516 * easily attach individual weights at the right place during the refinement
1517 * process, i.e., after all refine and coarsen flags have been set correctly
1518 * for hp-adaptation and right before repartitioning for load balancing is
1519 * about to happen. Functions can be registered that will attach weights in
1520 * the form that @f$a (n_\text{dofs})^b@f$ with a provided pair of parameters
1521 * @f$(a,b)@f$. We register such a function in the following.
1522 *
1523
1524 *
1525 * For load balancing, efficient solvers like the one we use should scale
1526 * linearly with the number of degrees of freedom owned. We set the
1527 * parameters for cell weighting correspondingly: A weighting factor of @f$1@f$
1528 * and an exponent of @f$1@f$ (see the definitions of the `weighting_factor` and
1529 * `weighting_exponent` above).
1530 *
1531 * @code
1532 * cell_weights = std::make_unique<parallel::CellWeights<dim>>(
1533 * dof_handler,
1535 * {prm.weighting_factor, prm.weighting_exponent}));
1536 *
1537 * @endcode
1538 *
1539 * In h-adaptive applications, we ensure a 2:1 mesh balance by limiting the
1540 * difference of refinement levels of neighboring cells to one. With the
1541 * second call in the following code snippet, we will ensure the same for
1542 * p-levels on neighboring cells: levels of future finite elements are not
1543 * allowed to differ by more than a specified difference. The function
1544 * hp::Refinement::limit_p_level_difference takes care of this, but needs to
1545 * be connected to a very specific signal in the parallel context. The issue
1546 * is that we need to know how the mesh will be actually refined to set
1547 * future FE indices accordingly. As we ask the p4est oracle to perform
1548 * refinement, we need to ensure that the Triangulation has been updated
1549 * with the adaptation flags of the oracle first. An instantiation of
1551 * that for the duration of its life. Thus, we will create an object of this
1552 * class right before limiting the p-level difference, and connect the
1553 * corresponding lambda function to the signal
1554 * Triangulation::Signals::post_p4est_refinement, which will be triggered
1555 * after the oracle got refined, but before the Triangulation is refined.
1556 * Furthermore, we specify that this function will be connected to the front
1557 * of the signal, to ensure that the modification is performed before any
1558 * other function connected to the same signal.
1559 *
1560 * @code
1561 * triangulation.signals.post_p4est_refinement.connect(
1562 * [&, min_fe_index]() {
1564 * refine_modifier(triangulation);
1566 * prm.max_p_level_difference,
1567 * /*contains=*/min_fe_index);
1568 * },
1569 * boost::signals2::at_front);
1570 * }
1571 *
1572 *
1573 *
1574 * @endcode
1575 *
1576 *
1577 * <a name="LaplaceProbleminitialize_grid"></a>
1578 * <h4>LaplaceProblem::initialize_grid</h4>
1579 *
1580
1581 *
1582 * For a L-shaped domain, we could use the function GridGenerator::hyper_L()
1583 * as demonstrated in @ref step_50 "step-50". However in the 2D case, that particular
1584 * function removes the first quadrant, while we need the fourth quadrant
1585 * removed in our scenario. Thus, we will use a different function
1586 * GridGenerator::subdivided_hyper_L() which gives us more options to create
1587 * the mesh. Furthermore, we formulate that function in a way that it also
1588 * generates a 3D mesh: the 2D L-shaped domain will basically elongated by 1
1589 * in the positive z-direction.
1590 *
1591
1592 *
1593 * We first pretend to build a GridGenerator::subdivided_hyper_rectangle().
1594 * The parameters that we need to provide are Point objects for the lower left
1595 * and top right corners, as well as the number of repetitions that the base
1596 * mesh will have in each direction. We provide them for the first two
1597 * dimensions and treat the higher third dimension separately.
1598 *
1599
1600 *
1601 * To create a L-shaped domain, we need to remove the excess cells. For this,
1602 * we specify the <code>cells_to_remove</code> accordingly. We would like to
1603 * remove one cell in every cell from the negative direction, but remove one
1604 * from the positive x-direction.
1605 *
1606
1607 *
1608 * On the coarse grid, we set the initial active FE indices and distribute the
1609 * degrees of freedom once. We do that in order to assign the hp::FECollection
1610 * to the DoFHandler, so that all cells know how many DoFs they are going to
1611 * have. This step is mandatory for the weighted load balancing algorithm,
1612 * which will be called implicitly in
1613 * parallel::distributed::Triangulation::refine_global().
1614 *
1615 * @code
1616 * template <int dim>
1617 * void LaplaceProblem<dim>::initialize_grid()
1618 * {
1619 * TimerOutput::Scope t(computing_timer, "initialize grid");
1620 *
1621 * std::vector<unsigned int> repetitions(dim);
1622 * Point<dim> bottom_left, top_right;
1623 * for (unsigned int d = 0; d < dim; ++d)
1624 * if (d < 2)
1625 * {
1626 * repetitions[d] = 2;
1627 * bottom_left[d] = -1.;
1628 * top_right[d] = 1.;
1629 * }
1630 * else
1631 * {
1632 * repetitions[d] = 1;
1633 * bottom_left[d] = 0.;
1634 * top_right[d] = 1.;
1635 * }
1636 *
1637 * std::vector<int> cells_to_remove(dim, 1);
1638 * cells_to_remove[0] = -1;
1639 *
1641 * triangulation, repetitions, bottom_left, top_right, cells_to_remove);
1642 *
1643 * const unsigned int min_fe_index = prm.min_p_degree - 1;
1644 * for (const auto &cell : dof_handler.active_cell_iterators())
1645 * if (cell->is_locally_owned())
1646 * cell->set_active_fe_index(min_fe_index);
1647 *
1648 * dof_handler.distribute_dofs(fe_collection);
1649 *
1650 * triangulation.refine_global(prm.min_h_level);
1651 * }
1652 *
1653 *
1654 *
1655 * @endcode
1656 *
1657 *
1658 * <a name="LaplaceProblemsetup_system"></a>
1659 * <h4>LaplaceProblem::setup_system</h4>
1660 *
1661
1662 *
1663 * This function looks exactly the same to the one of @ref step_40 "step-40", but you will
1664 * notice the absence of the system matrix as well as the scaffold that
1665 * surrounds it. Instead, we will initialize the MatrixFree formulation of the
1666 * <code>laplace_operator</code> here. For boundary conditions, we will use
1667 * the Solution class introduced earlier in this tutorial.
1668 *
1669 * @code
1670 * template <int dim>
1671 * void LaplaceProblem<dim>::setup_system()
1672 * {
1673 * TimerOutput::Scope t(computing_timer, "setup system");
1674 *
1675 * dof_handler.distribute_dofs(fe_collection);
1676 *
1677 * locally_owned_dofs = dof_handler.locally_owned_dofs();
1678 * locally_relevant_dofs =
1680 *
1681 * locally_relevant_solution.reinit(locally_owned_dofs,
1682 * locally_relevant_dofs,
1683 * mpi_communicator);
1684 * system_rhs.reinit(locally_owned_dofs, mpi_communicator);
1685 *
1686 * constraints.clear();
1687 * constraints.reinit(locally_relevant_dofs);
1688 * DoFTools::make_hanging_node_constraints(dof_handler, constraints);
1690 * mapping_collection, dof_handler, 0, Solution<dim>(), constraints);
1691 * constraints.close();
1692 *
1693 * laplace_operator.reinit(mapping_collection,
1694 * dof_handler,
1695 * quadrature_collection,
1696 * constraints,
1697 * system_rhs);
1698 * }
1699 *
1700 *
1701 *
1702 * @endcode
1703 *
1704 *
1705 * <a name="LaplaceProblemprint_diagnostics"></a>
1706 * <h4>LaplaceProblem::print_diagnostics</h4>
1707 *
1708
1709 *
1710 * This is a function that prints additional diagnostics about the equation
1711 * system and its partitioning. In addition to the usual global number of
1712 * active cells and degrees of freedom, we also output their local
1713 * equivalents. For a regulated output, we will communicate the local
1714 * quantities with a Utilities::MPI::gather operation to the first process
1715 * which will then output all information. Output of local quantities is
1716 * limited to the first 8 processes to avoid cluttering the terminal.
1717 *
1718
1719 *
1720 * Furthermore, we would like to print the frequencies of the polynomial
1721 * degrees in the numerical discretization. Since this information is only
1722 * stored locally, we will count the finite elements on locally owned cells
1723 * and later communicate them via Utilities::MPI::sum.
1724 *
1725 * @code
1726 * template <int dim>
1727 * void LaplaceProblem<dim>::print_diagnostics()
1728 * {
1729 * const unsigned int first_n_processes =
1730 * std::min<unsigned int>(8,
1731 * Utilities::MPI::n_mpi_processes(mpi_communicator));
1732 * const bool output_cropped =
1733 * first_n_processes < Utilities::MPI::n_mpi_processes(mpi_communicator);
1734 *
1735 * {
1736 * pcout << " Number of active cells: "
1737 * << triangulation.n_global_active_cells() << std::endl
1738 * << " by partition: ";
1739 *
1740 * std::vector<unsigned int> n_active_cells_per_subdomain =
1741 * Utilities::MPI::gather(mpi_communicator,
1742 * triangulation.n_locally_owned_active_cells());
1743 * for (unsigned int i = 0; i < first_n_processes; ++i)
1744 * pcout << ' ' << n_active_cells_per_subdomain[i];
1745 * if (output_cropped)
1746 * pcout << " ...";
1747 * pcout << std::endl;
1748 * }
1749 *
1750 * {
1751 * pcout << " Number of degrees of freedom: " << dof_handler.n_dofs()
1752 * << std::endl
1753 * << " by partition: ";
1754 *
1755 * std::vector<types::global_dof_index> n_dofs_per_subdomain =
1756 * Utilities::MPI::gather(mpi_communicator,
1757 * dof_handler.n_locally_owned_dofs());
1758 * for (unsigned int i = 0; i < first_n_processes; ++i)
1759 * pcout << ' ' << n_dofs_per_subdomain[i];
1760 * if (output_cropped)
1761 * pcout << " ...";
1762 * pcout << std::endl;
1763 * }
1764 *
1765 * {
1766 * std::vector<types::global_dof_index> n_constraints_per_subdomain =
1767 * Utilities::MPI::gather(mpi_communicator, constraints.n_constraints());
1768 *
1769 * pcout << " Number of constraints: "
1770 * << std::accumulate(n_constraints_per_subdomain.begin(),
1771 * n_constraints_per_subdomain.end(),
1772 * 0)
1773 * << std::endl
1774 * << " by partition: ";
1775 * for (unsigned int i = 0; i < first_n_processes; ++i)
1776 * pcout << ' ' << n_constraints_per_subdomain[i];
1777 * if (output_cropped)
1778 * pcout << " ...";
1779 * pcout << std::endl;
1780 * }
1781 *
1782 * {
1783 * std::vector<unsigned int> n_fe_indices(fe_collection.size(), 0);
1784 * for (const auto &cell : dof_handler.active_cell_iterators())
1785 * if (cell->is_locally_owned())
1786 * n_fe_indices[cell->active_fe_index()]++;
1787 *
1788 * Utilities::MPI::sum(n_fe_indices, mpi_communicator, n_fe_indices);
1789 *
1790 * pcout << " Frequencies of poly. degrees:";
1791 * for (unsigned int i = 0; i < fe_collection.size(); ++i)
1792 * if (n_fe_indices[i] > 0)
1793 * pcout << ' ' << fe_collection[i].degree << ':' << n_fe_indices[i];
1794 * pcout << std::endl;
1795 * }
1796 * }
1797 *
1798 *
1799 *
1800 * @endcode
1801 *
1802 *
1803 * <a name="LaplaceProblemsolve_system"></a>
1804 * <h4>LaplaceProblem::solve_system</h4>
1805 *
1806
1807 *
1808 * The scaffold around the solution is similar to the one of @ref step_40 "step-40". We
1809 * prepare a vector that matches the requirements of MatrixFree and collect
1810 * the locally-relevant degrees of freedoms we solved the equation system. The
1811 * solution happens with the function introduced earlier.
1812 *
1813 * @code
1814 * template <int dim>
1815 * void LaplaceProblem<dim>::solve_system()
1816 * {
1817 * TimerOutput::Scope t(computing_timer, "solve system");
1818 *
1819 * LinearAlgebra::distributed::Vector<double> completely_distributed_solution;
1820 * laplace_operator.initialize_dof_vector(completely_distributed_solution);
1821 *
1822 * SolverControl solver_control(system_rhs.size(),
1823 * prm.tolerance_factor * system_rhs.l2_norm());
1824 *
1825 * solve_with_gmg(solver_control,
1826 * laplace_operator,
1827 * completely_distributed_solution,
1828 * system_rhs,
1829 * prm.mg_data,
1830 * mapping_collection,
1831 * dof_handler,
1832 * quadrature_collection);
1833 *
1834 * pcout << " Solved in " << solver_control.last_step() << " iterations."
1835 * << std::endl;
1836 *
1837 * constraints.distribute(completely_distributed_solution);
1838 *
1839 * locally_relevant_solution.copy_locally_owned_data_from(
1840 * completely_distributed_solution);
1841 * locally_relevant_solution.update_ghost_values();
1842 * }
1843 *
1844 *
1845 *
1846 * @endcode
1847 *
1848 *
1849 * <a name="LaplaceProblemcompute_indicators"></a>
1850 * <h4>LaplaceProblem::compute_indicators</h4>
1851 *
1852
1853 *
1854 * This function contains only a part of the typical <code>refine_grid</code>
1855 * function from other tutorials and is new in that sense. Here, we will only
1856 * calculate all indicators for adaptation with actually refining the grid. We
1857 * do this for the purpose of writing all indicators to the file system, so we
1858 * store them for later.
1859 *
1860
1861 *
1862 * Since we are dealing the an elliptic problem, we will make use of the
1863 * KellyErrorEstimator again, but with a slight difference. Modifying the
1864 * scaling factor of the underlying face integrals to be dependent on the
1865 * actual polynomial degree of the neighboring elements is favorable in
1866 * hp-adaptive applications @cite davydov2017hp. We can do this by specifying
1867 * the very last parameter from the additional ones you notices. The others
1868 * are actually just the defaults.
1869 *
1870
1871 *
1872 * For the purpose of hp-adaptation, we will calculate smoothness estimates
1873 * with the strategy presented in the tutorial introduction and use the
1874 * implementation in SmoothnessEstimator::Legendre. In the Parameters struct,
1875 * we set the minimal polynomial degree to 2 as it seems that the smoothness
1876 * estimation algorithms have trouble with linear elements.
1877 *
1878 * @code
1879 * template <int dim>
1880 * void LaplaceProblem<dim>::compute_indicators()
1881 * {
1882 * TimerOutput::Scope t(computing_timer, "compute indicators");
1883 *
1884 * estimated_error_per_cell.grow_or_shrink(triangulation.n_active_cells());
1886 * dof_handler,
1887 * face_quadrature_collection,
1888 * std::map<types::boundary_id, const Function<dim> *>(),
1889 * locally_relevant_solution,
1890 * estimated_error_per_cell,
1891 * /*component_mask=*/ComponentMask(),
1892 * /*coefficients=*/nullptr,
1893 * /*n_threads=*/numbers::invalid_unsigned_int,
1894 * /*subdomain_id=*/numbers::invalid_subdomain_id,
1895 * /*material_id=*/numbers::invalid_material_id,
1896 * /*strategy=*/
1898 *
1899 * hp_decision_indicators.grow_or_shrink(triangulation.n_active_cells());
1901 * dof_handler,
1902 * locally_relevant_solution,
1903 * hp_decision_indicators);
1904 * }
1905 *
1906 *
1907 *
1908 * @endcode
1909 *
1910 *
1911 * <a name="LaplaceProblemadapt_resolution"></a>
1912 * <h4>LaplaceProblem::adapt_resolution</h4>
1913 *
1914
1915 *
1916 * With the previously calculated indicators, we will finally flag all cells
1917 * for adaptation and also execute refinement in this function. As in previous
1918 * tutorials, we will use the "fixed number" strategy, but now for
1919 * hp-adaptation.
1920 *
1921 * @code
1922 * template <int dim>
1923 * void LaplaceProblem<dim>::adapt_resolution()
1924 * {
1925 * TimerOutput::Scope t(computing_timer, "adapt resolution");
1926 *
1927 * @endcode
1928 *
1929 * First, we will set refine and coarsen flags based on the error estimates
1930 * on each cell. There is nothing new here.
1931 *
1932
1933 *
1934 * We will use general refine and coarsen fractions that have been
1935 * elaborated in the other deal.II tutorials: using the fixed number
1936 * strategy, we will flag 30% of all cells for refinement and 3% for
1937 * coarsening, as provided in the Parameters struct.
1938 *
1939 * @code
1941 * triangulation,
1942 * estimated_error_per_cell,
1943 * prm.refine_fraction,
1944 * prm.coarsen_fraction);
1945 *
1946 * @endcode
1947 *
1948 * Next, we will make all adjustments for hp-adaptation. We want to refine
1949 * and coarsen those cells flagged in the previous step, but need to decide
1950 * if we would like to do it by adjusting the grid resolution or the
1951 * polynomial degree.
1952 *
1953
1954 *
1955 * The next function call sets future FE indices according to the previously
1956 * calculated smoothness indicators as p-adaptation indicators. These
1957 * indices will only be set on those cells that have refine or coarsen flags
1958 * assigned.
1959 *
1960
1961 *
1962 * For the p-adaptation fractions, we will take an educated guess. Since we
1963 * only expect a single singularity in our scenario, i.e., in the origin of
1964 * the domain, and a smooth solution anywhere else, we would like to
1965 * strongly prefer to use p-adaptation over h-adaptation. This reflects in
1966 * our choice of a fraction of 90% for both p-refinement and p-coarsening.
1967 *
1968 * @code
1970 * hp_decision_indicators,
1971 * prm.p_refine_fraction,
1972 * prm.p_coarsen_fraction);
1973 *
1974 * @endcode
1975 *
1976 * At this stage, we have both the future FE indices and the classic refine
1977 * and coarsen flags set, from which the latter will be interpreted by
1979 * We would like to only impose one type of adaptation on cells, which is
1980 * what the next function will sort out for us. In short, on cells which
1981 * have both types of indicators assigned, we will favor the p-adaptation
1982 * one and remove the h-adaptation one.
1983 *
1984 * @code
1985 * hp::Refinement::choose_p_over_h(dof_handler);
1986 *
1987 * @endcode
1988 *
1989 * After setting all indicators, we will remove those that exceed the
1990 * specified limits of the provided level ranges in the Parameters struct.
1991 * This limitation naturally arises for p-adaptation as the number of
1992 * supplied finite elements is limited. In addition, we registered a custom
1993 * hierarchy for p-adaptation in the constructor. Now, we need to do this
1994 * manually in the h-adaptive context like in @ref step_31 "step-31".
1995 *
1996
1997 *
1998 * We will iterate over all cells on the designated min and max levels and
1999 * remove the corresponding flags. As an alternative, we could also flag
2000 * these cells for p-adaptation by setting future FE indices accordingly
2001 * instead of simply clearing the refine and coarsen flags.
2002 *
2003 * @code
2004 * Assert(triangulation.n_levels() >= prm.min_h_level + 1 &&
2005 * triangulation.n_levels() <= prm.max_h_level + 1,
2006 * ExcInternalError());
2007 *
2008 * if (triangulation.n_levels() > prm.max_h_level)
2009 * for (const auto &cell :
2010 * triangulation.active_cell_iterators_on_level(prm.max_h_level))
2011 * cell->clear_refine_flag();
2012 *
2013 * for (const auto &cell :
2014 * triangulation.active_cell_iterators_on_level(prm.min_h_level))
2015 * cell->clear_coarsen_flag();
2016 *
2017 * @endcode
2018 *
2019 * In the end, we are left to execute coarsening and refinement. Here, not
2020 * only the grid will be updated, but also all previous future FE indices
2021 * will become active.
2022 *
2023
2024 *
2025 * Remember that we have attached functions to triangulation signals in the
2026 * constructor, will be triggered in this function call. So there is even
2027 * more happening: weighted repartitioning will be performed to ensure load
2028 * balancing, as well as we will limit the difference of p-levels between
2029 * neighboring cells.
2030 *
2031 * @code
2032 * triangulation.execute_coarsening_and_refinement();
2033 * }
2034 *
2035 *
2036 *
2037 * @endcode
2038 *
2039 *
2040 * <a name="LaplaceProblemoutput_results"></a>
2041 * <h4>LaplaceProblem::output_results</h4>
2042 *
2043
2044 *
2045 * Writing results to the file system in parallel applications works exactly
2046 * like in @ref step_40 "step-40". In addition to the data containers that we prepared
2047 * throughout the tutorial, we would also like to write out the polynomial
2048 * degree of each finite element on the grid as well as the subdomain each
2049 * cell belongs to. We prepare necessary containers for this in the scope of
2050 * this function.
2051 *
2052 * @code
2053 * template <int dim>
2054 * void LaplaceProblem<dim>::output_results(const unsigned int cycle)
2055 * {
2056 * TimerOutput::Scope t(computing_timer, "output results");
2057 *
2058 * Vector<float> fe_degrees(triangulation.n_active_cells());
2059 * for (const auto &cell : dof_handler.active_cell_iterators())
2060 * if (cell->is_locally_owned())
2061 * fe_degrees(cell->active_cell_index()) = cell->get_fe().degree;
2062 *
2063 * Vector<float> subdomain(triangulation.n_active_cells());
2064 * for (auto &subd : subdomain)
2065 * subd = triangulation.locally_owned_subdomain();
2066 *
2067 * DataOut<dim> data_out;
2068 * data_out.attach_dof_handler(dof_handler);
2069 * data_out.add_data_vector(locally_relevant_solution, "solution");
2070 * data_out.add_data_vector(fe_degrees, "fe_degree");
2071 * data_out.add_data_vector(subdomain, "subdomain");
2072 * data_out.add_data_vector(estimated_error_per_cell, "error");
2073 * data_out.add_data_vector(hp_decision_indicators, "hp_indicator");
2074 * data_out.build_patches(mapping_collection);
2075 *
2076 * data_out.write_vtu_with_pvtu_record(
2077 * "./", "solution", cycle, mpi_communicator, 2, 1);
2078 * }
2079 *
2080 *
2081 *
2082 * @endcode
2083 *
2084 *
2085 * <a name="LaplaceProblemrun"></a>
2086 * <h4>LaplaceProblem::run</h4>
2087 *
2088
2089 *
2090 * The actual run function again looks very familiar to @ref step_40 "step-40". The only
2091 * addition is the bracketed section that precedes the actual cycle loop.
2092 * Here, we will pre-calculate the Legendre transformation matrices. In
2093 * general, these will be calculated on the fly via lazy allocation whenever a
2094 * certain matrix is needed. For timing purposes however, we would like to
2095 * calculate them all at once before the actual time measurement begins. We
2096 * will thus designate their calculation to their own scope.
2097 *
2098 * @code
2099 * template <int dim>
2100 * void LaplaceProblem<dim>::run()
2101 * {
2102 * pcout << "Running with Trilinos on "
2103 * << Utilities::MPI::n_mpi_processes(mpi_communicator)
2104 * << " MPI rank(s)..." << std::endl;
2105 *
2106 * {
2107 * pcout << "Calculating transformation matrices..." << std::endl;
2108 * TimerOutput::Scope t(computing_timer, "calculate transformation");
2109 * legendre->precalculate_all_transformation_matrices();
2110 * }
2111 *
2112 * for (unsigned int cycle = 0; cycle < prm.n_cycles; ++cycle)
2113 * {
2114 * pcout << "Cycle " << cycle << ':' << std::endl;
2115 *
2116 * if (cycle == 0)
2117 * initialize_grid();
2118 * else
2119 * adapt_resolution();
2120 *
2121 * setup_system();
2122 *
2123 * print_diagnostics();
2124 *
2125 * solve_system();
2126 *
2127 * compute_indicators();
2128 *
2129 * if (Utilities::MPI::n_mpi_processes(mpi_communicator) <= 32)
2130 * output_results(cycle);
2131 *
2132 * computing_timer.print_summary();
2133 * computing_timer.reset();
2134 *
2135 * pcout << std::endl;
2136 * }
2137 * }
2138 * } // namespace Step75
2139 *
2140 *
2141 *
2142 * @endcode
2143 *
2144 *
2145 * <a name="main"></a>
2146 * <h4>main()</h4>
2147 *
2148
2149 *
2150 * The final function is the <code>main</code> function that will ultimately
2151 * create and run a LaplaceOperator instantiation. Its structure is similar to
2152 * most other tutorial programs.
2153 *
2154 * @code
2155 * int main(int argc, char *argv[])
2156 * {
2157 * try
2158 * {
2159 * using namespace dealii;
2160 * using namespace Step75;
2161 *
2162 * Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
2163 *
2164 * Parameters prm;
2165 * LaplaceProblem<2> laplace_problem(prm);
2166 * laplace_problem.run();
2167 * }
2168 * catch (std::exception &exc)
2169 * {
2170 * std::cerr << std::endl
2171 * << std::endl
2172 * << "----------------------------------------------------"
2173 * << std::endl;
2174 * std::cerr << "Exception on processing: " << std::endl
2175 * << exc.what() << std::endl
2176 * << "Aborting!" << std::endl
2177 * << "----------------------------------------------------"
2178 * << std::endl;
2179 *
2180 * return 1;
2181 * }
2182 * catch (...)
2183 * {
2184 * std::cerr << std::endl
2185 * << std::endl
2186 * << "----------------------------------------------------"
2187 * << std::endl;
2188 * std::cerr << "Unknown exception!" << std::endl
2189 * << "Aborting!" << std::endl
2190 * << "----------------------------------------------------"
2191 * << std::endl;
2192 * return 1;
2193 * }
2194 *
2195 * return 0;
2196 * }
2197 * @endcode
2198<a name="Results"></a><h1>Results</h1>
2199
2200
2201When you run the program with the given parameters on four processes in
2202release mode, your terminal output should look like this:
2203@code
2204Running with Trilinos on 4 MPI rank(s)...
2205Calculating transformation matrices...
2206Cycle 0:
2207 Number of active cells: 3072
2208 by partition: 768 768 768 768
2209 Number of degrees of freedom: 12545
2210 by partition: 3201 3104 3136 3104
2211 Number of constraints: 542
2212 by partition: 165 74 138 165
2213 Frequencies of poly. degrees: 2:3072
2214 Solved in 7 iterations.
2215
2216
2217+---------------------------------------------+------------+------------+
2218| Total wallclock time elapsed since start | 0.598s | |
2219| | | |
2220| Section | no. calls | wall time | % of total |
2221+---------------------------------+-----------+------------+------------+
2222| calculate transformation | 1 | 0.0533s | 8.9% |
2223| compute indicators | 1 | 0.0177s | 3% |
2224| initialize grid | 1 | 0.0397s | 6.6% |
2225| output results | 1 | 0.0844s | 14% |
2226| setup system | 1 | 0.0351s | 5.9% |
2227| solve system | 1 | 0.362s | 61% |
2228+---------------------------------+-----------+------------+------------+
2229
2230
2231Cycle 1:
2232 Number of active cells: 3351
2233 by partition: 875 761 843 872
2234 Number of degrees of freedom: 18223
2235 by partition: 4535 4735 4543 4410
2236 Number of constraints: 1202
2237 by partition: 303 290 326 283
2238 Frequencies of poly. degrees: 2:2523 3:828
2239 Solved in 7 iterations.
2240
2241
2242+---------------------------------------------+------------+------------+
2243| Total wallclock time elapsed since start | 0.442s | |
2244| | | |
2245| Section | no. calls | wall time | % of total |
2246+---------------------------------+-----------+------------+------------+
2247| adapt resolution | 1 | 0.0189s | 4.3% |
2248| compute indicators | 1 | 0.0135s | 3% |
2249| output results | 1 | 0.064s | 14% |
2250| setup system | 1 | 0.0232s | 5.2% |
2251| solve system | 1 | 0.322s | 73% |
2252+---------------------------------+-----------+------------+------------+
2253
2254
2255...
2256
2257
2258Cycle 7:
2259 Number of active cells: 5610
2260 by partition: 1324 1483 1482 1321
2261 Number of degrees of freedom: 82062
2262 by partition: 21116 19951 20113 20882
2263 Number of constraints: 14383
2264 by partition: 3825 3225 3557 3776
2265 Frequencies of poly. degrees: 2:1130 3:1283 4:2727 5:465 6:5
2266 Solved in 7 iterations.
2267
2268
2269+---------------------------------------------+------------+------------+
2270| Total wallclock time elapsed since start | 0.932s | |
2271| | | |
2272| Section | no. calls | wall time | % of total |
2273+---------------------------------+-----------+------------+------------+
2274| adapt resolution | 1 | 0.0182s | 1.9% |
2275| compute indicators | 1 | 0.0173s | 1.9% |
2276| output results | 1 | 0.0572s | 6.1% |
2277| setup system | 1 | 0.0252s | 2.7% |
2278| solve system | 1 | 0.813s | 87% |
2279+---------------------------------+-----------+------------+------------+
2280@endcode
2281
2282When running the code with more processes, you will notice slight
2283differences in the number of active cells and degrees of freedom. This
2284is due to the fact that solver and preconditioner depend on the
2285partitioning of the problem, which might yield to slight differences of
2286the solution in the last digits and ultimately yields to different
2287adaptation behavior.
2288
2289Furthermore, the number of iterations for the solver stays about the
2290same in all cycles despite hp-adaptation, indicating the robustness of
2291the proposed algorithms and promising good scalability for even larger
2292problem sizes and on more processes.
2293
2294Let us have a look at the graphical output of the program. After all
2295refinement cycles in the given parameter configuration, the actual
2296discretized function space looks like the following with its
2297partitioning on twelve processes on the left and the polynomial degrees
2298of finite elements on the right. In the left picture, each color
2299represents a unique subdomain. In the right picture, the lightest color
2300corresponds to the polynomial degree two and the darkest one corresponds
2301to degree six:
2302
2303<div class="twocolumn" style="width: 80%; text-align: center;">
2304 <div>
2305 <img src="https://www.dealii.org/images/steps/developer/step-75.subdomains-07.svg"
2306 alt="Partitioning after seven refinements.">
2307 </div>
2308 <div>
2309 <img src="https://www.dealii.org/images/steps/developer/step-75.fedegrees-07.svg"
2310 alt="Local approximation degrees after seven refinements.">
2311 </div>
2312</div>
2313
2314
2315
2316<a name="extensions"></a>
2317<a name="Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
2318
2319
2320This tutorial shows only one particular way how to use parallel
2321hp-adaptive finite element methods. In the following paragraphs, you
2322will get to know which alternatives are possible. Most of these
2323extensions are already part of https://github.com/marcfehling/hpbox/,
2324which provides you with implementation examples that you can play
2325around with.
2326
2327
2328<a name="Differenthpdecisionstrategies"></a><h4>Different hp-decision strategies</h4>
2329
2330
2331The deal.II library offers multiple strategies to decide which type of
2332adaptation to impose on cells: either adjust the grid resolution or
2333change the polynomial degree. We only presented the <i>Legendre
2334coefficient decay</i> strategy in this tutorial, while @ref step_27 "step-27"
2335demonstrated the <i>Fourier</i> equivalent of the same idea.
2336
2337See the "possibilities for extensions" section of @ref step_27 "step-27" for an
2338overview over these strategies, or the corresponding documentation
2339for a detailed description.
2340
2341There, another strategy is mentioned that has not been shown in any
2342tutorial so far: the strategy based on <i>refinement history</i>. The
2343usage of this method for parallel distributed applications is more
2344tricky than the others, so we will highlight the challenges that come
2345along with it. We need information about the final state of refinement
2346flags, and we need to transfer the solution across refined meshes. For
2347the former, we need to attach the hp::Refinement::predict_error()
2348function to the Triangulation::Signals::post_p4est_refinement signal in
2349a way that it will be called <i>after</i> the
2350hp::Refinement::limit_p_level_difference() function. At this stage, all
2351refinement flags and future FE indices are terminally set and a reliable
2352prediction of the error is possible. The predicted error then needs to
2353be transferred across refined meshes with the aid of
2354parallel::distributed::CellDataTransfer.
2355
2356Try implementing one of these strategies into this tutorial and observe
2357the subtle changes to the results. You will notice that all strategies
2358are capable of identifying the singularities near the reentrant corners
2359and will perform @f$h@f$-refinement in these regions, while preferring
2360@f$p@f$-refinement in the bulk domain. A detailed comparison of these
2361strategies is presented in @cite fehling2020 .
2362
2363
2364<a name="Solvewithmatrixbasedmethods"></a><h4>Solve with matrix-based methods</h4>
2365
2366
2367This tutorial focuses solely on matrix-free strategies. All hp-adaptive
2368algorithms however also work with matrix-based approaches in the
2369parallel distributed context.
2370
2371To create a system matrix, you can either use the
2372LaplaceOperator::get_system_matrix() function, or use an
2373<code>assemble_system()</code> function similar to the one of @ref step_27 "step-27".
2374You can then pass the system matrix to the solver as usual.
2375
2376You can time the results of both matrix-based and matrix-free
2377implementations, quantify the speed-up, and convince yourself which
2378variant is faster.
2379
2380
2381<a name="Multigridvariants"></a><h4>Multigrid variants</h4>
2382
2383
2384For sake of simplicity, we have restricted ourselves to a single type of
2385coarse-grid solver (CG with AMG), smoother (Chebyshev smoother with
2386point Jacobi preconditioner), and geometric-coarsening scheme (global
2387coarsening) within the multigrid algorithm. Feel free to try out
2388alternatives and investigate their performance and robustness.
2389 *
2390 *
2391<a name="PlainProg"></a>
2392<h1> The plain program</h1>
2393@include "step-75.cc"
2394*/
void reinit(const IndexSet &local_constraints=IndexSet())
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)
Definition: data_out.cc:1064
void reinit(const Triangulation< dim, spacedim > &tria)
void evaluate(const EvaluationFlags::EvaluationFlags evaluation_flag)
Definition: fe_q.h:549
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, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), 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 resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&...args)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename PreconditionerType::AdditionalData &additional_data=typename PreconditionerType::AdditionalData())
const DoFHandler< dim > & get_dof_handler(const unsigned int dof_handler_index=0) const
void initialize_dof_vector(VectorType &vec, const unsigned int dof_handler_index=0) const
void cell_loop(const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false) const
void reinit(const MappingType &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< number2 > &constraint, const QuadratureType &quad, const AdditionalData &additional_data=AdditionalData())
Definition: point.h:111
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
@ wall_times
Definition: timer.h:653
void coarsen_global(const unsigned int times=1)
virtual void execute_coarsening_and_refinement()
Definition: vector.h:109
unsigned int size() const
Definition: collection.h:264
static WeightingFunction ndofs_weighting(const std::pair< float, float > &coefficients)
Definition: cell_weights.cc:77
@ update_gradients
Shape function gradients.
Point< 2 > second
Definition: grid_out.cc:4604
Point< 2 > first
Definition: grid_out.cc:4603
unsigned int level
Definition: grid_out.cc:4606
__global__ void set(Number *val, const Number s, const size_type N)
std::string write_vtu_with_pvtu_record(const std::string &directory, const std::string &filename_without_extension, const unsigned int counter, const MPI_Comm &mpi_communicator, const unsigned int n_digits_for_counter=numbers::invalid_unsigned_int, const unsigned int n_groups=0) const
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:439
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1144
std::array< double, dim > to_spherical(const Point< dim > &point)
void hyper_L(Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false)
void subdivided_hyper_L(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &repetitions, const Point< dim > &bottom_left, const Point< dim > &top_right, const std::vector< int > &n_cells_to_remove)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2084
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ diagonal
Matrix is diagonal.
@ general
No special properties.
std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > create_geometric_coarsening_sequence(const Triangulation< dim, spacedim > &tria)
unsigned int create_next_polynomial_coarsening_degree(const unsigned int degree, const PolynomialCoarseningSequenceType &p_sequence)
std::vector< unsigned int > create_polynomial_coarsening_sequence(const unsigned int max_degree, const PolynomialCoarseningSequenceType &p_sequence)
void compute_diagonal(const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, VectorType &diagonal_global, const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &local_vmult, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
void compute_matrix(const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const AffineConstraints< Number > &constraints, MatrixType &matrix, const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &local_vmult, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
FESeries::Legendre< dim, spacedim > default_fe_series(const hp::FECollection< dim, spacedim > &fe_collection, const unsigned int component=numbers::invalid_unsigned_int)
void coefficient_decay(FESeries::Legendre< dim, spacedim > &fe_legendre, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &smoothness_indicators, const VectorTools::NormType regression_strategy=VectorTools::Linfty_norm, const double smallest_abs_coefficient=1e-10, const bool only_flagged_cells=false)
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
void call(const std::function< RT()> &function, internal::return_value< RT > &ret_val)
void free(T *&pointer)
Definition: cuda.h:97
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:151
std::vector< T > gather(const MPI_Comm &comm, const T &object_to_send, const unsigned int root_process=0)
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:140
T max(const T &t, const MPI_Comm &mpi_communicator)
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=ComponentMask())
void run(const Iterator &begin, const typename identity< Iterator >::type &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
Definition: work_stream.h:474
bool limit_p_level_difference(const ::DoFHandler< dim, spacedim > &dof_handler, const unsigned int max_difference=1, const unsigned int contains_fe_index=0)
Definition: refinement.cc:803
void p_adaptivity_fixed_number(const ::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const double p_refine_fraction=0.5, const double p_coarsen_fraction=0.5, const ComparisonFunction< typename identity< Number >::type > &compare_refine=std::greater_equal< Number >(), const ComparisonFunction< typename identity< Number >::type > &compare_coarsen=std::less_equal< Number >())
Definition: refinement.cc:248
void predict_error(const ::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &error_indicators, Vector< Number > &predicted_errors, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
Definition: refinement.cc:550
Definition: hp.h:118
int(&) functions(const void *v1, const void *v2)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
Definition: mg.h:82
const types::material_id invalid_material_id
Definition: types.h:233
const types::subdomain_id invalid_subdomain_id
Definition: types.h:281
static const unsigned int invalid_unsigned_int
Definition: types.h:201
void refine_and_coarsen_fixed_number(parallel::distributed::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())
double legendre(unsigned int l, double x)
Definition: cmath.h:75
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
Definition: types.h:32
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
UpdateFlags mapping_update_flags
Definition: matrix_free.h:367
boost::signals2::signal< void()> post_p4est_refinement
Definition: tria.h:2388