Reference documentation for deal.II version 9.6.0
\(\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-59.h
Go to the documentation of this file.
1
812true,
815 *   }
816 *  
817 *  
818 *  
819 * @endcode
820 *
821 * Since the Laplacian is symmetric, the `Tvmult()` (needed by the multigrid
822 * smoother interfaces) operation is simply forwarded to the `vmult()` case.
823 *
824
825 *
826 *
827 * @code
828 *   template <int dim, int fe_degree, typename number>
829 *   void LaplaceOperator<dim, fe_degree, number>::Tvmult(
831 *   const LinearAlgebra::distributed::Vector<number> &src) const
832 *   {
833 *   vmult(dst, src);
834 *   }
835 *  
836 *  
837 *  
838 * @endcode
839 *
840 * The cell operation is very similar to @ref step_37 "step-37". We do not use a
841 * coefficient here, though. The second difference is that we replaced the
842 * two steps of FEEvaluation::read_dof_values() followed by
843 * FEEvaluation::evaluate() by a single function call
844 * FEEvaluation::gather_evaluate() which internally calls the sequence of
845 * the two individual methods. Likewise, FEEvaluation::integrate_scatter()
846 * implements the sequence of FEEvaluation::integrate() followed by
847 * FEEvaluation::distribute_local_to_global(). In this case, these new
848 * functions merely save two lines of code. However, we use them for the
849 * analogy with FEFaceEvaluation where they are more important as
850 * explained below.
851 *
852
853 *
854 *
855 * @code
856 *   template <int dim, int fe_degree, typename number>
857 *   void LaplaceOperator<dim, fe_degree, number>::apply_cell(
858 *   const MatrixFree<dim, number> &data,
859 *   LinearAlgebra::distributed::Vector<number> &dst,
860 *   const LinearAlgebra::distributed::Vector<number> &src,
861 *   const std::pair<unsigned int, unsigned int> &cell_range) const
862 *   {
864 *   for (unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
865 *   {
866 *   phi.reinit(cell);
867 *   phi.gather_evaluate(src, EvaluationFlags::gradients);
868 *   for (const unsigned int q : phi.quadrature_point_indices())
869 *   phi.submit_gradient(phi.get_gradient(q), q);
870 *   phi.integrate_scatter(EvaluationFlags::gradients, dst);
871 *   }
872 *   }
873 *  
874 *  
875 *  
876 * @endcode
877 *
878 * The face operation implements the terms of the interior penalty method in
879 * analogy to @ref step_39 "step-39", as explained in the introduction. We need two
880 * evaluator objects for this task, one for handling the solution that comes
881 * from the cell on one of the two sides of an interior face, and one for
882 * handling the solution from the other side. The evaluators for face
883 * integrals are called FEFaceEvaluation and take a boolean argument in the
884 * second slot of the constructor to indicate which of the two sides the
885 * evaluator should belong two. In FEFaceEvaluation and MatrixFree, we call
886 * one of the two sides the `interior` one and the other the `exterior`
887 * one. The name `exterior` refers to the fact that the evaluator from both
888 * sides will return the same normal vector. For the `interior` side, the
889 * normal vector points outwards, whereas it points inwards on the other
890 * side, and is opposed to the outer normal vector of that cell. Apart from
891 * the new class name, we again get a range of items to work with in
892 * analogy to what was discussed in @ref step_37 "step-37", but for the interior faces in
893 * this case. Note that the data structure of MatrixFree forms batches of
894 * faces that are analogous to the batches of cells for the cell
895 * integrals. All faces within a batch involve different cell numbers but
896 * have the face number within the reference cell, have the same refinement
897 * configuration (no refinement or the same subface), and the same
898 * orientation, to keep SIMD operations simple and efficient.
899 *
900
901 *
902 * Note that there is no implied meaning in interior versus exterior except
903 * the logic decision of the orientation of the normal, which is pretty
904 * random internally. One can in no way rely on a certain pattern of
905 * assigning interior versus exterior flags, as the decision is made for the
906 * sake of access regularity and uniformity in the MatrixFree setup
907 * routines. Since most sane DG methods are conservative, i.e., fluxes look
908 * the same from both sides of an interface, the mathematics are unaltered
909 * if the interior/exterior flags are switched and normal vectors get the
910 * opposite sign.
911 *
912
913 *
914 *
915 * @code
916 *   template <int dim, int fe_degree, typename number>
917 *   void LaplaceOperator<dim, fe_degree, number>::apply_face(
918 *   const MatrixFree<dim, number> &data,
921 *   const std::pair<unsigned int, unsigned int> &face_range) const
922 *   {
924 *   true);
926 *   false);
927 *   for (unsigned int face = face_range.first; face < face_range.second; ++face)
928 *   {
929 * @endcode
930 *
931 * On a given batch of faces, we first update the pointers to the
932 * current face and then access the vector. As mentioned above, we
933 * combine the vector access with the evaluation. In the case of face
934 * integrals, the data access into the vector can be reduced for the
935 * special case of an FE_DGQHermite basis as explained for the data
936 * exchange above: Since only @f$2(k+1)^{d-1}@f$ out of the @f$(k+1)^d@f$ cell
937 * degrees of freedom get multiplied by a non-zero value or derivative
938 * of a shape function, this structure can be utilized for the
939 * evaluation, significantly reducing the data access. The reduction
940 * of the data access is not only beneficial because it reduces the
941 * data in flight and thus helps caching, but also because the data
942 * access to faces is often more irregular than for cell integrals when
943 * gathering values from cells that are farther apart in the index
944 * list of cells.
945 *
946 * @code
947 *   phi_inner.reinit(face);
948 *   phi_inner.gather_evaluate(src,
951 *   phi_outer.reinit(face);
952 *   phi_outer.gather_evaluate(src,
955 *  
956 * @endcode
957 *
958 * The next two statements compute the penalty parameter for the
959 * interior penalty method. As explained in the introduction, we would
960 * like to have a scaling like @f$\frac{1}{h_\text{i}}@f$ of the length
961 * @f$h_\text{i}@f$ normal to the face. For a general non-Cartesian mesh,
962 * this length must be computed by the product of the inverse Jacobian
963 * times the normal vector in real coordinates. From this vector of
964 * `dim` components, we must finally pick the component that is
965 * oriented normal to the reference cell. In the geometry data stored
966 * in MatrixFree, a permutation of the components in the Jacobian is
967 * applied such that this latter direction is always the last
968 * component `dim-1` (this is beneficial because reference-cell
969 * derivative sorting can be made agnostic of the direction of the
970 * face). This means that we can simply access the last entry `dim-1`
971 * and must not look up the local face number in
972 * `data.get_face_info(face).interior_face_no` and
973 * `data.get_face_info(face).exterior_face_no`. Finally, we must also
974 * take the absolute value of these factors as the normal could point
975 * into either positive or negative direction.
976 *
977 * @code
978 *   const VectorizedArray<number> inverse_length_normal_to_face =
979 *   0.5 * (std::abs((phi_inner.normal_vector(0) *
980 *   phi_inner.inverse_jacobian(0))[dim - 1]) +
981 *   std::abs((phi_outer.normal_vector(0) *
982 *   phi_outer.inverse_jacobian(0))[dim - 1]));
983 *   const VectorizedArray<number> sigma =
984 *   inverse_length_normal_to_face * get_penalty_factor();
985 *  
986 * @endcode
987 *
988 * In the loop over the quadrature points, we eventually compute all
989 * contributions to the interior penalty scheme. According to the
990 * formulas in the introduction, the value of the test function gets
991 * multiplied by the difference of the jump in the solution times the
992 * penalty parameter and the average of the normal derivative in real
993 * space. Since the two evaluators for interior and exterior sides get
994 * different signs due to the jump, we pass the result with a
995 * different sign here. The normal derivative of the test function
996 * gets multiplied by the negative jump in the solution between the
997 * interior and exterior side. This term, coined adjoint consistency
998 * term, must also include the factor of @f$\frac{1}{2}@f$ in the code in
999 * accordance with its relation to the primal consistency term that
1000 * gets the factor of one half due to the average in the test function
1001 * slot.
1002 *
1003 * @code
1004 *   for (const unsigned int q : phi_inner.quadrature_point_indices())
1005 *   {
1006 *   const VectorizedArray<number> solution_jump =
1007 *   (phi_inner.get_value(q) - phi_outer.get_value(q));
1008 *   const VectorizedArray<number> average_normal_derivative =
1009 *   (phi_inner.get_normal_derivative(q) +
1010 *   phi_outer.get_normal_derivative(q)) *
1011 *   number(0.5);
1012 *   const VectorizedArray<number> test_by_value =
1013 *   solution_jump * sigma - average_normal_derivative;
1014 *  
1015 *   phi_inner.submit_value(test_by_value, q);
1016 *   phi_outer.submit_value(-test_by_value, q);
1017 *  
1018 *   phi_inner.submit_normal_derivative(-solution_jump * number(0.5), q);
1019 *   phi_outer.submit_normal_derivative(-solution_jump * number(0.5), q);
1020 *   }
1021 *  
1022 * @endcode
1023 *
1024 * Once we are done with the loop over quadrature points, we can do
1025 * the sum factorization operations for the integration loops on faces
1026 * and sum the results into the result vector, using the
1027 * `integrate_scatter` function. The name `scatter` reflects the
1028 * distribution of the vector data into scattered positions in the
1029 * vector using the same pattern as in `gather_evaluate`. Like before,
1030 * the combined integrate + write operation allows us to reduce the
1031 * data access.
1032 *
1033 * @code
1034 *   phi_inner.integrate_scatter(EvaluationFlags::values |
1036 *   dst);
1037 *   phi_outer.integrate_scatter(EvaluationFlags::values |
1039 *   dst);
1040 *   }
1041 *   }
1042 *  
1043 *  
1044 *  
1045 * @endcode
1046 *
1047 * The boundary face function follows by and large the interior face
1048 * function. The only difference is the fact that we do not have a separate
1049 * FEFaceEvaluation object that provides us with exterior values @f$u^+@f$, but
1050 * we must define them from the boundary conditions and interior values
1051 * @f$u^-@f$. As explained in the introduction, we use @f$u^+ = -u^- + 2
1052 * g_\text{D}@f$ and @f$\mathbf{n}^-\cdot \nabla u^+ = \mathbf{n}^-\cdot \nabla
1053 * u^-@f$ on Dirichlet boundaries and @f$u^+=u^-@f$ and @f$\mathbf{n}^-\cdot \nabla
1054 * u^+ = -\mathbf{n}^-\cdot \nabla u^- + 2 g_\text{N}@f$ on Neumann
1055 * boundaries. Since this operation implements the homogeneous part, i.e.,
1056 * the matrix-vector product, we must neglect the boundary functions
1057 * @f$g_\text{D}@f$ and @f$g_\text{N}@f$ here, and added them to the right hand side
1058 * in `LaplaceProblem::compute_rhs()`. Note that due to extension of the
1059 * solution @f$u^-@f$ to the exterior via @f$u^+@f$, we can keep all factors @f$0.5@f$
1060 * the same as in the inner face function, see also the discussion in
1061 * @ref step_39 "step-39".
1062 *
1063
1064 *
1065 * There is one catch at this point: The implementation below uses a boolean
1066 * variable `is_dirichlet` to switch between the Dirichlet and the Neumann
1067 * cases. However, we solve a problem where we also want to impose periodic
1068 * boundary conditions on some boundaries, namely along those in the @f$x@f$
1069 * direction. One might wonder how those conditions should be handled
1070 * here. The answer is that MatrixFree automatically treats periodic
1071 * boundaries as what they are technically, namely an inner face where the
1072 * solution values of two adjacent cells meet and must be treated by proper
1073 * numerical fluxes. Thus, all the faces on the periodic boundaries will
1074 * appear in the `apply_face()` function and not in this one.
1075 *
1076
1077 *
1078 *
1079 * @code
1080 *   template <int dim, int fe_degree, typename number>
1081 *   void LaplaceOperator<dim, fe_degree, number>::apply_boundary(
1082 *   const MatrixFree<dim, number> &data,
1085 *   const std::pair<unsigned int, unsigned int> &face_range) const
1086 *   {
1088 *   true);
1089 *   for (unsigned int face = face_range.first; face < face_range.second; ++face)
1090 *   {
1091 *   phi_inner.reinit(face);
1092 *   phi_inner.gather_evaluate(src,
1095 *  
1096 *   const VectorizedArray<number> inverse_length_normal_to_face = std::abs((
1097 *   phi_inner.normal_vector(0) * phi_inner.inverse_jacobian(0))[dim - 1]);
1098 *   const VectorizedArray<number> sigma =
1099 *   inverse_length_normal_to_face * get_penalty_factor();
1100 *  
1101 *   const bool is_dirichlet = (data.get_boundary_id(face) == 0);
1102 *  
1103 *   for (const unsigned int q : phi_inner.quadrature_point_indices())
1104 *   {
1105 *   const VectorizedArray<number> u_inner = phi_inner.get_value(q);
1106 *   const VectorizedArray<number> u_outer =
1107 *   is_dirichlet ? -u_inner : u_inner;
1108 *   const VectorizedArray<number> normal_derivative_inner =
1109 *   phi_inner.get_normal_derivative(q);
1110 *   const VectorizedArray<number> normal_derivative_outer =
1111 *   is_dirichlet ? normal_derivative_inner : -normal_derivative_inner;
1112 *   const VectorizedArray<number> solution_jump = (u_inner - u_outer);
1113 *   const VectorizedArray<number> average_normal_derivative =
1114 *   (normal_derivative_inner + normal_derivative_outer) * number(0.5);
1115 *   const VectorizedArray<number> test_by_value =
1116 *   solution_jump * sigma - average_normal_derivative;
1117 *   phi_inner.submit_normal_derivative(-solution_jump * number(0.5), q);
1118 *   phi_inner.submit_value(test_by_value, q);
1119 *   }
1120 *   phi_inner.integrate_scatter(EvaluationFlags::values |
1122 *   dst);
1123 *   }
1124 *   }
1125 *  
1126 *  
1127 *  
1128 * @endcode
1129 *
1130 * Next we turn to the preconditioner initialization. As explained in the
1131 * introduction, we want to construct an (approximate) inverse of the cell
1132 * matrices from a product of 1d mass and Laplace matrices. Our first task
1133 * is to compute the 1d matrices, which we do by first creating a 1d finite
1134 * element. Instead of anticipating FE_DGQHermite<1> here, we get the finite
1135 * element's name from DoFHandler, replace the @p dim argument (2 or 3) by 1
1136 * to create a 1d name, and construct the 1d element by using FETools.
1137 *
1138
1139 *
1140 *
1141 * @code
1142 *   template <int dim, int fe_degree, typename number>
1143 *   void PreconditionBlockJacobi<dim, fe_degree, number>::initialize(
1144 *   const LaplaceOperator<dim, fe_degree, number> &op)
1145 *   {
1146 *   data = op.get_matrix_free();
1147 *  
1148 *   std::string name = data->get_dof_handler().get_fe().get_name();
1149 *   name.replace(name.find('<') + 1, 1, "1");
1150 *   std::unique_ptr<FiniteElement<1>> fe_1d = FETools::get_fe_by_name<1>(name);
1151 *  
1152 * @endcode
1153 *
1154 * As for computing the 1d matrices on the unit element, we simply write
1155 * down what a typical assembly procedure over rows and columns of the
1156 * matrix as well as the quadrature points would do. We select the same
1157 * Laplace matrices once and for all using the coefficients 0.5 for
1158 * interior faces (but possibly scaled differently in different directions
1159 * as a result of the mesh). Thus, we make a slight mistake at the
1160 * Dirichlet boundary (where the correct factor would be 1 for the
1161 * derivative terms and 2 for the penalty term, see @ref step_39 "step-39") or at the
1162 * Neumann boundary where the factor should be zero. Since we only use
1163 * this class as a smoother inside a multigrid scheme, this error is not
1164 * going to have any significant effect and merely affects smoothing
1165 * quality.
1166 *
1167 * @code
1168 *   const unsigned int N = fe_degree + 1;
1169 *   FullMatrix<double> laplace_unscaled(N, N);
1170 *   std::array<Table<2, VectorizedArray<number>>, dim> mass_matrices;
1171 *   std::array<Table<2, VectorizedArray<number>>, dim> laplace_matrices;
1172 *   for (unsigned int d = 0; d < dim; ++d)
1173 *   {
1174 *   mass_matrices[d].reinit(N, N);
1175 *   laplace_matrices[d].reinit(N, N);
1176 *   }
1177 *  
1178 *   const QGauss<1> quadrature(N);
1179 *   for (unsigned int i = 0; i < N; ++i)
1180 *   for (unsigned int j = 0; j < N; ++j)
1181 *   {
1182 *   double sum_mass = 0, sum_laplace = 0;
1183 *   for (unsigned int q = 0; q < quadrature.size(); ++q)
1184 *   {
1185 *   sum_mass += (fe_1d->shape_value(i, quadrature.point(q)) *
1186 *   fe_1d->shape_value(j, quadrature.point(q))) *
1187 *   quadrature.weight(q);
1188 *   sum_laplace += (fe_1d->shape_grad(i, quadrature.point(q))[0] *
1189 *   fe_1d->shape_grad(j, quadrature.point(q))[0]) *
1190 *   quadrature.weight(q);
1191 *   }
1192 *   for (unsigned int d = 0; d < dim; ++d)
1193 *   mass_matrices[d](i, j) = sum_mass;
1194 *  
1195 * @endcode
1196 *
1197 * The left and right boundary terms assembled by the next two
1198 * statements appear to have somewhat arbitrary signs, but those are
1199 * correct as can be verified by looking at @ref step_39 "step-39" and inserting
1200 * the value -1 and 1 for the normal vector in the 1d case.
1201 *
1202 * @code
1203 *   sum_laplace +=
1204 *   (1. * fe_1d->shape_value(i, Point<1>()) *
1205 *   fe_1d->shape_value(j, Point<1>()) * op.get_penalty_factor() +
1206 *   0.5 * fe_1d->shape_grad(i, Point<1>())[0] *
1207 *   fe_1d->shape_value(j, Point<1>()) +
1208 *   0.5 * fe_1d->shape_grad(j, Point<1>())[0] *
1209 *   fe_1d->shape_value(i, Point<1>()));
1210 *  
1211 *   sum_laplace +=
1212 *   (1. * fe_1d->shape_value(i, Point<1>(1.0)) *
1213 *   fe_1d->shape_value(j, Point<1>(1.0)) * op.get_penalty_factor() -
1214 *   0.5 * fe_1d->shape_grad(i, Point<1>(1.0))[0] *
1215 *   fe_1d->shape_value(j, Point<1>(1.0)) -
1216 *   0.5 * fe_1d->shape_grad(j, Point<1>(1.0))[0] *
1217 *   fe_1d->shape_value(i, Point<1>(1.0)));
1218 *  
1219 *   laplace_unscaled(i, j) = sum_laplace;
1220 *   }
1221 *  
1222 * @endcode
1223 *
1224 * Next, we go through the cells and pass the scaled matrices to
1225 * TensorProductMatrixSymmetricSum to actually compute the generalized
1226 * eigenvalue problem for representing the inverse: Since the matrix
1227 * approximation is constructed as @f$A\otimes M + M\otimes A@f$ and the
1228 * weights are constant for each element, we can apply all weights on the
1229 * Laplace matrix and simply keep the mass matrices unscaled. In the loop
1230 * over cells, we want to make use of the geometry compression provided by
1231 * the MatrixFree class and check if the current geometry is the same as
1232 * on the last cell batch, in which case there is nothing to do. This
1233 * compression can be accessed by
1234 * FEEvaluation::get_mapping_data_index_offset() once `reinit()` has been
1235 * called.
1236 *
1237
1238 *
1239 * Once we have accessed the inverse Jacobian through the FEEvaluation
1240 * access function (we take the one for the zeroth quadrature point as
1241 * they should be the same on all quadrature points for a Cartesian cell),
1242 * we check that it is diagonal and then extract the determinant of the
1243 * original Jacobian, i.e., the inverse of the determinant of the inverse
1244 * Jacobian, and set the weight as @f$\text{det}(J) / h_d^2@f$ according to
1245 * the 1d Laplacian times @f$d-1@f$ copies of the @ref GlossMassMatrix "mass matrix".
1246 *
1247 * @code
1248 *   cell_matrices.clear();
1249 *   FEEvaluation<dim, fe_degree, fe_degree + 1, 1, number> phi(*data);
1250 *   unsigned int old_mapping_data_index = numbers::invalid_unsigned_int;
1251 *   for (unsigned int cell = 0; cell < data->n_cell_batches(); ++cell)
1252 *   {
1253 *   phi.reinit(cell);
1254 *  
1255 *   if (phi.get_mapping_data_index_offset() == old_mapping_data_index)
1256 *   continue;
1257 *  
1258 *   Tensor<2, dim, VectorizedArray<number>> inverse_jacobian =
1259 *   phi.inverse_jacobian(0);
1260 *  
1261 *   for (unsigned int d = 0; d < dim; ++d)
1262 *   for (unsigned int e = 0; e < dim; ++e)
1263 *   if (d != e)
1264 *   for (unsigned int v = 0; v < VectorizedArray<number>::size(); ++v)
1265 *   AssertThrow(inverse_jacobian[d][e][v] == 0.,
1266 *   ExcNotImplemented());
1267 *  
1268 *   VectorizedArray<number> jacobian_determinant = inverse_jacobian[0][0];
1269 *   for (unsigned int e = 1; e < dim; ++e)
1270 *   jacobian_determinant *= inverse_jacobian[e][e];
1271 *   jacobian_determinant = 1. / jacobian_determinant;
1272 *  
1273 *   for (unsigned int d = 0; d < dim; ++d)
1274 *   {
1275 *   const VectorizedArray<number> scaling_factor =
1276 *   inverse_jacobian[d][d] * inverse_jacobian[d][d] *
1277 *   jacobian_determinant;
1278 *  
1279 * @endcode
1280 *
1281 * Once we know the factor by which we should scale the Laplace
1282 * matrix, we apply this weight to the unscaled DG Laplace matrix
1283 * and send the array to the class TensorProductMatrixSymmetricSum
1284 * for computing the generalized eigenvalue problem mentioned in
1285 * the introduction.
1286 *
1287
1288 *
1289 *
1290 * @code
1291 *   for (unsigned int i = 0; i < N; ++i)
1292 *   for (unsigned int j = 0; j < N; ++j)
1293 *   laplace_matrices[d](i, j) =
1294 *   scaling_factor * laplace_unscaled(i, j);
1295 *   }
1296 *   if (cell_matrices.size() <= phi.get_mapping_data_index_offset())
1297 *   cell_matrices.resize(phi.get_mapping_data_index_offset() + 1);
1298 *   cell_matrices[phi.get_mapping_data_index_offset()].reinit(
1299 *   mass_matrices, laplace_matrices);
1300 *   }
1301 *   }
1302 *  
1303 *  
1304 *  
1305 * @endcode
1306 *
1307 * The vmult function for the approximate block-Jacobi preconditioner is
1308 * very simple in the DG context: We simply need to read the values of the
1309 * current cell batch, apply the inverse for the given entry in the array of
1310 * tensor product matrix, and write the result back. In this loop, we
1311 * overwrite the content in `dst` rather than first setting the entries to
1312 * zero. This is legitimate for a DG method because every cell has
1313 * independent degrees of freedom. Furthermore, we manually write out the
1314 * loop over all cell batches, rather than going through
1315 * MatrixFree::cell_loop(). We do this because we know that we are not going
1316 * to need data exchange over the MPI network here as all computations are
1317 * done on the cells held locally on each processor.
1318 *
1319
1320 *
1321 *
1322 * @code
1323 *   template <int dim, int fe_degree, typename number>
1324 *   void PreconditionBlockJacobi<dim, fe_degree, number>::vmult(
1325 *   LinearAlgebra::distributed::Vector<number> &dst,
1326 *   const LinearAlgebra::distributed::Vector<number> &src) const
1327 *   {
1328 *   adjust_ghost_range_if_necessary(*data, dst);
1329 *   adjust_ghost_range_if_necessary(*data, src);
1330 *  
1331 *   FEEvaluation<dim, fe_degree, fe_degree + 1, 1, number> phi(*data);
1332 *   for (unsigned int cell = 0; cell < data->n_cell_batches(); ++cell)
1333 *   {
1334 *   phi.reinit(cell);
1335 *   phi.read_dof_values(src);
1336 *   cell_matrices[phi.get_mapping_data_index_offset()].apply_inverse(
1337 *   ArrayView<VectorizedArray<number>>(phi.begin_dof_values(),
1338 *   phi.dofs_per_cell),
1339 *   ArrayView<const VectorizedArray<number>>(phi.begin_dof_values(),
1340 *   phi.dofs_per_cell));
1341 *   phi.set_dof_values(dst);
1342 *   }
1343 *   }
1344 *  
1345 *  
1346 *  
1347 * @endcode
1348 *
1349 * The definition of the LaplaceProblem class is very similar to
1350 * @ref step_37 "step-37". One difference is the fact that we add the element degree as a
1351 * template argument to the class, which would allow us to more easily
1352 * include more than one degree in the same program by creating different
1353 * instances in the `main()` function. The second difference is the
1354 * selection of the element, FE_DGQHermite, which is specialized for this
1355 * kind of equations.
1356 *
1357
1358 *
1359 *
1360 * @code
1361 *   template <int dim, int fe_degree>
1362 *   class LaplaceProblem
1363 *   {
1364 *   public:
1365 *   LaplaceProblem();
1366 *   void run();
1367 *  
1368 *   private:
1369 *   void setup_system();
1370 *   void compute_rhs();
1371 *   void solve();
1372 *   void analyze_results() const;
1373 *  
1374 *   #ifdef DEAL_II_WITH_P4EST
1375 *   parallel::distributed::Triangulation<dim> triangulation;
1376 *   #else
1377 *   Triangulation<dim> triangulation;
1378 *   #endif
1379 *  
1380 *   const FE_DGQHermite<dim> fe;
1381 *   DoFHandler<dim> dof_handler;
1382 *  
1383 *   MappingQ1<dim> mapping;
1384 *  
1385 *   using SystemMatrixType = LaplaceOperator<dim, fe_degree, double>;
1386 *   SystemMatrixType system_matrix;
1387 *  
1388 *   using LevelMatrixType = LaplaceOperator<dim, fe_degree, float>;
1389 *   MGLevelObject<LevelMatrixType> mg_matrices;
1390 *  
1391 *   LinearAlgebra::distributed::Vector<double> solution;
1392 *   LinearAlgebra::distributed::Vector<double> system_rhs;
1393 *  
1394 *   double setup_time;
1395 *   ConditionalOStream pcout;
1396 *   ConditionalOStream time_details;
1397 *   };
1398 *  
1399 *  
1400 *  
1401 *   template <int dim, int fe_degree>
1402 *   LaplaceProblem<dim, fe_degree>::LaplaceProblem()
1403 *   #ifdef DEAL_II_WITH_P4EST
1404 *   : triangulation(MPI_COMM_WORLD,
1405 *   Triangulation<dim>::limit_level_difference_at_vertices,
1406 *   parallel::distributed::Triangulation<
1407 *   dim>::construct_multigrid_hierarchy)
1408 *   #else
1409 *   : triangulation(Triangulation<dim>::limit_level_difference_at_vertices)
1410 *   #endif
1411 *   , fe(fe_degree)
1412 *   , dof_handler(triangulation)
1413 *   , setup_time(0.)
1414 *   , pcout(std::cout, Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
1415 *   , time_details(std::cout,
1416 *   false &&
1417 *   Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
1418 *   {}
1419 *  
1420 *  
1421 *  
1422 * @endcode
1423 *
1424 * The setup function differs in two aspects from @ref step_37 "step-37". The first is that
1425 * we do not need to interpolate any constraints for the discontinuous
1426 * ansatz space, and simply pass a dummy AffineConstraints object into
1427 * Matrixfree::reinit(). The second change arises because we need to tell
1428 * MatrixFree to also initialize the data structures for faces. We do this
1429 * by setting update flags for the inner and boundary faces,
1430 * respectively. On the boundary faces, we need both the function values,
1431 * their gradients, JxW values (for integration), the normal vectors, and
1432 * quadrature points (for the evaluation of the boundary conditions),
1433 * whereas we only need shape function values, gradients, JxW values, and
1434 * normal vectors for interior faces. The face data structures in MatrixFree
1435 * are always built as soon as one of `mapping_update_flags_inner_faces` or
1436 * `mapping_update_flags_boundary_faces` are different from the default
1437 * value `update_default` of UpdateFlags.
1438 *
1439
1440 *
1441 *
1442 * @code
1443 *   template <int dim, int fe_degree>
1444 *   void LaplaceProblem<dim, fe_degree>::setup_system()
1445 *   {
1446 *   Timer time;
1447 *   setup_time = 0;
1448 *  
1449 *   system_matrix.clear();
1450 *   mg_matrices.clear_elements();
1451 *  
1452 *   dof_handler.distribute_dofs(fe);
1453 *   dof_handler.distribute_mg_dofs();
1454 *  
1455 *   pcout << "Number of degrees of freedom: " << dof_handler.n_dofs()
1456 *   << std::endl;
1457 *  
1458 *   setup_time += time.wall_time();
1459 *   time_details << "Distribute DoFs " << time.wall_time() << " s"
1460 *   << std::endl;
1461 *   time.restart();
1462 *  
1463 *   AffineConstraints<double> dummy;
1464 *   dummy.close();
1465 *  
1466 *   {
1467 *   typename MatrixFree<dim, double>::AdditionalData additional_data;
1468 *   additional_data.tasks_parallel_scheme =
1469 *   MatrixFree<dim, double>::AdditionalData::none;
1470 *   additional_data.mapping_update_flags =
1471 *   (update_gradients | update_JxW_values | update_quadrature_points);
1472 *   additional_data.mapping_update_flags_inner_faces =
1473 *   (update_gradients | update_JxW_values | update_normal_vectors);
1474 *   additional_data.mapping_update_flags_boundary_faces =
1475 *   (update_gradients | update_JxW_values | update_normal_vectors |
1476 *   update_quadrature_points);
1477 *   const auto system_mf_storage =
1478 *   std::make_shared<MatrixFree<dim, double>>();
1479 *   system_mf_storage->reinit(
1480 *   mapping, dof_handler, dummy, QGauss<1>(fe.degree + 1), additional_data);
1481 *   system_matrix.initialize(system_mf_storage);
1482 *   }
1483 *  
1484 *   system_matrix.initialize_dof_vector(solution);
1485 *   system_matrix.initialize_dof_vector(system_rhs);
1486 *  
1487 *   setup_time += time.wall_time();
1488 *   time_details << "Setup matrix-free system " << time.wall_time() << " s"
1489 *   << std::endl;
1490 *   time.restart();
1491 *  
1492 *   const unsigned int nlevels = triangulation.n_global_levels();
1493 *   mg_matrices.resize(0, nlevels - 1);
1494 *  
1495 *   for (unsigned int level = 0; level < nlevels; ++level)
1496 *   {
1497 *   typename MatrixFree<dim, float>::AdditionalData additional_data;
1498 *   additional_data.tasks_parallel_scheme =
1499 *   MatrixFree<dim, float>::AdditionalData::none;
1500 *   additional_data.mapping_update_flags =
1501 *   (update_gradients | update_JxW_values);
1502 *   additional_data.mapping_update_flags_inner_faces =
1503 *   (update_gradients | update_JxW_values);
1504 *   additional_data.mapping_update_flags_boundary_faces =
1505 *   (update_gradients | update_JxW_values);
1506 *   additional_data.mg_level = level;
1507 *   const auto mg_mf_storage_level =
1508 *   std::make_shared<MatrixFree<dim, float>>();
1509 *   mg_mf_storage_level->reinit(mapping,
1510 *   dof_handler,
1511 *   dummy,
1512 *   QGauss<1>(fe.degree + 1),
1513 *   additional_data);
1514 *  
1515 *   mg_matrices[level].initialize(mg_mf_storage_level);
1516 *   }
1517 *   setup_time += time.wall_time();
1518 *   time_details << "Setup matrix-free levels " << time.wall_time() << " s"
1519 *   << std::endl;
1520 *   }
1521 *  
1522 *  
1523 *  
1524 * @endcode
1525 *
1526 * The computation of the right hand side is a bit more complicated than in
1527 * @ref step_37 "step-37". The cell term now consists of the negative Laplacian of the
1528 * analytical solution, `RightHandSide`, for which we need to first split up
1529 * the Point of VectorizedArray fields, i.e., a batch of points, into a
1530 * single point by evaluating all lanes in the VectorizedArray
1531 * separately. Remember that the number of lanes depends on the hardware; it
1532 * could be 1 for systems that do not offer vectorization (or where deal.II
1533 * does not have intrinsics), but it could also be 8 or 16 on AVX-512 of
1534 * recent Intel architectures.
1535 *
1536 * @code
1537 *   template <int dim, int fe_degree>
1538 *   void LaplaceProblem<dim, fe_degree>::compute_rhs()
1539 *   {
1540 *   Timer time;
1541 *   system_rhs = 0;
1542 *   const MatrixFree<dim, double> &data = *system_matrix.get_matrix_free();
1543 *   FEEvaluation<dim, fe_degree> phi(data);
1544 *   RightHandSide<dim> rhs_func;
1545 *   Solution<dim> exact_solution;
1546 *   for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1547 *   {
1548 *   phi.reinit(cell);
1549 *   for (const unsigned int q : phi.quadrature_point_indices())
1550 *   {
1551 *   VectorizedArray<double> rhs_val = VectorizedArray<double>();
1552 *   Point<dim, VectorizedArray<double>> point_batch =
1553 *   phi.quadrature_point(q);
1554 *   for (unsigned int v = 0; v < VectorizedArray<double>::size(); ++v)
1555 *   {
1556 *   Point<dim> single_point;
1557 *   for (unsigned int d = 0; d < dim; ++d)
1558 *   single_point[d] = point_batch[d][v];
1559 *   rhs_val[v] = rhs_func.value(single_point);
1560 *   }
1561 *   phi.submit_value(rhs_val, q);
1562 *   }
1563 *   phi.integrate_scatter(EvaluationFlags::values, system_rhs);
1564 *   }
1565 *  
1566 * @endcode
1567 *
1568 * Secondly, we also need to apply the Dirichlet and Neumann boundary
1569 * conditions. This function is the missing part of to the function
1570 * `LaplaceOperator::apply_boundary()` function once the exterior solution
1571 * values @f$u^+ = -u^- + 2 g_\text{D}@f$ and @f$\mathbf{n}^-\cdot \nabla u^+ =
1572 * \mathbf{n}^-\cdot \nabla u^-@f$ on Dirichlet boundaries and @f$u^+=u^-@f$ and
1573 * @f$\mathbf{n}^-\cdot \nabla u^+ = -\mathbf{n}^-\cdot \nabla u^- + 2
1574 * g_\text{N}@f$ on Neumann boundaries are inserted and expanded in terms of
1575 * the boundary functions @f$g_\text{D}@f$ and @f$g_\text{N}@f$. One thing to
1576 * remember is that we move the boundary conditions to the right hand
1577 * side, so the sign is the opposite from what we imposed on the solution
1578 * part.
1579 *
1580
1581 *
1582 * We could have issued both the cell and the boundary part through a
1583 * MatrixFree::loop part, but we choose to manually write the full loop
1584 * over all faces to learn how the index layout of face indices is set up
1585 * in MatrixFree: Both the inner faces and the boundary faces share the
1586 * index range, and all batches of inner faces have lower numbers than the
1587 * batches of boundary cells. A single index for both variants allows us
1588 * to easily use the same data structure FEFaceEvaluation for both cases
1589 * that attaches to the same data field, just at different positions. The
1590 * number of inner face batches (where a batch is due to the combination
1591 * of several faces into one for vectorization) is given by
1592 * MatrixFree::n_inner_face_batches(), whereas the number of boundary face
1593 * batches is given by MatrixFree::n_boundary_face_batches().
1594 *
1595 * @code
1596 *   FEFaceEvaluation<dim, fe_degree> phi_face(data, true);
1597 *   for (unsigned int face = data.n_inner_face_batches();
1598 *   face < data.n_inner_face_batches() + data.n_boundary_face_batches();
1599 *   ++face)
1600 *   {
1601 *   phi_face.reinit(face);
1602 *  
1603 *   const VectorizedArray<double> inverse_length_normal_to_face = std::abs(
1604 *   (phi_face.normal_vector(0) * phi_face.inverse_jacobian(0))[dim - 1]);
1605 *   const VectorizedArray<double> sigma =
1606 *   inverse_length_normal_to_face * system_matrix.get_penalty_factor();
1607 *  
1608 *   for (const unsigned int q : phi_face.quadrature_point_indices())
1609 *   {
1610 *   VectorizedArray<double> test_value = VectorizedArray<double>(),
1611 *   test_normal_derivative =
1612 *   VectorizedArray<double>();
1613 *   Point<dim, VectorizedArray<double>> point_batch =
1614 *   phi_face.quadrature_point(q);
1615 *  
1616 *   for (unsigned int v = 0; v < VectorizedArray<double>::size(); ++v)
1617 *   {
1618 *   Point<dim> single_point;
1619 *   for (unsigned int d = 0; d < dim; ++d)
1620 *   single_point[d] = point_batch[d][v];
1621 *  
1622 * @endcode
1623 *
1624 * The MatrixFree class lets us query the boundary_id of the
1625 * current face batch. Remember that MatrixFree sets up the
1626 * batches for vectorization such that all faces within a
1627 * batch have the same properties, which includes their
1628 * `boundary_id`. Thus, we can query that id here for the
1629 * current face index `face` and either impose the Dirichlet
1630 * case (where we add something to the function value) or the
1631 * Neumann case (where we add something to the normal
1632 * derivative).
1633 *
1634 * @code
1635 *   if (data.get_boundary_id(face) == 0)
1636 *   test_value[v] = 2.0 * exact_solution.value(single_point);
1637 *   else
1638 *   {
1639 *   Tensor<1, dim> normal;
1640 *   for (unsigned int d = 0; d < dim; ++d)
1641 *   normal[d] = phi_face.normal_vector(q)[d][v];
1642 *   test_normal_derivative[v] =
1643 *   -normal * exact_solution.gradient(single_point);
1644 *   }
1645 *   }
1646 *   phi_face.submit_value(test_value * sigma - test_normal_derivative,
1647 *   q);
1648 *   phi_face.submit_normal_derivative(-0.5 * test_value, q);
1649 *   }
1650 *   phi_face.integrate_scatter(EvaluationFlags::values |
1651 *   EvaluationFlags::gradients,
1652 *   system_rhs);
1653 *   }
1654 *  
1655 * @endcode
1656 *
1657 * Since we have manually run the loop over cells rather than using
1658 * MatrixFree::loop(), we must not forget to perform the data exchange
1659 * with MPI - or actually, we would not need that for DG elements here
1660 * because each cell carries its own degrees of freedom and cell and
1661 * boundary integrals only evaluate quantities on the locally owned
1662 * cells. The coupling to neighboring subdomain only comes in by the inner
1663 * face integrals, which we have not done here. That said, it does not
1664 * hurt to call this function here, so we do it as a reminder of what
1665 * happens inside MatrixFree::loop().
1666 *
1667 * @code
1668 *   system_rhs.compress(VectorOperation::add);
1669 *   setup_time += time.wall_time();
1670 *   time_details << "Compute right hand side " << time.wall_time()
1671 *   << " s\n";
1672 *   }
1673 *  
1674 *  
1675 *  
1676 * @endcode
1677 *
1678 * The `solve()` function is copied almost verbatim from @ref step_37 "step-37". We set up
1679 * the same multigrid ingredients, namely the level transfer, a smoother,
1680 * and a coarse grid solver. The only difference is the fact that we do not
1681 * use the diagonal of the Laplacian for the preconditioner of the Chebyshev
1682 * iteration used for smoothing, but instead our newly resolved class
1683 * `%PreconditionBlockJacobi`. The mechanisms are the same, though.
1684 *
1685 * @code
1686 *   template <int dim, int fe_degree>
1687 *   void LaplaceProblem<dim, fe_degree>::solve()
1688 *   {
1689 *   Timer time;
1690 *   MGTransferMatrixFree<dim, float> mg_transfer;
1691 *   mg_transfer.build(dof_handler);
1692 *   setup_time += time.wall_time();
1693 *   time_details << "MG build transfer time " << time.wall_time()
1694 *   << " s\n";
1695 *   time.restart();
1696 *  
1697 *   using SmootherType =
1698 *   PreconditionChebyshev<LevelMatrixType,
1699 *   LinearAlgebra::distributed::Vector<float>,
1700 *   PreconditionBlockJacobi<dim, fe_degree, float>>;
1701 *   mg::SmootherRelaxation<SmootherType,
1702 *   LinearAlgebra::distributed::Vector<float>>
1703 *   mg_smoother;
1704 *   MGLevelObject<typename SmootherType::AdditionalData> smoother_data;
1705 *   smoother_data.resize(0, triangulation.n_global_levels() - 1);
1706 *   for (unsigned int level = 0; level < triangulation.n_global_levels();
1707 *   ++level)
1708 *   {
1709 *   if (level > 0)
1710 *   {
1711 *   smoother_data[level].smoothing_range = 15.;
1712 *   smoother_data[level].degree = 3;
1713 *   smoother_data[level].eig_cg_n_iterations = 10;
1714 *   }
1715 *   else
1716 *   {
1717 *   smoother_data[0].smoothing_range = 2e-2;
1718 *   smoother_data[0].degree = numbers::invalid_unsigned_int;
1719 *   smoother_data[0].eig_cg_n_iterations = mg_matrices[0].m();
1720 *   }
1721 *   smoother_data[level].preconditioner =
1722 *   std::make_shared<PreconditionBlockJacobi<dim, fe_degree, float>>();
1723 *   smoother_data[level].preconditioner->initialize(mg_matrices[level]);
1724 *   }
1725 *   mg_smoother.initialize(mg_matrices, smoother_data);
1726 *  
1727 *   MGCoarseGridApplySmoother<LinearAlgebra::distributed::Vector<float>>
1728 *   mg_coarse;
1729 *   mg_coarse.initialize(mg_smoother);
1730 *  
1731 *   mg::Matrix<LinearAlgebra::distributed::Vector<float>> mg_matrix(
1732 *   mg_matrices);
1733 *  
1734 *   Multigrid<LinearAlgebra::distributed::Vector<float>> mg(
1735 *   mg_matrix, mg_coarse, mg_transfer, mg_smoother, mg_smoother);
1736 *  
1737 *   PreconditionMG<dim,
1738 *   LinearAlgebra::distributed::Vector<float>,
1739 *   MGTransferMatrixFree<dim, float>>
1740 *   preconditioner(dof_handler, mg, mg_transfer);
1741 *  
1742 *   SolverControl solver_control(10000, 1e-12 * system_rhs.l2_norm());
1743 *   SolverCG<LinearAlgebra::distributed::Vector<double>> cg(solver_control);
1744 *   setup_time += time.wall_time();
1745 *   time_details << "MG build smoother time " << time.wall_time()
1746 *   << "s\n";
1747 *   pcout << "Total setup time " << setup_time << " s\n";
1748 *  
1749 *   time.reset();
1750 *   time.start();
1751 *   cg.solve(system_matrix, solution, system_rhs, preconditioner);
1752 *  
1753 *   pcout << "Time solve (" << solver_control.last_step() << " iterations) "
1754 *   << time.wall_time() << " s" << std::endl;
1755 *   }
1756 *  
1757 *  
1758 *  
1759 * @endcode
1760 *
1761 * Since we have solved a problem with analytic solution, we want to verify
1762 * the correctness of our implementation by computing the L2 error of the
1763 * numerical result against the analytic solution.
1764 *
1765
1766 *
1767 *
1768 * @code
1769 *   template <int dim, int fe_degree>
1770 *   void LaplaceProblem<dim, fe_degree>::analyze_results() const
1771 *   {
1772 *   Vector<float> error_per_cell(triangulation.n_active_cells());
1773 *   VectorTools::integrate_difference(mapping,
1774 *   dof_handler,
1775 *   solution,
1776 *   Solution<dim>(),
1777 *   error_per_cell,
1778 *   QGauss<dim>(fe.degree + 2),
1779 *   VectorTools::L2_norm);
1780 *   pcout << "Verification via L2 error: "
1781 *   << std::sqrt(
1782 *   Utilities::MPI::sum(error_per_cell.norm_sqr(), MPI_COMM_WORLD))
1783 *   << std::endl;
1784 *   }
1785 *  
1786 *  
1787 *  
1788 * @endcode
1789 *
1790 * The `run()` function sets up the initial grid and then runs the multigrid
1791 * program in the usual way. As a domain, we choose a rectangle with
1792 * periodic boundary conditions in the @f$x@f$-direction, a Dirichlet condition
1793 * on the front face in @f$y@f$ direction (i.e., the face with index number 2,
1794 * with boundary id equal to 0), and Neumann conditions on the back face as
1795 * well as the two faces in @f$z@f$ direction for the 3d case (with boundary id
1796 * equal to 1). The extent of the domain is a bit different in the @f$x@f$
1797 * direction (where we want to achieve a periodic solution given the
1798 * definition of `Solution`) as compared to the @f$y@f$ and @f$z@f$ directions.
1799 *
1800
1801 *
1802 *
1803 * @code
1804 *   template <int dim, int fe_degree>
1805 *   void LaplaceProblem<dim, fe_degree>::run()
1806 *   {
1807 *   const unsigned int n_ranks =
1808 *   Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
1809 *   pcout << "Running with " << n_ranks << " MPI process"
1810 *   << (n_ranks > 1 ? "es" : "") << ", element " << fe.get_name()
1811 *   << std::endl
1812 *   << std::endl;
1813 *   for (unsigned int cycle = 0; cycle < 9 - dim; ++cycle)
1814 *   {
1815 *   pcout << "Cycle " << cycle << std::endl;
1816 *  
1817 *   if (cycle == 0)
1818 *   {
1819 *   Point<dim> upper_right;
1820 *   upper_right[0] = 2.5;
1821 *   for (unsigned int d = 1; d < dim; ++d)
1822 *   upper_right[d] = 2.8;
1823 *   GridGenerator::hyper_rectangle(triangulation,
1824 *   Point<dim>(),
1825 *   upper_right);
1826 *   triangulation.begin_active()->face(0)->set_boundary_id(10);
1827 *   triangulation.begin_active()->face(1)->set_boundary_id(11);
1828 *   triangulation.begin_active()->face(2)->set_boundary_id(0);
1829 *   for (unsigned int f = 3;
1830 *   f < triangulation.begin_active()->n_faces();
1831 *   ++f)
1832 *   triangulation.begin_active()->face(f)->set_boundary_id(1);
1833 *  
1834 *   std::vector<GridTools::PeriodicFacePair<
1835 *   typename Triangulation<dim>::cell_iterator>>
1836 *   periodic_faces;
1837 *   GridTools::collect_periodic_faces(
1838 *   triangulation, 10, 11, 0, periodic_faces);
1839 *   triangulation.add_periodicity(periodic_faces);
1840 *  
1841 *   triangulation.refine_global(6 - 2 * dim);
1842 *   }
1843 *   triangulation.refine_global(1);
1844 *   setup_system();
1845 *   compute_rhs();
1846 *   solve();
1847 *   analyze_results();
1848 *   pcout << std::endl;
1849 *   };
1850 *   }
1851 *   } // namespace Step59
1852 *  
1853 *  
1854 *  
1855 * @endcode
1856 *
1857 * There is nothing unexpected in the `main()` function. We call `MPI_Init()`
1858 * through the `MPI_InitFinalize` class, pass on the two parameters on the
1859 * dimension and the degree set at the top of the file, and run the Laplace
1860 * problem.
1861 *
1862
1863 *
1864 *
1865 * @code
1866 *   int main(int argc, char *argv[])
1867 *   {
1868 *   try
1869 *   {
1870 *   using namespace Step59;
1871 *  
1872 *   Utilities::MPI::MPI_InitFinalize mpi_init(argc, argv, 1);
1873 *  
1874 *   LaplaceProblem<dimension, degree_finite_element> laplace_problem;
1875 *   laplace_problem.run();
1876 *   }
1877 *   catch (std::exception &exc)
1878 *   {
1879 *   std::cerr << std::endl
1880 *   << std::endl
1881 *   << "----------------------------------------------------"
1882 *   << std::endl;
1883 *   std::cerr << "Exception on processing: " << std::endl
1884 *   << exc.what() << std::endl
1885 *   << "Aborting!" << std::endl
1886 *   << "----------------------------------------------------"
1887 *   << std::endl;
1888 *   return 1;
1889 *   }
1890 *   catch (...)
1891 *   {
1892 *   std::cerr << std::endl
1893 *   << std::endl
1894 *   << "----------------------------------------------------"
1895 *   << std::endl;
1896 *   std::cerr << "Unknown exception!" << std::endl
1897 *   << "Aborting!" << std::endl
1898 *   << "----------------------------------------------------"
1899 *   << std::endl;
1900 *   return 1;
1901 *   }
1902 *  
1903 *   return 0;
1904 *   }
1905 * @endcode
1906<a name="step_59-Results"></a><h1>Results</h1>
1907
1908
1909<a name="step_59-Programoutput"></a><h3>Program output</h3>
1910
1911
1912Like in @ref step_37 "step-37", we evaluate the multigrid solver in terms of run time. In
1913two space dimensions with elements of degree 8, a possible output could look
1914as follows:
1915@code
1916Running with 12 MPI processes, element FE_DGQHermite<2>(8)
1917
1918Cycle 0
1919Number of degrees of freedom: 5184
1920Total setup time 0.0282445 s
1921Time solve (14 iterations) 0.0110712 s
1922Verification via L2 error: 1.66232e-07
1923
1924Cycle 1
1925Number of degrees of freedom: 20736
1926Total setup time 0.0126282 s
1927Time solve (14 iterations) 0.0157021 s
1928Verification via L2 error: 2.91505e-10
1929
1930Cycle 2
1931Number of degrees of freedom: 82944
1932Total setup time 0.0227573 s
1933Time solve (14 iterations) 0.026568 s
1934Verification via L2 error: 6.64514e-13
1935
1936Cycle 3
1937Number of degrees of freedom: 331776
1938Total setup time 0.0604685 s
1939Time solve (14 iterations) 0.0628356 s
1940Verification via L2 error: 5.57513e-13
1941
1942Cycle 4
1943Number of degrees of freedom: 1327104
1944Total setup time 0.154359 s
1945Time solve (13 iterations) 0.219555 s
1946Verification via L2 error: 3.08139e-12
1947
1948Cycle 5
1949Number of degrees of freedom: 5308416
1950Total setup time 0.467764 s
1951Time solve (13 iterations) 1.1821 s
1952Verification via L2 error: 3.90334e-12
1953
1954Cycle 6
1955Number of degrees of freedom: 21233664
1956Total setup time 1.73263 s
1957Time solve (13 iterations) 5.21054 s
1958Verification via L2 error: 4.94543e-12
1959@endcode
1960
1961Like in @ref step_37 "step-37", the number of CG iterations remains constant with increasing
1962problem size. The iteration counts are a bit higher, which is because we use a
1963lower degree of the Chebyshev polynomial (2 vs 5 in @ref step_37 "step-37") and because the
1964interior penalty discretization has a somewhat larger spread in
1965eigenvalues. Nonetheless, 13 iterations to reduce the residual by 12 orders of
1966magnitude, or almost a factor of 9 per iteration, indicates an overall very
1967efficient method. In particular, we can solve a system with 21 million degrees
1968of freedom in 5 seconds when using 12 cores, which is a very good
1969efficiency. Of course, in 2D we are well inside the regime of roundoff for a
1970polynomial degree of 8 &ndash; as a matter of fact, around 83k DoFs or 0.025s
1971would have been enough to fully converge this (simple) analytic solution
1972here.
1973
1974Not much changes if we run the program in three spatial dimensions, except for
1975the fact that we now use do something more useful with the higher polynomial
1976degree and increasing mesh sizes, as the roundoff errors are only obtained at
1977the finest mesh. Still, it is remarkable that we can solve a 3D Laplace
1978problem with a wave of three periods to roundoff accuracy on a twelve-core
1979machine pretty easily - using about 3.5 GB of memory in total for the second
1980to largest case with 24m DoFs, taking not more than eight seconds. The largest
1981case uses 30GB of memory with 191m DoFs.
1982
1983@code
1984Running with 12 MPI processes, element FE_DGQHermite<3>(8)
1985
1986Cycle 0
1987Number of degrees of freedom: 5832
1988Total setup time 0.0210681 s
1989Time solve (15 iterations) 0.0956945 s
1990Verification via L2 error: 0.0297194
1991
1992Cycle 1
1993Number of degrees of freedom: 46656
1994Total setup time 0.0452428 s
1995Time solve (15 iterations) 0.113827 s
1996Verification via L2 error: 9.55733e-05
1997
1998Cycle 2
1999Number of degrees of freedom: 373248
2000Total setup time 0.190423 s
2001Time solve (15 iterations) 0.218309 s
2002Verification via L2 error: 2.6868e-07
2003
2004Cycle 3
2005Number of degrees of freedom: 2985984
2006Total setup time 0.627914 s
2007Time solve (15 iterations) 1.0595 s
2008Verification via L2 error: 4.6918e-10
2009
2010Cycle 4
2011Number of degrees of freedom: 23887872
2012Total setup time 2.85215 s
2013Time solve (15 iterations) 8.30576 s
2014Verification via L2 error: 9.38583e-13
2015
2016Cycle 5
2017Number of degrees of freedom: 191102976
2018Total setup time 16.1324 s
2019Time solve (15 iterations) 65.57 s
2020Verification via L2 error: 3.17875e-13
2021@endcode
2022
2023<a name="step_59-Comparisonofefficiencyatdifferentpolynomialdegrees"></a><h3>Comparison of efficiency at different polynomial degrees</h3>
2024
2025
2026In the introduction and in-code comments, it was mentioned several times that
2027high orders are treated very efficiently with the FEEvaluation and
2028FEFaceEvaluation evaluators. Now, we want to substantiate these claims by
2029looking at the throughput of the 3D multigrid solver for various polynomial
2030degrees. We collect the times as follows: We first run a solver at problem
2031size close to ten million, indicated in the first four table rows, and record
2032the timings. Then, we normalize the throughput by recording the number of
2033million degrees of freedom solved per second (MDoFs/s) to be able to compare
2034the efficiency of the different degrees, which is computed by dividing the
2035number of degrees of freedom by the solver time.
2036
2037<table align="center" class="doxtable">
2038 <tr>
2039 <th>degree</th>
2040 <th>1</th>
2041 <th>2</th>
2042 <th>3</th>
2043 <th>4</th>
2044 <th>5</th>
2045 <th>6</th>
2046 <th>7</th>
2047 <th>8</th>
2048 <th>9</th>
2049 <th>10</th>
2050 <th>11</th>
2051 <th>12</th>
2052 </tr>
2053 <tr>
2054 <th>Number of DoFs</th>
2055 <td>2097152</td>
2056 <td>7077888</td>
2057 <td>16777216</td>
2058 <td>32768000</td>
2059 <td>7077888</td>
2060 <td>11239424</td>
2061 <td>16777216</td>
2062 <td>23887872</td>
2063 <td>32768000</td>
2064 <td>43614208</td>
2065 <td>7077888</td>
2066 <td>8998912</td>
2067 </tr>
2068 <tr>
2069 <th>Number of iterations</th>
2070 <td>13</td>
2071 <td>12</td>
2072 <td>12</td>
2073 <td>12</td>
2074 <td>13</td>
2075 <td>13</td>
2076 <td>15</td>
2077 <td>15</td>
2078 <td>17</td>
2079 <td>19</td>
2080 <td>18</td>
2081 <td>18</td>
2082 </tr>
2083 <tr>
2084 <th>Solver time [s]</th>
2085 <td>0.713</td>
2086 <td>2.150</td>
2087 <td>4.638</td>
2088 <td>8.803</td>
2089 <td>2.041</td>
2090 <td>3.295</td>
2091 <td>5.723</td>
2092 <td>8.306</td>
2093 <td>12.75</td>
2094 <td>19.25</td>
2095 <td>3.530</td>
2096 <td>4.814</td>
2097 </tr>
2098 <tr>
2099 <th>MDoFs/s</th>
2100 <td>2.94</td>
2101 <td>3.29</td>
2102 <td>3.62</td>
2103 <td>3.72</td>
2104 <td>3.47</td>
2105 <td>3.41</td>
2106 <td>2.93</td>
2107 <td>2.88</td>
2108 <td>2.57</td>
2109 <td>2.27</td>
2110 <td>2.01</td>
2111 <td>1.87</td>
2112 </tr>
2113</table>
2114
2115We clearly see how the efficiency per DoF initially improves until it reaches
2116a maximum for the polynomial degree @f$k=4@f$. This effect is surprising, not only
2117because higher polynomial degrees often yield a vastly better solution, but
2118especially also when having matrix-based schemes in mind where the denser
2119coupling at higher degree leads to a monotonously decreasing throughput (and a
2120drastic one in 3D, with @f$k=4@f$ being more than ten times slower than
2121@f$k=1@f$!). For higher degrees, the throughput decreases a bit, which is both due
2122to an increase in the number of iterations (going from 12 at @f$k=2,3,4@f$ to 19
2123at @f$k=10@f$) and due to the @f$\mathcal O(k)@f$ complexity of operator
2124evaluation. Nonetheless, efficiency as the time to solution would be still
2125better for higher polynomial degrees because they have better convergence rates (at least
2126for problems as simple as this one): For @f$k=12@f$, we reach roundoff accuracy
2127already with 1 million DoFs (solver time less than a second), whereas for @f$k=8@f$
2128we need 24 million DoFs and 8 seconds. For @f$k=5@f$, the error is around
2129@f$10^{-9}@f$ with 57m DoFs and thus still far away from roundoff, despite taking 16
2130seconds.
2131
2132Note that the above numbers are a bit pessimistic because they include the
2133time it takes the Chebyshev smoother to compute an eigenvalue estimate, which
2134is around 10 percent of the solver time. If the system is solved several times
2135(as e.g. common in fluid dynamics), this eigenvalue cost is only paid once and
2136faster times become available.
2137
2138<a name="step_59-Evaluationofefficiencyofingredients"></a><h3>Evaluation of efficiency of ingredients</h3>
2139
2140
2141Finally, we take a look at some of the special ingredients presented in this
2142tutorial program, namely the FE_DGQHermite basis in particular and the
2143specification of MatrixFree::DataAccessOnFaces. In the following table, the
2144third row shows the optimized solver above, the fourth row shows the timings
2145with only the MatrixFree::DataAccessOnFaces set to `unspecified` rather than
2146the optimal `gradients`, and the last one with replacing FE_DGQHermite by the
2147basic FE_DGQ elements where both the MPI exchange are more expensive and the
2148operations done by FEFaceEvaluation::gather_evaluate() and
2149FEFaceEvaluation::integrate_scatter().
2150
2151<table align="center" class="doxtable">
2152 <tr>
2153 <th>degree</th>
2154 <th>1</th>
2155 <th>2</th>
2156 <th>3</th>
2157 <th>4</th>
2158 <th>5</th>
2159 <th>6</th>
2160 <th>7</th>
2161 <th>8</th>
2162 <th>9</th>
2163 <th>10</th>
2164 <th>11</th>
2165 <th>12</th>
2166 </tr>
2167 <tr>
2168 <th>Number of DoFs</th>
2169 <td>2097152</td>
2170 <td>7077888</td>
2171 <td>16777216</td>
2172 <td>32768000</td>
2173 <td>7077888</td>
2174 <td>11239424</td>
2175 <td>16777216</td>
2176 <td>23887872</td>
2177 <td>32768000</td>
2178 <td>43614208</td>
2179 <td>7077888</td>
2180 <td>8998912</td>
2181 </tr>
2182 <tr>
2183 <th>Solver time optimized as in tutorial [s]</th>
2184 <td>0.713</td>
2185 <td>2.150</td>
2186 <td>4.638</td>
2187 <td>8.803</td>
2188 <td>2.041</td>
2189 <td>3.295</td>
2190 <td>5.723</td>
2191 <td>8.306</td>
2192 <td>12.75</td>
2193 <td>19.25</td>
2194 <td>3.530</td>
2195 <td>4.814</td>
2196 </tr>
2197 <tr>
2198 <th>Solver time MatrixFree::DataAccessOnFaces::unspecified [s]</th>
2199 <td>0.711</td>
2200 <td>2.151</td>
2201 <td>4.675</td>
2202 <td>8.968</td>
2203 <td>2.243</td>
2204 <td>3.655</td>
2205 <td>6.277</td>
2206 <td>9.082</td>
2207 <td>13.50</td>
2208 <td>20.05</td>
2209 <td>3.817</td>
2210 <td>5.178</td>
2211 </tr>
2212 <tr>
2213 <th>Solver time FE_DGQ [s]</th>
2214 <td>0.712</td>
2215 <td>2.041</td>
2216 <td>5.066</td>
2217 <td>9.335</td>
2218 <td>2.379</td>
2219 <td>3.802</td>
2220 <td>6.564</td>
2221 <td>9.714</td>
2222 <td>14.54</td>
2223 <td>22.76</td>
2224 <td>4.148</td>
2225 <td>5.857</td>
2226 </tr>
2227</table>
2228
2229The data in the table shows that not using MatrixFree::DataAccessOnFaces
2230increases costs by around 10% for higher polynomial degrees. For lower
2231degrees, the difference is obviously less pronounced because the
2232volume-to-surface ratio is more beneficial and less data needs to be
2233exchanged. The difference is larger when looking at the matrix-vector product
2234only, rather than the full multigrid solver shown here, with around 20% worse
2235timings just because of the MPI communication.
2236
2237For @f$k=1@f$ and @f$k=2@f$, the Hermite-like basis functions do obviously not really
2238pay off (indeed, for @f$k=1@f$ the polynomials are exactly the same as for FE_DGQ)
2239and the results are similar as with the FE_DGQ basis. However, for degrees
2240starting at three, we see an increasing advantage for FE_DGQHermite, showing
2241the effectiveness of these basis functions.
2242
2243<a name="step_59-Possibilitiesforextension"></a><h3>Possibilities for extension</h3>
2244
2245
2246As mentioned in the introduction, the fast diagonalization method as realized
2247here is tied to a Cartesian mesh with constant coefficients. When dealing with
2248meshes that contain deformed cells or with variable coefficients, it is common
2249to determine a nearby Cartesian mesh cell as an approximation. This can be
2250done with the class TensorProductMatrixSymmetricSumCollection. Here, one can
2251insert cell matrices similarly to the PreconditionBlockJacobi::initialize()
2252function of this tutorial program. The benefit of the collection class is that
2253cells on which the coefficient of the PDE has the same value can re-use the
2254same Laplacian matrix, which reduces the memory consumption for the inverse
2255matrices. As compared to the algorithm implemented in this tutorial program,
2256one would define the length scales as the distances between opposing
2257faces. For continuous elements, the code project <a
2258href=https://github.com/peterrum/dealii-dd-and-schwarz">Cache-optimized and
2259low-overhead implementations of multigrid smoothers for high-order FEM
2260computations</a> presents the computation for continuous elements. There is
2261currently no infrastructure in deal.II to automatically generate the 1D
2262matrices for discontinuous elements with SIP-DG discretization, as opposed to
2263continuous elements, where we provide
2264TensorProductMatrixCreator::create_laplace_tensor_product_matrix().
2265
2266Another way of extending the program would be to include support for adaptive
2267meshes. While the classical approach of defining interface operations at edges
2268of different refinement level, as discussed in @ref step_39 "step-39", is one possibility,
2269for Poisson-type problems another option is typically more beneficial. Using
2270the class MGTransferGlobalCoarsening, which is explained in the @ref step_75 "step-75"
2271tutorial program, one can deal with meshes of hanging nodes on all levels. An
2272algorithmic improvement can be obtained by combining the discontinuous
2273function space with the auxiliary continuous finite element space of the same
2274polynomial degree. This idea, introduced by Antonietti et al.
2275@cite antonietti2016uniform in 2016, allows making the multigrid convergence
2276independent of the penalty parameter. As demonstrated by Fehn et al.
2277@cite fehn2020hybrid, this also gives considerably lower iteration counts than
2278a multigrid solver directly working on levels with discontinuous function
2279spaces. The latter work also proposes p-multigrid techniques and combination
2280with algebraic multigrid coarse spaces as a means to efficiently solve Poisson
2281problems with high-order discontinuous Galerkin discretizations on complicated
2282geometries, representing the current state-of-the-art for simple Poisson-type
2283problems. The class MGTransferGlobalCoarsening provides features for each of
2284these three coarsening variants, the discontinuous-continuous auxiliary
2285function concept, p-multigrid, and traditional h-multigrid. The main
2286ingredient is to define an appropriate MGTwoLevelTransfer object and call
2287MGTwoLevelTransfer::reinit_geometric_transfer() or
2288MGTwoLevelTransfer::reinit_polynomial_transfer(), respectively.
2289 *
2290 *
2291<a name="step_59-PlainProg"></a>
2292<h1> The plain program</h1>
2293@include "step-59.cc"
2294*/
gradient_type get_gradient(const unsigned int q_point) const
void submit_gradient(const gradient_type grad_in, const unsigned int q_point)
void read_dof_values(const VectorType &src, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip())
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
__global__ void reduction(Number *result, const Number *v, const size_type N)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
Expression sign(const Expression &x)
void random(DoFHandler< dim, spacedim > &dof_handler)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ general
No special properties.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm mpi_communicator)
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
T scatter(const MPI_Comm comm, const std::vector< T > &objects_to_send, const unsigned int root_process=0)
int(& functions)(const void *v1, const void *v2)
STL namespace.
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)