deal.II version GIT relicensing-1721-g8100761196 2024-08-31 12:30:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
step-31.h
Go to the documentation of this file.
1,
1268 *   const unsigned int /*component*/ = 0) const override
1269 *   {
1270 *   return 0;
1271 *   }
1272 *  
1273 *   virtual void vector_value(const Point<dim> &p,
1274 *   Vector<double> &value) const override
1275 *   {
1276 *   for (unsigned int c = 0; c < this->n_components; ++c)
1277 *   value(c) = TemperatureInitialValues<dim>::value(p, c);
1278 *   }
1279 *   };
1280 *  
1281 *  
1282 *  
1283 *   template <int dim>
1284 *   class TemperatureRightHandSide : public Function<dim>
1285 *   {
1286 *   public:
1287 *   TemperatureRightHandSide()
1288 *   : Function<dim>(1)
1289 *   {}
1290 *  
1291 *   virtual double value(const Point<dim> &p,
1292 *   const unsigned int component = 0) const override
1293 *   {
1294 *   (void)component;
1295 *   Assert(component == 0,
1296 *   ExcMessage("Invalid operation for a scalar function."));
1297 *  
1298 *   Assert((dim == 2) || (dim == 3), ExcNotImplemented());
1299 *  
1300 *   static const Point<dim> source_centers[3] = {
1301 *   (dim == 2 ? Point<dim>(.3, .1) : Point<dim>(.3, .5, .1)),
1302 *   (dim == 2 ? Point<dim>(.45, .1) : Point<dim>(.45, .5, .1)),
1303 *   (dim == 2 ? Point<dim>(.75, .1) : Point<dim>(.75, .5, .1))};
1304 *   static const double source_radius = (dim == 2 ? 1. / 32 : 1. / 8);
1305 *  
1306 *   return ((source_centers[0].distance(p) < source_radius) ||
1307 *   (source_centers[1].distance(p) < source_radius) ||
1308 *   (source_centers[2].distance(p) < source_radius) ?
1309 *   1 :
1310 *   0);
1311 *   }
1312 *  
1313 *   virtual void vector_value(const Point<dim> &p,
1314 *   Vector<double> &value) const override
1315 *   {
1316 *   for (unsigned int c = 0; c < this->n_components; ++c)
1317 *   value(c) = TemperatureRightHandSide<dim>::value(p, c);
1318 *   }
1319 *   };
1320 *   } // namespace EquationData
1321 *  
1322 *  
1323 *  
1324 * @endcode
1325 *
1326 *
1327 * <a name="step_31-Linearsolversandpreconditioners"></a>
1328 * <h3>Linear solvers and preconditioners</h3>
1329 *
1330
1331 *
1332 * This section introduces some objects that are used for the solution of
1333 * the linear equations of the Stokes system that we need to solve in each
1334 * time step. Many of the ideas used here are the same as in @ref step_20 "step-20", where
1335 * Schur complement based preconditioners and solvers have been introduced,
1336 * with the actual interface taken from @ref step_22 "step-22" (in particular the
1337 * discussion in the "Results" section of @ref step_22 "step-22", in which we introduce
1338 * alternatives to the direct Schur complement approach). Note, however,
1339 * that here we don't use the Schur complement to solve the Stokes
1340 * equations, though an approximate Schur complement (the mass matrix on the
1341 * pressure space) appears in the preconditioner.
1342 *
1343 * @code
1344 *   namespace LinearSolvers
1345 *   {
1346 * @endcode
1347 *
1348 *
1349 * <a name="step_31-ThecodeInverseMatrixcodeclasstemplate"></a>
1350 * <h4>The <code>InverseMatrix</code> class template</h4>
1351 *
1352
1353 *
1354 * This class is an interface to calculate the action of an "inverted"
1355 * matrix on a vector (using the <code>vmult</code> operation) in the same
1356 * way as the corresponding class in @ref step_22 "step-22": when the product of an
1357 * object of this class is requested, we solve a linear equation system
1358 * with that matrix using the CG method, accelerated by a preconditioner
1359 * of (templated) class <code>PreconditionerType</code>.
1360 *
1361
1362 *
1363 * In a minor deviation from the implementation of the same class in
1364 * @ref step_22 "step-22", we make the <code>vmult</code> function take any
1365 * kind of vector type (it will yield compiler errors, however, if the
1366 * matrix does not allow a matrix-vector product with this kind of
1367 * vector).
1368 *
1369
1370 *
1371 * Secondly, we catch any exceptions that the solver may have thrown. The
1372 * reason is as follows: When debugging a program like this one
1373 * occasionally makes a mistake of passing an indefinite or nonsymmetric
1374 * matrix or preconditioner to the current class. The solver will, in that
1375 * case, not converge and throw a run-time exception. If not caught here
1376 * it will propagate up the call stack and may end up in
1377 * <code>main()</code> where we output an error message that will say that
1378 * the CG solver failed. The question then becomes: Which CG solver? The
1379 * one that inverted the mass matrix? The one that inverted the top left
1380 * block with the Laplace operator? Or a CG solver in one of the several
1381 * other nested places where we use linear solvers in the current code? No
1382 * indication about this is present in a run-time exception because it
1383 * doesn't store the stack of calls through which we got to the place
1384 * where the exception was generated.
1385 *
1386
1387 *
1388 * So rather than letting the exception propagate freely up to
1389 * <code>main()</code>, we acknowledge that in the current context
1390 * there is little that an outer function can do if the inner
1391 * solver fails. As a consequence, we deal with the situation by
1392 * catching the exception and letting the program fail by
1393 * triggering an assertion with a `false` condition (which of
1394 * course always fails) and that uses the error message associated
1395 * with the caught exception (returned by calling `e.what()`) as
1396 * error text. In other words, instead of letting the error
1397 * message be produced by `main()`, we rather report it here where
1398 * we can abort the program preserving information about where the
1399 * problem happened.
1400 *
1401 * @code
1402 *   template <class MatrixType, class PreconditionerType>
1403 *   class InverseMatrix : public Subscriptor
1404 *   {
1405 *   public:
1406 *   InverseMatrix(const MatrixType &m,
1407 *   const PreconditionerType &preconditioner);
1408 *  
1409 *  
1410 *   template <typename VectorType>
1411 *   void vmult(VectorType &dst, const VectorType &src) const;
1412 *  
1413 *   private:
1415 *   const PreconditionerType &preconditioner;
1416 *   };
1417 *  
1418 *  
1419 *   template <class MatrixType, class PreconditionerType>
1420 *   InverseMatrix<MatrixType, PreconditionerType>::InverseMatrix(
1421 *   const MatrixType &m,
1422 *   const PreconditionerType &preconditioner)
1423 *   : matrix(&m)
1424 *   , preconditioner(preconditioner)
1425 *   {}
1426 *  
1427 *  
1428 *  
1429 *   template <class MatrixType, class PreconditionerType>
1430 *   template <typename VectorType>
1431 *   void InverseMatrix<MatrixType, PreconditionerType>::vmult(
1432 *   VectorType &dst,
1433 *   const VectorType &src) const
1434 *   {
1435 *   SolverControl solver_control(src.size(), 1e-7 * src.l2_norm());
1436 *   SolverCG<VectorType> cg(solver_control);
1437 *  
1438 *   dst = 0;
1439 *  
1440 *   try
1441 *   {
1442 *   cg.solve(*matrix, dst, src, preconditioner);
1443 *   }
1444 *   catch (std::exception &e)
1445 *   {
1446 *   Assert(false, ExcMessage(e.what()));
1447 *   }
1448 *   }
1449 *  
1450 * @endcode
1451 *
1452 *
1453 * <a name="step_31-Schurcomplementpreconditioner"></a>
1454 * <h4>Schur complement preconditioner</h4>
1455 *
1456
1457 *
1458 * This is the implementation of the Schur complement preconditioner as
1459 * described in detail in the introduction. As opposed to @ref step_20 "step-20" and
1460 * @ref step_22 "step-22", we solve the block system all-at-once using GMRES, and use the
1461 * Schur complement of the block structured matrix to build a good
1462 * preconditioner instead.
1463 *
1464
1465 *
1466 * Let's have a look at the ideal preconditioner matrix
1467 * @f$P=\left(\begin{array}{cc} A & 0 \\ B & -S \end{array}\right)@f$
1468 * described in the introduction. If we apply this matrix in the solution
1469 * of a linear system, convergence of an iterative GMRES solver will be
1470 * governed by the matrix @f{eqnarray*} P^{-1}\left(\begin{array}{cc} A &
1471 * B^T \\ B & 0 \end{array}\right) = \left(\begin{array}{cc} I & A^{-1}
1472 * B^T \\ 0 & I \end{array}\right), @f} which indeed is very simple. A
1473 * GMRES solver based on exact matrices would converge in one iteration,
1474 * since all eigenvalues are equal (any Krylov method takes at most as
1475 * many iterations as there are distinct eigenvalues). Such a
1476 * preconditioner for the blocked Stokes system has been proposed by
1477 * Silvester and Wathen ("Fast iterative solution of stabilised Stokes
1478 * systems part II. Using general block preconditioners", SIAM
1479 * J. Numer. Anal., 31 (1994), pp. 1352-1367).
1480 *
1481
1482 *
1483 * Replacing @f$P@f$ by @f$\tilde{P}@f$ keeps that spirit alive: the product
1484 * @f$P^{-1} A@f$ will still be close to a matrix with eigenvalues 1 with a
1485 * distribution that does not depend on the problem size. This lets us
1486 * hope to be able to get a number of GMRES iterations that is
1487 * problem-size independent.
1488 *
1489
1490 *
1491 * The deal.II users who have already gone through the @ref step_20 "step-20" and @ref step_22 "step-22"
1492 * tutorials can certainly imagine how we're going to implement this. We
1493 * replace the exact inverse matrices in @f$P^{-1}@f$ by some approximate
1494 * inverses built from the InverseMatrix class, and the inverse Schur
1495 * complement will be approximated by the pressure mass matrix @f$M_p@f$
1496 * (weighted by @f$\eta^{-1}@f$ as mentioned in the introduction). As pointed
1497 * out in the results section of @ref step_22 "step-22", we can replace the exact inverse
1498 * of @f$A@f$ by just the application of a preconditioner, in this case
1499 * on a vector Laplace matrix as was explained in the introduction. This
1500 * does increase the number of (outer) GMRES iterations, but is still
1501 * significantly cheaper than an exact inverse, which would require
1502 * between 20 and 35 CG iterations for <em>each</em> outer solver step
1503 * (using the AMG preconditioner).
1504 *
1505
1506 *
1507 * Having the above explanations in mind, we define a preconditioner class
1508 * with a <code>vmult</code> functionality, which is all we need for the
1509 * interaction with the usual solver functions further below in the
1510 * program code.
1511 *
1512
1513 *
1514 * First the declarations. These are similar to the definition of the
1515 * Schur complement in @ref step_20 "step-20", with the difference that we need some more
1516 * preconditioners in the constructor and that the matrices we use here
1517 * are built upon Trilinos:
1518 *
1519 * @code
1520 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1521 *   class BlockSchurPreconditioner : public Subscriptor
1522 *   {
1523 *   public:
1524 *   BlockSchurPreconditioner(
1525 *   const TrilinosWrappers::BlockSparseMatrix &S,
1526 *   const InverseMatrix<TrilinosWrappers::SparseMatrix,
1527 *   PreconditionerTypeMp> &Mpinv,
1528 *   const PreconditionerTypeA &Apreconditioner);
1529 *  
1530 *   void vmult(TrilinosWrappers::MPI::BlockVector &dst,
1531 *   const TrilinosWrappers::MPI::BlockVector &src) const;
1532 *  
1533 *   private:
1535 *   stokes_matrix;
1536 *   const SmartPointer<const InverseMatrix<TrilinosWrappers::SparseMatrix,
1537 *   PreconditionerTypeMp>>
1538 *   m_inverse;
1539 *   const PreconditionerTypeA &a_preconditioner;
1540 *  
1541 *   mutable TrilinosWrappers::MPI::Vector tmp;
1542 *   };
1543 *  
1544 *  
1545 *  
1546 * @endcode
1547 *
1548 * When using a TrilinosWrappers::MPI::Vector or a
1549 * TrilinosWrappers::MPI::BlockVector, the Vector is initialized using an
1550 * IndexSet. IndexSet is used not only to resize the
1551 * TrilinosWrappers::MPI::Vector but it also associates an index in the
1552 * TrilinosWrappers::MPI::Vector with a degree of freedom (see @ref step_40 "step-40" for
1553 * a more detailed explanation). The function complete_index_set() creates
1554 * an IndexSet where every valid index is part of the set. Note that this
1555 * program can only be run sequentially and will throw an exception if used
1556 * in parallel.
1557 *
1558 * @code
1559 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1560 *   BlockSchurPreconditioner<PreconditionerTypeA, PreconditionerTypeMp>::
1561 *   BlockSchurPreconditioner(
1562 *   const TrilinosWrappers::BlockSparseMatrix &S,
1563 *   const InverseMatrix<TrilinosWrappers::SparseMatrix,
1564 *   PreconditionerTypeMp> &Mpinv,
1565 *   const PreconditionerTypeA &Apreconditioner)
1566 *   : stokes_matrix(&S)
1567 *   , m_inverse(&Mpinv)
1568 *   , a_preconditioner(Apreconditioner)
1569 *   , tmp(complete_index_set(stokes_matrix->block(1, 1).m()))
1570 *   {}
1571 *  
1572 *  
1573 * @endcode
1574 *
1575 * Next is the <code>vmult</code> function. We implement the action of
1576 * @f$P^{-1}@f$ as described above in three successive steps. In formulas, we
1577 * want to compute @f$Y=P^{-1}X@f$ where @f$X,Y@f$ are both vectors with two block
1578 * components.
1579 *
1580
1581 *
1582 * The first step multiplies the velocity part of the vector by a
1583 * preconditioner of the matrix @f$A@f$, i.e., we compute @f$Y_0={\tilde
1584 * A}^{-1}X_0@f$. The resulting velocity vector is then multiplied by @f$B@f$
1585 * and subtracted from the pressure, i.e., we want to compute @f$X_1-BY_0@f$.
1586 * This second step only acts on the pressure vector and is accomplished
1587 * by the residual function of our matrix classes, except that the sign is
1588 * wrong. Consequently, we change the sign in the temporary pressure
1589 * vector and finally multiply by the inverse pressure mass matrix to get
1590 * the final pressure vector, completing our work on the Stokes
1591 * preconditioner:
1592 *
1593 * @code
1594 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1595 *   void
1596 *   BlockSchurPreconditioner<PreconditionerTypeA, PreconditionerTypeMp>::vmult(
1598 *   const TrilinosWrappers::MPI::BlockVector &src) const
1599 *   {
1600 *   a_preconditioner.vmult(dst.block(0), src.block(0));
1601 *   stokes_matrix->block(1, 0).residual(tmp, dst.block(0), src.block(1));
1602 *   tmp *= -1;
1603 *   m_inverse->vmult(dst.block(1), tmp);
1604 *   }
1605 *   } // namespace LinearSolvers
1606 *  
1607 *  
1608 *  
1609 * @endcode
1610 *
1611 *
1612 * <a name="step_31-ThecodeBoussinesqFlowProblemcodeclasstemplate"></a>
1613 * <h3>The <code>BoussinesqFlowProblem</code> class template</h3>
1614 *
1615
1616 *
1617 * The definition of the class that defines the top-level logic of solving
1618 * the time-dependent Boussinesq problem is mainly based on the @ref step_22 "step-22"
1619 * tutorial program. The main differences are that now we also have to solve
1620 * for the temperature equation, which forces us to have a second DoFHandler
1621 * object for the temperature variable as well as matrices, right hand
1622 * sides, and solution vectors for the current and previous time steps. As
1623 * mentioned in the introduction, all linear algebra objects are going to
1624 * use wrappers of the corresponding Trilinos functionality.
1625 *
1626
1627 *
1628 * The member functions of this class are reminiscent of @ref step_21 "step-21", where we
1629 * also used a staggered scheme that first solve the flow equations (here
1630 * the Stokes equations, in @ref step_21 "step-21" Darcy flow) and then update the advected
1631 * quantity (here the temperature, there the saturation). The functions that
1632 * are new are mainly concerned with determining the time step, as well as
1633 * the proper size of the artificial viscosity stabilization.
1634 *
1635
1636 *
1637 * The last three variables indicate whether the various matrices or
1638 * preconditioners need to be rebuilt the next time the corresponding build
1639 * functions are called. This allows us to move the corresponding
1640 * <code>if</code> into the respective function and thereby keeping our main
1641 * <code>run()</code> function clean and easy to read.
1642 *
1643 * @code
1644 *   template <int dim>
1645 *   class BoussinesqFlowProblem
1646 *   {
1647 *   public:
1648 *   BoussinesqFlowProblem();
1649 *   void run();
1650 *  
1651 *   private:
1652 *   void setup_dofs();
1653 *   void assemble_stokes_preconditioner();
1654 *   void build_stokes_preconditioner();
1655 *   void assemble_stokes_system();
1656 *   void assemble_temperature_system(const double maximal_velocity);
1657 *   void assemble_temperature_matrix();
1658 *   double get_maximal_velocity() const;
1659 *   std::pair<double, double> get_extrapolated_temperature_range() const;
1660 *   void solve();
1661 *   void output_results() const;
1662 *   void refine_mesh(const unsigned int max_grid_level);
1663 *  
1664 *   double compute_viscosity(
1665 *   const std::vector<double> &old_temperature,
1666 *   const std::vector<double> &old_old_temperature,
1667 *   const std::vector<Tensor<1, dim>> &old_temperature_grads,
1668 *   const std::vector<Tensor<1, dim>> &old_old_temperature_grads,
1669 *   const std::vector<double> &old_temperature_laplacians,
1670 *   const std::vector<double> &old_old_temperature_laplacians,
1671 *   const std::vector<Tensor<1, dim>> &old_velocity_values,
1672 *   const std::vector<Tensor<1, dim>> &old_old_velocity_values,
1673 *   const std::vector<double> &gamma_values,
1674 *   const double global_u_infty,
1675 *   const double global_T_variation,
1676 *   const double cell_diameter) const;
1677 *  
1678 *  
1680 *   double global_Omega_diameter;
1681 *  
1682 *   const unsigned int stokes_degree;
1683 *   const FESystem<dim> stokes_fe;
1684 *   DoFHandler<dim> stokes_dof_handler;
1685 *   AffineConstraints<double> stokes_constraints;
1686 *  
1687 *   std::vector<IndexSet> stokes_partitioning;
1688 *   TrilinosWrappers::BlockSparseMatrix stokes_matrix;
1689 *   TrilinosWrappers::BlockSparseMatrix stokes_preconditioner_matrix;
1690 *  
1691 *   TrilinosWrappers::MPI::BlockVector stokes_solution;
1692 *   TrilinosWrappers::MPI::BlockVector old_stokes_solution;
1693 *   TrilinosWrappers::MPI::BlockVector stokes_rhs;
1694 *  
1695 *  
1696 *   const unsigned int temperature_degree;
1697 *   const FE_Q<dim> temperature_fe;
1698 *   DoFHandler<dim> temperature_dof_handler;
1699 *   AffineConstraints<double> temperature_constraints;
1700 *  
1701 *   TrilinosWrappers::SparseMatrix temperature_mass_matrix;
1702 *   TrilinosWrappers::SparseMatrix temperature_stiffness_matrix;
1703 *   TrilinosWrappers::SparseMatrix temperature_matrix;
1704 *  
1705 *   TrilinosWrappers::MPI::Vector temperature_solution;
1706 *   TrilinosWrappers::MPI::Vector old_temperature_solution;
1707 *   TrilinosWrappers::MPI::Vector old_old_temperature_solution;
1708 *   TrilinosWrappers::MPI::Vector temperature_rhs;
1709 *  
1710 *  
1711 *   double time_step;
1712 *   double old_time_step;
1713 *   unsigned int timestep_number;
1714 *  
1715 *   std::shared_ptr<TrilinosWrappers::PreconditionAMG> Amg_preconditioner;
1716 *   std::shared_ptr<TrilinosWrappers::PreconditionIC> Mp_preconditioner;
1717 *  
1718 *   bool rebuild_stokes_matrix;
1719 *   bool rebuild_temperature_matrices;
1720 *   bool rebuild_stokes_preconditioner;
1721 *   };
1722 *  
1723 *  
1724 * @endcode
1725 *
1726 *
1727 * <a name="step_31-BoussinesqFlowProblemclassimplementation"></a>
1728 * <h3>BoussinesqFlowProblem class implementation</h3>
1729 *
1730
1731 *
1732 *
1733 * <a name="step_31-BoussinesqFlowProblemBoussinesqFlowProblem"></a>
1734 * <h4>BoussinesqFlowProblem::BoussinesqFlowProblem</h4>
1735 *
1736
1737 *
1738 * The constructor of this class is an extension of the constructor in
1739 * @ref step_22 "step-22". We need to add the various variables that concern the
1740 * temperature. As discussed in the introduction, we are going to use
1741 * @f$Q_2^d\times Q_1@f$ (Taylor-Hood) elements again for the Stokes part, and
1742 * @f$Q_2@f$ elements for the temperature. However, by using variables that
1743 * store the polynomial degree of the Stokes and temperature finite
1744 * elements, it is easy to consistently modify the degree of the elements as
1745 * well as all quadrature formulas used on them downstream. Moreover, we
1746 * initialize the time stepping as well as the options for matrix assembly
1747 * and preconditioning:
1748 *
1749 * @code
1750 *   template <int dim>
1751 *   BoussinesqFlowProblem<dim>::BoussinesqFlowProblem()
1753 *   , global_Omega_diameter(std::numeric_limits<double>::quiet_NaN())
1754 *   , stokes_degree(1)
1755 *   , stokes_fe(FE_Q<dim>(stokes_degree + 1) ^ dim, FE_Q<dim>(stokes_degree))
1756 *   , stokes_dof_handler(triangulation)
1757 *   ,
1758 *  
1759 *   temperature_degree(2)
1760 *   , temperature_fe(temperature_degree)
1761 *   , temperature_dof_handler(triangulation)
1762 *   ,
1763 *  
1764 *   time_step(0)
1765 *   , old_time_step(0)
1766 *   , timestep_number(0)
1767 *   , rebuild_stokes_matrix(true)
1768 *   , rebuild_temperature_matrices(true)
1769 *   , rebuild_stokes_preconditioner(true)
1770 *   {}
1771 *  
1772 *  
1773 *  
1774 * @endcode
1775 *
1776 *
1777 * <a name="step_31-BoussinesqFlowProblemget_maximal_velocity"></a>
1778 * <h4>BoussinesqFlowProblem::get_maximal_velocity</h4>
1779 *
1780
1781 *
1782 * Starting the real functionality of this class is a helper function that
1783 * determines the maximum (@f$L_\infty@f$) velocity in the domain (at the
1784 * quadrature points, in fact). How it works should be relatively obvious to
1785 * all who have gotten to this point of the tutorial. Note that since we are
1786 * only interested in the velocity, rather than using
1787 * <code>stokes_fe_values.get_function_values</code> to get the values of
1788 * the entire Stokes solution (velocities and pressures) we use
1789 * <code>stokes_fe_values[velocities].get_function_values</code> to extract
1790 * only the velocities part. This has the additional benefit that we get it
1791 * as a Tensor<1,dim>, rather than some components in a Vector<double>,
1792 * allowing us to process it right away using the <code>norm()</code>
1793 * function to get the magnitude of the velocity.
1794 *
1795
1796 *
1797 * The only point worth thinking about a bit is how to choose the quadrature
1798 * points we use here. Since the goal of this function is to find the
1799 * maximal velocity over a domain by looking at quadrature points on each
1800 * cell. So we should ask how we should best choose these quadrature points
1801 * on each cell. To this end, recall that if we had a single @f$Q_1@f$ field
1802 * (rather than the vector-valued field of higher order) then the maximum
1803 * would be attained at a vertex of the mesh. In other words, we should use
1804 * the QTrapezoid class that has quadrature points only at the vertices of
1805 * cells.
1806 *
1807
1808 *
1809 * For higher order shape functions, the situation is more complicated: the
1810 * maxima and minima may be attained at points between the support points of
1811 * shape functions (for the usual @f$Q_p@f$ elements the support points are the
1812 * equidistant Lagrange interpolation points); furthermore, since we are
1813 * looking for the maximum magnitude of a vector-valued quantity, we can
1814 * even less say with certainty where the set of potential maximal points
1815 * are. Nevertheless, intuitively if not provably, the Lagrange
1816 * interpolation points appear to be a better choice than the Gauss points.
1817 *
1818
1819 *
1820 * There are now different methods to produce a quadrature formula with
1821 * quadrature points equal to the interpolation points of the finite
1822 * element. One option would be to use the
1823 * FiniteElement::get_unit_support_points() function, reduce the output to a
1824 * unique set of points to avoid duplicate function evaluations, and create
1825 * a Quadrature object using these points. Another option, chosen here, is
1826 * to use the QTrapezoid class and combine it with the QIterated class that
1827 * repeats the QTrapezoid formula on a number of sub-cells in each coordinate
1828 * direction. To cover all support points, we need to iterate it
1829 * <code>stokes_degree+1</code> times since this is the polynomial degree of
1830 * the Stokes element in use:
1831 *
1832 * @code
1833 *   template <int dim>
1834 *   double BoussinesqFlowProblem<dim>::get_maximal_velocity() const
1835 *   {
1836 *   const QIterated<dim> quadrature_formula(QTrapezoid<1>(), stokes_degree + 1);
1837 *   const unsigned int n_q_points = quadrature_formula.size();
1838 *  
1839 *   FEValues<dim> fe_values(stokes_fe, quadrature_formula, update_values);
1840 *   std::vector<Tensor<1, dim>> velocity_values(n_q_points);
1841 *   double max_velocity = 0;
1842 *  
1843 *   const FEValuesExtractors::Vector velocities(0);
1844 *  
1845 *   for (const auto &cell : stokes_dof_handler.active_cell_iterators())
1846 *   {
1847 *   fe_values.reinit(cell);
1848 *   fe_values[velocities].get_function_values(stokes_solution,
1849 *   velocity_values);
1850 *  
1851 *   for (unsigned int q = 0; q < n_q_points; ++q)
1852 *   max_velocity = std::max(max_velocity, velocity_values[q].norm());
1853 *   }
1854 *  
1855 *   return max_velocity;
1856 *   }
1857 *  
1858 *  
1859 *  
1860 * @endcode
1861 *
1862 *
1863 * <a name="step_31-BoussinesqFlowProblemget_extrapolated_temperature_range"></a>
1864 * <h4>BoussinesqFlowProblem::get_extrapolated_temperature_range</h4>
1865 *
1866
1867 *
1868 * Next a function that determines the minimum and maximum temperature at
1869 * quadrature points inside @f$\Omega@f$ when extrapolated from the two previous
1870 * time steps to the current one. We need this information in the
1871 * computation of the artificial viscosity parameter @f$\nu@f$ as discussed in
1872 * the introduction.
1873 *
1874
1875 *
1876 * The formula for the extrapolated temperature is
1877 * @f$\left(1+\frac{k_n}{k_{n-1}} \right)T^{n-1} + \frac{k_n}{k_{n-1}}
1878 * T^{n-2}@f$. The way to compute it is to loop over all quadrature points and
1879 * update the maximum and minimum value if the current value is
1880 * bigger/smaller than the previous one. We initialize the variables that
1881 * store the max and min before the loop over all quadrature points by the
1882 * smallest and the largest number representable as a double. Then we know
1883 * for a fact that it is larger/smaller than the minimum/maximum and that
1884 * the loop over all quadrature points is ultimately going to update the
1885 * initial value with the correct one.
1886 *
1887
1888 *
1889 * The only other complication worth mentioning here is that in the first
1890 * time step, @f$T^{k-2}@f$ is not yet available of course. In that case, we can
1891 * only use @f$T^{k-1}@f$ which we have from the initial temperature. As
1892 * quadrature points, we use the same choice as in the previous function
1893 * though with the difference that now the number of repetitions is
1894 * determined by the polynomial degree of the temperature field.
1895 *
1896 * @code
1897 *   template <int dim>
1898 *   std::pair<double, double>
1899 *   BoussinesqFlowProblem<dim>::get_extrapolated_temperature_range() const
1900 *   {
1901 *   const QIterated<dim> quadrature_formula(QTrapezoid<1>(),
1902 *   temperature_degree);
1903 *   const unsigned int n_q_points = quadrature_formula.size();
1904 *  
1905 *   FEValues<dim> fe_values(temperature_fe, quadrature_formula, update_values);
1906 *   std::vector<double> old_temperature_values(n_q_points);
1907 *   std::vector<double> old_old_temperature_values(n_q_points);
1908 *  
1909 *   if (timestep_number != 0)
1910 *   {
1911 *   double min_temperature = std::numeric_limits<double>::max(),
1912 *   max_temperature = -std::numeric_limits<double>::max();
1913 *  
1914 *   for (const auto &cell : temperature_dof_handler.active_cell_iterators())
1915 *   {
1916 *   fe_values.reinit(cell);
1917 *   fe_values.get_function_values(old_temperature_solution,
1918 *   old_temperature_values);
1919 *   fe_values.get_function_values(old_old_temperature_solution,
1920 *   old_old_temperature_values);
1921 *  
1922 *   for (unsigned int q = 0; q < n_q_points; ++q)
1923 *   {
1924 *   const double temperature =
1925 *   (1. + time_step / old_time_step) * old_temperature_values[q] -
1926 *   time_step / old_time_step * old_old_temperature_values[q];
1927 *  
1928 *   min_temperature = std::min(min_temperature, temperature);
1929 *   max_temperature = std::max(max_temperature, temperature);
1930 *   }
1931 *   }
1932 *  
1933 *   return std::make_pair(min_temperature, max_temperature);
1934 *   }
1935 *   else
1936 *   {
1937 *   double min_temperature = std::numeric_limits<double>::max(),
1938 *   max_temperature = -std::numeric_limits<double>::max();
1939 *  
1940 *   for (const auto &cell : temperature_dof_handler.active_cell_iterators())
1941 *   {
1942 *   fe_values.reinit(cell);
1943 *   fe_values.get_function_values(old_temperature_solution,
1944 *   old_temperature_values);
1945 *  
1946 *   for (unsigned int q = 0; q < n_q_points; ++q)
1947 *   {
1948 *   const double temperature = old_temperature_values[q];
1949 *  
1950 *   min_temperature = std::min(min_temperature, temperature);
1951 *   max_temperature = std::max(max_temperature, temperature);
1952 *   }
1953 *   }
1954 *  
1955 *   return std::make_pair(min_temperature, max_temperature);
1956 *   }
1957 *   }
1958 *  
1959 *  
1960 *  
1961 * @endcode
1962 *
1963 *
1964 * <a name="step_31-BoussinesqFlowProblemcompute_viscosity"></a>
1965 * <h4>BoussinesqFlowProblem::compute_viscosity</h4>
1966 *
1967
1968 *
1969 * The last of the tool functions computes the artificial viscosity
1970 * parameter @f$\nu|_K@f$ on a cell @f$K@f$ as a function of the extrapolated
1971 * temperature, its gradient and Hessian (second derivatives), the velocity,
1972 * the right hand side @f$\gamma@f$ all on the quadrature points of the current
1973 * cell, and various other parameters as described in detail in the
1974 * introduction.
1975 *
1976
1977 *
1978 * There are some universal constants worth mentioning here. First, we need
1979 * to fix @f$\beta@f$; we choose @f$\beta=0.017\cdot dim@f$, a choice discussed in
1980 * detail in the results section of this tutorial program. The second is the
1981 * exponent @f$\alpha@f$; @f$\alpha=1@f$ appears to work fine for the current
1982 * program, even though some additional benefit might be expected from
1983 * choosing @f$\alpha = 2@f$. Finally, there is one thing that requires special
1984 * casing: In the first time step, the velocity equals zero, and the formula
1985 * for @f$\nu|_K@f$ is not defined. In that case, we return @f$\nu|_K=5\cdot 10^3
1986 * \cdot h_K@f$, a choice admittedly more motivated by heuristics than
1987 * anything else (it is in the same order of magnitude, however, as the
1988 * value returned for most cells on the second time step).
1989 *
1990
1991 *
1992 * The rest of the function should be mostly obvious based on the material
1993 * discussed in the introduction:
1994 *
1995 * @code
1996 *   template <int dim>
1997 *   double BoussinesqFlowProblem<dim>::compute_viscosity(
1998 *   const std::vector<double> &old_temperature,
1999 *   const std::vector<double> &old_old_temperature,
2000 *   const std::vector<Tensor<1, dim>> &old_temperature_grads,
2001 *   const std::vector<Tensor<1, dim>> &old_old_temperature_grads,
2002 *   const std::vector<double> &old_temperature_laplacians,
2003 *   const std::vector<double> &old_old_temperature_laplacians,
2004 *   const std::vector<Tensor<1, dim>> &old_velocity_values,
2005 *   const std::vector<Tensor<1, dim>> &old_old_velocity_values,
2006 *   const std::vector<double> &gamma_values,
2007 *   const double global_u_infty,
2008 *   const double global_T_variation,
2009 *   const double cell_diameter) const
2010 *   {
2011 *   constexpr double beta = 0.017 * dim;
2012 *   constexpr double alpha = 1.0;
2013 *  
2014 *   if (global_u_infty == 0)
2015 *   return 5e-3 * cell_diameter;
2016 *  
2017 *   const unsigned int n_q_points = old_temperature.size();
2018 *  
2019 *   double max_residual = 0;
2020 *   double max_velocity = 0;
2021 *  
2022 *   for (unsigned int q = 0; q < n_q_points; ++q)
2023 *   {
2024 *   const Tensor<1, dim> u =
2025 *   (old_velocity_values[q] + old_old_velocity_values[q]) / 2;
2026 *  
2027 *   const double dT_dt =
2028 *   (old_temperature[q] - old_old_temperature[q]) / old_time_step;
2029 *   const double u_grad_T =
2030 *   u * (old_temperature_grads[q] + old_old_temperature_grads[q]) / 2;
2031 *  
2032 *   const double kappa_Delta_T =
2033 *   EquationData::kappa *
2034 *   (old_temperature_laplacians[q] + old_old_temperature_laplacians[q]) /
2035 *   2;
2036 *  
2037 *   const double residual =
2038 *   std::abs((dT_dt + u_grad_T - kappa_Delta_T - gamma_values[q]) *
2039 *   std::pow((old_temperature[q] + old_old_temperature[q]) / 2,
2040 *   alpha - 1.));
2041 *  
2042 *   max_residual = std::max(residual, max_residual);
2043 *   max_velocity = std::max(std::sqrt(u * u), max_velocity);
2044 *   }
2045 *  
2046 *   const double c_R = std::pow(2., (4. - 2 * alpha) / dim);
2047 *   const double global_scaling = c_R * global_u_infty * global_T_variation *
2048 *   std::pow(global_Omega_diameter, alpha - 2.);
2049 *  
2050 *   return (
2051 *   beta * max_velocity *
2052 *   std::min(cell_diameter,
2053 *   std::pow(cell_diameter, alpha) * max_residual / global_scaling));
2054 *   }
2055 *  
2056 *  
2057 *  
2058 * @endcode
2059 *
2060 *
2061 * <a name="step_31-BoussinesqFlowProblemsetup_dofs"></a>
2062 * <h4>BoussinesqFlowProblem::setup_dofs</h4>
2063 *
2064
2065 *
2066 * This is the function that sets up the DoFHandler objects we have here
2067 * (one for the Stokes part and one for the temperature part) as well as set
2068 * to the right sizes the various objects required for the linear algebra in
2069 * this program. Its basic operations are similar to what we do in @ref step_22 "step-22".
2070 *
2071
2072 *
2073 * The body of the function first enumerates all degrees of freedom for the
2074 * Stokes and temperature systems. For the Stokes part, degrees of freedom
2075 * are then sorted to ensure that velocities precede pressure DoFs so that
2076 * we can partition the Stokes matrix into a @f$2\times 2@f$ matrix. As a
2077 * difference to @ref step_22 "step-22", we do not perform any additional DoF
2078 * renumbering. In that program, it paid off since our solver was heavily
2079 * dependent on ILU's, whereas we use AMG here which is not sensitive to the
2080 * DoF numbering. The IC preconditioner for the inversion of the pressure
2081 * mass matrix would of course take advantage of a Cuthill-McKee like
2082 * renumbering, but its costs are low compared to the velocity portion, so
2083 * the additional work does not pay off.
2084 *
2085
2086 *
2087 * We then proceed with the generation of the hanging node constraints that
2088 * arise from adaptive grid refinement for both DoFHandler objects. For the
2089 * velocity, we impose no-flux boundary conditions @f$\mathbf{u}\cdot
2090 * \mathbf{n}=0@f$ by adding constraints to the object that already stores the
2091 * hanging node constraints matrix. The second parameter in the function
2092 * describes the first of the velocity components in the total dof vector,
2093 * which is zero here. The variable <code>no_normal_flux_boundaries</code>
2094 * denotes the boundary indicators for which to set the no flux boundary
2095 * conditions; here, this is boundary indicator zero.
2096 *
2097
2098 *
2099 * After having done so, we count the number of degrees of freedom in the
2100 * various blocks:
2101 *
2102 * @code
2103 *   template <int dim>
2104 *   void BoussinesqFlowProblem<dim>::setup_dofs()
2105 *   {
2106 *   std::vector<unsigned int> stokes_sub_blocks(dim + 1, 0);
2107 *   stokes_sub_blocks[dim] = 1;
2108 *  
2109 *   {
2110 *   stokes_dof_handler.distribute_dofs(stokes_fe);
2111 *   DoFRenumbering::component_wise(stokes_dof_handler, stokes_sub_blocks);
2112 *  
2113 *   stokes_constraints.clear();
2114 *   DoFTools::make_hanging_node_constraints(stokes_dof_handler,
2115 *   stokes_constraints);
2116 *   const std::set<types::boundary_id> no_normal_flux_boundaries = {0};
2117 *   VectorTools::compute_no_normal_flux_constraints(stokes_dof_handler,
2118 *   0,
2119 *   no_normal_flux_boundaries,
2120 *   stokes_constraints);
2121 *   stokes_constraints.close();
2122 *   }
2123 *   {
2124 *   temperature_dof_handler.distribute_dofs(temperature_fe);
2125 *  
2126 *   temperature_constraints.clear();
2127 *   DoFTools::make_hanging_node_constraints(temperature_dof_handler,
2128 *   temperature_constraints);
2129 *   temperature_constraints.close();
2130 *   }
2131 *  
2132 *   const std::vector<types::global_dof_index> stokes_dofs_per_block =
2133 *   DoFTools::count_dofs_per_fe_block(stokes_dof_handler, stokes_sub_blocks);
2134 *  
2135 *   const types::global_dof_index n_u = stokes_dofs_per_block[0],
2136 *   n_p = stokes_dofs_per_block[1],
2137 *   n_T = temperature_dof_handler.n_dofs();
2138 *  
2139 *   std::cout << "Number of active cells: " << triangulation.n_active_cells()
2140 *   << " (on " << triangulation.n_levels() << " levels)" << std::endl
2141 *   << "Number of degrees of freedom: " << n_u + n_p + n_T << " ("
2142 *   << n_u << '+' << n_p << '+' << n_T << ')' << std::endl
2143 *   << std::endl;
2144 *  
2145 * @endcode
2146 *
2147 * The next step is to create the sparsity pattern for the Stokes and
2148 * temperature system matrices as well as the preconditioner matrix from
2149 * which we build the Stokes preconditioner. As in @ref step_22 "step-22", we choose to
2150 * create the pattern by
2151 * using the blocked version of DynamicSparsityPattern.
2152 *
2153
2154 *
2155 * So, we first release the memory stored in the matrices, then set up an
2156 * object of type BlockDynamicSparsityPattern consisting of
2157 * @f$2\times 2@f$ blocks (for the Stokes system matrix and preconditioner) or
2158 * DynamicSparsityPattern (for the temperature part). We then
2159 * fill these objects with the nonzero pattern, taking into account that
2160 * for the Stokes system matrix, there are no entries in the
2161 * pressure-pressure block (but all velocity vector components couple with
2162 * each other and with the pressure). Similarly, in the Stokes
2163 * preconditioner matrix, only the diagonal blocks are nonzero, since we
2164 * use the vector Laplacian as discussed in the introduction. This
2165 * operator only couples each vector component of the Laplacian with
2166 * itself, but not with the other vector components. (Application of the
2167 * constraints resulting from the no-flux boundary conditions will couple
2168 * vector components at the boundary again, however.)
2169 *
2170
2171 *
2172 * When generating the sparsity pattern, we directly apply the constraints
2173 * from hanging nodes and no-flux boundary conditions. This approach was
2174 * already used in @ref step_27 "step-27", but is different from the one in early
2175 * tutorial programs where we first built the original sparsity pattern
2176 * and only then added the entries resulting from constraints. The reason
2177 * for doing so is that later during assembly we are going to distribute
2178 * the constraints immediately when transferring local to global
2179 * dofs. Consequently, there will be no data written at positions of
2180 * constrained degrees of freedom, so we can let the
2181 * DoFTools::make_sparsity_pattern function omit these entries by setting
2182 * the last Boolean flag to <code>false</code>. Once the sparsity pattern
2183 * is ready, we can use it to initialize the Trilinos matrices. Since the
2184 * Trilinos matrices store the sparsity pattern internally, there is no
2185 * need to keep the sparsity pattern around after the initialization of
2186 * the matrix.
2187 *
2188 * @code
2189 *   stokes_partitioning.resize(2);
2190 *   stokes_partitioning[0] = complete_index_set(n_u);
2191 *   stokes_partitioning[1] = complete_index_set(n_p);
2192 *   {
2193 *   stokes_matrix.clear();
2194 *  
2195 *   BlockDynamicSparsityPattern dsp(stokes_dofs_per_block,
2196 *   stokes_dofs_per_block);
2197 *  
2198 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
2199 *  
2200 *   for (unsigned int c = 0; c < dim + 1; ++c)
2201 *   for (unsigned int d = 0; d < dim + 1; ++d)
2202 *   if (!((c == dim) && (d == dim)))
2203 *   coupling[c][d] = DoFTools::always;
2204 *   else
2205 *   coupling[c][d] = DoFTools::none;
2206 *  
2207 *   DoFTools::make_sparsity_pattern(
2208 *   stokes_dof_handler, coupling, dsp, stokes_constraints, false);
2209 *  
2210 *   stokes_matrix.reinit(dsp);
2211 *   }
2212 *  
2213 *   {
2214 *   Amg_preconditioner.reset();
2215 *   Mp_preconditioner.reset();
2216 *   stokes_preconditioner_matrix.clear();
2217 *  
2218 *   BlockDynamicSparsityPattern dsp(stokes_dofs_per_block,
2219 *   stokes_dofs_per_block);
2220 *  
2221 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
2222 *   for (unsigned int c = 0; c < dim + 1; ++c)
2223 *   for (unsigned int d = 0; d < dim + 1; ++d)
2224 *   if (c == d)
2225 *   coupling[c][d] = DoFTools::always;
2226 *   else
2227 *   coupling[c][d] = DoFTools::none;
2228 *  
2229 *   DoFTools::make_sparsity_pattern(
2230 *   stokes_dof_handler, coupling, dsp, stokes_constraints, false);
2231 *  
2232 *   stokes_preconditioner_matrix.reinit(dsp);
2233 *   }
2234 *  
2235 * @endcode
2236 *
2237 * The creation of the temperature matrix (or, rather, matrices, since we
2238 * provide a temperature mass matrix and a temperature @ref GlossStiffnessMatrix "stiffness matrix",
2239 * that will be added together for time discretization) follows the
2240 * generation of the Stokes matrix &ndash; except that it is much easier
2241 * here since we do not need to take care of any blocks or coupling
2242 * between components. Note how we initialize the three temperature
2243 * matrices: We only use the sparsity pattern for reinitialization of the
2244 * first matrix, whereas we use the previously generated matrix for the
2245 * two remaining reinits. The reason for doing so is that reinitialization
2246 * from an already generated matrix allows Trilinos to reuse the sparsity
2247 * pattern instead of generating a new one for each copy. This saves both
2248 * some time and memory.
2249 *
2250 * @code
2251 *   {
2252 *   temperature_mass_matrix.clear();
2253 *   temperature_stiffness_matrix.clear();
2254 *   temperature_matrix.clear();
2255 *  
2256 *   DynamicSparsityPattern dsp(n_T, n_T);
2257 *   DoFTools::make_sparsity_pattern(temperature_dof_handler,
2258 *   dsp,
2259 *   temperature_constraints,
2260 *   false);
2261 *  
2262 *   temperature_matrix.reinit(dsp);
2263 *   temperature_mass_matrix.reinit(temperature_matrix);
2264 *   temperature_stiffness_matrix.reinit(temperature_matrix);
2265 *   }
2266 *  
2267 * @endcode
2268 *
2269 * Lastly, we set the vectors for the Stokes solutions @f$\mathbf u^{n-1}@f$
2270 * and @f$\mathbf u^{n-2}@f$, as well as for the temperatures @f$T^{n}@f$,
2271 * @f$T^{n-1}@f$ and @f$T^{n-2}@f$ (required for time stepping) and all the system
2272 * right hand sides to their correct sizes and block structure:
2273 *
2274 * @code
2275 *   IndexSet temperature_partitioning = complete_index_set(n_T);
2276 *   stokes_solution.reinit(stokes_partitioning, MPI_COMM_WORLD);
2277 *   old_stokes_solution.reinit(stokes_partitioning, MPI_COMM_WORLD);
2278 *   stokes_rhs.reinit(stokes_partitioning, MPI_COMM_WORLD);
2279 *  
2280 *   temperature_solution.reinit(temperature_partitioning, MPI_COMM_WORLD);
2281 *   old_temperature_solution.reinit(temperature_partitioning, MPI_COMM_WORLD);
2282 *   old_old_temperature_solution.reinit(temperature_partitioning,
2283 *   MPI_COMM_WORLD);
2284 *  
2285 *   temperature_rhs.reinit(temperature_partitioning, MPI_COMM_WORLD);
2286 *   }
2287 *  
2288 *  
2289 *  
2290 * @endcode
2291 *
2292 *
2293 * <a name="step_31-BoussinesqFlowProblemassemble_stokes_preconditioner"></a>
2294 * <h4>BoussinesqFlowProblem::assemble_stokes_preconditioner</h4>
2295 *
2296
2297 *
2298 * This function assembles the matrix we use for preconditioning the Stokes
2299 * system. What we need are a vector Laplace matrix on the velocity
2300 * components and a mass matrix weighted by @f$\eta^{-1}@f$ on the pressure
2301 * component. We start by generating a quadrature object of appropriate
2302 * order, the FEValues object that can give values and gradients at the
2303 * quadrature points (together with quadrature weights). Next we create data
2304 * structures for the cell matrix and the relation between local and global
2305 * DoFs. The vectors <code>grad_phi_u</code> and <code>phi_p</code> are
2306 * going to hold the values of the basis functions in order to faster build
2307 * up the local matrices, as was already done in @ref step_22 "step-22". Before we start
2308 * the loop over all active cells, we have to specify which components are
2309 * pressure and which are velocity.
2310 *
2311 * @code
2312 *   template <int dim>
2313 *   void BoussinesqFlowProblem<dim>::assemble_stokes_preconditioner()
2314 *   {
2315 *   stokes_preconditioner_matrix = 0;
2316 *  
2317 *   const QGauss<dim> quadrature_formula(stokes_degree + 2);
2318 *   FEValues<dim> stokes_fe_values(stokes_fe,
2319 *   quadrature_formula,
2320 *   update_JxW_values | update_values |
2321 *   update_gradients);
2322 *  
2323 *   const unsigned int dofs_per_cell = stokes_fe.n_dofs_per_cell();
2324 *   const unsigned int n_q_points = quadrature_formula.size();
2325 *  
2326 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
2327 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
2328 *  
2329 *   std::vector<Tensor<2, dim>> grad_phi_u(dofs_per_cell);
2330 *   std::vector<double> phi_p(dofs_per_cell);
2331 *  
2332 *   const FEValuesExtractors::Vector velocities(0);
2333 *   const FEValuesExtractors::Scalar pressure(dim);
2334 *  
2335 *   for (const auto &cell : stokes_dof_handler.active_cell_iterators())
2336 *   {
2337 *   stokes_fe_values.reinit(cell);
2338 *   local_matrix = 0;
2339 *  
2340 * @endcode
2341 *
2342 * The creation of the local matrix is rather simple. There are only a
2343 * Laplace term (on the velocity) and a mass matrix weighted by
2344 * @f$\eta^{-1}@f$ to be generated, so the creation of the local matrix is
2345 * done in two lines. Once the local matrix is ready (loop over rows
2346 * and columns in the local matrix on each quadrature point), we get
2347 * the local DoF indices and write the local information into the
2348 * global matrix. We do this as in @ref step_27 "step-27", i.e., we directly apply the
2349 * constraints from hanging nodes locally. By doing so, we don't have
2350 * to do that afterwards, and we don't also write into entries of the
2351 * matrix that will actually be set to zero again later when
2352 * eliminating constraints.
2353 *
2354 * @code
2355 *   for (unsigned int q = 0; q < n_q_points; ++q)
2356 *   {
2357 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
2358 *   {
2359 *   grad_phi_u[k] = stokes_fe_values[velocities].gradient(k, q);
2360 *   phi_p[k] = stokes_fe_values[pressure].value(k, q);
2361 *   }
2362 *  
2363 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2364 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
2365 *   local_matrix(i, j) +=
2366 *   (EquationData::eta *
2367 *   scalar_product(grad_phi_u[i], grad_phi_u[j]) +
2368 *   (1. / EquationData::eta) * phi_p[i] * phi_p[j]) *
2369 *   stokes_fe_values.JxW(q);
2370 *   }
2371 *  
2372 *   cell->get_dof_indices(local_dof_indices);
2373 *   stokes_constraints.distribute_local_to_global(
2374 *   local_matrix, local_dof_indices, stokes_preconditioner_matrix);
2375 *   }
2376 *   }
2377 *  
2378 *  
2379 *  
2380 * @endcode
2381 *
2382 *
2383 * <a name="step_31-BoussinesqFlowProblembuild_stokes_preconditioner"></a>
2384 * <h4>BoussinesqFlowProblem::build_stokes_preconditioner</h4>
2385 *
2386
2387 *
2388 * This function generates the inner preconditioners that are going to be
2389 * used for the Schur complement block preconditioner. Since the
2390 * preconditioners need only to be regenerated when the matrices change,
2391 * this function does not have to do anything in case the matrices have not
2392 * changed (i.e., the flag <code>rebuild_stokes_preconditioner</code> has
2393 * the value <code>false</code>). Otherwise its first task is to call
2394 * <code>assemble_stokes_preconditioner</code> to generate the
2395 * preconditioner matrices.
2396 *
2397
2398 *
2399 * Next, we set up the preconditioner for the velocity-velocity matrix
2400 * @f$A@f$. As explained in the introduction, we are going to use an AMG
2401 * preconditioner based on a vector Laplace matrix @f$\hat{A}@f$ (which is
2402 * spectrally close to the Stokes matrix @f$A@f$). Usually, the
2403 * TrilinosWrappers::PreconditionAMG class can be seen as a good black-box
2404 * preconditioner which does not need any special knowledge. In this case,
2405 * however, we have to be careful: since we build an AMG for a vector
2406 * problem, we have to tell the preconditioner setup which dofs belong to
2407 * which vector component. We do this using the function
2408 * DoFTools::extract_constant_modes, a function that generates a set of
2409 * <code>dim</code> vectors, where each one has ones in the respective
2410 * component of the vector problem and zeros elsewhere. Hence, these are the
2411 * constant modes on each component, which explains the name of the
2412 * variable.
2413 *
2414 * @code
2415 *   template <int dim>
2416 *   void BoussinesqFlowProblem<dim>::build_stokes_preconditioner()
2417 *   {
2418 *   if (rebuild_stokes_preconditioner == false)
2419 *   return;
2420 *  
2421 *   std::cout << " Rebuilding Stokes preconditioner..." << std::flush;
2422 *  
2423 *   assemble_stokes_preconditioner();
2424 *  
2425 *   Amg_preconditioner = std::make_shared<TrilinosWrappers::PreconditionAMG>();
2426 *  
2427 *   const FEValuesExtractors::Vector velocity_components(0);
2428 *   const std::vector<std::vector<bool>> constant_modes =
2429 *   DoFTools::extract_constant_modes(
2430 *   stokes_dof_handler, stokes_fe.component_mask(velocity_components));
2431 *   TrilinosWrappers::PreconditionAMG::AdditionalData amg_data;
2432 *   amg_data.constant_modes = constant_modes;
2433 *  
2434 * @endcode
2435 *
2436 * Next, we set some more options of the AMG preconditioner. In
2437 * particular, we need to tell the AMG setup that we use quadratic basis
2438 * functions for the velocity matrix (this implies more nonzero elements
2439 * in the matrix, so that a more robust algorithm needs to be chosen
2440 * internally). Moreover, we want to be able to control how the coarsening
2441 * structure is build up. The way the Trilinos smoothed aggregation AMG
2442 * does this is to look which matrix entries are of similar size as the
2443 * diagonal entry in order to algebraically build a coarse-grid
2444 * structure. By setting the parameter <code>aggregation_threshold</code>
2445 * to 0.02, we specify that all entries that are more than two percent of
2446 * size of some diagonal pivots in that row should form one coarse grid
2447 * point. This parameter is rather ad hoc, and some fine-tuning of it can
2448 * influence the performance of the preconditioner. As a rule of thumb,
2449 * larger values of <code>aggregation_threshold</code> will decrease the
2450 * number of iterations, but increase the costs per iteration. A look at
2451 * the Trilinos documentation will provide more information on these
2452 * parameters. With this data set, we then initialize the preconditioner
2453 * with the matrix we want it to apply to.
2454 *
2455
2456 *
2457 * Finally, we also initialize the preconditioner for the inversion of the
2458 * pressure mass matrix. This matrix is symmetric and well-behaved, so we
2459 * can chose a simple preconditioner. We stick with an incomplete Cholesky
2460 * (IC) factorization preconditioner, which is designed for symmetric
2461 * matrices. We could have also chosen an SSOR preconditioner with
2462 * relaxation factor around 1.2, but IC is cheaper for our example. We
2463 * wrap the preconditioners into a <code>std::shared_ptr</code>
2464 * pointer, which makes it easier to recreate the preconditioner next time
2465 * around since we do not have to care about destroying the previously
2466 * used object.
2467 *
2468 * @code
2469 *   amg_data.elliptic = true;
2470 *   amg_data.higher_order_elements = true;
2471 *   amg_data.smoother_sweeps = 2;
2472 *   amg_data.aggregation_threshold = 0.02;
2473 *   Amg_preconditioner->initialize(stokes_preconditioner_matrix.block(0, 0),
2474 *   amg_data);
2475 *  
2476 *   Mp_preconditioner = std::make_shared<TrilinosWrappers::PreconditionIC>();
2477 *   Mp_preconditioner->initialize(stokes_preconditioner_matrix.block(1, 1));
2478 *  
2479 *   std::cout << std::endl;
2480 *  
2481 *   rebuild_stokes_preconditioner = false;
2482 *   }
2483 *  
2484 *  
2485 *  
2486 * @endcode
2487 *
2488 *
2489 * <a name="step_31-BoussinesqFlowProblemassemble_stokes_system"></a>
2490 * <h4>BoussinesqFlowProblem::assemble_stokes_system</h4>
2491 *
2492
2493 *
2494 * The time lag scheme we use for advancing the coupled Stokes-temperature
2495 * system forces us to split up the assembly (and the solution of linear
2496 * systems) into two step. The first one is to create the Stokes system
2497 * matrix and right hand side, and the second is to create matrix and right
2498 * hand sides for the temperature dofs, which depends on the result of the
2499 * linear system for the velocity.
2500 *
2501
2502 *
2503 * This function is called at the beginning of each time step. In the first
2504 * time step or if the mesh has changed, indicated by the
2505 * <code>rebuild_stokes_matrix</code>, we need to assemble the Stokes
2506 * matrix; on the other hand, if the mesh hasn't changed and the matrix is
2507 * already available, this is not necessary and all we need to do is
2508 * assemble the right hand side vector which changes in each time step.
2509 *
2510
2511 *
2512 * Regarding the technical details of implementation, not much has changed
2513 * from @ref step_22 "step-22". We reset matrix and vector, create a quadrature formula on
2514 * the cells, and then create the respective FEValues object. For the update
2515 * flags, we require basis function derivatives only in case of a full
2516 * assembly, since they are not needed for the right hand side; as always,
2517 * choosing the minimal set of flags depending on what is currently needed
2518 * makes the call to FEValues::reinit further down in the program more
2519 * efficient.
2520 *
2521
2522 *
2523 * There is one thing that needs to be commented &ndash; since we have a
2524 * separate finite element and DoFHandler for the temperature, we need to
2525 * generate a second FEValues object for the proper evaluation of the
2526 * temperature solution. This isn't too complicated to realize here: just
2527 * use the temperature structures and set an update flag for the basis
2528 * function values which we need for evaluation of the temperature
2529 * solution. The only important part to remember here is that the same
2530 * quadrature formula is used for both FEValues objects to ensure that we
2531 * get matching information when we loop over the quadrature points of the
2532 * two objects.
2533 *
2534
2535 *
2536 * The declarations proceed with some shortcuts for array sizes, the
2537 * creation of the local matrix and right hand side as well as the vector
2538 * for the indices of the local dofs compared to the global system.
2539 *
2540 * @code
2541 *   template <int dim>
2542 *   void BoussinesqFlowProblem<dim>::assemble_stokes_system()
2543 *   {
2544 *   std::cout << " Assembling..." << std::flush;
2545 *  
2546 *   if (rebuild_stokes_matrix == true)
2547 *   stokes_matrix = 0;
2548 *  
2549 *   stokes_rhs = 0;
2550 *  
2551 *   const QGauss<dim> quadrature_formula(stokes_degree + 2);
2552 *   FEValues<dim> stokes_fe_values(
2553 *   stokes_fe,
2554 *   quadrature_formula,
2555 *   update_values | update_quadrature_points | update_JxW_values |
2556 *   (rebuild_stokes_matrix == true ? update_gradients : UpdateFlags(0)));
2557 *  
2558 *   FEValues<dim> temperature_fe_values(temperature_fe,
2559 *   quadrature_formula,
2560 *   update_values);
2561 *  
2562 *   const unsigned int dofs_per_cell = stokes_fe.n_dofs_per_cell();
2563 *   const unsigned int n_q_points = quadrature_formula.size();
2564 *  
2565 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
2566 *   Vector<double> local_rhs(dofs_per_cell);
2567 *  
2568 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
2569 *  
2570 * @endcode
2571 *
2572 * Next we need a vector that will contain the values of the temperature
2573 * solution at the previous time level at the quadrature points to
2574 * assemble the source term in the right hand side of the momentum
2575 * equation. Let's call this vector <code>old_solution_values</code>.
2576 *
2577
2578 *
2579 * The set of vectors we create next hold the evaluations of the basis
2580 * functions as well as their gradients and symmetrized gradients that
2581 * will be used for creating the matrices. Putting these into their own
2582 * arrays rather than asking the FEValues object for this information each
2583 * time it is needed is an optimization to accelerate the assembly
2584 * process, see @ref step_22 "step-22" for details.
2585 *
2586
2587 *
2588 * The last two declarations are used to extract the individual blocks
2589 * (velocity, pressure, temperature) from the total FE system.
2590 *
2591 * @code
2592 *   std::vector<double> old_temperature_values(n_q_points);
2593 *  
2594 *   std::vector<Tensor<1, dim>> phi_u(dofs_per_cell);
2595 *   std::vector<SymmetricTensor<2, dim>> grads_phi_u(dofs_per_cell);
2596 *   std::vector<double> div_phi_u(dofs_per_cell);
2597 *   std::vector<double> phi_p(dofs_per_cell);
2598 *  
2599 *   const FEValuesExtractors::Vector velocities(0);
2600 *   const FEValuesExtractors::Scalar pressure(dim);
2601 *  
2602 * @endcode
2603 *
2604 * Now start the loop over all cells in the problem. We are working on two
2605 * different DoFHandlers for this assembly routine, so we must have two
2606 * different cell iterators for the two objects in use. This might seem a
2607 * bit peculiar, since both the Stokes system and the temperature system
2608 * use the same grid, but that's the only way to keep degrees of freedom
2609 * in sync. The first statements within the loop are again all very
2610 * familiar, doing the update of the finite element data as specified by
2611 * the update flags, zeroing out the local arrays and getting the values
2612 * of the old solution at the quadrature points. Then we are ready to loop
2613 * over the quadrature points on the cell.
2614 *
2615 * @code
2616 *   auto cell = stokes_dof_handler.begin_active();
2617 *   const auto endc = stokes_dof_handler.end();
2618 *   auto temperature_cell = temperature_dof_handler.begin_active();
2619 *  
2620 *   for (; cell != endc; ++cell, ++temperature_cell)
2621 *   {
2622 *   stokes_fe_values.reinit(cell);
2623 *   temperature_fe_values.reinit(temperature_cell);
2624 *  
2625 *   local_matrix = 0;
2626 *   local_rhs = 0;
2627 *  
2628 *   temperature_fe_values.get_function_values(old_temperature_solution,
2629 *   old_temperature_values);
2630 *  
2631 *   for (unsigned int q = 0; q < n_q_points; ++q)
2632 *   {
2633 *   const double old_temperature = old_temperature_values[q];
2634 *  
2635 * @endcode
2636 *
2637 * Next we extract the values and gradients of basis functions
2638 * relevant to the terms in the inner products. As shown in
2639 * @ref step_22 "step-22" this helps accelerate assembly.
2640 *
2641
2642 *
2643 * Once this is done, we start the loop over the rows and columns
2644 * of the local matrix and feed the matrix with the relevant
2645 * products. The right hand side is filled with the forcing term
2646 * driven by temperature in direction of gravity (which is
2647 * vertical in our example). Note that the right hand side term
2648 * is always generated, whereas the matrix contributions are only
2649 * updated when it is requested by the
2650 * <code>rebuild_matrices</code> flag.
2651 *
2652 * @code
2653 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
2654 *   {
2655 *   phi_u[k] = stokes_fe_values[velocities].value(k, q);
2656 *   if (rebuild_stokes_matrix)
2657 *   {
2658 *   grads_phi_u[k] =
2659 *   stokes_fe_values[velocities].symmetric_gradient(k, q);
2660 *   div_phi_u[k] =
2661 *   stokes_fe_values[velocities].divergence(k, q);
2662 *   phi_p[k] = stokes_fe_values[pressure].value(k, q);
2663 *   }
2664 *   }
2665 *  
2666 *   if (rebuild_stokes_matrix)
2667 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2668 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
2669 *   local_matrix(i, j) +=
2670 *   (EquationData::eta * 2 * (grads_phi_u[i] * grads_phi_u[j]) -
2671 *   div_phi_u[i] * phi_p[j] - phi_p[i] * div_phi_u[j]) *
2672 *   stokes_fe_values.JxW(q);
2673 *  
2674 *   const Point<dim> gravity =
2675 *   -((dim == 2) ? (Point<dim>(0, 1)) : (Point<dim>(0, 0, 1)));
2676 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2677 *   local_rhs(i) += (-EquationData::density * EquationData::beta *
2678 *   gravity * phi_u[i] * old_temperature) *
2679 *   stokes_fe_values.JxW(q);
2680 *   }
2681 *  
2682 * @endcode
2683 *
2684 * The last step in the loop over all cells is to enter the local
2685 * contributions into the global matrix and vector structures to the
2686 * positions specified in <code>local_dof_indices</code>. Again, we
2687 * let the AffineConstraints class do the insertion of the cell
2688 * matrix elements to the global matrix, which already condenses the
2689 * hanging node constraints.
2690 *
2691 * @code
2692 *   cell->get_dof_indices(local_dof_indices);
2693 *  
2694 *   if (rebuild_stokes_matrix == true)
2695 *   stokes_constraints.distribute_local_to_global(local_matrix,
2696 *   local_rhs,
2697 *   local_dof_indices,
2698 *   stokes_matrix,
2699 *   stokes_rhs);
2700 *   else
2701 *   stokes_constraints.distribute_local_to_global(local_rhs,
2702 *   local_dof_indices,
2703 *   stokes_rhs);
2704 *   }
2705 *  
2706 *   rebuild_stokes_matrix = false;
2707 *  
2708 *   std::cout << std::endl;
2709 *   }
2710 *  
2711 *  
2712 *  
2713 * @endcode
2714 *
2715 *
2716 * <a name="step_31-BoussinesqFlowProblemassemble_temperature_matrix"></a>
2717 * <h4>BoussinesqFlowProblem::assemble_temperature_matrix</h4>
2718 *
2719
2720 *
2721 * This function assembles the matrix in the temperature equation. The
2722 * temperature matrix consists of two parts, a mass matrix and the time step
2723 * size times a stiffness matrix given by a Laplace term times the amount of
2724 * diffusion. Since the matrix depends on the time step size (which varies
2725 * from one step to another), the temperature matrix needs to be updated
2726 * every time step. We could simply regenerate the matrices in every time
2727 * step, but this is not really efficient since mass and Laplace matrix do
2728 * only change when we change the mesh. Hence, we do this more efficiently
2729 * by generating two separate matrices in this function, one for the mass
2730 * matrix and one for the stiffness (diffusion) matrix. We will then sum up
2731 * the matrix plus the stiffness matrix times the time step size once we
2732 * know the actual time step.
2733 *
2734
2735 *
2736 * So the details for this first step are very simple. In case we need to
2737 * rebuild the matrix (i.e., the mesh has changed), we zero the data
2738 * structures, get a quadrature formula and a FEValues object, and create
2739 * local matrices, local dof indices and evaluation structures for the basis
2740 * functions.
2741 *
2742 * @code
2743 *   template <int dim>
2744 *   void BoussinesqFlowProblem<dim>::assemble_temperature_matrix()
2745 *   {
2746 *   if (rebuild_temperature_matrices == false)
2747 *   return;
2748 *  
2749 *   temperature_mass_matrix = 0;
2750 *   temperature_stiffness_matrix = 0;
2751 *  
2752 *   const QGauss<dim> quadrature_formula(temperature_degree + 2);
2753 *   FEValues<dim> temperature_fe_values(temperature_fe,
2754 *   quadrature_formula,
2755 *   update_values | update_gradients |
2756 *   update_JxW_values);
2757 *  
2758 *   const unsigned int dofs_per_cell = temperature_fe.n_dofs_per_cell();
2759 *   const unsigned int n_q_points = quadrature_formula.size();
2760 *  
2761 *   FullMatrix<double> local_mass_matrix(dofs_per_cell, dofs_per_cell);
2762 *   FullMatrix<double> local_stiffness_matrix(dofs_per_cell, dofs_per_cell);
2763 *  
2764 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
2765 *  
2766 *   std::vector<double> phi_T(dofs_per_cell);
2767 *   std::vector<Tensor<1, dim>> grad_phi_T(dofs_per_cell);
2768 *  
2769 * @endcode
2770 *
2771 * Now, let's start the loop over all cells in the triangulation. We need
2772 * to zero out the local matrices, update the finite element evaluations,
2773 * and then loop over the rows and columns of the matrices on each
2774 * quadrature point, where we then create the mass matrix and the
2775 * stiffness matrix (Laplace terms times the diffusion
2776 * <code>EquationData::kappa</code>. Finally, we let the constraints
2777 * object insert these values into the global matrix, and directly
2778 * condense the constraints into the matrix.
2779 *
2780 * @code
2781 *   for (const auto &cell : temperature_dof_handler.active_cell_iterators())
2782 *   {
2783 *   local_mass_matrix = 0;
2784 *   local_stiffness_matrix = 0;
2785 *  
2786 *   temperature_fe_values.reinit(cell);
2787 *  
2788 *   for (unsigned int q = 0; q < n_q_points; ++q)
2789 *   {
2790 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
2791 *   {
2792 *   grad_phi_T[k] = temperature_fe_values.shape_grad(k, q);
2793 *   phi_T[k] = temperature_fe_values.shape_value(k, q);
2794 *   }
2795 *  
2796 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2797 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
2798 *   {
2799 *   local_mass_matrix(i, j) +=
2800 *   (phi_T[i] * phi_T[j] * temperature_fe_values.JxW(q));
2801 *   local_stiffness_matrix(i, j) +=
2802 *   (EquationData::kappa * grad_phi_T[i] * grad_phi_T[j] *
2803 *   temperature_fe_values.JxW(q));
2804 *   }
2805 *   }
2806 *  
2807 *   cell->get_dof_indices(local_dof_indices);
2808 *  
2809 *   temperature_constraints.distribute_local_to_global(
2810 *   local_mass_matrix, local_dof_indices, temperature_mass_matrix);
2811 *   temperature_constraints.distribute_local_to_global(
2812 *   local_stiffness_matrix,
2813 *   local_dof_indices,
2814 *   temperature_stiffness_matrix);
2815 *   }
2816 *  
2817 *   rebuild_temperature_matrices = false;
2818 *   }
2819 *  
2820 *  
2821 *  
2822 * @endcode
2823 *
2824 *
2825 * <a name="step_31-BoussinesqFlowProblemassemble_temperature_system"></a>
2826 * <h4>BoussinesqFlowProblem::assemble_temperature_system</h4>
2827 *
2828
2829 *
2830 * This function does the second part of the assembly work on the
2831 * temperature matrix, the actual addition of pressure mass and stiffness
2832 * matrix (where the time step size comes into play), as well as the
2833 * creation of the velocity-dependent right hand side. The declarations for
2834 * the right hand side assembly in this function are pretty much the same as
2835 * the ones used in the other assembly routines, except that we restrict
2836 * ourselves to vectors this time. We are going to calculate residuals on
2837 * the temperature system, which means that we have to evaluate second
2838 * derivatives, specified by the update flag <code>update_hessians</code>.
2839 *
2840
2841 *
2842 * The temperature equation is coupled to the Stokes system by means of the
2843 * fluid velocity. These two parts of the solution are associated with
2844 * different DoFHandlers, so we again need to create a second FEValues
2845 * object for the evaluation of the velocity at the quadrature points.
2846 *
2847 * @code
2848 *   template <int dim>
2849 *   void BoussinesqFlowProblem<dim>::assemble_temperature_system(
2850 *   const double maximal_velocity)
2851 *   {
2852 *   const bool use_bdf2_scheme = (timestep_number != 0);
2853 *  
2854 *   if (use_bdf2_scheme == true)
2855 *   {
2856 *   temperature_matrix.copy_from(temperature_mass_matrix);
2857 *   temperature_matrix *=
2858 *   (2 * time_step + old_time_step) / (time_step + old_time_step);
2859 *   temperature_matrix.add(time_step, temperature_stiffness_matrix);
2860 *   }
2861 *   else
2862 *   {
2863 *   temperature_matrix.copy_from(temperature_mass_matrix);
2864 *   temperature_matrix.add(time_step, temperature_stiffness_matrix);
2865 *   }
2866 *  
2867 *   temperature_rhs = 0;
2868 *  
2869 *   const QGauss<dim> quadrature_formula(temperature_degree + 2);
2870 *   FEValues<dim> temperature_fe_values(temperature_fe,
2871 *   quadrature_formula,
2873 *   update_hessians |
2875 *   update_JxW_values);
2876 *   FEValues<dim> stokes_fe_values(stokes_fe,
2877 *   quadrature_formula,
2878 *   update_values);
2879 *  
2880 *   const unsigned int dofs_per_cell = temperature_fe.n_dofs_per_cell();
2881 *   const unsigned int n_q_points = quadrature_formula.size();
2882 *  
2883 *   Vector<double> local_rhs(dofs_per_cell);
2884 *  
2885 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
2886 *  
2887 * @endcode
2888 *
2889 * Next comes the declaration of vectors to hold the old and older
2890 * solution values (as a notation for time levels @f$n-1@f$ and
2891 * @f$n-2@f$, respectively) and gradients at quadrature points of the
2892 * current cell. We also declare an object to hold the temperature right
2893 * hand side values (<code>gamma_values</code>), and we again use
2894 * shortcuts for the temperature basis functions. Eventually, we need to
2895 * find the temperature extrema and the diameter of the computational
2896 * domain which will be used for the definition of the stabilization
2897 * parameter (we got the maximal velocity as an input to this function).
2898 *
2899 * @code
2900 *   std::vector<Tensor<1, dim>> old_velocity_values(n_q_points);
2901 *   std::vector<Tensor<1, dim>> old_old_velocity_values(n_q_points);
2902 *   std::vector<double> old_temperature_values(n_q_points);
2903 *   std::vector<double> old_old_temperature_values(n_q_points);
2904 *   std::vector<Tensor<1, dim>> old_temperature_grads(n_q_points);
2905 *   std::vector<Tensor<1, dim>> old_old_temperature_grads(n_q_points);
2906 *   std::vector<double> old_temperature_laplacians(n_q_points);
2907 *   std::vector<double> old_old_temperature_laplacians(n_q_points);
2908 *  
2909 *   EquationData::TemperatureRightHandSide<dim> temperature_right_hand_side;
2910 *   std::vector<double> gamma_values(n_q_points);
2911 *  
2912 *   std::vector<double> phi_T(dofs_per_cell);
2913 *   std::vector<Tensor<1, dim>> grad_phi_T(dofs_per_cell);
2914 *  
2915 *   const std::pair<double, double> global_T_range =
2916 *   get_extrapolated_temperature_range();
2917 *  
2918 *   const FEValuesExtractors::Vector velocities(0);
2919 *  
2920 * @endcode
2921 *
2922 * Now, let's start the loop over all cells in the triangulation. Again,
2923 * we need two cell iterators that walk in parallel through the cells of
2924 * the two involved DoFHandler objects for the Stokes and temperature
2925 * part. Within the loop, we first set the local rhs to zero, and then get
2926 * the values and derivatives of the old solution functions at the
2927 * quadrature points, since they are going to be needed for the definition
2928 * of the stabilization parameters and as coefficients in the equation,
2929 * respectively. Note that since the temperature has its own DoFHandler
2930 * and FEValues object we get the entire solution at the quadrature point
2931 * (which is the scalar temperature field only anyway) whereas for the
2932 * Stokes part we restrict ourselves to extracting the velocity part (and
2933 * ignoring the pressure part) by using
2934 * <code>stokes_fe_values[velocities].get_function_values</code>.
2935 *
2936 * @code
2937 *   auto cell = temperature_dof_handler.begin_active();
2938 *   const auto endc = temperature_dof_handler.end();
2939 *   auto stokes_cell = stokes_dof_handler.begin_active();
2940 *  
2941 *   for (; cell != endc; ++cell, ++stokes_cell)
2942 *   {
2943 *   local_rhs = 0;
2944 *  
2945 *   temperature_fe_values.reinit(cell);
2946 *   stokes_fe_values.reinit(stokes_cell);
2947 *  
2948 *   temperature_fe_values.get_function_values(old_temperature_solution,
2949 *   old_temperature_values);
2950 *   temperature_fe_values.get_function_values(old_old_temperature_solution,
2951 *   old_old_temperature_values);
2952 *  
2953 *   temperature_fe_values.get_function_gradients(old_temperature_solution,
2954 *   old_temperature_grads);
2955 *   temperature_fe_values.get_function_gradients(
2956 *   old_old_temperature_solution, old_old_temperature_grads);
2957 *  
2958 *   temperature_fe_values.get_function_laplacians(
2959 *   old_temperature_solution, old_temperature_laplacians);
2960 *   temperature_fe_values.get_function_laplacians(
2961 *   old_old_temperature_solution, old_old_temperature_laplacians);
2962 *  
2963 *   temperature_right_hand_side.value_list(
2964 *   temperature_fe_values.get_quadrature_points(), gamma_values);
2965 *  
2966 *   stokes_fe_values[velocities].get_function_values(stokes_solution,
2967 *   old_velocity_values);
2968 *   stokes_fe_values[velocities].get_function_values(
2969 *   old_stokes_solution, old_old_velocity_values);
2970 *  
2971 * @endcode
2972 *
2973 * Next, we calculate the artificial viscosity for stabilization
2974 * according to the discussion in the introduction using the dedicated
2975 * function. With that at hand, we can get into the loop over
2976 * quadrature points and local rhs vector components. The terms here
2977 * are quite lengthy, but their definition follows the time-discrete
2978 * system developed in the introduction of this program. The BDF-2
2979 * scheme needs one more term from the old time step (and involves
2980 * more complicated factors) than the backward Euler scheme that is
2981 * used for the first time step. When all this is done, we distribute
2982 * the local vector into the global one (including hanging node
2983 * constraints).
2984 *
2985 * @code
2986 *   const double nu =
2987 *   compute_viscosity(old_temperature_values,
2988 *   old_old_temperature_values,
2989 *   old_temperature_grads,
2990 *   old_old_temperature_grads,
2991 *   old_temperature_laplacians,
2992 *   old_old_temperature_laplacians,
2993 *   old_velocity_values,
2994 *   old_old_velocity_values,
2995 *   gamma_values,
2996 *   maximal_velocity,
2997 *   global_T_range.second - global_T_range.first,
2998 *   cell->diameter());
2999 *  
3000 *   for (unsigned int q = 0; q < n_q_points; ++q)
3001 *   {
3002 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
3003 *   {
3004 *   grad_phi_T[k] = temperature_fe_values.shape_grad(k, q);
3005 *   phi_T[k] = temperature_fe_values.shape_value(k, q);
3006 *   }
3007 *  
3008 *   const double T_term_for_rhs =
3009 *   (use_bdf2_scheme ?
3010 *   (old_temperature_values[q] * (1 + time_step / old_time_step) -
3011 *   old_old_temperature_values[q] * (time_step * time_step) /
3012 *   (old_time_step * (time_step + old_time_step))) :
3013 *   old_temperature_values[q]);
3014 *  
3015 *   const Tensor<1, dim> ext_grad_T =
3016 *   (use_bdf2_scheme ?
3017 *   (old_temperature_grads[q] * (1 + time_step / old_time_step) -
3018 *   old_old_temperature_grads[q] * time_step / old_time_step) :
3019 *   old_temperature_grads[q]);
3020 *  
3021 *   const Tensor<1, dim> extrapolated_u =
3022 *   (use_bdf2_scheme ?
3023 *   (old_velocity_values[q] * (1 + time_step / old_time_step) -
3024 *   old_old_velocity_values[q] * time_step / old_time_step) :
3025 *   old_velocity_values[q]);
3026 *  
3027 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3028 *   local_rhs(i) +=
3029 *   (T_term_for_rhs * phi_T[i] -
3030 *   time_step * extrapolated_u * ext_grad_T * phi_T[i] -
3031 *   time_step * nu * ext_grad_T * grad_phi_T[i] +
3032 *   time_step * gamma_values[q] * phi_T[i]) *
3033 *   temperature_fe_values.JxW(q);
3034 *   }
3035 *  
3036 *   cell->get_dof_indices(local_dof_indices);
3037 *   temperature_constraints.distribute_local_to_global(local_rhs,
3038 *   local_dof_indices,
3039 *   temperature_rhs);
3040 *   }
3041 *   }
3042 *  
3043 *  
3044 *  
3045 * @endcode
3046 *
3047 *
3048 * <a name="step_31-BoussinesqFlowProblemsolve"></a>
3049 * <h4>BoussinesqFlowProblem::solve</h4>
3050 *
3051
3052 *
3053 * This function solves the linear systems of equations. Following the
3054 * introduction, we start with the Stokes system, where we need to generate
3055 * our block Schur preconditioner. Since all the relevant actions are
3056 * implemented in the class <code>BlockSchurPreconditioner</code>, all we
3057 * have to do is to initialize the class appropriately. What we need to pass
3058 * down is an <code>InverseMatrix</code> object for the pressure mass
3059 * matrix, which we set up using the respective class together with the IC
3060 * preconditioner we already generated, and the AMG preconditioner for the
3061 * velocity-velocity matrix. Note that both <code>Mp_preconditioner</code>
3062 * and <code>Amg_preconditioner</code> are only pointers, so we use
3063 * <code>*</code> to pass down the actual preconditioner objects.
3064 *
3065
3066 *
3067 * Once the preconditioner is ready, we create a GMRES solver for the block
3068 * system. Since we are working with Trilinos data structures, we have to
3069 * set the respective template argument in the solver. GMRES needs to
3070 * internally store temporary vectors for each iteration (see the discussion
3071 * in the results section of @ref step_22 "step-22") &ndash; the more vectors it can use,
3072 * the better it will generally perform. To keep memory demands in check, we
3073 * set the number of vectors to 100. This means that up to 100 solver
3074 * iterations, every temporary vector can be stored. If the solver needs to
3075 * iterate more often to get the specified tolerance, it will work on a
3076 * reduced set of vectors by restarting at every 100 iterations.
3077 *
3078
3079 *
3080 * With this all set up, we solve the system and distribute the constraints
3081 * in the Stokes system, i.e., hanging nodes and no-flux boundary condition,
3082 * in order to have the appropriate solution values even at constrained
3083 * dofs. Finally, we write the number of iterations to the screen.
3084 *
3085 * @code
3086 *   template <int dim>
3087 *   void BoussinesqFlowProblem<dim>::solve()
3088 *   {
3089 *   std::cout << " Solving..." << std::endl;
3090 *  
3091 *   {
3092 *   const LinearSolvers::InverseMatrix<TrilinosWrappers::SparseMatrix,
3093 *   TrilinosWrappers::PreconditionIC>
3094 *   mp_inverse(stokes_preconditioner_matrix.block(1, 1),
3095 *   *Mp_preconditioner);
3096 *  
3097 *   const LinearSolvers::BlockSchurPreconditioner<
3098 *   TrilinosWrappers::PreconditionAMG,
3099 *   TrilinosWrappers::PreconditionIC>
3100 *   preconditioner(stokes_matrix, mp_inverse, *Amg_preconditioner);
3101 *  
3102 *   SolverControl solver_control(stokes_matrix.m(),
3103 *   1e-6 * stokes_rhs.l2_norm());
3104 *  
3105 *   SolverGMRES<TrilinosWrappers::MPI::BlockVector> gmres(
3106 *   solver_control,
3107 *   SolverGMRES<TrilinosWrappers::MPI::BlockVector>::AdditionalData(100));
3108 *  
3109 *   for (unsigned int i = 0; i < stokes_solution.size(); ++i)
3110 *   if (stokes_constraints.is_constrained(i))
3111 *   stokes_solution(i) = 0;
3112 *  
3113 *   gmres.solve(stokes_matrix, stokes_solution, stokes_rhs, preconditioner);
3114 *  
3115 *   stokes_constraints.distribute(stokes_solution);
3116 *  
3117 *   std::cout << " " << solver_control.last_step()
3118 *   << " GMRES iterations for Stokes subsystem." << std::endl;
3119 *   }
3120 *  
3121 * @endcode
3122 *
3123 * Once we know the Stokes solution, we can determine the new time step
3124 * from the maximal velocity. We have to do this to satisfy the CFL
3125 * condition since convection terms are treated explicitly in the
3126 * temperature equation, as discussed in the introduction. The exact form
3127 * of the formula used here for the time step is discussed in the results
3128 * section of this program.
3129 *
3130
3131 *
3132 * There is a snatch here. The formula contains a division by the maximum
3133 * value of the velocity. However, at the start of the computation, we
3134 * have a constant temperature field (we start with a constant
3135 * temperature, and it will be nonconstant only after the first time step
3136 * during which the source acts). Constant temperature means that no
3137 * buoyancy acts, and so the velocity is zero. Dividing by it will not
3138 * likely lead to anything good.
3139 *
3140
3141 *
3142 * To avoid the resulting infinite time step, we ask whether the maximal
3143 * velocity is very small (in particular smaller than the values we
3144 * encounter during any of the following time steps) and if so rather than
3145 * dividing by zero we just divide by a small value, resulting in a large
3146 * but finite time step.
3147 *
3148 * @code
3149 *   old_time_step = time_step;
3150 *   const double maximal_velocity = get_maximal_velocity();
3151 *  
3152 *   if (maximal_velocity >= 0.01)
3153 *   time_step = 1. / (1.7 * dim * std::sqrt(1. * dim)) / temperature_degree *
3154 *   GridTools::minimal_cell_diameter(triangulation) /
3155 *   maximal_velocity;
3156 *   else
3157 *   time_step = 1. / (1.7 * dim * std::sqrt(1. * dim)) / temperature_degree *
3158 *   GridTools::minimal_cell_diameter(triangulation) / .01;
3159 *  
3160 *   std::cout << " "
3161 *   << "Time step: " << time_step << std::endl;
3162 *  
3163 *   temperature_solution = old_temperature_solution;
3164 *  
3165 * @endcode
3166 *
3167 * Next we set up the temperature system and the right hand side using the
3168 * function <code>assemble_temperature_system()</code>. Knowing the
3169 * matrix and right hand side of the temperature equation, we set up a
3170 * preconditioner and a solver. The temperature matrix is a mass matrix
3171 * (with eigenvalues around one) plus a Laplace matrix (with eigenvalues
3172 * between zero and @f$ch^{-2}@f$) times a small number proportional to the
3173 * time step @f$k_n@f$. Hence, the resulting symmetric and positive definite
3174 * matrix has eigenvalues in the range @f$[1,1+k_nh^{-2}]@f$ (up to
3175 * constants). This matrix is only moderately ill conditioned even for
3176 * small mesh sizes and we get a reasonably good preconditioner by simple
3177 * means, for example with an incomplete Cholesky decomposition
3178 * preconditioner (IC) as we also use for preconditioning the pressure
3179 * mass matrix solver. As a solver, we choose the conjugate gradient
3180 * method CG. As before, we tell the solver to use Trilinos vectors via
3181 * the template argument <code>TrilinosWrappers::MPI::Vector</code>.
3182 * Finally, we solve, distribute the hanging node constraints and write out
3183 * the number of iterations.
3184 *
3185 * @code
3186 *   assemble_temperature_system(maximal_velocity);
3187 *   {
3188 *   SolverControl solver_control(temperature_matrix.m(),
3189 *   1e-8 * temperature_rhs.l2_norm());
3190 *   SolverCG<TrilinosWrappers::MPI::Vector> cg(solver_control);
3191 *  
3192 *   TrilinosWrappers::PreconditionIC preconditioner;
3193 *   preconditioner.initialize(temperature_matrix);
3194 *  
3195 *   cg.solve(temperature_matrix,
3196 *   temperature_solution,
3197 *   temperature_rhs,
3198 *   preconditioner);
3199 *  
3200 *   temperature_constraints.distribute(temperature_solution);
3201 *  
3202 *   std::cout << " " << solver_control.last_step()
3203 *   << " CG iterations for temperature." << std::endl;
3204 *  
3205 * @endcode
3206 *
3207 * At the end of this function, we step through the vector and read out
3208 * the maximum and minimum temperature value, which we also want to
3209 * output. This will come in handy when determining the correct constant
3210 * in the choice of time step as discuss in the results section of this
3211 * program.
3212 *
3213 * @code
3214 *   double min_temperature = temperature_solution(0),
3215 *   max_temperature = temperature_solution(0);
3216 *   for (unsigned int i = 0; i < temperature_solution.size(); ++i)
3217 *   {
3218 *   min_temperature =
3219 *   std::min<double>(min_temperature, temperature_solution(i));
3220 *   max_temperature =
3221 *   std::max<double>(max_temperature, temperature_solution(i));
3222 *   }
3223 *  
3224 *   std::cout << " Temperature range: " << min_temperature << ' '
3225 *   << max_temperature << std::endl;
3226 *   }
3227 *   }
3228 *  
3229 *  
3230 *  
3231 * @endcode
3232 *
3233 *
3234 * <a name="step_31-BoussinesqFlowProblemoutput_results"></a>
3235 * <h4>BoussinesqFlowProblem::output_results</h4>
3236 *
3237
3238 *
3239 * This function writes the solution to a VTK output file for visualization,
3240 * which is done every tenth time step. This is usually quite a simple task,
3241 * since the deal.II library provides functions that do almost all the job
3242 * for us. There is one new function compared to previous examples: We want
3243 * to visualize both the Stokes solution and the temperature as one data
3244 * set, but we have done all the calculations based on two different
3245 * DoFHandler objects. Luckily, the DataOut class is prepared to deal with
3246 * it. All we have to do is to not attach one single DoFHandler at the
3247 * beginning and then use that for all added vector, but specify the
3248 * DoFHandler to each vector separately. The rest is done as in @ref step_22 "step-22". We
3249 * create solution names (that are going to appear in the visualization
3250 * program for the individual components). The first <code>dim</code>
3251 * components are the vector velocity, and then we have pressure for the
3252 * Stokes part, whereas temperature is scalar. This information is read out
3253 * using the DataComponentInterpretation helper class. Next, we actually
3254 * attach the data vectors with their DoFHandler objects, build patches
3255 * according to the degree of freedom, which are (sub-) elements that
3256 * describe the data for visualization programs. Finally, we open a file
3257 * (that includes the time step number) and write the vtk data into it.
3258 *
3259 * @code
3260 *   template <int dim>
3261 *   void BoussinesqFlowProblem<dim>::output_results() const
3262 *   {
3263 *   if (timestep_number % 10 != 0)
3264 *   return;
3265 *  
3266 *   std::vector<std::string> stokes_names(dim, "velocity");
3267 *   stokes_names.emplace_back("p");
3268 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
3269 *   stokes_component_interpretation(
3270 *   dim + 1, DataComponentInterpretation::component_is_scalar);
3271 *   for (unsigned int i = 0; i < dim; ++i)
3272 *   stokes_component_interpretation[i] =
3273 *   DataComponentInterpretation::component_is_part_of_vector;
3274 *  
3275 *   DataOut<dim> data_out;
3276 *   data_out.add_data_vector(stokes_dof_handler,
3277 *   stokes_solution,
3278 *   stokes_names,
3279 *   stokes_component_interpretation);
3280 *   data_out.add_data_vector(temperature_dof_handler,
3281 *   temperature_solution,
3282 *   "T");
3283 *   data_out.build_patches(std::min(stokes_degree, temperature_degree));
3284 *  
3285 *   std::ofstream output("solution-" +
3286 *   Utilities::int_to_string(timestep_number, 4) + ".vtk");
3287 *   data_out.write_vtk(output);
3288 *   }
3289 *  
3290 *  
3291 *  
3292 * @endcode
3293 *
3294 *
3295 * <a name="step_31-BoussinesqFlowProblemrefine_mesh"></a>
3296 * <h4>BoussinesqFlowProblem::refine_mesh</h4>
3297 *
3298
3299 *
3300 * This function takes care of the adaptive mesh refinement. The three tasks
3301 * this function performs is to first find out which cells to
3302 * refine/coarsen, then to actually do the refinement and eventually
3303 * transfer the solution vectors between the two different grids. The first
3304 * task is simply achieved by using the well-established Kelly error
3305 * estimator on the temperature (it is the temperature we're mainly
3306 * interested in for this program, and we need to be accurate in regions of
3307 * high temperature gradients, also to not have too much numerical
3308 * diffusion). The second task is to actually do the remeshing. That
3309 * involves only basic functions as well, such as the
3310 * <code>refine_and_coarsen_fixed_fraction</code> that refines those cells
3311 * with the largest estimated error that together make up 80 per cent of the
3312 * error, and coarsens those cells with the smallest error that make up for
3313 * a combined 10 per cent of the error.
3314 *
3315
3316 *
3317 * If implemented like this, we would get a program that will not make much
3318 * progress: Remember that we expect temperature fields that are nearly
3319 * discontinuous (the diffusivity @f$\kappa@f$ is very small after all) and
3320 * consequently we can expect that a freely adapted mesh will refine further
3321 * and further into the areas of large gradients. This decrease in mesh size
3322 * will then be accompanied by a decrease in time step, requiring an
3323 * exceedingly large number of time steps to solve to a given final time. It
3324 * will also lead to meshes that are much better at resolving
3325 * discontinuities after several mesh refinement cycles than in the
3326 * beginning.
3327 *
3328
3329 *
3330 * In particular to prevent the decrease in time step size and the
3331 * correspondingly large number of time steps, we limit the maximal
3332 * refinement depth of the mesh. To this end, after the refinement indicator
3333 * has been applied to the cells, we simply loop over all cells on the
3334 * finest level and unselect them from refinement if they would result in
3335 * too high a mesh level.
3336 *
3337 * @code
3338 *   template <int dim>
3339 *   void
3340 *   BoussinesqFlowProblem<dim>::refine_mesh(const unsigned int max_grid_level)
3341 *   {
3342 *   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
3343 *  
3344 *   KellyErrorEstimator<dim>::estimate(temperature_dof_handler,
3345 *   QGauss<dim - 1>(temperature_degree + 1),
3346 *   {},
3347 *   temperature_solution,
3348 *   estimated_error_per_cell);
3349 *  
3351 *   estimated_error_per_cell,
3352 *   0.8,
3353 *   0.1);
3354 *   if (triangulation.n_levels() > max_grid_level)
3355 *   for (auto &cell :
3356 *   triangulation.active_cell_iterators_on_level(max_grid_level))
3357 *   cell->clear_refine_flag();
3358 *  
3359 * @endcode
3360 *
3361 * As part of mesh refinement we need to transfer the solution vectors
3362 * from the old mesh to the new one. To this end we use the
3363 * SolutionTransfer class and we have to prepare the solution vectors that
3364 * should be transferred to the new grid (we will lose the old grid once
3365 * we have done the refinement so the transfer has to happen concurrently
3366 * with refinement). What we definitely need are the current and the old
3367 * temperature (BDF-2 time stepping requires two old solutions). Since the
3368 * SolutionTransfer objects only support to transfer one object per dof
3369 * handler, we need to collect the two temperature solutions in one data
3370 * structure. Moreover, we choose to transfer the Stokes solution, too,
3371 * since we need the velocity at two previous time steps, of which only
3372 * one is calculated on the fly.
3373 *
3374
3375 *
3376 * Consequently, we initialize two SolutionTransfer objects for the Stokes
3377 * and temperature DoFHandler objects, by attaching them to the old dof
3378 * handlers. With this at place, we can prepare the triangulation and the
3379 * data vectors for refinement (in this order).
3380 *
3381 * @code
3382 *   const std::vector<TrilinosWrappers::MPI::Vector> x_temperature = {
3383 *   temperature_solution, old_temperature_solution};
3384 *   TrilinosWrappers::MPI::BlockVector x_stokes = stokes_solution;
3385 *  
3387 *   temperature_dof_handler);
3389 *   stokes_dof_handler);
3390 *  
3391 *   triangulation.prepare_coarsening_and_refinement();
3392 *   temperature_trans.prepare_for_coarsening_and_refinement(x_temperature);
3393 *   stokes_trans.prepare_for_coarsening_and_refinement(x_stokes);
3394 *  
3395 * @endcode
3396 *
3397 * Now everything is ready, so do the refinement and recreate the dof
3398 * structure on the new grid, and initialize the matrix structures and the
3399 * new vectors in the <code>setup_dofs</code> function. Next, we actually
3400 * perform the interpolation of the solutions between the grids. We create
3401 * another copy of temporary vectors for temperature (now corresponding to
3402 * the new grid), and let the interpolate function do the job. Then, the
3403 * resulting array of vectors is written into the respective vector member
3404 * variables.
3405 *
3406
3407 *
3408 * Remember that the set of constraints will be updated for the new
3409 * triangulation in the setup_dofs() call.
3410 *
3411 * @code
3412 *   triangulation.execute_coarsening_and_refinement();
3413 *   setup_dofs();
3414 *  
3415 *   std::vector<TrilinosWrappers::MPI::Vector> tmp = {
3416 *   TrilinosWrappers::MPI::Vector(temperature_solution),
3417 *   TrilinosWrappers::MPI::Vector(temperature_solution)};
3418 *   temperature_trans.interpolate(tmp);
3419 *  
3420 *   temperature_solution = tmp[0];
3421 *   old_temperature_solution = tmp[1];
3422 *  
3423 * @endcode
3424 *
3425 * After the solution has been transferred we then enforce the constraints
3426 * on the transferred solution.
3427 *
3428 * @code
3429 *   temperature_constraints.distribute(temperature_solution);
3430 *   temperature_constraints.distribute(old_temperature_solution);
3431 *  
3432 * @endcode
3433 *
3434 * For the Stokes vector, everything is just the same &ndash; except that
3435 * we do not need another temporary vector since we just interpolate a
3436 * single vector. In the end, we have to tell the program that the matrices
3437 * and preconditioners need to be regenerated, since the mesh has changed.
3438 *
3439 * @code
3440 *   stokes_trans.interpolate(stokes_solution);
3441 *  
3442 *   stokes_constraints.distribute(stokes_solution);
3443 *  
3444 *   rebuild_stokes_matrix = true;
3445 *   rebuild_temperature_matrices = true;
3446 *   rebuild_stokes_preconditioner = true;
3447 *   }
3448 *  
3449 *  
3450 *  
3451 * @endcode
3452 *
3453 *
3454 * <a name="step_31-BoussinesqFlowProblemrun"></a>
3455 * <h4>BoussinesqFlowProblem::run</h4>
3456 *
3457
3458 *
3459 * This function performs all the essential steps in the Boussinesq
3460 * program. It starts by setting up a grid (depending on the spatial
3461 * dimension, we choose some different level of initial refinement and
3462 * additional adaptive refinement steps, and then create a cube in
3463 * <code>dim</code> dimensions and set up the dofs for the first time. Since
3464 * we want to start the time stepping already with an adaptively refined
3465 * grid, we perform some pre-refinement steps, consisting of all assembly,
3466 * solution and refinement, but without actually advancing in time. Rather,
3467 * we use the vilified <code>goto</code> statement to jump out of the time
3468 * loop right after mesh refinement to start all over again on the new mesh
3469 * beginning at the <code>start_time_iteration</code> label. (The use of the
3470 * <code>goto</code> is discussed in @ref step_26 "step-26".)
3471 *
3472
3473 *
3474 * Before we start, we project the initial values to the grid and obtain the
3475 * first data for the <code>old_temperature_solution</code> vector. Then, we
3476 * initialize time step number and time step and start the time loop.
3477 *
3478 * @code
3479 *   template <int dim>
3480 *   void BoussinesqFlowProblem<dim>::run()
3481 *   {
3482 *   const unsigned int initial_refinement = (dim == 2 ? 4 : 2);
3483 *   const unsigned int n_pre_refinement_steps = (dim == 2 ? 4 : 3);
3484 *  
3485 *  
3486 *   GridGenerator::hyper_cube(triangulation);
3487 *   global_Omega_diameter = GridTools::diameter(triangulation);
3488 *  
3489 *   triangulation.refine_global(initial_refinement);
3490 *  
3491 *   setup_dofs();
3492 *  
3493 *   unsigned int pre_refinement_step = 0;
3494 *  
3495 *   start_time_iteration:
3496 *  
3497 *   VectorTools::project(temperature_dof_handler,
3498 *   temperature_constraints,
3499 *   QGauss<dim>(temperature_degree + 2),
3500 *   EquationData::TemperatureInitialValues<dim>(),
3501 *   old_temperature_solution);
3502 *  
3503 *   timestep_number = 0;
3504 *   time_step = old_time_step = 0;
3505 *  
3506 *   double time = 0;
3507 *  
3508 *   do
3509 *   {
3510 *   std::cout << "Timestep " << timestep_number << ": t=" << time
3511 *   << std::endl;
3512 *  
3513 * @endcode
3514 *
3515 * The first steps in the time loop are all obvious &ndash; we
3516 * assemble the Stokes system, the preconditioner, the temperature
3517 * matrix (matrices and preconditioner do actually only change in case
3518 * we've remeshed before), and then do the solve. Before going on with
3519 * the next time step, we have to check whether we should first finish
3520 * the pre-refinement steps or if we should remesh (every fifth time
3521 * step), refining up to a level that is consistent with initial
3522 * refinement and pre-refinement steps. Last in the loop is to advance
3523 * the solutions, i.e., to copy the solutions to the next "older" time
3524 * level.
3525 *
3526 * @code
3527 *   assemble_stokes_system();
3528 *   build_stokes_preconditioner();
3529 *   assemble_temperature_matrix();
3530 *  
3531 *   solve();
3532 *  
3533 *   output_results();
3534 *  
3535 *   std::cout << std::endl;
3536 *  
3537 *   if ((timestep_number == 0) &&
3538 *   (pre_refinement_step < n_pre_refinement_steps))
3539 *   {
3540 *   refine_mesh(initial_refinement + n_pre_refinement_steps);
3541 *   ++pre_refinement_step;
3542 *   goto start_time_iteration;
3543 *   }
3544 *   else if ((timestep_number > 0) && (timestep_number % 5 == 0))
3545 *   refine_mesh(initial_refinement + n_pre_refinement_steps);
3546 *  
3547 *   time += time_step;
3548 *   ++timestep_number;
3549 *  
3550 *   old_stokes_solution = stokes_solution;
3551 *   old_old_temperature_solution = old_temperature_solution;
3552 *   old_temperature_solution = temperature_solution;
3553 *   }
3554 * @endcode
3555 *
3556 * Do all the above until we arrive at time 100.
3557 *
3558 * @code
3559 *   while (time <= 100);
3560 *   }
3561 *   } // namespace Step31
3562 *  
3563 *  
3564 *  
3565 * @endcode
3566 *
3567 *
3568 * <a name="step_31-Thecodemaincodefunction"></a>
3569 * <h3>The <code>main</code> function</h3>
3570 *
3571
3572 *
3573 * The main function looks almost the same as in all other programs.
3574 *
3575
3576 *
3577 * There is one difference we have to be careful about. This program uses
3578 * Trilinos and, typically, Trilinos is configured so that it can run in
3579 * %parallel using MPI. This doesn't mean that it <i>has</i> to run in
3580 * %parallel, and in fact this program (unlike @ref step_32 "step-32") makes no attempt at
3581 * all to do anything in %parallel using MPI. Nevertheless, Trilinos wants the
3582 * MPI system to be initialized. We do that be creating an object of type
3583 * Utilities::MPI::MPI_InitFinalize that initializes MPI (if available) using
3584 * the arguments given to main() (i.e., <code>argc</code> and
3585 * <code>argv</code>) and de-initializes it again when the object goes out of
3586 * scope.
3587 *
3588 * @code
3589 *   int main(int argc, char *argv[])
3590 *   {
3591 *   try
3592 *   {
3593 *   using namespace dealii;
3594 *   using namespace Step31;
3595 *  
3596 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(
3597 *   argc, argv, numbers::invalid_unsigned_int);
3598 *  
3599 * @endcode
3600 *
3601 * This program can only be run in serial. Otherwise, throw an exception.
3602 *
3603 * @code
3604 *   AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) == 1,
3605 *   ExcMessage(
3606 *   "This program can only be run in serial, use ./step-31"));
3607 *  
3608 *   BoussinesqFlowProblem<2> flow_problem;
3609 *   flow_problem.run();
3610 *   }
3611 *   catch (std::exception &exc)
3612 *   {
3613 *   std::cerr << std::endl
3614 *   << std::endl
3615 *   << "----------------------------------------------------"
3616 *   << std::endl;
3617 *   std::cerr << "Exception on processing: " << std::endl
3618 *   << exc.what() << std::endl
3619 *   << "Aborting!" << std::endl
3620 *   << "----------------------------------------------------"
3621 *   << std::endl;
3622 *  
3623 *   return 1;
3624 *   }
3625 *   catch (...)
3626 *   {
3627 *   std::cerr << std::endl
3628 *   << std::endl
3629 *   << "----------------------------------------------------"
3630 *   << std::endl;
3631 *   std::cerr << "Unknown exception!" << std::endl
3632 *   << "Aborting!" << std::endl
3633 *   << "----------------------------------------------------"
3634 *   << std::endl;
3635 *   return 1;
3636 *   }
3637 *  
3638 *   return 0;
3639 *   }
3640 * @endcode
3641<a name="step_31-Results"></a><h1>Results</h1>
3642
3643
3644<a name="step_31-Resultsin2d"></a><h3> Results in 2d </h3>
3645
3646
3647When you run the program in 2d, the output will look something like
3648this:
3649<code>
3650<pre>
3651Number of active cells: 256 (on 5 levels)
3652Number of degrees of freedom: 3556 (2178+289+1089)
3653
3654Timestep 0: t=0
3655 Assembling...
3656 Rebuilding Stokes preconditioner...
3657 Solving...
3658 0 GMRES iterations for Stokes subsystem.
3659 Time step: 0.919118
3660 9 CG iterations for temperature.
3661 Temperature range: -0.16687 1.30011
3662
3663Number of active cells: 280 (on 6 levels)
3664Number of degrees of freedom: 4062 (2490+327+1245)
3665
3666Timestep 0: t=0
3667 Assembling...
3668 Rebuilding Stokes preconditioner...
3669 Solving...
3670 0 GMRES iterations for Stokes subsystem.
3671 Time step: 0.459559
3672 9 CG iterations for temperature.
3673 Temperature range: -0.0982971 0.598503
3674
3675Number of active cells: 520 (on 7 levels)
3676Number of degrees of freedom: 7432 (4562+589+2281)
3677
3678Timestep 0: t=0
3679 Assembling...
3680 Rebuilding Stokes preconditioner...
3681 Solving...
3682 0 GMRES iterations for Stokes subsystem.
3683 Time step: 0.229779
3684 9 CG iterations for temperature.
3685 Temperature range: -0.0551098 0.294493
3686
3687Number of active cells: 1072 (on 8 levels)
3688Number of degrees of freedom: 15294 (9398+1197+4699)
3689
3690Timestep 0: t=0
3691 Assembling...
3692 Rebuilding Stokes preconditioner...
3693 Solving...
3694 0 GMRES iterations for Stokes subsystem.
3695 Time step: 0.11489
3696 9 CG iterations for temperature.
3697 Temperature range: -0.0273524 0.156861
3698
3699Number of active cells: 2116 (on 9 levels)
3700Number of degrees of freedom: 30114 (18518+2337+9259)
3701
3702Timestep 0: t=0
3703 Assembling...
3704 Rebuilding Stokes preconditioner...
3705 Solving...
3706 0 GMRES iterations for Stokes subsystem.
3707 Time step: 0.0574449
3708 9 CG iterations for temperature.
3709 Temperature range: -0.014993 0.0738328
3710
3711Timestep 1: t=0.0574449
3712 Assembling...
3713 Solving...
3714 56 GMRES iterations for Stokes subsystem.
3715 Time step: 0.0574449
3716 9 CG iterations for temperature.
3717 Temperature range: -0.0273934 0.14488
3718
3719...
3720</pre>
3721</code>
3722
3723In the beginning we refine the mesh several times adaptively and
3724always return to time step zero to restart on the newly refined
3725mesh. Only then do we start the actual time iteration.
3726
3727The program runs for a while. The temperature field for time steps 0,
3728500, 1000, 1500, 2000, 3000, 4000, and 5000 looks like this (note that
3729the color scale used for the temperature is not always the same):
3730
3731<table align="center" class="doxtable">
3732 <tr>
3733 <td>
3734 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.00.png" alt="">
3735 </td>
3736 <td>
3737 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.01.png" alt="">
3738 </td>
3739 <td>
3740 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.02.png" alt="">
3741 </td>
3742 <td>
3743 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.03.png" alt="">
3744 </td>
3745 </tr>
3746 <tr>
3747 <td>
3748 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.04.png" alt="">
3749 </td>
3750 <td>
3751 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.05.png" alt="">
3752 </td>
3753 <td>
3754 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.06.png" alt="">
3755 </td>
3756 <td>
3757 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.solution.07.png" alt="">
3758 </td>
3759 </tr>
3760</table>
3761
3762The visualizations shown here were generated using a version of the example
3763which did not enforce the constraints after transferring the mesh.
3764
3765As can be seen, we have three heat sources that heat fluid and
3766therefore produce a buoyancy effect that lets hots pockets of fluid
3767rise up and swirl around. By a chimney effect, the three streams are
3768pressed together by fluid that comes from the outside and wants to
3769join the updraft party. Note that because the fluid is initially at
3770rest, those parts of the fluid that were initially over the sources
3771receive a longer heating time than that fluid that is later dragged
3772over the source by the fully developed flow field. It is therefore
3773hotter, a fact that can be seen in the red tips of the three
3774plumes. Note also the relatively fine features of the flow field, a
3775result of the sophisticated transport stabilization of the temperature
3776equation we have chosen.
3777
3778In addition to the pictures above, the following ones show the
3779adaptive mesh and the flow field at the same time steps:
3780
3781<table align="center" class="doxtable">
3782 <tr>
3783 <td>
3784 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.00.png" alt="">
3785 </td>
3786 <td>
3787 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.01.png" alt="">
3788 </td>
3789 <td>
3790 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.02.png" alt="">
3791 </td>
3792 <td>
3793 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.03.png" alt="">
3794 </td>
3795 </tr>
3796 <tr>
3797 <td>
3798 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.04.png" alt="">
3799 </td>
3800 <td>
3801 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.05.png" alt="">
3802 </td>
3803 <td>
3804 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.06.png" alt="">
3805 </td>
3806 <td>
3807 <img src="https://www.dealii.org/images/steps/developer/step-31.2d.grid.07.png" alt="">
3808 </td>
3809 </tr>
3810</table>
3811
3812
3813<a name="step_31-Resultsin3d"></a><h3> Results in 3d </h3>
3814
3815
3816The same thing can of course be done in 3d by changing the template
3817parameter to the BoussinesqFlowProblem object in <code>main()</code>
3818from 2 to 3, so that the output now looks like follows:
3819
3820<code>
3821<pre>
3822Number of active cells: 64 (on 3 levels)
3823Number of degrees of freedom: 3041 (2187+125+729)
3824
3825Timestep 0: t=0
3826 Assembling...
3827 Rebuilding Stokes preconditioner...
3828 Solving...
3829 0 GMRES iterations for Stokes subsystem.
3830 Time step: 2.45098
3831 9 CG iterations for temperature.
3832 Temperature range: -0.675683 4.94725
3833
3834Number of active cells: 288 (on 4 levels)
3835Number of degrees of freedom: 12379 (8943+455+2981)
3836
3837Timestep 0: t=0
3838 Assembling...
3839 Rebuilding Stokes preconditioner...
3840 Solving...
3841 0 GMRES iterations for Stokes subsystem.
3842 Time step: 1.22549
3843 9 CG iterations for temperature.
3844 Temperature range: -0.527701 2.25764
3845
3846Number of active cells: 1296 (on 5 levels)
3847Number of degrees of freedom: 51497 (37305+1757+12435)
3848
3849Timestep 0: t=0
3850 Assembling...
3851 Rebuilding Stokes preconditioner...
3852 Solving...
3853 0 GMRES iterations for Stokes subsystem.
3854 Time step: 0.612745
3855 10 CG iterations for temperature.
3856 Temperature range: -0.496942 0.847395
3857
3858Number of active cells: 5048 (on 6 levels)
3859Number of degrees of freedom: 192425 (139569+6333+46523)
3860
3861Timestep 0: t=0
3862 Assembling...
3863 Rebuilding Stokes preconditioner...
3864 Solving...
3865 0 GMRES iterations for Stokes subsystem.
3866 Time step: 0.306373
3867 10 CG iterations for temperature.
3868 Temperature range: -0.267683 0.497739
3869
3870Timestep 1: t=0.306373
3871 Assembling...
3872 Solving...
3873 27 GMRES iterations for Stokes subsystem.
3874 Time step: 0.306373
3875 10 CG iterations for temperature.
3876 Temperature range: -0.461787 0.958679
3877
3878...
3879</pre>
3880</code>
3881
3882Visualizing the temperature isocontours at time steps 0,
388350, 100, 150, 200, 300, 400, 500, 600, 700, and 800 yields the
3884following plots:
3885
3886<table align="center" class="doxtable">
3887 <tr>
3888 <td>
3889 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.00.png" alt="">
3890 </td>
3891 <td>
3892 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.01.png" alt="">
3893 </td>
3894 <td>
3895 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.02.png" alt="">
3896 </td>
3897 <td>
3898 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.03.png" alt="">
3899 </td>
3900 </tr>
3901 <tr>
3902 <td>
3903 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.04.png" alt="">
3904 </td>
3905 <td>
3906 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.05.png" alt="">
3907 </td>
3908 <td>
3909 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.06.png" alt="">
3910 </td>
3911 <td>
3912 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.07.png" alt="">
3913 </td>
3914 </tr>
3915 <tr>
3916 <td>
3917 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.08.png" alt="">
3918 </td>
3919 <td>
3920 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.09.png" alt="">
3921 </td>
3922 <td>
3923 <img src="https://www.dealii.org/images/steps/developer/step-31.3d.solution.10.png" alt="">
3924 </td>
3925 <td>
3926 </td>
3927 </tr>
3928</table>
3929
3930That the first picture looks like three hedgehogs stems from the fact that our
3931scheme essentially projects the source times the first time step size onto the
3932mesh to obtain the temperature field in the first time step. Since the source
3933function is discontinuous, we need to expect over- and undershoots from this
3934project. This is in fact what happens (it's easier to check this in 2d) and
3935leads to the crumpled appearance of the isosurfaces. The visualizations shown
3936here were generated using a version of the example which did not enforce the
3937constraints after transferring the mesh.
3938
3939
3940
3941<a name="step_31-Numericalexperimentstodetermineoptimalparameters"></a><h3> Numerical experiments to determine optimal parameters </h3>
3942
3943
3944The program as is has three parameters that we don't have much of a
3945theoretical handle on how to choose in an optimal way. These are:
3946<ul>
3947 <li>The time step must satisfy a CFL condition
3948 @f$k\le \min_K \frac{c_kh_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$. Here, @f$c_k@f$ is
3949 dimensionless, but what is the right value?
3950 <li>In the computation of the artificial viscosity,
3951@f{eqnarray*}{
3952 \nu_\alpha(T)|_K
3953 =
3954 \beta
3955 \|\mathbf{u}\|_{L^\infty(K)}
3956 \min\left\{
3957 h_K,
3958 h_K^\alpha
3959 \frac{\|R_\alpha(T)\|_{L^\infty(K)}}{c(\mathbf{u},T)}
3960 \right\},
3961@f}
3962 with @f$c(\mathbf{u},T) =
3963 c_R\ \|\mathbf{u}\|_{L^\infty(\Omega)} \ \mathrm{var}(T)
3964 \ |\mathrm{diam}(\Omega)|^{\alpha-2}@f$.
3965 Here, the choice of the dimensionless %numbers @f$\beta,c_R@f$ is of
3966 interest.
3967</ul>
3968In all of these cases, we will have to expect that the correct choice of each
3969value depends on that of the others, and most likely also on the space
3970dimension and polynomial degree of the finite element used for the
3971temperature. Below we'll discuss a few numerical experiments to choose
3972constants @f$c_k@f$ and @f$\beta@f$.
3973
3974Below, we will not discuss the choice of @f$c_R@f$. In the program, we set
3975it to @f$c_R=2^{\frac{4-2\alpha}{d}}@f$. The reason for this value is a
3976bit complicated and has more to do with the history of the program
3977than reasoning: while the correct formula for the global scaling
3978parameter @f$c(\mathbf{u},T)@f$ is shown above, the program (including the
3979version shipped with deal.II 6.2) initially had a bug in that we
3980computed
3981@f$c(\mathbf{u},T) =
3982 \|\mathbf{u}\|_{L^\infty(\Omega)} \ \mathrm{var}(T)
3983 \ \frac{1}{|\mathrm{diam}(\Omega)|^{\alpha-2}}@f$ instead, where
3984we had set the scaling parameter to one. Since we only computed on the
3985unit square/cube where @f$\mathrm{diam}(\Omega)=2^{1/d}@f$, this was
3986entirely equivalent to using the correct formula with
3987@f$c_R=\left(2^{1/d}\right)^{4-2\alpha}=2^{\frac{4-2\alpha}{d}}@f$. Since
3988this value for @f$c_R@f$ appears to work just fine for the current
3989program, we corrected the formula in the program and set @f$c_R@f$ to a
3990value that reproduces exactly the results we had before. We will,
3991however, revisit this issue again in @ref step_32 "step-32".
3992
3993Now, however, back to the discussion of what values of @f$c_k@f$ and
3994@f$\beta@f$ to choose:
3995
3996
3997<a name="step_31-Choosingicsubksubiandbeta"></a><h4> Choosing <i>c<sub>k</sub></i> and beta </h4>
3998
3999
4000These two constants are definitely linked in some way. The reason is easy to
4001see: In the case of a pure advection problem,
4002@f$\frac{\partial T}{\partial t} + \mathbf{u}\cdot\nabla T = \gamma@f$, any
4003explicit scheme has to satisfy a CFL condition of the form
4004@f$k\le \min_K \frac{c_k^a h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$. On the other hand,
4005for a pure diffusion problem,
4006@f$\frac{\partial T}{\partial t} + \nu \Delta T = \gamma@f$,
4007explicit schemes need to satisfy a condition
4008@f$k\le \min_K \frac{c_k^d h_K^2}{\nu}@f$. So given the form of @f$\nu@f$ above, an
4009advection diffusion problem like the one we have to solve here will result in
4010a condition of the form
4011@f$
4012k\le \min_K \min \left\{
4013 \frac{c_k^a h_K}{\|\mathbf{u}\|_{L^\infty(K)}},
4014 \frac{c_k^d h_K^2}{\beta \|\mathbf{u}\|_{L^\infty(K)} h_K}\right\}
4015 =
4016 \min_K \left( \min \left\{
4017 c_k^a,
4018 \frac{c_k^d}{\beta}\right\}
4019 \frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}} \right)
4020@f$.
4021It follows that we have to face the fact that we might want to choose @f$\beta@f$
4022larger to improve the stability of the numerical scheme (by increasing the
4023amount of artificial diffusion), but we have to pay a price in the form of
4024smaller, and consequently more time steps. In practice, one would therefore
4025like to choose @f$\beta@f$ as small as possible to keep the transport problem
4026sufficiently stabilized while at the same time trying to choose the time step
4027as large as possible to reduce the overall amount of work.
4028
4029The find the right balance, the only way is to do a few computational
4030experiments. Here's what we did: We modified the program slightly to allow
4031less mesh refinement (so we don't always have to wait that long) and to choose
4032@f$
4033 \nu(T)|_K
4034 =
4035 \beta
4036 \|\mathbf{u}\|_{L^\infty(K)} h_K
4037@f$ to eliminate the effect of the constant @f$c_R@f$ (we know that
4038solutions are stable by using this version of @f$\nu(T)@f$ as an artificial
4039viscosity, but that we can improve things -- i.e. make the solution
4040sharper -- by using the more complicated formula for this artificial
4041viscosity). We then run the program
4042for different values @f$c_k,\beta@f$ and observe maximal and minimal temperatures
4043in the domain. What we expect to see is this: If we choose the time step too
4044big (i.e. choose a @f$c_k@f$ bigger than theoretically allowed) then we will get
4045exponential growth of the temperature. If we choose @f$\beta@f$ too small, then
4046the transport stabilization becomes insufficient and the solution will show
4047significant oscillations but not exponential growth.
4048
4049
4050<a name="step_31-ResultsforQsub1subelements"></a><h5>Results for Q<sub>1</sub> elements</h5>
4051
4052
4053Here is what we get for
4054@f$\beta=0.01, \beta=0.1@f$, and @f$\beta=0.5@f$, different choices of @f$c_k@f$, and
4055bilinear elements (<code>temperature_degree=1</code>) in 2d:
4056
4057<table align="center" class="doxtable">
4058 <tr>
4059 <td>
4060 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q1.beta=0.01.png" alt="">
4061 </td>
4062 <td>
4063 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q1.beta=0.03.png" alt="">
4064 </td>
4065 </tr>
4066 <tr>
4067 <td>
4068 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q1.beta=0.1.png" alt="">
4069 </td>
4070 <td>
4071 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q1.beta=0.5.png" alt="">
4072 </td>
4073 </tr>
4074</table>
4075
4076The way to interpret these graphs goes like this: for @f$\beta=0.01@f$ and
4077@f$c_k=\frac 12,\frac 14@f$, we see exponential growth or at least large
4078variations, but if we choose
4079@f$k=\frac 18\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$
4080or smaller, then the scheme is
4081stable though a bit wobbly. For more artificial diffusion, we can choose
4082@f$k=\frac 14\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$
4083or smaller for @f$\beta=0.03@f$,
4084@f$k=\frac 13\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$
4085or smaller for @f$\beta=0.1@f$, and again need
4086@f$k=\frac 1{15}\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$
4087for @f$\beta=0.5@f$ (this time because much diffusion requires a small time
4088step).
4089
4090So how to choose? If we were simply interested in a large time step, then we
4091would go with @f$\beta=0.1@f$ and
4092@f$k=\frac 13\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$.
4093On the other hand, we're also interested in accuracy and here it may be of
4094interest to actually investigate what these curves show. To this end note that
4095we start with a zero temperature and that our sources are positive &mdash; so
4096we would intuitively expect that the temperature can never drop below
4097zero. But it does, a consequence of Gibb's phenomenon when using continuous
4098elements to approximate a discontinuous solution. We can therefore see that
4099choosing @f$\beta@f$ too small is bad: too little artificial diffusion leads to
4100over- and undershoots that aren't diffused away. On the other hand, for large
4101@f$\beta@f$, the minimum temperature drops below zero at the beginning but then
4102quickly diffuses back to zero.
4103
4104On the other hand, let's also look at the maximum temperature. Watching the
4105movie of the solution, we see that initially the fluid is at rest. The source
4106keeps heating the same volume of fluid whose temperature increases linearly at
4107the beginning until its buoyancy is able to move it upwards. The hottest part
4108of the fluid is therefore transported away from the solution and fluid taking
4109its place is heated for only a short time before being moved out of the source
4110region, therefore remaining cooler than the initial bubble. If @f$\kappa=0@f$
4111(in the program it is nonzero but very small) then the hottest part of the
4112fluid should be advected along with the flow with its temperature
4113constant. That's what we can see in the graphs with the smallest @f$\beta@f$: Once
4114the maximum temperature is reached, it hardly changes any more. On the other
4115hand, the larger the artificial diffusion, the more the hot spot is
4116diffused. Note that for this criterion, the time step size does not play a
4117significant role.
4118
4119So to sum up, likely the best choice would appear to be @f$\beta=0.03@f$
4120and @f$k=\frac 14\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$. The curve is
4121a bit wobbly, but overall pictures looks pretty reasonable with the
4122exception of some over and undershoots close to the start time due to
4123Gibb's phenomenon.
4124
4125
4126<a name="step_31-ResultsforQsub2subelements"></a><h5>Results for Q<sub>2</sub> elements</h5>
4127
4128
4129One can repeat the same sequence of experiments for higher order
4130elements as well. Here are the graphs for bi-quadratic shape functions
4131(<code>temperature_degree=2</code>) for the temperature, while we
4132retain the @f$Q_2/Q_1@f$ stable Taylor-Hood element for the Stokes system:
4133
4134<table align="center" class="doxtable">
4135 <tr>
4136 <td>
4137 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q2.beta=0.01.png" alt="">
4138 </td>
4139 <td>
4140 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q2.beta=0.03.png" alt="">
4141 </td>
4142 </tr>
4143 <tr>
4144 <td>
4145 <img src="https://www.dealii.org/images/steps/developer/step-31.timestep.q2.beta=0.1.png" alt="">
4146 </td>
4147 </tr>
4148</table>
4149
4150Again, small values of @f$\beta@f$ lead to less diffusion but we have to
4151choose the time step very small to keep things under control. Too
4152large values of @f$\beta@f$ make for more diffusion, but again require
4153small time steps. The best value would appear to be @f$\beta=0.03@f$, as
4154for the @f$Q_1@f$ element, and then we have to choose
4155@f$k=\frac 18\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$ &mdash; exactly
4156half the size for the @f$Q_1@f$ element, a fact that may not be surprising
4157if we state the CFL condition as the requirement that the time step be
4158small enough so that the distance transport advects in each time step
4159is no longer than one <i>grid point</i> away (which for @f$Q_1@f$ elements
4160is @f$h_K@f$, but for @f$Q_2@f$ elements is @f$h_K/2@f$). It turns out that @f$\beta@f$
4161needs to be slightly larger for obtaining stable results also late in
4162the simulation at times larger than 60, so we actually choose it as
4163@f$\beta = 0.034@f$ in the code.
4164
4165
4166<a name="step_31-Resultsfor3d"></a><h5>Results for 3d</h5>
4167
4168
4169One can repeat these experiments in 3d and find the optimal time step
4170for each value of @f$\beta@f$ and find the best value of @f$\beta@f$. What one
4171finds is that for the same @f$\beta@f$ already used in 2d, the time steps
4172needs to be a bit smaller, by around a factor of 1.2 or so. This is
4173easily explained: the time step restriction is
4174@f$k=\min_K \frac{ch_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$ where @f$h_K@f$ is
4175the <i>diameter</i> of the cell. However, what is really needed is the
4176distance between mesh points, which is @f$\frac{h_K}{\sqrt{d}}@f$. So a
4177more appropriate form would be
4178@f$k=\min_K \frac{ch_K}{\|\mathbf{u}\|_{L^\infty(K)}\sqrt{d}}@f$.
4179
4180The second find is that one needs to choose @f$\beta@f$ slightly bigger
4181(about @f$\beta=0.05@f$ or so). This then again reduces the time step we
4182can take.
4183
4184
4185
4186
4187<a name="step_31-Conclusions"></a><h5>Conclusions</h5>
4188
4189
4190Concluding, from the simple computations above, @f$\beta=0.034@f$ appears to be a
4191good choice for the stabilization parameter in 2d, and @f$\beta=0.05@f$ in 3d. In
4192a dimension independent way, we can model this as @f$\beta=0.017d@f$. If one does
4193longer computations (several thousand time steps) on finer meshes, one
4194realizes that the time step size is not quite small enough and that for
4195stability one will have to reduce the above values a bit more (by about a
4196factor of @f$\frac 78@f$).
4197
4198As a consequence, a formula that reconciles 2d, 3d, and variable polynomial
4199degree and takes all factors in account reads as follows:
4200@f{eqnarray*}{
4201 k =
4202 \frac 1{2 \cdot 1.7} \frac 1{\sqrt{d}}
4203 \frac 2d
4204 \frac 1{q_T}
4205 \frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}
4206 =
4207 \frac 1{1.7 d\sqrt{d}}
4208 \frac 1{q_T}
4209 \frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}.
4210@f}
4211In the first form (in the center of the equation), @f$\frac
42121{2 \cdot 1.7}@f$ is a universal constant, @f$\frac 1{\sqrt{d}}@f$
4213is the factor that accounts for the difference between cell diameter
4214and grid point separation,
4215@f$\frac 2d@f$ accounts for the increase in @f$\beta@f$ with space dimension,
4216@f$\frac 1{q_T}@f$ accounts for the distance between grid points for
4217higher order elements, and @f$\frac{h_K}{\|\mathbf{u}\|_{L^\infty(K)}}@f$
4218for the local speed of transport relative to the cell size. This is
4219the formula that we use in the program.
4220
4221As for the question of whether to use @f$Q_1@f$ or @f$Q_2@f$ elements for the
4222temperature, the following considerations may be useful: First,
4223solving the temperature equation is hardly a factor in the overall
4224scheme since almost the entire compute time goes into solving the
4225Stokes system in each time step. Higher order elements for the
4226temperature equation are therefore not a significant drawback. On the
4227other hand, if one compares the size of the over- and undershoots the
4228solution produces due to the discontinuous source description, one
4229notices that for the choice of @f$\beta@f$ and @f$k@f$ as above, the @f$Q_1@f$
4230solution dips down to around @f$-0.47@f$, whereas the @f$Q_2@f$ solution only
4231goes to @f$-0.13@f$ (remember that the exact solution should never become
4232negative at all. This means that the @f$Q_2@f$ solution is significantly
4233more accurate; the program therefore uses these higher order elements,
4234despite the penalty we pay in terms of smaller time steps.
4235
4236
4237<a name="step_31-Possibilitiesforextensions"></a><h3> Possibilities for extensions </h3>
4238
4239
4240There are various ways to extend the current program. Of particular interest
4241is, of course, to make it faster and/or increase the resolution of the
4242program, in particular in 3d. This is the topic of the @ref step_32 "step-32"
4243tutorial program which will implement strategies to solve this problem in
4244%parallel on a cluster. It is also the basis of the much larger open
4245source code ASPECT (see https://aspect.geodynamics.org/ ) that can solve realistic
4246problems and that constitutes the further development of @ref step_32 "step-32".
4247
4248Another direction would be to make the fluid flow more realistic. The program
4249was initially written to simulate various cases simulating the convection of
4250material in the earth's mantle, i.e. the zone between the outer earth core and
4251the solid earth crust: there, material is heated from below and cooled from
4252above, leading to thermal convection. The physics of this fluid are much more
4253complicated than shown in this program, however: The viscosity of mantle
4254material is strongly dependent on the temperature, i.e. @f$\eta=\eta(T)@f$, with
4255the dependency frequently modeled as a viscosity that is reduced exponentially
4256with rising temperature. Secondly, much of the dynamics of the mantle is
4257determined by chemical reactions, primarily phase changes of the various
4258crystals that make up the mantle; the buoyancy term on the right hand side of
4259the Stokes equations then depends not only on the temperature, but also on the
4260chemical composition at a given location which is advected by the flow field
4261but also changes as a function of pressure and temperature. We will
4262investigate some of these effects in later tutorial programs as well.
4263 *
4264 *
4265<a name="step_31-PlainProg"></a>
4266<h1> The plain program</h1>
4267@include "step-31.cc"
4268*/
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
Definition fe_q.h:554
const std::vector< Point< dim > > & get_unit_support_points() const
const unsigned int n_components
Definition function.h:163
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
Definition point.h:111
float depth
Point< 3 > vertices[4]
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1204
const Event initial
Definition event.cc:64
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)
Expression sign(const Expression &x)
void downstream(DoFHandler< dim, spacedim > &dof_handler, const Tensor< 1, spacedim > &direction, const bool dof_wise_renumbering=false)
void interpolate(const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max(), const VectorTools::NormType norm_type=VectorTools::L1_norm)
double diameter(const Triangulation< dim, spacedim > &tria)
@ matrix
Contents is actually a matrix.
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
constexpr ReturnType< rank, T >::value_type & extract(T &t, const ArrayType &indices)
VectorType::value_type * end(VectorType &V)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
unsigned int n_active_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition tria.cc:14889
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)
Definition loop.h:70
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
STL namespace.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation