733 *
const unsigned int )
const
741 *
class RightHandSide :
public Function<dim>
745 *
const unsigned int component = 0)
const override;
751 *
double RightHandSide<dim>::value(
const Point<dim> &p,
752 *
const unsigned int )
const
762 * The
class that
implements the exact pressure solution has an
763 * oddity in that we implement it as a vector-valued one with two
764 * components. (We say that it has two components in the constructor
765 * where we call the constructor of the base Function class.) In the
766 * `
value()` function, we do not test for the
value of the
767 * `component` argument, which implies that we return the same
value
768 * for both components of the vector-valued function. We do this
769 * because we describe the finite element in use in this program as
770 * a vector-valued system that contains the interior and the
771 * interface pressures, and when we compute errors, we will want to
772 * use the same pressure solution to test both of these components.
776 * class ExactPressure :
public Function<dim>
784 *
const unsigned int component)
const override;
790 *
double ExactPressure<dim>::value(
const Point<dim> &p,
791 *
const unsigned int )
const
819 *
return return_value;
827 * <a name=
"step_61-WGDarcyEquationclassimplementation"></a>
828 * <h3>WGDarcyEquation
class implementation</h3>
833 * <a name=
"step_61-WGDarcyEquationWGDarcyEquation"></a>
834 * <h4>WGDarcyEquation::WGDarcyEquation</h4>
838 * In
this constructor, we create a finite element space
for vector valued
839 *
functions, which will here include the ones used
for the interior and
840 *
interface pressures, @f$p^\circ@f$ and @f$p^\partial@f$.
844 * WGDarcyEquation<dim>::WGDarcyEquation(
const unsigned int degree)
856 * <a name=
"step_61-WGDarcyEquationmake_grid"></a>
857 * <h4>WGDarcyEquation::make_grid</h4>
861 * We generate a mesh on the unit square domain and
refine it.
865 *
void WGDarcyEquation<dim>::make_grid()
881 * <a name=
"step_61-WGDarcyEquationsetup_system"></a>
882 * <h4>WGDarcyEquation::setup_system</h4>
886 * After we have created the mesh above, we distribute degrees of
887 * freedom and resize matrices and vectors. The only piece of
888 * interest in
this function is how we
interpolate the boundary
889 * values
for the pressure. Since the pressure consists of interior
890 * and
interface components, we need to make sure that we only
891 *
interpolate onto that component of the vector-valued solution
892 * space that corresponds to the
interface pressures (as these are
893 * the only ones that are defined on the boundary of the domain). We
894 *
do this via a component mask
object for only the interface
899 *
void WGDarcyEquation<dim>::setup_system()
901 * dof_handler.distribute_dofs(fe);
902 * dof_handler_dgrt.distribute_dofs(fe_dgrt);
904 * std::cout <<
" Number of pressure degrees of freedom: "
905 * << dof_handler.n_dofs() << std::endl;
907 * solution.reinit(dof_handler.n_dofs());
908 * system_rhs.reinit(dof_handler.n_dofs());
912 * constraints.clear();
918 * BoundaryValues<dim>(),
920 * interface_pressure_mask);
921 * constraints.close();
927 * In the bilinear form, there is no integration term over faces
928 * between two neighboring cells, so we can just use
935 * sparsity_pattern.copy_from(dsp);
937 * system_matrix.reinit(sparsity_pattern);
945 * <a name=
"step_61-WGDarcyEquationassemble_system"></a>
946 * <h4>WGDarcyEquation::assemble_system</h4>
950 * This function is more interesting. As detailed in the
951 * introduction, the assembly of the linear system
requires us to
953 * element in the Raviart-Thomas space. As a consequence, we need to
954 * define a Raviart-Thomas finite element object, and have
FEValues
955 * objects that evaluate it at quadrature points. We then need to
956 * compute the
matrix @f$C^
K@f$ on every cell @f$K@f$,
for which we need the
957 * matrices @f$M^
K@f$ and @f$G^
K@f$ mentioned in the introduction.
961 * A
point that may not be obvious is that in all previous tutorial
963 * iterator from a
DoFHandler. This is so that one can call
964 * functions such as
FEValuesBase::get_function_values() that
965 * extract the values of a finite element function (represented by a
966 * vector of DoF values) on the quadrature points of a cell. For
967 * this operation to work, one needs to know which vector elements
968 * correspond to the degrees of freedom on a given cell -- i.e.,
969 * exactly the kind of information and operation provided by the
974 * We could create a
DoFHandler object for the "broken" Raviart-Thomas space
976 * least in the current function, we have no need for any globally defined
977 * degrees of freedom associated with this broken space, but really only
978 * need to reference the shape functions of such a space on the current
979 * cell. As a consequence, we use the fact that one can call
982 * can of course only provide us with information that only
983 * references information about cells, rather than degrees of freedom
984 * enumerated on these cells. So we can't use
986 *
FEValues::shape_value() to obtain the values of shape functions
987 * at quadrature points on the current cell. It is this kind of
988 * functionality we will make use of below. The variable that will
989 * give us this information about the Raviart-Thomas functions below
990 * is then the `fe_values_rt` (and corresponding `fe_face_values_rt`)
995 * Given this introduction, the following declarations should be
1000 *
void WGDarcyEquation<dim>::assemble_system()
1002 *
const QGauss<dim> quadrature_formula(fe_dgrt.degree + 1);
1003 *
const QGauss<dim - 1> face_quadrature_formula(fe_dgrt.degree + 1);
1006 * quadrature_formula,
1010 * face_quadrature_formula,
1016 * quadrature_formula,
1021 * face_quadrature_formula,
1028 *
const unsigned int dofs_per_cell_dgrt = fe_dgrt.n_dofs_per_cell();
1030 *
const unsigned int n_q_points = fe_values.get_quadrature().size();
1031 *
const unsigned int n_q_points_dgrt = fe_values_dgrt.get_quadrature().size();
1033 *
const unsigned int n_face_q_points = fe_face_values.get_quadrature().size();
1035 * RightHandSide<dim> right_hand_side;
1036 * std::vector<double> right_hand_side_values(n_q_points);
1038 *
const Coefficient<dim> coefficient;
1039 * std::vector<Tensor<2, dim>> coefficient_values(n_q_points);
1041 * std::vector<types::global_dof_index> local_dof_indices(
dofs_per_cell);
1046 * Next, let us declare the various cell matrices discussed in the
1060 * @p face component of the shape
functions.
1069 * This
finally gets us in position to
loop over all cells. On
1070 * each cell, we will
first calculate the various cell matrices
1071 * used to construct the local
matrix -- as they depend on the
1072 * cell in question, they need to be re-computed on each cell. We
1073 * need shape
functions for the Raviart-Thomas space as well,
for
1074 * which we need to create
first an iterator to the cell of the
1075 *
triangulation, which we can obtain by assignment from the cell
1079 *
for (
const auto &cell : dof_handler.active_cell_iterators())
1081 * fe_values.
reinit(cell);
1085 * fe_values_dgrt.reinit(cell_dgrt);
1087 * right_hand_side.value_list(fe_values.get_quadrature_points(),
1088 * right_hand_side_values);
1089 * coefficient.value_list(fe_values.get_quadrature_points(),
1090 * coefficient_values);
1094 * The
first cell
matrix we will compute is the @ref GlossMassMatrix
"mass matrix"
1095 *
for the Raviart-Thomas space. Hence, we need to
loop over
1099 * cell_matrix_M = 0;
1100 *
for (
unsigned int q = 0; q < n_q_points_dgrt; ++q)
1101 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1103 *
const Tensor<1, dim> v_i = fe_values_dgrt[velocities].value(i, q);
1104 *
for (
unsigned int k = 0; k < dofs_per_cell_dgrt; ++k)
1107 * fe_values_dgrt[velocities].value(k, q);
1108 * cell_matrix_M(i, k) += (v_i * v_k * fe_values_dgrt.JxW(q));
1113 * Next we take the inverse of
this matrix by
using
1115 * the coefficient
matrix @f$C^
K@f$ later. It is worth recalling
1116 * later that `cell_matrix_M` actually contains the *inverse*
1117 * of @f$M^
K@f$ after
this call.
1120 * cell_matrix_M.gauss_jordan();
1124 * From the introduction, we know that the right hand side
1125 * @f$G^
K@f$ of the equation that defines @f$C^
K@f$ is the difference
1126 * between a face integral and a cell integral. Here, we
1128 * interior. Each component of
this matrix is the integral of
1129 * a product between a basis function of the polynomial space
1130 * and the divergence of a basis function of the
1131 * Raviart-Thomas space. These basis
functions are defined in
1135 * cell_matrix_G = 0;
1136 *
for (
unsigned int q = 0; q < n_q_points; ++q)
1137 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1139 *
const double div_v_i =
1140 * fe_values_dgrt[velocities].divergence(i, q);
1143 *
const double phi_j_interior =
1144 * fe_values[pressure_interior].value(j, q);
1146 * cell_matrix_G(i, j) -=
1147 * (div_v_i * phi_j_interior * fe_values.JxW(q));
1155 * Each component is the integral of a product between a basis function
1156 * of the polynomial space and the dot product of a basis function of
1157 * the Raviart-Thomas space and the normal vector. So we
loop over all
1158 * the faces of the element and obtain the normal vector.
1161 *
for (
const auto &face : cell->face_iterators())
1163 * fe_face_values.
reinit(cell, face);
1164 * fe_face_values_dgrt.reinit(cell_dgrt, face);
1166 *
for (
unsigned int q = 0; q < n_face_q_points; ++q)
1168 *
const Tensor<1, dim> &normal = fe_face_values.normal_vector(q);
1170 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1173 * fe_face_values_dgrt[velocities].value(i, q);
1176 *
const double phi_j_face =
1177 * fe_face_values[pressure_face].value(j, q);
1179 * cell_matrix_G(i, j) +=
1180 * ((v_i * normal) * phi_j_face * fe_face_values.JxW(q));
1188 * @p cell_matrix_C is then the
matrix product between the
1190 * (where
this inverse is stored in @p cell_matrix_M):
1193 * cell_matrix_G.Tmmult(cell_matrix_C, cell_matrix_M);
1197 * Finally we can compute the local
matrix @f$A^
K@f$. Element
1198 * @f$A^K_{ij}@f$ is given by @f$\int_{
E} \sum_{k,
l} C_{ik} C_{jl}
1199 * (\mathbf{
K} \mathbf{v}_k) \cdot \mathbf{v}_l
1200 * \mathrm{
d}x@f$. We have calculated the coefficients @f$C@f$ in
1201 * the previous step, and so obtain the following after
1202 * suitably re-arranging the loops:
1206 *
for (
unsigned int q = 0; q < n_q_points_dgrt; ++q)
1208 *
for (
unsigned int k = 0; k < dofs_per_cell_dgrt; ++k)
1211 * fe_values_dgrt[velocities].value(k, q);
1212 *
for (
unsigned int l = 0;
l < dofs_per_cell_dgrt; ++
l)
1215 * fe_values_dgrt[velocities].value(l, q);
1217 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1218 *
for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1219 * local_matrix(i, j) +=
1220 * (coefficient_values[q] * cell_matrix_C[i][k] * v_k) *
1221 * cell_matrix_C[j][
l] * v_l * fe_values_dgrt.JxW(q);
1228 * Next, we calculate the right hand side, @f$\int_{
K} f q \mathrm{
d}x@f$:
1232 *
for (
unsigned int q = 0; q < n_q_points; ++q)
1233 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1235 * cell_rhs(i) += (fe_values[pressure_interior].value(i, q) *
1236 * right_hand_side_values[q] * fe_values.JxW(q));
1241 * The last step is to distribute components of the local
1242 *
matrix into the system
matrix and transfer components of
1243 * the cell right hand side into the system right hand side:
1246 * cell->get_dof_indices(local_dof_indices);
1247 * constraints.distribute_local_to_global(
1248 * local_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs);
1257 * <a name=
"step_61-WGDarcyEquationdimsolve"></a>
1258 * <h4>WGDarcyEquation<dim>::solve</h4>
1262 * This step is rather trivial and the same as in many previous
1263 * tutorial programs:
1266 *
template <
int dim>
1267 *
void WGDarcyEquation<dim>::solve()
1269 *
SolverControl solver_control(1000, 1e-8 * system_rhs.l2_norm());
1272 * constraints.distribute(solution);
1279 * <a name=
"step_61-WGDarcyEquationdimcompute_postprocessed_velocity"></a>
1280 * <h4>WGDarcyEquation<dim>::compute_postprocessed_velocity</h4>
1284 * In
this function, compute the velocity field from the pressure
1285 * solution previously computed. The
1286 * velocity is defined as @f$\mathbf{u}_h = \mathbf{Q}_h \left(
1287 * -\mathbf{
K}\nabla_{
w,
d}p_h \right)@f$, which
requires us to compute
1288 * many of the same terms as in the assembly of the system
matrix.
1289 * There are also the matrices @f$E^K,D^K@f$ we need to
assemble (see
1290 * the introduction) but they really just follow the same kind of
1295 * Computing the same matrices here as we have already done in the
1296 * `assemble_system()` function is of course wasteful in terms of
1297 * CPU time. Likewise, we
copy some of the code from there to
this
1298 * function, and
this is also generally a poor idea. A better
1299 * implementation might provide
for a function that encapsulates
1300 *
this duplicated code. One could also think of
using the classic
1301 * trade-off between computing efficiency and memory efficiency to
1302 * only compute the @f$C^
K@f$ matrices once per cell during the
1303 * assembly, storing them somewhere on the side, and re-
using them
1304 * here. (This is what @ref step_51
"step-51" does,
for example, where the
1305 * `assemble_system()` function takes an argument that determines
1306 * whether the local matrices are recomputed, and a similar approach
1307 * -- maybe with storing local matrices elsewhere -- could be
1308 * adapted
for the current program.)
1311 *
template <
int dim>
1312 *
void WGDarcyEquation<dim>::compute_postprocessed_velocity()
1314 * darcy_velocity.reinit(dof_handler_dgrt.n_dofs());
1316 *
const QGauss<dim> quadrature_formula(fe_dgrt.degree + 1);
1317 *
const QGauss<dim - 1> face_quadrature_formula(fe_dgrt.degree + 1);
1320 * quadrature_formula,
1325 * face_quadrature_formula,
1331 * quadrature_formula,
1337 * face_quadrature_formula,
1343 *
const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1344 *
const unsigned int dofs_per_cell_dgrt = fe_dgrt.n_dofs_per_cell();
1346 *
const unsigned int n_q_points = fe_values.get_quadrature().size();
1347 *
const unsigned int n_q_points_dgrt = fe_values_dgrt.get_quadrature().size();
1349 *
const unsigned int n_face_q_points = fe_face_values.get_quadrature().size();
1352 * std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1353 * std::vector<types::global_dof_index> local_dof_indices_dgrt(
1354 * dofs_per_cell_dgrt);
1365 *
const Coefficient<dim> coefficient;
1366 * std::vector<Tensor<2, dim>> coefficient_values(n_q_points_dgrt);
1374 * In the introduction, we explained how to calculate the numerical velocity
1375 * on the cell. We need the pressure solution values on each cell,
1376 * coefficients of the Gram
matrix and coefficients of the @f$L_2@f$ projection.
1377 * We have already calculated the global solution, so we will
extract the
1378 * cell solution from the global solution. The coefficients of the Gram
1379 *
matrix have been calculated when we assembled the system
matrix for the
1380 * pressures. We will
do the same way here. For the coefficients of the
1381 * projection, we
do matrix multiplication, i.e., the inverse of the Gram
1382 *
matrix times the
matrix with @f$(\mathbf{
K} \mathbf{
w}, \mathbf{
w})@f$ as
1383 * components. Then, we multiply all these coefficients and call them beta.
1384 * The numerical velocity is the product of beta and the basis functions of
1385 * the Raviart-Thomas space.
1389 * cell = dof_handler.begin_active(),
1390 * endc = dof_handler.end(), cell_dgrt = dof_handler_dgrt.begin_active();
1391 *
for (; cell != endc; ++cell, ++cell_dgrt)
1393 * fe_values.reinit(cell);
1394 * fe_values_dgrt.reinit(cell_dgrt);
1396 * coefficient.value_list(fe_values_dgrt.get_quadrature_points(),
1397 * coefficient_values);
1401 * The component of this <code>cell_matrix_E</code> is the integral of
1402 * @f$(\mathbf{
K} \mathbf{
w}, \mathbf{
w})@f$. <code>cell_matrix_M</code> is
1406 * cell_matrix_M = 0;
1407 * cell_matrix_E = 0;
1408 *
for (
unsigned int q = 0; q < n_q_points_dgrt; ++q)
1409 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1411 *
const Tensor<1, dim> v_i = fe_values_dgrt[velocities].value(i, q);
1412 *
for (
unsigned int k = 0; k < dofs_per_cell_dgrt; ++k)
1415 * fe_values_dgrt[velocities].value(k, q);
1417 * cell_matrix_E(i, k) +=
1418 * (coefficient_values[q] * v_i * v_k * fe_values_dgrt.JxW(q));
1420 * cell_matrix_M(i, k) += (v_i * v_k * fe_values_dgrt.JxW(q));
1426 * To compute the
matrix @f$D@f$ mentioned in the introduction, we
1427 * then need to evaluate @f$D=M^{-1}
E@f$ as explained in the
1431 * cell_matrix_M.gauss_jordan();
1432 * cell_matrix_M.mmult(cell_matrix_D, cell_matrix_E);
1436 * Then we also need, again, to compute the
matrix @f$C@f$ that is
1437 * used to evaluate the weak discrete
gradient. This is the
1438 * exact same code as used in the assembly of the system
1442 * cell_matrix_G = 0;
1443 *
for (
unsigned int q = 0; q < n_q_points; ++q)
1444 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1446 *
const double div_v_i =
1447 * fe_values_dgrt[velocities].divergence(i, q);
1448 *
for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1450 *
const double phi_j_interior =
1451 * fe_values[pressure_interior].value(j, q);
1453 * cell_matrix_G(i, j) -=
1454 * (div_v_i * phi_j_interior * fe_values.JxW(q));
1458 *
for (
const auto &face : cell->face_iterators())
1460 * fe_face_values.
reinit(cell, face);
1461 * fe_face_values_dgrt.reinit(cell_dgrt, face);
1463 *
for (
unsigned int q = 0; q < n_face_q_points; ++q)
1465 *
const Tensor<1, dim> &normal = fe_face_values.normal_vector(q);
1467 *
for (
unsigned int i = 0; i < dofs_per_cell_dgrt; ++i)
1470 * fe_face_values_dgrt[velocities].value(i, q);
1471 *
for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1473 *
const double phi_j_face =
1474 * fe_face_values[pressure_face].value(j, q);
1476 * cell_matrix_G(i, j) +=
1477 * ((v_i * normal) * phi_j_face * fe_face_values.JxW(q));
1482 * cell_matrix_G.Tmmult(cell_matrix_C, cell_matrix_M);
1486 * Finally, we need to
extract the pressure unknowns that
1487 * correspond to the current cell:
1490 * cell->get_dof_values(solution, cell_solution);
1494 * We are now in a position to compute the local velocity
1495 * unknowns (with respect to the Raviart-Thomas space we are
1496 * projecting the term @f$-\mathbf K \nabla_{
w,
d} p_h@f$ into):
1499 * cell_velocity = 0;
1500 *
for (
unsigned int k = 0; k < dofs_per_cell_dgrt; ++k)
1501 *
for (
unsigned int j = 0; j < dofs_per_cell_dgrt; ++j)
1502 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1503 * cell_velocity(k) +=
1504 * -(cell_solution(i) * cell_matrix_C(i, j) * cell_matrix_D(k, j));
1508 * We compute Darcy velocity.
1509 * This is same as cell_velocity but used to graph Darcy velocity.
1512 * cell_dgrt->get_dof_indices(local_dof_indices_dgrt);
1513 *
for (
unsigned int k = 0; k < dofs_per_cell_dgrt; ++k)
1514 *
for (
unsigned int j = 0; j < dofs_per_cell_dgrt; ++j)
1515 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1516 * darcy_velocity(local_dof_indices_dgrt[k]) +=
1517 * -(cell_solution(i) * cell_matrix_C(i, j) * cell_matrix_D(k, j));
1526 * <a name=
"step_61-WGDarcyEquationdimcompute_pressure_error"></a>
1527 * <h4>WGDarcyEquation<dim>::compute_pressure_error</h4>
1531 * This part is to calculate the @f$L_2@f$ error of the pressure. We
1532 * define a vector that holds the
norm of the error on each cell.
1533 * Next, we use VectorTool::integrate_difference() to compute the
1534 * error in the @f$L_2@f$ norm on each cell. However, we really only
1535 * care about the error in the interior component of the solution
1536 * vector (we can't even evaluate the interface pressures at the
1537 * quadrature points because these are all located in the interior
1538 * of cells) and consequently have to use a weight function that
1539 * ensures that the interface component of the solution variable is
1541 * arguments indicate which component we want to select (component
1542 * zero, i.e., the interior pressures) and how many components there
1543 * are in total (two).
1546 * template <
int dim>
1547 *
void WGDarcyEquation<dim>::compute_pressure_error()
1553 * ExactPressure<dim>(),
1554 * difference_per_cell,
1557 * &select_interior_pressure);
1559 *
const double L2_error = difference_per_cell.l2_norm();
1560 * std::cout <<
"L2_error_pressure " << L2_error << std::endl;
1568 * <a name=
"step_61-WGDarcyEquationdimcompute_velocity_error"></a>
1569 * <h4>WGDarcyEquation<dim>::compute_velocity_error</h4>
1573 * In
this function, we evaluate @f$L_2@f$ errors
for the velocity on
1574 * each cell, and @f$L_2@f$ errors
for the flux on faces. The function
1575 * relies on the `compute_postprocessed_velocity()` function having
1576 * previous computed, which computes the velocity field based on the
1577 * pressure solution that has previously been computed.
1581 * We are going to evaluate velocities on each cell and calculate
1582 * the difference between numerical and exact velocities.
1585 *
template <
int dim>
1586 *
void WGDarcyEquation<dim>::compute_velocity_errors()
1588 *
const QGauss<dim> quadrature_formula(fe_dgrt.degree + 1);
1589 *
const QGauss<dim - 1> face_quadrature_formula(fe_dgrt.degree + 1);
1592 * quadrature_formula,
1598 * face_quadrature_formula,
1604 *
const unsigned int n_q_points_dgrt = fe_values_dgrt.get_quadrature().size();
1605 *
const unsigned int n_face_q_points_dgrt =
1606 * fe_face_values_dgrt.get_quadrature().size();
1608 * std::vector<Tensor<1, dim>> velocity_values(n_q_points_dgrt);
1609 * std::vector<Tensor<1, dim>> velocity_face_values(n_face_q_points_dgrt);
1613 *
const ExactVelocity<dim> exact_velocity;
1615 *
double L2_err_velocity_cell_sqr_global = 0;
1616 *
double L2_err_flux_sqr = 0;
1620 * Having previously computed the postprocessed velocity, we here
1621 * only have to
extract the corresponding values on each cell and
1622 * face and compare it to the exact values.
1625 *
for (
const auto &cell_dgrt : dof_handler_dgrt.active_cell_iterators())
1627 * fe_values_dgrt.
reinit(cell_dgrt);
1631 * First compute the @f$L_2@f$ error between the postprocessed velocity
1632 * field and the exact one:
1635 * fe_values_dgrt[velocities].get_function_values(darcy_velocity,
1637 *
double L2_err_velocity_cell_sqr_local = 0;
1638 *
for (
unsigned int q = 0; q < n_q_points_dgrt; ++q)
1642 * exact_velocity.value(fe_values_dgrt.quadrature_point(q));
1644 * L2_err_velocity_cell_sqr_local +=
1645 * ((velocity - true_velocity) * (velocity - true_velocity) *
1646 * fe_values_dgrt.JxW(q));
1648 * L2_err_velocity_cell_sqr_global += L2_err_velocity_cell_sqr_local;
1652 * For reconstructing the flux we need the size of cells and
1653 * faces. Since fluxes are calculated on faces, we have the
1654 *
loop over all four faces of each cell. To calculate the
1655 * face velocity, we
extract values at the quadrature points from the
1656 * `darcy_velocity` which we have computed previously. Then, we
1657 * calculate the squared velocity error in normal direction. Finally, we
1658 * calculate the @f$L_2@f$ flux error on the cell by appropriately scaling
1659 * with face and cell areas and add it to the global error.
1662 *
const double cell_area = cell_dgrt->measure();
1663 *
for (
const auto &face_dgrt : cell_dgrt->face_iterators())
1665 * const double face_length = face_dgrt->measure();
1666 * fe_face_values_dgrt.reinit(cell_dgrt, face_dgrt);
1667 * fe_face_values_dgrt[velocities].get_function_values(
1668 * darcy_velocity, velocity_face_values);
1670 *
double L2_err_flux_face_sqr_local = 0;
1671 *
for (
unsigned int q = 0; q < n_face_q_points_dgrt; ++q)
1675 * exact_velocity.value(fe_face_values_dgrt.quadrature_point(q));
1678 * fe_face_values_dgrt.normal_vector(q);
1680 * L2_err_flux_face_sqr_local +=
1681 * ((velocity * normal - true_velocity * normal) *
1682 * (velocity * normal - true_velocity * normal) *
1683 * fe_face_values_dgrt.JxW(q));
1685 *
const double err_flux_each_face =
1686 * L2_err_flux_face_sqr_local / face_length * cell_area;
1687 * L2_err_flux_sqr += err_flux_each_face;
1693 * After adding up errors over all cells and faces, we take the
1694 * square root and get the @f$L_2@f$ errors of velocity and
1695 * flux. These we output to screen.
1698 *
const double L2_err_velocity_cell =
1699 *
std::sqrt(L2_err_velocity_cell_sqr_global);
1700 *
const double L2_err_flux_face =
std::sqrt(L2_err_flux_sqr);
1702 * std::cout <<
"L2_error_vel: " << L2_err_velocity_cell << std::endl
1703 * <<
"L2_error_flux: " << L2_err_flux_face << std::endl;
1710 * <a name=
"step_61-WGDarcyEquationoutput_results"></a>
1711 * <h4>WGDarcyEquation::output_results</h4>
1715 * We have two sets of results to output: the interior solution and
1716 * the skeleton solution. We use <code>
DataOut</code> to visualize
1717 * interior results. The graphical output
for the skeleton results
1722 * In both of the output files, both the interior and the face
1723 * variables are stored. For the
interface output, the output file
1724 * simply contains the interpolation of the interior pressures onto
1725 * the faces, but because it is undefined which of the two interior
1726 * pressure variables you get from the two adjacent cells, it is
1727 * best to ignore the interior pressure in the
interface output
1728 * file. Conversely,
for the cell interior output file, it is of
1729 * course impossible to show any
interface pressures @f$p^\partial@f$,
1730 * because these are only available on interfaces and not cell
1731 * interiors. Consequently, you will see them shown as an
invalid
1732 *
value (such as an infinity).
1736 * For the cell interior output, we also want to output the velocity
1737 * variables. This is a bit tricky since it lives on the same mesh
1738 * but uses a different
DoFHandler object (the pressure variables live
1739 * on the `dof_handler`
object, the Darcy velocity on the `dof_handler_dgrt`
1740 *
object). Fortunately, there are variations of the
1742 *
DoFHandler a vector corresponds to, and consequently we can visualize
1743 * the data from both
DoFHandler objects within the same file.
1746 * template <
int dim>
1747 *
void WGDarcyEquation<dim>::output_results() const
1754 * First attach the pressure solution to the
DataOut object:
1757 *
const std::vector<std::string> solution_names = {
"interior_pressure",
1758 *
"interface_pressure"};
1763 * Then
do the same with the Darcy velocity field, and
continue
1764 * with writing everything out into a file.
1767 *
const std::vector<std::string> velocity_names(dim,
"velocity");
1768 *
const std::vector<
1770 * velocity_component_interpretation(
1772 * data_out.add_data_vector(dof_handler_dgrt,
1775 * velocity_component_interpretation);
1777 * data_out.build_patches(fe.degree);
1778 * std::ofstream output(
"solution_interior.vtu");
1779 * data_out.write_vtu(output);
1784 * data_out_faces.attach_dof_handler(dof_handler);
1785 * data_out_faces.add_data_vector(solution,
"Pressure_Face");
1786 * data_out_faces.build_patches(fe.degree);
1787 * std::ofstream face_output(
"solution_interface.vtu");
1788 * data_out_faces.write_vtu(face_output);
1796 * <a name=
"step_61-WGDarcyEquationrun"></a>
1797 * <h4>WGDarcyEquation::run</h4>
1801 * This is the
final function of the main
class. It calls the other
functions
1805 *
template <
int dim>
1806 *
void WGDarcyEquation<dim>::run()
1808 * std::cout <<
"Solving problem in " << dim <<
" space dimensions."
1812 * assemble_system();
1814 * compute_postprocessed_velocity();
1815 * compute_pressure_error();
1816 * compute_velocity_errors();
1826 * <a name=
"step_61-Thecodemaincodefunction"></a>
1827 * <h3>The <code>main</code> function</h3>
1831 * This is the main function. We can change the dimension here to
run in 3
d.
1838 * Step61::WGDarcyEquation<2> wg_darcy(0);
1841 *
catch (std::exception &exc)
1843 * std::cerr << std::endl
1845 * <<
"----------------------------------------------------"
1847 * std::cerr <<
"Exception on processing: " << std::endl
1848 * << exc.what() << std::endl
1849 * <<
"Aborting!" << std::endl
1850 * <<
"----------------------------------------------------"
1856 * std::cerr << std::endl
1858 * <<
"----------------------------------------------------"
1860 * std::cerr <<
"Unknown exception!" << std::endl
1861 * <<
"Aborting!" << std::endl
1862 * <<
"----------------------------------------------------"
1870<a name=
"step_61-Results"></a><h1>Results</h1>
1873We
run the program with a right hand side that will produce the
1874solution @f$p = \sin(\pi x) \sin(\pi y)@f$ and with homogeneous Dirichlet
1875boundary conditions in the domain @f$\Omega = (0,1)^2@f$. In addition, we
1876choose the coefficient matrix in the differential
operator
1878@f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$, @f$\mbox{WG}(Q_1,Q_1;RT_{[1]})@f$ and
1879@f$\mbox{WG}(Q_2,Q_2;RT_{[2]})@f$ element combinations, which one can
1880select by
using the appropriate constructor argument
for the
1881`WGDarcyEquation`
object in `main()`. We will then visualize pressure
1882values in interiors of cells and on faces. We want to see that the
1883pressure maximum is around 1 and the minimum is around 0. With mesh
1884refinement, the convergence rates of pressure, velocity and flux
1885should then be around 1
for @f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$ , 2
for
1886@f$\mbox{WG}(Q_1,Q_1;RT_{[1]})@f$, and 3
for
1887@f$\mbox{WG}(Q_2,Q_2;RT_{[2]})@f$.
1890<a name=
"step_61-TestresultsoniWGQsub0subQsub0subRTsub0subi"></a><h3>Test results on <i>WG(Q<sub>0</sub>,Q<sub>0</sub>;RT<sub>[0]</sub>)</i></h3>
1893The following figures show interior pressures and face pressures
using the
1894@f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$ element. The mesh is refined 2 times (top)
1895and 4 times (bottom), respectively. (This number can be adjusted in the
1896`make_grid()` function.) When the mesh is coarse, one can see
1897the face pressures @f$p^\partial@f$ neatly between the values of the interior
1898pressures @f$p^\circ@f$ on the two adjacent cells.
1900<table align=
"center">
1902 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg000_2d_2.png" alt=
""></td>
1903 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg000_3d_2.png" alt=
""></td>
1906 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg000_2d_4.png" alt=
""></td>
1907 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg000_3d_4.png" alt=
""></td>
1911From the figures, we can see that with the mesh refinement, the maximum and
1912minimum pressure values are approaching the values we expect.
1913Since the mesh is a rectangular mesh and
numbers of cells in each direction is even, we
1914have
symmetric solutions. From the 3
d figures on the right,
1915we can see that on @f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$, the pressure is a constant
1916in the interior of the cell, as expected.
1918<a name=
"step_61-Convergencetableforik0i"></a><h4>Convergence table for <i>k=0</i></h4>
1921We run the code with differently refined meshes (chosen in the `make_grid()` function)
1922and get the following convergence rates of pressure,
1923velocity, and flux (as defined in the introduction).
1925<table align=
"center" class=
"doxtable">
1927 <th>number of refinements </th><th> @f$\|p-p_h^\circ\|@f$ </th><th> @f$\|\mathbf{u}-\mathbf{u}_h\|@f$ </th><th> @f$\|(\mathbf{u}-\mathbf{u}_h) \cdot \mathbf{n}\|@f$ </th>
1930 <td> 2 </td><td> 1.587e-01 </td><td> 5.113e-01 </td><td> 7.062e-01 </td>
1933 <td> 3 </td><td> 8.000e-02 </td><td> 2.529e-01 </td><td> 3.554e-01 </td>
1936 <td> 4 </td><td> 4.006e-02 </td><td> 1.260e-01 </td><td> 1.780e-01 </td>
1939 <td> 5 </td><td> 2.004e-02 </td><td> 6.297e-02 </td><td> 8.902e-02 </td>
1942 <th>Conv.rate </th><th> 1.00 </th><th> 1.00 </th><th> 1.00 </th>
1946We can see that the convergence rates of @f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$ are around 1.
1947This, of course, matches our theoretical expectations.
1950<a name=
"step_61-TestresultsoniWGQsub1subQsub1subRTsub1subi"></a><h3>Test results on <i>WG(Q<sub>1</sub>,Q<sub>1</sub>;RT<sub>[1]</sub>)</i></h3>
1953We can repeat the experiment from above
using the next higher polynomial
1955The following figures are interior pressures and face pressures implemented using
1956@f$\mbox{WG}(Q_1,Q_1;RT_{[1]})@f$. The mesh is refined 4 times. Compared to the
1957previous figures
using
1958@f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$, on each cell, the solution is no longer constant
1959on each cell, as we now use bilinear polynomials to
do the approximation.
1960Consequently, there are 4 pressure values in one interior, 2 pressure values on
1963<table align=
"center">
1965 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg111_2d_4.png" alt=
""></td>
1966 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg111_3d_4.png" alt=
""></td>
1970Compared to the corresponding image
for the @f$\mbox{WG}(Q_0,Q_0;RT_{[0]})@f$
1971combination, the solution is now substantially more accurate and, in
1972particular so close to being continuous at the interfaces that we can
1973no longer distinguish the interface pressures @f$p^\partial@f$ from the
1974interior pressures @f$p^\circ@f$ on the adjacent cells.
1976<a name=
"step_61-Convergencetableforik1i"></a><h4>Convergence table for <i>k=1</i></h4>
1979The following are the convergence rates of pressure, velocity, and flux
1980we obtain from
using the @f$\mbox{WG}(Q_1,Q_1;RT_{[1]})@f$ element combination:
1982<table align=
"center" class=
"doxtable">
1984 <th>number of refinements </th><th> @f$\|p-p_h^\circ\|@f$ </th><th> @f$\|\mathbf{u}-\mathbf{u}_h\|@f$ </th><th> @f$\|(\mathbf{u}-\mathbf{u}_h) \cdot \mathbf{n}\|@f$ </th>
1987 <td> 2 </td><td> 1.613e-02 </td><td> 5.093e-02 </td><td> 7.167e-02 </td>
1990 <td> 3 </td><td> 4.056e-03 </td><td> 1.276e-02 </td><td> 1.802e-02 </td>
1993 <td> 4 </td><td> 1.015e-03 </td><td> 3.191e-03 </td><td> 4.512e-03 </td>
1996 <td> 5 </td><td> 2.540e-04 </td><td> 7.979e-04 </td><td> 1.128e-03 </td>
1999 <th>Conv.rate </th><th> 2.00 </th><th> 2.00 </th><th> 2.00 </th>
2003The convergence rates of @f$WG(Q_1,Q_1;RT_{[1]})@f$ are around 2, as expected.
2007<a name=
"step_61-TestresultsoniWGQsub2subQsub2subRTsub2subi"></a><h3>Test results on <i>WG(Q<sub>2</sub>,Q<sub>2</sub>;RT<sub>[2]</sub>)</i></h3>
2010Let us go one polynomial degree higher.
2011The following are interior pressures and face pressures implemented
using
2012@f$WG(Q_2,Q_2;RT_{[2]})@f$, with mesh size @f$h = 1/32@f$ (i.e., 5 global mesh
2013refinement steps). In the program, we use
2014`data_out_face.build_patches(fe.degree)` when generating graphical output
2016we divide each 2
d cell interior into 4 subcells in order to provide a better
2017visualization of the quadratic polynomials.
2018<table align=
"center">
2020 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg222_2d_5.png" alt=
""></td>
2021 <td><img src=
"https://www.dealii.org/images/steps/developer/step-61.wg222_3d_5.png" alt=
""></td>
2026<a name=
"step_61-Convergencetableforik2i"></a><h4>Convergence table for <i>k=2</i></h4>
2029As before, we can generate convergence data
for the
2030@f$L_2@f$ errors of pressure, velocity, and flux
2031using the @f$\mbox{WG}(Q_2,Q_2;RT_{[2]})@f$ combination:
2033<table align=
"center" class=
"doxtable">
2035 <th>number of refinements </th><th> @f$\|p-p_h^\circ\|@f$ </th><th> @f$\|\mathbf{u}-\mathbf{u}_h\|@f$ </th><th> @f$\|(\mathbf{u}-\mathbf{u}_h) \cdot \mathbf{n}\|@f$ </th>
2038 <td> 2 </td><td> 1.072e-03 </td><td> 3.375e-03 </td><td> 4.762e-03 </td>
2041 <td> 3 </td><td> 1.347e-04 </td><td> 4.233e-04 </td><td> 5.982e-04 </td>
2044 <td> 4 </td><td> 1.685e-05 </td><td> 5.295e-05 </td><td> 7.487e-05 </td>
2047 <td> 5 </td><td> 2.107e-06 </td><td> 6.620e-06 </td><td> 9.362e-06 </td>
2050 <th>Conv.rate </th><th> 3.00 </th><th> 3.00 </th><th> 3.00 </th>
2054Once more, the convergence rates of @f$\mbox{WG}(Q_2,Q_2;RT_{[2]})@f$ is
2055as expected, with values around 3.
2058<a name=
"step_61-PlainProg"></a>
2059<h1> The plain program</h1>
2060@include
"step-61.cc"
std::vector< bool > component_mask
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
const unsigned int dofs_per_cell
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
const Quadrature< dim > quadrature
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
virtual value_type value(const Point< dim > &p) const
TensorFunction(const time_type initial_time=time_type(0.0))
unsigned int n_active_cells() const
void refine_global(const unsigned int times=1)
unsigned int n_cells() const
typename ActiveSelector::active_cell_iterator active_cell_iterator
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())
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
DataComponentInterpretation
@ component_is_part_of_vector
void approximate(const SynchronousIterators< std::tuple< typename DoFHandler< dim, spacedim >::active_cell_iterator, Vector< float >::iterator > > &cell, const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
SymmetricTensor< 2, dim, Number > E(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr ReturnType< rank, T >::value_type & extract(T &t, const ArrayType &indices)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
void copy(const T *begin, const T *end, U *dest)
int(& functions)(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static constexpr double PI
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation