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