deal.II version GIT relicensing-1721-g8100761196 2024-08-31 12:30:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
step-44.h
Go to the documentation of this file.
1)
1884 *   {}
1885 *  
1886 * @endcode
1887 *
1888 * Solve for the displacement using a Newton-Raphson method. We break this
1889 * function into the nonlinear loop and the function that solves the
1890 * linearized Newton-Raphson step:
1891 *
1892 * @code
1893 *   void solve_nonlinear_timestep(BlockVector<double> &solution_delta);
1894 *  
1895 *   std::pair<unsigned int, double>
1896 *   solve_linear_system(BlockVector<double> &newton_update);
1897 *  
1898 * @endcode
1899 *
1900 * Solution retrieval as well as post-processing and writing data to file :
1901 *
1902 * @code
1904 *   get_total_solution(const BlockVector<double> &solution_delta) const;
1905 *  
1906 *   void output_results() const;
1907 *  
1908 * @endcode
1909 *
1910 * Finally, some member variables that describe the current state: A
1911 * collection of the parameters used to describe the problem setup...
1912 *
1913 * @code
1914 *   Parameters::AllParameters parameters;
1915 *  
1916 * @endcode
1917 *
1918 * ...the volume of the reference configuration...
1919 *
1920 * @code
1921 *   double vol_reference;
1922 *  
1923 * @endcode
1924 *
1925 * ...and description of the geometry on which the problem is solved:
1926 *
1927 * @code
1929 *  
1930 * @endcode
1931 *
1932 * Also, keep track of the current time and the time spent evaluating
1933 * certain functions
1934 *
1935 * @code
1936 *   Time time;
1937 *   mutable TimerOutput timer;
1938 *  
1939 * @endcode
1940 *
1941 * A storage object for quadrature point information. As opposed to
1942 * @ref step_18 "step-18", deal.II's native quadrature point data manager is employed
1943 * here.
1944 *
1945 * @code
1946 *   CellDataStorage<typename Triangulation<dim>::cell_iterator,
1947 *   PointHistory<dim>>
1948 *   quadrature_point_history;
1949 *  
1950 * @endcode
1951 *
1952 * A description of the finite-element system including the displacement
1953 * polynomial degree, the degree-of-freedom handler, number of DoFs per
1954 * cell and the extractor objects used to retrieve information from the
1955 * solution vectors:
1956 *
1957 * @code
1958 *   const unsigned int degree;
1959 *   const FESystem<dim> fe;
1960 *   DoFHandler<dim> dof_handler;
1961 *   const unsigned int dofs_per_cell;
1962 *   const FEValuesExtractors::Vector u_fe;
1963 *   const FEValuesExtractors::Scalar p_fe;
1964 *   const FEValuesExtractors::Scalar J_fe;
1965 *  
1966 * @endcode
1967 *
1968 * Description of how the block-system is arranged. There are 3 blocks,
1969 * the first contains a vector DOF @f$\mathbf{u}@f$ while the other two
1970 * describe scalar DOFs, @f$\widetilde{p}@f$ and @f$\widetilde{J}@f$.
1971 *
1972 * @code
1973 *   static const unsigned int n_blocks = 3;
1974 *   static const unsigned int n_components = dim + 2;
1975 *   static const unsigned int first_u_component = 0;
1976 *   static const unsigned int p_component = dim;
1977 *   static const unsigned int J_component = dim + 1;
1978 *  
1979 *   enum
1980 *   {
1981 *   u_dof = 0,
1982 *   p_dof = 1,
1983 *   J_dof = 2
1984 *   };
1985 *  
1986 *   std::vector<types::global_dof_index> dofs_per_block;
1987 *   std::vector<types::global_dof_index> element_indices_u;
1988 *   std::vector<types::global_dof_index> element_indices_p;
1989 *   std::vector<types::global_dof_index> element_indices_J;
1990 *  
1991 * @endcode
1992 *
1993 * Rules for Gauss-quadrature on both the cell and faces. The number of
1994 * quadrature points on both cells and faces is recorded.
1995 *
1996 * @code
1997 *   const QGauss<dim> qf_cell;
1998 *   const QGauss<dim - 1> qf_face;
1999 *   const unsigned int n_q_points;
2000 *   const unsigned int n_q_points_f;
2001 *  
2002 * @endcode
2003 *
2004 * Objects that store the converged solution and right-hand side vectors,
2005 * as well as the tangent matrix. There is an AffineConstraints object used
2006 * to keep track of constraints. We make use of a sparsity pattern
2007 * designed for a block system.
2008 *
2009 * @code
2010 *   AffineConstraints<double> constraints;
2011 *   BlockSparsityPattern sparsity_pattern;
2012 *   BlockSparseMatrix<double> tangent_matrix;
2013 *   BlockVector<double> system_rhs;
2014 *   BlockVector<double> solution_n;
2015 *  
2016 * @endcode
2017 *
2018 * Then define a number of variables to store norms and update norms and
2019 * normalization factors.
2020 *
2021 * @code
2022 *   struct Errors
2023 *   {
2024 *   Errors()
2025 *   : norm(1.0)
2026 *   , u(1.0)
2027 *   , p(1.0)
2028 *   , J(1.0)
2029 *   {}
2030 *  
2031 *   void reset()
2032 *   {
2033 *   norm = 1.0;
2034 *   u = 1.0;
2035 *   p = 1.0;
2036 *   J = 1.0;
2037 *   }
2038 *   void normalize(const Errors &rhs)
2039 *   {
2040 *   if (rhs.norm != 0.0)
2041 *   norm /= rhs.norm;
2042 *   if (rhs.u != 0.0)
2043 *   u /= rhs.u;
2044 *   if (rhs.p != 0.0)
2045 *   p /= rhs.p;
2046 *   if (rhs.J != 0.0)
2047 *   J /= rhs.J;
2048 *   }
2049 *  
2050 *   double norm, u, p, J;
2051 *   };
2052 *  
2053 *   Errors error_residual, error_residual_0, error_residual_norm, error_update,
2054 *   error_update_0, error_update_norm;
2055 *  
2056 * @endcode
2057 *
2058 * Methods to calculate error measures
2059 *
2060 * @code
2061 *   void get_error_residual(Errors &error_residual);
2062 *  
2063 *   void get_error_update(const BlockVector<double> &newton_update,
2064 *   Errors &error_update);
2065 *  
2066 *   std::pair<double, double> get_error_dilation() const;
2067 *  
2068 * @endcode
2069 *
2070 * Compute the volume in the spatial configuration
2071 *
2072 * @code
2073 *   double compute_vol_current() const;
2074 *  
2075 * @endcode
2076 *
2077 * Print information to screen in a pleasing way...
2078 *
2079 * @code
2080 *   static void print_conv_header();
2081 *  
2082 *   void print_conv_footer();
2083 *   };
2084 *  
2085 * @endcode
2086 *
2087 *
2088 * <a name="step_44-ImplementationofthecodeSolidcodeclass"></a>
2089 * <h3>Implementation of the <code>Solid</code> class</h3>
2090 *
2091
2092 *
2093 *
2094 * <a name="step_44-Publicinterface"></a>
2095 * <h4>Public interface</h4>
2096 *
2097
2098 *
2099 * We initialize the Solid class using data extracted from the parameter file.
2100 *
2101 * @code
2102 *   template <int dim>
2103 *   Solid<dim>::Solid(const std::string &input_file)
2104 *   : parameters(input_file)
2105 *   , vol_reference(0.)
2106 *   , triangulation(Triangulation<dim>::maximum_smoothing)
2107 *   , time(parameters.end_time, parameters.delta_t)
2108 *   , timer(std::cout, TimerOutput::summary, TimerOutput::wall_times)
2109 *   , degree(parameters.poly_degree)
2110 *   ,
2111 * @endcode
2112 *
2113 * The Finite Element System is composed of dim continuous displacement
2114 * DOFs, and discontinuous pressure and dilatation DOFs. In an attempt to
2115 * satisfy the Babuska-Brezzi or LBB stability conditions (see Hughes
2116 * (2000)), we set up a @f$Q_n \times DGP_{n-1} \times DGP_{n-1}@f$
2117 * system. @f$Q_2 \times DGP_1 \times DGP_1@f$ elements satisfy this
2118 * condition, while @f$Q_1 \times DGP_0 \times DGP_0@f$ elements do
2119 * not. However, it has been shown that the latter demonstrate good
2120 * convergence characteristics nonetheless.
2121 *
2122 * @code
2123 *   fe(FE_Q<dim>(parameters.poly_degree) ^ dim, // displacement
2124 *   FE_DGP<dim>(parameters.poly_degree - 1), // pressure
2125 *   FE_DGP<dim>(parameters.poly_degree - 1)) // dilatation
2126 *   , dof_handler(triangulation)
2127 *   , dofs_per_cell(fe.n_dofs_per_cell())
2128 *   , u_fe(first_u_component)
2129 *   , p_fe(p_component)
2130 *   , J_fe(J_component)
2131 *   , dofs_per_block(n_blocks)
2132 *   , qf_cell(parameters.quad_order)
2133 *   , qf_face(parameters.quad_order)
2134 *   , n_q_points(qf_cell.size())
2135 *   , n_q_points_f(qf_face.size())
2136 *   {
2137 *   Assert(dim == 2 || dim == 3,
2138 *   ExcMessage("This problem only works in 2 or 3 space dimensions."));
2139 *   determine_component_extractors();
2140 *   }
2141 *  
2142 *  
2143 * @endcode
2144 *
2145 * In solving the quasi-static problem, the time becomes a loading parameter,
2146 * i.e. we increasing the loading linearly with time, making the two concepts
2147 * interchangeable. We choose to increment time linearly using a constant time
2148 * step size.
2149 *
2150
2151 *
2152 * We start the function with preprocessing, setting the initial dilatation
2153 * values, and then output the initial grid before starting the simulation
2154 * proper with the first time (and loading)
2155 * increment.
2156 *
2157
2158 *
2159 * Care must be taken (or at least some thought given) when imposing the
2160 * constraint @f$\widetilde{J}=1@f$ on the initial solution field. The constraint
2161 * corresponds to the determinant of the deformation gradient in the
2162 * undeformed configuration, which is the identity tensor. We use
2163 * FE_DGP bases to interpolate the dilatation field, thus we can't
2164 * simply set the corresponding dof to unity as they correspond to the
2165 * coefficients of a truncated Legendre polynomial.
2166 * Thus we use the VectorTools::project function to do the work for us.
2167 * The VectorTools::project function requires an argument
2168 * indicating the hanging node constraints. We have none in this program
2169 * So we have to create a constraint object. In its original state, constraint
2170 * objects are unsorted, and have to be sorted (using the
2171 * AffineConstraints::close function) before they can be used. Have a look at
2172 * @ref step_21 "step-21" for more information. We only need to enforce the initial condition
2173 * on the dilatation. In order to do this, we make use of a
2174 * ComponentSelectFunction which acts as a mask and sets the J_component of
2175 * n_components to 1. This is exactly what we want. Have a look at its usage
2176 * in @ref step_20 "step-20" for more information.
2177 *
2178 * @code
2179 *   template <int dim>
2180 *   void Solid<dim>::run()
2181 *   {
2182 *   make_grid();
2183 *   system_setup();
2184 *   {
2185 *   AffineConstraints<double> constraints;
2186 *   constraints.close();
2187 *  
2188 *   const ComponentSelectFunction<dim> J_mask(J_component, n_components);
2189 *  
2191 *   dof_handler, constraints, QGauss<dim>(degree + 2), J_mask, solution_n);
2192 *   }
2193 *   output_results();
2194 *   time.increment();
2195 *  
2196 * @endcode
2197 *
2198 * We then declare the incremental solution update @f$\varDelta
2199 * \mathbf{\Xi} \dealcoloneq \{\varDelta \mathbf{u},\varDelta \widetilde{p},
2200 * \varDelta \widetilde{J} \}@f$ and start the loop over the time domain.
2201 *
2202
2203 *
2204 * At the beginning, we reset the solution update for this time step...
2205 *
2206 * @code
2207 *   BlockVector<double> solution_delta(dofs_per_block);
2208 *   while (time.current() < time.end())
2209 *   {
2210 *   solution_delta = 0.0;
2211 *  
2212 * @endcode
2213 *
2214 * ...solve the current time step and update total solution vector
2215 * @f$\mathbf{\Xi}_{\textrm{n}} = \mathbf{\Xi}_{\textrm{n-1}} +
2216 * \varDelta \mathbf{\Xi}@f$...
2217 *
2218 * @code
2219 *   solve_nonlinear_timestep(solution_delta);
2220 *   solution_n += solution_delta;
2221 *  
2222 * @endcode
2223 *
2224 * ...and plot the results before moving on happily to the next time
2225 * step:
2226 *
2227 * @code
2228 *   output_results();
2229 *   time.increment();
2230 *   }
2231 *   }
2232 *  
2233 *  
2234 * @endcode
2235 *
2236 *
2237 * <a name="step_44-Privateinterface"></a>
2238 * <h3>Private interface</h3>
2239 *
2240
2241 *
2242 *
2243 * <a name="step_44-Threadingbuildingblocksstructures"></a>
2244 * <h4>Threading-building-blocks structures</h4>
2245 *
2246
2247 *
2248 * The first group of private member functions is related to parallelization.
2249 * We use the Threading Building Blocks library (TBB) to perform as many
2250 * computationally intensive distributed tasks as possible. In particular, we
2251 * assemble the tangent matrix and right hand side vector, the static
2252 * condensation contributions, and update data stored at the quadrature points
2253 * using TBB. Our main tool for this is the WorkStream class (see the @ref
2254 * threads topic for more information).
2255 *
2256
2257 *
2258 * Firstly we deal with the tangent matrix and right-hand side assembly
2259 * structures. The PerTaskData object stores local contributions to the global
2260 * system.
2261 *
2262 * @code
2263 *   template <int dim>
2264 *   struct Solid<dim>::PerTaskData_ASM
2265 *   {
2266 *   FullMatrix<double> cell_matrix;
2267 *   Vector<double> cell_rhs;
2268 *   std::vector<types::global_dof_index> local_dof_indices;
2269 *  
2270 *   PerTaskData_ASM(const unsigned int dofs_per_cell)
2271 *   : cell_matrix(dofs_per_cell, dofs_per_cell)
2272 *   , cell_rhs(dofs_per_cell)
2273 *   , local_dof_indices(dofs_per_cell)
2274 *   {}
2275 *  
2276 *   void reset()
2277 *   {
2278 *   cell_matrix = 0.0;
2279 *   cell_rhs = 0.0;
2280 *   }
2281 *   };
2282 *  
2283 *  
2284 * @endcode
2285 *
2286 * On the other hand, the ScratchData object stores the larger objects such as
2287 * the shape-function values array (<code>Nx</code>) and a shape function
2288 * gradient and symmetric gradient vector which we will use during the
2289 * assembly.
2290 *
2291 * @code
2292 *   template <int dim>
2293 *   struct Solid<dim>::ScratchData_ASM
2294 *   {
2295 *   FEValues<dim> fe_values;
2296 *   FEFaceValues<dim> fe_face_values;
2297 *  
2298 *   std::vector<std::vector<double>> Nx;
2299 *   std::vector<std::vector<Tensor<2, dim>>> grad_Nx;
2300 *   std::vector<std::vector<SymmetricTensor<2, dim>>> symm_grad_Nx;
2301 *  
2302 *   ScratchData_ASM(const FiniteElement<dim> &fe_cell,
2303 *   const QGauss<dim> &qf_cell,
2304 *   const UpdateFlags uf_cell,
2305 *   const QGauss<dim - 1> &qf_face,
2306 *   const UpdateFlags uf_face)
2307 *   : fe_values(fe_cell, qf_cell, uf_cell)
2308 *   , fe_face_values(fe_cell, qf_face, uf_face)
2309 *   , Nx(qf_cell.size(), std::vector<double>(fe_cell.n_dofs_per_cell()))
2310 *   , grad_Nx(qf_cell.size(),
2311 *   std::vector<Tensor<2, dim>>(fe_cell.n_dofs_per_cell()))
2312 *   , symm_grad_Nx(qf_cell.size(),
2313 *   std::vector<SymmetricTensor<2, dim>>(
2314 *   fe_cell.n_dofs_per_cell()))
2315 *   {}
2316 *  
2317 *   ScratchData_ASM(const ScratchData_ASM &rhs)
2318 *   : fe_values(rhs.fe_values.get_fe(),
2319 *   rhs.fe_values.get_quadrature(),
2320 *   rhs.fe_values.get_update_flags())
2321 *   , fe_face_values(rhs.fe_face_values.get_fe(),
2322 *   rhs.fe_face_values.get_quadrature(),
2323 *   rhs.fe_face_values.get_update_flags())
2324 *   , Nx(rhs.Nx)
2325 *   , grad_Nx(rhs.grad_Nx)
2326 *   , symm_grad_Nx(rhs.symm_grad_Nx)
2327 *   {}
2328 *  
2329 *   void reset()
2330 *   {
2331 *   const unsigned int n_q_points = Nx.size();
2332 *   const unsigned int n_dofs_per_cell = Nx[0].size();
2333 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2334 *   {
2335 *   Assert(Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2336 *   Assert(grad_Nx[q_point].size() == n_dofs_per_cell,
2337 *   ExcInternalError());
2338 *   Assert(symm_grad_Nx[q_point].size() == n_dofs_per_cell,
2339 *   ExcInternalError());
2340 *   for (unsigned int k = 0; k < n_dofs_per_cell; ++k)
2341 *   {
2342 *   Nx[q_point][k] = 0.0;
2343 *   grad_Nx[q_point][k] = 0.0;
2344 *   symm_grad_Nx[q_point][k] = 0.0;
2345 *   }
2346 *   }
2347 *   }
2348 *   };
2349 *  
2350 *  
2351 * @endcode
2352 *
2353 * Then we define structures to assemble the statically condensed tangent
2354 * matrix. Recall that we wish to solve for a displacement-based formulation.
2355 * We do the condensation at the element level as the @f$\widetilde{p}@f$ and
2356 * @f$\widetilde{J}@f$ fields are element-wise discontinuous. As these operations
2357 * are matrix-based, we need to set up a number of matrices to store the local
2358 * contributions from a number of the tangent matrix sub-blocks. We place
2359 * these in the PerTaskData struct.
2360 *
2361
2362 *
2363 * We choose not to reset any data in the <code>reset()</code> function as the
2364 * matrix extraction and replacement tools will take care of this
2365 *
2366 * @code
2367 *   template <int dim>
2368 *   struct Solid<dim>::PerTaskData_SC
2369 *   {
2371 *   std::vector<types::global_dof_index> local_dof_indices;
2372 *  
2373 *   FullMatrix<double> k_orig;
2374 *   FullMatrix<double> k_pu;
2375 *   FullMatrix<double> k_pJ;
2376 *   FullMatrix<double> k_JJ;
2377 *   FullMatrix<double> k_pJ_inv;
2378 *   FullMatrix<double> k_bbar;
2379 *   FullMatrix<double> A;
2380 *   FullMatrix<double> B;
2381 *   FullMatrix<double> C;
2382 *  
2383 *   PerTaskData_SC(const unsigned int dofs_per_cell,
2384 *   const unsigned int n_u,
2385 *   const unsigned int n_p,
2386 *   const unsigned int n_J)
2387 *   : cell_matrix(dofs_per_cell, dofs_per_cell)
2388 *   , local_dof_indices(dofs_per_cell)
2389 *   , k_orig(dofs_per_cell, dofs_per_cell)
2390 *   , k_pu(n_p, n_u)
2391 *   , k_pJ(n_p, n_J)
2392 *   , k_JJ(n_J, n_J)
2393 *   , k_pJ_inv(n_p, n_J)
2394 *   , k_bbar(n_u, n_u)
2395 *   , A(n_J, n_u)
2396 *   , B(n_J, n_u)
2397 *   , C(n_p, n_u)
2398 *   {}
2399 *  
2400 *   void reset()
2401 *   {}
2402 *   };
2403 *  
2404 *  
2405 * @endcode
2406 *
2407 * The ScratchData object for the operations we wish to perform here is empty
2408 * since we need no temporary data, but it still needs to be defined for the
2409 * current implementation of TBB in deal.II. So we create a dummy struct for
2410 * this purpose.
2411 *
2412 * @code
2413 *   template <int dim>
2414 *   struct Solid<dim>::ScratchData_SC
2415 *   {
2416 *   void reset()
2417 *   {}
2418 *   };
2419 *  
2420 *  
2421 * @endcode
2422 *
2423 * And finally we define the structures to assist with updating the quadrature
2424 * point information. Similar to the SC assembly process, we do not need the
2425 * PerTaskData object (since there is nothing to store here) but must define
2426 * one nonetheless. Note that this is because for the operation that we have
2427 * here -- updating the data on quadrature points -- the operation is purely
2428 * local: the things we do on every cell get consumed on every cell, without
2429 * any global aggregation operation as is usually the case when using the
2430 * WorkStream class. The fact that we still have to define a per-task data
2431 * structure points to the fact that the WorkStream class may be ill-suited to
2432 * this operation (we could, in principle simply create a new task using
2433 * Threads::new_task for each cell) but there is not much harm done to doing
2434 * it this way anyway.
2435 * Furthermore, should there be different material models associated with a
2436 * quadrature point, requiring varying levels of computational expense, then
2437 * the method used here could be advantageous.
2438 *
2439 * @code
2440 *   template <int dim>
2441 *   struct Solid<dim>::PerTaskData_UQPH
2442 *   {
2443 *   void reset()
2444 *   {}
2445 *   };
2446 *  
2447 *  
2448 * @endcode
2449 *
2450 * The ScratchData object will be used to store an alias for the solution
2451 * vector so that we don't have to copy this large data structure. We then
2452 * define a number of vectors to extract the solution values and gradients at
2453 * the quadrature points.
2454 *
2455 * @code
2456 *   template <int dim>
2457 *   struct Solid<dim>::ScratchData_UQPH
2458 *   {
2459 *   const BlockVector<double> &solution_total;
2460 *  
2461 *   std::vector<Tensor<2, dim>> solution_grads_u_total;
2462 *   std::vector<double> solution_values_p_total;
2463 *   std::vector<double> solution_values_J_total;
2464 *  
2465 *   FEValues<dim> fe_values;
2466 *  
2467 *   ScratchData_UQPH(const FiniteElement<dim> &fe_cell,
2468 *   const QGauss<dim> &qf_cell,
2469 *   const UpdateFlags uf_cell,
2470 *   const BlockVector<double> &solution_total)
2471 *   : solution_total(solution_total)
2472 *   , solution_grads_u_total(qf_cell.size())
2473 *   , solution_values_p_total(qf_cell.size())
2474 *   , solution_values_J_total(qf_cell.size())
2475 *   , fe_values(fe_cell, qf_cell, uf_cell)
2476 *   {}
2477 *  
2478 *   ScratchData_UQPH(const ScratchData_UQPH &rhs)
2479 *   : solution_total(rhs.solution_total)
2480 *   , solution_grads_u_total(rhs.solution_grads_u_total)
2481 *   , solution_values_p_total(rhs.solution_values_p_total)
2482 *   , solution_values_J_total(rhs.solution_values_J_total)
2483 *   , fe_values(rhs.fe_values.get_fe(),
2484 *   rhs.fe_values.get_quadrature(),
2485 *   rhs.fe_values.get_update_flags())
2486 *   {}
2487 *  
2488 *   void reset()
2489 *   {
2490 *   const unsigned int n_q_points = solution_grads_u_total.size();
2491 *   for (unsigned int q = 0; q < n_q_points; ++q)
2492 *   {
2493 *   solution_grads_u_total[q] = 0.0;
2494 *   solution_values_p_total[q] = 0.0;
2495 *   solution_values_J_total[q] = 0.0;
2496 *   }
2497 *   }
2498 *   };
2499 *  
2500 *  
2501 * @endcode
2502 *
2503 *
2504 * <a name="step_44-Solidmake_grid"></a>
2505 * <h4>Solid::make_grid</h4>
2506 *
2507
2508 *
2509 * On to the first of the private member functions. Here we create the
2510 * triangulation of the domain, for which we choose the scaled cube with each
2511 * face given a boundary ID number. The grid must be refined at least once
2512 * for the indentation problem.
2513 *
2514
2515 *
2516 * We then determine the volume of the reference configuration and print it
2517 * for comparison:
2518 *
2519 * @code
2520 *   template <int dim>
2521 *   void Solid<dim>::make_grid()
2522 *   {
2523 *   GridGenerator::hyper_rectangle(
2524 *   triangulation,
2525 *   (dim == 3 ? Point<dim>(0.0, 0.0, 0.0) : Point<dim>(0.0, 0.0)),
2526 *   (dim == 3 ? Point<dim>(1.0, 1.0, 1.0) : Point<dim>(1.0, 1.0)),
2527 *   true);
2528 *   GridTools::scale(parameters.scale, triangulation);
2529 *   triangulation.refine_global(std::max(1U, parameters.global_refinement));
2530 *  
2531 *   vol_reference = GridTools::volume(triangulation);
2532 *   std::cout << "Grid:\n\t Reference volume: " << vol_reference << std::endl;
2533 *  
2534 * @endcode
2535 *
2536 * Since we wish to apply a Neumann BC to a patch on the top surface, we
2537 * must find the cell faces in this part of the domain and mark them with
2538 * a distinct boundary ID number. The faces we are looking for are on the
2539 * +y surface and will get boundary ID 6 (zero through five are already
2540 * used when creating the six faces of the cube domain):
2541 *
2542 * @code
2543 *   for (const auto &cell : triangulation.active_cell_iterators())
2544 *   for (const auto &face : cell->face_iterators())
2545 *   {
2546 *   if (face->at_boundary() == true &&
2547 *   face->center()[1] == 1.0 * parameters.scale)
2548 *   {
2549 *   if (dim == 3)
2550 *   {
2551 *   if (face->center()[0] < 0.5 * parameters.scale &&
2552 *   face->center()[2] < 0.5 * parameters.scale)
2553 *   face->set_boundary_id(6);
2554 *   }
2555 *   else
2556 *   {
2557 *   if (face->center()[0] < 0.5 * parameters.scale)
2558 *   face->set_boundary_id(6);
2559 *   }
2560 *   }
2561 *   }
2562 *   }
2563 *  
2564 *  
2565 * @endcode
2566 *
2567 *
2568 * <a name="step_44-Solidsystem_setup"></a>
2569 * <h4>Solid::system_setup</h4>
2570 *
2571
2572 *
2573 * Next we describe how the FE system is setup. We first determine the number
2574 * of components per block. Since the displacement is a vector component, the
2575 * first dim components belong to it, while the next two describe scalar
2576 * pressure and dilatation DOFs.
2577 *
2578 * @code
2579 *   template <int dim>
2580 *   void Solid<dim>::system_setup()
2581 *   {
2582 *   timer.enter_subsection("Setup system");
2583 *  
2584 *   std::vector<unsigned int> block_component(n_components,
2585 *   u_dof); // Displacement
2586 *   block_component[p_component] = p_dof; // Pressure
2587 *   block_component[J_component] = J_dof; // Dilatation
2588 *  
2589 * @endcode
2590 *
2591 * The DOF handler is then initialized and we renumber the grid in an
2592 * efficient manner. We also record the number of DOFs per block.
2593 *
2594 * @code
2595 *   dof_handler.distribute_dofs(fe);
2596 *   DoFRenumbering::Cuthill_McKee(dof_handler);
2597 *   DoFRenumbering::component_wise(dof_handler, block_component);
2598 *  
2599 *   dofs_per_block =
2600 *   DoFTools::count_dofs_per_fe_block(dof_handler, block_component);
2601 *  
2602 *   std::cout << "Triangulation:"
2603 *   << "\n\t Number of active cells: "
2604 *   << triangulation.n_active_cells()
2605 *   << "\n\t Number of degrees of freedom: " << dof_handler.n_dofs()
2606 *   << std::endl;
2607 *  
2608 * @endcode
2609 *
2610 * Setup the sparsity pattern and tangent matrix
2611 *
2612 * @code
2613 *   tangent_matrix.clear();
2614 *   {
2615 *   BlockDynamicSparsityPattern dsp(dofs_per_block, dofs_per_block);
2616 *  
2617 * @endcode
2618 *
2619 * The global system matrix initially has the following structure
2620 * @f{align*}{
2621 * \underbrace{\begin{bmatrix}
2622 * \mathsf{\mathbf{K}}_{uu} & \mathsf{\mathbf{K}}_{u\widetilde{p}} &
2623 * \mathbf{0}
2624 * \\ \mathsf{\mathbf{K}}_{\widetilde{p}u} & \mathbf{0} &
2625 * \mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}
2626 * \\ \mathbf{0} & \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}} &
2627 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
2628 * \end{bmatrix}}_{\mathsf{\mathbf{K}}(\mathbf{\Xi}_{\textrm{i}})}
2629 * \underbrace{\begin{bmatrix}
2630 * d \mathsf{u}
2631 * \\ d \widetilde{\mathsf{\mathbf{p}}}
2632 * \\ d \widetilde{\mathsf{\mathbf{J}}}
2633 * \end{bmatrix}}_{d \mathbf{\Xi}}
2634 * =
2635 * \underbrace{\begin{bmatrix}
2636 * \mathsf{\mathbf{F}}_{u}(\mathbf{u}_{\textrm{i}})
2637 * \\ \mathsf{\mathbf{F}}_{\widetilde{p}}(\widetilde{p}_{\textrm{i}})
2638 * \\ \mathsf{\mathbf{F}}_{\widetilde{J}}(\widetilde{J}_{\textrm{i}})
2639 * \end{bmatrix}}_{ \mathsf{\mathbf{F}}(\mathbf{\Xi}_{\textrm{i}}) } \, .
2640 * @f}
2641 * We optimize the sparsity pattern to reflect this structure
2642 * and prevent unnecessary data creation for the right-diagonal
2643 * block components.
2644 *
2645 * @code
2646 *   Table<2, DoFTools::Coupling> coupling(n_components, n_components);
2647 *   for (unsigned int ii = 0; ii < n_components; ++ii)
2648 *   for (unsigned int jj = 0; jj < n_components; ++jj)
2649 *   if (((ii < p_component) && (jj == J_component)) ||
2650 *   ((ii == J_component) && (jj < p_component)) ||
2651 *   ((ii == p_component) && (jj == p_component)))
2652 *   coupling[ii][jj] = DoFTools::none;
2653 *   else
2654 *   coupling[ii][jj] = DoFTools::always;
2655 *   DoFTools::make_sparsity_pattern(
2656 *   dof_handler, coupling, dsp, constraints, false);
2657 *   sparsity_pattern.copy_from(dsp);
2658 *   }
2659 *  
2660 *   tangent_matrix.reinit(sparsity_pattern);
2661 *  
2662 * @endcode
2663 *
2664 * We then set up storage vectors
2665 *
2666 * @code
2667 *   system_rhs.reinit(dofs_per_block);
2668 *   solution_n.reinit(dofs_per_block);
2669 *  
2670 * @endcode
2671 *
2672 * ...and finally set up the quadrature
2673 * point history:
2674 *
2675 * @code
2676 *   setup_qph();
2677 *  
2678 *   timer.leave_subsection();
2679 *   }
2680 *  
2681 *  
2682 * @endcode
2683 *
2684 *
2685 * <a name="step_44-Soliddetermine_component_extractors"></a>
2686 * <h4>Solid::determine_component_extractors</h4>
2687 * Next we compute some information from the FE system that describes which
2688 * local element DOFs are attached to which block component. This is used
2689 * later to extract sub-blocks from the global matrix.
2690 *
2691
2692 *
2693 * In essence, all we need is for the FESystem object to indicate to which
2694 * block component a DOF on the reference cell is attached to. Currently, the
2695 * interpolation fields are setup such that 0 indicates a displacement DOF, 1
2696 * a pressure DOF and 2 a dilatation DOF.
2697 *
2698 * @code
2699 *   template <int dim>
2700 *   void Solid<dim>::determine_component_extractors()
2701 *   {
2702 *   element_indices_u.clear();
2703 *   element_indices_p.clear();
2704 *   element_indices_J.clear();
2705 *  
2706 *   for (unsigned int k = 0; k < fe.n_dofs_per_cell(); ++k)
2707 *   {
2708 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
2709 *   if (k_group == u_dof)
2710 *   element_indices_u.push_back(k);
2711 *   else if (k_group == p_dof)
2712 *   element_indices_p.push_back(k);
2713 *   else if (k_group == J_dof)
2714 *   element_indices_J.push_back(k);
2715 *   else
2716 *   {
2717 *   Assert(k_group <= J_dof, ExcInternalError());
2718 *   }
2719 *   }
2720 *   }
2721 *  
2722 * @endcode
2723 *
2724 *
2725 * <a name="step_44-Solidsetup_qph"></a>
2726 * <h4>Solid::setup_qph</h4>
2727 * The method used to store quadrature information is already described in
2728 * @ref step_18 "step-18". Here we implement a similar setup for a SMP machine.
2729 *
2730
2731 *
2732 * Firstly the actual QPH data objects are created. This must be done only
2733 * once the grid is refined to its finest level.
2734 *
2735 * @code
2736 *   template <int dim>
2737 *   void Solid<dim>::setup_qph()
2738 *   {
2739 *   std::cout << " Setting up quadrature point data..." << std::endl;
2740 *  
2741 *   quadrature_point_history.initialize(triangulation.begin_active(),
2742 *   triangulation.end(),
2743 *   n_q_points);
2744 *  
2745 * @endcode
2746 *
2747 * Next we set up the initial quadrature point data.
2748 * Note that when the quadrature point data is retrieved,
2749 * it is returned as a vector of smart pointers.
2750 *
2751 * @code
2752 *   for (const auto &cell : triangulation.active_cell_iterators())
2753 *   {
2754 *   const std::vector<std::shared_ptr<PointHistory<dim>>> lqph =
2755 *   quadrature_point_history.get_data(cell);
2756 *   Assert(lqph.size() == n_q_points, ExcInternalError());
2757 *  
2758 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2759 *   lqph[q_point]->setup_lqp(parameters);
2760 *   }
2761 *   }
2762 *  
2763 * @endcode
2764 *
2765 *
2766 * <a name="step_44-Solidupdate_qph_incremental"></a>
2767 * <h4>Solid::update_qph_incremental</h4>
2768 * As the update of QP information occurs frequently and involves a number of
2769 * expensive operations, we define a multithreaded approach to distributing
2770 * the task across a number of CPU cores.
2771 *
2772
2773 *
2774 * To start this, we first we need to obtain the total solution as it stands
2775 * at this Newton increment and then create the initial copy of the scratch
2776 * and copy data objects:
2777 *
2778 * @code
2779 *   template <int dim>
2780 *   void
2781 *   Solid<dim>::update_qph_incremental(const BlockVector<double> &solution_delta)
2782 *   {
2783 *   timer.enter_subsection("Update QPH data");
2784 *   std::cout << " UQPH " << std::flush;
2785 *  
2786 *   const BlockVector<double> solution_total(
2787 *   get_total_solution(solution_delta));
2788 *  
2789 *   const UpdateFlags uf_UQPH(update_values | update_gradients);
2790 *   PerTaskData_UQPH per_task_data_UQPH;
2791 *   ScratchData_UQPH scratch_data_UQPH(fe, qf_cell, uf_UQPH, solution_total);
2792 *  
2793 * @endcode
2794 *
2795 * We then pass them and the one-cell update function to the WorkStream to
2796 * be processed:
2797 *
2798 * @code
2799 *   WorkStream::run(dof_handler.active_cell_iterators(),
2800 *   *this,
2801 *   &Solid::update_qph_incremental_one_cell,
2802 *   &Solid::copy_local_to_global_UQPH,
2803 *   scratch_data_UQPH,
2804 *   per_task_data_UQPH);
2805 *  
2806 *   timer.leave_subsection();
2807 *   }
2808 *  
2809 *  
2810 * @endcode
2811 *
2812 * Now we describe how we extract data from the solution vector and pass it
2813 * along to each QP storage object for processing.
2814 *
2815 * @code
2816 *   template <int dim>
2817 *   void Solid<dim>::update_qph_incremental_one_cell(
2818 *   const typename DoFHandler<dim>::active_cell_iterator &cell,
2819 *   ScratchData_UQPH &scratch,
2820 *   PerTaskData_UQPH & /*data*/)
2821 *   {
2822 *   const std::vector<std::shared_ptr<PointHistory<dim>>> lqph =
2823 *   quadrature_point_history.get_data(cell);
2824 *   Assert(lqph.size() == n_q_points, ExcInternalError());
2825 *  
2826 *   Assert(scratch.solution_grads_u_total.size() == n_q_points,
2827 *   ExcInternalError());
2828 *   Assert(scratch.solution_values_p_total.size() == n_q_points,
2829 *   ExcInternalError());
2830 *   Assert(scratch.solution_values_J_total.size() == n_q_points,
2831 *   ExcInternalError());
2832 *  
2833 *   scratch.reset();
2834 *  
2835 * @endcode
2836 *
2837 * We first need to find the values and gradients at quadrature points
2838 * inside the current cell and then we update each local QP using the
2839 * displacement gradient and total pressure and dilatation solution
2840 * values:
2841 *
2842 * @code
2843 *   scratch.fe_values.reinit(cell);
2844 *   scratch.fe_values[u_fe].get_function_gradients(
2845 *   scratch.solution_total, scratch.solution_grads_u_total);
2846 *   scratch.fe_values[p_fe].get_function_values(
2847 *   scratch.solution_total, scratch.solution_values_p_total);
2848 *   scratch.fe_values[J_fe].get_function_values(
2849 *   scratch.solution_total, scratch.solution_values_J_total);
2850 *  
2851 *   for (const unsigned int q_point :
2852 *   scratch.fe_values.quadrature_point_indices())
2853 *   lqph[q_point]->update_values(scratch.solution_grads_u_total[q_point],
2854 *   scratch.solution_values_p_total[q_point],
2855 *   scratch.solution_values_J_total[q_point]);
2856 *   }
2857 *  
2858 *  
2859 * @endcode
2860 *
2861 *
2862 * <a name="step_44-Solidsolve_nonlinear_timestep"></a>
2863 * <h4>Solid::solve_nonlinear_timestep</h4>
2864 *
2865
2866 *
2867 * The next function is the driver method for the Newton-Raphson scheme. At
2868 * its top we create a new vector to store the current Newton update step,
2869 * reset the error storage objects and print solver header.
2870 *
2871 * @code
2872 *   template <int dim>
2873 *   void Solid<dim>::solve_nonlinear_timestep(BlockVector<double> &solution_delta)
2874 *   {
2875 *   std::cout << std::endl
2876 *   << "Timestep " << time.get_timestep() << " @ " << time.current()
2877 *   << 's' << std::endl;
2878 *  
2879 *   BlockVector<double> newton_update(dofs_per_block);
2880 *  
2881 *   error_residual.reset();
2882 *   error_residual_0.reset();
2883 *   error_residual_norm.reset();
2884 *   error_update.reset();
2885 *   error_update_0.reset();
2886 *   error_update_norm.reset();
2887 *  
2888 *   print_conv_header();
2889 *  
2890 * @endcode
2891 *
2892 * We now perform a number of Newton iterations to iteratively solve the
2893 * nonlinear problem. Since the problem is fully nonlinear and we are
2894 * using a full Newton method, the data stored in the tangent matrix and
2895 * right-hand side vector is not reusable and must be cleared at each
2896 * Newton step. We then initially build the linear system and
2897 * check for convergence (and store this value in the first iteration).
2898 * The unconstrained DOFs of the rhs vector hold the out-of-balance
2899 * forces, and collectively determine whether or not the equilibrium
2900 * solution has been attained.
2901 *
2902
2903 *
2904 * Although for this particular problem we could potentially construct the
2905 * RHS vector before assembling the system matrix, for the sake of
2906 * extensibility we choose not to do so. The benefit to assembling the RHS
2907 * vector and system matrix separately is that the latter is an expensive
2908 * operation and we can potentially avoid an extra assembly process by not
2909 * assembling the tangent matrix when convergence is attained. However, this
2910 * makes parallelizing the code using MPI more difficult. Furthermore, when
2911 * extending the problem to the transient case additional contributions to
2912 * the RHS may result from the time discretization and application of
2913 * constraints for the velocity and acceleration fields.
2914 *
2915 * @code
2916 *   unsigned int newton_iteration = 0;
2917 *   for (; newton_iteration < parameters.max_iterations_NR; ++newton_iteration)
2918 *   {
2919 *   std::cout << ' ' << std::setw(2) << newton_iteration << ' '
2920 *   << std::flush;
2921 *  
2922 * @endcode
2923 *
2924 * We construct the linear system, but hold off on solving it
2925 * (a step that should be significantly more expensive than assembly):
2926 *
2927 * @code
2928 *   make_constraints(newton_iteration);
2929 *   assemble_system();
2930 *  
2931 * @endcode
2932 *
2933 * We can now determine the normalized residual error and check for
2934 * solution convergence:
2935 *
2936 * @code
2937 *   get_error_residual(error_residual);
2938 *   if (newton_iteration == 0)
2939 *   error_residual_0 = error_residual;
2940 *  
2941 *   error_residual_norm = error_residual;
2942 *   error_residual_norm.normalize(error_residual_0);
2943 *  
2944 *   if (newton_iteration > 0 && error_update_norm.u <= parameters.tol_u &&
2945 *   error_residual_norm.u <= parameters.tol_f)
2946 *   {
2947 *   std::cout << " CONVERGED! " << std::endl;
2948 *   print_conv_footer();
2949 *  
2950 *   break;
2951 *   }
2952 *  
2953 * @endcode
2954 *
2955 * If we have decided that we want to continue with the iteration, we
2956 * solve the linearized system:
2957 *
2958 * @code
2959 *   const std::pair<unsigned int, double> lin_solver_output =
2960 *   solve_linear_system(newton_update);
2961 *  
2962 * @endcode
2963 *
2964 * We can now determine the normalized Newton update error:
2965 *
2966 * @code
2967 *   get_error_update(newton_update, error_update);
2968 *   if (newton_iteration == 0)
2969 *   error_update_0 = error_update;
2970 *  
2971 *   error_update_norm = error_update;
2972 *   error_update_norm.normalize(error_update_0);
2973 *  
2974 * @endcode
2975 *
2976 * Lastly, since we implicitly accept the solution step we can perform
2977 * the actual update of the solution increment for the current time
2978 * step, update all quadrature point information pertaining to
2979 * this new displacement and stress state and continue iterating:
2980 *
2981 * @code
2982 *   solution_delta += newton_update;
2983 *   update_qph_incremental(solution_delta);
2984 *  
2985 *   std::cout << " | " << std::fixed << std::setprecision(3) << std::setw(7)
2986 *   << std::scientific << lin_solver_output.first << " "
2987 *   << lin_solver_output.second << " "
2988 *   << error_residual_norm.norm << " " << error_residual_norm.u
2989 *   << " " << error_residual_norm.p << " "
2990 *   << error_residual_norm.J << " " << error_update_norm.norm
2991 *   << " " << error_update_norm.u << " " << error_update_norm.p
2992 *   << " " << error_update_norm.J << " " << std::endl;
2993 *   }
2994 *  
2995 * @endcode
2996 *
2997 * At the end, if it turns out that we have in fact done more iterations
2998 * than the parameter file allowed, we raise an exception that can be
2999 * caught in the main() function. The call <code>AssertThrow(condition,
3000 * exc_object)</code> is in essence equivalent to <code>if (!cond) throw
3001 * exc_object;</code> but the former form fills certain fields in the
3002 * exception object that identify the location (filename and line number)
3003 * where the exception was raised to make it simpler to identify where the
3004 * problem happened.
3005 *
3006 * @code
3007 *   AssertThrow(newton_iteration < parameters.max_iterations_NR,
3008 *   ExcMessage("No convergence in nonlinear solver!"));
3009 *   }
3010 *  
3011 *  
3012 * @endcode
3013 *
3014 *
3015 * <a name="step_44-Solidprint_conv_headerandSolidprint_conv_footer"></a>
3016 * <h4>Solid::print_conv_header and Solid::print_conv_footer</h4>
3017 *
3018
3019 *
3020 * This program prints out data in a nice table that is updated
3021 * on a per-iteration basis. The next two functions set up the table
3022 * header and footer:
3023 *
3024 * @code
3025 *   template <int dim>
3026 *   void Solid<dim>::print_conv_header()
3027 *   {
3028 *   static const unsigned int l_width = 150;
3029 *  
3030 *   for (unsigned int i = 0; i < l_width; ++i)
3031 *   std::cout << '_';
3032 *   std::cout << std::endl;
3033 *  
3034 *   std::cout << " SOLVER STEP "
3035 *   << " | LIN_IT LIN_RES RES_NORM "
3036 *   << " RES_U RES_P RES_J NU_NORM "
3037 *   << " NU_U NU_P NU_J " << std::endl;
3038 *  
3039 *   for (unsigned int i = 0; i < l_width; ++i)
3040 *   std::cout << '_';
3041 *   std::cout << std::endl;
3042 *   }
3043 *  
3044 *  
3045 *  
3046 *   template <int dim>
3047 *   void Solid<dim>::print_conv_footer()
3048 *   {
3049 *   static const unsigned int l_width = 150;
3050 *  
3051 *   for (unsigned int i = 0; i < l_width; ++i)
3052 *   std::cout << '_';
3053 *   std::cout << std::endl;
3054 *  
3055 *   const std::pair<double, double> error_dil = get_error_dilation();
3056 *  
3057 *   std::cout << "Relative errors:" << std::endl
3058 *   << "Displacement:\t" << error_update.u / error_update_0.u
3059 *   << std::endl
3060 *   << "Force: \t\t" << error_residual.u / error_residual_0.u
3061 *   << std::endl
3062 *   << "Dilatation:\t" << error_dil.first << std::endl
3063 *   << "v / V_0:\t" << error_dil.second * vol_reference << " / "
3064 *   << vol_reference << " = " << error_dil.second << std::endl;
3065 *   }
3066 *  
3067 *  
3068 * @endcode
3069 *
3070 *
3071 * <a name="step_44-Solidget_error_dilation"></a>
3072 * <h4>Solid::get_error_dilation</h4>
3073 *
3074
3075 *
3076 * Calculate the volume of the domain in the spatial configuration
3077 *
3078 * @code
3079 *   template <int dim>
3080 *   double Solid<dim>::compute_vol_current() const
3081 *   {
3082 *   double vol_current = 0.0;
3083 *  
3084 *   FEValues<dim> fe_values(fe, qf_cell, update_JxW_values);
3085 *  
3086 *   for (const auto &cell : triangulation.active_cell_iterators())
3087 *   {
3088 *   fe_values.reinit(cell);
3089 *  
3090 * @endcode
3091 *
3092 * In contrast to that which was previously called for,
3093 * in this instance the quadrature point data is specifically
3094 * non-modifiable since we will only be accessing data.
3095 * We ensure that the right get_data function is called by
3096 * marking this update function as constant.
3097 *
3098 * @code
3099 *   const std::vector<std::shared_ptr<const PointHistory<dim>>> lqph =
3100 *   quadrature_point_history.get_data(cell);
3101 *   Assert(lqph.size() == n_q_points, ExcInternalError());
3102 *  
3103 *   for (const unsigned int q_point : fe_values.quadrature_point_indices())
3104 *   {
3105 *   const double det_F_qp = lqph[q_point]->get_det_F();
3106 *   const double JxW = fe_values.JxW(q_point);
3107 *  
3108 *   vol_current += det_F_qp * JxW;
3109 *   }
3110 *   }
3111 *   Assert(vol_current > 0.0, ExcInternalError());
3112 *   return vol_current;
3113 *   }
3114 *  
3115 * @endcode
3116 *
3117 * Calculate how well the dilatation @f$\widetilde{J}@f$ agrees with @f$J
3118 * \dealcoloneq \textrm{det}\ \mathbf{F}@f$ from the @f$L^2@f$ error @f$ \bigl[
3119 * \int_{\Omega_0} {[ J - \widetilde{J}]}^{2}\textrm{d}V \bigr]^{1/2}@f$.
3120 * We also return the ratio of the current volume of the
3121 * domain to the reference volume. This is of interest for incompressible
3122 * media where we want to check how well the isochoric constraint has been
3123 * enforced.
3124 *
3125 * @code
3126 *   template <int dim>
3127 *   std::pair<double, double> Solid<dim>::get_error_dilation() const
3128 *   {
3129 *   double dil_L2_error = 0.0;
3130 *  
3131 *   FEValues<dim> fe_values(fe, qf_cell, update_JxW_values);
3132 *  
3133 *   for (const auto &cell : triangulation.active_cell_iterators())
3134 *   {
3135 *   fe_values.reinit(cell);
3136 *  
3137 *   const std::vector<std::shared_ptr<const PointHistory<dim>>> lqph =
3138 *   quadrature_point_history.get_data(cell);
3139 *   Assert(lqph.size() == n_q_points, ExcInternalError());
3140 *  
3141 *   for (const unsigned int q_point : fe_values.quadrature_point_indices())
3142 *   {
3143 *   const double det_F_qp = lqph[q_point]->get_det_F();
3144 *   const double J_tilde_qp = lqph[q_point]->get_J_tilde();
3145 *   const double the_error_qp_squared =
3146 *   Utilities::fixed_power<2>((det_F_qp - J_tilde_qp));
3147 *   const double JxW = fe_values.JxW(q_point);
3148 *  
3149 *   dil_L2_error += the_error_qp_squared * JxW;
3150 *   }
3151 *   }
3152 *  
3153 *   return std::make_pair(std::sqrt(dil_L2_error),
3154 *   compute_vol_current() / vol_reference);
3155 *   }
3156 *  
3157 *  
3158 * @endcode
3159 *
3160 *
3161 * <a name="step_44-Solidget_error_residual"></a>
3162 * <h4>Solid::get_error_residual</h4>
3163 *
3164
3165 *
3166 * Determine the true residual error for the problem. That is, determine the
3167 * error in the residual for the unconstrained degrees of freedom. Note that
3168 * to do so, we need to ignore constrained DOFs by setting the residual in
3169 * these vector components to zero.
3170 *
3171 * @code
3172 *   template <int dim>
3173 *   void Solid<dim>::get_error_residual(Errors &error_residual)
3174 *   {
3175 *   BlockVector<double> error_res(dofs_per_block);
3176 *  
3177 *   for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
3178 *   if (!constraints.is_constrained(i))
3179 *   error_res(i) = system_rhs(i);
3180 *  
3181 *   error_residual.norm = error_res.l2_norm();
3182 *   error_residual.u = error_res.block(u_dof).l2_norm();
3183 *   error_residual.p = error_res.block(p_dof).l2_norm();
3184 *   error_residual.J = error_res.block(J_dof).l2_norm();
3185 *   }
3186 *  
3187 *  
3188 * @endcode
3189 *
3190 *
3191 * <a name="step_44-Solidget_error_update"></a>
3192 * <h4>Solid::get_error_update</h4>
3193 *
3194
3195 *
3196 * Determine the true Newton update error for the problem
3197 *
3198 * @code
3199 *   template <int dim>
3200 *   void Solid<dim>::get_error_update(const BlockVector<double> &newton_update,
3201 *   Errors &error_update)
3202 *   {
3203 *   BlockVector<double> error_ud(dofs_per_block);
3204 *   for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
3205 *   if (!constraints.is_constrained(i))
3206 *   error_ud(i) = newton_update(i);
3207 *  
3208 *   error_update.norm = error_ud.l2_norm();
3209 *   error_update.u = error_ud.block(u_dof).l2_norm();
3210 *   error_update.p = error_ud.block(p_dof).l2_norm();
3211 *   error_update.J = error_ud.block(J_dof).l2_norm();
3212 *   }
3213 *  
3214 *  
3215 *  
3216 * @endcode
3217 *
3218 *
3219 * <a name="step_44-Solidget_total_solution"></a>
3220 * <h4>Solid::get_total_solution</h4>
3221 *
3222
3223 *
3224 * This function provides the total solution, which is valid at any Newton
3225 * step. This is required as, to reduce computational error, the total
3226 * solution is only updated at the end of the timestep.
3227 *
3228 * @code
3229 *   template <int dim>
3230 *   BlockVector<double> Solid<dim>::get_total_solution(
3231 *   const BlockVector<double> &solution_delta) const
3232 *   {
3233 *   BlockVector<double> solution_total(solution_n);
3234 *   solution_total += solution_delta;
3235 *   return solution_total;
3236 *   }
3237 *  
3238 *  
3239 * @endcode
3240 *
3241 *
3242 * <a name="step_44-Solidassemble_system"></a>
3243 * <h4>Solid::assemble_system</h4>
3244 *
3245
3246 *
3247 * Since we use TBB for assembly, we simply setup a copy of the
3248 * data structures required for the process and pass them, along
3249 * with the assembly functions to the WorkStream object for processing. Note
3250 * that we must ensure that the matrix and RHS vector are reset before any
3251 * assembly operations can occur. Furthermore, since we are describing a
3252 * problem with Neumann BCs, we will need the face normals and so must specify
3253 * this in the face update flags.
3254 *
3255 * @code
3256 *   template <int dim>
3257 *   void Solid<dim>::assemble_system()
3258 *   {
3259 *   timer.enter_subsection("Assemble system");
3260 *   std::cout << " ASM_SYS " << std::flush;
3261 *  
3262 *   tangent_matrix = 0.0;
3263 *   system_rhs = 0.0;
3264 *  
3265 *   const UpdateFlags uf_cell(update_values | update_gradients |
3266 *   update_JxW_values);
3267 *   const UpdateFlags uf_face(update_values | update_normal_vectors |
3268 *   update_JxW_values);
3269 *  
3270 *   PerTaskData_ASM per_task_data(dofs_per_cell);
3271 *   ScratchData_ASM scratch_data(fe, qf_cell, uf_cell, qf_face, uf_face);
3272 *  
3273 * @endcode
3274 *
3275 * The syntax used here to pass data to the WorkStream class
3276 * is discussed in @ref step_13 "step-13".
3277 *
3278 * @code
3279 *   WorkStream::run(
3280 *   dof_handler.active_cell_iterators(),
3281 *   [this](const typename DoFHandler<dim>::active_cell_iterator &cell,
3282 *   ScratchData_ASM &scratch,
3283 *   PerTaskData_ASM &data) {
3284 *   this->assemble_system_one_cell(cell, scratch, data);
3285 *   },
3286 *   [this](const PerTaskData_ASM &data) {
3287 *   this->constraints.distribute_local_to_global(data.cell_matrix,
3288 *   data.cell_rhs,
3289 *   data.local_dof_indices,
3290 *   tangent_matrix,
3291 *   system_rhs);
3292 *   },
3293 *   scratch_data,
3294 *   per_task_data);
3295 *  
3296 *   timer.leave_subsection();
3297 *   }
3298 *  
3299 * @endcode
3300 *
3301 * Of course, we still have to define how we assemble the tangent matrix
3302 * contribution for a single cell. We first need to reset and initialize some
3303 * of the scratch data structures and retrieve some basic information
3304 * regarding the DOF numbering on this cell. We can precalculate the cell
3305 * shape function values and gradients. Note that the shape function gradients
3306 * are defined with regard to the current configuration. That is
3307 * @f$\textrm{grad}\ \boldsymbol{\varphi} = \textrm{Grad}\ \boldsymbol{\varphi}
3308 * \ \mathbf{F}^{-1}@f$.
3309 *
3310 * @code
3311 *   template <int dim>
3312 *   void Solid<dim>::assemble_system_one_cell(
3313 *   const typename DoFHandler<dim>::active_cell_iterator &cell,
3314 *   ScratchData_ASM &scratch,
3315 *   PerTaskData_ASM &data) const
3316 *   {
3317 *   data.reset();
3318 *   scratch.reset();
3319 *   scratch.fe_values.reinit(cell);
3320 *   cell->get_dof_indices(data.local_dof_indices);
3321 *  
3322 *   const std::vector<std::shared_ptr<const PointHistory<dim>>> lqph =
3323 *   quadrature_point_history.get_data(cell);
3324 *   Assert(lqph.size() == n_q_points, ExcInternalError());
3325 *  
3326 *   for (const unsigned int q_point :
3327 *   scratch.fe_values.quadrature_point_indices())
3328 *   {
3329 *   const Tensor<2, dim> F_inv = lqph[q_point]->get_F_inv();
3330 *   for (const unsigned int k : scratch.fe_values.dof_indices())
3331 *   {
3332 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
3333 *  
3334 *   if (k_group == u_dof)
3335 *   {
3336 *   scratch.grad_Nx[q_point][k] =
3337 *   scratch.fe_values[u_fe].gradient(k, q_point) * F_inv;
3338 *   scratch.symm_grad_Nx[q_point][k] =
3339 *   symmetrize(scratch.grad_Nx[q_point][k]);
3340 *   }
3341 *   else if (k_group == p_dof)
3342 *   scratch.Nx[q_point][k] =
3343 *   scratch.fe_values[p_fe].value(k, q_point);
3344 *   else if (k_group == J_dof)
3345 *   scratch.Nx[q_point][k] =
3346 *   scratch.fe_values[J_fe].value(k, q_point);
3347 *   else
3348 *   Assert(k_group <= J_dof, ExcInternalError());
3349 *   }
3350 *   }
3351 *  
3352 * @endcode
3353 *
3354 * Now we build the local cell @ref GlossStiffnessMatrix "stiffness matrix" and RHS vector. Since the
3355 * global and local system matrices are symmetric, we can exploit this
3356 * property by building only the lower half of the local matrix and copying
3357 * the values to the upper half. So we only assemble half of the
3358 * @f$\mathsf{\mathbf{k}}_{uu}@f$, @f$\mathsf{\mathbf{k}}_{\widetilde{p}
3359 * \widetilde{p}} = \mathbf{0}@f$, @f$\mathsf{\mathbf{k}}_{\widetilde{J}
3360 * \widetilde{J}}@f$ blocks, while the whole
3361 * @f$\mathsf{\mathbf{k}}_{\widetilde{p} \widetilde{J}}@f$,
3362 * @f$\mathsf{\mathbf{k}}_{u \widetilde{J}} = \mathbf{0}@f$,
3363 * @f$\mathsf{\mathbf{k}}_{u \widetilde{p}}@f$ blocks are built.
3364 *
3365
3366 *
3367 * In doing so, we first extract some configuration dependent variables
3368 * from our quadrature history objects for the current quadrature point.
3369 *
3370 * @code
3371 *   for (const unsigned int q_point :
3372 *   scratch.fe_values.quadrature_point_indices())
3373 *   {
3374 *   const SymmetricTensor<2, dim> tau = lqph[q_point]->get_tau();
3375 *   const Tensor<2, dim> tau_ns = lqph[q_point]->get_tau();
3376 *   const SymmetricTensor<4, dim> Jc = lqph[q_point]->get_Jc();
3377 *   const double det_F = lqph[q_point]->get_det_F();
3378 *   const double p_tilde = lqph[q_point]->get_p_tilde();
3379 *   const double J_tilde = lqph[q_point]->get_J_tilde();
3380 *   const double dPsi_vol_dJ = lqph[q_point]->get_dPsi_vol_dJ();
3381 *   const double d2Psi_vol_dJ2 = lqph[q_point]->get_d2Psi_vol_dJ2();
3382 *   const SymmetricTensor<2, dim> &I =
3383 *   Physics::Elasticity::StandardTensors<dim>::I;
3384 *  
3385 * @endcode
3386 *
3387 * These two tensors store some precomputed data. Their use will
3388 * explained shortly.
3389 *
3390 * @code
3391 *   SymmetricTensor<2, dim> symm_grad_Nx_i_x_Jc;
3392 *   Tensor<1, dim> grad_Nx_i_comp_i_x_tau;
3393 *  
3394 * @endcode
3395 *
3396 * Next we define some aliases to make the assembly process easier to
3397 * follow.
3398 *
3399 * @code
3400 *   const std::vector<double> &N = scratch.Nx[q_point];
3401 *   const std::vector<SymmetricTensor<2, dim>> &symm_grad_Nx =
3402 *   scratch.symm_grad_Nx[q_point];
3403 *   const std::vector<Tensor<2, dim>> &grad_Nx = scratch.grad_Nx[q_point];
3404 *   const double JxW = scratch.fe_values.JxW(q_point);
3405 *  
3406 *   for (const unsigned int i : scratch.fe_values.dof_indices())
3407 *   {
3408 *   const unsigned int component_i =
3409 *   fe.system_to_component_index(i).first;
3410 *   const unsigned int i_group = fe.system_to_base_index(i).first.first;
3411 *  
3412 * @endcode
3413 *
3414 * We first compute the contributions
3415 * from the internal forces. Note, by
3416 * definition of the rhs as the negative
3417 * of the residual, these contributions
3418 * are subtracted.
3419 *
3420 * @code
3421 *   if (i_group == u_dof)
3422 *   data.cell_rhs(i) -= (symm_grad_Nx[i] * tau) * JxW;
3423 *   else if (i_group == p_dof)
3424 *   data.cell_rhs(i) -= N[i] * (det_F - J_tilde) * JxW;
3425 *   else if (i_group == J_dof)
3426 *   data.cell_rhs(i) -= N[i] * (dPsi_vol_dJ - p_tilde) * JxW;
3427 *   else
3428 *   Assert(i_group <= J_dof, ExcInternalError());
3429 *  
3430 * @endcode
3431 *
3432 * Before we go into the inner loop, we have one final chance to
3433 * introduce some optimizations. We've already taken into account
3434 * the symmetry of the system, and we can now precompute some
3435 * common terms that are repeatedly applied in the inner loop.
3436 * We won't be excessive here, but will rather focus on expensive
3437 * operations, namely those involving the rank-4 material stiffness
3438 * tensor and the rank-2 stress tensor.
3439 *
3440
3441 *
3442 * What we may observe is that both of these tensors are contracted
3443 * with shape function gradients indexed on the "i" DoF. This
3444 * implies that this particular operation remains constant as we
3445 * loop over the "j" DoF. For that reason, we can extract this from
3446 * the inner loop and save the many operations that, for each
3447 * quadrature point and DoF index "i" and repeated over index "j"
3448 * are required to double contract a rank-2 symmetric tensor with a
3449 * rank-4 symmetric tensor, and a rank-1 tensor with a rank-2
3450 * tensor.
3451 *
3452
3453 *
3454 * At the loss of some readability, this small change will reduce
3455 * the assembly time of the symmetrized system by about half when
3456 * using the simulation default parameters, and becomes more
3457 * significant as the h-refinement level increases.
3458 *
3459 * @code
3460 *   if (i_group == u_dof)
3461 *   {
3462 *   symm_grad_Nx_i_x_Jc = symm_grad_Nx[i] * Jc;
3463 *   grad_Nx_i_comp_i_x_tau = grad_Nx[i][component_i] * tau_ns;
3464 *   }
3465 *  
3466 * @endcode
3467 *
3468 * Now we're prepared to compute the tangent matrix contributions:
3469 *
3470 * @code
3471 *   for (const unsigned int j :
3472 *   scratch.fe_values.dof_indices_ending_at(i))
3473 *   {
3474 *   const unsigned int component_j =
3475 *   fe.system_to_component_index(j).first;
3476 *   const unsigned int j_group =
3477 *   fe.system_to_base_index(j).first.first;
3478 *  
3479 * @endcode
3480 *
3481 * This is the @f$\mathsf{\mathbf{k}}_{uu}@f$
3482 * contribution. It comprises a material contribution, and a
3483 * geometrical stress contribution which is only added along
3484 * the local matrix diagonals:
3485 *
3486 * @code
3487 *   if ((i_group == j_group) && (i_group == u_dof))
3488 *   {
3489 * @endcode
3490 *
3491 * The material contribution:
3492 *
3493 * @code
3494 *   data.cell_matrix(i, j) += symm_grad_Nx_i_x_Jc *
3495 *   symm_grad_Nx[j] * JxW;
3496 *  
3497 * @endcode
3498 *
3499 * The geometrical stress contribution:
3500 *
3501 * @code
3502 *   if (component_i == component_j)
3503 *   data.cell_matrix(i, j) +=
3504 *   grad_Nx_i_comp_i_x_tau * grad_Nx[j][component_j] * JxW;
3505 *   }
3506 * @endcode
3507 *
3508 * Next is the @f$\mathsf{\mathbf{k}}_{ \widetilde{p} u}@f$
3509 * contribution
3510 *
3511 * @code
3512 *   else if ((i_group == p_dof) && (j_group == u_dof))
3513 *   {
3514 *   data.cell_matrix(i, j) += N[i] * det_F *
3515 *   (symm_grad_Nx[j] * I) * JxW;
3516 *   }
3517 * @endcode
3518 *
3519 * and lastly the @f$\mathsf{\mathbf{k}}_{ \widetilde{J}
3520 * \widetilde{p}}@f$ and @f$\mathsf{\mathbf{k}}_{ \widetilde{J}
3521 * \widetilde{J}}@f$ contributions:
3522 *
3523 * @code
3524 *   else if ((i_group == J_dof) && (j_group == p_dof))
3525 *   data.cell_matrix(i, j) -= N[i] * N[j] * JxW;
3526 *   else if ((i_group == j_group) && (i_group == J_dof))
3527 *   data.cell_matrix(i, j) += N[i] * d2Psi_vol_dJ2 * N[j] * JxW;
3528 *   else
3529 *   Assert((i_group <= J_dof) && (j_group <= J_dof),
3530 *   ExcInternalError());
3531 *   }
3532 *   }
3533 *   }
3534 *  
3535 * @endcode
3536 *
3537 * Next we assemble the Neumann contribution. We first check to see it the
3538 * cell face exists on a boundary on which a traction is applied and add
3539 * the contribution if this is the case.
3540 *
3541 * @code
3542 *   for (const auto &face : cell->face_iterators())
3543 *   if (face->at_boundary() && face->boundary_id() == 6)
3544 *   {
3545 *   scratch.fe_face_values.reinit(cell, face);
3546 *  
3547 *   for (const unsigned int f_q_point :
3548 *   scratch.fe_face_values.quadrature_point_indices())
3549 *   {
3550 *   const Tensor<1, dim> &N =
3551 *   scratch.fe_face_values.normal_vector(f_q_point);
3552 *  
3553 * @endcode
3554 *
3555 * Using the face normal at this quadrature point we specify the
3556 * traction in reference configuration. For this problem, a
3557 * defined pressure is applied in the reference configuration.
3558 * The direction of the applied traction is assumed not to
3559 * evolve with the deformation of the domain. The traction is
3560 * defined using the first Piola-Kirchhoff stress is simply
3561 * @f$\mathbf{t} = \mathbf{P}\mathbf{N} = [p_0 \mathbf{I}]
3562 * \mathbf{N} = p_0 \mathbf{N}@f$ We use the time variable to
3563 * linearly ramp up the pressure load.
3564 *
3565
3566 *
3567 * Note that the contributions to the right hand side vector we
3568 * compute here only exist in the displacement components of the
3569 * vector.
3570 *
3571 * @code
3572 *   static const double p0 =
3573 *   -4.0 / (parameters.scale * parameters.scale);
3574 *   const double time_ramp = (time.current() / time.end());
3575 *   const double pressure = p0 * parameters.p_p0 * time_ramp;
3576 *   const Tensor<1, dim> traction = pressure * N;
3577 *  
3578 *   for (const unsigned int i : scratch.fe_values.dof_indices())
3579 *   {
3580 *   const unsigned int i_group =
3581 *   fe.system_to_base_index(i).first.first;
3582 *  
3583 *   if (i_group == u_dof)
3584 *   {
3585 *   const unsigned int component_i =
3586 *   fe.system_to_component_index(i).first;
3587 *   const double Ni =
3588 *   scratch.fe_face_values.shape_value(i, f_q_point);
3589 *   const double JxW = scratch.fe_face_values.JxW(f_q_point);
3590 *  
3591 *   data.cell_rhs(i) += (Ni * traction[component_i]) * JxW;
3592 *   }
3593 *   }
3594 *   }
3595 *   }
3596 *  
3597 * @endcode
3598 *
3599 * Finally, we need to copy the lower half of the local matrix into the
3600 * upper half:
3601 *
3602 * @code
3603 *   for (const unsigned int i : scratch.fe_values.dof_indices())
3604 *   for (const unsigned int j :
3605 *   scratch.fe_values.dof_indices_starting_at(i + 1))
3606 *   data.cell_matrix(i, j) = data.cell_matrix(j, i);
3607 *   }
3608 *  
3609 *  
3610 *  
3611 * @endcode
3612 *
3613 *
3614 * <a name="step_44-Solidmake_constraints"></a>
3615 * <h4>Solid::make_constraints</h4>
3616 * The constraints for this problem are simple to describe.
3617 * In this particular example, the boundary values will be calculated for
3618 * the two first iterations of Newton's algorithm. In general, one would
3619 * build non-homogeneous constraints in the zeroth iteration (that is, when
3620 * `apply_dirichlet_bc == true` in the code block that follows) and build
3621 * only the corresponding homogeneous constraints in the following step. While
3622 * the current example has only homogeneous constraints, previous experiences
3623 * have shown that a common error is forgetting to add the extra condition
3624 * when refactoring the code to specific uses. This could lead to errors that
3625 * are hard to debug. In this spirit, we choose to make the code more verbose
3626 * in terms of what operations are performed at each Newton step.
3627 *
3628 * @code
3629 *   template <int dim>
3630 *   void Solid<dim>::make_constraints(const int it_nr)
3631 *   {
3632 * @endcode
3633 *
3634 * Since we (a) are dealing with an iterative Newton method, (b) are using
3635 * an incremental formulation for the displacement, and (c) apply the
3636 * constraints to the incremental displacement field, any non-homogeneous
3637 * constraints on the displacement update should only be specified at the
3638 * zeroth iteration. No subsequent contributions are to be made since the
3639 * constraints will be exactly satisfied after that iteration.
3640 *
3641 * @code
3642 *   const bool apply_dirichlet_bc = (it_nr == 0);
3643 *  
3644 * @endcode
3645 *
3646 * Furthermore, after the first Newton iteration within a timestep, the
3647 * constraints remain the same and we do not need to modify or rebuild them
3648 * so long as we do not clear the @p constraints object.
3649 *
3650 * @code
3651 *   if (it_nr > 1)
3652 *   {
3653 *   std::cout << " --- " << std::flush;
3654 *   return;
3655 *   }
3656 *  
3657 *   std::cout << " CST " << std::flush;
3658 *  
3659 *   if (apply_dirichlet_bc)
3660 *   {
3661 * @endcode
3662 *
3663 * At the zeroth Newton iteration we wish to apply the full set of
3664 * non-homogeneous and homogeneous constraints that represent the
3665 * boundary conditions on the displacement increment. Since in general
3666 * the constraints may be different at each time step, we need to clear
3667 * the constraints matrix and completely rebuild it. An example case
3668 * would be if a surface is accelerating; in such a scenario the change
3669 * in displacement is non-constant between each time step.
3670 *
3671 * @code
3672 *   constraints.clear();
3673 *  
3674 * @endcode
3675 *
3676 * The boundary conditions for the indentation problem in 3d are as
3677 * follows: On the -x, -y and -z faces (IDs 0,2,4) we set up a symmetry
3678 * condition to allow only planar movement while the +x and +z faces
3679 * (IDs 1,5) are traction free. In this contrived problem, part of the
3680 * +y face (ID 3) is set to have no motion in the x- and z-component.
3681 * Finally, as described earlier, the other part of the +y face has an
3682 * the applied pressure but is also constrained in the x- and
3683 * z-directions.
3684 *
3685
3686 *
3687 * In the following, we will have to tell the function interpolation
3688 * boundary values which components of the solution vector should be
3689 * constrained (i.e., whether it's the x-, y-, z-displacements or
3690 * combinations thereof). This is done using ComponentMask objects (see
3691 * @ref GlossComponentMask) which we can get from the finite element if we
3692 * provide it with an extractor object for the component we wish to
3693 * select. To this end we first set up such extractor objects and later
3694 * use it when generating the relevant component masks:
3695 *
3696 * @code
3697 *   const FEValuesExtractors::Scalar x_displacement(0);
3698 *   const FEValuesExtractors::Scalar y_displacement(1);
3699 *  
3700 *   {
3701 *   const int boundary_id = 0;
3702 *  
3704 *   dof_handler,
3705 *   boundary_id,
3706 *   Functions::ZeroFunction<dim>(n_components),
3707 *   constraints,
3708 *   fe.component_mask(x_displacement));
3709 *   }
3710 *   {
3711 *   const int boundary_id = 2;
3712 *  
3714 *   dof_handler,
3715 *   boundary_id,
3716 *   Functions::ZeroFunction<dim>(n_components),
3717 *   constraints,
3718 *   fe.component_mask(y_displacement));
3719 *   }
3720 *  
3721 *   if (dim == 3)
3722 *   {
3723 *   const FEValuesExtractors::Scalar z_displacement(2);
3724 *  
3725 *   {
3726 *   const int boundary_id = 3;
3727 *  
3729 *   dof_handler,
3730 *   boundary_id,
3731 *   Functions::ZeroFunction<dim>(n_components),
3732 *   constraints,
3733 *   (fe.component_mask(x_displacement) |
3734 *   fe.component_mask(z_displacement)));
3735 *   }
3736 *   {
3737 *   const int boundary_id = 4;
3738 *  
3740 *   dof_handler,
3741 *   boundary_id,
3742 *   Functions::ZeroFunction<dim>(n_components),
3743 *   constraints,
3744 *   fe.component_mask(z_displacement));
3745 *   }
3746 *  
3747 *   {
3748 *   const int boundary_id = 6;
3749 *  
3751 *   dof_handler,
3752 *   boundary_id,
3753 *   Functions::ZeroFunction<dim>(n_components),
3754 *   constraints,
3755 *   (fe.component_mask(x_displacement) |
3756 *   fe.component_mask(z_displacement)));
3757 *   }
3758 *   }
3759 *   else
3760 *   {
3761 *   {
3762 *   const int boundary_id = 3;
3763 *  
3765 *   dof_handler,
3766 *   boundary_id,
3767 *   Functions::ZeroFunction<dim>(n_components),
3768 *   constraints,
3769 *   (fe.component_mask(x_displacement)));
3770 *   }
3771 *   {
3772 *   const int boundary_id = 6;
3773 *  
3775 *   dof_handler,
3776 *   boundary_id,
3777 *   Functions::ZeroFunction<dim>(n_components),
3778 *   constraints,
3779 *   (fe.component_mask(x_displacement)));
3780 *   }
3781 *   }
3782 *   }
3783 *   else
3784 *   {
3785 * @endcode
3786 *
3787 * As all Dirichlet constraints are fulfilled exactly after the zeroth
3788 * Newton iteration, we want to ensure that no further modification are
3789 * made to those entries. This implies that we want to convert
3790 * all non-homogeneous Dirichlet constraints into homogeneous ones.
3791 *
3792
3793 *
3794 * In this example the procedure to do this is quite straightforward,
3795 * and in fact we can (and will) circumvent any unnecessary operations
3796 * when only homogeneous boundary conditions are applied.
3797 * In a more general problem one should be mindful of hanging node
3798 * and periodic constraints, which may also introduce some
3799 * inhomogeneities. It might then be advantageous to keep disparate
3800 * objects for the different types of constraints, and merge them
3801 * together once the homogeneous Dirichlet constraints have been
3802 * constructed.
3803 *
3804 * @code
3805 *   if (constraints.has_inhomogeneities())
3806 *   {
3807 * @endcode
3808 *
3809 * Since the affine constraints were finalized at the previous
3810 * Newton iteration, they may not be modified directly. So
3811 * we need to copy them to another temporary object and make
3812 * modification there. Once we're done, we'll transfer them
3813 * back to the main @p constraints object.
3814 *
3815 * @code
3816 *   AffineConstraints<double> homogeneous_constraints(constraints);
3817 *   for (unsigned int dof = 0; dof != dof_handler.n_dofs(); ++dof)
3818 *   if (homogeneous_constraints.is_inhomogeneously_constrained(dof))
3819 *   homogeneous_constraints.set_inhomogeneity(dof, 0.0);
3820 *  
3821 *   constraints.clear();
3822 *   constraints.copy_from(homogeneous_constraints);
3823 *   }
3824 *   }
3825 *  
3826 *   constraints.close();
3827 *   }
3828 *  
3829 * @endcode
3830 *
3831 *
3832 * <a name="step_44-Solidassemble_sc"></a>
3833 * <h4>Solid::assemble_sc</h4>
3834 * Solving the entire block system is a bit problematic as there are no
3835 * contributions to the @f$\mathsf{\mathbf{K}}_{ \widetilde{J} \widetilde{J}}@f$
3836 * block, rendering it noninvertible (when using an iterative solver).
3837 * Since the pressure and dilatation variables DOFs are discontinuous, we can
3838 * condense them out to form a smaller displacement-only system which
3839 * we will then solve and subsequently post-process to retrieve the
3840 * pressure and dilatation solutions.
3841 *
3842
3843 *
3844 * The static condensation process could be performed at a global level but we
3845 * need the inverse of one of the blocks. However, since the pressure and
3846 * dilatation variables are discontinuous, the static condensation (SC)
3847 * operation can also be done on a per-cell basis and we can produce the
3848 * inverse of the block-diagonal
3849 * @f$\mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}@f$ block by inverting the
3850 * local blocks. We can again use TBB to do this since each operation will be
3851 * independent of one another.
3852 *
3853
3854 *
3855 * Using the TBB via the WorkStream class, we assemble the contributions to
3856 * form
3857 * @f$
3858 * \mathsf{\mathbf{K}}_{\textrm{con}}
3859 * = \bigl[ \mathsf{\mathbf{K}}_{uu} +
3860 * \overline{\overline{\mathsf{\mathbf{K}}}}~ \bigr]
3861 * @f$
3862 * from each element's contributions. These contributions are then added to
3863 * the global stiffness matrix. Given this description, the following two
3864 * functions should be clear:
3865 *
3866 * @code
3867 *   template <int dim>
3868 *   void Solid<dim>::assemble_sc()
3869 *   {
3870 *   timer.enter_subsection("Perform static condensation");
3871 *   std::cout << " ASM_SC " << std::flush;
3872 *  
3873 *   PerTaskData_SC per_task_data(dofs_per_cell,
3874 *   element_indices_u.size(),
3875 *   element_indices_p.size(),
3876 *   element_indices_J.size());
3877 *   ScratchData_SC scratch_data;
3878 *  
3879 *   WorkStream::run(dof_handler.active_cell_iterators(),
3880 *   *this,
3881 *   &Solid::assemble_sc_one_cell,
3882 *   &Solid::copy_local_to_global_sc,
3883 *   scratch_data,
3884 *   per_task_data);
3885 *  
3886 *   timer.leave_subsection();
3887 *   }
3888 *  
3889 *  
3890 *   template <int dim>
3891 *   void Solid<dim>::copy_local_to_global_sc(const PerTaskData_SC &data)
3892 *   {
3893 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3894 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
3895 *   tangent_matrix.add(data.local_dof_indices[i],
3896 *   data.local_dof_indices[j],
3897 *   data.cell_matrix(i, j));
3898 *   }
3899 *  
3900 *  
3901 * @endcode
3902 *
3903 * Now we describe the static condensation process. As per usual, we must
3904 * first find out which global numbers the degrees of freedom on this cell
3905 * have and reset some data structures:
3906 *
3907 * @code
3908 *   template <int dim>
3909 *   void Solid<dim>::assemble_sc_one_cell(
3910 *   const typename DoFHandler<dim>::active_cell_iterator &cell,
3911 *   ScratchData_SC &scratch,
3912 *   PerTaskData_SC &data)
3913 *   {
3914 *   data.reset();
3915 *   scratch.reset();
3916 *   cell->get_dof_indices(data.local_dof_indices);
3917 *  
3918 * @endcode
3919 *
3920 * We now extract the contribution of the dofs associated with the current
3921 * cell to the global stiffness matrix. The discontinuous nature of the
3922 * @f$\widetilde{p}@f$ and @f$\widetilde{J}@f$ interpolations mean that their is
3923 * no coupling of the local contributions at the global level. This is not
3924 * the case with the @f$\mathbf{u}@f$ dof. In other words,
3925 * @f$\mathsf{\mathbf{k}}_{\widetilde{J} \widetilde{p}}@f$,
3926 * @f$\mathsf{\mathbf{k}}_{\widetilde{p} \widetilde{p}}@f$ and
3927 * @f$\mathsf{\mathbf{k}}_{\widetilde{J} \widetilde{p}}@f$, when extracted
3928 * from the global stiffness matrix are the element contributions. This
3929 * is not the case for @f$\mathsf{\mathbf{k}}_{uu}@f$.
3930 *
3931
3932 *
3933 * Note: A lower-case symbol is used to denote element stiffness matrices.
3934 *
3935
3936 *
3937 * Currently the matrix corresponding to
3938 * the dof associated with the current element
3939 * (denoted somewhat loosely as @f$\mathsf{\mathbf{k}}@f$)
3940 * is of the form:
3941 * @f{align*}{
3942 * \begin{bmatrix}
3943 * \mathsf{\mathbf{k}}_{uu} & \mathsf{\mathbf{k}}_{u\widetilde{p}}
3944 * & \mathbf{0}
3945 * \\ \mathsf{\mathbf{k}}_{\widetilde{p}u} & \mathbf{0} &
3946 * \mathsf{\mathbf{k}}_{\widetilde{p}\widetilde{J}}
3947 * \\ \mathbf{0} & \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{p}} &
3948 * \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{J}} \end{bmatrix}
3949 * @f}
3950 *
3951
3952 *
3953 * We now need to modify it such that it appear as
3954 * @f{align*}{
3955 * \begin{bmatrix}
3956 * \mathsf{\mathbf{k}}_{\textrm{con}} &
3957 * \mathsf{\mathbf{k}}_{u\widetilde{p}} & \mathbf{0}
3958 * \\ \mathsf{\mathbf{k}}_{\widetilde{p}u} & \mathbf{0} &
3959 * \mathsf{\mathbf{k}}_{\widetilde{p}\widetilde{J}}^{-1}
3960 * \\ \mathbf{0} & \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{p}} &
3961 * \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{J}} \end{bmatrix}
3962 * @f}
3963 * with @f$\mathsf{\mathbf{k}}_{\textrm{con}} = \bigl[
3964 * \mathsf{\mathbf{k}}_{uu} +\overline{\overline{\mathsf{\mathbf{k}}}}~
3965 * \bigr]@f$ where @f$ \overline{\overline{\mathsf{\mathbf{k}}}}
3966 * \dealcoloneq \mathsf{\mathbf{k}}_{u\widetilde{p}}
3967 * \overline{\mathsf{\mathbf{k}}} \mathsf{\mathbf{k}}_{\widetilde{p}u}
3968 * @f$
3969 * and
3970 * @f$
3971 * \overline{\mathsf{\mathbf{k}}} =
3972 * \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{p}}^{-1}
3973 * \mathsf{\mathbf{k}}_{\widetilde{J}\widetilde{J}}
3974 * \mathsf{\mathbf{k}}_{\widetilde{p}\widetilde{J}}^{-1}
3975 * @f$.
3976 *
3977
3978 *
3979 * At this point, we need to take note of
3980 * the fact that global data already exists
3981 * in the @f$\mathsf{\mathbf{K}}_{uu}@f$,
3982 * @f$\mathsf{\mathbf{K}}_{\widetilde{p} \widetilde{J}}@f$
3983 * and
3984 * @f$\mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{p}}@f$
3985 * sub-blocks. So if we are to modify them, we must account for the data
3986 * that is already there (i.e. simply add to it or remove it if
3987 * necessary). Since the copy_local_to_global operation is a "+="
3988 * operation, we need to take this into account
3989 *
3990
3991 *
3992 * For the @f$\mathsf{\mathbf{K}}_{uu}@f$ block in particular, this means that
3993 * contributions have been added from the surrounding cells, so we need to
3994 * be careful when we manipulate this block. We can't just erase the
3995 * sub-blocks.
3996 *
3997
3998 *
3999 * This is the strategy we will employ to get the sub-blocks we want:
4000 *
4001
4002 *
4003 * - @f$ {\mathsf{\mathbf{k}}}_{\textrm{store}}@f$:
4004 * Since we don't have access to @f$\mathsf{\mathbf{k}}_{uu}@f$,
4005 * but we know its contribution is added to
4006 * the global @f$\mathsf{\mathbf{K}}_{uu}@f$ matrix, we just want
4007 * to add the element wise
4008 * static-condensation @f$\overline{\overline{\mathsf{\mathbf{k}}}}@f$.
4009 *
4010
4011 *
4012 * - @f$\mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}@f$:
4013 * Similarly, @f$\mathsf{\mathbf{k}}_{\widetilde{p}
4014 * \widetilde{J}}@f$ exists in
4015 * the subblock. Since the copy
4016 * operation is a += operation, we
4017 * need to subtract the existing
4018 * @f$\mathsf{\mathbf{k}}_{\widetilde{p} \widetilde{J}}@f$
4019 * submatrix in addition to
4020 * "adding" that which we wish to
4021 * replace it with.
4022 *
4023
4024 *
4025 * - @f$\mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{p}}@f$:
4026 * Since the global matrix
4027 * is symmetric, this block is the
4028 * same as the one above and we
4029 * can simply use
4030 * @f$\mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}@f$
4031 * as a substitute for this one.
4032 *
4033
4034 *
4035 * We first extract element data from the
4036 * system matrix. So first we get the
4037 * entire subblock for the cell, then
4038 * extract @f$\mathsf{\mathbf{k}}@f$
4039 * for the dofs associated with
4040 * the current element
4041 *
4042 * @code
4043 *   data.k_orig.extract_submatrix_from(tangent_matrix,
4044 *   data.local_dof_indices,
4045 *   data.local_dof_indices);
4046 * @endcode
4047 *
4048 * and next the local matrices for
4049 * @f$\mathsf{\mathbf{k}}_{ \widetilde{p} u}@f$
4050 * @f$\mathsf{\mathbf{k}}_{ \widetilde{p} \widetilde{J}}@f$
4051 * and
4052 * @f$\mathsf{\mathbf{k}}_{ \widetilde{J} \widetilde{J}}@f$:
4053 *
4054 * @code
4055 *   data.k_pu.extract_submatrix_from(data.k_orig,
4056 *   element_indices_p,
4057 *   element_indices_u);
4058 *   data.k_pJ.extract_submatrix_from(data.k_orig,
4059 *   element_indices_p,
4060 *   element_indices_J);
4061 *   data.k_JJ.extract_submatrix_from(data.k_orig,
4062 *   element_indices_J,
4063 *   element_indices_J);
4064 *  
4065 * @endcode
4066 *
4067 * To get the inverse of @f$\mathsf{\mathbf{k}}_{\widetilde{p}
4068 * \widetilde{J}}@f$, we invert it directly. This operation is relatively
4069 * inexpensive since @f$\mathsf{\mathbf{k}}_{\widetilde{p} \widetilde{J}}@f$
4070 * since block-diagonal.
4071 *
4072 * @code
4073 *   data.k_pJ_inv.invert(data.k_pJ);
4074 *  
4075 * @endcode
4076 *
4077 * Now we can make condensation terms to
4078 * add to the @f$\mathsf{\mathbf{k}}_{uu}@f$
4079 * block and put them in
4080 * the cell local matrix
4081 * @f$
4082 * \mathsf{\mathbf{A}}
4083 * =
4084 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}
4085 * \mathsf{\mathbf{k}}_{\widetilde{p} u}
4086 * @f$:
4087 *
4088 * @code
4089 *   data.k_pJ_inv.mmult(data.A, data.k_pu);
4090 * @endcode
4091 *
4092 * @f$
4093 * \mathsf{\mathbf{B}}
4094 * =
4095 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{J}}
4096 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}
4097 * \mathsf{\mathbf{k}}_{\widetilde{p} u}
4098 * @f$
4099 *
4100 * @code
4101 *   data.k_JJ.mmult(data.B, data.A);
4102 * @endcode
4103 *
4104 * @f$
4105 * \mathsf{\mathbf{C}}
4106 * =
4107 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{p}}
4108 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{J}}
4109 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}
4110 * \mathsf{\mathbf{k}}_{\widetilde{p} u}
4111 * @f$
4112 *
4113 * @code
4114 *   data.k_pJ_inv.Tmmult(data.C, data.B);
4115 * @endcode
4116 *
4117 * @f$
4118 * \overline{\overline{\mathsf{\mathbf{k}}}}
4119 * =
4120 * \mathsf{\mathbf{k}}_{u \widetilde{p}}
4121 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{p}}
4122 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{J} \widetilde{J}}
4123 * \mathsf{\mathbf{k}}^{-1}_{\widetilde{p} \widetilde{J}}
4124 * \mathsf{\mathbf{k}}_{\widetilde{p} u}
4125 * @f$
4126 *
4127 * @code
4128 *   data.k_pu.Tmmult(data.k_bbar, data.C);
4129 *   data.k_bbar.scatter_matrix_to(element_indices_u,
4130 *   element_indices_u,
4131 *   data.cell_matrix);
4132 *  
4133 * @endcode
4134 *
4135 * Next we place
4136 * @f$\mathsf{\mathbf{k}}^{-1}_{ \widetilde{p} \widetilde{J}}@f$
4137 * in the
4138 * @f$\mathsf{\mathbf{k}}_{ \widetilde{p} \widetilde{J}}@f$
4139 * block for post-processing. Note again
4140 * that we need to remove the
4141 * contribution that already exists there.
4142 *
4143 * @code
4144 *   data.k_pJ_inv.add(-1.0, data.k_pJ);
4145 *   data.k_pJ_inv.scatter_matrix_to(element_indices_p,
4146 *   element_indices_J,
4147 *   data.cell_matrix);
4148 *   }
4149 *  
4150 * @endcode
4151 *
4152 *
4153 * <a name="step_44-Solidsolve_linear_system"></a>
4154 * <h4>Solid::solve_linear_system</h4>
4155 * We now have all of the necessary components to use one of two possible
4156 * methods to solve the linearised system. The first is to perform static
4157 * condensation on an element level, which requires some alterations
4158 * to the tangent matrix and RHS vector. Alternatively, the full block
4159 * system can be solved by performing condensation on a global level.
4160 * Below we implement both approaches.
4161 *
4162 * @code
4163 *   template <int dim>
4164 *   std::pair<unsigned int, double>
4165 *   Solid<dim>::solve_linear_system(BlockVector<double> &newton_update)
4166 *   {
4167 *   unsigned int lin_it = 0;
4168 *   double lin_res = 0.0;
4169 *  
4170 *   if (parameters.use_static_condensation == true)
4171 *   {
4172 * @endcode
4173 *
4174 * Firstly, here is the approach using the (permanent) augmentation of
4175 * the tangent matrix. For the following, recall that
4176 * @f{align*}{
4177 * \mathsf{\mathbf{K}}_{\textrm{store}}
4178 * \dealcoloneq
4179 * \begin{bmatrix}
4180 * \mathsf{\mathbf{K}}_{\textrm{con}} &
4181 * \mathsf{\mathbf{K}}_{u\widetilde{p}} & \mathbf{0}
4182 * \\ \mathsf{\mathbf{K}}_{\widetilde{p}u} & \mathbf{0} &
4183 * \mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}^{-1}
4184 * \\ \mathbf{0} &
4185 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}} &
4186 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}} \end{bmatrix} \, .
4187 * @f}
4188 * and
4189 * @f{align*}{
4190 * d \widetilde{\mathsf{\mathbf{p}}}
4191 * & =
4192 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4193 * \bigl[
4194 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4195 * -
4196 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4197 * d \widetilde{\mathsf{\mathbf{J}}} \bigr]
4198 * \\ d \widetilde{\mathsf{\mathbf{J}}}
4199 * & =
4200 * \mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}^{-1}
4201 * \bigl[
4202 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4203 * - \mathsf{\mathbf{K}}_{\widetilde{p}u} d
4204 * \mathsf{\mathbf{u}} \bigr]
4205 * \\ \Rightarrow d \widetilde{\mathsf{\mathbf{p}}}
4206 * &= \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4207 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4208 * -
4209 * \underbrace{\bigl[\mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4210 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4211 * \mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}^{-1}\bigr]}_{\overline{\mathsf{\mathbf{K}}}}\bigl[
4212 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4213 * - \mathsf{\mathbf{K}}_{\widetilde{p}u} d
4214 * \mathsf{\mathbf{u}} \bigr]
4215 * @f}
4216 * and thus
4217 * @f[
4218 * \underbrace{\bigl[ \mathsf{\mathbf{K}}_{uu} +
4219 * \overline{\overline{\mathsf{\mathbf{K}}}}~ \bigr]
4220 * }_{\mathsf{\mathbf{K}}_{\textrm{con}}} d
4221 * \mathsf{\mathbf{u}}
4222 * =
4223 * \underbrace{
4224 * \Bigl[
4225 * \mathsf{\mathbf{F}}_{u}
4226 * - \mathsf{\mathbf{K}}_{u\widetilde{p}} \bigl[
4227 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4228 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4229 * -
4230 * \overline{\mathsf{\mathbf{K}}}\mathsf{\mathbf{F}}_{\widetilde{p}}
4231 * \bigr]
4232 * \Bigr]}_{\mathsf{\mathbf{F}}_{\textrm{con}}}
4233 * @f]
4234 * where
4235 * @f[
4236 * \overline{\overline{\mathsf{\mathbf{K}}}} \dealcoloneq
4237 * \mathsf{\mathbf{K}}_{u\widetilde{p}}
4238 * \overline{\mathsf{\mathbf{K}}}
4239 * \mathsf{\mathbf{K}}_{\widetilde{p}u} \, .
4240 * @f]
4241 *
4242
4243 *
4244 * At the top, we allocate two temporary vectors to help with the
4245 * static condensation, and variables to store the number of
4246 * linear solver iterations and the (hopefully converged) residual.
4247 *
4248 * @code
4249 *   BlockVector<double> A(dofs_per_block);
4250 *   BlockVector<double> B(dofs_per_block);
4251 *  
4252 *  
4253 * @endcode
4254 *
4255 * In the first step of this function, we solve for the incremental
4256 * displacement @f$d\mathbf{u}@f$. To this end, we perform static
4257 * condensation to make
4258 * @f$\mathsf{\mathbf{K}}_{\textrm{con}}
4259 * = \bigl[ \mathsf{\mathbf{K}}_{uu} +
4260 * \overline{\overline{\mathsf{\mathbf{K}}}}~ \bigr]@f$
4261 * and put
4262 * @f$\mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}@f$
4263 * in the original @f$\mathsf{\mathbf{K}}_{\widetilde{p} \widetilde{J}}@f$
4264 * block. That is, we make @f$\mathsf{\mathbf{K}}_{\textrm{store}}@f$.
4265 *
4266 * @code
4267 *   {
4268 *   assemble_sc();
4269 *  
4270 * @endcode
4271 *
4272 * @f$
4273 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4274 * =
4275 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4276 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4277 * @f$
4278 *
4279 * @code
4280 *   tangent_matrix.block(p_dof, J_dof)
4281 *   .vmult(A.block(J_dof), system_rhs.block(p_dof));
4282 * @endcode
4283 *
4284 * @f$
4285 * \mathsf{\mathbf{B}}_{\widetilde{J}}
4286 * =
4287 * \mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{J}}
4288 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4289 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4290 * @f$
4291 *
4292 * @code
4293 *   tangent_matrix.block(J_dof, J_dof)
4294 *   .vmult(B.block(J_dof), A.block(J_dof));
4295 * @endcode
4296 *
4297 * @f$
4298 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4299 * =
4300 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4301 * -
4302 * \mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{J}}
4303 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4304 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4305 * @f$
4306 *
4307 * @code
4308 *   A.block(J_dof) = system_rhs.block(J_dof);
4309 *   A.block(J_dof) -= B.block(J_dof);
4310 * @endcode
4311 *
4312 * @f$
4313 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4314 * =
4315 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{J} \widetilde{p}}
4316 * [
4317 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4318 * -
4319 * \mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{J}}
4320 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4321 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4322 * ]
4323 * @f$
4324 *
4325 * @code
4326 *   tangent_matrix.block(p_dof, J_dof)
4327 *   .Tvmult(A.block(p_dof), A.block(J_dof));
4328 * @endcode
4329 *
4330 * @f$
4331 * \mathsf{\mathbf{A}}_{u}
4332 * =
4333 * \mathsf{\mathbf{K}}_{u \widetilde{p}}
4334 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{J} \widetilde{p}}
4335 * [
4336 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4337 * -
4338 * \mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{J}}
4339 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4340 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4341 * ]
4342 * @f$
4343 *
4344 * @code
4345 *   tangent_matrix.block(u_dof, p_dof)
4346 *   .vmult(A.block(u_dof), A.block(p_dof));
4347 * @endcode
4348 *
4349 * @f$
4350 * \mathsf{\mathbf{F}}_{\text{con}}
4351 * =
4352 * \mathsf{\mathbf{F}}_{u}
4353 * -
4354 * \mathsf{\mathbf{K}}_{u \widetilde{p}}
4355 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{J} \widetilde{p}}
4356 * [
4357 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4358 * -
4359 * \mathsf{\mathbf{K}}_{\widetilde{J} \widetilde{J}}
4360 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p} \widetilde{J}}
4361 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4362 * ]
4363 * @f$
4364 *
4365 * @code
4366 *   system_rhs.block(u_dof) -= A.block(u_dof);
4367 *  
4368 *   timer.enter_subsection("Linear solver");
4369 *   std::cout << " SLV " << std::flush;
4370 *   if (parameters.type_lin == "CG")
4371 *   {
4372 *   const auto solver_its = static_cast<unsigned int>(
4373 *   tangent_matrix.block(u_dof, u_dof).m() *
4374 *   parameters.max_iterations_lin);
4375 *   const double tol_sol =
4376 *   parameters.tol_lin * system_rhs.block(u_dof).l2_norm();
4377 *  
4378 *   SolverControl solver_control(solver_its, tol_sol);
4379 *  
4380 *   GrowingVectorMemory<Vector<double>> GVM;
4381 *   SolverCG<Vector<double>> solver_CG(solver_control, GVM);
4382 *  
4383 * @endcode
4384 *
4385 * We've chosen by default a SSOR preconditioner as it appears to
4386 * provide the fastest solver convergence characteristics for this
4387 * problem on a single-thread machine. However, this might not be
4388 * true for different problem sizes.
4389 *
4390 * @code
4392 *   preconditioner(parameters.preconditioner_type,
4393 *   parameters.preconditioner_relaxation);
4394 *   preconditioner.use_matrix(tangent_matrix.block(u_dof, u_dof));
4395 *  
4396 *   solver_CG.solve(tangent_matrix.block(u_dof, u_dof),
4397 *   newton_update.block(u_dof),
4398 *   system_rhs.block(u_dof),
4399 *   preconditioner);
4400 *  
4401 *   lin_it = solver_control.last_step();
4402 *   lin_res = solver_control.last_value();
4403 *   }
4404 *   else if (parameters.type_lin == "Direct")
4405 *   {
4406 * @endcode
4407 *
4408 * Otherwise if the problem is small
4409 * enough, a direct solver can be
4410 * utilised.
4411 *
4412 * @code
4413 *   SparseDirectUMFPACK A_direct;
4414 *   A_direct.initialize(tangent_matrix.block(u_dof, u_dof));
4415 *   A_direct.vmult(newton_update.block(u_dof),
4416 *   system_rhs.block(u_dof));
4417 *  
4418 *   lin_it = 1;
4419 *   lin_res = 0.0;
4420 *   }
4421 *   else
4422 *   Assert(false, ExcMessage("Linear solver type not implemented"));
4423 *  
4424 *   timer.leave_subsection();
4425 *   }
4426 *  
4427 * @endcode
4428 *
4429 * Now that we have the displacement update, distribute the constraints
4430 * back to the Newton update:
4431 *
4432 * @code
4433 *   constraints.distribute(newton_update);
4434 *  
4435 *   timer.enter_subsection("Linear solver postprocessing");
4436 *   std::cout << " PP " << std::flush;
4437 *  
4438 * @endcode
4439 *
4440 * The next step after solving the displacement
4441 * problem is to post-process to get the
4442 * dilatation solution from the
4443 * substitution:
4444 * @f$
4445 * d \widetilde{\mathsf{\mathbf{J}}}
4446 * = \mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}^{-1} \bigl[
4447 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4448 * - \mathsf{\mathbf{K}}_{\widetilde{p}u} d \mathsf{\mathbf{u}}
4449 * \bigr]
4450 * @f$
4451 *
4452 * @code
4453 *   {
4454 * @endcode
4455 *
4456 * @f$
4457 * \mathsf{\mathbf{A}}_{\widetilde{p}}
4458 * =
4459 * \mathsf{\mathbf{K}}_{\widetilde{p}u} d \mathsf{\mathbf{u}}
4460 * @f$
4461 *
4462 * @code
4463 *   tangent_matrix.block(p_dof, u_dof)
4464 *   .vmult(A.block(p_dof), newton_update.block(u_dof));
4465 * @endcode
4466 *
4467 * @f$
4468 * \mathsf{\mathbf{A}}_{\widetilde{p}}
4469 * =
4470 * -\mathsf{\mathbf{K}}_{\widetilde{p}u} d \mathsf{\mathbf{u}}
4471 * @f$
4472 *
4473 * @code
4474 *   A.block(p_dof) *= -1.0;
4475 * @endcode
4476 *
4477 * @f$
4478 * \mathsf{\mathbf{A}}_{\widetilde{p}}
4479 * =
4480 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4481 * -\mathsf{\mathbf{K}}_{\widetilde{p}u} d \mathsf{\mathbf{u}}
4482 * @f$
4483 *
4484 * @code
4485 *   A.block(p_dof) += system_rhs.block(p_dof);
4486 * @endcode
4487 *
4488 * @f$
4489 * d\mathsf{\mathbf{\widetilde{J}}}
4490 * =
4491 * \mathsf{\mathbf{K}}^{-1}_{\widetilde{p}\widetilde{J}}
4492 * [
4493 * \mathsf{\mathbf{F}}_{\widetilde{p}}
4494 * -\mathsf{\mathbf{K}}_{\widetilde{p}u} d \mathsf{\mathbf{u}}
4495 * ]
4496 * @f$
4497 *
4498 * @code
4499 *   tangent_matrix.block(p_dof, J_dof)
4500 *   .vmult(newton_update.block(J_dof), A.block(p_dof));
4501 *   }
4502 *  
4503 * @endcode
4504 *
4505 * we ensure here that any Dirichlet constraints
4506 * are distributed on the updated solution:
4507 *
4508 * @code
4509 *   constraints.distribute(newton_update);
4510 *  
4511 * @endcode
4512 *
4513 * Finally we solve for the pressure
4514 * update with the substitution:
4515 * @f$
4516 * d \widetilde{\mathsf{\mathbf{p}}}
4517 * =
4518 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4519 * \bigl[
4520 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4521 * - \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4522 * d \widetilde{\mathsf{\mathbf{J}}}
4523 * \bigr]
4524 * @f$
4525 *
4526 * @code
4527 *   {
4528 * @endcode
4529 *
4530 * @f$
4531 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4532 * =
4533 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4534 * d \widetilde{\mathsf{\mathbf{J}}}
4535 * @f$
4536 *
4537 * @code
4538 *   tangent_matrix.block(J_dof, J_dof)
4539 *   .vmult(A.block(J_dof), newton_update.block(J_dof));
4540 * @endcode
4541 *
4542 * @f$
4543 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4544 * =
4545 * -\mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4546 * d \widetilde{\mathsf{\mathbf{J}}}
4547 * @f$
4548 *
4549 * @code
4550 *   A.block(J_dof) *= -1.0;
4551 * @endcode
4552 *
4553 * @f$
4554 * \mathsf{\mathbf{A}}_{\widetilde{J}}
4555 * =
4556 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4557 * -
4558 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4559 * d \widetilde{\mathsf{\mathbf{J}}}
4560 * @f$
4561 *
4562 * @code
4563 *   A.block(J_dof) += system_rhs.block(J_dof);
4564 * @endcode
4565 *
4566 * and finally....
4567 * @f$
4568 * d \widetilde{\mathsf{\mathbf{p}}}
4569 * =
4570 * \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}
4571 * \bigl[
4572 * \mathsf{\mathbf{F}}_{\widetilde{J}}
4573 * - \mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{J}}
4574 * d \widetilde{\mathsf{\mathbf{J}}}
4575 * \bigr]
4576 * @f$
4577 *
4578 * @code
4579 *   tangent_matrix.block(p_dof, J_dof)
4580 *   .Tvmult(newton_update.block(p_dof), A.block(J_dof));
4581 *   }
4582 *  
4583 * @endcode
4584 *
4585 * We are now at the end, so we distribute all
4586 * constrained dofs back to the Newton
4587 * update:
4588 *
4589 * @code
4590 *   constraints.distribute(newton_update);
4591 *  
4592 *   timer.leave_subsection();
4593 *   }
4594 *   else
4595 *   {
4596 *   std::cout << " ------ " << std::flush;
4597 *  
4598 *   timer.enter_subsection("Linear solver");
4599 *   std::cout << " SLV " << std::flush;
4600 *  
4601 *   if (parameters.type_lin == "CG")
4602 *   {
4603 * @endcode
4604 *
4605 * Manual condensation of the dilatation and pressure fields on
4606 * a local level, and subsequent post-processing, took quite a
4607 * bit of effort to achieve. To recap, we had to produce the
4608 * inverse matrix
4609 * @f$\mathsf{\mathbf{K}}_{\widetilde{p}\widetilde{J}}^{-1}@f$, which
4610 * was permanently written into the global tangent matrix. We then
4611 * permanently modified @f$\mathsf{\mathbf{K}}_{uu}@f$ to produce
4612 * @f$\mathsf{\mathbf{K}}_{\textrm{con}}@f$. This involved the
4613 * extraction and manipulation of local sub-blocks of the tangent
4614 * matrix. After solving for the displacement, the individual
4615 * matrix-vector operations required to solve for dilatation and
4616 * pressure were carefully implemented. Contrast these many sequence
4617 * of steps to the much simpler and transparent implementation using
4618 * functionality provided by the LinearOperator class.
4619 *
4620
4621 *
4622 * For ease of later use, we define some aliases for
4623 * blocks in the RHS vector
4624 *
4625 * @code
4626 *   const Vector<double> &f_u = system_rhs.block(u_dof);
4627 *   const Vector<double> &f_p = system_rhs.block(p_dof);
4628 *   const Vector<double> &f_J = system_rhs.block(J_dof);
4629 *  
4630 * @endcode
4631 *
4632 * ... and for blocks in the Newton update vector.
4633 *
4634 * @code
4635 *   Vector<double> &d_u = newton_update.block(u_dof);
4636 *   Vector<double> &d_p = newton_update.block(p_dof);
4637 *   Vector<double> &d_J = newton_update.block(J_dof);
4638 *  
4639 * @endcode
4640 *
4641 * We next define some linear operators for the tangent matrix
4642 * sub-blocks We will exploit the symmetry of the system, so not all
4643 * blocks are required.
4644 *
4645 * @code
4646 *   const auto K_uu =
4647 *   linear_operator(tangent_matrix.block(u_dof, u_dof));
4648 *   const auto K_up =
4649 *   linear_operator(tangent_matrix.block(u_dof, p_dof));
4650 *   const auto K_pu =
4651 *   linear_operator(tangent_matrix.block(p_dof, u_dof));
4652 *   const auto K_Jp =
4653 *   linear_operator(tangent_matrix.block(J_dof, p_dof));
4654 *   const auto K_JJ =
4655 *   linear_operator(tangent_matrix.block(J_dof, J_dof));
4656 *  
4657 * @endcode
4658 *
4659 * We then construct a LinearOperator that represents the inverse of
4660 * (square block)
4661 * @f$\mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}@f$. Since it is
4662 * diagonal (or, when a higher order ansatz it used, nearly
4663 * diagonal), a Jacobi preconditioner is suitable.
4664 *
4665 * @code
4667 *   preconditioner_K_Jp_inv("jacobi");
4668 *   preconditioner_K_Jp_inv.use_matrix(
4669 *   tangent_matrix.block(J_dof, p_dof));
4670 *   ReductionControl solver_control_K_Jp_inv(
4671 *   static_cast<unsigned int>(tangent_matrix.block(J_dof, p_dof).m() *
4672 *   parameters.max_iterations_lin),
4673 *   1.0e-30,
4674 *   parameters.tol_lin);
4675 *   SolverSelector<Vector<double>> solver_K_Jp_inv;
4676 *   solver_K_Jp_inv.select("cg");
4677 *   solver_K_Jp_inv.set_control(solver_control_K_Jp_inv);
4678 *   const auto K_Jp_inv =
4679 *   inverse_operator(K_Jp, solver_K_Jp_inv, preconditioner_K_Jp_inv);
4680 *  
4681 * @endcode
4682 *
4683 * Now we can construct that transpose of
4684 * @f$\mathsf{\mathbf{K}}_{\widetilde{J}\widetilde{p}}^{-1}@f$ and a
4685 * linear operator that represents the condensed operations
4686 * @f$\overline{\mathsf{\mathbf{K}}}@f$ and
4687 * @f$\overline{\overline{\mathsf{\mathbf{K}}}}@f$ and the final
4688 * augmented matrix
4689 * @f$\mathsf{\mathbf{K}}_{\textrm{con}}@f$.
4690 * Note that the schur_complement() operator could also be of use
4691 * here, but for clarity and the purpose of demonstrating the
4692 * similarities between the formulation and implementation of the
4693 * linear solution scheme, we will perform these operations
4694 * manually.
4695 *
4696 * @code
4697 *   const auto K_pJ_inv = transpose_operator(K_Jp_inv);
4698 *   const auto K_pp_bar = K_Jp_inv * K_JJ * K_pJ_inv;
4699 *   const auto K_uu_bar_bar = K_up * K_pp_bar * K_pu;
4700 *   const auto K_uu_con = K_uu + K_uu_bar_bar;
4701 *  
4702 * @endcode
4703 *
4704 * Lastly, we define an operator for inverse of augmented stiffness
4705 * matrix, namely @f$\mathsf{\mathbf{K}}_{\textrm{con}}^{-1}@f$. Note
4706 * that the preconditioner for the augmented stiffness matrix is
4707 * different to the case when we use static condensation. In this
4708 * instance, the preconditioner is based on a non-modified
4709 * @f$\mathsf{\mathbf{K}}_{uu}@f$, while with the first approach we
4710 * actually modified the entries of this sub-block. However, since
4711 * @f$\mathsf{\mathbf{K}}_{\textrm{con}}@f$ and
4712 * @f$\mathsf{\mathbf{K}}_{uu}@f$ operate on the same space, it remains
4713 * adequate for this problem.
4714 *
4715 * @code
4717 *   preconditioner_K_con_inv(parameters.preconditioner_type,
4718 *   parameters.preconditioner_relaxation);
4719 *   preconditioner_K_con_inv.use_matrix(
4720 *   tangent_matrix.block(u_dof, u_dof));
4721 *   ReductionControl solver_control_K_con_inv(
4722 *   static_cast<unsigned int>(tangent_matrix.block(u_dof, u_dof).m() *
4723 *   parameters.max_iterations_lin),
4724 *   1.0e-30,
4725 *   parameters.tol_lin);
4726 *   SolverSelector<Vector<double>> solver_K_con_inv;
4727 *   solver_K_con_inv.select("cg");
4728 *   solver_K_con_inv.set_control(solver_control_K_con_inv);
4729 *   const auto K_uu_con_inv =
4730 *   inverse_operator(K_uu_con,
4731 *   solver_K_con_inv,
4732 *   preconditioner_K_con_inv);
4733 *  
4734 * @endcode
4735 *
4736 * Now we are in a position to solve for the displacement field.
4737 * We can nest the linear operations, and the result is immediately
4738 * written to the Newton update vector.
4739 * It is clear that the implementation closely mimics the derivation
4740 * stated in the introduction.
4741 *
4742 * @code
4743 *   d_u =
4744 *   K_uu_con_inv * (f_u - K_up * (K_Jp_inv * f_J - K_pp_bar * f_p));
4745 *  
4746 *   timer.leave_subsection();
4747 *  
4748 * @endcode
4749 *
4750 * The operations need to post-process for the dilatation and
4751 * pressure fields are just as easy to express.
4752 *
4753 * @code
4754 *   timer.enter_subsection("Linear solver postprocessing");
4755 *   std::cout << " PP " << std::flush;
4756 *  
4757 *   d_J = K_pJ_inv * (f_p - K_pu * d_u);
4758 *   d_p = K_Jp_inv * (f_J - K_JJ * d_J);
4759 *  
4760 *   lin_it = solver_control_K_con_inv.last_step();
4761 *   lin_res = solver_control_K_con_inv.last_value();
4762 *   }
4763 *   else if (parameters.type_lin == "Direct")
4764 *   {
4765 * @endcode
4766 *
4767 * Solve the full block system with
4768 * a direct solver. As it is relatively
4769 * robust, it may be immune to problem
4770 * arising from the presence of the zero
4771 * @f$\mathsf{\mathbf{K}}_{ \widetilde{J} \widetilde{J}}@f$
4772 * block.
4773 *
4774 * @code
4775 *   SparseDirectUMFPACK A_direct;
4776 *   A_direct.initialize(tangent_matrix);
4777 *   A_direct.vmult(newton_update, system_rhs);
4778 *  
4779 *   lin_it = 1;
4780 *   lin_res = 0.0;
4781 *  
4782 *   std::cout << " -- " << std::flush;
4783 *   }
4784 *   else
4785 *   Assert(false, ExcMessage("Linear solver type not implemented"));
4786 *  
4787 *   timer.leave_subsection();
4788 *  
4789 * @endcode
4790 *
4791 * Finally, we again ensure here that any Dirichlet
4792 * constraints are distributed on the updated solution:
4793 *
4794 * @code
4795 *   constraints.distribute(newton_update);
4796 *   }
4797 *  
4798 *   return std::make_pair(lin_it, lin_res);
4799 *   }
4800 *  
4801 * @endcode
4802 *
4803 *
4804 * <a name="step_44-Solidoutput_results"></a>
4805 * <h4>Solid::output_results</h4>
4806 * Here we present how the results are written to file to be viewed
4807 * using ParaView or VisIt. The method is similar to that shown in previous
4808 * tutorials so will not be discussed in detail.
4809 *
4810
4811 *
4812 * @note As of 2023, Visit 3.3.3 can still not deal with higher-order cells.
4813 * Rather, it simply reports that there is no data to show. To view the
4814 * results of this program with Visit, you will want to comment out the
4815 * line that sets `output_flags.write_higher_order_cells = true;`. On the
4816 * other hand, Paraview is able to understand VTU files with higher order
4817 * cells just fine.
4818 *
4819 * @code
4820 *   template <int dim>
4821 *   void Solid<dim>::output_results() const
4822 *   {
4823 *   DataOut<dim> data_out;
4824 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
4825 *   data_component_interpretation(
4827 *   data_component_interpretation.push_back(
4829 *   data_component_interpretation.push_back(
4831 *  
4832 *   std::vector<std::string> solution_name(dim, "displacement");
4833 *   solution_name.emplace_back("pressure");
4834 *   solution_name.emplace_back("dilatation");
4835 *  
4836 *   DataOutBase::VtkFlags output_flags;
4837 *   output_flags.write_higher_order_cells = true;
4838 *   output_flags.physical_units["displacement"] = "m";
4839 *   data_out.set_flags(output_flags);
4840 *  
4841 *   data_out.attach_dof_handler(dof_handler);
4842 *   data_out.add_data_vector(solution_n,
4843 *   solution_name,
4845 *   data_component_interpretation);
4846 *  
4847 * @endcode
4848 *
4849 * Since we are dealing with a large deformation problem, it would be nice
4850 * to display the result on a displaced grid! The MappingQEulerian class
4851 * linked with the DataOut class provides an interface through which this
4852 * can be achieved without physically moving the grid points in the
4853 * Triangulation object ourselves. We first need to copy the solution to
4854 * a temporary vector and then create the Eulerian mapping. We also
4855 * specify the polynomial degree to the DataOut object in order to produce
4856 * a more refined output data set when higher order polynomials are used.
4857 *
4858 * @code
4859 *   Vector<double> soln(solution_n.size());
4860 *   for (unsigned int i = 0; i < soln.size(); ++i)
4861 *   soln(i) = solution_n(i);
4862 *   const MappingQEulerian<dim> q_mapping(degree, dof_handler, soln);
4863 *   data_out.build_patches(q_mapping, degree);
4864 *  
4865 *   std::ofstream output("solution-" + std::to_string(dim) + "d-" +
4866 *   std::to_string(time.get_timestep()) + ".vtu");
4867 *   data_out.write_vtu(output);
4868 *   }
4869 *  
4870 *   } // namespace Step44
4871 *  
4872 *  
4873 * @endcode
4874 *
4875 *
4876 * <a name="step_44-Mainfunction"></a>
4877 * <h3>Main function</h3>
4878 * Lastly we provide the main driver function which appears
4879 * no different to the other tutorials.
4880 *
4881 * @code
4882 *   int main()
4883 *   {
4884 *   using namespace Step44;
4885 *  
4886 *   try
4887 *   {
4888 *   const unsigned int dim = 3;
4889 *   Solid<dim> solid("parameters.prm");
4890 *   solid.run();
4891 *   }
4892 *   catch (std::exception &exc)
4893 *   {
4894 *   std::cerr << std::endl
4895 *   << std::endl
4896 *   << "----------------------------------------------------"
4897 *   << std::endl;
4898 *   std::cerr << "Exception on processing: " << std::endl
4899 *   << exc.what() << std::endl
4900 *   << "Aborting!" << std::endl
4901 *   << "----------------------------------------------------"
4902 *   << std::endl;
4903 *  
4904 *   return 1;
4905 *   }
4906 *   catch (...)
4907 *   {
4908 *   std::cerr << std::endl
4909 *   << std::endl
4910 *   << "----------------------------------------------------"
4911 *   << std::endl;
4912 *   std::cerr << "Unknown exception!" << std::endl
4913 *   << "Aborting!" << std::endl
4914 *   << "----------------------------------------------------"
4915 *   << std::endl;
4916 *   return 1;
4917 *   }
4918 *  
4919 *   return 0;
4920 *   }
4921 * @endcode
4922<a name="step_44-Results"></a><h1>Results</h1>
4923
4924
4925Firstly, we present a comparison of a series of 3-d results with those
4926in the literature (see Reese et al (2000)) to demonstrate that the program works as expected.
4927
4928We begin with a comparison of the convergence with mesh refinement for the @f$Q_1-DGPM_0-DGPM_0@f$ and
4929@f$Q_2-DGPM_1-DGPM_1@f$ formulations, as summarised in the figure below.
4930The vertical displacement of the midpoint of the upper surface of the block is used to assess convergence.
4931Both schemes demonstrate good convergence properties for varying values of the load parameter @f$p/p_0@f$.
4932The results agree with those in the literature.
4933The lower-order formulation typically overestimates the displacement for low levels of refinement,
4934while the higher-order interpolation scheme underestimates it, but be a lesser degree.
4935This benchmark, and a series of others not shown here, give us confidence that the code is working
4936as it should.
4937
4938<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
4939 <tr>
4940 <td align="center">
4941 <img src="https://www.dealii.org/images/steps/developer/step-44.Q1-P0_convergence.png" alt="">
4942 <p align="center">
4943 Convergence of the @f$Q_1-DGPM_0-DGPM_0@f$ formulation in 3-d.
4944 </p>
4945 </td>
4946 <td align="center">
4947 <img src="https://www.dealii.org/images/steps/developer/step-44.Q2-P1_convergence.png" alt="">
4948 <p align="center">
4949 Convergence of the @f$Q_2-DGPM_1-DGPM_1@f$ formulation in 3-d.
4950 </p>
4951 </td>
4952 </tr>
4953</table>
4954
4955
4956A typical screen output generated by running the problem is shown below.
4957The particular case demonstrated is that of the @f$Q_2-DGPM_1-DGPM_1@f$ formulation.
4958It is clear that, using the Newton-Raphson method, quadratic convergence of the solution is obtained.
4959Solution convergence is achieved within 5 Newton increments for all time-steps.
4960The converged displacement's @f$L_2@f$-norm is several orders of magnitude less than the geometry scale.
4961
4962@code
4963Grid:
4964 Reference volume: 1e-09
4966 Number of active cells: 64
4967 Number of degrees of freedom: 2699
4968 Setting up quadrature point data...
4969
4970Timestep 1 @ 0.1s
4971___________________________________________________________________________________________________________________________________________________________
4972 SOLVER STEP | LIN_IT LIN_RES RES_NORM RES_U RES_P RES_J NU_NORM NU_U NU_P NU_J
4973___________________________________________________________________________________________________________________________________________________________
4974 0 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 786 2.118e-06 1.000e+00 1.000e+00 0.000e+00 0.000e+00 1.000e+00 1.000e+00 1.000e+00 1.000e+00
4975 1 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 552 1.031e-03 8.563e-02 8.563e-02 9.200e-13 3.929e-08 1.060e-01 3.816e-02 1.060e-01 1.060e-01
4976 2 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 667 5.602e-06 2.482e-03 2.482e-03 3.373e-15 2.982e-10 2.936e-03 2.053e-04 2.936e-03 2.936e-03
4977 3 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 856 6.469e-10 2.129e-06 2.129e-06 2.245e-19 1.244e-13 1.887e-06 7.289e-07 1.887e-06 1.887e-06
4978 4 ASM_R CONVERGED!
4979___________________________________________________________________________________________________________________________________________________________
4980Relative errors:
4981Displacement: 7.289e-07
4982Force: 2.451e-10
4983Dilatation: 1.353e-07
4984v / V_0: 1.000e-09 / 1.000e-09 = 1.000e+00
4985
4986
4987[...]
4988
4989Timestep 10 @ 1.000e+00s
4990___________________________________________________________________________________________________________________________________________________________
4991 SOLVER STEP | LIN_IT LIN_RES RES_NORM RES_U RES_P RES_J NU_NORM NU_U NU_P NU_J
4992___________________________________________________________________________________________________________________________________________________________
4993 0 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 874 2.358e-06 1.000e+00 1.000e+00 1.000e+00 1.000e+00 1.000e+00 1.000e+00 1.000e+00 1.000e+00
4994 1 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 658 2.942e-04 1.544e-01 1.544e-01 1.208e+13 1.855e+06 6.014e-02 7.398e-02 6.014e-02 6.014e-02
4995 2 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 790 2.206e-06 2.908e-03 2.908e-03 7.302e+10 2.067e+03 2.716e-03 1.433e-03 2.716e-03 2.717e-03
4996 3 ASM_R ASM_K CST ASM_SC SLV PP UQPH | 893 2.374e-09 1.919e-06 1.919e-06 4.527e+07 4.100e+00 1.672e-06 6.842e-07 1.672e-06 1.672e-06
4997 4 ASM_R CONVERGED!
4998___________________________________________________________________________________________________________________________________________________________
4999Relative errors:
5000Displacement: 6.842e-07
5001Force: 8.995e-10
5002Dilatation: 1.528e-06
5003v / V_0: 1.000e-09 / 1.000e-09 = 1.000e+00
5004@endcode
5005
5006
5007
5008Using the Timer class, we can discern which parts of the code require the highest computational expense.
5009For a case with a large number of degrees-of-freedom (i.e. a high level of refinement), a typical output of the Timer is given below.
5010Much of the code in the tutorial has been developed based on the optimizations described,
5011discussed and demonstrated in @ref step_18 "step-18" and others.
5012With over 93% of the time being spent in the linear solver, it is obvious that it may be necessary
5013to invest in a better solver for large three-dimensional problems.
5014The SSOR preconditioner is not multithreaded but is effective for this class of solid problems.
5015It may be beneficial to investigate the use of another solver such as those available through the Trilinos library.
5016
5017
5018@code
5019+---------------------------------------------+------------+------------+
5020| Total wallclock time elapsed since start | 9.874e+02s | |
5021| | | |
5022| Section | no. calls | wall time | % of total |
5023+---------------------------------+-----------+------------+------------+
5024| Assemble system right-hand side | 53 | 1.727e+00s | 1.75e-01% |
5025| Assemble tangent matrix | 43 | 2.707e+01s | 2.74e+00% |
5026| Linear solver | 43 | 9.248e+02s | 9.37e+01% |
5027| Linear solver postprocessing | 43 | 2.743e-02s | 2.78e-03% |
5028| Perform static condensation | 43 | 1.437e+01s | 1.46e+00% |
5029| Setup system | 1 | 3.897e-01s | 3.95e-02% |
5030| Update QPH data | 43 | 5.770e-01s | 5.84e-02% |
5031+---------------------------------+-----------+------------+------------+
5032@endcode
5033
5034
5035We then used ParaView to visualize the results for two cases.
5036The first was for the coarsest grid and the lowest-order interpolation method: @f$Q_1-DGPM_0-DGPM_0@f$.
5037The second was on a refined grid using a @f$Q_2-DGPM_1-DGPM_1@f$ formulation.
5038The vertical component of the displacement, the pressure @f$\widetilde{p}@f$ and the dilatation @f$\widetilde{J}@f$ fields
5039are shown below.
5040
5041
5042For the first case it is clear that the coarse spatial discretization coupled with large displacements leads to a low quality solution
5043(the loading ratio is @f$p/p_0=80@f$).
5044Additionally, the pressure difference between elements is very large.
5045The constant pressure field on the element means that the large pressure gradient is not captured.
5046However, it should be noted that locking, which would be present in a standard @f$Q_1@f$ displacement formulation does not arise
5047even in this poorly discretised case.
5048The final vertical displacement of the tracked node on the top surface of the block is still within 12.5% of the converged solution.
5049The pressure solution is very coarse and has large jumps between adjacent cells.
5050It is clear that the volume nearest to the applied traction undergoes compression while the outer extents
5051of the domain are in a state of expansion.
5052The dilatation solution field and pressure field are clearly linked,
5053with positive dilatation indicating regions of positive pressure and negative showing regions placed in compression.
5054As discussed in the Introduction, a compressive pressure has a negative sign
5055while an expansive pressure takes a positive sign.
5056This stems from the definition of the volumetric strain energy function
5057and is opposite to the physically realistic interpretation of pressure.
5058
5059
5060<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
5061 <tr>
5062 <td align="center">
5063 <img src="https://www.dealii.org/images/steps/developer/step-44.Q1-P0_gr_1_p_ratio_80-displacement.png" alt="">
5064 <p align="center">
5065 Z-displacement solution for the 3-d problem.
5066 </p>
5067 </td>
5068 <td align="center">
5069 <img src="https://www.dealii.org/images/steps/developer/step-44.Q1-P0_gr_1_p_ratio_80-pressure.png" alt="">
5070 <p align="center">
5071 Discontinuous piece-wise constant pressure field.
5072 </p>
5073 </td>
5074 <td align="center">
5075 <img src="https://www.dealii.org/images/steps/developer/step-44.Q1-P0_gr_1_p_ratio_80-dilatation.png" alt="">
5076 <p align="center">
5077 Discontinuous piece-wise constant dilatation field.
5078 </p>
5079 </td>
5080 </tr>
5081</table>
5082
5083Combining spatial refinement and a higher-order interpolation scheme results in a high-quality solution.
5084Three grid refinements coupled with a @f$Q_2-DGPM_1-DGPM_1@f$ formulation produces
5085a result that clearly captures the mechanics of the problem.
5086The deformation of the traction surface is well resolved.
5087We can now observe the actual extent of the applied traction, with the maximum force being applied
5088at the central point of the surface causing the largest compression.
5089Even though very high strains are experienced in the domain,
5090especially at the boundary of the region of applied traction,
5091the solution remains accurate.
5092The pressure field is captured in far greater detail than before.
5093There is a clear distinction and transition between regions of compression and expansion,
5094and the linear approximation of the pressure field allows a refined visualization
5095of the pressure at the sub-element scale.
5096It should however be noted that the pressure field remains discontinuous
5097and could be smoothed on a continuous grid for the post-processing purposes.
5098
5099
5100
5101<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
5102 <tr>
5103 <td align="center">
5104 <img src="https://www.dealii.org/images/steps/developer/step-44.Q2-P1_gr_3_p_ratio_80-displacement.png" alt="">
5105 <p align="center">
5106 Z-displacement solution for the 3-d problem.
5107 </p>
5108 </td>
5109 <td align="center">
5110 <img src="https://www.dealii.org/images/steps/developer/step-44.Q2-P1_gr_3_p_ratio_80-pressure.png" alt="">
5111 <p align="center">
5112 Discontinuous linear pressure field.
5113 </p>
5114 </td>
5115 <td align="center">
5116 <img src="https://www.dealii.org/images/steps/developer/step-44.Q2-P1_gr_3_p_ratio_80-dilatation.png" alt="">
5117 <p align="center">
5118 Discontinuous linear dilatation field.
5119 </p>
5120 </td>
5121 </tr>
5122</table>
5123
5124This brief analysis of the results demonstrates that the three-field formulation is effective
5125in circumventing volumetric locking for highly-incompressible media.
5126The mixed formulation is able to accurately simulate the displacement of a
5127near-incompressible block under compression.
5128The command-line output indicates that the volumetric change under extreme compression resulted in
5129less than 0.01% volume change for a Poisson's ratio of 0.4999.
5130
5131In terms of run-time, the @f$Q_2-DGPM_1-DGPM_1@f$ formulation tends to be more computationally expensive
5132than the @f$Q_1-DGPM_0-DGPM_0@f$ for a similar number of degrees-of-freedom
5133(produced by adding an extra grid refinement level for the lower-order interpolation).
5134This is shown in the graph below for a batch of tests run consecutively on a single 4-core (8-thread) machine.
5135The increase in computational time for the higher-order method is likely due to
5136the increased band-width required for the higher-order elements.
5137As previously mentioned, the use of a better solver and preconditioner may mitigate the
5138expense of using a higher-order formulation.
5139It was observed that for the given problem using the multithreaded Jacobi preconditioner can reduce the
5140computational runtime by up to 72% (for the worst case being a higher-order formulation with a large number
5141of degrees-of-freedom) in comparison to the single-thread SSOR preconditioner.
5142However, it is the author's experience that the Jacobi method of preconditioning may not be suitable for
5143some finite-strain problems involving alternative constitutive models.
5144
5145
5146<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
5147 <tr>
5148 <td align="center">
5149 <img src="https://www.dealii.org/images/steps/developer/step-44.Normalised_runtime.png" alt="">
5150 <p align="center">
5151 Runtime on a 4-core machine, normalised against the lowest grid resolution @f$Q_1-DGPM_0-DGPM_0@f$ solution that utilised a SSOR preconditioner.
5152 </p>
5153 </td>
5154 </tr>
5155</table>
5156
5157
5158Lastly, results for the displacement solution for the 2-d problem are showcased below for
5159two different levels of grid refinement.
5160It is clear that due to the extra constraints imposed by simulating in 2-d that the resulting
5161displacement field, although qualitatively similar, is different to that of the 3-d case.
5162
5163
5164<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
5165 <tr>
5166 <td align="center">
5167 <img src="https://www.dealii.org/images/steps/developer/step-44.2d-gr_2.png" alt="">
5168 <p align="center">
5169 Y-displacement solution in 2-d for 2 global grid refinement levels.
5170 </p>
5171 </td>
5172 <td align="center">
5173 <img src="https://www.dealii.org/images/steps/developer/step-44.2d-gr_5.png" alt="">
5174 <p align="center">
5175 Y-displacement solution in 2-d for 5 global grid refinement levels.
5176 </p>
5177 </td>
5178 </tr>
5179</table>
5180
5181<a name="step-44-extensions"></a>
5182<a name="step_44-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
5183
5184
5185There are a number of obvious extensions for this work:
5186
5187- Firstly, an additional constraint could be added to the free-energy
5188 function in order to enforce a high degree of incompressibility in
5189 materials. An additional Lagrange multiplier would be introduced,
5190 but this could most easily be dealt with using the principle of
5191 augmented Lagrange multipliers. This is demonstrated in <em>Simo and
5192 Taylor (1991) </em>.
5193- The constitutive relationship used in this
5194 model is relatively basic. It may be beneficial to split the material
5195 class into two separate classes, one dealing with the volumetric
5196 response and the other the isochoric response, and produce a generic
5197 materials class (i.e. having abstract virtual functions that derived
5198 classes have to implement) that would allow for the addition of more complex
5199 material models. Such models could include other hyperelastic
5200 materials, plasticity and viscoelastic materials and others.
5201- The program has been developed for solving problems on single-node
5202 multicore machines. With a little effort, the program could be
5203 extended to a large-scale computing environment through the use of
5204 PETSc or Trilinos, using a similar technique to that demonstrated in
5205 @ref step_40 "step-40". This would mostly involve changes to the setup, assembly,
5206 <code>PointHistory</code> and linear solver routines.
5207- As this program assumes quasi-static equilibrium, extensions to
5208 include dynamic effects would be necessary to study problems where
5209 inertial effects are important, e.g. problems involving impact.
5210- Load and solution limiting procedures may be necessary for highly
5211 nonlinear problems. It is possible to add a linesearch algorithm to
5212 limit the step size within a Newton increment to ensure optimum
5213 convergence. It may also be necessary to use a load limiting method,
5214 such as the Riks method, to solve unstable problems involving
5215 geometric instability such as buckling and snap-through.
5216- Many physical problems involve contact. It is possible to include
5217 the effect of frictional or frictionless contact between objects
5218 into this program. This would involve the addition of an extra term
5219 in the free-energy functional and therefore an addition to the
5220 assembly routine. One would also need to manage the contact problem
5221 (detection and stress calculations) itself. An alternative to
5222 additional penalty terms in the free-energy functional would be to
5223 use active set methods such as the one used in @ref step_41 "step-41".
5224- The complete condensation procedure using LinearOperators has been
5225 coded into the linear solver routine. This could also have been
5226 achieved through the application of the schur_complement()
5227 operator to condense out one or more of the fields in a more
5228 automated manner.
5229- Finally, adaptive mesh refinement, as demonstrated in @ref step_6 "step-6" and
5230 @ref step_18 "step-18", could provide additional solution accuracy.
5231 *
5232 *
5233<a name="step_44-PlainProg"></a>
5234<h1> The plain program</h1>
5235@include "step-44.cc"
5236*/
void select(const std::string &name)
void initialize(const SparsityPattern &sparsity_pattern)
Definition timer.h:117
Point< 3 > center
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
#define Assert(cond, exc)
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
LinearOperator< Domain, Range, Payload > transpose_operator(const LinearOperator< Range, Domain, Payload > &op)
LinearOperator< Range_2, Domain_2, Payload > schur_complement(const LinearOperator< Domain_1, Range_1, Payload > &A_inv, const LinearOperator< Range_1, Domain_2, Payload > &B, const LinearOperator< Range_2, Domain_1, Payload > &C, const LinearOperator< Range_2, Domain_2, Payload > &D)
LinearOperator< Domain, Range, Payload > inverse_operator(const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const Preconditioner &preconditioner)
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
UpdateFlags
Task< RT > new_task(const std::function< RT()> &function)
std::vector< value_type > split(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
const Event initial
Definition event.cc:64
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
Expression sign(const Expression &x)
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
double volume(const Triangulation< dim, spacedim > &tria)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ diagonal
Matrix is diagonal.
@ general
No special properties.
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
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * end(VectorType &V)
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
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)
bool check(const ConstraintKinds kind_in, const unsigned int dim)
DEAL_II_HOST constexpr TableIndices< 2 > merge(const TableIndices< 2 > &previous_indices, const unsigned int new_index, const unsigned int position)
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)
Definition types.h:32
unsigned int boundary_id
Definition types.h:144
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation