Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Quasi_static_Finite_strain_Quasi_incompressible_ViscoElasticity.h
Go to the documentation of this file.
1
224 *  
225 *   #include <deal.II/base/mpi.h>
226 *   #include <deal.II/base/function.h>
227 *   #include <deal.II/base/parameter_handler.h>
228 *   #include <deal.II/base/point.h>
229 *   #include <deal.II/base/quadrature_lib.h>
230 *   #include <deal.II/base/symmetric_tensor.h>
231 *   #include <deal.II/base/tensor.h>
232 *   #include <deal.II/base/timer.h>
233 *   #include <deal.II/base/work_stream.h>
234 *   #include <deal.II/dofs/dof_renumbering.h>
235 *   #include <deal.II/dofs/dof_tools.h>
236 *   #include <deal.II/base/quadrature_point_data.h>
237 *   #include <deal.II/grid/filtered_iterator.h>
238 *   #include <deal.II/grid/grid_generator.h>
239 *   #include <deal.II/grid/grid_tools.h>
240 *   #include <deal.II/grid/grid_in.h>
241 *   #include <deal.II/grid/manifold_lib.h>
242 *   #include <deal.II/grid/tria.h>
243 *   #include <deal.II/fe/fe_dgp_monomial.h>
244 *   #include <deal.II/fe/fe_q.h>
245 *   #include <deal.II/fe/fe_system.h>
246 *   #include <deal.II/fe/fe_tools.h>
247 *   #include <deal.II/fe/fe_values.h>
248 *   #include <deal.II/fe/mapping_q_eulerian.h>
249 *   #include <deal.II/lac/block_sparsity_pattern.h>
250 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
251 *   #include <deal.II/lac/affine_constraints.h>
252 *   #include <deal.II/lac/full_matrix.h>
253 *   #include <deal.II/lac/solver_selector.h>
254 *   #include <deal.II/lac/trilinos_block_sparse_matrix.h>
255 *   #include <deal.II/lac/trilinos_precondition.h>
256 *   #include <deal.II/lac/trilinos_sparsity_pattern.h>
257 *   #include <deal.II/lac/trilinos_sparse_matrix.h>
258 *   #include <deal.II/lac/trilinos_vector.h>
259 *   #include <deal.II/lac/linear_operator.h>
260 *   #include <deal.II/lac/packaged_operation.h>
261 *   #include <deal.II/lac/trilinos_linear_operator.h>
262 *   #include <deal.II/numerics/data_out.h>
263 *   #include <deal.II/numerics/vector_tools.h>
264 *   #include <deal.II/physics/transformations.h>
265 *   #include <deal.II/physics/elasticity/kinematics.h>
266 *   #include <deal.II/physics/elasticity/standard_tensors.h>
267 *   #include <iostream>
268 *   #include <fstream>
269 *   #include <numeric>
270 *  
271 *   #include <deal.II/grid/grid_out.h>
272 *  
273 *   namespace ViscoElasStripHole
274 *   {
275 *   using namespace dealii;
276 *   namespace LA = TrilinosWrappers;
277 *   namespace Parameters
278 *   {
279 *   struct BoundaryConditions
280 *   {
281 *   BoundaryConditions();
282 *  
283 *   std::string driver;
284 *   double stretch;
285 *   double pressure;
286 *   double load_time;
287 *  
288 *   const types::boundary_id boundary_id_minus_X;
289 *   const types::boundary_id boundary_id_plus_X;
290 *   const types::boundary_id boundary_id_minus_Y;
291 *   const types::boundary_id boundary_id_plus_Y;
292 *   const types::boundary_id boundary_id_minus_Z;
293 *   const types::boundary_id boundary_id_plus_Z;
294 *   const types::boundary_id boundary_id_hole;
295 *   const types::manifold_id manifold_id_hole;
296 *  
297 *   static void
298 *   declare_parameters(ParameterHandler &prm);
299 *   void
300 *   parse_parameters(ParameterHandler &prm);
301 *   };
302 *   BoundaryConditions::BoundaryConditions()
303 *   :
304 *   driver ("Neumann"),
305 *   stretch (2.0),
306 *   pressure(0.0),
307 *   load_time(2.5),
308 *   boundary_id_minus_X (1),
309 *   boundary_id_plus_X (2),
310 *   boundary_id_minus_Y (3),
311 *   boundary_id_plus_Y (4),
312 *   boundary_id_minus_Z (5),
313 *   boundary_id_plus_Z (6),
314 *   boundary_id_hole (10),
315 *   manifold_id_hole (10)
316 *   { }
317 *   void BoundaryConditions::declare_parameters(ParameterHandler &prm)
318 *   {
319 *   prm.enter_subsection("Boundary conditions");
320 *   {
321 *   prm.declare_entry("Driver", "Dirichlet",
322 *   Patterns::Selection("Dirichlet|Neumann"),
323 *   "Driver boundary condition for the problem");
324 *   prm.declare_entry("Final stretch", "2.0",
325 *   Patterns::Double(1.0),
326 *   "Positive stretch applied length-ways to the strip");
327 *   prm.declare_entry("Applied pressure", "0.0",
328 *   Patterns::Double(-1e3,1e3),
329 *   "Hydrostatic pressure applied (in the referential configuration) to the interior surface of the hole");
330 *   prm.declare_entry("Load time", "2.5",
331 *   Patterns::Double(0.0),
332 *   "Total time over which the stretch/pressure is ramped up");
333 *   }
334 *   prm.leave_subsection();
335 *   }
336 *   void BoundaryConditions::parse_parameters(ParameterHandler &prm)
337 *   {
338 *   prm.enter_subsection("Boundary conditions");
339 *   {
340 *   driver = prm.get("Driver");
341 *   stretch = prm.get_double("Final stretch");
342 *   pressure = prm.get_double("Applied pressure");
343 *   load_time = prm.get_double("Load time");
344 *   }
345 *   prm.leave_subsection();
346 *   }
347 *   struct FESystem
348 *   {
349 *   unsigned int poly_degree;
350 *   unsigned int quad_order;
351 *   static void
352 *   declare_parameters(ParameterHandler &prm);
353 *   void
354 *   parse_parameters(ParameterHandler &prm);
355 *   };
356 *   void FESystem::declare_parameters(ParameterHandler &prm)
357 *   {
358 *   prm.enter_subsection("Finite element system");
359 *   {
360 *   prm.declare_entry("Polynomial degree", "2",
361 *   Patterns::Integer(0),
362 *   "Displacement system polynomial order");
363 *   prm.declare_entry("Quadrature order", "3",
364 *   Patterns::Integer(0),
365 *   "Gauss quadrature order");
366 *   }
367 *   prm.leave_subsection();
368 *   }
369 *   void FESystem::parse_parameters(ParameterHandler &prm)
370 *   {
371 *   prm.enter_subsection("Finite element system");
372 *   {
373 *   poly_degree = prm.get_integer("Polynomial degree");
374 *   quad_order = prm.get_integer("Quadrature order");
375 *   }
376 *   prm.leave_subsection();
377 *   }
378 *   struct Geometry
379 *   {
380 *   double length;
381 *   double width;
382 *   double thickness;
383 *   double hole_diameter;
384 *   double hole_division_fraction;
385 *   unsigned int n_repetitions_xy;
386 *   unsigned int n_repetitions_z;
387 *   unsigned int global_refinement;
388 *   double scale;
389 *   static void
390 *   declare_parameters(ParameterHandler &prm);
391 *   void
392 *   parse_parameters(ParameterHandler &prm);
393 *   };
394 *   void Geometry::declare_parameters(ParameterHandler &prm)
395 *   {
396 *   prm.enter_subsection("Geometry");
397 *   {
398 *   prm.declare_entry("Length", "100.0",
399 *   Patterns::Double(0.0),
400 *   "Total sample length");
401 *   prm.declare_entry("Width", "50.0",
402 *   Patterns::Double(0.0),
403 *   "Total sample width");
404 *   prm.declare_entry("Thickness", "5.0",
405 *   Patterns::Double(0.0),
406 *   "Total sample thickness");
407 *   prm.declare_entry("Hole diameter", "20.0",
408 *   Patterns::Double(0.0),
409 *   "Hole diameter");
410 *   prm.declare_entry("Hole division fraction", "0.5",
411 *   Patterns::Double(0.0,1.0),
412 *   "A geometric factor affecting the discretisation near the hole");
413 *   prm.declare_entry("Number of subdivisions in cross-section", "2",
414 *   Patterns::Integer(1.0),
415 *   "A factor defining the number of initial grid subdivisions in the cross-section");
416 *   prm.declare_entry("Number of subdivisions thickness", "6",
417 *   Patterns::Integer(1.0),
418 *   "A factor defining the number of initial grid subdivisions through the thickness");
419 *   prm.declare_entry("Global refinement", "2",
420 *   Patterns::Integer(0),
421 *   "Global refinement level");
422 *   prm.declare_entry("Grid scale", "1e-3",
423 *   Patterns::Double(0.0),
424 *   "Global grid scaling factor");
425 *   }
426 *   prm.leave_subsection();
427 *   }
428 *   void Geometry::parse_parameters(ParameterHandler &prm)
429 *   {
430 *   prm.enter_subsection("Geometry");
431 *   {
432 *   length = prm.get_double("Length");
433 *   width = prm.get_double("Width");
434 *   thickness = prm.get_double("Thickness");
435 *   hole_diameter = prm.get_double("Hole diameter");
436 *   hole_division_fraction = prm.get_double("Hole division fraction");
437 *   n_repetitions_xy = prm.get_integer("Number of subdivisions in cross-section");
438 *   n_repetitions_z = prm.get_integer("Number of subdivisions thickness");
439 *   global_refinement = prm.get_integer("Global refinement");
440 *   scale = prm.get_double("Grid scale");
441 *   }
442 *   prm.leave_subsection();
443 *   }
444 *   struct Materials
445 *   {
446 *   double nu_e;
447 *   double mu_e;
448 *   double mu_v;
449 *   double tau_v;
450 *   static void
451 *   declare_parameters(ParameterHandler &prm);
452 *   void
453 *   parse_parameters(ParameterHandler &prm);
454 *   };
455 *   void Materials::declare_parameters(ParameterHandler &prm)
456 *   {
457 *   prm.enter_subsection("Material properties");
458 *   {
459 *   prm.declare_entry("Poisson's ratio", "0.4999",
460 *   Patterns::Double(-1.0,0.5),
461 *   "Poisson's ratio");
462 *   prm.declare_entry("Elastic shear modulus", "80.194e6",
463 *   Patterns::Double(0.0),
464 *   "Elastic shear modulus");
465 *   prm.declare_entry("Viscous shear modulus", "80.194e6",
466 *   Patterns::Double(0.0),
467 *   "Viscous shear modulus");
468 *   prm.declare_entry("Viscous relaxation time", "2.0",
469 *   Patterns::Double(0.0),
470 *   "Viscous relaxation time");
471 *   }
472 *   prm.leave_subsection();
473 *   }
474 *   void Materials::parse_parameters(ParameterHandler &prm)
475 *   {
476 *   prm.enter_subsection("Material properties");
477 *   {
478 *   nu_e = prm.get_double("Poisson's ratio");
479 *   mu_e = prm.get_double("Elastic shear modulus");
480 *   mu_v = prm.get_double("Viscous shear modulus");
481 *   tau_v = prm.get_double("Viscous relaxation time");
482 *   }
483 *   prm.leave_subsection();
484 *   }
485 *   struct LinearSolver
486 *   {
487 *   std::string type_lin;
488 *   double tol_lin;
489 *   double max_iterations_lin;
490 *   static void
491 *   declare_parameters(ParameterHandler &prm);
492 *   void
493 *   parse_parameters(ParameterHandler &prm);
494 *   };
495 *   void LinearSolver::declare_parameters(ParameterHandler &prm)
496 *   {
497 *   prm.enter_subsection("Linear solver");
498 *   {
499 *   prm.declare_entry("Solver type", "cg",
501 *   "Type of solver used to solve the linear system");
502 *   prm.declare_entry("Residual", "1e-6",
503 *   Patterns::Double(0.0),
504 *   "Linear solver residual (scaled by residual norm)");
505 *   prm.declare_entry("Max iteration multiplier", "1",
506 *   Patterns::Double(1.0),
507 *   "Linear solver iterations (multiples of the system matrix size)");
508 *   }
509 *   prm.leave_subsection();
510 *   }
511 *   void LinearSolver::parse_parameters(ParameterHandler &prm)
512 *   {
513 *   prm.enter_subsection("Linear solver");
514 *   {
515 *   type_lin = prm.get("Solver type");
516 *   tol_lin = prm.get_double("Residual");
517 *   max_iterations_lin = prm.get_double("Max iteration multiplier");
518 *   }
519 *   prm.leave_subsection();
520 *   }
521 *   struct NonlinearSolver
522 *   {
523 *   unsigned int max_iterations_NR;
524 *   double tol_f;
525 *   double tol_u;
526 *   static void
527 *   declare_parameters(ParameterHandler &prm);
528 *   void
529 *   parse_parameters(ParameterHandler &prm);
530 *   };
531 *   void NonlinearSolver::declare_parameters(ParameterHandler &prm)
532 *   {
533 *   prm.enter_subsection("Nonlinear solver");
534 *   {
535 *   prm.declare_entry("Max iterations Newton-Raphson", "10",
536 *   Patterns::Integer(0),
537 *   "Number of Newton-Raphson iterations allowed");
538 *   prm.declare_entry("Tolerance displacement", "1.0e-6",
539 *   Patterns::Double(0.0),
540 *   "Displacement error tolerance");
541 *   prm.declare_entry("Tolerance force", "1.0e-9",
542 *   Patterns::Double(0.0),
543 *   "Force residual tolerance");
544 *   }
545 *   prm.leave_subsection();
546 *   }
547 *   void NonlinearSolver::parse_parameters(ParameterHandler &prm)
548 *   {
549 *   prm.enter_subsection("Nonlinear solver");
550 *   {
551 *   max_iterations_NR = prm.get_integer("Max iterations Newton-Raphson");
552 *   tol_f = prm.get_double("Tolerance force");
553 *   tol_u = prm.get_double("Tolerance displacement");
554 *   }
555 *   prm.leave_subsection();
556 *   }
557 *   struct Time
558 *   {
559 *   double delta_t;
560 *   double end_time;
561 *   static void
562 *   declare_parameters(ParameterHandler &prm);
563 *   void
564 *   parse_parameters(ParameterHandler &prm);
565 *   };
566 *   void Time::declare_parameters(ParameterHandler &prm)
567 *   {
568 *   prm.enter_subsection("Time");
569 *   {
570 *   prm.declare_entry("End time", "1",
571 *   Patterns::Double(),
572 *   "End time");
573 *   prm.declare_entry("Time step size", "0.1",
574 *   Patterns::Double(),
575 *   "Time step size");
576 *   }
577 *   prm.leave_subsection();
578 *   }
579 *   void Time::parse_parameters(ParameterHandler &prm)
580 *   {
581 *   prm.enter_subsection("Time");
582 *   {
583 *   end_time = prm.get_double("End time");
584 *   delta_t = prm.get_double("Time step size");
585 *   }
586 *   prm.leave_subsection();
587 *   }
588 *   struct AllParameters
589 *   : public BoundaryConditions,
590 *   public FESystem,
591 *   public Geometry,
592 *   public Materials,
593 *   public LinearSolver,
594 *   public NonlinearSolver,
595 *   public Time
596 *   {
597 *   AllParameters(const std::string &input_file);
598 *   static void
599 *   declare_parameters(ParameterHandler &prm);
600 *   void
601 *   parse_parameters(ParameterHandler &prm);
602 *   };
603 *   AllParameters::AllParameters(const std::string &input_file)
604 *   {
605 *   ParameterHandler prm;
606 *   declare_parameters(prm);
607 *   prm.parse_input(input_file);
608 *   parse_parameters(prm);
609 *   }
610 *   void AllParameters::declare_parameters(ParameterHandler &prm)
611 *   {
612 *   BoundaryConditions::declare_parameters(prm);
613 *   FESystem::declare_parameters(prm);
614 *   Geometry::declare_parameters(prm);
615 *   Materials::declare_parameters(prm);
616 *   LinearSolver::declare_parameters(prm);
617 *   NonlinearSolver::declare_parameters(prm);
618 *   Time::declare_parameters(prm);
619 *   }
620 *   void AllParameters::parse_parameters(ParameterHandler &prm)
621 *   {
622 *   BoundaryConditions::parse_parameters(prm);
623 *   FESystem::parse_parameters(prm);
624 *   Geometry::parse_parameters(prm);
625 *   Materials::parse_parameters(prm);
626 *   LinearSolver::parse_parameters(prm);
627 *   NonlinearSolver::parse_parameters(prm);
628 *   Time::parse_parameters(prm);
629 *   }
630 *   }
631 *   class Time
632 *   {
633 *   public:
634 *   Time (const double time_end,
635 *   const double delta_t)
636 *   :
637 *   timestep(0),
638 *   time_current(0.0),
639 *   time_end(time_end),
640 *   delta_t(delta_t)
641 *   {}
642 *   virtual ~Time()
643 *   {}
644 *   double current() const
645 *   {
646 *   return time_current;
647 *   }
648 *   double end() const
649 *   {
650 *   return time_end;
651 *   }
652 *   double get_delta_t() const
653 *   {
654 *   return delta_t;
655 *   }
656 *   unsigned int get_timestep() const
657 *   {
658 *   return timestep;
659 *   }
660 *   void increment()
661 *   {
662 *   time_current += delta_t;
663 *   ++timestep;
664 *   }
665 *   private:
666 *   unsigned int timestep;
667 *   double time_current;
668 *   const double time_end;
669 *   const double delta_t;
670 *   };
671 *   template <int dim>
672 *   class Material_Compressible_Three_Field_Linear_Viscoelastic
673 *   {
674 *   public:
675 *   Material_Compressible_Three_Field_Linear_Viscoelastic(const double mu_e,
676 *   const double nu_e,
677 *   const double mu_v,
678 *   const double tau_v,
679 *   const Time &time)
680 *   :
681 *   kappa((2.0 * mu_e * (1.0 + nu_e)) / (3.0 * (1.0 - 2.0 * nu_e))),
682 *   mu_e(mu_e),
683 *   mu_v(mu_v),
684 *   tau_v(tau_v),
685 *   time(time),
688 *   {
689 *   Assert(kappa > 0, ExcInternalError());
690 *   }
691 *   ~Material_Compressible_Three_Field_Linear_Viscoelastic()
692 *   {}
693 *  
695 *   get_tau(const Tensor<2,dim> &F,
696 *   const double &p_tilde) const
697 *   {
698 *   return get_tau_iso(F) + get_tau_vol(F,p_tilde);
699 *   }
700 *   SymmetricTensor<4,dim> get_Jc(const Tensor<2,dim> &F,
701 *   const double &p_tilde) const
702 *   {
703 *   return get_Jc_iso(F) + get_Jc_vol(F,p_tilde);
704 *   }
705 *   double
706 *   get_dPsi_vol_dJ(const double &J_tilde) const
707 *   {
708 *   return (kappa / 2.0) * (J_tilde - 1.0 / J_tilde);
709 *   }
710 *   double
711 *   get_d2Psi_vol_dJ2(const double &J_tilde) const
712 *   {
713 *   return ( (kappa / 2.0) * (1.0 + 1.0 / (J_tilde * J_tilde)));
714 *   }
715 *   void
716 *   update_internal_equilibrium(const Tensor<2, dim> &F,
717 *   const double &/*p_tilde*/,
718 *   const double &/*J_tilde*/)
719 *   {
720 *   const double det_F = determinant(F);
721 *   const SymmetricTensor<2,dim> C_bar = std::pow(det_F, -2.0 / dim) * Physics::Elasticity::Kinematics::C(F);
722 * @endcode
723 *
724 * Linder2011 eq 54
725 * Assumes first-oder backward Euler time discretisation
726 *
727 * @code
728 *   Q_n_t = (1.0/(1.0 + time.get_delta_t()/tau_v))*(Q_t1 + (time.get_delta_t()/tau_v)*invert(C_bar));
729 *   }
730 *   void
731 *   update_end_timestep()
732 *   {
733 *   Q_t1 = Q_n_t;
734 *   }
735 *  
736 *   protected:
737 *   const double kappa;
738 *   const double mu_e;
739 *   const double mu_v;
740 *   const double tau_v;
741 *   const Time &time;
742 *   SymmetricTensor<2,dim> Q_n_t; // Value of internal variable at this Newton step and timestep
743 *   SymmetricTensor<2,dim> Q_t1; // Value of internal variable at the previous timestep
744 *  
746 *   get_tau_vol(const Tensor<2,dim> &F,
747 *   const double &p_tilde) const
748 *   {
749 *   const double det_F = determinant(F);
750 *  
751 *   return p_tilde * det_F * Physics::Elasticity::StandardTensors<dim>::I;
752 *   }
754 *   get_tau_iso(const Tensor<2,dim> &F) const
755 *   {
756 *   return Physics::Elasticity::StandardTensors<dim>::dev_P * get_tau_bar(F);
757 *   }
759 *   get_tau_bar(const Tensor<2,dim> &F) const
760 *   {
761 *   const double det_F = determinant(F);
762 *   const Tensor<2,dim> F_bar = std::pow(det_F, -1.0 / dim) * F;
763 *   const SymmetricTensor<2,dim> b_bar = std::pow(det_F, -2.0 / dim) * symmetrize(F * transpose(F));
764 * @endcode
765 *
766 * Elastic Neo-Hookean + Linder2011 eq 47
767 *
768 * @code
769 *   return mu_e * b_bar
770 *   + mu_v * symmetrize(F_bar*static_cast<Tensor<2,dim> >(Q_n_t)*transpose(F_bar));
771 *   }
772 *   SymmetricTensor<4, dim> get_Jc_vol(const Tensor<2,dim> &F,
773 *   const double &p_tilde) const
774 *   {
775 *   const double det_F = determinant(F);
776 *   return p_tilde * det_F
779 *   }
780 *   SymmetricTensor<4, dim> get_Jc_iso(const Tensor<2,dim> &F) const
781 *   {
782 *   const SymmetricTensor<2, dim> tau_bar = get_tau_bar(F);
783 *   const SymmetricTensor<2, dim> tau_iso = get_tau_iso(F);
784 *   const SymmetricTensor<4, dim> tau_iso_x_I
785 *   = outer_product(tau_iso,
787 *   const SymmetricTensor<4, dim> I_x_tau_iso
789 *   tau_iso);
790 *   const SymmetricTensor<4, dim> c_bar = get_c_bar(F);
791 *   return (2.0 / dim) * trace(tau_bar)
793 *   - (2.0 / dim) * (tau_iso_x_I + I_x_tau_iso)
796 *   }
797 *   SymmetricTensor<4, dim> get_c_bar(const Tensor<2,dim> &/*F*/) const
798 *   {
799 * @endcode
800 *
801 * Elastic Neo-Hookean + Linder2011 eq 56
802 *
803 * @code
804 *   return -2.0*mu_v*((time.get_delta_t()/tau_v)/(1.0 + time.get_delta_t()/tau_v))*Physics::Elasticity::StandardTensors<dim>::S;
805 *   }
806 *   };
807 *   template <int dim>
808 *   class PointHistory
809 *   {
810 *   public:
811 *   PointHistory()
812 *   {}
813 *   virtual ~PointHistory()
814 *   {}
815 *   void
816 *   setup_lqp (const Parameters::AllParameters &parameters,
817 *   const Time &time)
818 *   {
819 *   material.reset(new Material_Compressible_Three_Field_Linear_Viscoelastic<dim>(
820 *   parameters.mu_e, parameters.nu_e,
821 *   parameters.mu_v, parameters.tau_v,
822 *   time));
823 *   }
824 *  
826 *   get_tau(const Tensor<2, dim> &F,
827 *   const double &p_tilde) const
828 *   {
829 *   return material->get_tau(F, p_tilde);
830 *   }
832 *   get_Jc(const Tensor<2, dim> &F,
833 *   const double &p_tilde) const
834 *   {
835 *   return material->get_Jc(F, p_tilde);
836 *   }
837 *   double
838 *   get_dPsi_vol_dJ(const double &J_tilde) const
839 *   {
840 *   return material->get_dPsi_vol_dJ(J_tilde);
841 *   }
842 *   double
843 *   get_d2Psi_vol_dJ2(const double &J_tilde) const
844 *   {
845 *   return material->get_d2Psi_vol_dJ2(J_tilde);
846 *   }
847 *   void
848 *   update_internal_equilibrium(const Tensor<2, dim> &F,
849 *   const double &p_tilde,
850 *   const double &J_tilde)
851 *   {
852 *   material->update_internal_equilibrium(F,p_tilde,J_tilde);
853 *   }
854 *   void
855 *   update_end_timestep()
856 *   {
857 *   material->update_end_timestep();
858 *   }
859 *   private:
860 *   std::shared_ptr< Material_Compressible_Three_Field_Linear_Viscoelastic<dim> > material;
861 *   };
862 *   template <int dim>
863 *   class Solid
864 *   {
865 *   public:
866 *   Solid(const std::string &input_file);
867 *   virtual
868 *   ~Solid();
869 *   void
870 *   run();
871 *   private:
872 *   struct PerTaskData_ASM;
873 *   struct ScratchData_ASM;
874 *   void
875 *   make_grid();
876 *   void
877 *   make_2d_quarter_plate_with_hole(Triangulation<2> &tria_2d,
878 *   const double half_length,
879 *   const double half_width,
880 *   const double hole_radius,
881 *   const unsigned int n_repetitions_xy = 1,
882 *   const double hole_division_fraction = 0.25);
883 *   void
884 *   setup_system(LA::MPI::BlockVector &solution_delta);
885 *   void
886 *   determine_component_extractors();
887 *   void
888 *   assemble_system(const LA::MPI::BlockVector &solution_delta);
889 *   void
890 *   assemble_system_one_cell(const typename DoFHandler<dim>::active_cell_iterator &cell,
891 *   ScratchData_ASM &scratch,
892 *   PerTaskData_ASM &data) const;
893 *   void
894 *   copy_local_to_global_system(const PerTaskData_ASM &data);
895 *   void
896 *   make_constraints(const int &it_nr);
897 *   void
898 *   setup_qph();
899 *   void
900 *   solve_nonlinear_timestep(LA::MPI::BlockVector &solution_delta);
901 *   std::pair<unsigned int, double>
902 *   solve_linear_system(LA::MPI::BlockVector &newton_update);
903 *   LA::MPI::BlockVector
904 *   get_solution_total(const LA::MPI::BlockVector &solution_delta) const;
905 *   void
906 *   update_end_timestep();
907 *   void
908 *   output_results(const unsigned int timestep,
909 *   const double current_time) const;
910 *   void
911 *   compute_vertex_positions(std::vector<double> &real_time,
912 *   std::vector<std::vector<Point<dim> > > &tracked_vertices,
913 *   const LA::MPI::BlockVector &solution_total) const;
914 *  
915 * @endcode
916 *
917 * Parallel communication
918 *
919 * @code
920 *   MPI_Comm mpi_communicator;
921 *   const unsigned int n_mpi_processes;
922 *   const unsigned int this_mpi_process;
923 *   mutable ConditionalOStream pcout;
924 *  
925 *   Parameters::AllParameters parameters;
927 *   Time time;
928 *   mutable TimerOutput timer;
930 *   PointHistory<dim> > quadrature_point_history;
931 *   const unsigned int degree;
932 *   const FESystem<dim> fe;
933 *   DoFHandler<dim> dof_handler;
934 *   const unsigned int dofs_per_cell;
935 *   const FEValuesExtractors::Vector u_fe;
936 *   const FEValuesExtractors::Scalar p_fe;
937 *   const FEValuesExtractors::Scalar J_fe;
938 *   static const unsigned int n_blocks = 3;
939 *   static const unsigned int n_components = dim + 2;
940 *   static const unsigned int first_u_component = 0;
941 *   static const unsigned int p_component = dim;
942 *   static const unsigned int J_component = dim + 1;
943 *   enum
944 *   {
945 *   u_block = 0,
946 *   p_block = 1,
947 *   J_block = 2
948 *   };
949 * @endcode
950 *
951 * Block data
952 *
953 * @code
954 *   std::vector<unsigned int> block_component;
955 *  
956 * @endcode
957 *
958 * DoF index data
959 *
960 * @code
961 *   std::vector<IndexSet> all_locally_owned_dofs;
962 *   IndexSet locally_owned_dofs;
963 *   IndexSet locally_relevant_dofs;
964 *   std::vector<IndexSet> locally_owned_partitioning;
965 *   std::vector<IndexSet> locally_relevant_partitioning;
966 *   std::vector<types::global_dof_index> dofs_per_block;
967 *   std::vector<types::global_dof_index> element_indices_u;
968 *   std::vector<types::global_dof_index> element_indices_p;
969 *   std::vector<types::global_dof_index> element_indices_J;
970 *   const QGauss<dim> qf_cell;
971 *   const QGauss<dim - 1> qf_face;
972 *   const unsigned int n_q_points;
973 *   const unsigned int n_q_points_f;
974 *   AffineConstraints<double> constraints;
975 *   LA::BlockSparseMatrix tangent_matrix;
976 *   LA::MPI::BlockVector system_rhs;
977 *   LA::MPI::BlockVector solution_n;
978 *   struct Errors
979 *   {
980 *   Errors()
981 *   :
982 *   norm(1.0), u(1.0), p(1.0), J(1.0)
983 *   {}
984 *   void reset()
985 *   {
986 *   norm = 1.0;
987 *   u = 1.0;
988 *   p = 1.0;
989 *   J = 1.0;
990 *   }
991 *   void normalise(const Errors &rhs)
992 *   {
993 *   if (rhs.norm != 0.0)
994 *   norm /= rhs.norm;
995 *   if (rhs.u != 0.0)
996 *   u /= rhs.u;
997 *   if (rhs.p != 0.0)
998 *   p /= rhs.p;
999 *   if (rhs.J != 0.0)
1000 *   J /= rhs.J;
1001 *   }
1002 *   double norm, u, p, J;
1003 *   };
1004 *   Errors error_residual, error_residual_0, error_residual_norm, error_update,
1005 *   error_update_0, error_update_norm;
1006 *   void
1007 *   get_error_residual(Errors &error_residual);
1008 *   void
1009 *   get_error_update(const LA::MPI::BlockVector &newton_update,
1010 *   Errors &error_update);
1011 *   std::pair<double, std::pair<double,double> >
1012 *   get_error_dilation(const LA::MPI::BlockVector &solution_total) const;
1013 *   void
1014 *   print_conv_header();
1015 *   void
1016 *   print_conv_footer(const LA::MPI::BlockVector &solution_delta);
1017 *   };
1018 *   template <int dim>
1019 *   Solid<dim>::Solid(const std::string &input_file)
1020 *   :
1021 *   mpi_communicator(MPI_COMM_WORLD),
1022 *   n_mpi_processes (Utilities::MPI::n_mpi_processes(mpi_communicator)),
1023 *   this_mpi_process (Utilities::MPI::this_mpi_process(mpi_communicator)),
1024 *   pcout(std::cout, this_mpi_process == 0),
1025 *   parameters(input_file),
1027 *   time(parameters.end_time, parameters.delta_t),
1028 *   timer(mpi_communicator,
1029 *   pcout,
1032 *   degree(parameters.poly_degree),
1033 *   fe(FE_Q<dim>(parameters.poly_degree), dim, // displacement
1034 *   FE_DGPMonomial<dim>(parameters.poly_degree - 1), 1, // pressure
1035 *   FE_DGPMonomial<dim>(parameters.poly_degree - 1), 1), // dilatation
1036 *   dof_handler(triangulation),
1037 *   dofs_per_cell (fe.dofs_per_cell),
1038 *   u_fe(first_u_component),
1039 *   p_fe(p_component),
1040 *   J_fe(J_component),
1041 *   dofs_per_block(n_blocks),
1042 *   qf_cell(parameters.quad_order),
1043 *   qf_face(parameters.quad_order),
1044 *   n_q_points (qf_cell.size()),
1045 *   n_q_points_f (qf_face.size())
1046 *   {
1047 *   Assert(dim==2 || dim==3, ExcMessage("This problem only works in 2 or 3 space dimensions."));
1048 *   determine_component_extractors();
1049 *   }
1050 *   template <int dim>
1051 *   Solid<dim>::~Solid()
1052 *   {
1053 *   dof_handler.clear();
1054 *   }
1055 *   template <int dim>
1056 *   void Solid<dim>::run()
1057 *   {
1058 *   LA::MPI::BlockVector solution_delta;
1059 *  
1060 *   make_grid();
1061 *   setup_system(solution_delta);
1062 *   {
1063 *   AffineConstraints<double> constraints;
1064 *   constraints.close();
1066 *   J_mask (J_component, n_components);
1067 *   VectorTools::project (dof_handler,
1068 *   constraints,
1069 *   QGauss<dim>(degree+2),
1070 *   J_mask,
1071 *   solution_n);
1072 *   }
1073 *   output_results(time.get_timestep(), time.current());
1074 *   time.increment();
1075 *  
1076 * @endcode
1077 *
1078 * Some points for post-processing
1079 *
1080 * @code
1081 *   std::vector<double> real_time;
1082 *   real_time.push_back(0);
1083 *   std::vector<std::vector<Point<dim> > > tracked_vertices (4);
1084 *   {
1085 *   Point<dim> p;
1086 *   p[1] = parameters.length/2.0;
1087 *   tracked_vertices[0].push_back(p*parameters.scale);
1088 *   }
1089 *   {
1090 *   Point<dim> p;
1091 *   p[1] = parameters.hole_diameter/2.0;
1092 *   tracked_vertices[1].push_back(p*parameters.scale);
1093 *   }
1094 *   {
1095 *   Point<dim> p;
1096 *   p[0] = parameters.hole_diameter/2.0;
1097 *   tracked_vertices[2].push_back(p*parameters.scale);
1098 *   }
1099 *   {
1100 *   Point<dim> p;
1101 *   p[0] = parameters.width/2.0;
1102 *   tracked_vertices[3].push_back(p*parameters.scale);
1103 *   }
1104 *  
1105 *   while (time.current() < time.end()+0.01*time.get_delta_t())
1106 *   {
1107 *   solve_nonlinear_timestep(solution_delta);
1108 *   solution_n += solution_delta;
1109 *   solution_delta = 0.0;
1110 *   output_results(time.get_timestep(), time.current());
1111 *   compute_vertex_positions(real_time,
1112 *   tracked_vertices,
1113 *   get_solution_total(solution_delta));
1114 *   update_end_timestep();
1115 *   time.increment();
1116 *   }
1117 *  
1118 *   pcout << "\n\n*** Spatial position history for tracked vertices ***" << std::endl;
1119 *   for (unsigned int t=0; t<real_time.size(); ++t)
1120 *   {
1121 *   if (t == 0)
1122 *   {
1123 *   pcout << "Time,";
1124 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
1125 *   {
1126 *   for (unsigned int d=0; d<dim; ++d)
1127 *   {
1128 *   pcout << "Point " << p << " [" << d << "]";
1129 *   if (!(p == tracked_vertices.size()-1 && d == dim-1))
1130 *   pcout << ",";
1131 *   }
1132 *   }
1133 *   pcout << std::endl;
1134 *   }
1135 *  
1136 *   pcout << std::setprecision(6);
1137 *   pcout << real_time[t] << ",";
1138 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
1139 *   {
1140 *   Assert(tracked_vertices[p].size() == real_time.size(),
1141 *   ExcMessage("Vertex not tracked at each timestep"));
1142 *   for (unsigned int d=0; d<dim; ++d)
1143 *   {
1144 *   pcout << tracked_vertices[p][t][d];
1145 *   if (!(p == tracked_vertices.size()-1 && d == dim-1))
1146 *   pcout << ",";
1147 *   }
1148 *   }
1149 *   pcout << std::endl;
1150 *   }
1151 *   }
1152 *   template <int dim>
1153 *   struct Solid<dim>::PerTaskData_ASM
1154 *   {
1156 *   Vector<double> cell_rhs;
1157 *   std::vector<types::global_dof_index> local_dof_indices;
1158 *   PerTaskData_ASM(const unsigned int dofs_per_cell)
1159 *   :
1160 *   cell_matrix(dofs_per_cell, dofs_per_cell),
1161 *   cell_rhs(dofs_per_cell),
1162 *   local_dof_indices(dofs_per_cell)
1163 *   {}
1164 *   void reset()
1165 *   {
1166 *   cell_matrix = 0.0;
1167 *   cell_rhs = 0.0;
1168 *   }
1169 *   };
1170 *   template <int dim>
1171 *   struct Solid<dim>::ScratchData_ASM
1172 *   {
1173 *   const LA::MPI::BlockVector &solution_total;
1174 *  
1175 * @endcode
1176 *
1177 * Integration helper
1178 *
1179 * @code
1180 *   FEValues<dim> fe_values_ref;
1181 *   FEFaceValues<dim> fe_face_values_ref;
1182 *  
1183 * @endcode
1184 *
1185 * Quadrature point solution
1186 *
1187 * @code
1188 *   std::vector<Tensor<2, dim> > solution_grads_u_total;
1189 *   std::vector<double> solution_values_p_total;
1190 *   std::vector<double> solution_values_J_total;
1191 *  
1192 * @endcode
1193 *
1194 * Shape function values and gradients
1195 *
1196 * @code
1197 *   std::vector<std::vector<double> > Nx;
1198 *   std::vector<std::vector<Tensor<2, dim> > > grad_Nx;
1199 *   std::vector<std::vector<SymmetricTensor<2, dim> > > symm_grad_Nx;
1200 *  
1201 *   ScratchData_ASM(const FiniteElement<dim> &fe_cell,
1202 *   const QGauss<dim> &qf_cell, const UpdateFlags uf_cell,
1203 *   const QGauss<dim - 1> & qf_face, const UpdateFlags uf_face,
1204 *   const LA::MPI::BlockVector &solution_total)
1205 *   :
1206 *   solution_total (solution_total),
1207 *   fe_values_ref(fe_cell, qf_cell, uf_cell),
1208 *   fe_face_values_ref(fe_cell, qf_face, uf_face),
1209 *   solution_grads_u_total(qf_cell.size()),
1210 *   solution_values_p_total(qf_cell.size()),
1211 *   solution_values_J_total(qf_cell.size()),
1212 *   Nx(qf_cell.size(),
1213 *   std::vector<double>(fe_cell.dofs_per_cell)),
1214 *   grad_Nx(qf_cell.size(),
1215 *   std::vector<Tensor<2, dim> >(fe_cell.dofs_per_cell)),
1216 *   symm_grad_Nx(qf_cell.size(),
1217 *   std::vector<SymmetricTensor<2, dim> >
1218 *   (fe_cell.dofs_per_cell))
1219 *   {}
1220 *   ScratchData_ASM(const ScratchData_ASM &rhs)
1221 *   :
1222 *   solution_total (rhs.solution_total),
1223 *   fe_values_ref(rhs.fe_values_ref.get_fe(),
1224 *   rhs.fe_values_ref.get_quadrature(),
1225 *   rhs.fe_values_ref.get_update_flags()),
1226 *   fe_face_values_ref(rhs.fe_face_values_ref.get_fe(),
1227 *   rhs.fe_face_values_ref.get_quadrature(),
1228 *   rhs.fe_face_values_ref.get_update_flags()),
1229 *   solution_grads_u_total(rhs.solution_grads_u_total),
1230 *   solution_values_p_total(rhs.solution_values_p_total),
1231 *   solution_values_J_total(rhs.solution_values_J_total),
1232 *   Nx(rhs.Nx),
1233 *   grad_Nx(rhs.grad_Nx),
1234 *   symm_grad_Nx(rhs.symm_grad_Nx)
1235 *   {}
1236 *   void reset()
1237 *   {
1238 *   const unsigned int n_q_points = solution_grads_u_total.size();
1239 *   const unsigned int n_dofs_per_cell = Nx[0].size();
1240 *  
1241 *   Assert(solution_grads_u_total.size() == n_q_points,
1242 *   ExcInternalError());
1243 *   Assert(solution_values_p_total.size() == n_q_points,
1244 *   ExcInternalError());
1245 *   Assert(solution_values_J_total.size() == n_q_points,
1246 *   ExcInternalError());
1247 *   Assert(Nx.size() == n_q_points,
1248 *   ExcInternalError());
1249 *   Assert(grad_Nx.size() == n_q_points,
1250 *   ExcInternalError());
1251 *   Assert(symm_grad_Nx.size() == n_q_points,
1252 *   ExcInternalError());
1253 *  
1254 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1255 *   {
1256 *   Assert( Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
1257 *   Assert( grad_Nx[q_point].size() == n_dofs_per_cell,
1258 *   ExcInternalError());
1259 *   Assert( symm_grad_Nx[q_point].size() == n_dofs_per_cell,
1260 *   ExcInternalError());
1261 *  
1262 *   solution_grads_u_total[q_point] = 0.0;
1263 *   solution_values_p_total[q_point] = 0.0;
1264 *   solution_values_J_total[q_point] = 0.0;
1265 *   for (unsigned int k = 0; k < n_dofs_per_cell; ++k)
1266 *   {
1267 *   Nx[q_point][k] = 0.0;
1268 *   grad_Nx[q_point][k] = 0.0;
1269 *   symm_grad_Nx[q_point][k] = 0.0;
1270 *   }
1271 *   }
1272 *   }
1273 *   };
1274 *   template <>
1275 *   void Solid<2>::make_grid()
1276 *   {
1277 *   const int dim = 2;
1278 *   const double tol = 1e-12;
1279 *   make_2d_quarter_plate_with_hole(triangulation,
1280 *   parameters.length/2.0,
1281 *   parameters.width/2.0,
1282 *   parameters.hole_diameter/2.0,
1283 *   parameters.n_repetitions_xy,
1284 *   parameters.hole_division_fraction);
1285 *  
1286 * @endcode
1287 *
1288 * Clear boundary ID's
1289 *
1290 * @code
1291 *   for (typename Triangulation<dim>::active_cell_iterator
1292 *   cell = triangulation.begin_active();
1293 *   cell != triangulation.end(); ++cell)
1294 *   {
1295 *   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1296 *   if (cell->face(face)->at_boundary())
1297 *   {
1298 *   cell->face(face)->set_all_boundary_ids(0);
1299 *   }
1300 *   }
1301 *  
1302 * @endcode
1303 *
1304 * Set boundary IDs and and manifolds
1305 *
1306 * @code
1307 *   const Point<dim> centre (0,0);
1308 *   for (typename Triangulation<dim>::active_cell_iterator
1309 *   cell = triangulation.begin_active();
1310 *   cell != triangulation.end(); ++cell)
1311 *   {
1312 *   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1313 *   if (cell->face(face)->at_boundary())
1314 *   {
1315 * @endcode
1316 *
1317 * Set boundary IDs
1318 *
1319 * @code
1320 *   if (std::abs(cell->face(face)->center()[0] - 0.0) < tol)
1321 *   {
1322 *   cell->face(face)->set_boundary_id(parameters.boundary_id_minus_X);
1323 *   }
1324 *   else if (std::abs(cell->face(face)->center()[0] - parameters.width/2.0) < tol)
1325 *   {
1326 *   cell->face(face)->set_boundary_id(parameters.boundary_id_plus_X);
1327 *   }
1328 *   else if (std::abs(cell->face(face)->center()[1] - 0.0) < tol)
1329 *   {
1330 *   cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Y);
1331 *   }
1332 *   else if (std::abs(cell->face(face)->center()[1] - parameters.length/2.0) < tol)
1333 *   {
1334 *   cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Y);
1335 *   }
1336 *   else
1337 *   {
1338 *   for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1339 *   if (std::abs(cell->vertex(vertex).distance(centre) - parameters.hole_diameter/2.0) < tol)
1340 *   {
1341 *   cell->face(face)->set_boundary_id(parameters.boundary_id_hole);
1342 *   break;
1343 *   }
1344 *   }
1345 *  
1346 * @endcode
1347 *
1348 * Set manifold IDs
1349 *
1350 * @code
1351 *   for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1352 *   if (std::abs(cell->vertex(vertex).distance(centre) - parameters.hole_diameter/2.0) < tol)
1353 *   {
1354 *   cell->face(face)->set_manifold_id(parameters.manifold_id_hole);
1355 *   break;
1356 *   }
1357 *   }
1358 *   }
1359 *   static SphericalManifold<dim> spherical_manifold (centre);
1360 *   triangulation.set_manifold(parameters.manifold_id_hole,spherical_manifold);
1361 *   triangulation.refine_global(parameters.global_refinement);
1362 *   GridTools::scale(parameters.scale,triangulation);
1363 *   }
1364 *   template <>
1365 *   void Solid<3>::make_grid()
1366 *   {
1367 *   const int dim = 3;
1368 *   const double tol = 1e-12;
1369 *   Triangulation<2> tria_2d;
1370 *   make_2d_quarter_plate_with_hole(tria_2d,
1371 *   parameters.length/2.0,
1372 *   parameters.width/2.0,
1373 *   parameters.hole_diameter/2.0,
1374 *   parameters.n_repetitions_xy,
1375 *   parameters.hole_division_fraction);
1376 *   GridGenerator::extrude_triangulation(tria_2d,
1377 *   parameters.n_repetitions_z+1,
1378 *   parameters.thickness/2.0,
1379 *   triangulation);
1380 *  
1381 * @endcode
1382 *
1383 * Clear boundary ID's
1384 *
1385 * @code
1387 *   cell = triangulation.begin_active();
1388 *   cell != triangulation.end(); ++cell)
1389 *   {
1390 *   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1391 *   if (cell->face(face)->at_boundary())
1392 *   {
1393 *   cell->face(face)->set_all_boundary_ids(0);
1394 *   }
1395 *   }
1396 *  
1397 * @endcode
1398 *
1399 * Set boundary IDs and and manifolds
1400 *
1401 * @code
1402 *   const Point<dim> direction (0,0,1);
1403 *   const Point<dim> centre (0,0,0);
1405 *   cell = triangulation.begin_active();
1406 *   cell != triangulation.end(); ++cell)
1407 *   {
1408 *   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1409 *   if (cell->face(face)->at_boundary())
1410 *   {
1411 * @endcode
1412 *
1413 * Set boundary IDs
1414 *
1415 * @code
1416 *   if (std::abs(cell->face(face)->center()[0] - 0.0) < tol)
1417 *   {
1418 *   cell->face(face)->set_boundary_id(parameters.boundary_id_minus_X);
1419 *   }
1420 *   else if (std::abs(cell->face(face)->center()[0] - parameters.width/2.0) < tol)
1421 *   {
1422 *   cell->face(face)->set_boundary_id(parameters.boundary_id_plus_X);
1423 *   }
1424 *   else if (std::abs(cell->face(face)->center()[1] - 0.0) < tol)
1425 *   {
1426 *   cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Y);
1427 *   }
1428 *   else if (std::abs(cell->face(face)->center()[1] - parameters.length/2.0) < tol)
1429 *   {
1430 *   cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Y);
1431 *   }
1432 *   else if (std::abs(cell->face(face)->center()[2] - 0.0) < tol)
1433 *   {
1434 *   cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Z);
1435 *   }
1436 *   else if (std::abs(cell->face(face)->center()[2] - parameters.thickness/2.0) < tol)
1437 *   {
1438 *   cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Z);
1439 *   }
1440 *   else
1441 *   {
1442 *   for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1443 *   {
1444 * @endcode
1445 *
1446 * Project the cell vertex to the XY plane and
1447 * test the distance from the cylinder axis
1448 *
1449 * @code
1450 *   Point<dim> vertex_proj = cell->vertex(vertex);
1451 *   vertex_proj[2] = 0.0;
1452 *   if (std::abs(vertex_proj.distance(centre) - parameters.hole_diameter/2.0) < tol)
1453 *   {
1454 *   cell->face(face)->set_boundary_id(parameters.boundary_id_hole);
1455 *   break;
1456 *   }
1457 *   }
1458 *   }
1459 *  
1460 * @endcode
1461 *
1462 * Set manifold IDs
1463 *
1464 * @code
1465 *   for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1466 *   {
1467 * @endcode
1468 *
1469 * Project the cell vertex to the XY plane and
1470 * test the distance from the cylinder axis
1471 *
1472 * @code
1473 *   Point<dim> vertex_proj = cell->vertex(vertex);
1474 *   vertex_proj[2] = 0.0;
1475 *   if (std::abs(vertex_proj.distance(centre) - parameters.hole_diameter/2.0) < 1e-12)
1476 *   {
1477 * @endcode
1478 *
1479 * Set manifold ID on face and edges
1480 *
1481 * @code
1482 *   cell->face(face)->set_all_manifold_ids(parameters.manifold_id_hole);
1483 *   break;
1484 *   }
1485 *   }
1486 *   }
1487 *   }
1488 *   static CylindricalManifold<dim> cylindrical_manifold (direction,centre);
1489 *   triangulation.set_manifold(parameters.manifold_id_hole,cylindrical_manifold);
1490 *   triangulation.refine_global(parameters.global_refinement);
1491 *   GridTools::scale(parameters.scale,triangulation);
1492 *   }
1493 *   template <int dim>
1494 *   void Solid<dim>::make_2d_quarter_plate_with_hole(Triangulation<2> &tria_2d,
1495 *   const double half_length,
1496 *   const double half_width,
1497 *   const double hole_radius,
1498 *   const unsigned int n_repetitions_xy,
1499 *   const double hole_division_fraction)
1500 *   {
1501 *   const double length = 2.0*half_length;
1502 *   const double width = 2.0*half_width;
1503 *   const double hole_diameter = 2.0*hole_radius;
1504 *  
1505 *   const double internal_width = hole_diameter + hole_division_fraction*(width - hole_diameter);
1506 *   Triangulation<2> tria_quarter_plate_hole;
1507 *   {
1508 *   Triangulation<2> tria_plate_hole;
1510 *   hole_diameter/2.0,
1511 *   internal_width/2.0);
1512 *  
1513 *   std::set<typename Triangulation<2>::active_cell_iterator > cells_to_remove;
1515 *   cell = tria_plate_hole.begin_active();
1516 *   cell != tria_plate_hole.end(); ++cell)
1517 *   {
1518 * @endcode
1519 *
1520 * Remove all cells that are not in the first quadrant
1521 *
1522 * @code
1523 *   if (cell->center()[0] < 0.0 || cell->center()[1] < 0.0)
1524 *   cells_to_remove.insert(cell);
1525 *   }
1526 *   Assert(cells_to_remove.size() > 0, ExcInternalError());
1527 *   Assert(cells_to_remove.size() != tria_plate_hole.n_active_cells(), ExcInternalError());
1528 *   GridGenerator::create_triangulation_with_removed_cells(tria_plate_hole,cells_to_remove,tria_quarter_plate_hole);
1529 *   }
1530 *  
1531 *   Triangulation<2> tria_cut_plate;
1532 *   {
1533 *   Triangulation<2> tria_plate;
1534 * @endcode
1535 *
1536 * Subdivide the plate so that we're left one
1537 * cell to remove (we'll replace this with the
1538 * plate with the hole) and then make the
1539 * rest of the subdivisions so that we're left
1540 * with cells with a decent aspect ratio
1541 *
1542 * @code
1543 *   std::vector<std::vector<double> > step_sizes;
1544 *   {
1545 *   std::vector<double> subdivision_width;
1546 *   subdivision_width.push_back(internal_width/2.0);
1547 *   const double width_remaining = (width - internal_width)/2.0;
1548 *   const unsigned int n_subs = static_cast<unsigned int>(std::max(1.0,std::ceil(width_remaining/(internal_width/2.0))));
1549 *   Assert(n_subs>0, ExcInternalError());
1550 *   for (unsigned int s=0; s<n_subs; ++s)
1551 *   subdivision_width.push_back(width_remaining/n_subs);
1552 *   step_sizes.push_back(subdivision_width);
1553 *  
1554 *   const double sum_half_width = std::accumulate(subdivision_width.begin(), subdivision_width.end(), 0.0);
1555 *   (void)sum_half_width;
1556 *   Assert(std::abs(sum_half_width-width/2.0) < 1e-12, ExcInternalError());
1557 *   }
1558 *   {
1559 *   std::vector<double> subdivision_length;
1560 *   subdivision_length.push_back(internal_width/2.0);
1561 *   const double length_remaining = (length - internal_width)/2.0;
1562 *   const unsigned int n_subs = static_cast<unsigned int>(std::max(1.0,std::ceil(length_remaining/(internal_width/2.0))));
1563 *   Assert(n_subs>0, ExcInternalError());
1564 *   for (unsigned int s=0; s<n_subs; ++s)
1565 *   subdivision_length.push_back(length_remaining/n_subs);
1566 *   step_sizes.push_back(subdivision_length);
1567 *  
1568 *   const double sum_half_length = std::accumulate(subdivision_length.begin(), subdivision_length.end(), 0.0);
1569 *   (void)sum_half_length;
1570 *   Assert(std::abs(sum_half_length-length/2.0) < 1e-12, ExcInternalError());
1571 *   }
1572 *  
1573 *   GridGenerator::subdivided_hyper_rectangle(tria_plate,
1574 *   step_sizes,
1575 *   Point<2>(0.0, 0.0),
1576 *   Point<2>(width/2.0, length/2.0));
1577 *  
1578 *   std::set<typename Triangulation<2>::active_cell_iterator > cells_to_remove;
1579 *   for (typename Triangulation<2>::active_cell_iterator
1580 *   cell = tria_plate.begin_active();
1581 *   cell != tria_plate.end(); ++cell)
1582 *   {
1583 * @endcode
1584 *
1585 * Remove all cells that are in the first quadrant
1586 *
1587 * @code
1588 *   if (cell->center()[0] < internal_width/2.0 && cell->center()[1] < internal_width/2.0)
1589 *   cells_to_remove.insert(cell);
1590 *   }
1591 *   Assert(cells_to_remove.size() > 0, ExcInternalError());
1592 *   Assert(cells_to_remove.size() != tria_plate.n_active_cells(), ExcInternalError());
1593 *   GridGenerator::create_triangulation_with_removed_cells(tria_plate,cells_to_remove,tria_cut_plate);
1594 *   }
1595 *  
1596 *   Triangulation<2> tria_2d_not_flat;
1597 *   GridGenerator::merge_triangulations(tria_quarter_plate_hole,
1598 *   tria_cut_plate,
1599 *   tria_2d_not_flat);
1600 *  
1601 * @endcode
1602 *
1603 * Attach a manifold to the curved boundary and refine
1604 * Note: We can only guarantee that the vertices sit on
1605 * the curve, so we must test with their position instead
1606 * of the cell centre.
1607 *
1608 * @code
1609 *   const Point<2> centre_2d (0,0);
1610 *   for (typename Triangulation<2>::active_cell_iterator
1611 *   cell = tria_2d_not_flat.begin_active();
1612 *   cell != tria_2d_not_flat.end(); ++cell)
1613 *   {
1614 *   for (unsigned int face=0; face<GeometryInfo<2>::faces_per_cell; ++face)
1615 *   if (cell->face(face)->at_boundary())
1616 *   for (unsigned int vertex=0; vertex<GeometryInfo<2>::vertices_per_face; ++vertex)
1617 *   if (std::abs(cell->vertex(vertex).distance(centre_2d) - hole_diameter/2.0) < 1e-12)
1618 *   {
1619 *   cell->face(face)->set_manifold_id(10);
1620 *   break;
1621 *   }
1622 *   }
1623 *   SphericalManifold<2> spherical_manifold_2d (centre_2d);
1624 *   tria_2d_not_flat.set_manifold(10,spherical_manifold_2d);
1625 *   tria_2d_not_flat.refine_global(std::max (1U, n_repetitions_xy));
1626 *   tria_2d_not_flat.reset_manifold(10); // Clear manifold
1627 *  
1628 *   GridGenerator::flatten_triangulation(tria_2d_not_flat,tria_2d);
1629 *   }
1630 *   template <int dim>
1631 *   void Solid<dim>::setup_system(LA::MPI::BlockVector &solution_delta)
1632 *   {
1633 *   timer.enter_subsection("Setup system");
1634 *   pcout << "Setting up linear system..." << std::endl;
1635 *  
1636 * @endcode
1637 *
1638 * Partition triangulation
1639 *
1640 * @code
1641 *   GridTools::partition_triangulation (n_mpi_processes,
1642 *   triangulation);
1643 *  
1644 *   block_component = std::vector<unsigned int> (n_components, u_block); // Displacement
1645 *   block_component[p_component] = p_block; // Pressure
1646 *   block_component[J_component] = J_block; // Dilatation
1647 *   dof_handler.distribute_dofs(fe);
1648 *   DoFRenumbering::Cuthill_McKee(dof_handler);
1649 *   DoFRenumbering::component_wise(dof_handler, block_component);
1650 *  
1651 * @endcode
1652 *
1653 * Count DoFs in each block
1654 *
1655 * @code
1656 *   dofs_per_block = DoFTools::count_dofs_per_fe_block(dof_handler, block_component);
1657 *  
1658 *   all_locally_owned_dofs = DoFTools::locally_owned_dofs_per_subdomain (dof_handler);
1659 *   std::vector<IndexSet> all_locally_relevant_dofs
1660 *   = DoFTools::locally_relevant_dofs_per_subdomain (dof_handler);
1661 *  
1662 *   locally_owned_dofs.clear();
1663 *   locally_owned_partitioning.clear();
1664 *   Assert(all_locally_owned_dofs.size() > this_mpi_process, ExcInternalError());
1665 *   locally_owned_dofs = all_locally_owned_dofs[this_mpi_process];
1666 *  
1667 *   locally_relevant_dofs.clear();
1668 *   locally_relevant_partitioning.clear();
1669 *   Assert(all_locally_relevant_dofs.size() > this_mpi_process, ExcInternalError());
1670 *   locally_relevant_dofs = all_locally_relevant_dofs[this_mpi_process];
1671 *  
1672 *   locally_owned_partitioning.reserve(n_blocks);
1673 *   locally_relevant_partitioning.reserve(n_blocks);
1674 *   for (unsigned int b=0; b<n_blocks; ++b)
1675 *   {
1676 *   const types::global_dof_index idx_begin
1677 *   = std::accumulate(dofs_per_block.begin(),
1678 *   std::next(dofs_per_block.begin(),b), 0);
1679 *   const types::global_dof_index idx_end
1680 *   = std::accumulate(dofs_per_block.begin(),
1681 *   std::next(dofs_per_block.begin(),b+1), 0);
1682 *   locally_owned_partitioning.push_back(locally_owned_dofs.get_view(idx_begin, idx_end));
1683 *   locally_relevant_partitioning.push_back(locally_relevant_dofs.get_view(idx_begin, idx_end));
1684 *   }
1685 *  
1686 *   pcout
1687 *   << " Number of active cells: " << triangulation.n_active_cells()
1688 *   << " (by partition:";
1689 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
1690 *   pcout
1691 *   << (p==0 ? ' ' : '+')
1692 *   << (GridTools::count_cells_with_subdomain_association (triangulation,p));
1693 *   pcout << ")" << std::endl;
1694 *  
1695 *   pcout
1696 *   << " Number of degrees of freedom: " << dof_handler.n_dofs()
1697 *   << " (by partition:";
1698 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
1699 *   pcout
1700 *   << (p==0 ? ' ' : '+')
1701 *   << (DoFTools::count_dofs_with_subdomain_association (dof_handler,p));
1702 *   pcout << ")" << std::endl;
1703 *   pcout
1704 *   << " Number of degrees of freedom per block: "
1705 *   << "[n_u, n_p, n_J] = ["
1706 *   << dofs_per_block[u_block] << ", "
1707 *   << dofs_per_block[p_block] << ", "
1708 *   << dofs_per_block[J_block] << "]"
1709 *   << std::endl;
1710 *  
1711 *  
1712 *   Table<2, DoFTools::Coupling> coupling(n_components, n_components);
1713 *   for (unsigned int ii = 0; ii < n_components; ++ii)
1714 *   for (unsigned int jj = 0; jj < n_components; ++jj)
1715 *   if (((ii < p_component) && (jj == J_component))
1716 *   || ((ii == J_component) && (jj < p_component))
1717 *   || ((ii == p_component) && (jj == p_component)))
1718 *   coupling[ii][jj] = DoFTools::none;
1719 *   else
1720 *   coupling[ii][jj] = DoFTools::always;
1721 *  
1722 *   TrilinosWrappers::BlockSparsityPattern bsp (locally_owned_partitioning,
1723 *   locally_owned_partitioning,
1724 *   locally_relevant_partitioning,
1725 *   mpi_communicator);
1726 *   DoFTools::make_sparsity_pattern (dof_handler, bsp,
1727 *   constraints, false,
1728 *   this_mpi_process);
1729 *   bsp.compress();
1730 *   tangent_matrix.reinit (bsp);
1731 *  
1732 * @endcode
1733 *
1734 * We then set up storage vectors
1735 *
1736 * @code
1737 *   system_rhs.reinit(locally_owned_partitioning,
1738 *   mpi_communicator);
1739 *   solution_n.reinit(locally_owned_partitioning,
1740 *   mpi_communicator);
1741 *   solution_delta.reinit(locally_owned_partitioning,
1742 *   mpi_communicator);
1743 *   setup_qph();
1744 *   timer.leave_subsection();
1745 *   }
1746 *   template <int dim>
1747 *   void
1748 *   Solid<dim>::determine_component_extractors()
1749 *   {
1750 *   element_indices_u.clear();
1751 *   element_indices_p.clear();
1752 *   element_indices_J.clear();
1753 *   for (unsigned int k = 0; k < fe.dofs_per_cell; ++k)
1754 *   {
1755 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
1756 *   if (k_group == u_block)
1757 *   element_indices_u.push_back(k);
1758 *   else if (k_group == p_block)
1759 *   element_indices_p.push_back(k);
1760 *   else if (k_group == J_block)
1761 *   element_indices_J.push_back(k);
1762 *   else
1763 *   {
1764 *   Assert(k_group <= J_block, ExcInternalError());
1765 *   }
1766 *   }
1767 *   }
1768 *   template <int dim>
1769 *   void Solid<dim>::setup_qph()
1770 *   {
1771 *   pcout << "Setting up quadrature point data..." << std::endl;
1772 *   quadrature_point_history.initialize(triangulation.begin_active(),
1773 *   triangulation.end(),
1774 *   n_q_points);
1775 *   FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1776 *   cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1777 *   dof_handler.begin_active()),
1778 *   endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1779 *   dof_handler.end());
1780 *   for (; cell!=endc; ++cell)
1781 *   {
1782 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1783 *   const std::vector<std::shared_ptr<PointHistory<dim> > > lqph =
1784 *   quadrature_point_history.get_data(cell);
1785 *   Assert(lqph.size() == n_q_points, ExcInternalError());
1786 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1787 *   lqph[q_point]->setup_lqp(parameters, time);
1788 *   }
1789 *   }
1790 *   template <int dim>
1791 *   void
1792 *   Solid<dim>::solve_nonlinear_timestep(LA::MPI::BlockVector &solution_delta)
1793 *   {
1794 *   pcout << std::endl
1795 *   << "Timestep " << time.get_timestep() << " @ "
1796 *   << time.current() << "s of "
1797 *   << time.end() << "s" << std::endl;
1798 *   LA::MPI::BlockVector newton_update(locally_owned_partitioning,
1799 *   mpi_communicator);
1800 *   error_residual.reset();
1801 *   error_residual_0.reset();
1802 *   error_residual_norm.reset();
1803 *   error_update.reset();
1804 *   error_update_0.reset();
1805 *   error_update_norm.reset();
1806 *   print_conv_header();
1807 *   unsigned int newton_iteration = 0;
1808 *   for (; newton_iteration < parameters.max_iterations_NR;
1809 *   ++newton_iteration)
1810 *   {
1811 *   pcout << " " << std::setw(2) << newton_iteration << " " << std::flush;
1812 *   make_constraints(newton_iteration);
1813 *   assemble_system(solution_delta);
1814 *   get_error_residual(error_residual);
1815 *   if (newton_iteration == 0)
1816 *   error_residual_0 = error_residual;
1817 *   error_residual_norm = error_residual;
1818 *   error_residual_norm.normalise(error_residual_0);
1819 *   if (newton_iteration > 0 &&
1820 *   (error_update_norm.u <= parameters.tol_u &&
1821 *   error_residual_norm.u <= parameters.tol_f) )
1822 *   {
1823 *   pcout << " CONVERGED! " << std::endl;
1824 *   print_conv_footer(solution_delta);
1825 *   break;
1826 *   }
1827 *   const std::pair<unsigned int, double>
1828 *   lin_solver_output = solve_linear_system(newton_update);
1829 *   get_error_update(newton_update, error_update);
1830 *   if (newton_iteration == 0)
1831 *   error_update_0 = error_update;
1832 *   error_update_norm = error_update;
1833 *   error_update_norm.normalise(error_update_0);
1834 *   solution_delta += newton_update;
1835 *   newton_update = 0.0;
1836 *   pcout << " | " << std::fixed << std::setprecision(3) << std::setw(7)
1837 *   << std::scientific << lin_solver_output.first << " "
1838 *   << lin_solver_output.second << " " << error_residual_norm.norm
1839 *   << " " << error_residual_norm.u << " "
1840 *   << error_residual_norm.p << " " << error_residual_norm.J
1841 *   << " " << error_update_norm.norm << " " << error_update_norm.u
1842 *   << " " << error_update_norm.p << " " << error_update_norm.J
1843 *   << " " << std::endl;
1844 *   }
1845 *   AssertThrow (newton_iteration <= parameters.max_iterations_NR,
1846 *   ExcMessage("No convergence in nonlinear solver!"));
1847 *   }
1848 *   template <int dim>
1849 *   void Solid<dim>::print_conv_header()
1850 *   {
1851 *   pcout << std::string(132,'_') << std::endl;
1852 *   pcout << " SOLVER STEP "
1853 *   << " | LIN_IT LIN_RES RES_NORM "
1854 *   << " RES_U RES_P RES_J NU_NORM "
1855 *   << " NU_U NU_P NU_J " << std::endl;
1856 *   pcout << std::string(132,'_') << std::endl;
1857 *   }
1858 *   template <int dim>
1859 *   void Solid<dim>::print_conv_footer(const LA::MPI::BlockVector &solution_delta)
1860 *   {
1861 *   pcout << std::string(132,'_') << std::endl;
1862 *   const std::pair<double,std::pair<double,double> > error_dil = get_error_dilation(get_solution_total(solution_delta));
1863 *   pcout << "Relative errors:" << std::endl
1864 *   << "Displacement:\t" << error_update.u / error_update_0.u << std::endl
1865 *   << "Force: \t\t" << error_residual.u / error_residual_0.u << std::endl
1866 *   << "Dilatation:\t" << error_dil.first << std::endl
1867 *   << "v / V_0:\t" << error_dil.second.second << " / " << error_dil.second.first
1868 *   << " = " << (error_dil.second.second/error_dil.second.first) << std::endl;
1869 *   }
1870 *   template <int dim>
1871 *   std::pair<double,std::pair<double,double> >
1872 *   Solid<dim>::get_error_dilation(const LA::MPI::BlockVector &solution_total) const
1873 *   {
1874 *   double vol_reference = 0.0;
1875 *   double vol_current = 0.0;
1876 *   double dil_L2_error = 0.0;
1877 *   FEValues<dim> fe_values_ref(fe, qf_cell,
1878 *   update_values | update_gradients | update_JxW_values);
1879 *   std::vector<Tensor<2, dim> > solution_grads_u_total (qf_cell.size());
1880 *   std::vector<double> solution_values_J_total (qf_cell.size());
1881 *   FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1882 *   cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1883 *   dof_handler.begin_active()),
1884 *   endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1885 *   dof_handler.end());
1886 *   for (; cell != endc; ++cell)
1887 *   {
1888 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1889 *   fe_values_ref.reinit(cell);
1890 *   fe_values_ref[u_fe].get_function_gradients(solution_total,
1891 *   solution_grads_u_total);
1892 *   fe_values_ref[J_fe].get_function_values(solution_total,
1893 *   solution_values_J_total);
1894 *   const std::vector<std::shared_ptr<const PointHistory<dim> > > lqph =
1895 *   quadrature_point_history.get_data(cell);
1896 *   Assert(lqph.size() == n_q_points, ExcInternalError());
1897 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1898 *   {
1899 *   const double det_F_qp = determinant(Physics::Elasticity::Kinematics::F(solution_grads_u_total[q_point]));
1900 *   const double J_tilde_qp = solution_values_J_total[q_point];
1901 *   const double the_error_qp_squared = std::pow((det_F_qp - J_tilde_qp),
1902 *   2);
1903 *   const double JxW = fe_values_ref.JxW(q_point);
1904 *   dil_L2_error += the_error_qp_squared * JxW;
1905 *   vol_reference += JxW;
1906 *   vol_current += det_F_qp * JxW;
1907 *   }
1908 *   }
1909 *   Assert(vol_current > 0.0, ExcInternalError());
1910 * @endcode
1911 *
1912 * Sum across all processors
1913 *
1914 * @code
1915 *   dil_L2_error = Utilities::MPI::sum(dil_L2_error,mpi_communicator);
1916 *   vol_reference = Utilities::MPI::sum(vol_reference,mpi_communicator);
1917 *   vol_current = Utilities::MPI::sum(vol_current,mpi_communicator);
1918 *  
1919 *   return std::make_pair(std::sqrt(dil_L2_error),
1920 *   std::make_pair(vol_reference,vol_current));
1921 *   }
1922 *   template <int dim>
1923 *   void Solid<dim>::get_error_residual(Errors &error_residual)
1924 *   {
1925 * @endcode
1926 *
1927 * Construct a residual vector that has the values for all of its
1928 * constrained DoFs set to zero.
1929 *
1930 * @code
1931 *   LA::MPI::BlockVector error_res (system_rhs);
1932 *   constraints.set_zero(error_res);
1933 *   error_residual.norm = error_res.l2_norm();
1934 *   error_residual.u = error_res.block(u_block).l2_norm();
1935 *   error_residual.p = error_res.block(p_block).l2_norm();
1936 *   error_residual.J = error_res.block(J_block).l2_norm();
1937 *   }
1938 *   template <int dim>
1939 *   void Solid<dim>::get_error_update(const LA::MPI::BlockVector &newton_update,
1940 *   Errors &error_update)
1941 *   {
1942 * @endcode
1943 *
1944 * Construct a update vector that has the values for all of its
1945 * constrained DoFs set to zero.
1946 *
1947 * @code
1948 *   LA::MPI::BlockVector error_ud (newton_update);
1949 *   constraints.set_zero(error_ud);
1950 *   error_update.norm = error_ud.l2_norm();
1951 *   error_update.u = error_ud.block(u_block).l2_norm();
1952 *   error_update.p = error_ud.block(p_block).l2_norm();
1953 *   error_update.J = error_ud.block(J_block).l2_norm();
1954 *   }
1955 *   template <int dim>
1956 *   LA::MPI::BlockVector
1957 *   Solid<dim>::get_solution_total(const LA::MPI::BlockVector &solution_delta) const
1958 *   {
1959 * @endcode
1960 *
1961 * Cell interpolation -> Ghosted vector
1962 *
1963 * @code
1964 *   LA::MPI::BlockVector solution_total (locally_owned_partitioning,
1965 *   locally_relevant_partitioning,
1966 *   mpi_communicator,
1967 *   /*vector_writable = */ false);
1968 *   LA::MPI::BlockVector tmp (solution_total);
1969 *   solution_total = solution_n;
1970 *   tmp = solution_delta;
1971 *   solution_total += tmp;
1972 *   return solution_total;
1973 *   }
1974 *   template <int dim>
1975 *   void Solid<dim>::assemble_system(const LA::MPI::BlockVector &solution_delta)
1976 *   {
1977 *   timer.enter_subsection("Assemble system");
1978 *   pcout << " ASM_SYS " << std::flush;
1979 *   tangent_matrix = 0.0;
1980 *   system_rhs = 0.0;
1981 *   const LA::MPI::BlockVector solution_total(get_solution_total(solution_delta));
1982 *   const UpdateFlags uf_cell(update_values |
1983 *   update_gradients |
1984 *   update_JxW_values);
1985 *   const UpdateFlags uf_face(update_values |
1986 *   update_normal_vectors |
1987 *   update_JxW_values);
1988 *   PerTaskData_ASM per_task_data(dofs_per_cell);
1989 *   ScratchData_ASM scratch_data(fe, qf_cell, uf_cell, qf_face, uf_face, solution_total);
1990 *  
1991 *   FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1992 *   cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1993 *   dof_handler.begin_active()),
1994 *   endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1995 *   dof_handler.end());
1996 *   for (; cell != endc; ++cell)
1997 *   {
1998 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1999 *   assemble_system_one_cell(cell, scratch_data, per_task_data);
2000 *   copy_local_to_global_system(per_task_data);
2001 *   }
2002 *   tangent_matrix.compress(VectorOperation::add);
2003 *   system_rhs.compress(VectorOperation::add);
2004 *   timer.leave_subsection();
2005 *   }
2006 *   template <int dim>
2007 *   void Solid<dim>::copy_local_to_global_system(const PerTaskData_ASM &data)
2008 *   {
2009 *   constraints.distribute_local_to_global(data.cell_matrix, data.cell_rhs,
2010 *   data.local_dof_indices,
2011 *   tangent_matrix, system_rhs);
2012 *   }
2013 *   template <int dim>
2014 *   void
2015 *   Solid<dim>::assemble_system_one_cell(const typename DoFHandler<dim>::active_cell_iterator &cell,
2016 *   ScratchData_ASM &scratch,
2017 *   PerTaskData_ASM &data) const
2018 *   {
2019 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2020 *  
2021 *   data.reset();
2022 *   scratch.reset();
2023 *   scratch.fe_values_ref.reinit(cell);
2024 *   cell->get_dof_indices(data.local_dof_indices);
2025 *   const std::vector<std::shared_ptr<const PointHistory<dim> > > lqph =
2026 *   quadrature_point_history.get_data(cell);
2027 *   Assert(lqph.size() == n_q_points, ExcInternalError());
2028 *  
2029 * @endcode
2030 *
2031 * Update quadrature point solution
2032 *
2033 * @code
2034 *   scratch.fe_values_ref[u_fe].get_function_gradients(scratch.solution_total,
2035 *   scratch.solution_grads_u_total);
2036 *   scratch.fe_values_ref[p_fe].get_function_values(scratch.solution_total,
2037 *   scratch.solution_values_p_total);
2038 *   scratch.fe_values_ref[J_fe].get_function_values(scratch.solution_total,
2039 *   scratch.solution_values_J_total);
2040 *  
2041 * @endcode
2042 *
2043 * Update shape functions and their gradients (push-forward)
2044 *
2045 * @code
2046 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2047 *   {
2048 *   const Tensor<2, dim> F = Physics::Elasticity::Kinematics::F(scratch.solution_grads_u_total[q_point]);
2049 *   const Tensor<2, dim> F_inv = invert(F);
2050 *  
2051 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
2052 *   {
2053 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
2054 *   if (k_group == u_block)
2055 *   {
2056 *   scratch.grad_Nx[q_point][k] = scratch.fe_values_ref[u_fe].gradient(k, q_point)
2057 *   * F_inv;
2058 *   scratch.symm_grad_Nx[q_point][k] = symmetrize(scratch.grad_Nx[q_point][k]);
2059 *   }
2060 *   else if (k_group == p_block)
2061 *   scratch.Nx[q_point][k] = scratch.fe_values_ref[p_fe].value(k,
2062 *   q_point);
2063 *   else if (k_group == J_block)
2064 *   scratch.Nx[q_point][k] = scratch.fe_values_ref[J_fe].value(k,
2065 *   q_point);
2066 *   else
2067 *   Assert(k_group <= J_block, ExcInternalError());
2068 *   }
2069 *   }
2070 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2071 *   {
2072 *   const SymmetricTensor<2, dim> &I = Physics::Elasticity::StandardTensors<dim>::I;
2073 *   const Tensor<2, dim> F = Physics::Elasticity::Kinematics::F(scratch.solution_grads_u_total[q_point]);
2074 *   const double det_F = determinant(F);
2075 *   const double &p_tilde = scratch.solution_values_p_total[q_point];
2076 *   const double &J_tilde = scratch.solution_values_J_total[q_point];
2077 *   Assert(det_F > 0, ExcInternalError());
2078 *  
2079 *   {
2080 *   PointHistory<dim> *lqph_q_point_nc = const_cast<PointHistory<dim>*>(lqph[q_point].get());
2081 *   lqph_q_point_nc->update_internal_equilibrium(F,p_tilde,J_tilde);
2082 *   }
2083 *  
2084 *   const SymmetricTensor<2, dim> tau = lqph[q_point]->get_tau(F,p_tilde);
2085 *   const Tensor<2, dim> tau_ns (tau);
2086 *   const SymmetricTensor<4, dim> Jc = lqph[q_point]->get_Jc(F,p_tilde);
2087 *   const double dPsi_vol_dJ = lqph[q_point]->get_dPsi_vol_dJ(J_tilde);
2088 *   const double d2Psi_vol_dJ2 = lqph[q_point]->get_d2Psi_vol_dJ2(J_tilde);
2089 *  
2090 *   const std::vector<double> &Nx = scratch.Nx[q_point];
2091 *   const std::vector<Tensor<2, dim> > &grad_Nx = scratch.grad_Nx[q_point];
2092 *   const std::vector<SymmetricTensor<2, dim> > &symm_grad_Nx = scratch.symm_grad_Nx[q_point];
2093 *   const double JxW = scratch.fe_values_ref.JxW(q_point);
2094 *  
2095 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2096 *   {
2097 *   const unsigned int component_i = fe.system_to_component_index(i).first;
2098 *   const unsigned int i_group = fe.system_to_base_index(i).first.first;
2099 *   if (i_group == u_block)
2100 *   data.cell_rhs(i) -= (symm_grad_Nx[i] * tau) * JxW;
2101 *   else if (i_group == p_block)
2102 *   data.cell_rhs(i) -= Nx[i] * (det_F - J_tilde) * JxW;
2103 *   else if (i_group == J_block)
2104 *   data.cell_rhs(i) -= Nx[i] * (dPsi_vol_dJ - p_tilde) * JxW;
2105 *   else
2106 *   Assert(i_group <= J_block, ExcInternalError());
2107 *  
2108 *   for (unsigned int j = 0; j <= i; ++j)
2109 *   {
2110 *   const unsigned int component_j = fe.system_to_component_index(j).first;
2111 *   const unsigned int j_group = fe.system_to_base_index(j).first.first;
2112 *   if ((i_group == u_block) && (j_group == u_block))
2113 *   {
2114 *   data.cell_matrix(i, j) += symm_grad_Nx[i] * Jc // The material contribution:
2115 *   * symm_grad_Nx[j] * JxW;
2116 *   if (component_i == component_j) // geometrical stress contribution
2117 *   data.cell_matrix(i, j) += grad_Nx[i][component_i] * tau_ns
2118 *   * grad_Nx[j][component_j] * JxW;
2119 *   }
2120 *   else if ((i_group == u_block) && (j_group == p_block))
2121 *   {
2122 *   data.cell_matrix(i, j) += (symm_grad_Nx[i] * I)
2123 *   * Nx[j] * det_F
2124 *   * JxW;
2125 *   }
2126 *   else if ((i_group == p_block) && (j_group == u_block))
2127 *   {
2128 *   data.cell_matrix(i, j) += Nx[i] * det_F
2129 *   * (symm_grad_Nx[j] * I)
2130 *   * JxW;
2131 *   }
2132 *   else if ((i_group == p_block) && (j_group == J_block))
2133 *   data.cell_matrix(i, j) -= Nx[i] * Nx[j] * JxW;
2134 *   else if ((i_group == J_block) && (j_group == p_block))
2135 *   data.cell_matrix(i, j) -= Nx[i] * Nx[j] * JxW;
2136 *   else if ((i_group == J_block) && (j_group == J_block))
2137 *   data.cell_matrix(i, j) += Nx[i] * d2Psi_vol_dJ2 * Nx[j] * JxW;
2138 *   else
2139 *   Assert((i_group <= J_block) && (j_group <= J_block),
2140 *   ExcInternalError());
2141 *   }
2142 *   }
2143 *   }
2144 *  
2145 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2146 *   for (unsigned int j = i + 1; j < dofs_per_cell; ++j)
2147 *   data.cell_matrix(i, j) = data.cell_matrix(j, i);
2148 *  
2149 *   if (parameters.driver == "Neumann")
2150 *   for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
2151 *   ++face)
2152 *   if (cell->face(face)->at_boundary() == true
2153 *   && cell->face(face)->boundary_id() == parameters.boundary_id_plus_Y)
2154 *   {
2155 *   scratch.fe_face_values_ref.reinit(cell, face);
2156 *   for (unsigned int f_q_point = 0; f_q_point < n_q_points_f;
2157 *   ++f_q_point)
2158 *   {
2159 *   const Tensor<1, dim> &N =
2160 *   scratch.fe_face_values_ref.normal_vector(f_q_point);
2161 *   static const double pressure_nom = parameters.pressure
2162 *   / (parameters.scale * parameters.scale);
2163 *   const double time_ramp = (time.current() < parameters.load_time ?
2164 *   time.current() / parameters.load_time : 1.0);
2165 *   const double pressure = -pressure_nom * time_ramp;
2166 *   const Tensor<1, dim> traction = pressure * N;
2167 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2168 *   {
2169 *   const unsigned int i_group =
2170 *   fe.system_to_base_index(i).first.first;
2171 *   if (i_group == u_block)
2172 *   {
2173 *   const unsigned int component_i =
2174 *   fe.system_to_component_index(i).first;
2175 *   const double Ni =
2176 *   scratch.fe_face_values_ref.shape_value(i,
2177 *   f_q_point);
2178 *   const double JxW = scratch.fe_face_values_ref.JxW(
2179 *   f_q_point);
2180 *   data.cell_rhs(i) += (Ni * traction[component_i])
2181 *   * JxW;
2182 *   }
2183 *   }
2184 *   }
2185 *   }
2186 *   }
2187 *   template <int dim>
2188 *   void Solid<dim>::make_constraints(const int &it_nr)
2189 *   {
2190 *   pcout << " CST " << std::flush;
2191 *   if (it_nr > 1)
2192 *   return;
2193 *   constraints.clear();
2194 *   const bool apply_dirichlet_bc = (it_nr == 0);
2195 *   const FEValuesExtractors::Scalar x_displacement(0);
2196 *   const FEValuesExtractors::Scalar y_displacement(1);
2197 *   {
2198 *   const int boundary_id = parameters.boundary_id_minus_X;
2199 *   if (apply_dirichlet_bc == true)
2200 *   VectorTools::interpolate_boundary_values(dof_handler,
2201 *   boundary_id,
2202 *   Functions::ZeroFunction<dim>(n_components),
2203 *   constraints,
2204 *   fe.component_mask(x_displacement));
2205 *   else
2206 *   VectorTools::interpolate_boundary_values(dof_handler,
2207 *   boundary_id,
2208 *   Functions::ZeroFunction<dim>(n_components),
2209 *   constraints,
2210 *   fe.component_mask(x_displacement));
2211 *   }
2212 *   {
2213 *   const int boundary_id = parameters.boundary_id_minus_Y;
2214 *   if (apply_dirichlet_bc == true)
2215 *   VectorTools::interpolate_boundary_values(dof_handler,
2216 *   boundary_id,
2217 *   Functions::ZeroFunction<dim>(n_components),
2218 *   constraints,
2219 *   fe.component_mask(y_displacement));
2220 *   else
2221 *   VectorTools::interpolate_boundary_values(dof_handler,
2222 *   boundary_id,
2223 *   Functions::ZeroFunction<dim>(n_components),
2224 *   constraints,
2225 *   fe.component_mask(y_displacement));
2226 *   }
2227 *   if (dim==3)
2228 *   {
2229 *   const FEValuesExtractors::Scalar z_displacement(2);
2230 *   {
2231 *   const int boundary_id = parameters.boundary_id_minus_Z;
2232 *   if (apply_dirichlet_bc == true)
2233 *   VectorTools::interpolate_boundary_values(dof_handler,
2234 *   boundary_id,
2235 *   Functions::ZeroFunction<dim>(n_components),
2236 *   constraints,
2237 *   fe.component_mask(z_displacement));
2238 *   else
2239 *   VectorTools::interpolate_boundary_values(dof_handler,
2240 *   boundary_id,
2241 *   Functions::ZeroFunction<dim>(n_components),
2242 *   constraints,
2243 *   fe.component_mask(z_displacement));
2244 *   }
2245 *   {
2246 *   const int boundary_id = parameters.boundary_id_plus_Z;
2247 *   if (apply_dirichlet_bc == true)
2248 *   VectorTools::interpolate_boundary_values(dof_handler,
2249 *   boundary_id,
2250 *   Functions::ZeroFunction<dim>(n_components),
2251 *   constraints,
2252 *   fe.component_mask(z_displacement));
2253 *   else
2254 *   VectorTools::interpolate_boundary_values(dof_handler,
2255 *   boundary_id,
2256 *   Functions::ZeroFunction<dim>(n_components),
2257 *   constraints,
2258 *   fe.component_mask(z_displacement));
2259 *   }
2260 *   }
2261 *   if (parameters.driver == "Dirichlet")
2262 *   {
2263 *   const int boundary_id = parameters.boundary_id_plus_Y;
2264 *   if (apply_dirichlet_bc == true)
2265 *   {
2266 *  
2267 *   if (time.current() < parameters.load_time+0.01*time.get_delta_t())
2268 *   {
2269 *   const double delta_length = parameters.length*(parameters.stretch - 1.0)*parameters.scale;
2270 *   const unsigned int n_stretch_steps = static_cast<unsigned int>(parameters.load_time/time.get_delta_t());
2271 *   const double delta_u_y = delta_length/2.0/n_stretch_steps;
2272 *   VectorTools::interpolate_boundary_values(dof_handler,
2273 *   boundary_id,
2274 *   Functions::ConstantFunction<dim>(delta_u_y,n_components),
2275 *   constraints,
2276 *   fe.component_mask(y_displacement));
2277 *   }
2278 *   else
2279 *   VectorTools::interpolate_boundary_values(dof_handler,
2280 *   boundary_id,
2281 *   Functions::ZeroFunction<dim>(n_components),
2282 *   constraints,
2283 *   fe.component_mask(y_displacement));
2284 *   }
2285 *   else
2286 *   VectorTools::interpolate_boundary_values(dof_handler,
2287 *   boundary_id,
2288 *   Functions::ZeroFunction<dim>(n_components),
2289 *   constraints,
2290 *   fe.component_mask(y_displacement));
2291 *   }
2292 *   constraints.close();
2293 *   }
2294 *   template <int dim>
2295 *   std::pair<unsigned int, double>
2296 *   Solid<dim>::solve_linear_system(LA::MPI::BlockVector &newton_update)
2297 *   {
2298 *   unsigned int lin_it = 0;
2299 *   double lin_res = 0.0;
2300 *  
2301 *   timer.enter_subsection("Linear solver");
2302 *   pcout << " SLV " << std::flush;
2303 *  
2304 *   const LA::MPI::Vector &f_u = system_rhs.block(u_block);
2305 *   const LA::MPI::Vector &f_p = system_rhs.block(p_block);
2306 *   const LA::MPI::Vector &f_J = system_rhs.block(J_block);
2307 *   LA::MPI::Vector &d_u = newton_update.block(u_block);
2308 *   LA::MPI::Vector &d_p = newton_update.block(p_block);
2309 *   LA::MPI::Vector &d_J = newton_update.block(J_block);
2310 *   const auto K_uu = linear_operator<LA::MPI::Vector>(tangent_matrix.block(u_block, u_block));
2311 *   const auto K_up = linear_operator<LA::MPI::Vector>(tangent_matrix.block(u_block, p_block));
2312 *   const auto K_pu = linear_operator<LA::MPI::Vector>(tangent_matrix.block(p_block, u_block));
2313 *   const auto K_Jp = linear_operator<LA::MPI::Vector>(tangent_matrix.block(J_block, p_block));
2314 *   const auto K_JJ = linear_operator<LA::MPI::Vector>(tangent_matrix.block(J_block, J_block));
2315 *  
2316 *   LA::PreconditionJacobi preconditioner_K_Jp_inv;
2317 *   preconditioner_K_Jp_inv.initialize(
2318 *   tangent_matrix.block(J_block, p_block),
2319 *   LA::PreconditionJacobi::AdditionalData());
2320 *   ReductionControl solver_control_K_Jp_inv (
2321 *   static_cast<unsigned int>(tangent_matrix.block(J_block, p_block).m()
2322 *   * parameters.max_iterations_lin),
2323 *   1.0e-30, 1e-6);
2324 *   ::SolverCG<LA::MPI::Vector> solver_K_Jp_inv (solver_control_K_Jp_inv);
2325 *  
2326 *   const auto K_Jp_inv = inverse_operator(K_Jp,
2327 *   solver_K_Jp_inv,
2328 *   preconditioner_K_Jp_inv);
2329 *   const auto K_pJ_inv = transpose_operator(K_Jp_inv);
2330 *   const auto K_pp_bar = K_Jp_inv * K_JJ * K_pJ_inv;
2331 *   const auto K_uu_bar_bar = K_up * K_pp_bar * K_pu;
2332 *   const auto K_uu_con = K_uu + K_uu_bar_bar;
2333 *  
2334 *   LA::PreconditionAMG preconditioner_K_con_inv;
2335 *   preconditioner_K_con_inv.initialize(
2336 *   tangent_matrix.block(u_block, u_block),
2337 *   LA::PreconditionAMG::AdditionalData(
2338 *   true /*elliptic*/,
2339 *   (parameters.poly_degree > 1 /*higher_order_elements*/)) );
2340 *   ReductionControl solver_control_K_con_inv (
2341 *   static_cast<unsigned int>(tangent_matrix.block(u_block, u_block).m()
2342 *   * parameters.max_iterations_lin),
2343 *   1.0e-30, parameters.tol_lin);
2344 *   ::SolverSelector<LA::MPI::Vector> solver_K_con_inv;
2345 *   solver_K_con_inv.select(parameters.type_lin);
2346 *   solver_K_con_inv.set_control(solver_control_K_con_inv);
2347 *   const auto K_uu_con_inv = inverse_operator(K_uu_con,
2348 *   solver_K_con_inv,
2349 *   preconditioner_K_con_inv);
2350 *  
2351 *   d_u = K_uu_con_inv*(f_u - K_up*(K_Jp_inv*f_J - K_pp_bar*f_p));
2352 *   lin_it = solver_control_K_con_inv.last_step();
2353 *   lin_res = solver_control_K_con_inv.last_value();
2354 *   timer.leave_subsection();
2355 *  
2356 *   timer.enter_subsection("Linear solver postprocessing");
2357 *   d_J = K_pJ_inv*(f_p - K_pu*d_u);
2358 *   d_p = K_Jp_inv*(f_J - K_JJ*d_J);
2359 *   timer.leave_subsection();
2360 *  
2361 *   constraints.distribute(newton_update);
2362 *   return std::make_pair(lin_it, lin_res);
2363 *   }
2364 *   template <int dim>
2365 *   void
2366 *   Solid<dim>::update_end_timestep ()
2367 *   {
2368 *   FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
2369 *   cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2370 *   dof_handler.begin_active()),
2371 *   endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2372 *   dof_handler.end());
2373 *   for (; cell != endc; ++cell)
2374 *   {
2375 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2376 *   const std::vector<std::shared_ptr<PointHistory<dim> > > lqph =
2377 *   quadrature_point_history.get_data(cell);
2378 *   Assert(lqph.size() == n_q_points, ExcInternalError());
2379 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2380 *   lqph[q_point]->update_end_timestep();
2381 *   }
2382 *   }
2383 *  
2384 *   template<int dim>
2385 *   class FilteredDataOut : public DataOut<dim>
2386 *   {
2387 *   public:
2388 *   FilteredDataOut (const unsigned int subdomain_id)
2389 *   :
2390 *   subdomain_id (subdomain_id)
2391 *   {}
2392 *  
2393 *   virtual ~FilteredDataOut() {}
2394 *  
2395 *   virtual typename DataOut<dim>::cell_iterator
2396 *   first_cell ()
2397 *   {
2398 *   auto cell = this->dofs->begin_active();
2399 *   while ((cell != this->dofs->end()) &&
2400 *   (cell->subdomain_id() != subdomain_id))
2401 *   ++cell;
2402 *   return cell;
2403 *   }
2404 *  
2405 *   virtual typename DataOut<dim>::cell_iterator
2406 *   next_cell (const typename DataOut<dim>::cell_iterator &old_cell)
2407 *   {
2408 *   if (old_cell != this->dofs->end())
2409 *   {
2410 *   const IteratorFilters::SubdomainEqualTo predicate(subdomain_id);
2411 *   return
2412 *   ++(FilteredIterator
2413 *   <typename DataOut<dim>::cell_iterator>
2414 *   (predicate,old_cell));
2415 *   }
2416 *   else
2417 *   return old_cell;
2418 *   }
2419 *  
2420 *   private:
2421 *   const unsigned int subdomain_id;
2422 *   };
2423 *  
2424 *   template <int dim>
2425 *   void Solid<dim>::output_results(const unsigned int timestep,
2426 *   const double current_time) const
2427 *   {
2428 * @endcode
2429 *
2430 * Output -> Ghosted vector
2431 *
2432 * @code
2433 *   LA::MPI::BlockVector solution_total (locally_owned_partitioning,
2434 *   locally_relevant_partitioning,
2435 *   mpi_communicator,
2436 *   /*vector_writable = */ false);
2437 *   LA::MPI::BlockVector residual (locally_owned_partitioning,
2438 *   locally_relevant_partitioning,
2439 *   mpi_communicator,
2440 *   /*vector_writable = */ false);
2441 *   solution_total = solution_n;
2442 *   residual = system_rhs;
2443 *   residual *= -1.0;
2444 *  
2445 * @endcode
2446 *
2447 * --- Additional data ---
2448 *
2449 * @code
2450 *   Vector<double> material_id;
2451 *   Vector<double> polynomial_order;
2452 *   material_id.reinit(triangulation.n_active_cells());
2453 *   polynomial_order.reinit(triangulation.n_active_cells());
2454 *   std::vector<types::subdomain_id> partition_int (triangulation.n_active_cells());
2455 *  
2456 *   FilteredDataOut<dim> data_out(this_mpi_process);
2457 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
2458 *   data_component_interpretation(dim,
2459 *   DataComponentInterpretation::component_is_part_of_vector);
2460 *   data_component_interpretation.push_back(DataComponentInterpretation::component_is_scalar);
2461 *   data_component_interpretation.push_back(DataComponentInterpretation::component_is_scalar);
2462 *  
2463 *   GridTools::get_subdomain_association (triangulation, partition_int);
2464 *  
2465 * @endcode
2466 *
2467 * Can't use filtered iterators here because the cell
2468 * count "c" is incorrect for the parallel case
2469 *
2470 * @code
2471 *   unsigned int c = 0;
2473 *   cell = dof_handler.begin_active(),
2474 *   endc = dof_handler.end();
2475 *   for (; cell!=endc; ++cell, ++c)
2476 *   {
2477 *   if (cell->subdomain_id() != this_mpi_process) continue;
2478 *  
2479 *   material_id(c) = static_cast<int>(cell->material_id());
2480 *   }
2481 *  
2482 *   std::vector<std::string> solution_name(n_components, "solution_");
2483 *   std::vector<std::string> residual_name(n_components, "residual_");
2484 *   for (unsigned int c=0; c<n_components; ++c)
2485 *   {
2486 *   if (block_component[c] == u_block)
2487 *   {
2488 *   solution_name[c] += "u";
2489 *   residual_name[c] += "u";
2490 *   }
2491 *   else if (block_component[c] == p_block)
2492 *   {
2493 *   solution_name[c] += "p";
2494 *   residual_name[c] += "p";
2495 *   }
2496 *   else if (block_component[c] == J_block)
2497 *   {
2498 *   solution_name[c] += "J";
2499 *   residual_name[c] += "J";
2500 *   }
2501 *   else
2502 *   {
2503 *   Assert(c <= J_block, ExcInternalError());
2504 *   }
2505 *   }
2506 *  
2507 *   data_out.attach_dof_handler(dof_handler);
2508 *   data_out.add_data_vector(solution_total,
2509 *   solution_name,
2511 *   data_component_interpretation);
2512 *   data_out.add_data_vector(residual,
2513 *   residual_name,
2515 *   data_component_interpretation);
2516 *   const Vector<double> partitioning(partition_int.begin(),
2517 *   partition_int.end());
2518 *   data_out.add_data_vector (material_id, "material_id");
2519 *   data_out.add_data_vector (partitioning, "partitioning");
2520 *   data_out.build_patches(degree);
2521 *  
2522 *   struct Filename
2523 *   {
2524 *   static std::string get_filename_vtu (unsigned int process,
2525 *   unsigned int timestep,
2526 *   const unsigned int n_digits = 4)
2527 *   {
2528 *   std::ostringstream filename_vtu;
2529 *   filename_vtu
2530 *   << "solution-"
2531 *   << (std::to_string(dim) + "d")
2532 *   << "."
2533 *   << Utilities::int_to_string (process, n_digits)
2534 *   << "."
2535 *   << Utilities::int_to_string(timestep, n_digits)
2536 *   << ".vtu";
2537 *   return filename_vtu.str();
2538 *   }
2539 *  
2540 *   static std::string get_filename_pvtu (unsigned int timestep,
2541 *   const unsigned int n_digits = 4)
2542 *   {
2543 *   std::ostringstream filename_vtu;
2544 *   filename_vtu
2545 *   << "solution-"
2546 *   << (std::to_string(dim) + "d")
2547 *   << "."
2548 *   << Utilities::int_to_string(timestep, n_digits)
2549 *   << ".pvtu";
2550 *   return filename_vtu.str();
2551 *   }
2552 *  
2553 *   static std::string get_filename_pvd (void)
2554 *   {
2555 *   std::ostringstream filename_vtu;
2556 *   filename_vtu
2557 *   << "solution-"
2558 *   << (std::to_string(dim) + "d")
2559 *   << ".pvd";
2560 *   return filename_vtu.str();
2561 *   }
2562 *   };
2563 *  
2564 * @endcode
2565 *
2566 * Write out main data file
2567 *
2568 * @code
2569 *   const std::string filename_vtu = Filename::get_filename_vtu(this_mpi_process, timestep);
2570 *   std::ofstream output(filename_vtu.c_str());
2571 *   data_out.write_vtu(output);
2572 *  
2573 * @endcode
2574 *
2575 * Collection of files written in parallel
2576 * This next set of steps should only be performed
2577 * by master process
2578 *
2579 * @code
2580 *   if (this_mpi_process == 0)
2581 *   {
2582 * @endcode
2583 *
2584 * List of all files written out at this timestep by all processors
2585 *
2586 * @code
2587 *   std::vector<std::string> parallel_filenames_vtu;
2588 *   for (unsigned int p=0; p < n_mpi_processes; ++p)
2589 *   {
2590 *   parallel_filenames_vtu.push_back(Filename::get_filename_vtu(p, timestep));
2591 *   }
2592 *  
2593 *   const std::string filename_pvtu (Filename::get_filename_pvtu(timestep));
2594 *   std::ofstream pvtu_master(filename_pvtu.c_str());
2595 *   data_out.write_pvtu_record(pvtu_master,
2596 *   parallel_filenames_vtu);
2597 *  
2598 * @endcode
2599 *
2600 * Time dependent data master file
2601 *
2602 * @code
2603 *   static std::vector<std::pair<double,std::string> > time_and_name_history;
2604 *   time_and_name_history.push_back (std::make_pair (current_time,
2605 *   filename_pvtu));
2606 *   const std::string filename_pvd (Filename::get_filename_pvd());
2607 *   std::ofstream pvd_output (filename_pvd.c_str());
2608 *   DataOutBase::write_pvd_record (pvd_output, time_and_name_history);
2609 *   }
2610 *   }
2611 *   template <int dim>
2612 *   void Solid<dim>::compute_vertex_positions(std::vector<double> &real_time,
2613 *   std::vector<std::vector<Point<dim> > > &tracked_vertices,
2614 *   const LA::MPI::BlockVector &solution_total) const
2615 *   {
2616 *   real_time.push_back(time.current());
2617 *  
2618 *   std::vector<bool> vertex_found (tracked_vertices.size(), false);
2619 *   std::vector<Tensor<1,dim> > vertex_update (tracked_vertices.size());
2620 *  
2622 *   cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2623 *   dof_handler.begin_active()),
2624 *   endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2625 *   dof_handler.end());
2626 *   for (; cell != endc; ++cell)
2627 *   {
2628 *   Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2629 *   for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2630 *   {
2631 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
2632 *   {
2633 *   if (vertex_found[p] == true) continue;
2634 *  
2635 *   const Point<dim> pt_ref = tracked_vertices[p][0];
2636 *   if (cell->vertex(v).distance(pt_ref) < 1e-6*parameters.scale)
2637 *   {
2638 *   for (unsigned int d=0; d<dim; ++d)
2639 *   vertex_update[p][d] = solution_total(cell->vertex_dof_index(v,u_block+d));
2640 *  
2641 *   vertex_found[p] = true;
2642 *   }
2643 *   }
2644 *   }
2645 *   }
2646 *  
2647 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
2648 *   {
2649 *   const int found_on_n_processes = Utilities::MPI::sum(int(vertex_found[p]), mpi_communicator);
2650 *   Assert(found_on_n_processes>0, ExcMessage("Vertex not found on any processor"));
2651 *   Tensor<1,dim> update;
2652 *   for (unsigned int d=0; d<dim; ++d)
2653 *   update[d] = Utilities::MPI::sum(vertex_update[p][d], mpi_communicator);
2654 *   update /= found_on_n_processes;
2655 *   tracked_vertices[p].push_back(tracked_vertices[p][0] + update);
2656 *   }
2657 *  
2658 *   }
2659 *   }
2660 *   int main (int argc, char *argv[])
2661 *   {
2662 *   using namespace dealii;
2663 *   using namespace ViscoElasStripHole;
2664 *  
2665 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
2666 *  
2667 *   try
2668 *   {
2669 *   const unsigned int dim = 2; // Works in both 2d and 3d
2670 *   Solid<dim> solid("parameters.prm");
2671 *   solid.run();
2672 *   }
2673 *   catch (std::exception &exc)
2674 *   {
2675 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
2676 *   {
2677 *   std::cerr << std::endl << std::endl
2678 *   << "----------------------------------------------------"
2679 *   << std::endl;
2680 *   std::cerr << "Exception on processing: " << std::endl << exc.what()
2681 *   << std::endl << "Aborting!" << std::endl
2682 *   << "----------------------------------------------------"
2683 *   << std::endl;
2684 *   return 1;
2685 *   }
2686 *   }
2687 *   catch (...)
2688 *   {
2689 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
2690 *   {
2691 *   std::cerr << std::endl << std::endl
2692 *   << "----------------------------------------------------"
2693 *   << std::endl;
2694 *   std::cerr << "Unknown exception!" << std::endl << "Aborting!"
2695 *   << std::endl
2696 *   << "----------------------------------------------------"
2697 *   << std::endl;
2698 *   return 1;
2699 *   }
2700 *   }
2701 *   return 0;
2702 *   }
2703 * @endcode
2704
2705
2706*/
Definition fe_q.h:554
Definition point.h:111
@ wall_times
Definition timer.h:651
void refine_global(const unsigned int times=1)
cell_iterator end() const
active_cell_iterator begin_active(const unsigned int level=0) const
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Point< 2 > first
Definition grid_out.cc:4623
__global__ void set(Number *val, const Number s, const size_type N)
#define Assert(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
UpdateFlags
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
void write_pvd_record(std::ostream &out, const std::vector< std::pair< double, std::string > > &times_and_names)
void hyper_cube_with_cylindrical_hole(Triangulation< dim, spacedim > &triangulation, const double inner_radius=.25, const double outer_radius=.5, const double L=.5, const unsigned int repetitions=1, const bool colorize=false)
void cylinder(Triangulation< dim > &tria, const double radius=1., const double half_length=1.)
void create_triangulation_with_removed_cells(const Triangulation< dim, spacedim > &input_triangulation, const std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > &cells_to_remove, Triangulation< dim, spacedim > &result)
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Definition operators.h:49
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * end(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:92
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
void project(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const bool enforce_zero_boundary=false, const Quadrature< dim - 1 > &q_boundary=(dim > 1 ? QGauss< dim - 1 >(2) :Quadrature< dim - 1 >()), const bool project_to_boundary_first=false)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int material_id
Definition types.h:167
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)