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