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