Reference documentation for deal.II version 9.3.2
\(\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\}}\)
step-47.h
Go to the documentation of this file.
1  = 0) const override
670  * {
671  * return std::sin(PI * p[0]) * std::sin(PI * p[1]);
672  * }
673  *
674  * virtual Tensor<1, dim>
675  * gradient(const Point<dim> &p,
676  * const unsigned int /*component*/ = 0) const override
677  * {
678  * Tensor<1, dim> r;
679  * r[0] = PI * std::cos(PI * p[0]) * std::sin(PI * p[1]);
680  * r[1] = PI * std::cos(PI * p[1]) * std::sin(PI * p[0]);
681  * return r;
682  * }
683  *
684  * virtual void
685  * hessian_list(const std::vector<Point<dim>> & points,
686  * std::vector<SymmetricTensor<2, dim>> &hessians,
687  * const unsigned int /*component*/ = 0) const override
688  * {
689  * for (unsigned i = 0; i < points.size(); ++i)
690  * {
691  * const double x = points[i][0];
692  * const double y = points[i][1];
693  *
694  * hessians[i][0][0] = -PI * PI * std::sin(PI * x) * std::sin(PI * y);
695  * hessians[i][0][1] = PI * PI * std::cos(PI * x) * std::cos(PI * y);
696  * hessians[i][1][1] = -PI * PI * std::sin(PI * x) * std::sin(PI * y);
697  * }
698  * }
699  * };
700  *
701  *
702  * template <int dim>
703  * class RightHandSide : public Function<dim>
704  * {
705  * public:
706  * static_assert(dim == 2, "Only dim==2 is implemented");
707  *
708  * virtual double value(const Point<dim> &p,
709  * const unsigned int /*component*/ = 0) const override
710  *
711  * {
712  * return 4 * std::pow(PI, 4.0) * std::sin(PI * p[0]) *
713  * std::sin(PI * p[1]);
714  * }
715  * };
716  * } // namespace ExactSolution
717  *
718  *
719  *
720  * @endcode
721  *
722  *
723  * <a name="Themainclass"></a>
724  * <h3>The main class</h3>
725  *
726 
727  *
728  * The following is the principal class of this tutorial program. It has
729  * the structure of many of the other tutorial programs and there should
730  * really be nothing particularly surprising about its contents or
731  * the constructor that follows it.
732  *
733  * @code
734  * template <int dim>
735  * class BiharmonicProblem
736  * {
737  * public:
738  * BiharmonicProblem(const unsigned int fe_degree);
739  *
740  * void run();
741  *
742  * private:
743  * void make_grid();
744  * void setup_system();
745  * void assemble_system();
746  * void solve();
747  * void compute_errors();
748  * void output_results(const unsigned int iteration) const;
749  *
751  *
752  * MappingQ<dim> mapping;
753  *
754  * FE_Q<dim> fe;
755  * DoFHandler<dim> dof_handler;
756  * AffineConstraints<double> constraints;
757  *
758  * SparsityPattern sparsity_pattern;
759  * SparseMatrix<double> system_matrix;
760  *
761  * Vector<double> solution;
762  * Vector<double> system_rhs;
763  * };
764  *
765  *
766  *
767  * template <int dim>
768  * BiharmonicProblem<dim>::BiharmonicProblem(const unsigned int fe_degree)
769  * : mapping(1)
770  * , fe(fe_degree)
771  * , dof_handler(triangulation)
772  * {}
773  *
774  *
775  *
776  * @endcode
777  *
778  * Next up are the functions that create the initial mesh (a once refined
779  * unit square) and set up the constraints, vectors, and matrices on
780  * each mesh. Again, both of these are essentially unchanged from many
781  * previous tutorial programs.
782  *
783  * @code
784  * template <int dim>
785  * void BiharmonicProblem<dim>::make_grid()
786  * {
788  * triangulation.refine_global(1);
789  *
790  * std::cout << "Number of active cells: " << triangulation.n_active_cells()
791  * << std::endl
792  * << "Total number of cells: " << triangulation.n_cells()
793  * << std::endl;
794  * }
795  *
796  *
797  *
798  * template <int dim>
799  * void BiharmonicProblem<dim>::setup_system()
800  * {
801  * dof_handler.distribute_dofs(fe);
802  *
803  * std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
804  * << std::endl;
805  *
806  * constraints.clear();
807  * DoFTools::make_hanging_node_constraints(dof_handler, constraints);
808  *
810  * 0,
811  * ExactSolution::Solution<dim>(),
812  * constraints);
813  * constraints.close();
814  *
815  *
816  * DynamicSparsityPattern dsp(dof_handler.n_dofs());
817  * DoFTools::make_flux_sparsity_pattern(dof_handler, dsp, constraints, true);
818  * sparsity_pattern.copy_from(dsp);
819  * system_matrix.reinit(sparsity_pattern);
820  *
821  * solution.reinit(dof_handler.n_dofs());
822  * system_rhs.reinit(dof_handler.n_dofs());
823  * }
824  *
825  *
826  *
827  * @endcode
828  *
829  *
830  * <a name="Assemblingthelinearsystem"></a>
831  * <h4>Assembling the linear system</h4>
832  *
833 
834  *
835  * The following pieces of code are more interesting. They all relate to the
836  * assembly of the linear system. While assembling the cell-interior terms
837  * is not of great difficulty -- that works in essence like the assembly
838  * of the corresponding terms of the Laplace equation, and you have seen
839  * how this works in @ref step_4 "step-4" or @ref step_6 "step-6", for example -- the difficulty
840  * is with the penalty terms in the formulation. These require the evaluation
841  * of gradients of shape functions at interfaces of cells. At the least,
842  * one would therefore need to use two FEFaceValues objects, but if one of the
843  * two sides is adaptively refined, then one actually needs an FEFaceValues
844  * and one FESubfaceValues objects; one also needs to keep track which
845  * shape functions live where, and finally we need to ensure that every
846  * face is visited only once. All of this is a substantial overhead to the
847  * logic we really want to implement (namely the penalty terms in the
848  * bilinear form). As a consequence, we will make use of the
849  * FEInterfaceValues class -- a helper class in deal.II that allows us
850  * to abstract away the two FEFaceValues or FESubfaceValues objects and
851  * directly access what we really care about: jumps, averages, etc.
852  *
853 
854  *
855  * But this doesn't yet solve our problem of having to keep track of
856  * which faces we have already visited when we loop over all cells and
857  * all of their faces. To make this process simpler, we use the
858  * MeshWorker::mesh_loop() function that provides a simple interface
859  * for this task: Based on the ideas outlined in the WorkStream
860  * namespace documentation, MeshWorker::mesh_loop() requires three
861  * functions that do work on cells, interior faces, and boundary
862  * faces. These functions work on scratch objects for intermediate
863  * results, and then copy the result of their computations into
864  * copy data objects from where a copier function copies them into
865  * the global matrix and right hand side objects.
866  *
867 
868  *
869  * The following structures then provide the scratch and copy objects
870  * that are necessary for this approach. You may look up the WorkStream
871  * namespace as well as the
872  * @ref threads "Parallel computing with multiple processors"
873  * module for more information on how they typically work.
874  *
875  * @code
876  * template <int dim>
877  * struct ScratchData
878  * {
879  * ScratchData(const Mapping<dim> & mapping,
880  * const FiniteElement<dim> &fe,
881  * const unsigned int quadrature_degree,
882  * const UpdateFlags update_flags,
883  * const UpdateFlags interface_update_flags)
884  * : fe_values(mapping, fe, QGauss<dim>(quadrature_degree), update_flags)
885  * , fe_interface_values(mapping,
886  * fe,
887  * QGauss<dim - 1>(quadrature_degree),
888  * interface_update_flags)
889  * {}
890  *
891  *
892  * ScratchData(const ScratchData<dim> &scratch_data)
893  * : fe_values(scratch_data.fe_values.get_mapping(),
894  * scratch_data.fe_values.get_fe(),
895  * scratch_data.fe_values.get_quadrature(),
896  * scratch_data.fe_values.get_update_flags())
897  * , fe_interface_values(scratch_data.fe_values.get_mapping(),
898  * scratch_data.fe_values.get_fe(),
899  * scratch_data.fe_interface_values.get_quadrature(),
900  * scratch_data.fe_interface_values.get_update_flags())
901  * {}
902  *
903  * FEValues<dim> fe_values;
904  * FEInterfaceValues<dim> fe_interface_values;
905  * };
906  *
907  *
908  *
909  * struct CopyData
910  * {
911  * CopyData(const unsigned int dofs_per_cell)
912  * : cell_matrix(dofs_per_cell, dofs_per_cell)
913  * , cell_rhs(dofs_per_cell)
914  * , local_dof_indices(dofs_per_cell)
915  * {}
916  *
917  *
918  * CopyData(const CopyData &) = default;
919  *
920  *
921  * CopyData(CopyData &&) = default;
922  *
923  *
924  * ~CopyData() = default;
925  *
926  *
927  * CopyData &operator=(const CopyData &) = default;
928  *
929  *
930  * CopyData &operator=(CopyData &&) = default;
931  *
932  *
933  * struct FaceData
934  * {
935  * FullMatrix<double> cell_matrix;
936  * std::vector<types::global_dof_index> joint_dof_indices;
937  * };
938  *
939  * FullMatrix<double> cell_matrix;
940  * Vector<double> cell_rhs;
941  * std::vector<types::global_dof_index> local_dof_indices;
942  * std::vector<FaceData> face_data;
943  * };
944  *
945  *
946  *
947  * @endcode
948  *
949  * The more interesting part is where we actually assemble the linear system.
950  * Fundamentally, this function has five parts:
951  * - The definition of the `cell_worker` lambda function, a small
952  * function that is defined within the `assemble_system()`
953  * function and that will be responsible for computing the local
954  * integrals on an individual cell. It will work on a copy of the
955  * `ScratchData` class and put its results into the corresponding
956  * `CopyData` object.
957  * - The definition of the `face_worker` lambda function that does
958  * the integration of all terms that live on the interfaces between
959  * cells.
960  * - The definition of the `boundary_worker` function that does the
961  * same but for cell faces located on the boundary of the domain.
962  * - The definition of the `copier` function that is responsible
963  * for copying all of the data the previous three functions have
964  * put into copy objects for a single cell, into the global matrix
965  * and right hand side.
966  *
967 
968  *
969  * The fifth part is the one where we bring all of this together.
970  *
971 
972  *
973  * Let us go through each of these pieces necessary for the assembly
974  * in turns.
975  *
976  * @code
977  * template <int dim>
978  * void BiharmonicProblem<dim>::assemble_system()
979  * {
980  * using Iterator = typename DoFHandler<dim>::active_cell_iterator;
981  *
982  * @endcode
983  *
984  * The first piece is the `cell_worker` that does the assembly
985  * on the cell interiors. It is a (lambda) function that takes
986  * a cell (input), a scratch object, and a copy object (output)
987  * as arguments. It looks like the assembly functions of many
988  * other of the tutorial programs, or at least the body of the
989  * loop over all cells.
990  *
991 
992  *
993  * The terms we integrate here are the cell contribution
994  * @f{align*}{
995  * A^K_{ij} = \int_K \nabla^2\varphi_i(x) : \nabla^2\varphi_j(x) dx
996  * @f}
997  * to the global matrix, and
998  * @f{align*}{
999  * f^K_i = \int_K \varphi_i(x) f(x) dx
1000  * @f}
1001  * to the right hand side vector.
1002  *
1003 
1004  *
1005  * We use the same technique as used in the assembly of @ref step_22 "step-22"
1006  * to accelerate the function: Instead of calling
1007  * `fe_values.shape_hessian(i, qpoint)` in the innermost loop,
1008  * we create a variable `hessian_i` that evaluates this
1009  * value once in the loop over `i` and re-use the so-evaluated
1010  * value in the loop over `j`. For symmetry, we do the same with a
1011  * variable `hessian_j`, although it is indeed only used once and
1012  * we could have left the call to `fe_values.shape_hessian(j,qpoint)`
1013  * in the instruction that computes the scalar product between
1014  * the two terms.
1015  *
1016  * @code
1017  * auto cell_worker = [&](const Iterator & cell,
1018  * ScratchData<dim> &scratch_data,
1019  * CopyData & copy_data) {
1020  * copy_data.cell_matrix = 0;
1021  * copy_data.cell_rhs = 0;
1022  *
1023  * FEValues<dim> &fe_values = scratch_data.fe_values;
1024  * fe_values.reinit(cell);
1025  *
1026  * cell->get_dof_indices(copy_data.local_dof_indices);
1027  *
1028  * const ExactSolution::RightHandSide<dim> right_hand_side;
1029  *
1030  * const unsigned int dofs_per_cell =
1031  * scratch_data.fe_values.get_fe().n_dofs_per_cell();
1032  *
1033  * for (unsigned int qpoint = 0; qpoint < fe_values.n_quadrature_points;
1034  * ++qpoint)
1035  * {
1036  * for (unsigned int i = 0; i < dofs_per_cell; ++i)
1037  * {
1038  * const Tensor<2, dim> &hessian_i =
1039  * fe_values.shape_hessian(i, qpoint);
1040  *
1041  * for (unsigned int j = 0; j < dofs_per_cell; ++j)
1042  * {
1043  * const Tensor<2, dim> &hessian_j =
1044  * fe_values.shape_hessian(j, qpoint);
1045  *
1046  * copy_data.cell_matrix(i, j) +=
1047  * scalar_product(hessian_i, // nabla^2 phi_i(x)
1048  * hessian_j) * // nabla^2 phi_j(x)
1049  * fe_values.JxW(qpoint); // dx
1050  * }
1051  *
1052  * copy_data.cell_rhs(i) +=
1053  * fe_values.shape_value(i, qpoint) * // phi_i(x)
1054  * right_hand_side.value(
1055  * fe_values.quadrature_point(qpoint)) * // f(x)
1056  * fe_values.JxW(qpoint); // dx
1057  * }
1058  * }
1059  * };
1060  *
1061  *
1062  * @endcode
1063  *
1064  * The next building block is the one that assembles penalty terms on each
1065  * of the interior faces of the mesh. As described in the documentation of
1066  * MeshWorker::mesh_loop(), this function receives arguments that denote
1067  * a cell and its neighboring cell, as well as (for each of the two
1068  * cells) the face (and potentially sub-face) we have to integrate
1069  * over. Again, we also get a scratch object, and a copy object
1070  * for putting the results in.
1071  *
1072 
1073  *
1074  * The function has three parts itself. At the top, we initialize
1075  * the FEInterfaceValues object and create a new `CopyData::FaceData`
1076  * object to store our input in. This gets pushed to the end of the
1077  * `copy_data.face_data` variable. We need to do this because
1078  * the number of faces (or subfaces) over which we integrate for a
1079  * given cell differs from cell to cell, and the sizes of these
1080  * matrices also differ, depending on what degrees of freedom
1081  * are adjacent to the face or subface. As discussed in the documentation
1082  * of MeshWorker::mesh_loop(), the copy object is reset every time a new
1083  * cell is visited, so that what we push to the end of
1084  * `copy_data.face_data()` is really all that the later `copier` function
1085  * gets to see when it copies the contributions of each cell to the global
1086  * matrix and right hand side objects.
1087  *
1088  * @code
1089  * auto face_worker = [&](const Iterator & cell,
1090  * const unsigned int &f,
1091  * const unsigned int &sf,
1092  * const Iterator & ncell,
1093  * const unsigned int &nf,
1094  * const unsigned int &nsf,
1095  * ScratchData<dim> & scratch_data,
1096  * CopyData & copy_data) {
1097  * FEInterfaceValues<dim> &fe_interface_values =
1098  * scratch_data.fe_interface_values;
1099  * fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
1100  *
1101  * copy_data.face_data.emplace_back();
1102  * CopyData::FaceData &copy_data_face = copy_data.face_data.back();
1103  *
1104  * copy_data_face.joint_dof_indices =
1105  * fe_interface_values.get_interface_dof_indices();
1106  *
1107  * const unsigned int n_interface_dofs =
1108  * fe_interface_values.n_current_interface_dofs();
1109  * copy_data_face.cell_matrix.reinit(n_interface_dofs, n_interface_dofs);
1110  *
1111  * @endcode
1112  *
1113  * The second part deals with determining what the penalty
1114  * parameter should be. By looking at the units of the various
1115  * terms in the bilinear form, it is clear that the penalty has
1116  * to have the form @f$\frac{\gamma}{h_K}@f$ (i.e., one over length
1117  * scale), but it is not a priori obvious how one should choose
1118  * the dimension-less number @f$\gamma@f$. From the discontinuous
1119  * Galerkin theory for the Laplace equation, one might
1120  * conjecture that the right choice is @f$\gamma=p(p+1)@f$ is the
1121  * right choice, where @f$p@f$ is the polynomial degree of the
1122  * finite element used. We will discuss this choice in a bit
1123  * more detail in the results section of this program.
1124  *
1125 
1126  *
1127  * In the formula above, @f$h_K@f$ is the size of cell @f$K@f$. But this
1128  * is not quite so straightforward either: If one uses highly
1129  * stretched cells, then a more involved theory says that @f$h@f$
1130  * should be replaced by the diameter of cell @f$K@f$ normal to the
1131  * direction of the edge in question. It turns out that there
1132  * is a function in deal.II for that. Secondly, @f$h_K@f$ may be
1133  * different when viewed from the two different sides of a face.
1134  *
1135 
1136  *
1137  * To stay on the safe side, we take the maximum of the two values.
1138  * We will note that it is possible that this computation has to be
1139  * further adjusted if one were to use hanging nodes resulting from
1140  * adaptive mesh refinement.
1141  *
1142  * @code
1143  * const unsigned int p = fe.degree;
1144  * const double gamma_over_h =
1145  * std::max((1.0 * p * (p + 1) /
1146  * cell->extent_in_direction(
1147  * GeometryInfo<dim>::unit_normal_direction[f])),
1148  * (1.0 * p * (p + 1) /
1149  * ncell->extent_in_direction(
1150  * GeometryInfo<dim>::unit_normal_direction[nf])));
1151  *
1152  * @endcode
1153  *
1154  * Finally, and as usual, we loop over the quadrature points and
1155  * indices `i` and `j` to add up the contributions of this face
1156  * or sub-face. These are then stored in the
1157  * `copy_data.face_data` object created above. As for the cell
1158  * worker, we pull the evaluation of averages and jumps out of
1159  * the loops if possible, introducing local variables that store
1160  * these results. The assembly then only needs to use these
1161  * local variables in the innermost loop. Regarding the concrete
1162  * formula this code implements, recall that the interface terms
1163  * of the bilinear form were as follows:
1164  * @f{align*}{
1165  * -\sum_{e \in \mathbb{F}} \int_{e}
1166  * \jump{ \frac{\partial v_h}{\partial \mathbf n}}
1167  * \average{\frac{\partial^2 u_h}{\partial \mathbf n^2}} \ ds
1168  * -\sum_{e \in \mathbb{F}} \int_{e}
1169  * \average{\frac{\partial^2 v_h}{\partial \mathbf n^2}}
1170  * \jump{\frac{\partial u_h}{\partial \mathbf n}} \ ds
1171  * + \sum_{e \in \mathbb{F}}
1172  * \frac{\gamma}{h_e}
1173  * \int_e
1174  * \jump{\frac{\partial v_h}{\partial \mathbf n}}
1175  * \jump{\frac{\partial u_h}{\partial \mathbf n}} \ ds.
1176  * @f}
1177  *
1178  * @code
1179  * for (unsigned int qpoint = 0;
1180  * qpoint < fe_interface_values.n_quadrature_points;
1181  * ++qpoint)
1182  * {
1183  * const auto &n = fe_interface_values.normal(qpoint);
1184  *
1185  * for (unsigned int i = 0; i < n_interface_dofs; ++i)
1186  * {
1187  * const double av_hessian_i_dot_n_dot_n =
1188  * (fe_interface_values.average_hessian(i, qpoint) * n * n);
1189  * const double jump_grad_i_dot_n =
1190  * (fe_interface_values.jump_gradient(i, qpoint) * n);
1191  *
1192  * for (unsigned int j = 0; j < n_interface_dofs; ++j)
1193  * {
1194  * const double av_hessian_j_dot_n_dot_n =
1195  * (fe_interface_values.average_hessian(j, qpoint) * n * n);
1196  * const double jump_grad_j_dot_n =
1197  * (fe_interface_values.jump_gradient(j, qpoint) * n);
1198  *
1199  * copy_data_face.cell_matrix(i, j) +=
1200  * (-av_hessian_i_dot_n_dot_n // - {grad^2 v n n }
1201  * * jump_grad_j_dot_n // [grad u n]
1202  * - av_hessian_j_dot_n_dot_n // - {grad^2 u n n }
1203  * * jump_grad_i_dot_n // [grad v n]
1204  * + // +
1205  * gamma_over_h * // gamma/h
1206  * jump_grad_i_dot_n * // [grad v n]
1207  * jump_grad_j_dot_n) * // [grad u n]
1208  * fe_interface_values.JxW(qpoint); // dx
1209  * }
1210  * }
1211  * }
1212  * };
1213  *
1214  *
1215  * @endcode
1216  *
1217  * The third piece is to do the same kind of assembly for faces that
1218  * are at the boundary. The idea is the same as above, of course,
1219  * with only the difference that there are now penalty terms that
1220  * also go into the right hand side.
1221  *
1222 
1223  *
1224  * As before, the first part of the function simply sets up some
1225  * helper objects:
1226  *
1227  * @code
1228  * auto boundary_worker = [&](const Iterator & cell,
1229  * const unsigned int &face_no,
1230  * ScratchData<dim> & scratch_data,
1231  * CopyData & copy_data) {
1232  * FEInterfaceValues<dim> &fe_interface_values =
1233  * scratch_data.fe_interface_values;
1234  * fe_interface_values.reinit(cell, face_no);
1235  * const auto &q_points = fe_interface_values.get_quadrature_points();
1236  *
1237  * copy_data.face_data.emplace_back();
1238  * CopyData::FaceData &copy_data_face = copy_data.face_data.back();
1239  *
1240  * const unsigned int n_dofs =
1241  * fe_interface_values.n_current_interface_dofs();
1242  * copy_data_face.joint_dof_indices =
1243  * fe_interface_values.get_interface_dof_indices();
1244  *
1245  * copy_data_face.cell_matrix.reinit(n_dofs, n_dofs);
1246  *
1247  * const std::vector<double> &JxW = fe_interface_values.get_JxW_values();
1248  * const std::vector<Tensor<1, dim>> &normals =
1249  * fe_interface_values.get_normal_vectors();
1250  *
1251  *
1252  * const ExactSolution::Solution<dim> exact_solution;
1253  * std::vector<Tensor<1, dim>> exact_gradients(q_points.size());
1254  * exact_solution.gradient_list(q_points, exact_gradients);
1255  *
1256  *
1257  * @endcode
1258  *
1259  * Positively, because we now only deal with one cell adjacent to the
1260  * face (as we are on the boundary), the computation of the penalty
1261  * factor @f$\gamma@f$ is substantially simpler:
1262  *
1263  * @code
1264  * const unsigned int p = fe.degree;
1265  * const double gamma_over_h =
1266  * (1.0 * p * (p + 1) /
1267  * cell->extent_in_direction(
1268  * GeometryInfo<dim>::unit_normal_direction[face_no]));
1269  *
1270  * @endcode
1271  *
1272  * The third piece is the assembly of terms. This is now
1273  * slightly more involved since these contains both terms for
1274  * the matrix and for the right hand side. The former is exactly
1275  * the same as for the interior faces stated above if one just
1276  * defines the jump and average appropriately (which is what the
1277  * FEInterfaceValues class does). The latter requires us to
1278  * evaluate the boundary conditions @f$j(\mathbf x)@f$, which in the
1279  * current case (where we know the exact solution) we compute
1280  * from @f$j(\mathbf x) = \frac{\partial u(\mathbf x)}{\partial
1281  * {\mathbf n}}@f$. The term to be added to the right hand side
1282  * vector is then
1283  * @f$\frac{\gamma}{h_e}\int_e
1284  * \jump{\frac{\partial v_h}{\partial \mathbf n}} j \ ds@f$.
1285  *
1286  * @code
1287  * for (unsigned int qpoint = 0; qpoint < q_points.size(); ++qpoint)
1288  * {
1289  * const auto &n = normals[qpoint];
1290  *
1291  * for (unsigned int i = 0; i < n_dofs; ++i)
1292  * {
1293  * const double av_hessian_i_dot_n_dot_n =
1294  * (fe_interface_values.average_hessian(i, qpoint) * n * n);
1295  * const double jump_grad_i_dot_n =
1296  * (fe_interface_values.jump_gradient(i, qpoint) * n);
1297  *
1298  * for (unsigned int j = 0; j < n_dofs; ++j)
1299  * {
1300  * const double av_hessian_j_dot_n_dot_n =
1301  * (fe_interface_values.average_hessian(j, qpoint) * n * n);
1302  * const double jump_grad_j_dot_n =
1303  * (fe_interface_values.jump_gradient(j, qpoint) * n);
1304  *
1305  * copy_data_face.cell_matrix(i, j) +=
1306  * (-av_hessian_i_dot_n_dot_n // - {grad^2 v n n}
1307  * * jump_grad_j_dot_n // [grad u n]
1308  *
1309  * - av_hessian_j_dot_n_dot_n // - {grad^2 u n n}
1310  * * jump_grad_i_dot_n // [grad v n]
1311  *
1312  * + gamma_over_h // gamma/h
1313  * * jump_grad_i_dot_n // [grad v n]
1314  * * jump_grad_j_dot_n // [grad u n]
1315  * ) *
1316  * JxW[qpoint]; // dx
1317  * }
1318  *
1319  * copy_data.cell_rhs(i) +=
1320  * (-av_hessian_i_dot_n_dot_n * // - {grad^2 v n n }
1321  * (exact_gradients[qpoint] * n) // (grad u_exact . n)
1322  * + // +
1323  * gamma_over_h // gamma/h
1324  * * jump_grad_i_dot_n // [grad v n]
1325  * * (exact_gradients[qpoint] * n) // (grad u_exact . n)
1326  * ) *
1327  * JxW[qpoint]; // dx
1328  * }
1329  * }
1330  * };
1331  *
1332  * @endcode
1333  *
1334  * Part 4 is a small function that copies the data produced by the
1335  * cell, interior, and boundary face assemblers above into the
1336  * global matrix and right hand side vector. There really is not
1337  * very much to do here: We distribute the cell matrix and right
1338  * hand side contributions as we have done in almost all of the
1339  * other tutorial programs using the constraints objects. We then
1340  * also have to do the same for the face matrix contributions
1341  * that have gained content for the faces (interior and boundary)
1342  * and that the `face_worker` and `boundary_worker` have added
1343  * to the `copy_data.face_data` array.
1344  *
1345  * @code
1346  * auto copier = [&](const CopyData &copy_data) {
1347  * constraints.distribute_local_to_global(copy_data.cell_matrix,
1348  * copy_data.cell_rhs,
1349  * copy_data.local_dof_indices,
1350  * system_matrix,
1351  * system_rhs);
1352  *
1353  * for (auto &cdf : copy_data.face_data)
1354  * {
1355  * constraints.distribute_local_to_global(cdf.cell_matrix,
1356  * cdf.joint_dof_indices,
1357  * system_matrix);
1358  * }
1359  * };
1360  *
1361  *
1362  * @endcode
1363  *
1364  * Having set all of this up, what remains is to just create a scratch
1365  * and copy data object and call the MeshWorker::mesh_loop() function
1366  * that then goes over all cells and faces, calls the respective workers
1367  * on them, and then the copier function that puts things into the
1368  * global matrix and right hand side. As an additional benefit,
1369  * MeshWorker::mesh_loop() does all of this in parallel, using
1370  * as many processor cores as your machine happens to have.
1371  *
1372  * @code
1373  * const unsigned int n_gauss_points = dof_handler.get_fe().degree + 1;
1374  * ScratchData<dim> scratch_data(mapping,
1375  * fe,
1376  * n_gauss_points,
1377  * update_values | update_gradients |
1378  * update_hessians | update_quadrature_points |
1379  * update_JxW_values,
1380  * update_values | update_gradients |
1381  * update_hessians | update_quadrature_points |
1382  * update_JxW_values | update_normal_vectors);
1383  * CopyData copy_data(dof_handler.get_fe().n_dofs_per_cell());
1384  * MeshWorker::mesh_loop(dof_handler.begin_active(),
1385  * dof_handler.end(),
1386  * cell_worker,
1387  * copier,
1388  * scratch_data,
1389  * copy_data,
1390  * MeshWorker::assemble_own_cells |
1391  * MeshWorker::assemble_boundary_faces |
1392  * MeshWorker::assemble_own_interior_faces_once,
1393  * boundary_worker,
1394  * face_worker);
1395  * }
1396  *
1397  *
1398  *
1399  * @endcode
1400  *
1401  *
1402  * <a name="Solvingthelinearsystemandpostprocessing"></a>
1403  * <h4>Solving the linear system and postprocessing</h4>
1404  *
1405 
1406  *
1407  * The show is essentially over at this point: The remaining functions are
1408  * not overly interesting or novel. The first one simply uses a direct
1409  * solver to solve the linear system (see also @ref step_29 "step-29"):
1410  *
1411  * @code
1412  * template <int dim>
1413  * void BiharmonicProblem<dim>::solve()
1414  * {
1415  * std::cout << " Solving system..." << std::endl;
1416  *
1417  * SparseDirectUMFPACK A_direct;
1418  * A_direct.initialize(system_matrix);
1419  * A_direct.vmult(solution, system_rhs);
1420  *
1421  * constraints.distribute(solution);
1422  * }
1423  *
1424  *
1425  *
1426  * @endcode
1427  *
1428  * The next function evaluates the error between the computed solution
1429  * and the exact solution (which is known here because we have chosen
1430  * the right hand side and boundary values in a way so that we know
1431  * the corresponding solution). In the first two code blocks below,
1432  * we compute the error in the @f$L_2@f$ norm and the @f$H^1@f$ semi-norm.
1433  *
1434  * @code
1435  * template <int dim>
1436  * void BiharmonicProblem<dim>::compute_errors()
1437  * {
1438  * {
1439  * Vector<float> norm_per_cell(triangulation.n_active_cells());
1440  * VectorTools::integrate_difference(mapping,
1441  * dof_handler,
1442  * solution,
1443  * ExactSolution::Solution<dim>(),
1444  * norm_per_cell,
1445  * QGauss<dim>(fe.degree + 2),
1446  * VectorTools::L2_norm);
1447  * const double error_norm =
1448  * VectorTools::compute_global_error(triangulation,
1449  * norm_per_cell,
1450  * VectorTools::L2_norm);
1451  * std::cout << " Error in the L2 norm : " << error_norm
1452  * << std::endl;
1453  * }
1454  *
1455  * {
1456  * Vector<float> norm_per_cell(triangulation.n_active_cells());
1457  * VectorTools::integrate_difference(mapping,
1458  * dof_handler,
1459  * solution,
1460  * ExactSolution::Solution<dim>(),
1461  * norm_per_cell,
1462  * QGauss<dim>(fe.degree + 2),
1463  * VectorTools::H1_seminorm);
1464  * const double error_norm =
1465  * VectorTools::compute_global_error(triangulation,
1466  * norm_per_cell,
1467  * VectorTools::H1_seminorm);
1468  * std::cout << " Error in the H1 seminorm : " << error_norm
1469  * << std::endl;
1470  * }
1471  *
1472  * @endcode
1473  *
1474  * Now also compute an approximation to the @f$H^2@f$ seminorm error. The actual
1475  * @f$H^2@f$ seminorm would require us to integrate second derivatives of the
1476  * solution @f$u_h@f$, but given the Lagrange shape functions we use, @f$u_h@f$ of
1477  * course has kinks at the interfaces between cells, and consequently second
1478  * derivatives are singular at interfaces. As a consequence, we really only
1479  * integrate over the interior of cells and ignore the interface
1480  * contributions. This is *not* an equivalent norm to the energy norm for
1481  * the problem, but still gives us an idea of how fast the error converges.
1482  *
1483 
1484  *
1485  * We note that one could address this issue by defining a norm that
1486  * is equivalent to the energy norm. This would involve adding up not
1487  * only the integrals over cell interiors as we do below, but also adding
1488  * penalty terms for the jump of the derivative of @f$u_h@f$ across interfaces,
1489  * with an appropriate scaling of the two kinds of terms. We will leave
1490  * this for later work.
1491  *
1492  * @code
1493  * {
1494  * const QGauss<dim> quadrature_formula(fe.degree + 2);
1495  * ExactSolution::Solution<dim> exact_solution;
1496  * Vector<double> error_per_cell(triangulation.n_active_cells());
1497  *
1498  * FEValues<dim> fe_values(mapping,
1499  * fe,
1500  * quadrature_formula,
1501  * update_values | update_hessians |
1502  * update_quadrature_points | update_JxW_values);
1503  *
1504  * FEValuesExtractors::Scalar scalar(0);
1505  * const unsigned int n_q_points = quadrature_formula.size();
1506  *
1507  * std::vector<SymmetricTensor<2, dim>> exact_hessians(n_q_points);
1508  * std::vector<Tensor<2, dim>> hessians(n_q_points);
1509  * for (auto &cell : dof_handler.active_cell_iterators())
1510  * {
1511  * fe_values.reinit(cell);
1512  * fe_values[scalar].get_function_hessians(solution, hessians);
1513  * exact_solution.hessian_list(fe_values.get_quadrature_points(),
1514  * exact_hessians);
1515  *
1516  * double local_error = 0;
1517  * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1518  * {
1519  * local_error +=
1520  * ((exact_hessians[q_point] - hessians[q_point]).norm_square() *
1521  * fe_values.JxW(q_point));
1522  * }
1523  * error_per_cell[cell->active_cell_index()] = std::sqrt(local_error);
1524  * }
1525  *
1526  * const double error_norm = error_per_cell.l2_norm();
1527  * std::cout << " Error in the broken H2 seminorm: " << error_norm
1528  * << std::endl;
1529  * }
1530  * }
1531  *
1532  *
1533  *
1534  * @endcode
1535  *
1536  * Equally uninteresting is the function that generates graphical output.
1537  * It looks exactly like the one in @ref step_6 "step-6", for example.
1538  *
1539  * @code
1540  * template <int dim>
1541  * void
1542  * BiharmonicProblem<dim>::output_results(const unsigned int iteration) const
1543  * {
1544  * std::cout << " Writing graphical output..." << std::endl;
1545  *
1546  * DataOut<dim> data_out;
1547  *
1548  * data_out.attach_dof_handler(dof_handler);
1549  * data_out.add_data_vector(solution, "solution");
1550  * data_out.build_patches();
1551  *
1552  * const std::string filename =
1553  * ("output_" + Utilities::int_to_string(iteration, 6) + ".vtu");
1554  * std::ofstream output_vtu(filename);
1555  * data_out.write_vtu(output_vtu);
1556  * }
1557  *
1558  *
1559  *
1560  * @endcode
1561  *
1562  * The same is true for the `run()` function: Just like in previous
1563  * programs.
1564  *
1565  * @code
1566  * template <int dim>
1567  * void BiharmonicProblem<dim>::run()
1568  * {
1569  * make_grid();
1570  *
1571  * const unsigned int n_cycles = 4;
1572  * for (unsigned int cycle = 0; cycle < n_cycles; ++cycle)
1573  * {
1574  * std::cout << "Cycle " << cycle << " of " << n_cycles << std::endl;
1575  *
1576  * triangulation.refine_global(1);
1577  * setup_system();
1578  *
1579  * assemble_system();
1580  * solve();
1581  *
1582  * output_results(cycle);
1583  *
1584  * compute_errors();
1585  * std::cout << std::endl;
1586  * }
1587  * }
1588  * } // namespace Step47
1589  *
1590  *
1591  *
1592  * @endcode
1593  *
1594  *
1595  * <a name="Themainfunction"></a>
1596  * <h3>The main() function</h3>
1597  *
1598 
1599  *
1600  * Finally for the `main()` function. There is, again, not very much to see
1601  * here: It looks like the ones in previous tutorial programs. There
1602  * is a variable that allows selecting the polynomial degree of the element
1603  * we want to use for solving the equation. Because the C0IP formulation
1604  * we use requires the element degree to be at least two, we check with
1605  * an assertion that whatever one sets for the polynomial degree actually
1606  * makes sense.
1607  *
1608  * @code
1609  * int main()
1610  * {
1611  * try
1612  * {
1613  * using namespace dealii;
1614  * using namespace Step47;
1615  *
1616  * const unsigned int fe_degree = 2;
1617  * Assert(fe_degree >= 2,
1618  * ExcMessage("The C0IP formulation for the biharmonic problem "
1619  * "only works if one uses elements of polynomial "
1620  * "degree at least 2."));
1621  *
1622  * BiharmonicProblem<2> biharmonic_problem(fe_degree);
1623  * biharmonic_problem.run();
1624  * }
1625  * catch (std::exception &exc)
1626  * {
1627  * std::cerr << std::endl
1628  * << std::endl
1629  * << "----------------------------------------------------"
1630  * << std::endl;
1631  * std::cerr << "Exception on processing: " << std::endl
1632  * << exc.what() << std::endl
1633  * << "Aborting!" << std::endl
1634  * << "----------------------------------------------------"
1635  * << std::endl;
1636  *
1637  * return 1;
1638  * }
1639  * catch (...)
1640  * {
1641  * std::cerr << std::endl
1642  * << std::endl
1643  * << "----------------------------------------------------"
1644  * << std::endl;
1645  * std::cerr << "Unknown exception!" << std::endl
1646  * << "Aborting!" << std::endl
1647  * << "----------------------------------------------------"
1648  * << std::endl;
1649  * return 1;
1650  * }
1651  *
1652  * return 0;
1653  * }
1654  * @endcode
1655 <a name="Results"></a><h1>Results</h1>
1656 
1657 
1658 We run the program with right hand side and boundary values as
1659 discussed in the introduction. These will produce the
1660 solution @f$u = \sin(\pi x) \sin(\pi y)@f$ on the domain @f$\Omega = (0,1)^2@f$.
1661 We test this setup using @f$Q_2@f$, @f$Q_3@f$, and @f$Q_4@f$ elements, which one can
1662 change via the `fe_degree` variable in the `main()` function. With mesh
1663 refinement, the @f$L_2@f$ convergence rates, @f$H^1@f$-seminorm rate,
1664 and @f$H^2@f$-seminorm convergence of @f$u@f$
1665 should then be around 2, 2, 1 for @f$Q_2@f$ (with the @f$L_2@f$ norm
1666 sub-optimal as discussed in the introduction); 4, 3, 2 for
1667 @f$Q_3@f$; and 5, 4, 3 for @f$Q_4@f$, respectively.
1668 
1669 From the literature, it is not immediately clear what
1670 the penalty parameter @f$\gamma@f$ should be. For example,
1671 @cite Brenner2009 state that it needs to be larger than one, and
1672 choose @f$\gamma=5@f$. The FEniCS/Dolphin tutorial chooses it as
1673 @f$\gamma=8@f$, see
1674 https://fenicsproject.org/docs/dolfin/1.6.0/python/demo/documented/biharmonic/python/documentation.html
1675 . @cite Wells2007 uses a value for @f$\gamma@f$ larger than the
1676 number of edges belonging to an element for Kirchhoff plates (see
1677 their Section 4.2). This suggests that maybe
1678 @f$\gamma = 1@f$, @f$2@f$, are too small; on the other hand, a value
1679 @f$p(p+1)@f$ would be reasonable,
1680 where @f$p@f$ is the degree of polynomials. The last of these choices is
1681 the one one would expect to work by comparing
1682 to the discontinuous Galerkin formulations for the Laplace equation
1683 (see, for example, the discussions in @ref step_39 "step-39" and @ref step_74 "step-74"),
1684 and it will turn out to also work here.
1685 But we should check what value of @f$\gamma@f$ is right, and we will do so
1686 below; changing @f$\gamma@f$ is easy in the two `face_worker` and
1687 `boundary_worker` functions defined in `assemble_system()`.
1688 
1689 
1690 <a name="TestresultsoniQsub2subiwithigammapp1i"></a><h3>Test results on <i>Q<sub>2</sub></i> with <i>&gamma; = p(p+1)</i> </h3>
1691 
1692 
1693 We run the code with differently refined meshes
1694 and get the following convergence rates.
1695 
1696 <table align="center" class="doxtable">
1697  <tr>
1698  <th>Number of refinements </th><th> @f$\|u-u_h^\circ\|_{L_2}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^1}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^2}@f$ </th><th> Conv. rates </th>
1699  </tr>
1700  <tr>
1701  <td> 2 </td><td> 8.780e-03 </td><td> </td><td> 7.095e-02 </td><td> </td><td> 1.645 </td><td> </td>
1702  </tr>
1703  <tr>
1704  <td> 3 </td><td> 3.515e-03 </td><td> 1.32 </td><td> 2.174e-02 </td><td> 1.70 </td><td> 8.121e-01 </td><td> 1.018 </td>
1705  </tr>
1706  <tr>
1707  <td> 4 </td><td> 1.103e-03 </td><td> 1.67 </td><td> 6.106e-03 </td><td> 1.83 </td><td> 4.015e-01 </td><td> 1.016 </td>
1708  </tr>
1709  <tr>
1710  <td> 5 </td><td> 3.084e-04 </td><td> 1.83 </td><td> 1.622e-03 </td><td> 1.91 </td><td> 1.993e-01 </td><td> 1.010 </td>
1711  </tr>
1712 </table>
1713 We can see that the @f$L_2@f$ convergence rates are around 2,
1714 @f$H^1@f$-seminorm convergence rates are around 2,
1715 and @f$H^2@f$-seminorm convergence rates are around 1. The latter two
1716 match the theoretically expected rates; for the former, we have no
1717 theorem but are not surprised that it is sub-optimal given the remark
1718 in the introduction.
1719 
1720 
1721 <a name="TestresultsoniQsub3subiwithigammapp1i"></a><h3>Test results on <i>Q<sub>3</sub></i> with <i>&gamma; = p(p+1)</i> </h3>
1722 
1723 
1724 
1725 <table align="center" class="doxtable">
1726  <tr>
1727  <th>Number of refinements </th><th> @f$\|u-u_h^\circ\|_{L_2}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^1}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^2}@f$ </th><th> Conv. rates </th>
1728  </tr>
1729  <tr>
1730  <td> 2 </td><td> 2.045e-04 </td><td> </td><td> 4.402e-03 </td><td> </td><td> 1.641e-01 </td><td> </td>
1731  </tr>
1732  <tr>
1733  <td> 3 </td><td> 1.312e-05 </td><td> 3.96 </td><td> 5.537e-04 </td><td> 2.99 </td><td> 4.096e-02 </td><td> 2.00 </td>
1734  </tr>
1735  <tr>
1736  <td> 4 </td><td> 8.239e-07 </td><td> 3.99 </td><td> 6.904e-05 </td><td> 3.00 </td><td> 1.023e-02 </td><td> 2.00 </td>
1737  </tr>
1738  <tr>
1739  <td> 5 </td><td> 5.158e-08 </td><td> 3.99 </td><td> 8.621e-06 </td><td> 3.00 </td><td> 2.558e-03 </td><td> 2.00 </td>
1740  </tr>
1741 </table>
1742 We can see that the @f$L_2@f$ convergence rates are around 4,
1743 @f$H^1@f$-seminorm convergence rates are around 3,
1744 and @f$H^2@f$-seminorm convergence rates are around 2.
1745 This, of course, matches our theoretical expectations.
1746 
1747 
1748 <a name="TestresultsoniQsub4subiwithigammapp1i"></a><h3>Test results on <i>Q<sub>4</sub></i> with <i>&gamma; = p(p+1)</i> </h3>
1749 
1750 
1751 <table align="center" class="doxtable">
1752  <tr>
1753  <th>Number of refinements </th><th> @f$\|u-u_h^\circ\|_{L_2}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^1}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^2}@f$ </th><th> Conv. rates </th>
1754  </tr>
1755  <tr>
1756  <td> 2 </td><td> 6.510e-06 </td><td> </td><td> 2.215e-04 </td><td> </td><td> 1.275e-02 </td><td> </td>
1757  </tr>
1758  <tr>
1759  <td> 3 </td><td> 2.679e-07 </td><td> 4.60 </td><td> 1.569e-05 </td><td> 3.81 </td><td> 1.496e-03 </td><td> 3.09 </td>
1760  </tr>
1761  <tr>
1762  <td> 4 </td><td> 9.404e-09 </td><td> 4.83 </td><td> 1.040e-06 </td><td> 3.91 </td><td> 1.774e-04 </td><td> 3.07 </td>
1763  </tr>
1764  <tr>
1765  <td> 5 </td><td> 7.943e-10 </td><td> 3.56 </td><td> 6.693e-08 </td><td> 3.95 </td><td> 2.150e-05 </td><td> 3.04 </td>
1766  </tr>
1767 </table>
1768 We can see that the @f$L_2@f$ norm convergence rates are around 5,
1769 @f$H^1@f$-seminorm convergence rates are around 4,
1770 and @f$H^2@f$-seminorm convergence rates are around 3.
1771 On the finest mesh, the @f$L_2@f$ norm convergence rate
1772 is much smaller than our theoretical expectations
1773 because the linear solver becomes the limiting factor due
1774 to round-off. Of course the @f$L_2@f$ error is also very small already in
1775 that case.
1776 
1777 
1778 <a name="TestresultsoniQsub2subiwithigamma1i"></a><h3>Test results on <i>Q<sub>2</sub></i> with <i>&gamma; = 1</i> </h3>
1779 
1780 
1781 For comparison with the results above, let us now also consider the
1782 case where we simply choose @f$\gamma=1@f$:
1783 
1784 <table align="center" class="doxtable">
1785  <tr>
1786  <th>Number of refinements </th><th> @f$\|u-u_h^\circ\|_{L_2}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^1}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^2}@f$ </th><th> Conv. rates </th>
1787  </tr>
1788  <tr>
1789  <td> 2 </td><td> 7.350e-02 </td><td> </td><td> 7.323e-01 </td><td> </td><td> 10.343 </td><td> </td>
1790  </tr>
1791  <tr>
1792  <td> 3 </td><td> 6.798e-03 </td><td> 3.43 </td><td> 1.716e-01 </td><td> 2.09 </td><td>4.836 </td><td> 1.09 </td>
1793  </tr>
1794  <tr>
1795  <td> 4 </td><td> 9.669e-04 </td><td> 2.81 </td><td> 6.436e-02 </td><td> 1.41 </td><td> 3.590 </td><td> 0.430 </td>
1796  </tr>
1797  <tr>
1798  <td> 5 </td><td> 1.755e-04 </td><td> 2.46 </td><td> 2.831e-02 </td><td> 1.18 </td><td>3.144 </td><td> 0.19 </td>
1799  </tr>
1800 </table>
1801 Although @f$L_2@f$ norm convergence rates of @f$u@f$ more or less
1802 follows the theoretical expectations,
1803 the @f$H^1@f$-seminorm and @f$H^2@f$-seminorm do not seem to converge as expected.
1804 Comparing results from @f$\gamma = 1@f$ and @f$\gamma = p(p+1)@f$, it is clear that
1805 @f$\gamma = p(p+1)@f$ is a better penalty.
1806 Given that @f$\gamma=1@f$ is already too small for @f$Q_2@f$ elements, it may not be surprising that if one repeated the
1807 experiment with a @f$Q_3@f$ element, the results are even more disappointing: One again only obtains convergence
1808 rates of 2, 1, zero -- i.e., no better than for the @f$Q_2@f$ element (although the errors are smaller in magnitude).
1809 Maybe surprisingly, however, one obtains more or less the expected convergence orders when using @f$Q_4@f$
1810 elements. Regardless, this uncertainty suggests that @f$\gamma=1@f$ is at best a risky choice, and at worst an
1811 unreliable one and that we should choose @f$\gamma@f$ larger.
1812 
1813 
1814 <a name="TestresultsoniQsub2subiwithigamma2i"></a><h3>Test results on <i>Q<sub>2</sub></i> with <i>&gamma; = 2</i> </h3>
1815 
1816 
1817 Since @f$\gamma=1@f$ is clearly too small, one might conjecture that
1818 @f$\gamma=2@f$ might actually work better. Here is what one obtains in
1819 that case:
1820 
1821 <table align="center" class="doxtable">
1822  <tr>
1823  <th>Number of refinements </th><th> @f$\|u-u_h^\circ\|_{L_2}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^1}@f$ </th><th> Conv. rates </th><th> @f$|u-u_h|_{H^2}@f$ </th><th> Conv. rates </th>
1824  </tr>
1825  <tr>
1826  <td> 2 </td><td> 4.133e-02 </td><td> </td><td> 2.517e-01 </td><td> </td><td> 3.056 </td><td> </td>
1827  </tr>
1828  <tr>
1829  <td> 3 </td><td> 6.500e-03 </td><td>2.66 </td><td> 5.916e-02 </td><td> 2.08 </td><td>1.444 </td><td> 1.08 </td>
1830  </tr>
1831  <tr>
1832  <td> 4 </td><td> 6.780e-04 </td><td> 3.26 </td><td> 1.203e-02 </td><td> 2.296 </td><td> 6.151e-01 </td><td> 1.231 </td>
1833  </tr>
1834  <tr>
1835  <td> 5 </td><td> 1.622e-04 </td><td> 2.06 </td><td> 2.448e-03 </td><td> 2.297 </td><td> 2.618e-01 </td><td> 1.232 </td>
1836  </tr>
1837 </table>
1838 In this case, the convergence rates more or less follow the
1839 theoretical expectations, but, compared to the results from @f$\gamma =
1840 p(p+1)@f$, are more variable.
1841 Again, we could repeat this kind of experiment for @f$Q_3@f$ and @f$Q_4@f$ elements. In both cases, we will find that we
1842 obtain roughly the expected convergence rates. Of more interest may then be to compare the absolute
1843 size of the errors. While in the table above, for the @f$Q_2@f$ case, the errors on the finest grid are comparable between
1844 the @f$\gamma=p(p+1)@f$ and @f$\gamma=2@f$ case, for @f$Q_3@f$ the errors are substantially larger for @f$\gamma=2@f$ than for
1845 @f$\gamma=p(p+1)@f$. The same is true for the @f$Q_4@f$ case.
1846 
1847 
1848 <a name="Conclusionsforthechoiceofthepenaltyparameter"></a><h3> Conclusions for the choice of the penalty parameter </h3>
1849 
1850 
1851 The conclusions for which of the "reasonable" choices one should use for the penalty parameter
1852 is that @f$\gamma=p(p+1)@f$ yields the expected results. It is, consequently, what the code
1853 uses as currently written.
1854 
1855 
1856 <a name="Possibilitiesforextensions"></a><h3> Possibilities for extensions </h3>
1857 
1858 
1859 There are a number of obvious extensions to this program that would
1860 make sense:
1861 
1862 - The program uses a square domain and a uniform mesh. Real problems
1863  don't come this way, and one should verify convergence also on
1864  domains with other shapes and, in particular, curved boundaries. One
1865  may also be interested in resolving areas of less regularity by
1866  using adaptive mesh refinement.
1867 
1868 - From a more theoretical perspective, the convergence results above
1869  only used the "broken" @f$H^2@f$ seminorm @f$|\cdot|^\circ_{H^2}@f$ instead
1870  of the "equivalent" norm @f$|\cdot|_h@f$. This is good enough to
1871  convince ourselves that the program isn't fundamentally
1872  broken. However, it might be interesting to measure the error in the
1873  actual norm for which we have theoretical results. Implementing this
1874  addition should not be overly difficult using, for example, the
1875  FEInterfaceValues class combined with MeshWorker::mesh_loop() in the
1876  same spirit as we used for the assembly of the linear system.
1877 
1878 
1879 <a name="Derivationforthesimplysupportedplates"></a> <h4> Derivation for the simply supported plates </h4>
1880 
1881 
1882  Similar to the "clamped" boundary condition addressed in the implementation,
1883  we will derive the @f$C^0@f$ IP finite element scheme for simply supported plates:
1884  @f{align*}{
1885  \Delta^2 u(\mathbf x) &= f(\mathbf x)
1886  \qquad \qquad &&\forall \mathbf x \in \Omega,
1887  u(\mathbf x) &= g(\mathbf x) \qquad \qquad
1888  &&\forall \mathbf x \in \partial\Omega, \\
1889  \Delta u(\mathbf x) &= h(\mathbf x) \qquad \qquad
1890  &&\forall \mathbf x \in \partial\Omega.
1891  @f}
1892  We multiply the biharmonic equation by the test function @f$v_h@f$ and integrate over @f$ K @f$ and get:
1893  @f{align*}{
1894  \int_K v_h (\Delta^2 u_h)
1895  &= \int_K (D^2 v_h) : (D^2 u_h)
1896  + \int_{\partial K} v_h \frac{\partial (\Delta u_h)}{\partial \mathbf{n}}
1897  -\int_{\partial K} (\nabla v_h) \cdot (\frac{\partial \nabla u_h}{\partial \mathbf{n}}).
1898  @f}
1899 
1900  Summing up over all cells @f$K \in \mathbb{T}@f$,since normal directions of @f$\Delta u_h@f$ are pointing at
1901  opposite directions on each interior edge shared by two cells and @f$v_h = 0@f$ on @f$\partial \Omega@f$,
1902  @f{align*}{
1903  \sum_{K \in \mathbb{T}} \int_{\partial K} v_h \frac{\partial (\Delta u_h)}{\partial \mathbf{n}} = 0,
1904  @f}
1905  and by the definition of jump over cell interfaces,
1906  @f{align*}{
1907  -\sum_{K \in \mathbb{T}} \int_{\partial K} (\nabla v_h) \cdot (\frac{\partial \nabla u_h}{\partial \mathbf{n}}) = -\sum_{e \in \mathbb{F}} \int_{e} \jump{\frac{\partial v_h}{\partial \mathbf{n}}} (\frac{\partial^2 u_h}{\partial \mathbf{n^2}}).
1908  @f}
1909  We separate interior faces and boundary faces of the domain,
1910  @f{align*}{
1911  -\sum_{K \in \mathbb{T}} \int_{\partial K} (\nabla v_h) \cdot (\frac{\partial \nabla u_h}{\partial \mathbf{n}}) = -\sum_{e \in \mathbb{F}^i} \int_{e} \jump{\frac{\partial v_h}{\partial \mathbf{n}}} (\frac{\partial^2 u_h}{\partial \mathbf{n^2}})
1912  - \sum_{e \in \partial \Omega} \int_{e} \jump{\frac{\partial v_h}{\partial \mathbf{n}}} h,
1913  @f}
1914  where @f$\mathbb{F}^i@f$ is the set of interior faces.
1915  This leads us to
1916  @f{align*}{
1917  \sum_{K \in \mathbb{T}} \int_K (D^2 v_h) : (D^2 u_h) \ dx - \sum_{e \in \mathbb{F}^i} \int_{e} \jump{\frac{\partial v_h}{\partial \mathbf{n}}} (\frac{\partial^2 u_h}{\partial \mathbf{n^2}}) \ ds
1918  = \sum_{K \in \mathbb{T}}\int_{K} v_h f \ dx + \sum_{e\subset\partial\Omega} \int_{e} \jump{\frac{\partial v_h}{\partial \mathbf{n}}} h \ ds.
1919  @f}
1920 
1921  In order to symmetrize and stabilize the discrete problem,
1922  we add symmetrization and stabilization term.
1923  We finally get the @f$C^0@f$ IP finite element scheme for the biharmonic equation:
1924  find @f$u_h@f$ such that @f$u_h =g@f$ on @f$\partial \Omega@f$ and
1925  @f{align*}{
1926  \mathcal{A}(v_h,u_h)&=\mathcal{F}(v_h) \quad \text{holds for all test functions } v_h,
1927  @f}
1928  where
1929  @f{align*}{
1930  \mathcal{A}(v_h,u_h):=&\sum_{K \in \mathbb{T}}\int_K D^2v_h:D^2u_h \ dx
1931  \\
1932  &
1933  -\sum_{e \in \mathbb{F}^i} \int_{e}
1934  \jump{\frac{\partial v_h}{\partial \mathbf n}}
1935  \average{\frac{\partial^2 u_h}{\partial \mathbf n^2}} \ ds
1936  -\sum_{e \in \mathbb{F}^i} \int_{e}
1937  \average{\frac{\partial^2 v_h}{\partial \mathbf n^2}}
1938  \jump{\frac{\partial u_h}{\partial \mathbf n}} \ ds
1939  \\
1940  &+ \sum_{e \in \mathbb{F}^i}
1941  \frac{\gamma}{h_e}
1942  \int_e
1943  \jump{\frac{\partial v_h}{\partial \mathbf n}}
1944  \jump{\frac{\partial u_h}{\partial \mathbf n}} \ ds,
1945  @f}
1946  and
1947  @f{align*}{
1948  \mathcal{F}(v_h)&:=\sum_{K \in \mathbb{T}}\int_{K} v_h f \ dx
1949  +
1950  \sum_{e\subset\partial\Omega}
1951  \int_e \jump{\frac{\partial v_h}{\partial \mathbf n}} h \ ds.
1952  @f}
1953  The implementation of this boundary case is similar to the "clamped" version
1954  except that `boundary_worker` is no longer needed for system assembling
1955  and the right hand side is changed according to the formulation.
1956  *
1957  *
1958 <a name="PlainProg"></a>
1959 <h1> The plain program</h1>
1960 @include "step-47.cc"
1961 */
Definition: fe_q.h:549
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
Definition: point.h:111
__global__ void set(Number *val, const Number s, const size_type N)
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern)
const Event initial
Definition: event.cc:65
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
static const types::blas_int one
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim >>> &Du)
Definition: divergence.h:472
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask())
void run(const Iterator &begin, const typename identity< Iterator >::type &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
Definition: work_stream.h:691
int(&) functions(const void *v1, const void *v2)
static constexpr double PI
Definition: numbers.h:231
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation