deal.II version GIT relicensing-2173-gae8fc9d14b 2024-11-24 06:40:00+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
step-50.h
Go to the documentation of this file.
1,
482 *   const unsigned int /*component*/ = 0) const override
483 *   {
484 *   return 1.0;
485 *   }
486 *  
487 *  
488 *   template <typename number>
490 *   value(const Point<dim, VectorizedArray<number>> & /*p*/,
491 *   const unsigned int /*component*/ = 0) const
492 *   {
493 *   return VectorizedArray<number>(1.0);
494 *   }
495 *   };
496 *  
497 *  
498 * @endcode
499 *
500 * This next class represents the diffusion coefficient. We use a variable
501 * coefficient which is 100.0 at any point where at least one coordinate is
502 * less than -0.5, and 1.0 at all other points. As above, a separate value()
503 * returning a VectorizedArray is used for the matrix-free code. An @p
504 * average() function computes the arithmetic average for a set of points.
505 *
506 * @code
507 *   template <int dim>
508 *   class Coefficient : public Function<dim>
509 *   {
510 *   public:
511 *   virtual double value(const Point<dim> &p,
512 *   const unsigned int /*component*/ = 0) const override;
513 *  
514 *   template <typename number>
516 *   const unsigned int /*component*/ = 0) const;
517 *  
518 *   template <typename number>
519 *   number average_value(const std::vector<Point<dim, number>> &points) const;
520 *  
521 * @endcode
522 *
523 * When using a coefficient in the MatrixFree framework, we also
524 * need a function that creates a Table of coefficient values for a
525 * set of cells provided by the MatrixFree operator argument here.
526 *
527 * @code
528 *   template <typename number>
529 *   std::shared_ptr<Table<2, VectorizedArray<number>>> make_coefficient_table(
530 *   const MatrixFree<dim, number, VectorizedArray<number>> &mf_storage) const;
531 *   };
532 *  
533 *  
534 *  
535 *   template <int dim>
536 *   double Coefficient<dim>::value(const Point<dim> &p, const unsigned int) const
537 *   {
538 *   for (int d = 0; d < dim; ++d)
539 *   {
540 *   if (p[d] < -0.5)
541 *   return 100.0;
542 *   }
543 *   return 1.0;
544 *   }
545 *  
546 *  
547 *  
548 *   template <int dim>
549 *   template <typename number>
551 *   Coefficient<dim>::value(const Point<dim, VectorizedArray<number>> &p,
552 *   const unsigned int) const
553 *   {
554 *   VectorizedArray<number> return_value = VectorizedArray<number>(1.0);
555 *   for (unsigned int i = 0; i < VectorizedArray<number>::size(); ++i)
556 *   {
557 *   for (int d = 0; d < dim; ++d)
558 *   if (p[d][i] < -0.5)
559 *   {
560 *   return_value[i] = 100.0;
561 *   break;
562 *   }
563 *   }
564 *  
565 *   return return_value;
566 *   }
567 *  
568 *  
569 *  
570 *   template <int dim>
571 *   template <typename number>
572 *   number Coefficient<dim>::average_value(
573 *   const std::vector<Point<dim, number>> &points) const
574 *   {
575 *   number average(0);
576 *   for (unsigned int i = 0; i < points.size(); ++i)
577 *   average += value(points[i]);
578 *   average /= points.size();
579 *  
580 *   return average;
581 *   }
582 *  
583 *  
584 *  
585 *   template <int dim>
586 *   template <typename number>
587 *   std::shared_ptr<Table<2, VectorizedArray<number>>>
588 *   Coefficient<dim>::make_coefficient_table(
589 *   const MatrixFree<dim, number, VectorizedArray<number>> &mf_storage) const
590 *   {
591 *   auto coefficient_table =
592 *   std::make_shared<Table<2, VectorizedArray<number>>>();
593 *  
594 *   FEEvaluation<dim, -1, 0, 1, number> fe_eval(mf_storage);
595 *  
596 *   const unsigned int n_cells = mf_storage.n_cell_batches();
597 *  
598 *   coefficient_table->reinit(n_cells, 1);
599 *  
600 *   for (unsigned int cell = 0; cell < n_cells; ++cell)
601 *   {
602 *   fe_eval.reinit(cell);
603 *  
604 *   VectorizedArray<number> average_value = 0.;
605 *   for (const unsigned int q : fe_eval.quadrature_point_indices())
606 *   average_value += value(fe_eval.quadrature_point(q));
607 *   average_value /= fe_eval.n_q_points;
608 *  
609 *   (*coefficient_table)(cell, 0) = average_value;
610 *   }
611 *  
612 *   return coefficient_table;
613 *   }
614 *  
615 *  
616 *  
617 * @endcode
618 *
619 *
620 * <a name="step_50-Runtimeparameters"></a>
621 * <h3>Run time parameters</h3>
622 *
623
624 *
625 * We will use ParameterHandler to pass in parameters at runtime. The
626 * structure @p Settings parses and stores these parameters to be queried
627 * throughout the program.
628 *
629 * @code
630 *   struct Settings
631 *   {
632 *   bool try_parse(const std::string &prm_filename);
633 *  
634 *   enum SolverType
635 *   {
636 *   gmg_mb,
637 *   gmg_mf,
638 *   amg
639 *   };
640 *  
641 *   SolverType solver;
642 *  
643 *   int dimension;
644 *   double smoother_dampen;
645 *   unsigned int smoother_steps;
646 *   unsigned int n_steps;
647 *   bool output;
648 *   };
649 *  
650 *  
651 *  
652 *   bool Settings::try_parse(const std::string &prm_filename)
653 *   {
654 *   ParameterHandler prm;
655 *   prm.declare_entry("dim",
656 *   "2",
657 *   Patterns::Integer(),
658 *   "The problem dimension.");
659 *   prm.declare_entry("n_steps",
660 *   "10",
661 *   Patterns::Integer(0),
662 *   "Number of adaptive refinement steps.");
663 *   prm.declare_entry("smoother dampen",
664 *   "1.0",
665 *   Patterns::Double(0.0),
666 *   "Dampen factor for the smoother.");
667 *   prm.declare_entry("smoother steps",
668 *   "1",
669 *   Patterns::Integer(1),
670 *   "Number of smoother steps.");
671 *   prm.declare_entry("solver",
672 *   "MF",
673 *   Patterns::Selection("MF|MB|AMG"),
674 *   "Switch between matrix-free GMG, "
675 *   "matrix-based GMG, and AMG.");
676 *   prm.declare_entry("output",
677 *   "false",
678 *   Patterns::Bool(),
679 *   "Output graphical results.");
680 *  
681 *   if (prm_filename.empty())
682 *   {
683 *   std::cout
684 *   << "**** Error: No input file provided!\n"
685 *   << "**** Error: Call this program as './step-50 input.prm\n"
686 *   << '\n'
687 *   << "**** You may want to use one of the input files in this\n"
688 *   << "**** directory, or use the following default values\n"
689 *   << "**** to create an input file:\n";
690 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
691 *   prm.print_parameters(std::cout, ParameterHandler::PRM);
692 *   return false;
693 *   }
694 *  
695 *   try
696 *   {
697 *   prm.parse_input(prm_filename);
698 *   }
699 *   catch (std::exception &e)
700 *   {
701 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
702 *   std::cerr << e.what() << std::endl;
703 *   return false;
704 *   }
705 *  
706 *   if (prm.get("solver") == "MF")
707 *   this->solver = gmg_mf;
708 *   else if (prm.get("solver") == "MB")
709 *   this->solver = gmg_mb;
710 *   else if (prm.get("solver") == "AMG")
711 *   this->solver = amg;
712 *   else
713 *   AssertThrow(false, ExcNotImplemented());
714 *  
715 *   this->dimension = prm.get_integer("dim");
716 *   this->n_steps = prm.get_integer("n_steps");
717 *   this->smoother_dampen = prm.get_double("smoother dampen");
718 *   this->smoother_steps = prm.get_integer("smoother steps");
719 *   this->output = prm.get_bool("output");
720 *  
721 *   return true;
722 *   }
723 *  
724 *  
725 *  
726 * @endcode
727 *
728 *
729 * <a name="step_50-LaplaceProblemclass"></a>
730 * <h3>LaplaceProblem class</h3>
731 *
732
733 *
734 * This is the main class of the program. It looks very similar to
735 * @ref step_16 "step-16", @ref step_37 "step-37", and @ref step_40 "step-40". For the MatrixFree setup, we use the
736 * MatrixFreeOperators::LaplaceOperator class which defines `local_apply()`,
737 * `compute_diagonal()`, and `set_coefficient()` functions internally. Note
738 * that the polynomial degree is a template parameter of this class. This is
739 * necessary for the matrix-free code.
740 *
741 * @code
742 *   template <int dim, int degree>
743 *   class LaplaceProblem
744 *   {
745 *   public:
746 *   LaplaceProblem(const Settings &settings);
747 *   void run();
748 *  
749 *   private:
750 * @endcode
751 *
752 * We will use the following types throughout the program. First the
753 * matrix-based types, after that the matrix-free classes. For the
754 * matrix-free implementation, we use @p float for the level operators.
755 *
756 * @code
757 *   using MatrixType = LA::MPI::SparseMatrix;
758 *   using VectorType = LA::MPI::Vector;
759 *   using PreconditionAMG = LA::MPI::PreconditionAMG;
760 *  
761 *   using MatrixFreeLevelMatrix = MatrixFreeOperators::LaplaceOperator<
762 *   dim,
763 *   degree,
764 *   degree + 1,
765 *   1,
767 *   using MatrixFreeActiveMatrix = MatrixFreeOperators::LaplaceOperator<
768 *   dim,
769 *   degree,
770 *   degree + 1,
771 *   1,
773 *  
774 *   using MatrixFreeLevelVector = LinearAlgebra::distributed::Vector<float>;
775 *   using MatrixFreeActiveVector = LinearAlgebra::distributed::Vector<double>;
776 *  
777 *   void setup_system();
778 *   void setup_multigrid();
779 *   void assemble_system();
780 *   void assemble_multigrid();
781 *   void assemble_rhs();
782 *   void solve();
783 *   void estimate();
784 *   void refine_grid();
785 *   void output_results(const unsigned int cycle);
786 *  
787 *   Settings settings;
788 *  
789 *   MPI_Comm mpi_communicator;
790 *   ConditionalOStream pcout;
791 *  
793 *   const MappingQ1<dim> mapping;
794 *   const FE_Q<dim> fe;
795 *  
796 *   DoFHandler<dim> dof_handler;
797 *  
798 *   IndexSet locally_owned_dofs;
799 *   IndexSet locally_relevant_dofs;
800 *   AffineConstraints<double> constraints;
801 *  
802 *   MatrixType system_matrix;
803 *   MatrixFreeActiveMatrix mf_system_matrix;
804 *   VectorType solution;
805 *   VectorType right_hand_side;
806 *   Vector<double> estimated_error_square_per_cell;
807 *  
808 *   MGLevelObject<MatrixType> mg_matrix;
809 *   MGLevelObject<MatrixType> mg_interface_in;
810 *   MGConstrainedDoFs mg_constrained_dofs;
811 *  
812 *   MGLevelObject<MatrixFreeLevelMatrix> mf_mg_matrix;
813 *  
814 *   TimerOutput computing_timer;
815 *   };
816 *  
817 *  
818 * @endcode
819 *
820 * The only interesting part about the constructor is that we construct the
821 * multigrid hierarchy unless we use AMG. For that, we need to parse the
822 * run time parameters before this constructor completes.
823 *
824 * @code
825 *   template <int dim, int degree>
826 *   LaplaceProblem<dim, degree>::LaplaceProblem(const Settings &settings)
827 *   : settings(settings)
828 *   , mpi_communicator(MPI_COMM_WORLD)
829 *   , pcout(std::cout,
830 *   (Utilities::MPI::this_mpi_process(mpi_communicator) == 0))
831 *   , triangulation(
832 *   mpi_communicator,
834 *   (settings.solver == Settings::amg) ?
836 *   parallel::distributed::Triangulation<
838 *   , mapping()
839 *   , fe(degree)
840 *   , dof_handler(triangulation)
841 *   , computing_timer(pcout, TimerOutput::never, TimerOutput::wall_times)
842 *   {
843 *   GridGenerator::hyper_L(triangulation, -1., 1., /*colorize*/ false);
844 *   triangulation.refine_global(1);
845 *   }
846 *  
847 *  
848 *  
849 * @endcode
850 *
851 *
852 * <a name="step_50-LaplaceProblemsetup_system"></a>
853 * <h4>LaplaceProblem::setup_system()</h4>
854 *
855
856 *
857 * Unlike @ref step_16 "step-16" and @ref step_37 "step-37", we split the set up into two parts,
858 * setup_system() and setup_multigrid(). Here is the typical setup_system()
859 * function for the active mesh found in most tutorials. For matrix-free, the
860 * active mesh set up is similar to @ref step_37 "step-37"; for matrix-based (GMG and AMG
861 * solvers), the setup is similar to @ref step_40 "step-40".
862 *
863 * @code
864 *   template <int dim, int degree>
865 *   void LaplaceProblem<dim, degree>::setup_system()
866 *   {
867 *   TimerOutput::Scope timing(computing_timer, "Setup");
868 *  
869 *   dof_handler.distribute_dofs(fe);
870 *  
871 *   locally_relevant_dofs =
873 *   locally_owned_dofs = dof_handler.locally_owned_dofs();
874 *  
875 *   solution.reinit(locally_owned_dofs, mpi_communicator);
876 *   right_hand_side.reinit(locally_owned_dofs, mpi_communicator);
877 *   constraints.reinit(locally_owned_dofs, locally_relevant_dofs);
878 *   DoFTools::make_hanging_node_constraints(dof_handler, constraints);
879 *  
881 *   mapping, dof_handler, 0, Functions::ZeroFunction<dim>(), constraints);
882 *   constraints.close();
883 *  
884 *   switch (settings.solver)
885 *   {
886 *   case Settings::gmg_mf:
887 *   {
888 *   typename MatrixFree<dim, double>::AdditionalData additional_data;
889 *   additional_data.tasks_parallel_scheme =
891 *   additional_data.mapping_update_flags =
893 *   std::shared_ptr<MatrixFree<dim, double>> mf_storage =
894 *   std::make_shared<MatrixFree<dim, double>>();
895 *   mf_storage->reinit(mapping,
896 *   dof_handler,
897 *   constraints,
898 *   QGauss<1>(degree + 1),
899 *   additional_data);
900 *  
901 *   mf_system_matrix.initialize(mf_storage);
902 *  
903 *   const Coefficient<dim> coefficient;
904 *   mf_system_matrix.set_coefficient(
905 *   coefficient.make_coefficient_table(*mf_storage));
906 *  
907 *   break;
908 *   }
909 *  
910 *   case Settings::gmg_mb:
911 *   case Settings::amg:
912 *   {
913 *   #ifdef USE_PETSC_LA
914 *   DynamicSparsityPattern dsp(locally_relevant_dofs);
915 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints);
916 *  
918 *   locally_owned_dofs,
919 *   mpi_communicator,
920 *   locally_relevant_dofs);
921 *  
922 *   system_matrix.reinit(locally_owned_dofs,
923 *   locally_owned_dofs,
924 *   dsp,
925 *   mpi_communicator);
926 *   #else
927 *   TrilinosWrappers::SparsityPattern dsp(locally_owned_dofs,
928 *   locally_owned_dofs,
929 *   locally_relevant_dofs,
930 *   mpi_communicator);
931 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints);
932 *   dsp.compress();
933 *   system_matrix.reinit(dsp);
934 *   #endif
935 *  
936 *   break;
937 *   }
938 *  
939 *   default:
941 *   }
942 *   }
943 *  
944 * @endcode
945 *
946 *
947 * <a name="step_50-LaplaceProblemsetup_multigrid"></a>
948 * <h4>LaplaceProblem::setup_multigrid()</h4>
949 *
950
951 *
952 * This function does the multilevel setup for both matrix-free and
953 * matrix-based GMG. The matrix-free setup is similar to that of @ref step_37 "step-37", and
954 * the matrix-based is similar to @ref step_16 "step-16", except we must use appropriate
955 * distributed sparsity patterns.
956 *
957
958 *
959 * The function is not called for the AMG approach, but to err on the
960 * safe side, the main `switch` statement of this function
961 * nevertheless makes sure that the function only operates on known
962 * multigrid settings by throwing an assertion if the function were
963 * called for anything other than the two geometric multigrid methods.
964 *
965 * @code
966 *   template <int dim, int degree>
967 *   void LaplaceProblem<dim, degree>::setup_multigrid()
968 *   {
969 *   TimerOutput::Scope timing(computing_timer, "Setup multigrid");
970 *  
971 *   dof_handler.distribute_mg_dofs();
972 *  
973 *   mg_constrained_dofs.clear();
974 *   mg_constrained_dofs.initialize(dof_handler);
975 *  
976 *   const std::set<types::boundary_id> boundary_ids = {types::boundary_id(0)};
977 *   mg_constrained_dofs.make_zero_boundary_constraints(dof_handler,
978 *   boundary_ids);
979 *  
980 *   const unsigned int n_levels = triangulation.n_global_levels();
981 *  
982 *   switch (settings.solver)
983 *   {
984 *   case Settings::gmg_mf:
985 *   {
986 *   mf_mg_matrix.resize(0, n_levels - 1);
987 *  
988 *   for (unsigned int level = 0; level < n_levels; ++level)
989 *   {
990 *   AffineConstraints<double> level_constraints(
991 *   dof_handler.locally_owned_mg_dofs(level),
993 *   level));
994 *   for (const types::global_dof_index dof_index :
995 *   mg_constrained_dofs.get_boundary_indices(level))
996 *   level_constraints.constrain_dof_to_zero(dof_index);
997 *   level_constraints.close();
998 *  
999 *   typename MatrixFree<dim, float>::AdditionalData additional_data;
1000 *   additional_data.tasks_parallel_scheme =
1002 *   additional_data.mapping_update_flags =
1005 *   additional_data.mg_level = level;
1006 *   std::shared_ptr<MatrixFree<dim, float>> mf_storage_level(
1007 *   new MatrixFree<dim, float>());
1008 *   mf_storage_level->reinit(mapping,
1009 *   dof_handler,
1010 *   level_constraints,
1011 *   QGauss<1>(degree + 1),
1012 *   additional_data);
1013 *  
1014 *   mf_mg_matrix[level].initialize(mf_storage_level,
1015 *   mg_constrained_dofs,
1016 *   level);
1017 *  
1018 *   const Coefficient<dim> coefficient;
1019 *   mf_mg_matrix[level].set_coefficient(
1020 *   coefficient.make_coefficient_table(*mf_storage_level));
1021 *  
1022 *   mf_mg_matrix[level].compute_diagonal();
1023 *   }
1024 *  
1025 *   break;
1026 *   }
1027 *  
1028 *   case Settings::gmg_mb:
1029 *   {
1030 *   mg_matrix.resize(0, n_levels - 1);
1031 *   mg_matrix.clear_elements();
1032 *   mg_interface_in.resize(0, n_levels - 1);
1033 *   mg_interface_in.clear_elements();
1034 *  
1035 *   for (unsigned int level = 0; level < n_levels; ++level)
1036 *   {
1037 *   const IndexSet dof_set =
1039 *   level);
1040 *  
1041 *   {
1042 *   #ifdef USE_PETSC_LA
1043 *   DynamicSparsityPattern dsp(dof_set);
1044 *   MGTools::make_sparsity_pattern(dof_handler, dsp, level);
1045 *   dsp.compress();
1047 *   dsp,
1048 *   dof_handler.locally_owned_mg_dofs(level),
1049 *   mpi_communicator,
1050 *   dof_set);
1051 *  
1052 *   mg_matrix[level].reinit(
1053 *   dof_handler.locally_owned_mg_dofs(level),
1054 *   dof_handler.locally_owned_mg_dofs(level),
1055 *   dsp,
1056 *   mpi_communicator);
1057 *   #else
1059 *   dof_handler.locally_owned_mg_dofs(level),
1060 *   dof_handler.locally_owned_mg_dofs(level),
1061 *   dof_set,
1062 *   mpi_communicator);
1063 *   MGTools::make_sparsity_pattern(dof_handler, dsp, level);
1064 *  
1065 *   dsp.compress();
1066 *   mg_matrix[level].reinit(dsp);
1067 *   #endif
1068 *   }
1069 *  
1070 *   {
1071 *   #ifdef USE_PETSC_LA
1072 *   DynamicSparsityPattern dsp(dof_set);
1074 *   mg_constrained_dofs,
1075 *   dsp,
1076 *   level);
1077 *   dsp.compress();
1079 *   dsp,
1080 *   dof_handler.locally_owned_mg_dofs(level),
1081 *   mpi_communicator,
1082 *   dof_set);
1083 *  
1084 *   mg_interface_in[level].reinit(
1085 *   dof_handler.locally_owned_mg_dofs(level),
1086 *   dof_handler.locally_owned_mg_dofs(level),
1087 *   dsp,
1088 *   mpi_communicator);
1089 *   #else
1091 *   dof_handler.locally_owned_mg_dofs(level),
1092 *   dof_handler.locally_owned_mg_dofs(level),
1093 *   dof_set,
1094 *   mpi_communicator);
1095 *  
1097 *   mg_constrained_dofs,
1098 *   dsp,
1099 *   level);
1100 *   dsp.compress();
1101 *   mg_interface_in[level].reinit(dsp);
1102 *   #endif
1103 *   }
1104 *   }
1105 *   break;
1106 *   }
1107 *  
1108 *   default:
1110 *   }
1111 *   }
1112 *  
1113 *  
1114 * @endcode
1115 *
1116 *
1117 * <a name="step_50-LaplaceProblemassemble_system"></a>
1118 * <h4>LaplaceProblem::assemble_system()</h4>
1119 *
1120
1121 *
1122 * The assembly is split into three parts: `assemble_system()`,
1123 * `assemble_multigrid()`, and `assemble_rhs()`. The
1124 * `assemble_system()` function here assembles and stores the (global)
1125 * system matrix and the right-hand side for the matrix-based
1126 * methods. It is similar to the assembly in @ref step_40 "step-40".
1127 *
1128
1129 *
1130 * Note that the matrix-free method does not execute this function as it does
1131 * not need to assemble a matrix, and it will instead assemble the right-hand
1132 * side in assemble_rhs().
1133 *
1134 * @code
1135 *   template <int dim, int degree>
1136 *   void LaplaceProblem<dim, degree>::assemble_system()
1137 *   {
1138 *   TimerOutput::Scope timing(computing_timer, "Assemble");
1139 *  
1140 *   const QGauss<dim> quadrature_formula(degree + 1);
1141 *  
1142 *   FEValues<dim> fe_values(fe,
1143 *   quadrature_formula,
1146 *  
1147 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1148 *   const unsigned int n_q_points = quadrature_formula.size();
1149 *  
1150 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
1151 *   Vector<double> cell_rhs(dofs_per_cell);
1152 *  
1153 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1154 *  
1155 *   const Coefficient<dim> coefficient;
1156 *   RightHandSide<dim> rhs;
1157 *   std::vector<double> rhs_values(n_q_points);
1158 *  
1159 *   for (const auto &cell : dof_handler.active_cell_iterators())
1160 *   if (cell->is_locally_owned())
1161 *   {
1162 *   cell_matrix = 0;
1163 *   cell_rhs = 0;
1164 *  
1165 *   fe_values.reinit(cell);
1166 *  
1167 *   const double coefficient_value =
1168 *   coefficient.average_value(fe_values.get_quadrature_points());
1169 *   rhs.value_list(fe_values.get_quadrature_points(), rhs_values);
1170 *  
1171 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1172 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1173 *   {
1174 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1175 *   cell_matrix(i, j) +=
1176 *   coefficient_value * // epsilon(x)
1177 *   fe_values.shape_grad(i, q_point) * // * grad phi_i(x)
1178 *   fe_values.shape_grad(j, q_point) * // * grad phi_j(x)
1179 *   fe_values.JxW(q_point); // * dx
1180 *  
1181 *   cell_rhs(i) +=
1182 *   fe_values.shape_value(i, q_point) * // grad phi_i(x)
1183 *   rhs_values[q_point] * // * f(x)
1184 *   fe_values.JxW(q_point); // * dx
1185 *   }
1186 *  
1187 *   cell->get_dof_indices(local_dof_indices);
1188 *   constraints.distribute_local_to_global(cell_matrix,
1189 *   cell_rhs,
1190 *   local_dof_indices,
1191 *   system_matrix,
1192 *   right_hand_side);
1193 *   }
1194 *  
1195 *   system_matrix.compress(VectorOperation::add);
1196 *   right_hand_side.compress(VectorOperation::add);
1197 *   }
1198 *  
1199 *  
1200 * @endcode
1201 *
1202 *
1203 * <a name="step_50-LaplaceProblemassemble_multigrid"></a>
1204 * <h4>LaplaceProblem::assemble_multigrid()</h4>
1205 *
1206
1207 *
1208 * The following function assembles and stores the multilevel matrices for the
1209 * matrix-based GMG method. This function is similar to the one found in
1210 * @ref step_16 "step-16", only here it works for distributed meshes. This difference amounts
1211 * to adding a condition that we only assemble on locally owned level cells
1212 * and a call to compress() for each matrix that is built.
1213 *
1214 * @code
1215 *   template <int dim, int degree>
1216 *   void LaplaceProblem<dim, degree>::assemble_multigrid()
1217 *   {
1218 *   TimerOutput::Scope timing(computing_timer, "Assemble multigrid");
1219 *  
1220 *   const QGauss<dim> quadrature_formula(degree + 1);
1221 *  
1222 *   FEValues<dim> fe_values(fe,
1223 *   quadrature_formula,
1226 *  
1227 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1228 *   const unsigned int n_q_points = quadrature_formula.size();
1229 *  
1230 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
1231 *  
1232 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1233 *  
1234 *   const Coefficient<dim> coefficient;
1235 *  
1236 *   std::vector<AffineConstraints<double>> boundary_constraints(
1237 *   triangulation.n_global_levels());
1238 *   for (unsigned int level = 0; level < triangulation.n_global_levels();
1239 *   ++level)
1240 *   {
1241 *   boundary_constraints[level].reinit(
1242 *   dof_handler.locally_owned_mg_dofs(level),
1244 *  
1245 *   for (const types::global_dof_index dof_index :
1246 *   mg_constrained_dofs.get_refinement_edge_indices(level))
1247 *   boundary_constraints[level].constrain_dof_to_zero(dof_index);
1248 *   for (const types::global_dof_index dof_index :
1249 *   mg_constrained_dofs.get_boundary_indices(level))
1250 *   boundary_constraints[level].constrain_dof_to_zero(dof_index);
1251 *   boundary_constraints[level].close();
1252 *   }
1253 *  
1254 *   for (const auto &cell : dof_handler.cell_iterators())
1255 *   if (cell->level_subdomain_id() == triangulation.locally_owned_subdomain())
1256 *   {
1257 *   cell_matrix = 0;
1258 *   fe_values.reinit(cell);
1259 *  
1260 *   const double coefficient_value =
1261 *   coefficient.average_value(fe_values.get_quadrature_points());
1262 *  
1263 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1264 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1265 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1266 *   cell_matrix(i, j) +=
1267 *   coefficient_value * fe_values.shape_grad(i, q_point) *
1268 *   fe_values.shape_grad(j, q_point) * fe_values.JxW(q_point);
1269 *  
1270 *   cell->get_mg_dof_indices(local_dof_indices);
1271 *  
1272 *   boundary_constraints[cell->level()].distribute_local_to_global(
1273 *   cell_matrix, local_dof_indices, mg_matrix[cell->level()]);
1274 *  
1275 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1276 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1277 *   if (mg_constrained_dofs.is_interface_matrix_entry(
1278 *   cell->level(), local_dof_indices[i], local_dof_indices[j]))
1279 *   mg_interface_in[cell->level()].add(local_dof_indices[i],
1280 *   local_dof_indices[j],
1281 *   cell_matrix(i, j));
1282 *   }
1283 *  
1284 *   for (unsigned int i = 0; i < triangulation.n_global_levels(); ++i)
1285 *   {
1286 *   mg_matrix[i].compress(VectorOperation::add);
1287 *   mg_interface_in[i].compress(VectorOperation::add);
1288 *   }
1289 *   }
1290 *  
1291 *  
1292 *  
1293 * @endcode
1294 *
1295 *
1296 * <a name="step_50-LaplaceProblemassemble_rhs"></a>
1297 * <h4>LaplaceProblem::assemble_rhs()</h4>
1298 *
1299
1300 *
1301 * The final function in this triptych assembles the right-hand side
1302 * vector for the matrix-free method -- because in the matrix-free
1303 * framework, we don't have to assemble the matrix and can get away
1304 * with only assembling the right hand side. We could do this by extracting
1305 * the code from the `assemble_system()` function above that deals with the
1306 * right hand side, but we decide instead to go all in on the matrix-free
1307 * approach and do the assembly using that way as well.
1308 *
1309
1310 *
1311 * The result is a function that is similar
1312 * to the one found in the "Use FEEvaluation::read_dof_values_plain()
1313 * to avoid resolving constraints" subsection in the "Possibilities
1314 * for extensions" section of @ref step_37 "step-37".
1315 *
1316
1317 *
1318 * The reason for this function is that the MatrixFree operators do not take
1319 * into account non-homogeneous Dirichlet constraints, instead treating all
1320 * Dirichlet constraints as homogeneous. To account for this, the right-hand
1321 * side here is assembled as the residual @f$r_0 = f-Au_0@f$, where @f$u_0@f$ is a
1322 * zero vector except in the Dirichlet values. Then when solving, we have that
1323 * the solution is @f$u = u_0 + A^{-1}r_0@f$. This can be seen as a Newton
1324 * iteration on a linear system with initial guess @f$u_0@f$. The CG solve in the
1325 * `solve()` function below computes @f$A^{-1}r_0@f$ and the call to
1326 * `constraints.distribute()` (which directly follows) adds the @f$u_0@f$.
1327 *
1328
1329 *
1330 * Obviously, since we are considering a problem with zero Dirichlet boundary,
1331 * we could have taken a similar approach to @ref step_37 "step-37" `assemble_rhs()`, but
1332 * this additional work allows us to change the problem declaration if we so
1333 * choose.
1334 *
1335
1336 *
1337 * This function has two parts in the integration loop: applying the negative
1338 * of matrix @f$A@f$ to @f$u_0@f$ by submitting the negative of the gradient, and
1339 * adding the right-hand side contribution by submitting the value @f$f@f$. We
1340 * must be sure to use `read_dof_values_plain()` for evaluating @f$u_0@f$ as
1341 * `read_dof_values()` would set all Dirichlet values to zero.
1342 *
1343
1344 *
1345 * Finally, the system_rhs vector is of type LA::MPI::Vector, but the
1346 * MatrixFree class only work for
1347 * LinearAlgebra::distributed::Vector. Therefore we must
1348 * compute the right-hand side using MatrixFree functionality and then
1349 * use the functions in the `ChangeVectorType` namespace to copy it to
1350 * the correct type.
1351 *
1352 * @code
1353 *   template <int dim, int degree>
1354 *   void LaplaceProblem<dim, degree>::assemble_rhs()
1355 *   {
1356 *   TimerOutput::Scope timing(computing_timer, "Assemble right-hand side");
1357 *  
1358 *   MatrixFreeActiveVector solution_copy;
1359 *   MatrixFreeActiveVector right_hand_side_copy;
1360 *   mf_system_matrix.initialize_dof_vector(solution_copy);
1361 *   mf_system_matrix.initialize_dof_vector(right_hand_side_copy);
1362 *  
1363 *   solution_copy = 0.;
1364 *   constraints.distribute(solution_copy);
1365 *   solution_copy.update_ghost_values();
1366 *   right_hand_side_copy = 0;
1367 *   const Table<2, VectorizedArray<double>> &coefficient =
1368 *   *(mf_system_matrix.get_coefficient());
1369 *  
1370 *   RightHandSide<dim> right_hand_side_function;
1371 *  
1372 *   FEEvaluation<dim, degree, degree + 1, 1, double> phi(
1373 *   *mf_system_matrix.get_matrix_free());
1374 *  
1375 *   for (unsigned int cell = 0;
1376 *   cell < mf_system_matrix.get_matrix_free()->n_cell_batches();
1377 *   ++cell)
1378 *   {
1379 *   phi.reinit(cell);
1380 *   phi.read_dof_values_plain(solution_copy);
1381 *   phi.evaluate(EvaluationFlags::gradients);
1382 *  
1383 *   for (const unsigned int q : phi.quadrature_point_indices())
1384 *   {
1385 *   phi.submit_gradient(-1.0 *
1386 *   (coefficient(cell, 0) * phi.get_gradient(q)),
1387 *   q);
1388 *   phi.submit_value(
1389 *   right_hand_side_function.value(phi.quadrature_point(q)), q);
1390 *   }
1391 *  
1392 *   phi.integrate_scatter(EvaluationFlags::values |
1393 *   EvaluationFlags::gradients,
1394 *   right_hand_side_copy);
1395 *   }
1396 *  
1397 *   right_hand_side_copy.compress(VectorOperation::add);
1398 *  
1399 *   ChangeVectorTypes::copy(right_hand_side, right_hand_side_copy);
1400 *   }
1401 *  
1402 *  
1403 *  
1404 * @endcode
1405 *
1406 *
1407 * <a name="step_50-LaplaceProblemsolve"></a>
1408 * <h4>LaplaceProblem::solve()</h4>
1409 *
1410
1411 *
1412 * Here we set up the multigrid preconditioner, test the timing of a single
1413 * V-cycle, and solve the linear system. Unsurprisingly, this is one of the
1414 * places where the three methods differ the most.
1415 *
1416 * @code
1417 *   template <int dim, int degree>
1418 *   void LaplaceProblem<dim, degree>::solve()
1419 *   {
1420 *   TimerOutput::Scope timing(computing_timer, "Solve");
1421 *  
1422 *   SolverControl solver_control(1000, 1.e-10 * right_hand_side.l2_norm());
1423 *   solver_control.enable_history_data();
1424 *  
1425 *   solution = 0.;
1426 *  
1427 * @endcode
1428 *
1429 * The solver for the matrix-free GMG method is similar to @ref step_37 "step-37", apart
1430 * from adding some interface matrices in complete analogy to @ref step_16 "step-16".
1431 *
1432 * @code
1433 *   switch (settings.solver)
1434 *   {
1435 *   case Settings::gmg_mf:
1436 *   {
1437 *   computing_timer.enter_subsection("Solve: Preconditioner setup");
1438 *  
1439 *   MGTransferMatrixFree<dim, float> mg_transfer(mg_constrained_dofs);
1440 *   mg_transfer.build(dof_handler);
1441 *  
1442 *   SolverControl coarse_solver_control(1000, 1e-12, false, false);
1443 *   SolverCG<MatrixFreeLevelVector> coarse_solver(
1444 *   coarse_solver_control);
1445 *   PreconditionIdentity identity;
1446 *   MGCoarseGridIterativeSolver<MatrixFreeLevelVector,
1447 *   SolverCG<MatrixFreeLevelVector>,
1448 *   MatrixFreeLevelMatrix,
1449 *   PreconditionIdentity>
1450 *   coarse_grid_solver(coarse_solver, mf_mg_matrix[0], identity);
1451 *  
1452 *   using Smoother = PreconditionJacobi<MatrixFreeLevelMatrix>;
1453 *   MGSmootherPrecondition<MatrixFreeLevelMatrix,
1454 *   Smoother,
1455 *   MatrixFreeLevelVector>
1456 *   smoother;
1457 *   smoother.initialize(mf_mg_matrix,
1458 *   typename Smoother::AdditionalData(
1459 *   settings.smoother_dampen));
1460 *   smoother.set_steps(settings.smoother_steps);
1461 *  
1462 *   mg::Matrix<MatrixFreeLevelVector> mg_m(mf_mg_matrix);
1463 *  
1464 *   MGLevelObject<
1465 *   MatrixFreeOperators::MGInterfaceOperator<MatrixFreeLevelMatrix>>
1466 *   mg_interface_matrices;
1467 *   mg_interface_matrices.resize(0,
1468 *   triangulation.n_global_levels() - 1);
1469 *   for (unsigned int level = 0;
1470 *   level < triangulation.n_global_levels();
1471 *   ++level)
1472 *   mg_interface_matrices[level].initialize(mf_mg_matrix[level]);
1473 *   mg::Matrix<MatrixFreeLevelVector> mg_interface(
1474 *   mg_interface_matrices);
1475 *  
1476 *   Multigrid<MatrixFreeLevelVector> mg(
1477 *   mg_m, coarse_grid_solver, mg_transfer, smoother, smoother);
1478 *   mg.set_edge_matrices(mg_interface, mg_interface);
1479 *  
1480 *   PreconditionMG<dim,
1481 *   MatrixFreeLevelVector,
1482 *   MGTransferMatrixFree<dim, float>>
1483 *   preconditioner(dof_handler, mg, mg_transfer);
1484 *  
1485 * @endcode
1486 *
1487 * Copy the solution vector and right-hand side from LA::MPI::Vector
1488 * to LinearAlgebra::distributed::Vector so that we can solve.
1489 *
1490 * @code
1491 *   MatrixFreeActiveVector solution_copy;
1492 *   MatrixFreeActiveVector right_hand_side_copy;
1493 *   mf_system_matrix.initialize_dof_vector(solution_copy);
1494 *   mf_system_matrix.initialize_dof_vector(right_hand_side_copy);
1495 *  
1496 *   ChangeVectorTypes::copy(solution_copy, solution);
1497 *   ChangeVectorTypes::copy(right_hand_side_copy, right_hand_side);
1498 *   computing_timer.leave_subsection("Solve: Preconditioner setup");
1499 *  
1500 * @endcode
1501 *
1502 * Timing for 1 V-cycle.
1503 *
1504 * @code
1505 *   {
1506 *   TimerOutput::Scope timing(computing_timer,
1507 *   "Solve: 1 multigrid V-cycle");
1508 *   preconditioner.vmult(solution_copy, right_hand_side_copy);
1509 *   }
1510 *   solution_copy = 0.;
1511 *  
1512 * @endcode
1513 *
1514 * Solve the linear system, update the ghost values of the solution,
1515 * copy back to LA::MPI::Vector and distribute constraints.
1516 *
1517 * @code
1518 *   {
1519 *   SolverCG<MatrixFreeActiveVector> solver(solver_control);
1520 *  
1521 *   TimerOutput::Scope timing(computing_timer, "Solve: CG");
1522 *   solver.solve(mf_system_matrix,
1523 *   solution_copy,
1524 *   right_hand_side_copy,
1525 *   preconditioner);
1526 *   }
1527 *  
1528 *   solution_copy.update_ghost_values();
1529 *   ChangeVectorTypes::copy(solution, solution_copy);
1530 *   constraints.distribute(solution);
1531 *  
1532 *   break;
1533 *   }
1534 *  
1535 * @endcode
1536 *
1537 * Solver for the matrix-based GMG method, similar to @ref step_16 "step-16", only
1538 * using a Jacobi smoother instead of a SOR smoother (which is not
1539 * implemented in parallel).
1540 *
1541 * @code
1542 *   case Settings::gmg_mb:
1543 *   {
1544 *   computing_timer.enter_subsection("Solve: Preconditioner setup");
1545 *  
1546 *   MGTransferPrebuilt<VectorType> mg_transfer(mg_constrained_dofs);
1547 *   mg_transfer.build(dof_handler);
1548 *  
1549 *   SolverControl coarse_solver_control(1000, 1e-12, false, false);
1550 *   SolverCG<VectorType> coarse_solver(coarse_solver_control);
1551 *   PreconditionIdentity identity;
1552 *   MGCoarseGridIterativeSolver<VectorType,
1553 *   SolverCG<VectorType>,
1554 *   MatrixType,
1555 *   PreconditionIdentity>
1556 *   coarse_grid_solver(coarse_solver, mg_matrix[0], identity);
1557 *  
1558 *   using Smoother = LA::MPI::PreconditionJacobi;
1559 *   MGSmootherPrecondition<MatrixType, Smoother, VectorType> smoother;
1560 *  
1561 *   #ifdef USE_PETSC_LA
1562 *   smoother.initialize(mg_matrix);
1563 *   Assert(
1564 *   settings.smoother_dampen == 1.0,
1565 *   ExcNotImplemented(
1566 *   "PETSc's PreconditionJacobi has no support for a damping parameter."));
1567 *   #else
1568 *   smoother.initialize(mg_matrix, settings.smoother_dampen);
1569 *   #endif
1570 *  
1571 *   smoother.set_steps(settings.smoother_steps);
1572 *  
1573 *   mg::Matrix<VectorType> mg_m(mg_matrix);
1574 *   mg::Matrix<VectorType> mg_in(mg_interface_in);
1575 *   mg::Matrix<VectorType> mg_out(mg_interface_in);
1576 *  
1577 *   Multigrid<VectorType> mg(
1578 *   mg_m, coarse_grid_solver, mg_transfer, smoother, smoother);
1579 *   mg.set_edge_matrices(mg_out, mg_in);
1580 *  
1581 *  
1582 *   PreconditionMG<dim, VectorType, MGTransferPrebuilt<VectorType>>
1583 *   preconditioner(dof_handler, mg, mg_transfer);
1584 *  
1585 *   computing_timer.leave_subsection("Solve: Preconditioner setup");
1586 *  
1587 * @endcode
1588 *
1589 * Timing for 1 V-cycle.
1590 *
1591 * @code
1592 *   {
1593 *   TimerOutput::Scope timing(computing_timer,
1594 *   "Solve: 1 multigrid V-cycle");
1595 *   preconditioner.vmult(solution, right_hand_side);
1596 *   }
1597 *   solution = 0.;
1598 *  
1599 * @endcode
1600 *
1601 * Solve the linear system and distribute constraints.
1602 *
1603 * @code
1604 *   {
1605 *   SolverCG<VectorType> solver(solver_control);
1606 *  
1607 *   TimerOutput::Scope timing(computing_timer, "Solve: CG");
1608 *   solver.solve(system_matrix,
1609 *   solution,
1610 *   right_hand_side,
1611 *   preconditioner);
1612 *   }
1613 *  
1614 *   constraints.distribute(solution);
1615 *  
1616 *   break;
1617 *   }
1618 *  
1619 * @endcode
1620 *
1621 * Solver for the AMG method, similar to @ref step_40 "step-40".
1622 *
1623 * @code
1624 *   case Settings::amg:
1625 *   {
1626 *   computing_timer.enter_subsection("Solve: Preconditioner setup");
1627 *  
1628 *   PreconditionAMG preconditioner;
1629 *   PreconditionAMG::AdditionalData Amg_data;
1630 *  
1631 *   #ifdef USE_PETSC_LA
1632 *   Amg_data.symmetric_operator = true;
1633 *   #else
1634 *   Amg_data.elliptic = true;
1635 *   Amg_data.smoother_type = "Jacobi";
1636 *   Amg_data.higher_order_elements = true;
1637 *   Amg_data.smoother_sweeps = settings.smoother_steps;
1638 *   Amg_data.aggregation_threshold = 0.02;
1639 *   #endif
1640 *  
1641 *   Amg_data.output_details = false;
1642 *  
1643 *   preconditioner.initialize(system_matrix, Amg_data);
1644 *   computing_timer.leave_subsection("Solve: Preconditioner setup");
1645 *  
1646 * @endcode
1647 *
1648 * Timing for 1 V-cycle.
1649 *
1650 * @code
1651 *   {
1652 *   TimerOutput::Scope timing(computing_timer,
1653 *   "Solve: 1 multigrid V-cycle");
1654 *   preconditioner.vmult(solution, right_hand_side);
1655 *   }
1656 *   solution = 0.;
1657 *  
1658 * @endcode
1659 *
1660 * Solve the linear system and distribute constraints.
1661 *
1662 * @code
1663 *   {
1664 *   SolverCG<VectorType> solver(solver_control);
1665 *  
1666 *   TimerOutput::Scope timing(computing_timer, "Solve: CG");
1667 *   solver.solve(system_matrix,
1668 *   solution,
1669 *   right_hand_side,
1670 *   preconditioner);
1671 *   }
1672 *   constraints.distribute(solution);
1673 *  
1674 *   break;
1675 *   }
1676 *  
1677 *   default:
1678 *   DEAL_II_ASSERT_UNREACHABLE();
1679 *   }
1680 *  
1681 *   pcout << " Number of CG iterations: " << solver_control.last_step()
1682 *   << std::endl;
1683 *   }
1684 *  
1685 *  
1686 * @endcode
1687 *
1688 *
1689 * <a name="step_50-Theerrorestimator"></a>
1690 * <h3>The error estimator</h3>
1691 *
1692
1693 *
1694 * We use the FEInterfaceValues class to assemble an error estimator to decide
1695 * which cells to refine. See the exact definition of the cell and face
1696 * integrals in the introduction. To use the method, we define Scratch and
1697 * Copy objects for the MeshWorker::mesh_loop() with much of the following
1698 * code being in essence as was set up in @ref step_12 "step-12" already (or at least similar
1699 * in spirit).
1700 *
1701 * @code
1702 *   template <int dim>
1703 *   struct ScratchData
1704 *   {
1705 *   ScratchData(const Mapping<dim> &mapping,
1706 *   const FiniteElement<dim> &fe,
1707 *   const unsigned int quadrature_degree,
1708 *   const UpdateFlags update_flags,
1709 *   const UpdateFlags interface_update_flags)
1710 *   : fe_values(mapping, fe, QGauss<dim>(quadrature_degree), update_flags)
1711 *   , fe_interface_values(mapping,
1712 *   fe,
1713 *   QGauss<dim - 1>(quadrature_degree),
1714 *   interface_update_flags)
1715 *   {}
1716 *  
1717 *  
1718 *   ScratchData(const ScratchData<dim> &scratch_data)
1719 *   : fe_values(scratch_data.fe_values.get_mapping(),
1720 *   scratch_data.fe_values.get_fe(),
1721 *   scratch_data.fe_values.get_quadrature(),
1722 *   scratch_data.fe_values.get_update_flags())
1723 *   , fe_interface_values(scratch_data.fe_values.get_mapping(),
1724 *   scratch_data.fe_values.get_fe(),
1725 *   scratch_data.fe_interface_values.get_quadrature(),
1726 *   scratch_data.fe_interface_values.get_update_flags())
1727 *   {}
1728 *  
1729 *   FEValues<dim> fe_values;
1730 *   FEInterfaceValues<dim> fe_interface_values;
1731 *   };
1732 *  
1733 *  
1734 *  
1735 *   struct CopyData
1736 *   {
1737 *   CopyData()
1738 *   : cell_index(numbers::invalid_unsigned_int)
1739 *   , value(0.)
1740 *   {}
1741 *  
1742 *   struct FaceData
1743 *   {
1744 *   unsigned int cell_indices[2];
1745 *   double values[2];
1746 *   };
1747 *  
1748 *   unsigned int cell_index;
1749 *   double value;
1750 *   std::vector<FaceData> face_data;
1751 *   };
1752 *  
1753 *  
1754 *   template <int dim, int degree>
1755 *   void LaplaceProblem<dim, degree>::estimate()
1756 *   {
1757 *   TimerOutput::Scope timing(computing_timer, "Estimate");
1758 *  
1759 *   VectorType temp_solution;
1760 *   temp_solution.reinit(locally_owned_dofs,
1761 *   locally_relevant_dofs,
1762 *   mpi_communicator);
1763 *   temp_solution = solution;
1764 *  
1765 *   const Coefficient<dim> coefficient;
1766 *  
1767 *   estimated_error_square_per_cell.reinit(triangulation.n_active_cells());
1768 *  
1769 *   using Iterator = typename DoFHandler<dim>::active_cell_iterator;
1770 *  
1771 * @endcode
1772 *
1773 * Assembler for cell residual @f$h^2 \| f + \epsilon \triangle u \|_K^2@f$
1774 *
1775 * @code
1776 *   auto cell_worker = [&](const Iterator &cell,
1777 *   ScratchData<dim> &scratch_data,
1778 *   CopyData &copy_data) {
1779 *   FEValues<dim> &fe_values = scratch_data.fe_values;
1780 *   fe_values.reinit(cell);
1781 *  
1782 *   RightHandSide<dim> rhs;
1783 *   const double rhs_value = rhs.value(cell->center());
1784 *  
1785 *   const double nu = coefficient.value(cell->center());
1786 *  
1787 *   std::vector<Tensor<2, dim>> hessians(fe_values.n_quadrature_points);
1788 *   fe_values.get_function_hessians(temp_solution, hessians);
1789 *  
1790 *   copy_data.cell_index = cell->active_cell_index();
1791 *  
1792 *   double residual_norm_square = 0.;
1793 *   for (unsigned k = 0; k < fe_values.n_quadrature_points; ++k)
1794 *   {
1795 *   const double residual = (rhs_value + nu * trace(hessians[k]));
1796 *   residual_norm_square += residual * residual * fe_values.JxW(k);
1797 *   }
1798 *  
1799 *   copy_data.value =
1800 *   cell->diameter() * cell->diameter() * residual_norm_square;
1801 *   };
1802 *  
1803 * @endcode
1804 *
1805 * Assembler for face term @f$\sum_F h_F \| \jump{\epsilon \nabla u \cdot n}
1806 * \|_F^2@f$
1807 *
1808 * @code
1809 *   auto face_worker = [&](const Iterator &cell,
1810 *   const unsigned int &f,
1811 *   const unsigned int &sf,
1812 *   const Iterator &ncell,
1813 *   const unsigned int &nf,
1814 *   const unsigned int &nsf,
1815 *   ScratchData<dim> &scratch_data,
1816 *   CopyData &copy_data) {
1817 *   FEInterfaceValues<dim> &fe_interface_values =
1818 *   scratch_data.fe_interface_values;
1819 *   fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
1820 *  
1821 *   copy_data.face_data.emplace_back();
1822 *   CopyData::FaceData &copy_data_face = copy_data.face_data.back();
1823 *  
1824 *   copy_data_face.cell_indices[0] = cell->active_cell_index();
1825 *   copy_data_face.cell_indices[1] = ncell->active_cell_index();
1826 *  
1827 *   const double coeff1 = coefficient.value(cell->center());
1828 *   const double coeff2 = coefficient.value(ncell->center());
1829 *  
1830 *   std::vector<Tensor<1, dim>> grad_u[2];
1831 *  
1832 *   for (unsigned int i = 0; i < 2; ++i)
1833 *   {
1834 *   grad_u[i].resize(fe_interface_values.n_quadrature_points);
1835 *   fe_interface_values.get_fe_face_values(i).get_function_gradients(
1836 *   temp_solution, grad_u[i]);
1837 *   }
1838 *  
1839 *   double jump_norm_square = 0.;
1840 *  
1841 *   for (unsigned int qpoint = 0;
1842 *   qpoint < fe_interface_values.n_quadrature_points;
1843 *   ++qpoint)
1844 *   {
1845 *   const double jump =
1846 *   coeff1 * grad_u[0][qpoint] * fe_interface_values.normal(qpoint) -
1847 *   coeff2 * grad_u[1][qpoint] * fe_interface_values.normal(qpoint);
1848 *  
1849 *   jump_norm_square += jump * jump * fe_interface_values.JxW(qpoint);
1850 *   }
1851 *  
1852 *   const double h = cell->face(f)->measure();
1853 *   copy_data_face.values[0] = 0.5 * h * jump_norm_square;
1854 *   copy_data_face.values[1] = copy_data_face.values[0];
1855 *   };
1856 *  
1857 *   auto copier = [&](const CopyData &copy_data) {
1858 *   if (copy_data.cell_index != numbers::invalid_unsigned_int)
1859 *   estimated_error_square_per_cell[copy_data.cell_index] +=
1860 *   copy_data.value;
1861 *  
1862 *   for (const auto &cdf : copy_data.face_data)
1863 *   for (unsigned int j = 0; j < 2; ++j)
1864 *   estimated_error_square_per_cell[cdf.cell_indices[j]] += cdf.values[j];
1865 *   };
1866 *  
1867 *   const unsigned int n_gauss_points = degree + 1;
1868 *   ScratchData<dim> scratch_data(mapping,
1869 *   fe,
1870 *   n_gauss_points,
1871 *   update_hessians | update_quadrature_points |
1872 *   update_JxW_values,
1873 *   update_values | update_gradients |
1874 *   update_JxW_values | update_normal_vectors);
1875 *   CopyData copy_data;
1876 *  
1877 * @endcode
1878 *
1879 * We need to assemble each interior face once but we need to make sure that
1880 * both processes assemble the face term between a locally owned and a ghost
1881 * cell. This is achieved by setting the
1882 * MeshWorker::assemble_ghost_faces_both flag. We need to do this, because
1883 * we do not communicate the error estimator contributions here.
1884 *
1885 * @code
1886 *   MeshWorker::mesh_loop(dof_handler.begin_active(),
1887 *   dof_handler.end(),
1888 *   cell_worker,
1889 *   copier,
1890 *   scratch_data,
1891 *   copy_data,
1892 *   MeshWorker::assemble_own_cells |
1893 *   MeshWorker::assemble_ghost_faces_both |
1894 *   MeshWorker::assemble_own_interior_faces_once,
1895 *   /*boundary_worker=*/nullptr,
1896 *   face_worker);
1897 *  
1898 *   const double global_error_estimate =
1899 *   std::sqrt(Utilities::MPI::sum(estimated_error_square_per_cell.l1_norm(),
1900 *   mpi_communicator));
1901 *   pcout << " Global error estimate: " << global_error_estimate
1902 *   << std::endl;
1903 *   }
1904 *  
1905 *  
1906 * @endcode
1907 *
1908 *
1909 * <a name="step_50-LaplaceProblemrefine_grid"></a>
1910 * <h4>LaplaceProblem::refine_grid()</h4>
1911 *
1912
1913 *
1914 * We use the cell-wise estimator stored in the vector @p estimate_vector and
1915 * refine a fixed number of cells (chosen here to roughly double the number of
1916 * DoFs in each step).
1917 *
1918 * @code
1919 *   template <int dim, int degree>
1920 *   void LaplaceProblem<dim, degree>::refine_grid()
1921 *   {
1922 *   TimerOutput::Scope timing(computing_timer, "Refine grid");
1923 *  
1924 *   const double refinement_fraction = 1. / (std::pow(2.0, dim) - 1.);
1925 *   parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
1926 *   triangulation, estimated_error_square_per_cell, refinement_fraction, 0.0);
1927 *  
1928 *   triangulation.execute_coarsening_and_refinement();
1929 *   }
1930 *  
1931 *  
1932 * @endcode
1933 *
1934 *
1935 * <a name="step_50-LaplaceProblemoutput_results"></a>
1936 * <h4>LaplaceProblem::output_results()</h4>
1937 *
1938
1939 *
1940 * The output_results() function is similar to the ones found in many of the
1941 * tutorials (see @ref step_40 "step-40" for example).
1942 *
1943 * @code
1944 *   template <int dim, int degree>
1945 *   void LaplaceProblem<dim, degree>::output_results(const unsigned int cycle)
1946 *   {
1947 *   TimerOutput::Scope timing(computing_timer, "Output results");
1948 *  
1949 *   VectorType temp_solution;
1950 *   temp_solution.reinit(locally_owned_dofs,
1951 *   locally_relevant_dofs,
1952 *   mpi_communicator);
1953 *   temp_solution = solution;
1954 *  
1955 *   DataOut<dim> data_out;
1956 *   data_out.attach_dof_handler(dof_handler);
1957 *   data_out.add_data_vector(temp_solution, "solution");
1958 *  
1959 *   Vector<float> subdomain(triangulation.n_active_cells());
1960 *   for (unsigned int i = 0; i < subdomain.size(); ++i)
1961 *   subdomain(i) = triangulation.locally_owned_subdomain();
1962 *   data_out.add_data_vector(subdomain, "subdomain");
1963 *  
1964 *   Vector<float> level(triangulation.n_active_cells());
1965 *   for (const auto &cell : triangulation.active_cell_iterators())
1966 *   level(cell->active_cell_index()) = cell->level();
1967 *   data_out.add_data_vector(level, "level");
1968 *  
1969 *   if (estimated_error_square_per_cell.size() > 0)
1970 *   data_out.add_data_vector(estimated_error_square_per_cell,
1971 *   "estimated_error_square_per_cell");
1972 *  
1973 *   data_out.build_patches();
1974 *  
1975 *   const std::string pvtu_filename = data_out.write_vtu_with_pvtu_record(
1976 *   "", "solution", cycle, mpi_communicator, 2 /*n_digits*/, 1 /*n_groups*/);
1977 *  
1978 *   pcout << " Wrote " << pvtu_filename << std::endl;
1979 *   }
1980 *  
1981 *  
1982 * @endcode
1983 *
1984 *
1985 * <a name="step_50-LaplaceProblemrun"></a>
1986 * <h4>LaplaceProblem::run()</h4>
1987 *
1988
1989 *
1990 * As in most tutorials, this function calls the various functions defined
1991 * above to set up, assemble, solve, and output the results.
1992 *
1993 * @code
1994 *   template <int dim, int degree>
1995 *   void LaplaceProblem<dim, degree>::run()
1996 *   {
1997 *   for (unsigned int cycle = 0; cycle < settings.n_steps; ++cycle)
1998 *   {
1999 *   pcout << "Cycle " << cycle << ':' << std::endl;
2000 *   if (cycle > 0)
2001 *   refine_grid();
2002 *  
2003 *   pcout << " Number of active cells: "
2004 *   << triangulation.n_global_active_cells();
2005 *  
2006 * @endcode
2007 *
2008 * We only output level cell data for the GMG methods (same with DoF
2009 * data below). Note that the partition efficiency is irrelevant for AMG
2010 * since the level hierarchy is not distributed or used during the
2011 * computation.
2012 *
2013 * @code
2014 *   if (settings.solver == Settings::gmg_mf ||
2015 *   settings.solver == Settings::gmg_mb)
2016 *   pcout << " (" << triangulation.n_global_levels() << " global levels)"
2017 *   << std::endl
2018 *   << " Partition efficiency: "
2019 *   << 1.0 / MGTools::workload_imbalance(triangulation);
2020 *   pcout << std::endl;
2021 *  
2022 *   setup_system();
2023 *  
2024 * @endcode
2025 *
2026 * Only set up the multilevel hierarchy for GMG.
2027 *
2028 * @code
2029 *   if (settings.solver == Settings::gmg_mf ||
2030 *   settings.solver == Settings::gmg_mb)
2031 *   setup_multigrid();
2032 *  
2033 *   pcout << " Number of degrees of freedom: " << dof_handler.n_dofs();
2034 *   if (settings.solver == Settings::gmg_mf ||
2035 *   settings.solver == Settings::gmg_mb)
2036 *   {
2037 *   pcout << " (by level: ";
2038 *   for (unsigned int level = 0;
2039 *   level < triangulation.n_global_levels();
2040 *   ++level)
2041 *   pcout << dof_handler.n_dofs(level)
2042 *   << (level == triangulation.n_global_levels() - 1 ? ")" :
2043 *   ", ");
2044 *   }
2045 *   pcout << std::endl;
2046 *  
2047 * @endcode
2048 *
2049 * For the matrix-free method, we only assemble the right-hand side.
2050 * For both matrix-based methods, we assemble both active matrix and
2051 * right-hand side, and only assemble the multigrid matrices for
2052 * matrix-based GMG.
2053 *
2054 * @code
2055 *   if (settings.solver == Settings::gmg_mf)
2056 *   assemble_rhs();
2057 *   else /*gmg_mb or amg*/
2058 *   {
2059 *   assemble_system();
2060 *   if (settings.solver == Settings::gmg_mb)
2061 *   assemble_multigrid();
2062 *   }
2063 *  
2064 *   solve();
2065 *   estimate();
2066 *  
2067 *   if (settings.output)
2068 *   output_results(cycle);
2069 *  
2070 *   computing_timer.print_summary();
2071 *   computing_timer.reset();
2072 *   }
2073 *   }
2074 *   } // namespace Step50
2075 *  
2076 * @endcode
2077 *
2078 *
2079 * <a name="step_50-Themainfunction"></a>
2080 * <h3>The main() function</h3>
2081 *
2082
2083 *
2084 * This is a similar main function to @ref step_40 "step-40", with the exception that
2085 * we require the user to pass a .prm file as a sole command line
2086 * argument (see @ref step_29 "step-29" and the documentation of the ParameterHandler
2087 * class for a complete discussion of parameter files).
2088 *
2089 * @code
2090 *   int main(int argc, char *argv[])
2091 *   {
2092 *   using namespace dealii;
2093 *   using namespace Step50;
2094 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
2095 *  
2096 *   Settings settings;
2097 *   if (!settings.try_parse((argc > 1) ? (argv[1]) : ""))
2098 *   return 0;
2099 *  
2100 *   try
2101 *   {
2102 *   constexpr unsigned int fe_degree = 2;
2103 *  
2104 *   switch (settings.dimension)
2105 *   {
2106 *   case 2:
2107 *   {
2108 *   LaplaceProblem<2, fe_degree> test(settings);
2109 *   test.run();
2110 *  
2111 *   break;
2112 *   }
2113 *  
2114 *   case 3:
2115 *   {
2116 *   LaplaceProblem<3, fe_degree> test(settings);
2117 *   test.run();
2118 *  
2119 *   break;
2120 *   }
2121 *  
2122 *   default:
2123 *   DEAL_II_NOT_IMPLEMENTED();
2124 *   }
2125 *   }
2126 *   catch (std::exception &exc)
2127 *   {
2128 *   std::cerr << std::endl
2129 *   << std::endl
2130 *   << "----------------------------------------------------"
2131 *   << std::endl;
2132 *   std::cerr << "Exception on processing: " << std::endl
2133 *   << exc.what() << std::endl
2134 *   << "Aborting!" << std::endl
2135 *   << "----------------------------------------------------"
2136 *   << std::endl;
2137 *   MPI_Abort(MPI_COMM_WORLD, 1);
2138 *   return 1;
2139 *   }
2140 *   catch (...)
2141 *   {
2142 *   std::cerr << std::endl
2143 *   << std::endl
2144 *   << "----------------------------------------------------"
2145 *   << std::endl;
2146 *   std::cerr << "Unknown exception!" << std::endl
2147 *   << "Aborting!" << std::endl
2148 *   << "----------------------------------------------------"
2149 *   << std::endl;
2150 *   MPI_Abort(MPI_COMM_WORLD, 2);
2151 *   return 1;
2152 *   }
2153 *  
2154 *   return 0;
2155 *   }
2156 * @endcode
2157<a name="step_50-Results"></a><h1>Results</h1>
2158
2159
2160When you run the program using the following command
2161@code
2162mpirun -np 16 ./step-50 gmg_mf_2d.prm
2163@endcode
2164the screen output should look like the following:
2165@code
2166Cycle 0:
2167 Number of active cells: 12 (2 global levels)
2168 Partition efficiency: 0.1875
2169 Number of degrees of freedom: 65 (by level: 21, 65)
2170 Number of CG iterations: 10
2171 Global error estimate: 0.355373
2172 Wrote solution_00.pvtu
2173
2174
2175+---------------------------------------------+------------+------------+
2176| Total wallclock time elapsed since start | 0.0163s | |
2177| | | |
2178| Section | no. calls | wall time | % of total |
2179+---------------------------------+-----------+------------+------------+
2180| Assemble right-hand side | 1 | 0.000374s | 2.3% |
2181| Estimate | 1 | 0.000724s | 4.4% |
2182| Output results | 1 | 0.00277s | 17% |
2183| Setup | 1 | 0.00225s | 14% |
2184| Setup multigrid | 1 | 0.00181s | 11% |
2185| Solve | 1 | 0.00364s | 22% |
2186| Solve: 1 multigrid V-cycle | 1 | 0.000354s | 2.2% |
2187| Solve: CG | 1 | 0.00151s | 9.3% |
2188| Solve: Preconditioner setup | 1 | 0.00125s | 7.7% |
2189+---------------------------------+-----------+------------+------------+
2190
2191Cycle 1:
2192 Number of active cells: 24 (3 global levels)
2193 Partition efficiency: 0.276786
2194 Number of degrees of freedom: 139 (by level: 21, 65, 99)
2195 Number of CG iterations: 10
2196 Global error estimate: 0.216726
2197 Wrote solution_01.pvtu
2198
2199
2200+---------------------------------------------+------------+------------+
2201| Total wallclock time elapsed since start | 0.0169s | |
2202| | | |
2203| Section | no. calls | wall time | % of total |
2204+---------------------------------+-----------+------------+------------+
2205| Assemble right-hand side | 1 | 0.000309s | 1.8% |
2206| Estimate | 1 | 0.00156s | 9.2% |
2207| Output results | 1 | 0.00222s | 13% |
2208| Refine grid | 1 | 0.00278s | 16% |
2209| Setup | 1 | 0.00196s | 12% |
2210| Setup multigrid | 1 | 0.0023s | 14% |
2211| Solve | 1 | 0.00565s | 33% |
2212| Solve: 1 multigrid V-cycle | 1 | 0.000349s | 2.1% |
2213| Solve: CG | 1 | 0.00285s | 17% |
2214| Solve: Preconditioner setup | 1 | 0.00195s | 12% |
2215+---------------------------------+-----------+------------+------------+
2216
2217Cycle 2:
2218 Number of active cells: 51 (4 global levels)
2219 Partition efficiency: 0.41875
2220 Number of degrees of freedom: 245 (by level: 21, 65, 225, 25)
2221 Number of CG iterations: 11
2222 Global error estimate: 0.112098
2223 Wrote solution_02.pvtu
2224
2225
2226+---------------------------------------------+------------+------------+
2227| Total wallclock time elapsed since start | 0.0183s | |
2228| | | |
2229| Section | no. calls | wall time | % of total |
2230+---------------------------------+-----------+------------+------------+
2231| Assemble right-hand side | 1 | 0.000274s | 1.5% |
2232| Estimate | 1 | 0.00127s | 6.9% |
2233| Output results | 1 | 0.00227s | 12% |
2234| Refine grid | 1 | 0.0024s | 13% |
2235| Setup | 1 | 0.00191s | 10% |
2236| Setup multigrid | 1 | 0.00295s | 16% |
2237| Solve | 1 | 0.00702s | 38% |
2238| Solve: 1 multigrid V-cycle | 1 | 0.000398s | 2.2% |
2239| Solve: CG | 1 | 0.00376s | 21% |
2240| Solve: Preconditioner setup | 1 | 0.00238s | 13% |
2241+---------------------------------+-----------+------------+------------+
2242.
2243.
2244.
2245@endcode
2246Here, the timing of the `solve()` function is split up in 3 parts: setting
2247up the multigrid preconditioner, execution of a single multigrid V-cycle, and
2248the CG solver. The V-cycle that is timed is unnecessary for the overall solve
2249and only meant to give an insight at the different costs for AMG and GMG.
2250Also it should be noted that when using the AMG solver, "Workload imbalance"
2251is not included in the output since the hierarchy of coarse meshes is not
2252required.
2253
2254All results in this section are gathered on Intel Xeon Platinum 8280 (Cascade
2255Lake) nodes which have 56 cores and 192GB per node and support AVX-512 instructions,
2256allowing for vectorization over 8 doubles (vectorization used only in the matrix-free
2257computations). The code is compiled using gcc 7.1.0 with intel-mpi 17.0.3. Trilinos
225812.10.1 is used for the matrix-based GMG/AMG computations.
2259
2260We can then gather a variety of information by calling the program
2261with the input files that are provided in the directory in which
2262@ref step_50 "step-50" is located. Using these, and adjusting the number of mesh
2263refinement steps, we can produce information about how well the
2264program scales.
2265
2266The following table gives weak scaling timings for this program on up to 256M DoFs
2267and 7,168 processors. (Recall that weak scaling keeps the number of
2268degrees of freedom per processor constant while increasing the number of
2269processors; i.e., it considers larger and larger problems.)
2270Here, @f$\mathbb{E}@f$ is the partition efficiency from the
2271 introduction (also equal to 1.0/workload imbalance), "Setup" is a combination
2272of setup, setup multigrid, assemble, and assemble multigrid from the timing blocks,
2273and "Prec" is the preconditioner setup. Ideally all times would stay constant
2274over each problem size for the individual solvers, but since the partition
2275efficiency decreases from 0.371 to 0.161 from largest to smallest problem size,
2276we expect to see an approximately @f$0.371/0.161=2.3@f$ times increase in timings
2277for GMG. This is, in fact, pretty close to what we really get:
2278
2279<table align="center" class="doxtable">
2280<tr>
2281 <th colspan="4"></th>
2282 <th></th>
2283 <th colspan="4">MF-GMG</th>
2284 <th></th>
2285 <th colspan="4">MB-GMG</th>
2286 <th></th>
2287 <th colspan="4">AMG</th>
2288</tr>
2289<tr>
2290 <th align="right">Procs</th>
2291 <th align="right">Cycle</th>
2292 <th align="right">DoFs</th>
2293 <th align="right">@f$\mathbb{E}@f$</th>
2294 <th></th>
2295 <th align="right">Setup</th>
2296 <th align="right">Prec</th>
2297 <th align="right">Solve</th>
2298 <th align="right">Total</th>
2299 <th></th>
2300 <th align="right">Setup</th>
2301 <th align="right">Prec</th>
2302 <th align="right">Solve</th>
2303 <th align="right">Total</th>
2304 <th></th>
2305 <th align="right">Setup</th>
2306 <th align="right">Prec</th>
2307 <th align="right">Solve</th>
2308 <th align="right">Total</th>
2309</tr>
2310<tr>
2311 <td align="right">112</th>
2312 <td align="right">13</th>
2313 <td align="right">4M</th>
2314 <td align="right">0.37</th>
2315 <td></td>
2316 <td align="right">0.742</th>
2317 <td align="right">0.393</th>
2318 <td align="right">0.200</th>
2319 <td align="right">1.335</th>
2320 <td></td>
2321 <td align="right">1.714</th>
2322 <td align="right">2.934</th>
2323 <td align="right">0.716</th>
2324 <td align="right">5.364</th>
2325 <td></td>
2326 <td align="right">1.544</th>
2327 <td align="right">0.456</th>
2328 <td align="right">1.150</th>
2329 <td align="right">3.150</th>
2330</tr>
2331<tr>
2332 <td align="right">448</th>
2333 <td align="right">15</th>
2334 <td align="right">16M</th>
2335 <td align="right">0.29</th>
2336 <td></td>
2337 <td align="right">0.884</th>
2338 <td align="right">0.535</th>
2339 <td align="right">0.253</th>
2340 <td align="right">1.672</th>
2341 <td></td>
2342 <td align="right">1.927</th>
2343 <td align="right">3.776</th>
2344 <td align="right">1.190</th>
2345 <td align="right">6.893</th>
2346 <td></td>
2347 <td align="right">1.544</th>
2348 <td align="right">0.456</th>
2349 <td align="right">1.150</th>
2350 <td align="right">3.150</th>
2351</tr>
2352<tr>
2353 <td align="right">1,792</th>
2354 <td align="right">17</th>
2355 <td align="right">65M</th>
2356 <td align="right">0.22</th>
2357 <td></td>
2358 <td align="right">1.122</th>
2359 <td align="right">0.686</th>
2360 <td align="right">0.309</th>
2361 <td align="right">2.117</th>
2362 <td></td>
2363 <td align="right">2.171</th>
2364 <td align="right">4.862</th>
2365 <td align="right">1.660</th>
2366 <td align="right">8.693</th>
2367 <td></td>
2368 <td align="right">1.654</th>
2369 <td align="right">0.546</th>
2370 <td align="right">1.460</th>
2371 <td align="right">3.660</th>
2372</tr>
2373<tr>
2374 <td align="right">7,168</th>
2375 <td align="right">19</th>
2376 <td align="right">256M</th>
2377 <td align="right">0.16</th>
2378 <td></td>
2379 <td align="right">1.214</th>
2380 <td align="right">0.893</th>
2381 <td align="right">0.521</th>
2382 <td align="right">2.628</th>
2383 <td></td>
2384 <td align="right">2.386</th>
2385 <td align="right">7.260</th>
2386 <td align="right">2.560</th>
2387 <td align="right">12.206</th>
2388 <td></td>
2389 <td align="right">1.844</th>
2390 <td align="right">1.010</th>
2391 <td align="right">1.890</th>
2392 <td align="right">4.744</th>
2393</tr>
2394</table>
2395
2396On the other hand, the algebraic multigrid in the last set of columns
2397is relatively unaffected by the increasing imbalance of the mesh
2398hierarchy (because it doesn't use the mesh hierarchy) and the growth
2399in time is rather driven by other factors that are well documented in
2400the literature (most notably that the algorithmic complexity of
2401some parts of algebraic multigrid methods appears to be @f${\cal O}(N
2402\log N)@f$ instead of @f${\cal O}(N)@f$ for geometric multigrid).
2403
2404The upshort of the table above is that the matrix-free geometric multigrid
2405method appears to be the fastest approach to solving this equation if
2406not by a huge margin. Matrix-based methods, on the other hand, are
2407consistently the worst.
2408
2409The following figure provides strong scaling results for each method, i.e.,
2410we solve the same problem on more and more processors. Specifically,
2411we consider the problems after 16 mesh refinement cycles
2412(32M DoFs) and 19 cycles (256M DoFs), on between 56 to 28,672 processors:
2413
2414<img width="600px" src="https://www.dealii.org/images/steps/developer/step-50-strong-scaling.png" alt="">
2415
2416While the matrix-based GMG solver and AMG scale similarly and have a
2417similar time to solution (at least as long as there is a substantial
2418number of unknowns per processor -- say, several 10,000), the
2419matrix-free GMG solver scales much better and solves the finer problem
2420in roughly the same time as the AMG solver for the coarser mesh with
2421only an eighth of the number of processors. Conversely, it can solve the
2422same problem on the same number of processors in about one eighth the
2423time.
2424
2425
2426<a name="step_50-Possibilitiesforextensions"></a><h3> Possibilities for extensions </h3>
2427
2428
2429<a name="step_50-Testingconvergenceandhigherorderelements"></a><h4> Testing convergence and higher order elements </h4>
2430
2431
2432The finite element degree is currently hard-coded as 2, see the template
2433arguments of the main class. It is easy to change. To test, it would be
2434interesting to switch to a test problem with a reference solution. This way,
2435you can compare error rates.
2436
2437<a name="step_50-Coarsesolver"></a><h4> Coarse solver </h4>
2438
2439
2440A more interesting example would involve a more complicated coarse mesh (see
2441@ref step_49 "step-49" for inspiration). The issue in that case is that the coarsest
2442level of the mesh hierarchy is actually quite large, and one would
2443have to think about ways to solve the coarse level problem
2444efficiently. (This is not an issue for algebraic multigrid methods
2445because they would just continue to build coarser and coarser levels
2446of the matrix, regardless of their geometric origin.)
2447
2448In the program here, we simply solve the coarse level problem with a
2449Conjugate Gradient method without any preconditioner. That is acceptable
2450if the coarse problem is really small -- for example, if the coarse
2451mesh had a single cell, then the coarse mesh problems has a @f$9\times 9@f$
2452matrix in 2d, and a @f$27\times 27@f$ matrix in 3d; for the coarse mesh we
2453use on the @f$L@f$-shaped domain of the current program, these sizes are
2454@f$21\times 21@f$ in 2d and @f$117\times 117@f$ in 3d. But if the coarse mesh
2455consists of hundreds or thousands of cells, this approach will no
2456longer work and might start to dominate the overall run-time of each V-cycle.
2457A common approach is then to solve the coarse mesh problem using an
2458algebraic multigrid preconditioner; this would then, however, require
2459assembling the coarse matrix (even for the matrix-free version) as
2460input to the AMG implementation.
2461 *
2462 *
2463<a name="step_50-PlainProg"></a>
2464<h1> The plain program</h1>
2465@include "step-50.cc"
2466*/
Definition fe_q.h:554
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation="", const bool has_to_be_set=false)
Definition point.h:111
unsigned int level
Definition grid_out.cc:4626
#define AssertThrow(cond, exc)
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
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
#define DEAL_II_NOT_IMPLEMENTED()
std::vector< value_type > split(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
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 hyper_L(Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false)
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
@ matrix
Contents is actually a matrix.
PETScWrappers::PreconditionBoomerAMG PreconditionAMG
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
Tpetra::CrsMatrix< Number, LO, GO, NodeType< MemorySpace > > MatrixType
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
void make_interface_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs &mg_constrained_dofs, SparsityPatternBase &sparsity, const unsigned int level)
Definition mg_tools.cc:1013
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true)
Definition mg_tools.cc:575
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 > &)> &cell_operation, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0, const unsigned int first_vector_component=0)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
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 this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
std::string compress(const std::string &input)
Definition utilities.cc:389
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 run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &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)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition tria.cc:14871
int(&) functions(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition loop.h:70
Definition types.h:32
unsigned int boundary_id
Definition types.h:144
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
TasksParallelScheme tasks_parallel_scheme