deal.II version 9.7.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-97.h
Go to the documentation of this file.
1
1541 *  
1542 * @endcode
1543 *
1544 * These three get-functions are used to channel the mesh and the solution
1545 * to the next solver.
1546 *
1547 * @code
1548 *   const Triangulation<3> &get_tria() const
1549 *   {
1550 *   return triangulation;
1551 *   }
1552 *   const DoFHandler<3> &get_dof_handler() const
1553 *   {
1554 *   return dof_handler;
1555 *   }
1556 *   const Vector<double> &get_solution() const
1557 *   {
1558 *   return solution;
1559 *   }
1560 *  
1561 *   private:
1562 * @endcode
1563 *
1564 * The following data members are typical for all deal.II simulations:
1565 * triangulation, finite elements, dof handlers, etc. The constraints
1566 * are used to enforce the Dirichlet boundary conditions. The names of the
1567 * data members are self-explanatory.
1568 *
1569
1570 *
1571 *
1572 * @code
1573 *   Triangulation<3> triangulation;
1574 *  
1575 *   const FE_Nedelec<3> fe;
1576 *   DoFHandler<3> dof_handler;
1577 *   Vector<double> solution;
1578 *  
1579 *   SparseMatrix<double> system_matrix;
1580 *   Vector<double> system_rhs;
1581 *  
1582 *   AffineConstraints<double> constraints;
1583 *   SparsityPattern sparsity_pattern;
1584 *  
1585 *   SphericalManifold<3> sphere;
1586 *  
1587 *   const unsigned int refinement_parameter;
1588 *   const unsigned int mapping_degree;
1589 *   const double eta_squared;
1590 *   const std::string fname;
1591 *   TimerOutput timer;
1592 *  
1593 * @endcode
1594 *
1595 * The program utilizes the WorkStream technology. The @ref step_9 "step-9" tutorial
1596 * does a much better job of explaining the workings of WorkStream.
1597 * Reading the "WorkStream paper", see the glossary, is recommended.
1598 * The following structures and functions are related to WorkStream.
1599 *
1600 * @code
1601 *   struct AssemblyScratchData
1602 *   {
1603 *   AssemblyScratchData(const FiniteElement<3> &fe,
1604 *   const double eta_squared,
1605 *   const unsigned int mapping_degree);
1606 *  
1607 *   AssemblyScratchData(const AssemblyScratchData &scratch_data);
1608 *  
1609 *   const FreeCurrentDensity Jf;
1610 *  
1611 *   MappingQ<3> mapping;
1612 *   FEValues<3> fe_values;
1613 *  
1614 *   const unsigned int dofs_per_cell;
1615 *   const unsigned int n_q_points;
1616 *  
1617 *   std::vector<Tensor<1, 3>> Jf_list;
1618 *  
1619 *   const double eta_squared;
1620 *   const FEValuesExtractors::Vector ve;
1621 *   };
1622 *  
1623 *   struct AssemblyCopyData
1624 *   {
1626 *   Vector<double> cell_rhs;
1627 *   std::vector<types::global_dof_index> local_dof_indices;
1628 *   };
1629 *  
1630 *   void system_matrix_local(
1631 *   const typename DoFHandler<3>::active_cell_iterator &cell,
1632 *   AssemblyScratchData &scratch_data,
1633 *   AssemblyCopyData &copy_data);
1634 *  
1635 *   void copy_local_to_global(const AssemblyCopyData &copy_data);
1636 *   };
1637 *  
1638 *   Solver::Solver(const unsigned int p,
1639 *   const unsigned int r,
1640 *   const unsigned int mapping_degree,
1641 *   const double eta_squared,
1642 *   const std::string &fname)
1643 *   : fe(p)
1644 *   , refinement_parameter(r)
1645 *   , mapping_degree(mapping_degree)
1646 *   , eta_squared(eta_squared)
1647 *   , fname(fname)
1648 *   , timer(std::cout,
1649 *   (Settings::print_time_tables) ? TimerOutput::summary :
1650 *   TimerOutput::never,
1652 *   {}
1653 *  
1654 * @endcode
1655 *
1656 * The following function loads the mesh, assigns material IDs to all cells,
1657 * and attaches the spherical manifold to the mesh. The material IDs are
1658 * assigned on the basis of the distance from the center of a cell to the
1659 * origin. The spherical manifold is attached to a face if all vertices of
1660 * the face are at the same distance from the origin provided the cell is
1661 * outside the cube in the center of the mesh, see mesh description in the
1662 * introduction.
1663 *
1664 * @code
1665 *   void Solver::make_mesh()
1666 *   {
1667 *   GridIn<3> gridin;
1668 *  
1669 *   gridin.attach_triangulation(triangulation);
1670 *   std::ifstream ifs("sphere_r" + std::to_string(refinement_parameter) +
1671 *   ".msh");
1672 *   gridin.read_msh(ifs);
1673 *  
1674 *   triangulation.reset_all_manifolds();
1675 *  
1676 *   for (auto cell : triangulation.active_cell_iterators())
1677 *   {
1678 *   cell->set_material_id(
1679 *   Settings::material_id_free_space); // The cell is in free space.
1680 *  
1681 *   if ((cell->center().norm() > Settings::a1) &&
1682 *   (cell->center().norm() < Settings::b1))
1683 *   cell->set_material_id(
1684 *   Settings::material_id_core); // The cell is inside the core.
1685 *  
1686 *   if ((cell->center().norm() > Settings::a2) &&
1687 *   (cell->center().norm() < Settings::b2))
1688 *   cell->set_material_id(
1689 *   Settings::material_id_free_current); /* The cell is inside the Jf
1690 *   region. */
1691 *  
1692 *   for (unsigned int f = 0; f < cell->n_faces(); f++)
1693 *   {
1694 *   double dif_norm = 0.0;
1695 *   for (unsigned int v = 1; v < cell->face(f)->n_vertices(); v++)
1696 *   dif_norm += std::abs(cell->face(f)->vertex(0).norm() -
1697 *   cell->face(f)->vertex(v).norm());
1698 *  
1699 *   if ((dif_norm < Settings::eps) &&
1700 *   (cell->center().norm() > Settings::d1))
1701 *   cell->face(f)->set_all_manifold_ids(1);
1702 *   }
1703 *   }
1704 *  
1705 *   triangulation.set_manifold(1, sphere);
1706 *   }
1707 *  
1708 * @endcode
1709 *
1710 * This function initializes the dofs, applies the Dirichlet boundary
1711 * condition, and initializes the vectors and matrices.
1712 *
1713 * @code
1714 *   void Solver::setup()
1715 *   {
1716 *   dof_handler.reinit(triangulation);
1717 *   dof_handler.distribute_dofs(fe);
1718 *  
1719 * @endcode
1720 *
1721 * The following segment of the code applies the homogeneous Dirichlet
1722 * boundary condition. As discussed in the introduction, the Dirichlet
1723 * boundary condition is an essential condition and must be enforced
1724 * by constraining the system matrix. This segment of the code does the
1725 * constraining.
1726 *
1727 * @code
1728 *   constraints.clear();
1729 *  
1730 *   DoFTools::make_hanging_node_constraints(dof_handler, constraints);
1731 *  
1733 *   dof_handler,
1734 *   0, // The first vector component.
1736 *   Settings::boundary_id_infinity,
1737 *   constraints,
1738 *   MappingQ<3>(mapping_degree));
1739 *  
1740 *   constraints.close();
1741 *  
1742 * @endcode
1743 *
1744 * The rest of the function arranges the dofs in a sparsity pattern and
1745 * initializes the system matrices and vectors.
1746 *
1747 * @code
1748 *   DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
1749 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints, false);
1750 *   sparsity_pattern.copy_from(dsp);
1751 *  
1752 *   system_matrix.reinit(sparsity_pattern);
1753 *   solution.reinit(dof_handler.n_dofs());
1754 *   system_rhs.reinit(dof_handler.n_dofs());
1755 *   }
1756 *  
1757 * @endcode
1758 *
1759 * Formally, the following function assembles the system of linear equations.
1760 * In reality, however, it just spells all the magic words to get the
1761 * WorkStream going. The interesting part, i.e., the actual assembling of the
1762 * system matrix and the right-hand side, happens below in the function
1763 * Solver::system_matrix_local().
1764 *
1765 * @code
1766 *   void Solver::assemble()
1767 *   {
1768 *   WorkStream::run(dof_handler.begin_active(),
1769 *   dof_handler.end(),
1770 *   *this,
1771 *   &Solver::system_matrix_local,
1772 *   &Solver::copy_local_to_global,
1773 *   AssemblyScratchData(fe, eta_squared, mapping_degree),
1774 *   AssemblyCopyData());
1775 *   }
1776 *  
1777 * @endcode
1778 *
1779 * The following two constructors initialize scratch data from the input
1780 * parameters and from another object of the same type, i.e., a copy
1781 * constructor.
1782 *
1783 * @code
1784 *   Solver::AssemblyScratchData::AssemblyScratchData(
1785 *   const FiniteElement<3> &fe,
1786 *   const double eta_squared,
1787 *   const unsigned int mapping_degree)
1788 *   : Jf()
1789 *   , mapping(mapping_degree)
1790 *   , fe_values(mapping,
1791 *   fe,
1792 *   QGauss<3>(fe.degree + 2),
1794 *   update_JxW_values)
1795 *   , dofs_per_cell(fe_values.dofs_per_cell)
1796 *   , n_q_points(fe_values.get_quadrature().size())
1797 *   , Jf_list(n_q_points, Tensor<1, 3>())
1798 *   , eta_squared(eta_squared)
1799 *   , ve(0)
1800 *   {}
1801 *  
1802 *   Solver::AssemblyScratchData::AssemblyScratchData(
1803 *   const AssemblyScratchData &scratch_data)
1804 *   : Jf()
1805 *   , mapping(scratch_data.mapping.get_degree())
1806 *   , fe_values(mapping,
1807 *   scratch_data.fe_values.get_fe(),
1808 *   scratch_data.fe_values.get_quadrature(),
1810 *   update_JxW_values)
1811 *   , dofs_per_cell(fe_values.dofs_per_cell)
1812 *   , n_q_points(fe_values.get_quadrature().size())
1813 *   , Jf_list(n_q_points, Tensor<1, 3>())
1814 *   , eta_squared(scratch_data.eta_squared)
1815 *   , ve(0)
1816 *   {}
1817 *  
1818 * @endcode
1819 *
1820 * This function assembles a fraction of the system matrix and the system
1821 * right-hand side related to a single cell. These fractions are
1822 * `copy_data.cell_matrix` and `copy_data.cell_rhs`. They are copied into
1823 * the system matrix, @f$A_{ij}@f$, and the right-hand side, @f$b_i@f$, by the
1824 * function `Solver::copy_local_to_global()`.
1825 *
1826 * @code
1827 *   void Solver::system_matrix_local(
1828 *   const typename DoFHandler<3>::active_cell_iterator &cell,
1829 *   AssemblyScratchData &scratch_data,
1830 *   AssemblyCopyData &copy_data)
1831 *   {
1832 * @endcode
1833 *
1834 * First, we reinitialize the matrices and vectors related to the current
1835 * cell and compute the FE values.
1836 *
1837 * @code
1838 *   copy_data.cell_matrix.reinit(scratch_data.dofs_per_cell,
1839 *   scratch_data.dofs_per_cell);
1840 *  
1841 *   copy_data.cell_rhs.reinit(scratch_data.dofs_per_cell);
1842 *  
1843 *   copy_data.local_dof_indices.resize(scratch_data.dofs_per_cell);
1844 *  
1845 *   scratch_data.fe_values.reinit(cell);
1846 *  
1847 * @endcode
1848 *
1849 * Second, we compute the free-current density, @f$\vec{J}_f@f$, at the
1850 * quadrature points.
1851 *
1852 * @code
1853 *   scratch_data.Jf.value_list(scratch_data.fe_values.get_quadrature_points(),
1854 *   cell->material_id(),
1855 *   scratch_data.Jf_list);
1856 *  
1857 * @endcode
1858 *
1859 * Third, we compute the components of the cell matrix and cell right-hand
1860 * side. The labels of the integrals are the same as in the introduction to
1861 * this tutorial.
1862 *
1863 * @code
1864 *   for (unsigned int q_index = 0; q_index < scratch_data.n_q_points; ++q_index)
1865 *   {
1866 *   for (unsigned int i = 0; i < scratch_data.dofs_per_cell; ++i)
1867 *   {
1868 *   for (unsigned int j = 0; j < scratch_data.dofs_per_cell; ++j)
1869 *   {
1870 *   copy_data.cell_matrix(i, j) += // Integral I_a1+I_a3.
1871 *   (scratch_data.fe_values[scratch_data.ve].curl(
1872 *   i, q_index) * // curl N_i
1873 *   scratch_data.fe_values[scratch_data.ve].curl(
1874 *   j, q_index) // curl N_j
1875 *   + scratch_data.eta_squared * // eta^2
1876 *   scratch_data.fe_values[scratch_data.ve].value(
1877 *   i, q_index) * // N_i
1878 *   scratch_data.fe_values[scratch_data.ve].value(
1879 *   j, q_index) // N_j
1880 *   ) *
1881 *   scratch_data.fe_values.JxW(q_index); // dV
1882 *   }
1883 *   copy_data.cell_rhs(i) += // Integral I_b3-1.
1884 *   (scratch_data.Jf_list[q_index] *
1885 *   scratch_data.fe_values[scratch_data.ve].curl(i, q_index)) *
1886 *   scratch_data.fe_values.JxW(q_index); // J_f.(curl N_i)dV.
1887 *   }
1888 *   }
1889 *  
1890 * @endcode
1891 *
1892 * Finally, we query the dof indices on the current cell and store them
1893 * in the copy data structure, so we know to which locations of the system
1894 * matrix and right-hand side the components of the cell matrix and
1895 * cell right-hand side must be copied.
1896 *
1897 * @code
1898 *   cell->get_dof_indices(copy_data.local_dof_indices);
1899 *   }
1900 *  
1901 * @endcode
1902 *
1903 * This function copies the components of a cell matrix and a cell right-hand
1904 * side into the system matrix, @f$A_{i,j}@f$, and the system right-hand side,
1905 * @f$b_i@f$.
1906 *
1907 * @code
1908 *   void Solver::copy_local_to_global(const AssemblyCopyData &copy_data)
1909 *   {
1910 *   constraints.distribute_local_to_global(copy_data.cell_matrix,
1911 *   copy_data.cell_rhs,
1912 *   copy_data.local_dof_indices,
1913 *   system_matrix,
1914 *   system_rhs);
1915 *   }
1916 *  
1917 * @endcode
1918 *
1919 * This function solves the system of linear equations.
1920 * If `Settings::log_cg_convergence == true`, the convergence data is saved
1921 * into a file. In theory, a CG solver can solve an @f$m \times m@f$ system of
1922 * linear equations in at most @f$m@f$ steps. In practice, it can take more steps
1923 * to converge. The convergence of the algorithm depends on the spectral
1924 * properties of the system matrix. The best case is if the eigenvalues form a
1925 * compact cluster away from zero. In our case, however, the eigenvalues are
1926 * spread in between zero and the maximal eigenvalue. Consequently, we expect
1927 * a poor convergence and increase the maximal number of iteration steps by a
1928 * factor of 10, i.e., `10*system_rhs.size()`. The stopping condition is \f[
1929 * |\boldsymbol{b} - \boldsymbol{A}\boldsymbol{c}|
1930 * < 10^{-6} |\boldsymbol{b}|.
1931 * \f]
1932 * As soon as we use constraints, we must not forget to distribute them.
1933 *
1934 * @code
1935 *   void Solver::solve()
1936 *   {
1937 *   SolverControl control(10 * system_rhs.size(),
1938 *   1.0e-6 * system_rhs.l2_norm(),
1939 *   false,
1940 *   false);
1941 *  
1942 *   if (Settings::log_cg_convergence)
1943 *   control.enable_history_data();
1944 *  
1946 *   SolverCG<Vector<double>> cg(control, memory);
1947 *  
1948 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
1949 *   preconditioner.initialize(system_matrix, 1.2);
1950 *  
1951 *   cg.solve(system_matrix, solution, system_rhs, preconditioner);
1952 *  
1953 *   constraints.distribute(solution);
1954 *  
1955 *   if (Settings::log_cg_convergence)
1956 *   {
1957 *   const std::vector<double> history_data = control.get_history_data();
1958 *  
1959 *   std::ofstream ofs(fname + "_cg_convergence.csv");
1960 *  
1961 *   for (unsigned int i = 1; i < history_data.size(); i++)
1962 *   ofs << i << ", " << history_data[i] << "\n";
1963 *   }
1964 *   }
1965 *  
1966 * @endcode
1967 *
1968 * This function saves the computed current vector potential into a vtu file.
1969 *
1970 * @code
1971 *   void Solver::save() const
1972 *   {
1973 *   const std::vector<std::string> solution_names(3, "VectorField");
1974 *   const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1975 *   interpretation(3,
1977 *  
1978 *   DataOut<3> data_out;
1979 *  
1980 *   data_out.add_data_vector(dof_handler,
1981 *   solution,
1982 *   solution_names,
1983 *   interpretation);
1984 *  
1985 *   DataOutBase::VtkFlags flags;
1986 *   flags.write_higher_order_cells = true;
1987 *   data_out.set_flags(flags);
1988 *  
1989 *   const MappingQ<3> mapping(mapping_degree);
1990 *  
1991 *   data_out.build_patches(mapping,
1992 *   fe.degree + 2,
1994 *  
1995 *   std::ofstream ofs(fname + ".vtu");
1996 *   data_out.write_vtu(ofs);
1997 *   }
1998 *  
1999 *   void Solver::run()
2000 *   {
2001 *   {
2002 *   TimerOutput::Scope timer_section(timer, "Make mesh");
2003 *   make_mesh();
2004 *   }
2005 *   {
2006 *   TimerOutput::Scope timer_section(timer, "Setup");
2007 *   setup();
2008 *   }
2009 *   {
2010 *   TimerOutput::Scope timer_section(timer, "Assemble");
2011 *   assemble();
2012 *   }
2013 *   {
2014 *   TimerOutput::Scope timer_section(timer, "Solve");
2015 *   solve();
2016 *   }
2017 *   {
2018 *   TimerOutput::Scope timer_section(timer, "Save");
2019 *   save();
2020 *   }
2021 *   {
2022 *   TimerOutput::Scope timer_section(timer, "Clear");
2023 *   clear();
2024 *   }
2025 *   }
2026 *   } // namespace SolverT
2027 *  
2028 * @endcode
2029 *
2030 *
2031 * <a name="step_97-SolverA"></a>
2032 * <h3>Solver - A</h3>
2033 *
2034
2035 *
2036 * This name space contains all the code related to the computation of the
2037 * magnetic vector potential, @f$\vec{A}@f$. The main difference between this
2038 * solver and the solver for the current vector potential, @f$\vec{T}@f$, is in how
2039 * the information on the source is fed to respective solvers. The solver for
2040 * @f$\vec{T}@f$ is fed data sampled from the analytical closed-form expression for
2041 * @f$\vec{J}_f@f$. The solver for @f$\vec{A}@f$ is fed a field function, i.e., a
2042 * numerically computed current vector potential, @f$\vec{T}@f$.
2043 *
2044 * @code
2045 *   namespace SolverA
2046 *   {
2047 * @endcode
2048 *
2049 * This class describes the permeability in the entire problem domain. The
2050 * permeability is given by the definition of the problem, see the
2051 * introduction.
2052 *
2053 * @code
2054 *   class Permeability
2055 *   {
2056 *   public:
2057 *   void value_list(const types::material_id mid,
2058 *   std::vector<double> &values) const
2059 *   {
2060 *   if ((mid == Settings::material_id_free_space) ||
2061 *   (mid == Settings::material_id_free_current))
2062 *   std::fill(values.begin(), values.end(), Settings::mu_0);
2063 *  
2064 *   if (mid == Settings::material_id_core)
2065 *   std::fill(values.begin(), values.end(), Settings::mu_1);
2066 *   }
2067 *   };
2068 *  
2069 * @endcode
2070 *
2071 * This class describes the parameter @f$\gamma@f$ in the Robin boundary
2072 * condition. As soon as it is evaluated on the boundary, the permeability
2073 * equals to that of free space. Therefore, we evaluate the parameter gamma as
2074 * \f[
2075 * \gamma = \dfrac{1}{\mu_0 r}.
2076 * \f]
2077 *
2078 * @code
2079 *   class Gamma
2080 *   {
2081 *   public:
2082 *   void value_list(const std::vector<Point<3>> &r,
2083 *   std::vector<double> &values) const
2084 *   {
2085 *   Assert(r.size() == values.size(),
2086 *   ExcDimensionMismatch(r.size(), values.size()));
2087 *  
2088 *   for (unsigned int i = 0; i < values.size(); i++)
2089 *   values[i] = 1.0 / (Settings::mu_0 * r[i].norm());
2090 *   }
2091 *   };
2092 *  
2093 * @endcode
2094 *
2095 * This class implements the solver that minimizes the functional
2096 * @f$F(\vec{A})@f$. The numerically computed current vector potential, @f$\vec{T}@f$,
2097 * is fed to this solver by means of the input parameters `dof_handler_T` and
2098 * `solution_T`. Moreover, this solver reuses the mesh on which @f$\vec{T}@f$ has
2099 * been computed. The reference to the mesh is passed via the input parameter
2100 * `triangulation_T`.
2101 *
2102 * @code
2103 *   class Solver
2104 *   {
2105 *   public:
2106 *   Solver() = delete;
2107 *   Solver(const unsigned int p, // Degree of the Nedelec finite elements.
2108 *   const unsigned int mapping_degree,
2109 *   const Triangulation<3> &triangulation_T,
2110 *   const DoFHandler<3> &dof_handler_T,
2111 *   const Vector<double> &solution_T,
2112 *   const double eta_squared = 0.0,
2113 *   const std::string &fname = "data");
2114 *  
2115 *   void setup(); // Initializes dofs, vectors, matrices.
2116 *   void assemble(); // Assembles the system of linear equations.
2117 *   void solve(); // Solves the system of linear equations.
2118 *   void save() const; // Saves computed T into a vtu file.
2119 *   void clear() // Clears the memory for the next solver.
2120 *   {
2121 *   system_matrix.clear();
2122 *   system_rhs.reinit(0);
2123 *   }
2124 *   void run(); /* Executes the last five functions in the proper order
2125 *   and measures the execution time for each function. */
2126 *  
2127 *   const DoFHandler<3> &get_dof_handler() const
2128 *   {
2129 *   return dof_handler;
2130 *   }
2131 *   const Vector<double> &get_solution() const
2132 *   {
2133 *   return solution;
2134 *   }
2135 *  
2136 *   private:
2137 *   const Triangulation<3> &triangulation_T;
2138 *   const DoFHandler<3> &dof_handler_T;
2139 *   const Vector<double> &solution_T;
2140 *  
2141 * @endcode
2142 *
2143 * The following data members are typical for all deal.II simulations:
2144 * triangulation, finite elements, dof handlers, etc. The constraints
2145 * are used to enforce the Dirichlet boundary conditions. The names of the
2146 * data members are self-explanatory.
2147 *
2148 * @code
2149 *   const FE_Nedelec<3> fe;
2150 *   DoFHandler<3> dof_handler;
2151 *   Vector<double> solution;
2152 *  
2153 *   SparseMatrix<double> system_matrix;
2154 *   Vector<double> system_rhs;
2155 *  
2156 *   AffineConstraints<double> constraints;
2157 *   SparsityPattern sparsity_pattern;
2158 *  
2159 *   const unsigned int mapping_degree;
2160 *   const double eta_squared;
2161 *   const std::string fname;
2162 *   TimerOutput timer;
2163 *  
2164 * @endcode
2165 *
2166 * This time we have two dof handlers, `dof_handler_T` for @f$\vec{T}@f$ and
2167 * `dof_handler` for @f$\vec{A}@f$. The WorkStream needs to walk through
2168 * the two dof handlers synchronously. For this purpose we will pair two
2169 * active cell iterators (one from `dof_handler_T`, another from
2170 * `dof_handler`). For that we need the `IteratorPair` type.
2171 *
2172 * @code
2173 *   using IteratorTuple =
2174 *   std::tuple<typename DoFHandler<3>::active_cell_iterator,
2176 *  
2177 *   using IteratorPair = SynchronousIterators<IteratorTuple>;
2178 *  
2179 * @endcode
2180 *
2181 * The program utilizes the WorkStream technology. The @ref step_9 "step-9" tutorial
2182 * does a much better job of explaining the workings of WorkStream.
2183 * Reading the "WorkStream paper", see the glossary, is recommended.
2184 * The following structures and functions are related to WorkStream.
2185 *
2186 * @code
2187 *   struct AssemblyScratchData
2188 *   {
2189 *   AssemblyScratchData(const FiniteElement<3> &fe,
2190 *   const DoFHandler<3> &dof_hand_T,
2191 *   const Vector<double> &dofs_T,
2192 *   const unsigned int mapping_degree,
2193 *   const double eta_squared,
2194 *   const BoundaryConditionType boundary_condition_type);
2195 *  
2196 *   AssemblyScratchData(const AssemblyScratchData &scratch_data);
2197 *  
2198 *   const Permeability permeability;
2199 *   const Gamma gamma;
2200 *  
2201 *   MappingQ<3> mapping;
2202 *   FEValues<3> fe_values;
2203 *   FEFaceValues<3> fe_face_values;
2204 *  
2205 *   FEValues<3> fe_values_T;
2206 *  
2207 *   const unsigned int dofs_per_cell;
2208 *   const unsigned int n_q_points;
2209 *   const unsigned int n_q_points_face;
2210 *  
2211 *   std::vector<double> permeability_list;
2212 *   std::vector<double> gamma_list;
2213 *   std::vector<Tensor<1, 3>> T_values;
2214 *  
2215 *   const FEValuesExtractors::Vector ve;
2216 *  
2217 *   const DoFHandler<3> &dof_hand_T;
2218 *   const Vector<double> &dofs_T;
2219 *  
2220 *   const double eta_squared;
2221 *   const BoundaryConditionType boundary_condition_type;
2222 *   };
2223 *  
2224 *   struct AssemblyCopyData
2225 *   {
2227 *   Vector<double> cell_rhs;
2228 *   std::vector<types::global_dof_index> local_dof_indices;
2229 *   };
2230 *  
2231 *   void system_matrix_local(const IteratorPair &IP,
2232 *   AssemblyScratchData &scratch_data,
2233 *   AssemblyCopyData &copy_data);
2234 *  
2235 *   void copy_local_to_global(const AssemblyCopyData &copy_data);
2236 *   };
2237 *  
2238 *   Solver::Solver(const unsigned int p,
2239 *   const unsigned int mapping_degree,
2240 *   const Triangulation<3> &triangulation_T,
2241 *   const DoFHandler<3> &dof_handler_T,
2242 *   const Vector<double> &solution_T,
2243 *   const double eta_squared,
2244 *   const std::string &fname)
2245 *   : triangulation_T(triangulation_T)
2246 *   , dof_handler_T(dof_handler_T)
2247 *   , solution_T(solution_T)
2248 *   , fe(p)
2249 *   , mapping_degree(mapping_degree)
2250 *   , eta_squared(eta_squared)
2251 *   , fname(fname)
2252 *   , timer(std::cout,
2253 *   (Settings::print_time_tables) ? TimerOutput::summary :
2254 *   TimerOutput::never,
2256 *   {}
2257 *  
2258 * @endcode
2259 *
2260 * This function initializes the dofs, applies the Dirichlet boundary
2261 * condition, and initializes the vectors and matrices.
2262 *
2263 * @code
2264 *   void Solver::setup()
2265 *   {
2266 *   dof_handler.reinit(triangulation_T);
2267 *   dof_handler.distribute_dofs(fe);
2268 *  
2269 * @endcode
2270 *
2271 * The following segment of the code applies the homogeneous Dirichlet
2272 * boundary condition. As discussed in the introduction, the Dirichlet
2273 * boundary condition is an essential condition and must be enforced
2274 * by constraining the system matrix. This segment of code does the
2275 * constraining.
2276 *
2277 * @code
2278 *   constraints.clear();
2279 *  
2280 *   DoFTools::make_hanging_node_constraints(dof_handler, constraints);
2281 *  
2282 *   if (Settings::boundary_condition_type_A == Dirichlet)
2284 *   dof_handler,
2285 *   0, // The first vector component.
2287 *   Settings::boundary_id_infinity,
2288 *   constraints,
2289 *   MappingQ<3>(mapping_degree));
2290 *  
2291 *   constraints.close();
2292 *  
2293 * @endcode
2294 *
2295 * The rest of the function arranges the dofs in a sparsity pattern and
2296 * initializes the system matrix and the system vectors.
2297 *
2298 * @code
2299 *   DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
2300 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints, false);
2301 *  
2302 *   sparsity_pattern.copy_from(dsp);
2303 *   system_matrix.reinit(sparsity_pattern);
2304 *   solution.reinit(dof_handler.n_dofs());
2305 *   system_rhs.reinit(dof_handler.n_dofs());
2306 *   }
2307 *  
2308 * @endcode
2309 *
2310 * Formally, this function assembles the system of linear equations. In
2311 * reality, however, it just spells all the magic words to get the WorkStream
2312 * going. The interesting part, i.e., the actual assembling of the system
2313 * matrix and the right-hand side happens below in the
2314 * Solver::system_matrix_local function. Note that this time the first two
2315 * input parameters to `WorkStream::run` are pairs of iterators, not
2316 * iterators themselves as per usual. Note also the order in which we package
2317 * the iterators: first the iterator of `dof_handler`, then the iterator of
2318 * the `dof_handler_T`. We will extract them in the same order.
2319 *
2320 * @code
2321 *   void Solver::assemble()
2322 *   {
2323 *   WorkStream::run(IteratorPair(IteratorTuple(dof_handler.begin_active(),
2324 *   dof_handler_T.begin_active())),
2325 *   IteratorPair(
2326 *   IteratorTuple(dof_handler.end(), dof_handler_T.end())),
2327 *   *this,
2328 *   &Solver::system_matrix_local,
2329 *   &Solver::copy_local_to_global,
2330 *   AssemblyScratchData(fe,
2331 *   dof_handler_T,
2332 *   solution_T,
2333 *   mapping_degree,
2334 *   eta_squared,
2335 *   Settings::boundary_condition_type_A),
2336 *   AssemblyCopyData());
2337 *   }
2338 *  
2339 * @endcode
2340 *
2341 * The following two constructors initialize scratch data from the input
2342 * parameters and from another object of the same type, i.e., a copy
2343 * constructor.
2344 *
2345 * @code
2346 *   Solver::AssemblyScratchData::AssemblyScratchData(
2347 *   const FiniteElement<3> &fe,
2348 *   const DoFHandler<3> &dof_hand_T,
2349 *   const Vector<double> &dofs_T,
2350 *   const unsigned int mapping_degree,
2351 *   const double eta_squared,
2352 *   const BoundaryConditionType boundary_condition_type)
2353 *   : permeability()
2354 *   , gamma()
2355 *   , mapping(mapping_degree)
2356 *   , fe_values(mapping,
2357 *   fe,
2358 *   QGauss<3>(fe.degree + 2),
2360 *   , fe_face_values(mapping,
2361 *   fe,
2362 *   QGauss<2>(fe.degree + 2),
2365 *   , fe_values_T(mapping,
2366 *   dof_hand_T.get_fe(),
2367 *   QGauss<3>(fe.degree + 2),
2369 *   , dofs_per_cell(fe_values.dofs_per_cell)
2370 *   , n_q_points(fe_values.get_quadrature().size())
2371 *   , n_q_points_face(fe_face_values.get_quadrature().size())
2372 *   , permeability_list(n_q_points)
2373 *   , gamma_list(n_q_points_face)
2374 *   , T_values(n_q_points)
2375 *   , ve(0)
2376 *   , dof_hand_T(dof_hand_T)
2377 *   , dofs_T(dofs_T)
2378 *   , eta_squared(eta_squared)
2379 *   , boundary_condition_type(boundary_condition_type)
2380 *   {}
2381 *  
2382 *   Solver::AssemblyScratchData::AssemblyScratchData(
2383 *   const AssemblyScratchData &scratch_data)
2384 *   : permeability()
2385 *   , gamma()
2386 *   , mapping(scratch_data.mapping.get_degree())
2387 *   , fe_values(mapping,
2388 *   scratch_data.fe_values.get_fe(),
2389 *   scratch_data.fe_values.get_quadrature(),
2391 *   , fe_face_values(mapping,
2392 *   scratch_data.fe_face_values.get_fe(),
2393 *   scratch_data.fe_face_values.get_quadrature(),
2396 *   , fe_values_T(mapping,
2397 *   scratch_data.fe_values_T.get_fe(),
2398 *   scratch_data.fe_values_T.get_quadrature(),
2400 *   , dofs_per_cell(fe_values.dofs_per_cell)
2401 *   , n_q_points(fe_values.get_quadrature().size())
2402 *   , n_q_points_face(fe_face_values.get_quadrature().size())
2403 *   , permeability_list(n_q_points)
2404 *   , gamma_list(n_q_points_face)
2405 *   , T_values(n_q_points)
2406 *   , ve(0)
2407 *   , dof_hand_T(scratch_data.dof_hand_T)
2408 *   , dofs_T(scratch_data.dofs_T)
2409 *   , eta_squared(scratch_data.eta_squared)
2410 *   , boundary_condition_type(scratch_data.boundary_condition_type)
2411 *   {}
2412 *  
2413 * @endcode
2414 *
2415 * This function assembles a fraction of the system matrix and the system
2416 * right-hand side related to a single cell. These fractions are
2417 * `copy_data.cell_matrix` and `copy_data.cell_rhs`. They are copied into
2418 * to the system matrix, @f$A_{ij}@f$, and the right-hand side, @f$b_i@f$, by the
2419 * function `Solver::copy_local_to_global()`.
2420 *
2421 * @code
2422 *   void Solver::system_matrix_local(const IteratorPair &IP,
2423 *   AssemblyScratchData &scratch_data,
2424 *   AssemblyCopyData &copy_data)
2425 *   {
2426 * @endcode
2427 *
2428 * First we reinitialize the matrices and vectors related to the current
2429 * cell.
2430 *
2431 * @code
2432 *   copy_data.cell_matrix.reinit(scratch_data.dofs_per_cell,
2433 *   scratch_data.dofs_per_cell);
2434 *  
2435 *   copy_data.cell_rhs.reinit(scratch_data.dofs_per_cell);
2436 *  
2437 *   copy_data.local_dof_indices.resize(scratch_data.dofs_per_cell);
2438 *  
2439 * @endcode
2440 *
2441 * Second, we extract the cells from the pair. We extract them in the
2442 * correct order, see above.
2443 *
2444 * @code
2445 *   auto cell = std::get<0>(*IP);
2446 *   auto cell_T = std::get<1>(*IP);
2447 *  
2448 * @endcode
2449 *
2450 * Third, we compute the ordered FE values, the permeability, and the values
2451 * of the current vector potential, @f$\vec{T}@f$, on the cell.
2452 *
2453 * @code
2454 *   scratch_data.fe_values.reinit(cell);
2455 *   scratch_data.fe_values_T.reinit(cell_T);
2456 *  
2457 *   scratch_data.permeability.value_list(cell->material_id(),
2458 *   scratch_data.permeability_list);
2459 *  
2460 *   scratch_data.fe_values_T[scratch_data.ve].get_function_values(
2461 *   scratch_data.dofs_T, scratch_data.T_values);
2462 *  
2463 * @endcode
2464 *
2465 * Fourth, we compute the components of the cell matrix and cell right-hand
2466 * side. The labels of the integrals are the same as in the introduction to
2467 * this tutorial.
2468 *
2469 * @code
2470 *   for (unsigned int q_index = 0; q_index < scratch_data.n_q_points; ++q_index)
2471 *   {
2472 *   for (unsigned int i = 0; i < scratch_data.dofs_per_cell; ++i)
2473 *   {
2474 *   for (unsigned int j = 0; j < scratch_data.dofs_per_cell; ++j)
2475 *   {
2476 *   copy_data.cell_matrix(i, j) += // Integral I_a1+I_a3.
2477 *   (1.0 / scratch_data.permeability_list[q_index]) * // 1 / mu
2478 *   (scratch_data.fe_values[scratch_data.ve].curl(
2479 *   i, q_index) * // curl N_i
2480 *   scratch_data.fe_values[scratch_data.ve].curl(
2481 *   j, q_index) // curl N_j
2482 *   + scratch_data.eta_squared * // eta^2
2483 *   scratch_data.fe_values[scratch_data.ve].value(
2484 *   i, q_index) * // N_i
2485 *   scratch_data.fe_values[scratch_data.ve].value(
2486 *   j, q_index) // N_j
2487 *   ) *
2488 *   scratch_data.fe_values.JxW(q_index); // dV
2489 *   }
2490 *   copy_data.cell_rhs(i) += // Integral I_b3-1.
2491 *   (scratch_data.T_values[q_index] *
2492 *   scratch_data.fe_values[scratch_data.ve].curl(i, q_index)) *
2493 *   scratch_data.fe_values.JxW(q_index); // T.(curl N_i)dV
2494 *   }
2495 *   }
2496 *  
2497 * @endcode
2498 *
2499 * If the Robin boundary condition (first-order ABC) is ordered,
2500 * we compute an extra integral over the boundary.
2501 *
2502 * @code
2503 *   if (scratch_data.boundary_condition_type == BoundaryConditionType::Robin)
2504 *   {
2505 *   for (unsigned int f = 0; f < cell->n_faces(); ++f)
2506 *   {
2507 *   if (cell->face(f)->at_boundary())
2508 *   {
2509 *   scratch_data.fe_face_values.reinit(cell, f);
2510 *  
2511 *   for (unsigned int q_index_face = 0;
2512 *   q_index_face < scratch_data.n_q_points_face;
2513 *   ++q_index_face)
2514 *   {
2515 *   for (unsigned int i = 0; i < scratch_data.dofs_per_cell;
2516 *   ++i)
2517 *   {
2518 *   scratch_data.gamma.value_list(
2519 *   scratch_data.fe_face_values.get_quadrature_points(),
2520 *   scratch_data.gamma_list);
2521 *  
2522 *   for (unsigned int j = 0; j < scratch_data.dofs_per_cell;
2523 *   ++j)
2524 *   {
2525 *   copy_data.cell_matrix(i, j) += // Integral I_a2.
2526 *   scratch_data.gamma_list[q_index_face] * // gamma
2527 *   (cross_product_3d(
2528 *   scratch_data.fe_face_values.normal_vector(
2529 *   q_index_face),
2530 *   scratch_data.fe_face_values[scratch_data.ve]
2531 *   .value(i, q_index_face)) *
2532 *   cross_product_3d(
2533 *   scratch_data.fe_face_values.normal_vector(
2534 *   q_index_face),
2535 *   scratch_data.fe_face_values[scratch_data.ve]
2536 *   .value(j,
2537 *   q_index_face))) // (n x N_i).(n x N_j)
2538 *   * scratch_data.fe_face_values.JxW(
2539 *   q_index_face); // dS
2540 *   } // for j
2541 *   } // for i
2542 *   } // for q_index_face
2543 *   } // if (cell->face(f)->at_boundary())
2544 *   } // for f
2545 *   } /* if (scratch_data.boundary_condition_type ==
2546 *   BoundaryConditionType::Robin) */
2547 *  
2548 * @endcode
2549 *
2550 * Finally, we query the dof indices on the current cell and store them
2551 * in the copy data structure, so we know to which locations of the system
2552 * matrix and right-hand side the components of the cell matrix and
2553 * cell right-hand side must be copied.
2554 *
2555 * @code
2556 *   cell->get_dof_indices(copy_data.local_dof_indices);
2557 *   }
2558 *  
2559 * @endcode
2560 *
2561 * This function copies the components of a cell matrix and a cell right-hand
2562 * side into the system matrix, @f$A_{i,j}@f$, and the system right-hand side,
2563 * @f$b_i@f$.
2564 *
2565 * @code
2566 *   void Solver::copy_local_to_global(const AssemblyCopyData &copy_data)
2567 *   {
2568 *   constraints.distribute_local_to_global(copy_data.cell_matrix,
2569 *   copy_data.cell_rhs,
2570 *   copy_data.local_dof_indices,
2571 *   system_matrix,
2572 *   system_rhs);
2573 *   }
2574 *  
2575 * @endcode
2576 *
2577 * This function solves the system of linear equations.
2578 * If `Settings::log_cg_convergence == true`, the convergence data is saved
2579 * into a file. In theory, a CG solver can solve an @f$m \times m@f$ system of
2580 * linear equations in at most @f$m@f$ steps. In practice, it can take more steps
2581 * to converge. The convergence of the algorithm depends on the spectral
2582 * properties of the system matrix. The best case is if the eigenvalues form a
2583 * compact cluster away from zero. In our case, however, the eigenvalues are
2584 * spread in between zero and the maximal eigenvalue. Consequently, we expect
2585 * a poor convergence and increase the maximal number of iteration steps by a
2586 * factor of 10, i.e., `10*system_rhs.size()`. The stopping condition is \f[
2587 * |\boldsymbol{b} - \boldsymbol{A}\boldsymbol{c}|
2588 * < 10^{-6} |\boldsymbol{b}|.
2589 * \f]
2590 * As soon as we use constraints, we must not forget to distribute them.
2591 *
2592 * @code
2593 *   void Solver::solve()
2594 *   {
2595 *   SolverControl control(10 * system_rhs.size(),
2596 *   1.0e-6 * system_rhs.l2_norm(),
2597 *   false,
2598 *   false);
2599 *  
2600 *   if (Settings::log_cg_convergence)
2601 *   control.enable_history_data();
2602 *  
2604 *   SolverCG<Vector<double>> cg(control, memory);
2605 *  
2606 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
2607 *   preconditioner.initialize(system_matrix, 1.2);
2608 *  
2609 *   cg.solve(system_matrix, solution, system_rhs, preconditioner);
2610 *  
2611 *   constraints.distribute(solution);
2612 *  
2613 *   if (Settings::log_cg_convergence)
2614 *   {
2615 *   const std::vector<double> history_data = control.get_history_data();
2616 *  
2617 *   std::ofstream ofs(fname + "_cg_convergence.csv");
2618 *  
2619 *   for (unsigned int i = 1; i < history_data.size(); i++)
2620 *   ofs << i << ", " << history_data[i] << "\n";
2621 *   }
2622 *   }
2623 *  
2624 * @endcode
2625 *
2626 * This function saves the computed magnetic vector potential into a vtu file.
2627 *
2628 * @code
2629 *   void Solver::save() const
2630 *   {
2631 *   std::vector<std::string> solution_names(3, "VectorField");
2632 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
2633 *   interpretation(3,
2635 *  
2636 *   DataOut<3> data_out;
2637 *  
2638 *   data_out.add_data_vector(dof_handler,
2639 *   solution,
2640 *   solution_names,
2641 *   interpretation);
2642 *  
2643 *   DataOutBase::VtkFlags flags;
2644 *   flags.write_higher_order_cells = true;
2645 *   data_out.set_flags(flags);
2646 *  
2647 *   const MappingQ<3> mapping(mapping_degree);
2648 *  
2649 *   data_out.build_patches(mapping,
2650 *   fe.degree + 2,
2652 *  
2653 *   std::ofstream ofs(fname + ".vtu");
2654 *   data_out.write_vtu(ofs);
2655 *   }
2656 *  
2657 *   void Solver::run()
2658 *   {
2659 *   {
2660 *   TimerOutput::Scope timer_section(timer, "Setup");
2661 *   setup();
2662 *   }
2663 *   {
2664 *   TimerOutput::Scope timer_section(timer, "Assemble");
2665 *   assemble();
2666 *   }
2667 *   {
2668 *   TimerOutput::Scope timer_section(timer, "Solve");
2669 *   solve();
2670 *   }
2671 *   {
2672 *   TimerOutput::Scope timer_section(timer, "Save");
2673 *   save();
2674 *   }
2675 *   {
2676 *   TimerOutput::Scope timer_section(timer, "Clear");
2677 *   clear();
2678 *   }
2679 *   }
2680 *   } // namespace SolverA
2681 *  
2682 *  
2683 * @endcode
2684 *
2685 *
2686 * <a name="step_97-ProjectorfromHcurltoHdiv"></a>
2687 * <h3>Projector from H(curl) to H(div)</h3>
2688 * This name space contains all the code related to the conversion of the
2689 * magnetic vector potential, @f$\vec{A}@f$, into magnetic field, @f$\vec{B}@f$.
2690 * The magnetic vector potential is modeled by the FE_Nedelec finite elements,
2691 * while the magnetic field is modeled by the FE_RaviartThomas finite elements.
2692 * This code is also used for converting the current vector potential,
2693 * @f$\vec{T}@f$ into the free-current density, @f$\vec{J}_f@f$.
2694 *
2695 * @code
2696 *   namespace ProjectorHcurlToHdiv
2697 *   {
2698 *  
2699 * @endcode
2700 *
2701 * This class implements the solver that minimizes the functional @f$F(\vec{B})@f$
2702 * or @f$F(\vec{J}_f)@f$, see the introduction. The input vector field,
2703 * @f$\vec{A}@f$ or @f$\vec{T}@f$, is fed to the solver by means of the input
2704 * parameters `dof_handler_Hcurl` and `solution_Hcurl`. Moreover, this solver
2705 * reuses the mesh on which the input vector field has been computed. The
2706 * reference to the mesh is passed via the input parameter
2707 * `triangulation_Hcurl`. There are no constraints this time around as we are
2708 * not going to apply the Dirichlet boundary condition.
2709 *
2710 * @code
2711 *   class Solver
2712 *   {
2713 *   public:
2714 *   Solver() = delete;
2715 *   Solver(const unsigned int p, // Degree of the Raviart-Thomas finite elements
2716 *   const unsigned int mapping_degree,
2717 *   const Triangulation<3> &triangulation_Hcurl,
2718 *   const DoFHandler<3> &dof_handler_Hcurl,
2719 *   const Vector<double> &solution_Hcurl,
2720 *   const std::string &fname = "data",
2721 *   const Function<3> *exact_solution = nullptr);
2722 *  
2723 *   double get_L2_norm()
2724 *   {
2725 *   return L2_norm;
2726 *   };
2727 *  
2728 *   unsigned int get_n_cells() const
2729 *   {
2730 *   return triangulation_Hcurl.n_active_cells();
2731 *   }
2732 *  
2733 *   types::global_dof_index get_n_dofs() const
2734 *   {
2735 *   return dof_handler_Hdiv.n_dofs();
2736 *   }
2737 *  
2738 *   void setup(); // Initializes dofs, vectors, matrices.
2739 *   void assemble(); // Assembles the system of linear equations.
2740 *   void solve(); // Solves the system of linear equations.
2741 *   void save() const; // Saves computed T into a vtu file.
2742 *   void compute_error_norms(); // Computes L^2 error norm.
2743 *   void project_exact_solution_fcn(); // Projects exact solution.
2744 *   void clear()
2745 *   {
2746 *   system_matrix.clear();
2747 *   system_rhs.reinit(0);
2748 *   }
2749 *   void run(); /* Executes the last seven functions in the proper order
2750 *   and measures the execution time for each function. */
2751 *  
2752 *   private:
2753 *   const Triangulation<3> &triangulation_Hcurl;
2754 *   const DoFHandler<3> &dof_handler_Hcurl;
2755 *   const Vector<double> &solution_Hcurl;
2756 *  
2757 * @endcode
2758 *
2759 * The following data members are typical for all deal.II simulations:
2760 * triangulation, finite elements, dof handlers, etc. The constraints
2761 * are used to enforce the Dirichlet boundary conditions. The names of the
2762 * data members are self-explanatory.
2763 *
2764 * @code
2765 *   const FE_RaviartThomas<3> fe_Hdiv;
2766 *   DoFHandler<3> dof_handler_Hdiv;
2767 *  
2768 *   SparsityPattern sparsity_pattern;
2769 *   SparseMatrix<double> system_matrix;
2770 *  
2771 *   Vector<double> solution_Hdiv;
2772 *   Vector<double> system_rhs;
2773 *  
2774 *   Vector<double> projected_exact_solution;
2775 *  
2776 *   AffineConstraints<double> constraints;
2777 *  
2778 *   const Function<3> *exact_solution;
2779 *  
2780 *   const unsigned int mapping_degree;
2781 *  
2782 *   Vector<double> L2_per_cell;
2783 *   double L2_norm;
2784 *  
2785 *   const std::string fname;
2786 *   TimerOutput timer;
2787 *  
2788 * @endcode
2789 *
2790 * This time we have two dof handlers, `dof_handler_Hcurl` for the input
2791 * vector field and `dof_handler_Hdiv` for the output vector field. The
2792 * WorkStream needs to walk through the two dof handlers synchronously.
2793 * For this purpose we will pair two active cells iterators (one from
2794 * `dof_handler_Hcurl`, another from `dof_handler_Hdiv`) to be walked
2795 * through synchronously. For that we need the `IteratorPair` type.
2796 *
2797 * @code
2798 *   using IteratorTuple =
2799 *   std::tuple<typename DoFHandler<3>::active_cell_iterator,
2801 *  
2802 *   using IteratorPair = SynchronousIterators<IteratorTuple>;
2803 *  
2804 * @endcode
2805 *
2806 * The program utilizes the WorkStream technology. The @ref step_9 "step-9" tutorial
2807 * does a much better job of explaining the workings of WarkStream.
2808 * Reading the "WorkStream paper", see the glossary, is recommended.
2809 * The following structures and functions are related to WorkStream.
2810 *
2811 * @code
2812 *   struct AssemblyScratchData
2813 *   {
2814 *   AssemblyScratchData(const FiniteElement<3> &fe,
2815 *   const DoFHandler<3> &dof_handr_Hcurl,
2816 *   const Vector<double> &dofs_Hcurl,
2817 *   const unsigned int mapping_degree);
2818 *  
2819 *   AssemblyScratchData(const AssemblyScratchData &scratch_data);
2820 *  
2821 *   MappingQ<3> mapping;
2822 *   FEValues<3> fe_values_Hdiv;
2823 *   FEValues<3> fe_values_Hcurl;
2824 *  
2825 *   const unsigned int dofs_per_cell;
2826 *   const unsigned int n_q_points;
2827 *  
2828 *   std::vector<Tensor<1, 3>> curl_vec_in_Hcurl;
2829 *  
2830 *   const FEValuesExtractors::Vector ve;
2831 *  
2832 *   const DoFHandler<3> &dof_hand_Hcurl;
2833 *   const Vector<double> &dofs_Hcurl;
2834 *   };
2835 *  
2836 *   struct AssemblyCopyData
2837 *   {
2839 *   Vector<double> cell_rhs;
2840 *   std::vector<types::global_dof_index> local_dof_indices;
2841 *   };
2842 *  
2843 *   void system_matrix_local(const IteratorPair &IP,
2844 *   AssemblyScratchData &scratch_data,
2845 *   AssemblyCopyData &copy_data);
2846 *  
2847 *   void copy_local_to_global(const AssemblyCopyData &copy_data);
2848 *   };
2849 *  
2850 *   Solver::Solver(const unsigned int p,
2851 *   const unsigned int mapping_degree,
2852 *   const Triangulation<3> &triangulation_Hcurl,
2853 *   const DoFHandler<3> &dof_handler_Hcurl,
2854 *   const Vector<double> &solution_Hcurl,
2855 *   const std::string &fname,
2856 *   const Function<3> *exact_solution)
2857 *   : triangulation_Hcurl(triangulation_Hcurl)
2858 *   , dof_handler_Hcurl(dof_handler_Hcurl)
2859 *   , solution_Hcurl(solution_Hcurl)
2860 *   , fe_Hdiv(p)
2861 *   , exact_solution(exact_solution)
2862 *   , mapping_degree(mapping_degree)
2863 *   , fname(fname)
2864 *   , timer(std::cout,
2865 *   (Settings::print_time_tables) ? TimerOutput::summary :
2866 *   TimerOutput::never,
2868 *   {
2869 *   Assert(exact_solution != nullptr,
2870 *   ExcMessage("The exact solution is missing."));
2871 *   }
2872 *  
2873 * @endcode
2874 *
2875 * This function initializes the dofs, vectors and matrices. This time there
2876 * are no constraints as we do not apply Dirichlet boundary condition.
2877 *
2878 * @code
2879 *   void Solver::setup()
2880 *   {
2881 *   constraints.close();
2882 *  
2883 *   dof_handler_Hdiv.reinit(triangulation_Hcurl);
2884 *   dof_handler_Hdiv.distribute_dofs(fe_Hdiv);
2885 *  
2886 *   DynamicSparsityPattern dsp(dof_handler_Hdiv.n_dofs(),
2887 *   dof_handler_Hdiv.n_dofs());
2888 *   DoFTools::make_sparsity_pattern(dof_handler_Hdiv, dsp, constraints, false);
2889 *  
2890 *   sparsity_pattern.copy_from(dsp);
2891 *   system_matrix.reinit(sparsity_pattern);
2892 *   solution_Hdiv.reinit(dof_handler_Hdiv.n_dofs());
2893 *   system_rhs.reinit(dof_handler_Hdiv.n_dofs());
2894 *  
2895 *   if (Settings::project_exact_solution && exact_solution)
2896 *   projected_exact_solution.reinit(dof_handler_Hdiv.n_dofs());
2897 *  
2898 *   if (exact_solution)
2899 *   L2_per_cell.reinit(triangulation_Hcurl.n_active_cells());
2900 *   }
2901 *  
2902 * @endcode
2903 *
2904 * Formally, this function assembles the system of linear equations. In
2905 * reality, however, it just spells all the magic words to get the WorkStream
2906 * going. The interesting part, i.e., the actual assembling of the system
2907 * matrix and the right-hand side happens below in the
2908 * Solver::system_matrix_local function.
2909 *
2910 * @code
2911 *   void Solver::assemble()
2912 *   {
2913 *   WorkStream::run(IteratorPair(
2914 *   IteratorTuple(dof_handler_Hdiv.begin_active(),
2915 *   dof_handler_Hcurl.begin_active())),
2916 *   IteratorPair(IteratorTuple(dof_handler_Hdiv.end(),
2917 *   dof_handler_Hcurl.end())),
2918 *   *this,
2919 *   &Solver::system_matrix_local,
2920 *   &Solver::copy_local_to_global,
2921 *   AssemblyScratchData(fe_Hdiv,
2922 *   dof_handler_Hcurl,
2923 *   solution_Hcurl,
2924 *   mapping_degree),
2925 *   AssemblyCopyData());
2926 *   }
2927 *  
2928 * @endcode
2929 *
2930 * The following two constructors initialize scratch data from the input
2931 * parameters and from another object of the same type, i.e., a copy
2932 * constructor.
2933 *
2934 * @code
2935 *   Solver::AssemblyScratchData::AssemblyScratchData(
2936 *   const FiniteElement<3> &fe,
2937 *   const DoFHandler<3> &dof_hand_Hcurl,
2938 *   const Vector<double> &dofs_Hcurl,
2939 *   const unsigned int mapping_degree)
2940 *   : mapping(mapping_degree)
2941 *   , fe_values_Hdiv(mapping,
2942 *   fe,
2943 *   QGauss<3>(fe.degree + 2),
2945 *   , fe_values_Hcurl(mapping,
2946 *   dof_hand_Hcurl.get_fe(),
2947 *   QGauss<3>(fe.degree + 2),
2948 *   update_gradients)
2949 *   , dofs_per_cell(fe_values_Hdiv.dofs_per_cell)
2950 *   , n_q_points(fe_values_Hdiv.get_quadrature().size())
2951 *   , curl_vec_in_Hcurl(n_q_points)
2952 *   , ve(0)
2953 *   , dof_hand_Hcurl(dof_hand_Hcurl)
2954 *   , dofs_Hcurl(dofs_Hcurl)
2955 *   {}
2956 *  
2957 *   Solver::AssemblyScratchData::AssemblyScratchData(
2958 *   const AssemblyScratchData &scratch_data)
2959 *   : mapping(scratch_data.mapping.get_degree())
2960 *   , fe_values_Hdiv(mapping,
2961 *   scratch_data.fe_values_Hdiv.get_fe(),
2962 *   scratch_data.fe_values_Hdiv.get_quadrature(),
2964 *   , fe_values_Hcurl(mapping,
2965 *   scratch_data.fe_values_Hcurl.get_fe(),
2966 *   scratch_data.fe_values_Hcurl.get_quadrature(),
2967 *   update_gradients)
2968 *   , dofs_per_cell(fe_values_Hdiv.dofs_per_cell)
2969 *   , n_q_points(fe_values_Hdiv.get_quadrature().size())
2970 *   , curl_vec_in_Hcurl(scratch_data.n_q_points)
2971 *   , ve(0)
2972 *   , dof_hand_Hcurl(scratch_data.dof_hand_Hcurl)
2973 *   , dofs_Hcurl(scratch_data.dofs_Hcurl)
2974 *   {}
2975 *  
2976 * @endcode
2977 *
2978 * This function assembles a fraction of the system matrix and the system
2979 * right-hand side related to a single cell. These fractions are
2980 * `copy_data.cell_matrix` and `copy_data.cell_rhs`. They are copied into
2981 * to the system matrix, @f$A_{ij}@f$, and the right-hand side, @f$b_i@f$, by the
2982 * function `Solver::copy_local_to_global()`.
2983 *
2984 * @code
2985 *   void Solver::system_matrix_local(const IteratorPair &IP,
2986 *   AssemblyScratchData &scratch_data,
2987 *   AssemblyCopyData &copy_data)
2988 *   {
2989 * @endcode
2990 *
2991 * First we reinitialize the matrices and vectors related to the current
2992 * cell, update the FE values, and compute the curl of the input vector
2993 * field.
2994 *
2995 * @code
2996 *   copy_data.cell_matrix.reinit(scratch_data.dofs_per_cell,
2997 *   scratch_data.dofs_per_cell);
2998 *  
2999 *   copy_data.cell_rhs.reinit(scratch_data.dofs_per_cell);
3000 *  
3001 *   copy_data.local_dof_indices.resize(scratch_data.dofs_per_cell);
3002 *  
3003 *   scratch_data.fe_values_Hdiv.reinit(std::get<0>(*IP));
3004 *   scratch_data.fe_values_Hcurl.reinit(std::get<1>(*IP));
3005 *  
3006 * @endcode
3007 *
3008 * The variable `curl_vec_in_Hcurl` denotes the curl of the input vector
3009 * field, @f$\vec{\nabla} \times \vec{T}@f$ or @f$\vec{\nabla} \times \vec{A}@f$,
3010 * depending on the context.
3011 *
3012 * @code
3013 *   scratch_data.fe_values_Hcurl[scratch_data.ve].get_function_curls(
3014 *   scratch_data.dofs_Hcurl, scratch_data.curl_vec_in_Hcurl);
3015 *  
3016 * @endcode
3017 *
3018 * Second, we compute the components of the cell matrix and cell right-hand
3019 * side. The labels of the integrals are the same as in the introduction to
3020 * this tutorial.
3021 *
3022 * @code
3023 *   for (unsigned int q_index = 0; q_index < scratch_data.n_q_points; ++q_index)
3024 *   {
3025 *   for (unsigned int i = 0; i < scratch_data.dofs_per_cell; ++i)
3026 *   {
3027 *   for (unsigned int j = 0; j < scratch_data.dofs_per_cell; ++j)
3028 *   {
3029 *   copy_data.cell_matrix(i, j) += // Integral I_a
3030 *   scratch_data.fe_values_Hdiv[scratch_data.ve].value(i,
3031 *   q_index) *
3032 *   scratch_data.fe_values_Hdiv[scratch_data.ve].value(j,
3033 *   q_index) *
3034 *   scratch_data.fe_values_Hdiv.JxW(q_index);
3035 *   }
3036 *  
3037 *   copy_data.cell_rhs(i) += // Integral I_b
3038 *   scratch_data.curl_vec_in_Hcurl[q_index] *
3039 *   scratch_data.fe_values_Hdiv[scratch_data.ve].value(i, q_index) *
3040 *   scratch_data.fe_values_Hdiv.JxW(q_index);
3041 *   }
3042 *   }
3043 *  
3044 * @endcode
3045 *
3046 * Finally, we query the dof indices on the current cell and store them
3047 * in the copy data structure, so we know to which locations of the system
3048 * matrix and right-hand side the components of the cell matrix and
3049 * cell right-hand side must be copied.
3050 *
3051 * @code
3052 *   std::get<0>(*IP)->get_dof_indices(copy_data.local_dof_indices);
3053 *   }
3054 *  
3055 * @endcode
3056 *
3057 * This function copies the components of a cell matrix and a cell right-hand
3058 * side into the system matrix, @f$A_{i,j}@f$, and the system right-hand side,
3059 * @f$b_i@f$.
3060 *
3061 * @code
3062 *   void Solver::copy_local_to_global(const AssemblyCopyData &copy_data)
3063 *   {
3064 *   constraints.distribute_local_to_global(copy_data.cell_matrix,
3065 *   copy_data.cell_rhs,
3066 *   copy_data.local_dof_indices,
3067 *   system_matrix,
3068 *   system_rhs);
3069 *   }
3070 *  
3071 * @endcode
3072 *
3073 * The following two functions compute the error norms and project the exact
3074 * solution.
3075 *
3076 * @code
3077 *   void Solver::compute_error_norms()
3078 *   {
3079 *   const Weight weight;
3080 *   const Function<3, double> *mask = &weight;
3081 *  
3083 *   dof_handler_Hdiv,
3084 *   solution_Hdiv,
3085 *   *exact_solution,
3086 *   L2_per_cell,
3087 *   QGauss<3>(fe_Hdiv.degree + 4),
3089 *   mask);
3090 *  
3091 *   L2_norm = VectorTools::compute_global_error(triangulation_Hcurl,
3092 *   L2_per_cell,
3094 *   }
3095 *  
3096 *   void Solver::project_exact_solution_fcn()
3097 *   {
3098 *   AffineConstraints<double> constraints_empty;
3099 *   constraints_empty.close();
3100 *  
3101 *   VectorTools::project(MappingQ<3>(mapping_degree),
3102 *   dof_handler_Hdiv,
3103 *   constraints_empty,
3104 *   QGauss<3>(fe_Hdiv.degree + 2),
3105 *   *exact_solution,
3106 *   projected_exact_solution);
3107 *   }
3108 *  
3109 * @endcode
3110 *
3111 * This function solves the system of linear equations. This time we are
3112 * dealing with a mass matrix. It has good spectral properties. Consequently,
3113 * we do not use the factor of 10 as in preceding two solvers.
3114 * The stopping condition is
3115 * \f[
3116 * |\boldsymbol{b} - \boldsymbol{A}\boldsymbol{c}|
3117 * < 10^{-6} |\boldsymbol{b}|.
3118 * \f]
3119 *
3120 * @code
3121 *   void Solver::solve()
3122 *   {
3123 *   SolverControl control(system_rhs.size(),
3124 *   1.0e-6 * system_rhs.l2_norm(),
3125 *   false,
3126 *   false);
3127 *  
3128 *   if (Settings::log_cg_convergence)
3129 *   control.enable_history_data();
3130 *  
3132 *   SolverCG<Vector<double>> cg(control, memory);
3133 *  
3134 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
3135 *   preconditioner.initialize(system_matrix, 1.2);
3136 *  
3137 *   cg.solve(system_matrix, solution_Hdiv, system_rhs, preconditioner);
3138 *  
3139 *   if (Settings::log_cg_convergence)
3140 *   {
3141 *   const std::vector<double> history_data = control.get_history_data();
3142 *  
3143 *   std::ofstream ofs(fname + "_cg_convergence.csv");
3144 *  
3145 *   for (unsigned int i = 1; i < history_data.size(); i++)
3146 *   ofs << i << ", " << history_data[i] << "\n";
3147 *   }
3148 *   }
3149 *  
3150 * @endcode
3151 *
3152 * This function saves the computed fields into a vtu file. This time we also
3153 * save the projected exact solution and the @f$L^2@f$ error norm. The exact
3154 * solution is only saved if the `Settings::project_exact_solution = true`
3155 *
3156 * @code
3157 *   void Solver::save() const
3158 *   {
3159 *   std::vector<std::string> solution_names(3, "VectorField");
3160 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
3161 *   interpretation(3,
3163 *  
3164 *   DataOut<3> data_out;
3165 *  
3166 *   data_out.add_data_vector(dof_handler_Hdiv,
3167 *   solution_Hdiv,
3168 *   solution_names,
3169 *   interpretation);
3170 *  
3171 *   if (Settings::project_exact_solution)
3172 *   {
3173 *   std::vector<std::string> solution_names_ex(3, "VectorFieldExact");
3174 *  
3175 *   data_out.add_data_vector(dof_handler_Hdiv,
3176 *   projected_exact_solution,
3177 *   solution_names_ex,
3178 *   interpretation);
3179 *   }
3180 *  
3181 *   if (exact_solution)
3182 *   {
3183 *   data_out.add_data_vector(L2_per_cell, "L2norm");
3184 *   }
3185 *  
3186 *   DataOutBase::VtkFlags flags;
3187 *   flags.write_higher_order_cells = true;
3188 *   data_out.set_flags(flags);
3189 *  
3190 *   const MappingQ<3> mapping(mapping_degree);
3191 *  
3192 *   data_out.build_patches(mapping,
3193 *   fe_Hdiv.degree + 2,
3195 *  
3196 *   std::ofstream ofs(fname + ".vtu");
3197 *   data_out.write_vtu(ofs);
3198 *   }
3199 *  
3200 *   void Solver::run()
3201 *   {
3202 *   {
3203 *   TimerOutput::Scope timer_section(timer, "Setup");
3204 *   setup();
3205 *   }
3206 *   {
3207 *   TimerOutput::Scope timer_section(timer, "Assemble");
3208 *   assemble();
3209 *   }
3210 *   {
3211 *   TimerOutput::Scope timer_section(timer, "Solve");
3212 *   solve();
3213 *   }
3214 *  
3215 *   if (exact_solution)
3216 *   {
3217 *   {
3218 *   TimerOutput::Scope timer_section(timer, "Compute error norms");
3219 *   compute_error_norms();
3220 *   }
3221 *  
3222 *   if (Settings::project_exact_solution)
3223 *   {
3224 *   {
3225 *   TimerOutput::Scope timer_section(timer, "Project exact solution");
3226 *   project_exact_solution_fcn();
3227 *   }
3228 *   }
3229 *   }
3230 *   {
3231 *   TimerOutput::Scope timer_section(timer, "Save");
3232 *   save();
3233 *   }
3234 *   {
3235 *   TimerOutput::Scope timer_section(timer, "Clear");
3236 *   clear();
3237 *   }
3238 *   }
3239 *   } // namespace ProjectorHcurlToHdiv
3240 *  
3241 *  
3242 * @endcode
3243 *
3244 *
3245 * <a name="step_97-Themainloop"></a>
3246 * <h3>The main loop</h3>
3247 *
3248
3249 *
3250 * This class contains the main loop of the program.
3251 *
3252 * @code
3253 *   class MagneticProblem
3254 *   {
3255 *   public:
3256 *   void run()
3257 *   {
3258 *   if (Settings::n_threads_max)
3259 *   MultithreadInfo::set_thread_limit(Settings::n_threads_max);
3260 *  
3261 *   MainOutputTable table_Jf(3);
3262 *   MainOutputTable table_B(3);
3263 *  
3264 *   table_Jf.clear();
3265 *   table_B.clear();
3266 *  
3267 *   std::cout << "Solving for (p = " << Settings::fe_degree
3268 *   << "): " << std::flush;
3269 *  
3270 *   for (unsigned int r = 6; r < 10; r++) // Mesh refinement parameter.
3271 *   {
3272 *   table_Jf.add_value("r", r);
3273 *   table_Jf.add_value("p", Settings::fe_degree);
3274 *  
3275 *   table_B.add_value("r", r);
3276 *   table_B.add_value("p", Settings::fe_degree);
3277 *  
3278 * @endcode
3279 *
3280 * Stage 1. Computing @f$\vec{T}@f$.
3281 *
3282
3283 *
3284 *
3285 * @code
3286 *   std::cout << "T " << std::flush;
3287 *  
3288 *   SolverT::Solver T(Settings::fe_degree,
3289 *   r,
3290 *   Settings::mapping_degree,
3291 *   Settings::eta_squared_T,
3292 *   "T_p" + std::to_string(Settings::fe_degree) + "_r" +
3293 *   std::to_string(r));
3294 *  
3295 *   T.run();
3296 *  
3297 * @endcode
3298 *
3299 * Stage 2. Computing @f$\vec{J}_f@f$.
3300 *
3301
3302 *
3303 *
3304 * @code
3305 *   std::cout << "Jf " << std::flush;
3306 *  
3307 *   ExactSolutions::FreeCurrentDensity Jf_exact;
3308 *  
3309 *   ProjectorHcurlToHdiv::Solver Jf(Settings::fe_degree,
3310 *   Settings::mapping_degree,
3311 *   T.get_tria(),
3312 *   T.get_dof_handler(),
3313 *   T.get_solution(),
3314 *   "Jf_p" +
3315 *   std::to_string(Settings::fe_degree) +
3316 *   "_r" + std::to_string(r),
3317 *   &Jf_exact);
3318 *  
3319 *   Jf.run();
3320 *  
3321 *   table_Jf.add_value("ndofs", Jf.get_n_dofs());
3322 *   table_Jf.add_value("ncells", Jf.get_n_cells());
3323 *   table_Jf.add_value("L2", Jf.get_L2_norm());
3324 *  
3325 * @endcode
3326 *
3327 * Stage 3. Computing @f$\vec{A}@f$.
3328 *
3329
3330 *
3331 *
3332 * @code
3333 *   std::cout << "A " << std::flush;
3334 *  
3335 *   SolverA::Solver A(Settings::fe_degree,
3336 *   Settings::mapping_degree,
3337 *   T.get_tria(),
3338 *   T.get_dof_handler(),
3339 *   T.get_solution(),
3340 *   Settings::eta_squared_A,
3341 *   "A_p" + std::to_string(Settings::fe_degree) + "_r" +
3342 *   std::to_string(r));
3343 *  
3344 *   A.run();
3345 *  
3346 * @endcode
3347 *
3348 * Stage 4. Computing @f$\vec{B}@f$.
3349 *
3350
3351 *
3352 *
3353 * @code
3354 *   std::cout << "B " << std::flush;
3355 *  
3356 *   ExactSolutions::MagneticField B_exact;
3357 *  
3358 *   ProjectorHcurlToHdiv::Solver B(Settings::fe_degree,
3359 *   Settings::mapping_degree,
3360 *   T.get_tria(),
3361 *   A.get_dof_handler(),
3362 *   A.get_solution(),
3363 *   "B_p" +
3364 *   std::to_string(Settings::fe_degree) +
3365 *   "_r" + std::to_string(r),
3366 *   &B_exact);
3367 *   B.run();
3368 *  
3369 *   table_B.add_value("ndofs", B.get_n_dofs());
3370 *   table_B.add_value("ncells", B.get_n_cells());
3371 *   table_B.add_value("L2", B.get_L2_norm());
3372 * @endcode
3373 *
3374 * End stage 4.
3375 *
3376
3377 *
3378 *
3379 * @code
3380 *   table_Jf.save("table_Jf_p" + std::to_string(Settings::fe_degree));
3381 *   table_B.save("table_B_p" + std::to_string(Settings::fe_degree));
3382 *   }
3383 *   std::cout << std::endl;
3384 *   }
3385 *   };
3386 *  
3387 *   int main()
3388 *   {
3389 *   try
3390 *   {
3391 *   MagneticProblem problem;
3392 *   problem.run();
3393 *   }
3394 *   catch (std::exception &exc)
3395 *   {
3396 *   std::cerr << std::endl
3397 *   << std::endl
3398 *   << "----------------------------------------------------"
3399 *   << std::endl;
3400 *   std::cerr << "Exception on processing: " << std::endl
3401 *   << exc.what() << std::endl
3402 *   << "Aborting!" << std::endl
3403 *   << "----------------------------------------------------"
3404 *   << std::endl;
3405 *  
3406 *   return 1;
3407 *   }
3408 *   catch (...)
3409 *   {
3410 *   std::cerr << std::endl
3411 *   << std::endl
3412 *   << "----------------------------------------------------"
3413 *   << std::endl;
3414 *   std::cerr << "Unknown exception!" << std::endl
3415 *   << "Aborting!" << std::endl
3416 *   << "----------------------------------------------------"
3417 *   << std::endl;
3418 *   return 1;
3419 *   }
3420 *  
3421 *   return 0;
3422 *   }
3423 * @endcode
3424<a name="step_97-Results"></a><h1>Results</h1>
3425
3426
3427The program generates the following output in the command line interface by
3428default.
3429
3430@code
3431Solving for (p = 0): T Jf A B T Jf A B T Jf A B T Jf A B
3432@endcode
3433
3434The program assumes the finite elements of the lowermost degree, @f$p = 0@f$.
3435To change the degree of the finite elements, say @f$p = 2@f$, one needs to change
3436the setting `Settings::fe_degree = 2` and rebuild the program.
3437
3438The program also dumps a number of files in the current directory. In the default
3439configuration these files are:
3440- vtu files. They contain the computed vector fields. Recall that the spherical
3441 manifold is attached to many cell faces. Consequently, these cell faces are
3442 curved. They look more like patches of a sphere. Furthermore, the shape
3443 functions are mapped from the reference cell to the real mesh cells by the
3444 second-order mapping to accommodate the cells with curved faces. For these
3445 reasons, one needs to use a visualization software that can deal with curved
3446 faces and the higher-order mapping. A fresh version of ParaView is recommended.
3447 Visit will not do. The <a href="https://github.com/dealii/dealii/wiki/Notes-on-visualizing-high-order-output">
3448 Notes on visualizing high order output</a> provide more information on this topic.
3449- tex files. These files contain the convergence tables.
3450
3451The following provides examples of the convergence tables simulated with the
3452default settings for three different degrees of the finite elements,
3453@f$p = 0, 1, 2@f$.
3454
3455| p | r | cells | dofs |@f$\|e\|_{L^2}@f$|@f$\alpha_{L^2}@f$|
3456|-- |---|-------|---------|----------|------|
3457| 0 | 6 | 4625 | 13950 | 1.66e-01 | - |
3458| 0 | 7 | 7992 | 24084 | 1.38e-01 | 0.99 |
3459| 0 | 8 | 12691 | 38220 | 1.19e-01 | 0.99 |
3460| 0 | 9 | 18944 | 57024 | 1.04e-01 | 0.99 |
3461| 1 | 6 | 4625 | 111300 | 8.12e-04 | - |
3462| 1 | 7 | 7992 | 192240 | 4.97e-04 | 2.69 |
3463| 1 | 8 | 12691 | 305172 | 3.32e-04 | 2.61 |
3464| 1 | 9 | 18944 | 455424 | 2.37e-04 | 2.54 |
3465| 2 | 6 | 4625 | 375300 | 6.78e-04 | - |
3466| 2 | 7 | 7992 | 648324 | 3.94e-04 | 2.97 |
3467| 2 | 8 | 12691 | 1029294 | 2.49e-04 | 2.98 |
3468| 2 | 9 | 18944 | 1536192 | 1.67e-04 | 2.99 |
3469
3470**Table 1. Convergence table. Free-current density, @f$\vec{J}_f@f$.
3471
3472| p | r | cells | dofs |@f$\|e\|_{L^2}@f$|@f$\alpha_{L^2}@f$|
3473|-- |---|-------|---------|----------|------|
3474| 0 | 6 | 4625 | 13950 | 8.84e-08 | - |
3475| 0 | 7 | 7992 | 24084 | 7.36e-08 | 1.00 |
3476| 0 | 8 | 12691 | 38220 | 6.30e-08 | 1.01 |
3477| 0 | 9 | 18944 | 57024 | 5.51e-08 | 1.00 |
3478| 1 | 6 | 4625 | 111300 | 4.41e-09 | - |
3479| 1 | 7 | 7992 | 192240 | 3.11e-09 | 1.91 |
3480| 1 | 8 | 12691 | 305172 | 2.23e-09 | 2.18 |
3481| 1 | 9 | 18944 | 455424 | 1.71e-09 | 1.96 |
3482| 2 | 6 | 4625 | 375300 | 1.84e-10 | - |
3483| 2 | 7 | 7992 | 648324 | 1.03e-10 | 3.21 |
3484| 2 | 8 | 12691 | 1029294 | 6.08e-11 | 3.40 |
3485| 2 | 9 | 18944 | 1536192 | 4.04e-11 | 3.07 |
3486
3487**Table 2. Convergence table. Magnetic field, @f$\vec{B}@f$.
3488
3489The following notations were used in the headers of the tables:
3490
3491- p - the degree of the finite elements.
3492
3493- r - the mesh refinement parameter, i.e., the number of nodes on the transfinite
3494lines.
3495
3496- cells - the total amount of active cells.
3497
3498- dofs - the amount of degrees of freedom.
3499
3500-@f$\|e\|_{L^2}@f$ - the @f$L^2@f$ error norm.
3501
3502-@f$\alpha_{L^2}@f$ - the order of convergence of the @f$L^2@f$ error norm.
3503
3504If `Settings::log_cg_convergence = true`, the program saves the convergence data
3505of the CG solver into csv files.
3506
3507The vector representations of the calculated vector fields, @f$\vec{J}_f@f$ and
3508@f$\vec{B}@f$, are illustrated above by the first figure on this page. The figures
3509below illustrate slices of the magnitudes of these fields. The figures below
3510were simulated with @f$p = 2@f$ and @f$r = 9@f$. Visual inspection of the vector
3511potentials is not very informative as their conservative portions are unknown.
3512
3513@htmlonly
3514<p align="center">
3515 <img src="https://www.dealii.org/images/steps/developer/step-97-Jf.svg" alt="The result - free-current
3516 density" height="531">
3517</p>
3518@endhtmlonly
3519
3520@htmlonly
3521<p align="center">
3522 <img src="https://www.dealii.org/images/steps/developer/step-97-B.svg" alt="The result - magnetic field
3523 density" height="531">
3524</p>
3525@endhtmlonly
3526
3527<a name="step_97-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
3528
3529
3530Repeat the simulations for the three types of the boundary conditions,
3531Dirichlet, Neumann, and Robin. The Robin boundary condition is supposed to be
3532superior to the other two. Look at the simulated data to see that this is indeed
3533the case. You can save the projected exact solution next to the simulated
3534solutions into the vtu files, just set `Settings::project_exact_solution = true`.
3535ParaView has "Plot Over Line" filter. You can use this filter to visualize the
3536difference between the exact solution and a solution simulated with a particular
3537boundary condition. You can also draw conclusions by observing the convergence
3538tables. Keep in mind the @f$\eta^2@f$ parameter. Increase it if the CG solver chokes
3539while you are experimenting. Note that the benefits offered by the Robin
3540boundary condition are observed the best when higher-order finite elements are
3541used, i.e., @f$p = 1@f$ and @f$p = 2@f$.
3542
3543The Robin boundary condition as described above is also called the first-order
3544asymptotic boundary condition (ABC). There exist ABCs of higher orders
3545@cite gratkowski2010p. Implement and test the second-order ABC to see if it
3546performs any better. There exist improvised asymptotic boundary conditions, IABCs,
3547@cite meeker2013a. Try to implement the first order IABC.
3548 *
3549 *
3550<a name="step_97-PlainProg"></a>
3551<h1> The plain program</h1>
3552@include "step-97.cc"
3553*/
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
@ curved_inner_cells
Definition data_out.h:211
void attach_triangulation(Triangulation< dim, spacedim > &tria)
Definition grid_in.cc:153
static void set_thread_limit(const unsigned int max_threads=numbers::invalid_unsigned_int)
Definition point.h:113
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
@ cpu_and_wall_times_grouped
Definition timer.h:659
#define Assert(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
void project_boundary_values_curl_conforming_l2(const DoFHandler< dim, dim > &dof_handler, const unsigned int first_vector_component, const Function< dim, number > &boundary_function, const types::boundary_id boundary_component, AffineConstraints< number > &constraints, const Mapping< dim > &mapping)
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
@ matrix
Contents is actually a matrix.
constexpr char L
constexpr char T
constexpr types::blas_int zero
constexpr char A
constexpr types::blas_int one
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void apply(const Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type &team_member, const Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_data, const ViewTypeIn in, ViewTypeOut out)
constexpr ReturnType< rank, T >::value_type & extract(T &t, const ArrayType &indices)
double compute_global_error(const Triangulation< dim, spacedim > &tria, const InVector &cellwise_error, const NormType &norm, const double exponent=2.)
void integrate_difference(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
void project(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const bool enforce_zero_boundary=false, const Quadrature< dim - 1 > &q_boundary=(dim > 1 ? QGauss< dim - 1 >(2) :Quadrature< dim - 1 >()), const bool project_to_boundary_first=false)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
long double gamma(const unsigned int n)
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
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
Definition types.h:32
unsigned int material_id
Definition types.h:184
unsigned int global_dof_index
Definition types.h:94
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition tensor.h:2672