deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Nonlinear_PoroViscoelasticity.h
Go to the documentation of this file.
1
594 *  
595 * @endcode
596 *
597 * We start by including all the necessary deal.II header files and some C++
598 * related ones. They have been discussed in detail in previous tutorial
599 * programs, so you need only refer to past tutorials for details.
600 *
601
602 *
603 *
604 * @code
605 *   #include <deal.II/base/function.h>
606 *   #include <deal.II/base/parameter_handler.h>
607 *   #include <deal.II/base/point.h>
608 *   #include <deal.II/base/quadrature_lib.h>
609 *   #include <deal.II/base/symmetric_tensor.h>
610 *   #include <deal.II/base/tensor.h>
611 *   #include <deal.II/base/timer.h>
612 *   #include <deal.II/base/work_stream.h>
613 *   #include <deal.II/base/mpi.h>
614 *   #include <deal.II/base/quadrature_point_data.h>
615 *  
616 *   #include <deal.II/differentiation/ad.h>
617 *  
618 *   #include <deal.II/distributed/shared_tria.h>
619 *  
620 *   #include <deal.II/dofs/dof_renumbering.h>
621 *   #include <deal.II/dofs/dof_tools.h>
622 *   #include <deal.II/dofs/dof_accessor.h>
623 *  
624 *   #include <deal.II/grid/filtered_iterator.h>
625 *   #include <deal.II/grid/grid_generator.h>
626 *   #include <deal.II/grid/grid_tools.h>
627 *   #include <deal.II/grid/grid_in.h>
628 *   #include <deal.II/grid/grid_out.h>
629 *   #include <deal.II/grid/manifold_lib.h>
630 *   #include <deal.II/grid/tria_accessor.h>
631 *   #include <deal.II/grid/tria_iterator.h>
632 *  
633 *   #include <deal.II/fe/fe_dgp_monomial.h>
634 *   #include <deal.II/fe/fe_q.h>
635 *   #include <deal.II/fe/fe_system.h>
636 *   #include <deal.II/fe/fe_tools.h>
637 *   #include <deal.II/fe/fe_values.h>
638 *  
639 *   #include <deal.II/lac/block_sparsity_pattern.h>
640 *   #include <deal.II/lac/affine_constraints.h>
641 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
642 *   #include <deal.II/lac/full_matrix.h>
643 *   #include <deal.II/lac/linear_operator.h>
644 *   #include <deal.II/lac/packaged_operation.h>
645 *  
646 *   #include <deal.II/lac/trilinos_block_sparse_matrix.h>
647 *   #include <deal.II/lac/trilinos_linear_operator.h>
648 *   #include <deal.II/lac/trilinos_parallel_block_vector.h>
649 *   #include <deal.II/lac/trilinos_precondition.h>
650 *   #include <deal.II/lac/trilinos_sparse_matrix.h>
651 *   #include <deal.II/lac/trilinos_sparsity_pattern.h>
652 *   #include <deal.II/lac/trilinos_solver.h>
653 *   #include <deal.II/lac/trilinos_vector.h>
654 *  
655 *   #include <deal.II/lac/block_vector.h>
656 *   #include <deal.II/lac/vector.h>
657 *  
658 *   #include <deal.II/numerics/data_postprocessor.h>
659 *   #include <deal.II/numerics/data_out.h>
660 *   #include <deal.II/numerics/data_out_faces.h>
661 *   #include <deal.II/numerics/fe_field_function.h>
662 *   #include <deal.II/numerics/vector_tools.h>
663 *  
664 *   #include <deal.II/physics/transformations.h>
665 *   #include <deal.II/physics/elasticity/kinematics.h>
666 *   #include <deal.II/physics/elasticity/standard_tensors.h>
667 *  
668 *   #include <iostream>
669 *   #include <fstream>
670 *   #include <numeric>
671 *   #include <iomanip>
672 *  
673 *  
674 * @endcode
675 *
676 * We create a namespace for everything that relates to
677 * the nonlinear poro-viscoelastic formulation,
678 * and import all the deal.II function and class names into it:
679 *
680 * @code
681 *   namespace NonLinearPoroViscoElasticity
682 *   {
683 *   using namespace dealii;
684 *  
685 * @endcode
686 *
687 *
688 * <a name="nonlinear-poro-viscoelasticity.cc-Runtimeparameters"></a>
689 * <h3>Run-time parameters</h3>
690 *
691
692 *
693 * Set up a ParameterHandler object to read in the parameter choices at run-time
694 * introduced by the user through the file "parameters.prm"
695 *
696 * @code
697 *   namespace Parameters
698 *   {
699 * @endcode
700 *
701 *
702 * <a name="nonlinear-poro-viscoelasticity.cc-FiniteElementsystem"></a>
703 * <h4>Finite Element system</h4>
704 * Here we specify the polynomial order used to approximate the solution,
705 * both for the displacements and pressure unknowns.
706 * The quadrature order should be adjusted accordingly.
707 *
708 * @code
709 *   struct FESystem
710 *   {
711 *   unsigned int poly_degree_displ;
712 *   unsigned int poly_degree_pore;
713 *   unsigned int quad_order;
714 *  
715 *   static void
716 *   declare_parameters(ParameterHandler &prm);
717 *  
718 *   void
719 *   parse_parameters(ParameterHandler &prm);
720 *   };
721 *  
722 *   void FESystem::declare_parameters(ParameterHandler &prm)
723 *   {
724 *   prm.enter_subsection("Finite element system");
725 *   {
726 *   prm.declare_entry("Polynomial degree displ", "2",
727 *   Patterns::Integer(0),
728 *   "Displacement system polynomial order");
729 *  
730 *   prm.declare_entry("Polynomial degree pore", "1",
731 *   Patterns::Integer(0),
732 *   "Pore pressure system polynomial order");
733 *  
734 *   prm.declare_entry("Quadrature order", "3",
735 *   Patterns::Integer(0),
736 *   "Gauss quadrature order");
737 *   }
738 *   prm.leave_subsection();
739 *   }
740 *  
741 *   void FESystem::parse_parameters(ParameterHandler &prm)
742 *   {
743 *   prm.enter_subsection("Finite element system");
744 *   {
745 *   poly_degree_displ = prm.get_integer("Polynomial degree displ");
746 *   poly_degree_pore = prm.get_integer("Polynomial degree pore");
747 *   quad_order = prm.get_integer("Quadrature order");
748 *   }
749 *   prm.leave_subsection();
750 *   }
751 *  
752 * @endcode
753 *
754 *
755 * <a name="nonlinear-poro-viscoelasticity.cc-Geometry"></a>
756 * <h4>Geometry</h4>
757 * These parameters are related to the geometry definition and mesh generation.
758 * We select the type of problem to solve and introduce the desired load values.
759 *
760 * @code
761 *   struct Geometry
762 *   {
763 *   std::string geom_type;
764 *   unsigned int global_refinement;
765 *   double scale;
766 *   std::string load_type;
767 *   double load;
768 *   unsigned int num_cycle_sets;
769 *   double fluid_flow;
770 *   double drained_pressure;
771 *  
772 *   static void
773 *   declare_parameters(ParameterHandler &prm);
774 *  
775 *   void
776 *   parse_parameters(ParameterHandler &prm);
777 *   };
778 *  
779 *   void Geometry::declare_parameters(ParameterHandler &prm)
780 *   {
781 *   prm.enter_subsection("Geometry");
782 *   {
783 *   prm.declare_entry("Geometry type", "Ehlers_tube_step_load",
784 *   Patterns::Selection("Ehlers_tube_step_load"
785 *   "|Ehlers_tube_increase_load"
786 *   "|Ehlers_cube_consolidation"
787 *   "|Franceschini_consolidation"
788 *   "|Budday_cube_tension_compression"
789 *   "|Budday_cube_tension_compression_fully_fixed"
790 *   "|Budday_cube_shear_fully_fixed"),
791 *   "Type of geometry used. "
792 *   "For Ehlers verification examples see Ehlers and Eipper (1999). "
793 *   "For Franceschini brain consolidation see Franceschini et al. (2006)"
794 *   "For Budday brain examples see Budday et al. (2017)");
795 *  
796 *   prm.declare_entry("Global refinement", "1",
797 *   Patterns::Integer(0),
798 *   "Global refinement level");
799 *  
800 *   prm.declare_entry("Grid scale", "1.0",
801 *   Patterns::Double(0.0),
802 *   "Global grid scaling factor");
803 *  
804 *   prm.declare_entry("Load type", "pressure",
805 *   Patterns::Selection("pressure|displacement|none"),
806 *   "Type of loading");
807 *  
808 *   prm.declare_entry("Load value", "-7.5e+6",
809 *   Patterns::Double(),
810 *   "Loading value");
811 *  
812 *   prm.declare_entry("Number of cycle sets", "1",
813 *   Patterns::Integer(1,2),
814 *   "Number of times each set of 3 cycles is repeated, only for "
815 *   "Budday_cube_tension_compression and Budday_cube_tension_compression_fully_fixed. "
816 *   "Load value is doubled in second set, load rate is kept constant."
817 *   "Final time indicates end of second cycle set.");
818 *  
819 *   prm.declare_entry("Fluid flow value", "0.0",
820 *   Patterns::Double(),
821 *   "Prescribed fluid flow. Not implemented in any example yet.");
822 *  
823 *   prm.declare_entry("Drained pressure", "0.0",
824 *   Patterns::Double(),
825 *   "Increase of pressure value at drained boundary w.r.t the atmospheric pressure.");
826 *   }
827 *   prm.leave_subsection();
828 *   }
829 *  
830 *   void Geometry::parse_parameters(ParameterHandler &prm)
831 *   {
832 *   prm.enter_subsection("Geometry");
833 *   {
834 *   geom_type = prm.get("Geometry type");
835 *   global_refinement = prm.get_integer("Global refinement");
836 *   scale = prm.get_double("Grid scale");
837 *   load_type = prm.get("Load type");
838 *   load = prm.get_double("Load value");
839 *   num_cycle_sets = prm.get_integer("Number of cycle sets");
840 *   fluid_flow = prm.get_double("Fluid flow value");
841 *   drained_pressure = prm.get_double("Drained pressure");
842 *   }
843 *   prm.leave_subsection();
844 *   }
845 *  
846 * @endcode
847 *
848 *
849 * <a name="nonlinear-poro-viscoelasticity.cc-Materials"></a>
850 * <h4>Materials</h4>
851 *
852
853 *
854 * Here we select the type of material for the solid component
855 * and define the corresponding material parameters.
856 * Then we define he fluid data, including the type of
857 * seepage velocity definition to use.
858 *
859 * @code
860 *   struct Materials
861 *   {
862 *   std::string mat_type;
863 *   double lambda;
864 *   double mu;
865 *   double mu1_infty;
866 *   double mu2_infty;
867 *   double mu3_infty;
868 *   double alpha1_infty;
869 *   double alpha2_infty;
870 *   double alpha3_infty;
871 *   double mu1_mode_1;
872 *   double mu2_mode_1;
873 *   double mu3_mode_1;
874 *   double alpha1_mode_1;
875 *   double alpha2_mode_1;
876 *   double alpha3_mode_1;
877 *   double viscosity_mode_1;
878 *   std::string fluid_type;
879 *   double solid_vol_frac;
880 *   double kappa_darcy;
881 *   double init_intrinsic_perm;
882 *   double viscosity_FR;
883 *   double init_darcy_coef;
884 *   double weight_FR;
885 *   bool gravity_term;
886 *   int gravity_direction;
887 *   double gravity_value;
888 *   double density_FR;
889 *   double density_SR;
890 *   enum SymmetricTensorEigenvectorMethod eigen_solver;
891 *  
892 *   static void
893 *   declare_parameters(ParameterHandler &prm);
894 *  
895 *   void
896 *   parse_parameters(ParameterHandler &prm);
897 *   };
898 *  
899 *   void Materials::declare_parameters(ParameterHandler &prm)
900 *   {
901 *   prm.enter_subsection("Material properties");
902 *   {
903 *   prm.declare_entry("material", "Neo-Hooke",
904 *   Patterns::Selection("Neo-Hooke|Ogden|visco-Ogden"),
905 *   "Type of material used in the problem");
906 *  
907 *   prm.declare_entry("lambda", "8.375e6",
908 *   Patterns::Double(0,1e100),
909 *   "First Lamé parameter for extension function related to compactation point in solid material [Pa].");
910 *  
911 *   prm.declare_entry("shear modulus", "5.583e6",
912 *   Patterns::Double(0,1e100),
913 *   "shear modulus for Neo-Hooke materials [Pa].");
914 *  
915 *   prm.declare_entry("eigen solver", "QL Implicit Shifts",
916 *   Patterns::Selection("QL Implicit Shifts|Jacobi"),
917 *   "The type of eigen solver to be used for Ogden and visco-Ogden models.");
918 *  
919 *   prm.declare_entry("mu1", "0.0",
920 *   Patterns::Double(),
921 *   "Shear material parameter 'mu1' for Ogden material [Pa].");
922 *  
923 *   prm.declare_entry("mu2", "0.0",
924 *   Patterns::Double(),
925 *   "Shear material parameter 'mu2' for Ogden material [Pa].");
926 *  
927 *   prm.declare_entry("mu3", "0.0",
928 *   Patterns::Double(),
929 *   "Shear material parameter 'mu1' for Ogden material [Pa].");
930 *  
931 *   prm.declare_entry("alpha1", "1.0",
932 *   Patterns::Double(),
933 *   "Stiffness material parameter 'alpha1' for Ogden material [-].");
934 *  
935 *   prm.declare_entry("alpha2", "1.0",
936 *   Patterns::Double(),
937 *   "Stiffness material parameter 'alpha2' for Ogden material [-].");
938 *  
939 *   prm.declare_entry("alpha3", "1.0",
940 *   Patterns::Double(),
941 *   "Stiffness material parameter 'alpha3' for Ogden material [-].");
942 *  
943 *   prm.declare_entry("mu1_1", "0.0",
944 *   Patterns::Double(),
945 *   "Shear material parameter 'mu1' for first viscous mode in Ogden material [Pa].");
946 *  
947 *   prm.declare_entry("mu2_1", "0.0",
948 *   Patterns::Double(),
949 *   "Shear material parameter 'mu2' for first viscous mode in Ogden material [Pa].");
950 *  
951 *   prm.declare_entry("mu3_1", "0.0",
952 *   Patterns::Double(),
953 *   "Shear material parameter 'mu1' for first viscous mode in Ogden material [Pa].");
954 *  
955 *   prm.declare_entry("alpha1_1", "1.0",
956 *   Patterns::Double(),
957 *   "Stiffness material parameter 'alpha1' for first viscous mode in Ogden material [-].");
958 *  
959 *   prm.declare_entry("alpha2_1", "1.0",
960 *   Patterns::Double(),
961 *   "Stiffness material parameter 'alpha2' for first viscous mode in Ogden material [-].");
962 *  
963 *   prm.declare_entry("alpha3_1", "1.0",
964 *   Patterns::Double(),
965 *   "Stiffness material parameter 'alpha3' for first viscous mode in Ogden material [-].");
966 *  
967 *   prm.declare_entry("viscosity_1", "1e-10",
968 *   Patterns::Double(1e-10,1e100),
969 *   "Deformation-independent viscosity parameter 'eta_1' for first viscous mode in Ogden material [-].");
970 *  
971 *   prm.declare_entry("seepage definition", "Ehlers",
972 *   Patterns::Selection("Markert|Ehlers"),
973 *   "Type of formulation used to define the seepage velocity in the problem. "
974 *   "Choose between Markert formulation of deformation-dependent intrinsic permeability "
975 *   "and Ehlers formulation of deformation-dependent Darcy flow coefficient.");
976 *  
977 *   prm.declare_entry("initial solid volume fraction", "0.67",
978 *   Patterns::Double(0.001,0.999),
979 *   "Initial porosity (solid volume fraction, 0 < n_0s < 1)");
980 *  
981 *   prm.declare_entry("kappa", "0.0",
982 *   Patterns::Double(0,100),
983 *   "Deformation-dependency control parameter for specific permeability (kappa >= 0)");
984 *  
985 *   prm.declare_entry("initial intrinsic permeability", "0.0",
986 *   Patterns::Double(0,1e100),
987 *   "Initial intrinsic permeability parameter [m^2] (isotropic permeability). To be used with Markert formulation.");
988 *  
989 *   prm.declare_entry("fluid viscosity", "0.0",
990 *   Patterns::Double(0, 1e100),
991 *   "Effective shear viscosity parameter of the fluid [Pa·s, (N·s)/m^2]. To be used with Markert formulation.");
992 *  
993 *   prm.declare_entry("initial Darcy coefficient", "1.0e-4",
994 *   Patterns::Double(0,1e100),
995 *   "Initial Darcy flow coefficient [m/s] (isotropic permeability). To be used with Ehlers formulation.");
996 *  
997 *   prm.declare_entry("fluid weight", "1.0e4",
998 *   Patterns::Double(0, 1e100),
999 *   "Effective weight of the fluid [N/m^3]. To be used with Ehlers formulation.");
1000 *  
1001 *   prm.declare_entry("gravity term", "false",
1002 *   Patterns::Bool(),
1003 *   "Gravity term considered (true) or neglected (false)");
1004 *  
1005 *   prm.declare_entry("fluid density", "1.0",
1006 *   Patterns::Double(0,1e100),
1007 *   "Real (or effective) density of the fluid");
1008 *  
1009 *   prm.declare_entry("solid density", "1.0",
1010 *   Patterns::Double(0,1e100),
1011 *   "Real (or effective) density of the solid");
1012 *  
1013 *   prm.declare_entry("gravity direction", "2",
1014 *   Patterns::Integer(0,2),
1015 *   "Direction of gravity (unit vector 0 for x, 1 for y, 2 for z)");
1016 *  
1017 *   prm.declare_entry("gravity value", "-9.81",
1018 *   Patterns::Double(),
1019 *   "Value of gravity (be careful to have consistent units!)");
1020 *   }
1021 *   prm.leave_subsection();
1022 *   }
1023 *  
1024 *   void Materials::parse_parameters(ParameterHandler &prm)
1025 *   {
1026 *   prm.enter_subsection("Material properties");
1027 *   {
1028 * @endcode
1029 *
1030 * Solid
1031 *
1032 * @code
1033 *   mat_type = prm.get("material");
1034 *   lambda = prm.get_double("lambda");
1035 *   mu = prm.get_double("shear modulus");
1036 *   mu1_infty = prm.get_double("mu1");
1037 *   mu2_infty = prm.get_double("mu2");
1038 *   mu3_infty = prm.get_double("mu3");
1039 *   alpha1_infty = prm.get_double("alpha1");
1040 *   alpha2_infty = prm.get_double("alpha2");
1041 *   alpha3_infty = prm.get_double("alpha3");
1042 *   mu1_mode_1 = prm.get_double("mu1_1");
1043 *   mu2_mode_1 = prm.get_double("mu2_1");
1044 *   mu3_mode_1 = prm.get_double("mu3_1");
1045 *   alpha1_mode_1 = prm.get_double("alpha1_1");
1046 *   alpha2_mode_1 = prm.get_double("alpha2_1");
1047 *   alpha3_mode_1 = prm.get_double("alpha3_1");
1048 *   viscosity_mode_1 = prm.get_double("viscosity_1");
1049 * @endcode
1050 *
1051 * Fluid
1052 *
1053 * @code
1054 *   fluid_type = prm.get("seepage definition");
1055 *   solid_vol_frac = prm.get_double("initial solid volume fraction");
1056 *   kappa_darcy = prm.get_double("kappa");
1057 *   init_intrinsic_perm = prm.get_double("initial intrinsic permeability");
1058 *   viscosity_FR = prm.get_double("fluid viscosity");
1059 *   init_darcy_coef = prm.get_double("initial Darcy coefficient");
1060 *   weight_FR = prm.get_double("fluid weight");
1061 * @endcode
1062 *
1063 * Gravity effects
1064 *
1065 * @code
1066 *   gravity_term = prm.get_bool("gravity term");
1067 *   density_FR = prm.get_double("fluid density");
1068 *   density_SR = prm.get_double("solid density");
1069 *   gravity_direction = prm.get_integer("gravity direction");
1070 *   gravity_value = prm.get_double("gravity value");
1071 *  
1072 *   if ( (fluid_type == "Markert") && ((init_intrinsic_perm == 0.0) || (viscosity_FR == 0.0)) )
1073 *   AssertThrow(false, ExcMessage("Markert seepage velocity formulation requires the definition of "
1074 *   "'initial intrinsic permeability' and 'fluid viscosity' greater than 0.0."));
1075 *  
1076 *   if ( (fluid_type == "Ehlers") && ((init_darcy_coef == 0.0) || (weight_FR == 0.0)) )
1077 *   AssertThrow(false, ExcMessage("Ehler seepage velocity formulation requires the definition of "
1078 *   "'initial Darcy coefficient' and 'fluid weight' greater than 0.0."));
1079 *  
1080 *   const std::string eigen_solver_type = prm.get("eigen solver");
1081 *   if (eigen_solver_type == "QL Implicit Shifts")
1083 *   else if (eigen_solver_type == "Jacobi")
1084 *   eigen_solver = SymmetricTensorEigenvectorMethod::jacobi;
1085 *   else
1086 *   {
1087 *   AssertThrow(false, ExcMessage("Unknown eigen solver selected."));
1088 *   }
1089 *   }
1090 *   prm.leave_subsection();
1091 *   }
1092 *  
1093 * @endcode
1094 *
1095 *
1096 * <a name="nonlinear-poro-viscoelasticity.cc-Nonlinearsolver"></a>
1097 * <h4>Nonlinear solver</h4>
1098 *
1099
1100 *
1101 * We now define the tolerances and the maximum number of iterations for the
1102 * Newton-Raphson scheme used to solve the nonlinear system of governing equations.
1103 *
1104 * @code
1105 *   struct NonlinearSolver
1106 *   {
1107 *   unsigned int max_iterations_NR;
1108 *   double tol_f;
1109 *   double tol_u;
1110 *   double tol_p_fluid;
1111 *  
1112 *   static void
1113 *   declare_parameters(ParameterHandler &prm);
1114 *  
1115 *   void
1116 *   parse_parameters(ParameterHandler &prm);
1117 *   };
1118 *  
1119 *   void NonlinearSolver::declare_parameters(ParameterHandler &prm)
1120 *   {
1121 *   prm.enter_subsection("Nonlinear solver");
1122 *   {
1123 *   prm.declare_entry("Max iterations Newton-Raphson", "15",
1124 *   Patterns::Integer(0),
1125 *   "Number of Newton-Raphson iterations allowed");
1126 *  
1127 *   prm.declare_entry("Tolerance force", "1.0e-8",
1128 *   Patterns::Double(0.0),
1129 *   "Force residual tolerance");
1130 *  
1131 *   prm.declare_entry("Tolerance displacement", "1.0e-6",
1132 *   Patterns::Double(0.0),
1133 *   "Displacement error tolerance");
1134 *  
1135 *   prm.declare_entry("Tolerance pore pressure", "1.0e-6",
1136 *   Patterns::Double(0.0),
1137 *   "Pore pressure error tolerance");
1138 *   }
1139 *   prm.leave_subsection();
1140 *   }
1141 *  
1142 *   void NonlinearSolver::parse_parameters(ParameterHandler &prm)
1143 *   {
1144 *   prm.enter_subsection("Nonlinear solver");
1145 *   {
1146 *   max_iterations_NR = prm.get_integer("Max iterations Newton-Raphson");
1147 *   tol_f = prm.get_double("Tolerance force");
1148 *   tol_u = prm.get_double("Tolerance displacement");
1149 *   tol_p_fluid = prm.get_double("Tolerance pore pressure");
1150 *   }
1151 *   prm.leave_subsection();
1152 *   }
1153 *  
1154 * @endcode
1155 *
1156 *
1157 * <a name="nonlinear-poro-viscoelasticity.cc-Time"></a>
1158 * <h4>Time</h4>
1159 * Here we set the timestep size @f$ \varDelta t @f$ and the simulation end-time.
1160 *
1161 * @code
1162 *   struct Time
1163 *   {
1164 *   double end_time;
1165 *   double delta_t;
1166 *   static void
1167 *   declare_parameters(ParameterHandler &prm);
1168 *  
1169 *   void
1170 *   parse_parameters(ParameterHandler &prm);
1171 *   };
1172 *  
1173 *   void Time::declare_parameters(ParameterHandler &prm)
1174 *   {
1175 *   prm.enter_subsection("Time");
1176 *   {
1177 *   prm.declare_entry("End time", "10.0",
1178 *   Patterns::Double(),
1179 *   "End time");
1180 *  
1181 *   prm.declare_entry("Time step size", "0.002",
1182 *   Patterns::Double(1.0e-6),
1183 *   "Time step size. The value must be larger than the displacement error tolerance defined.");
1184 *   }
1185 *   prm.leave_subsection();
1186 *   }
1187 *  
1188 *   void Time::parse_parameters(ParameterHandler &prm)
1189 *   {
1190 *   prm.enter_subsection("Time");
1191 *   {
1192 *   end_time = prm.get_double("End time");
1193 *   delta_t = prm.get_double("Time step size");
1194 *   }
1195 *   prm.leave_subsection();
1196 *   }
1197 *  
1198 *  
1199 * @endcode
1200 *
1201 *
1202 * <a name="nonlinear-poro-viscoelasticity.cc-Output"></a>
1203 * <h4>Output</h4>
1204 * We can choose the frequency of the data for the output files.
1205 *
1206 * @code
1207 *   struct OutputParam
1208 *   {
1209 *  
1210 *   std::string outfiles_requested;
1211 *   unsigned int timestep_output;
1212 *   std::string outtype;
1213 *  
1214 *   static void
1215 *   declare_parameters(ParameterHandler &prm);
1216 *  
1217 *   void
1218 *   parse_parameters(ParameterHandler &prm);
1219 *   };
1220 *  
1221 *   void OutputParam::declare_parameters(ParameterHandler &prm)
1222 *   {
1223 *   prm.enter_subsection("Output parameters");
1224 *   {
1225 *   prm.declare_entry("Output files", "true",
1226 *   Patterns::Selection("true|false"),
1227 *   "Paraview output files to generate.");
1228 *   prm.declare_entry("Time step number output", "1",
1229 *   Patterns::Integer(0),
1230 *   "Output data for time steps multiple of the given "
1231 *   "integer value.");
1232 *   prm.declare_entry("Averaged results", "nodes",
1233 *   Patterns::Selection("elements|nodes"),
1234 *   "Output data associated with integration point values"
1235 *   " averaged on elements or on nodes.");
1236 *   }
1237 *   prm.leave_subsection();
1238 *   }
1239 *  
1240 *   void OutputParam::parse_parameters(ParameterHandler &prm)
1241 *   {
1242 *   prm.enter_subsection("Output parameters");
1243 *   {
1244 *   outfiles_requested = prm.get("Output files");
1245 *   timestep_output = prm.get_integer("Time step number output");
1246 *   outtype = prm.get("Averaged results");
1247 *   }
1248 *   prm.leave_subsection();
1249 *   }
1250 *  
1251 * @endcode
1252 *
1253 *
1254 * <a name="nonlinear-poro-viscoelasticity.cc-Allparameters"></a>
1255 * <h4>All parameters</h4>
1256 * We finally consolidate all of the above structures into a single container that holds all the run-time selections.
1257 *
1258 * @code
1259 *   struct AllParameters : public FESystem,
1260 *   public Geometry,
1261 *   public Materials,
1262 *   public NonlinearSolver,
1263 *   public Time,
1264 *   public OutputParam
1265 *   {
1266 *   AllParameters(const std::string &input_file);
1267 *  
1268 *   static void
1269 *   declare_parameters(ParameterHandler &prm);
1270 *  
1271 *   void
1272 *   parse_parameters(ParameterHandler &prm);
1273 *   };
1274 *  
1275 *   AllParameters::AllParameters(const std::string &input_file)
1276 *   {
1277 *   ParameterHandler prm;
1278 *   declare_parameters(prm);
1279 *   prm.parse_input(input_file);
1280 *   parse_parameters(prm);
1281 *   }
1282 *  
1283 *   void AllParameters::declare_parameters(ParameterHandler &prm)
1284 *   {
1285 *   FESystem::declare_parameters(prm);
1286 *   Geometry::declare_parameters(prm);
1287 *   Materials::declare_parameters(prm);
1288 *   NonlinearSolver::declare_parameters(prm);
1289 *   Time::declare_parameters(prm);
1290 *   OutputParam::declare_parameters(prm);
1291 *   }
1292 *  
1293 *   void AllParameters::parse_parameters(ParameterHandler &prm)
1294 *   {
1295 *   FESystem::parse_parameters(prm);
1296 *   Geometry::parse_parameters(prm);
1297 *   Materials::parse_parameters(prm);
1298 *   NonlinearSolver::parse_parameters(prm);
1299 *   Time::parse_parameters(prm);
1300 *   OutputParam::parse_parameters(prm);
1301 *   }
1302 *   }
1303 *  
1304 * @endcode
1305 *
1306 *
1307 * <a name="nonlinear-poro-viscoelasticity.cc-Timeclass"></a>
1308 * <h3>Time class</h3>
1309 * A simple class to store time data.
1310 * For simplicity we assume a constant time step size.
1311 *
1312 * @code
1313 *   class Time
1314 *   {
1315 *   public:
1316 *   Time (const double time_end,
1317 *   const double delta_t)
1318 *   :
1319 *   timestep(0),
1320 *   time_current(0.0),
1321 *   time_end(time_end),
1322 *   delta_t(delta_t)
1323 *   {}
1324 *  
1325 *   virtual ~Time()
1326 *   {}
1327 *  
1328 *   double get_current() const
1329 *   {
1330 *   return time_current;
1331 *   }
1332 *   double get_end() const
1333 *   {
1334 *   return time_end;
1335 *   }
1336 *   double get_delta_t() const
1337 *   {
1338 *   return delta_t;
1339 *   }
1340 *   unsigned int get_timestep() const
1341 *   {
1342 *   return timestep;
1343 *   }
1344 *   void increment_time ()
1345 *   {
1346 *   time_current += delta_t;
1347 *   ++timestep;
1348 *   }
1349 *  
1350 *   private:
1351 *   unsigned int timestep;
1352 *   double time_current;
1353 *   double time_end;
1354 *   const double delta_t;
1355 *   };
1356 *  
1357 * @endcode
1358 *
1359 *
1360 * <a name="nonlinear-poro-viscoelasticity.cc-Constitutiveequationforthesolidcomponentofthebiphasicmaterial"></a>
1361 * <h3>Constitutive equation for the solid component of the biphasic material</h3>
1362 *
1363
1364 *
1365 *
1366 * <a name="nonlinear-poro-viscoelasticity.cc-Baseclassgenerichyperelasticmaterial"></a>
1367 * <h4>Base class: generic hyperelastic material</h4>
1368 * The "extra" Kirchhoff stress in the solid component is the sum of isochoric
1369 * and a volumetric part:
1370 * @f$\mathbf{\tau} = \mathbf{\tau}_E^{(\bullet)} + \mathbf{\tau}^{\textrm{vol}}@f$.
1371 * The deviatoric part changes depending on the type of material model selected:
1372 * Neo-Hooken hyperelasticity, Ogden hyperelasticiy,
1373 * or a single-mode finite viscoelasticity based on the Ogden hyperelastic model.
1374 * In this base class we declare it as a virtual function,
1375 * and it will be defined for each model type in the corresponding derived class.
1376 * We define here the volumetric component, which depends on the
1377 * extension function @f$U(J_S)@f$ selected, and in this case is the same for all models.
1378 * We use the function proposed by
1379 * Ehlers & Eipper 1999 doi:10.1023/A:1006565509095.
1380 * We also define some public functions to access and update the internal variables.
1381 *
1382 * @code
1383 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> >
1384 *   class Material_Hyperelastic
1385 *   {
1386 *   public:
1387 *   Material_Hyperelastic(const Parameters::AllParameters &parameters,
1388 *   const Time &time)
1389 *   :
1390 *   n_OS (parameters.solid_vol_frac),
1391 *   lambda (parameters.lambda),
1392 *   time(time),
1393 *   det_F (1.0),
1394 *   det_F_converged (1.0),
1395 *   eigen_solver (parameters.eigen_solver)
1396 *   {}
1397 *   ~Material_Hyperelastic()
1398 *   {}
1399 *  
1401 *   get_tau_E(const Tensor<2,dim, NumberType> &F) const
1402 *   {
1403 *   return ( get_tau_E_base(F) + get_tau_E_ext_func(F) );
1404 *   }
1405 *  
1407 *   get_Cauchy_E(const Tensor<2, dim, NumberType> &F) const
1408 *   {
1409 *   const NumberType det_F = determinant(F);
1410 *   Assert(det_F > 0, ExcInternalError());
1411 *   return get_tau_E(F)*NumberType(1/det_F);
1412 *   }
1413 *  
1414 *   double
1415 *   get_converged_det_F() const
1416 *   {
1417 *   return det_F_converged;
1418 *   }
1419 *  
1420 *   virtual void
1421 *   update_end_timestep()
1422 *   {
1423 *   det_F_converged = det_F;
1424 *   }
1425 *  
1426 *   virtual void
1427 *   update_internal_equilibrium( const Tensor<2, dim, NumberType> &F )
1428 *   {
1429 *   det_F = Tensor<0,dim,double>(determinant(F));
1430 *   }
1431 *  
1432 *   virtual double
1433 *   get_viscous_dissipation( ) const = 0;
1434 *  
1435 *   const double n_OS;
1436 *   const double lambda;
1437 *   const Time &time;
1438 *   double det_F;
1439 *   double det_F_converged;
1440 *   const enum SymmetricTensorEigenvectorMethod eigen_solver;
1441 *  
1442 *   protected:
1444 *   get_tau_E_ext_func(const Tensor<2,dim, NumberType> &F) const
1445 *   {
1446 *   const NumberType det_F = determinant(F);
1447 *   Assert(det_F > 0, ExcInternalError());
1448 *  
1449 *   static const SymmetricTensor< 2, dim, double>
1451 *   return ( NumberType(lambda * (1.0-n_OS)*(1.0-n_OS)
1452 *   * (det_F/(1.0-n_OS) - det_F/(det_F-n_OS))) * I );
1453 *   }
1454 *  
1456 *   get_tau_E_base(const Tensor<2,dim, NumberType> &F) const = 0;
1457 *   };
1458 *  
1459 * @endcode
1460 *
1461 *
1462 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassNeoHookeanhyperelasticmaterial"></a>
1463 * <h4>Derived class: Neo-Hookean hyperelastic material</h4>
1464 *
1465 * @code
1466 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> >
1467 *   class NeoHooke : public Material_Hyperelastic < dim, NumberType >
1468 *   {
1469 *   public:
1470 *   NeoHooke(const Parameters::AllParameters &parameters,
1471 *   const Time &time)
1472 *   :
1473 *   Material_Hyperelastic< dim, NumberType > (parameters,time),
1474 *   mu(parameters.mu)
1475 *   {}
1476 *   virtual ~NeoHooke()
1477 *   {}
1478 *  
1479 *   double
1480 *   get_viscous_dissipation() const override
1481 *   {
1482 *   return 0.0;
1483 *   }
1484 *  
1485 *   protected:
1486 *   const double mu;
1487 *  
1488 *   SymmetricTensor<2, dim, NumberType>
1489 *   get_tau_E_base(const Tensor<2,dim, NumberType> &F) const override
1490 *   {
1491 *   static const SymmetricTensor< 2, dim, double>
1493 *  
1494 *   const bool use_standard_model = true;
1495 *  
1496 *   if (use_standard_model)
1497 *   {
1498 * @endcode
1499 *
1500 * Standard Neo-Hooke
1501 *
1502 * @code
1503 *   return ( mu * ( symmetrize(F * transpose(F)) - I ) );
1504 *   }
1505 *   else
1506 *   {
1507 * @endcode
1508 *
1509 * Neo-Hooke in terms of principal stretches
1510 *
1511 * @code
1512 *   const SymmetricTensor<2, dim, NumberType>
1513 *   B = symmetrize(F * transpose(F));
1514 *   const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1515 *   eigen_B = eigenvectors(B, this->eigen_solver);
1516 *  
1517 *   SymmetricTensor<2, dim, NumberType> B_ev;
1518 *   for (unsigned int d=0; d<dim; ++d)
1519 *   B_ev += eigen_B[d].first*symmetrize(outer_product(eigen_B[d].second,eigen_B[d].second));
1520 *  
1521 *   return ( mu*(B_ev-I) );
1522 *   }
1523 *   }
1524 *   };
1525 *  
1526 * @endcode
1527 *
1528 *
1529 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassOgdenhyperelasticmaterial"></a>
1530 * <h4>Derived class: Ogden hyperelastic material</h4>
1531 *
1532 * @code
1533 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> >
1534 *   class Ogden : public Material_Hyperelastic < dim, NumberType >
1535 *   {
1536 *   public:
1537 *   Ogden(const Parameters::AllParameters &parameters,
1538 *   const Time &time)
1539 *   :
1540 *   Material_Hyperelastic< dim, NumberType > (parameters,time),
1541 *   mu({parameters.mu1_infty,
1542 *   parameters.mu2_infty,
1543 *   parameters.mu3_infty}),
1544 *   alpha({parameters.alpha1_infty,
1545 *   parameters.alpha2_infty,
1546 *   parameters.alpha3_infty})
1547 *   {}
1548 *   virtual ~Ogden()
1549 *   {}
1550 *  
1551 *   double
1552 *   get_viscous_dissipation() const override
1553 *   {
1554 *   return 0.0;
1555 *   }
1556 *  
1557 *   protected:
1558 *   std::vector<double> mu;
1559 *   std::vector<double> alpha;
1560 *  
1561 *   SymmetricTensor<2, dim, NumberType>
1562 *   get_tau_E_base(const Tensor<2,dim, NumberType> &F) const override
1563 *   {
1564 *   const SymmetricTensor<2, dim, NumberType>
1565 *   B = symmetrize(F * transpose(F));
1566 *  
1567 *   const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1568 *   eigen_B = eigenvectors(B, this->eigen_solver);
1569 *  
1570 *   SymmetricTensor<2, dim, NumberType> tau;
1571 *   static const SymmetricTensor< 2, dim, double>
1573 *  
1574 *   for (unsigned int i = 0; i < 3; ++i)
1575 *   {
1576 *   for (unsigned int A = 0; A < dim; ++A)
1577 *   {
1578 *   SymmetricTensor<2, dim, NumberType> tau_aux1 = symmetrize(
1579 *   outer_product(eigen_B[A].second,eigen_B[A].second));
1580 *   tau_aux1 *= mu[i]*std::pow(eigen_B[A].first, (alpha[i]/2.) );
1581 *   tau += tau_aux1;
1582 *   }
1583 *   SymmetricTensor<2, dim, NumberType> tau_aux2 (I);
1584 *   tau_aux2 *= mu[i];
1585 *   tau -= tau_aux2;
1586 *   }
1587 *   return tau;
1588 *   }
1589 *   };
1590 *  
1591 * @endcode
1592 *
1593 *
1594 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassSinglemodeOgdenviscoelasticmaterial"></a>
1595 * <h4>Derived class: Single-mode Ogden viscoelastic material</h4>
1596 * We use the finite viscoelastic model described in
1597 * Reese & Govindjee (1998) doi:10.1016/S0020-7683(97)00217-5
1598 * The algorithm for the implicit exponential time integration is given in
1599 * Budday et al. (2017) doi: 10.1016/j.actbio.2017.06.024
1600 *
1601 * @code
1602 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> >
1603 *   class visco_Ogden : public Material_Hyperelastic < dim, NumberType >
1604 *   {
1605 *   public:
1606 *   visco_Ogden(const Parameters::AllParameters &parameters,
1607 *   const Time &time)
1608 *   :
1609 *   Material_Hyperelastic< dim, NumberType > (parameters,time),
1610 *   mu_infty({parameters.mu1_infty,
1611 *   parameters.mu2_infty,
1612 *   parameters.mu3_infty}),
1613 *   alpha_infty({parameters.alpha1_infty,
1614 *   parameters.alpha2_infty,
1615 *   parameters.alpha3_infty}),
1616 *   mu_mode_1({parameters.mu1_mode_1,
1617 *   parameters.mu2_mode_1,
1618 *   parameters.mu3_mode_1}),
1619 *   alpha_mode_1({parameters.alpha1_mode_1,
1620 *   parameters.alpha2_mode_1,
1621 *   parameters.alpha3_mode_1}),
1622 *   viscosity_mode_1(parameters.viscosity_mode_1),
1624 *   Cinv_v_1_converged(Physics::Elasticity::StandardTensors<dim>::I)
1625 *   {}
1626 *   virtual ~visco_Ogden()
1627 *   {}
1628 *  
1629 *   void
1630 *   update_internal_equilibrium( const Tensor<2, dim, NumberType> &F ) override
1631 *   {
1632 *   Material_Hyperelastic < dim, NumberType >::update_internal_equilibrium(F);
1633 *  
1634 *   this->Cinv_v_1 = this->Cinv_v_1_converged;
1635 *   SymmetricTensor<2, dim, NumberType> B_e_1_tr = symmetrize(F * this->Cinv_v_1 * transpose(F));
1636 *  
1637 *   const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1638 *   eigen_B_e_1_tr = eigenvectors(B_e_1_tr, this->eigen_solver);
1639 *  
1640 *   Tensor< 1, dim, NumberType > lambdas_e_1_tr;
1641 *   Tensor< 1, dim, NumberType > epsilon_e_1_tr;
1642 *   for (int a = 0; a < dim; ++a)
1643 *   {
1644 *   lambdas_e_1_tr[a] = std::sqrt(eigen_B_e_1_tr[a].first);
1645 *   epsilon_e_1_tr[a] = std::log(lambdas_e_1_tr[a]);
1646 *   }
1647 *  
1648 *   const double tolerance = 1e-8;
1649 *   double residual_check = tolerance*10.0;
1650 *   Tensor< 1, dim, NumberType > residual;
1651 *   Tensor< 2, dim, NumberType > tangent;
1653 *   NumberType J_e_1 = std::sqrt(determinant(B_e_1_tr));
1654 *  
1655 *   std::vector<NumberType> lambdas_e_1_iso(dim);
1657 *   int iteration = 0;
1658 *  
1659 *   Tensor< 1, dim, NumberType > lambdas_e_1;
1660 *   Tensor< 1, dim, NumberType > epsilon_e_1;
1661 *   epsilon_e_1 = epsilon_e_1_tr;
1662 *  
1663 *   while(residual_check > tolerance)
1664 *   {
1665 *   NumberType aux_J_e_1 = 1.0;
1666 *   for (unsigned int a = 0; a < dim; ++a)
1667 *   {
1668 *   lambdas_e_1[a] = std::exp(epsilon_e_1[a]);
1669 *   aux_J_e_1 *= lambdas_e_1[a];
1670 *   }
1671 *  
1672 *   J_e_1 = aux_J_e_1;
1673 *  
1674 *   for (unsigned int a = 0; a < dim; ++a)
1675 *   lambdas_e_1_iso[a] = lambdas_e_1[a]*std::pow(J_e_1,-1.0/dim);
1676 *  
1677 *   for (unsigned int a = 0; a < dim; ++a)
1678 *   {
1679 *   residual[a] = get_beta_mode_1(lambdas_e_1_iso, a);
1680 *   residual[a] *= this->time.get_delta_t()/(2.0*viscosity_mode_1);
1681 *   residual[a] += epsilon_e_1[a];
1682 *   residual[a] -= epsilon_e_1_tr[a];
1683 *  
1684 *   for (unsigned int b = 0; b < dim; ++b)
1685 *   {
1686 *   tangent[a][b] = get_gamma_mode_1(lambdas_e_1_iso, a, b);
1687 *   tangent[a][b] *= this->time.get_delta_t()/(2.0*viscosity_mode_1);
1688 *   tangent[a][b] += I[a][b];
1689 *   }
1690 *  
1691 *   }
1692 *   epsilon_e_1 -= invert(tangent)*residual;
1693 *  
1694 *   residual_check = 0.0;
1695 *   for (unsigned int a = 0; a < dim; ++a)
1696 *   {
1697 *   if ( std::abs(residual[a]) > residual_check)
1698 *   residual_check = std::abs(Tensor<0,dim,double>(residual[a]));
1699 *   }
1700 *   iteration += 1;
1701 *   if (iteration > 15 )
1702 *   AssertThrow(false, ExcMessage("No convergence in local Newton iteration for the "
1703 *   "viscoelastic exponential time integration algorithm."));
1704 *   }
1705 *  
1706 *   NumberType aux_J_e_1 = 1.0;
1707 *   for (unsigned int a = 0; a < dim; ++a)
1708 *   {
1709 *   lambdas_e_1[a] = std::exp(epsilon_e_1[a]);
1710 *   aux_J_e_1 *= lambdas_e_1[a];
1711 *   }
1712 *   J_e_1 = aux_J_e_1;
1713 *  
1714 *   for (unsigned int a = 0; a < dim; ++a)
1715 *   lambdas_e_1_iso[a] = lambdas_e_1[a]*std::pow(J_e_1,-1.0/dim);
1716 *  
1717 *   for (unsigned int a = 0; a < dim; ++a)
1718 *   {
1720 *   B_e_1_aux = symmetrize(outer_product(eigen_B_e_1_tr[a].second,eigen_B_e_1_tr[a].second));
1721 *   B_e_1_aux *= lambdas_e_1[a] * lambdas_e_1[a];
1722 *   B_e_1 += B_e_1_aux;
1723 *   }
1724 *  
1725 *   Tensor<2, dim, NumberType>Cinv_v_1_AD = symmetrize(invert(F) * B_e_1 * invert(transpose(F)));
1726 *  
1727 *   this->tau_neq_1 = 0;
1728 *   for (unsigned int a = 0; a < dim; ++a)
1729 *   {
1731 *   tau_neq_1_aux = symmetrize(outer_product(eigen_B_e_1_tr[a].second,eigen_B_e_1_tr[a].second));
1732 *   tau_neq_1_aux *= get_beta_mode_1(lambdas_e_1_iso, a);
1733 *   this->tau_neq_1 += tau_neq_1_aux;
1734 *   }
1735 *  
1736 * @endcode
1737 *
1738 * Store history
1739 *
1740 * @code
1741 *   for (unsigned int a = 0; a < dim; ++a)
1742 *   for (unsigned int b = 0; b < dim; ++b)
1743 *   this->Cinv_v_1[a][b]= Tensor<0,dim,double>(Cinv_v_1_AD[a][b]);
1744 *   }
1745 *  
1746 *   void update_end_timestep() override
1747 *   {
1748 *   Material_Hyperelastic < dim, NumberType >::update_end_timestep();
1749 *   this->Cinv_v_1_converged = this->Cinv_v_1;
1750 *   }
1751 *  
1752 *   double get_viscous_dissipation() const override
1753 *   {
1754 *   NumberType dissipation_term = get_tau_E_neq() * get_tau_E_neq(); //Double contract the two SymmetricTensor
1755 *   dissipation_term /= (2*viscosity_mode_1);
1756 *  
1757 *   return dissipation_term.val();
1758 *   }
1759 *  
1760 *   protected:
1761 *   std::vector<double> mu_infty;
1762 *   std::vector<double> alpha_infty;
1763 *   std::vector<double> mu_mode_1;
1764 *   std::vector<double> alpha_mode_1;
1765 *   double viscosity_mode_1;
1766 *   SymmetricTensor<2, dim, double> Cinv_v_1;
1767 *   SymmetricTensor<2, dim, double> Cinv_v_1_converged;
1769 *  
1771 *   get_tau_E_base(const Tensor<2,dim, NumberType> &F) const override
1772 *   {
1773 *   return ( get_tau_E_neq() + get_tau_E_eq(F) );
1774 *   }
1775 *  
1777 *   get_tau_E_eq(const Tensor<2,dim, NumberType> &F) const
1778 *   {
1780 *  
1781 *   std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim > eigen_B;
1782 *   eigen_B = eigenvectors(B, this->eigen_solver);
1783 *  
1785 *   static const SymmetricTensor< 2, dim, double>
1787 *  
1788 *   for (unsigned int i = 0; i < 3; ++i)
1789 *   {
1790 *   for (unsigned int A = 0; A < dim; ++A)
1791 *   {
1793 *   outer_product(eigen_B[A].second,eigen_B[A].second));
1794 *   tau_aux1 *= mu_infty[i]*std::pow(eigen_B[A].first, (alpha_infty[i]/2.) );
1795 *   tau += tau_aux1;
1796 *   }
1797 *   SymmetricTensor<2, dim, NumberType> tau_aux2 (I);
1798 *   tau_aux2 *= mu_infty[i];
1799 *   tau -= tau_aux2;
1800 *   }
1801 *   return tau;
1802 *   }
1803 *  
1805 *   get_tau_E_neq() const
1806 *   {
1807 *   return tau_neq_1;
1808 *   }
1809 *  
1810 *   NumberType
1811 *   get_beta_mode_1(std::vector< NumberType > &lambda, const int &A) const
1812 *   {
1813 *   NumberType beta = 0.0;
1814 *  
1815 *   for (unsigned int i = 0; i < 3; ++i) //3rd-order Ogden model
1816 *   {
1817 *  
1818 *   NumberType aux = 0.0;
1819 *   for (int p = 0; p < dim; ++p)
1820 *   aux += std::pow(lambda[p],alpha_mode_1[i]);
1821 *  
1822 *   aux *= -1.0/dim;
1823 *   aux += std::pow(lambda[A], alpha_mode_1[i]);
1824 *   aux *= mu_mode_1[i];
1825 *  
1826 *   beta += aux;
1827 *   }
1828 *   return beta;
1829 *   }
1830 *  
1831 *   NumberType
1832 *   get_gamma_mode_1(std::vector< NumberType > &lambda,
1833 *   const int &A,
1834 *   const int &B ) const
1835 *   {
1836 *   NumberType gamma = 0.0;
1837 *  
1838 *   if (A==B)
1839 *   {
1840 *   for (unsigned int i = 0; i < 3; ++i)
1841 *   {
1842 *   NumberType aux = 0.0;
1843 *   for (int p = 0; p < dim; ++p)
1844 *   aux += std::pow(lambda[p],alpha_mode_1[i]);
1845 *  
1846 *   aux *= 1.0/(dim*dim);
1847 *   aux += 1.0/dim * std::pow(lambda[A], alpha_mode_1[i]);
1848 *   aux *= mu_mode_1[i]*alpha_mode_1[i];
1849 *  
1850 *   gamma += aux;
1851 *   }
1852 *   }
1853 *   else
1854 *   {
1855 *   for (unsigned int i = 0; i < 3; ++i)
1856 *   {
1857 *   NumberType aux = 0.0;
1858 *   for (int p = 0; p < dim; ++p)
1859 *   aux += std::pow(lambda[p],alpha_mode_1[i]);
1860 *  
1861 *   aux *= 1.0/(dim*dim);
1862 *   aux -= 1.0/dim * std::pow(lambda[A], alpha_mode_1[i]);
1863 *   aux -= 1.0/dim * std::pow(lambda[B], alpha_mode_1[i]);
1864 *   aux *= mu_mode_1[i]*alpha_mode_1[i];
1865 *  
1866 *   gamma += aux;
1867 *   }
1868 *   }
1869 *  
1870 *   return gamma;
1871 *   }
1872 *   };
1873 *  
1874 *  
1875 * @endcode
1876 *
1877 *
1878 * <a name="nonlinear-poro-viscoelasticity.cc-Constitutiveequationforthefluidcomponentofthebiphasicmaterial"></a>
1879 * <h3>Constitutive equation for the fluid component of the biphasic material</h3>
1880 * We consider two slightly different definitions to define the seepage velocity with a Darcy-like law.
1881 * Ehlers & Eipper 1999, doi:10.1023/A:1006565509095
1882 * Markert 2007, doi:10.1007/s11242-007-9107-6
1883 * The selection of one or another is made by the user via the parameters file.
1884 *
1885 * @code
1886 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> >
1887 *   class Material_Darcy_Fluid
1888 *   {
1889 *   public:
1890 *   Material_Darcy_Fluid(const Parameters::AllParameters &parameters)
1891 *   :
1892 *   fluid_type(parameters.fluid_type),
1893 *   n_OS(parameters.solid_vol_frac),
1894 *   initial_intrinsic_permeability(parameters.init_intrinsic_perm),
1895 *   viscosity_FR(parameters.viscosity_FR),
1896 *   initial_darcy_coefficient(parameters.init_darcy_coef),
1897 *   weight_FR(parameters.weight_FR),
1898 *   kappa_darcy(parameters.kappa_darcy),
1899 *   gravity_term(parameters.gravity_term),
1900 *   density_FR(parameters.density_FR),
1901 *   gravity_direction(parameters.gravity_direction),
1902 *   gravity_value(parameters.gravity_value)
1903 *   {
1904 *   Assert(kappa_darcy >= 0, ExcInternalError());
1905 *   }
1906 *   ~Material_Darcy_Fluid()
1907 *   {}
1908 *  
1909 *   Tensor<1, dim, NumberType> get_seepage_velocity_current
1910 *   (const Tensor<2,dim, NumberType> &F,
1911 *   const Tensor<1,dim, NumberType> &grad_p_fluid) const
1912 *   {
1913 *   const NumberType det_F = determinant(F);
1914 *   Assert(det_F > 0.0, ExcInternalError());
1915 *  
1916 *   Tensor<2, dim, NumberType> permeability_term;
1917 *  
1918 *   if (fluid_type == "Markert")
1919 *   permeability_term = get_instrinsic_permeability_current(F) / viscosity_FR;
1920 *  
1921 *   else if (fluid_type == "Ehlers")
1922 *   permeability_term = get_darcy_flow_current(F) / weight_FR;
1923 *  
1924 *   else
1925 *   AssertThrow(false, ExcMessage(
1926 *   "Material_Darcy_Fluid --> Only Markert "
1927 *   "and Ehlers formulations have been implemented."));
1928 *  
1929 *   return ( -1.0 * permeability_term * det_F
1930 *   * (grad_p_fluid - get_body_force_FR_current()) );
1931 *   }
1932 *  
1933 *   double get_porous_dissipation(const Tensor<2,dim, NumberType> &F,
1934 *   const Tensor<1,dim, NumberType> &grad_p_fluid) const
1935 *   {
1936 *   NumberType dissipation_term;
1937 *   Tensor<1, dim, NumberType> seepage_velocity;
1938 *   Tensor<2, dim, NumberType> permeability_term;
1939 *  
1940 *   const NumberType det_F = determinant(F);
1941 *   Assert(det_F > 0.0, ExcInternalError());
1942 *  
1943 *   if (fluid_type == "Markert")
1944 *   {
1945 *   permeability_term = get_instrinsic_permeability_current(F) / viscosity_FR;
1946 *   seepage_velocity = get_seepage_velocity_current(F,grad_p_fluid);
1947 *   }
1948 *   else if (fluid_type == "Ehlers")
1949 *   {
1950 *   permeability_term = get_darcy_flow_current(F) / weight_FR;
1951 *   seepage_velocity = get_seepage_velocity_current(F,grad_p_fluid);
1952 *   }
1953 *   else
1954 *   AssertThrow(false, ExcMessage(
1955 *   "Material_Darcy_Fluid --> Only Markert and Ehlers "
1956 *   "formulations have been implemented."));
1957 *  
1958 *   dissipation_term = ( invert(permeability_term) * seepage_velocity ) * seepage_velocity;
1959 *   dissipation_term *= 1.0/(det_F*det_F);
1960 *   return Tensor<0,dim,double>(dissipation_term);
1961 *   }
1962 *  
1963 *   protected:
1964 *   const std::string fluid_type;
1965 *   const double n_OS;
1966 *   const double initial_intrinsic_permeability;
1967 *   const double viscosity_FR;
1968 *   const double initial_darcy_coefficient;
1969 *   const double weight_FR;
1970 *   const double kappa_darcy;
1971 *   const bool gravity_term;
1972 *   const double density_FR;
1973 *   const int gravity_direction;
1974 *   const double gravity_value;
1975 *  
1977 *   get_instrinsic_permeability_current(const Tensor<2,dim, NumberType> &F) const
1978 *   {
1979 *   static const SymmetricTensor< 2, dim, double>
1981 *   const Tensor<2, dim, NumberType> initial_instrinsic_permeability_tensor
1982 *   = Tensor<2, dim, double>(initial_intrinsic_permeability * I);
1983 *  
1984 *   const NumberType det_F = determinant(F);
1985 *   Assert(det_F > 0.0, ExcInternalError());
1986 *  
1987 *   const NumberType fraction = (det_F - n_OS)/(1 - n_OS);
1988 *   return ( NumberType (std::pow(fraction, kappa_darcy))
1989 *   * initial_instrinsic_permeability_tensor );
1990 *   }
1991 *  
1993 *   get_darcy_flow_current(const Tensor<2,dim, NumberType> &F) const
1994 *   {
1995 *   static const SymmetricTensor< 2, dim, double>
1997 *   const Tensor<2, dim, NumberType> initial_darcy_flow_tensor
1998 *   = Tensor<2, dim, double>(initial_darcy_coefficient * I);
1999 *  
2000 *   const NumberType det_F = determinant(F);
2001 *   Assert(det_F > 0.0, ExcInternalError());
2002 *  
2003 *   const NumberType fraction = (1.0 - (n_OS / det_F) )/(1.0 - n_OS);
2004 *   return ( NumberType (std::pow(fraction, kappa_darcy))
2005 *   * initial_darcy_flow_tensor);
2006 *   }
2007 *  
2009 *   get_body_force_FR_current() const
2010 *   {
2011 *   Tensor<1, dim, NumberType> body_force_FR_current;
2012 *  
2013 *   if (gravity_term == true)
2014 *   {
2015 *   Tensor<1, dim, NumberType> gravity_vector;
2016 *   gravity_vector[gravity_direction] = gravity_value;
2017 *   body_force_FR_current = density_FR * gravity_vector;
2018 *   }
2019 *   return body_force_FR_current;
2020 *   }
2021 *   };
2022 *  
2023 * @endcode
2024 *
2025 *
2026 * <a name="nonlinear-poro-viscoelasticity.cc-Quadraturepointhistory"></a>
2027 * <h3>Quadrature point history</h3>
2028 * As seen in @ref step_18 "step-18", the <code> PointHistory </code> class offers a method
2029 * for storing data at the quadrature points. Here each quadrature point
2030 * holds a pointer to a material description. Thus, different material models
2031 * can be used in different regions of the domain. Among other data, we
2032 * choose to store the "extra" Kirchhoff stress @f$\boldsymbol{\tau}_E@f$ and
2033 * the dissipation values @f$\mathcal{D}_p@f$ and @f$\mathcal{D}_v@f$.
2034 *
2035 * @code
2036 *   template <int dim, typename NumberType = Sacado::Fad::DFad<double> > //double>
2037 *   class PointHistory
2038 *   {
2039 *   public:
2040 *   PointHistory()
2041 *   {}
2042 *  
2043 *   virtual ~PointHistory()
2044 *   {}
2045 *  
2046 *   void setup_lqp (const Parameters::AllParameters &parameters,
2047 *   const Time &time)
2048 *   {
2049 *   if (parameters.mat_type == "Neo-Hooke")
2050 *   solid_material.reset(new NeoHooke<dim,NumberType>(parameters,time));
2051 *   else if (parameters.mat_type == "Ogden")
2052 *   solid_material.reset(new Ogden<dim,NumberType>(parameters,time));
2053 *   else if (parameters.mat_type == "visco-Ogden")
2054 *   solid_material.reset(new visco_Ogden<dim,NumberType>(parameters,time));
2055 *   else
2056 *   Assert (false, ExcMessage("Material type not implemented"));
2057 *  
2058 *   fluid_material.reset(new Material_Darcy_Fluid<dim,NumberType>(parameters));
2059 *   }
2060 *  
2062 *   get_tau_E(const Tensor<2, dim, NumberType> &F) const
2063 *   {
2064 *   return solid_material->get_tau_E(F);
2065 *   }
2066 *  
2068 *   get_Cauchy_E(const Tensor<2, dim, NumberType> &F) const
2069 *   {
2070 *   return solid_material->get_Cauchy_E(F);
2071 *   }
2072 *  
2073 *   double
2074 *   get_converged_det_F() const
2075 *   {
2076 *   return solid_material->get_converged_det_F();
2077 *   }
2078 *  
2079 *   void
2080 *   update_end_timestep()
2081 *   {
2082 *   solid_material->update_end_timestep();
2083 *   }
2084 *  
2085 *   void
2086 *   update_internal_equilibrium(const Tensor<2, dim, NumberType> &F )
2087 *   {
2088 *   solid_material->update_internal_equilibrium(F);
2089 *   }
2090 *  
2091 *   double
2092 *   get_viscous_dissipation() const
2093 *   {
2094 *   return solid_material->get_viscous_dissipation();
2095 *   }
2096 *  
2098 *   get_seepage_velocity_current (const Tensor<2,dim, NumberType> &F,
2099 *   const Tensor<1,dim, NumberType> &grad_p_fluid) const
2100 *   {
2101 *   return fluid_material->get_seepage_velocity_current(F, grad_p_fluid);
2102 *   }
2103 *  
2104 *   double
2105 *   get_porous_dissipation(const Tensor<2,dim, NumberType> &F,
2106 *   const Tensor<1,dim, NumberType> &grad_p_fluid) const
2107 *   {
2108 *   return fluid_material->get_porous_dissipation(F, grad_p_fluid);
2109 *   }
2110 *  
2112 *   get_overall_body_force (const Tensor<2,dim, NumberType> &F,
2113 *   const Parameters::AllParameters &parameters) const
2114 *   {
2115 *   Tensor<1, dim, NumberType> body_force;
2116 *  
2117 *   if (parameters.gravity_term == true)
2118 *   {
2119 *   const NumberType det_F_AD = determinant(F);
2120 *   Assert(det_F_AD > 0.0, ExcInternalError());
2121 *  
2122 *   const NumberType overall_density_ref
2123 *   = parameters.density_SR * parameters.solid_vol_frac
2124 *   + parameters.density_FR
2125 *   * (det_F_AD - parameters.solid_vol_frac);
2126 *  
2127 *   Tensor<1, dim, NumberType> gravity_vector;
2128 *   gravity_vector[parameters.gravity_direction] = parameters.gravity_value;
2129 *   body_force = overall_density_ref * gravity_vector;
2130 *   }
2131 *  
2132 *   return body_force;
2133 *   }
2134 *   private:
2135 *   std::shared_ptr< Material_Hyperelastic<dim, NumberType> > solid_material;
2136 *   std::shared_ptr< Material_Darcy_Fluid<dim, NumberType> > fluid_material;
2137 *   };
2138 *  
2139 * @endcode
2140 *
2141 *
2142 * <a name="nonlinear-poro-viscoelasticity.cc-Nonlinearporoviscoelasticsolid"></a>
2143 * <h3>Nonlinear poro-viscoelastic solid</h3>
2144 * The Solid class is the central class as it represents the problem at hand:
2145 * the nonlinear poro-viscoelastic solid
2146 *
2147 * @code
2148 *   template <int dim>
2149 *   class Solid
2150 *   {
2151 *   public:
2152 *   Solid(const Parameters::AllParameters &parameters);
2153 *   virtual ~Solid();
2154 *   void run();
2155 *  
2156 *   protected:
2157 *   using ADNumberType = Sacado::Fad::DFad<double>;
2158 *  
2159 *   std::ofstream outfile;
2160 *   std::ofstream pointfile;
2161 *  
2162 *   struct PerTaskData_ASM;
2163 *   template<typename NumberType = double> struct ScratchData_ASM;
2164 *  
2165 * @endcode
2166 *
2167 * Generate mesh
2168 *
2169 * @code
2170 *   virtual void make_grid() = 0;
2171 *  
2172 * @endcode
2173 *
2174 * Define points for post-processing
2175 *
2176 * @code
2177 *   virtual void define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) = 0;
2178 *  
2179 * @endcode
2180 *
2181 * Set up the finite element system to be solved:
2182 *
2183 * @code
2184 *   void system_setup(TrilinosWrappers::MPI::BlockVector &solution_delta_OUT);
2185 *  
2186 * @endcode
2187 *
2188 * Extract sub-blocks from the global matrix
2189 *
2190 * @code
2191 *   void determine_component_extractors();
2192 *  
2193 * @endcode
2194 *
2195 * Several functions to assemble the system and right hand side matrices using multithreading.
2196 *
2197 * @code
2198 *   void assemble_system
2199 *   (const TrilinosWrappers::MPI::BlockVector &solution_delta_OUT );
2200 *   void assemble_system_one_cell
2201 *   (const typename DoFHandler<dim>::active_cell_iterator &cell,
2202 *   ScratchData_ASM<ADNumberType> &scratch,
2203 *   PerTaskData_ASM &data) const;
2204 *   void copy_local_to_global_system(const PerTaskData_ASM &data);
2205 *  
2206 * @endcode
2207 *
2208 * Define boundary conditions
2209 *
2210 * @code
2211 *   virtual void make_constraints(const int &it_nr);
2212 *   virtual void make_dirichlet_constraints(AffineConstraints<double> &constraints) = 0;
2213 *   virtual Tensor<1,dim> get_neumann_traction
2214 *   (const types::boundary_id &boundary_id,
2215 *   const Point<dim> &pt,
2216 *   const Tensor<1,dim> &N) const = 0;
2217 *   virtual double get_prescribed_fluid_flow
2218 *   (const types::boundary_id &boundary_id,
2219 *   const Point<dim> &pt) const = 0;
2220 *   virtual types::boundary_id
2221 *   get_reaction_boundary_id_for_output () const = 0;
2222 *   virtual std::pair<types::boundary_id,types::boundary_id>
2223 *   get_drained_boundary_id_for_output () const = 0;
2224 *   virtual std::vector<double> get_dirichlet_load
2225 *   (const types::boundary_id &boundary_id,
2226 *   const int &direction) const = 0;
2227 *  
2228 * @endcode
2229 *
2230 * Create and update the quadrature points.
2231 *
2232 * @code
2233 *   void setup_qph();
2234 *  
2235 * @endcode
2236 *
2237 * Solve non-linear system using a Newton-Raphson scheme
2238 *
2239 * @code
2240 *   void solve_nonlinear_timestep(TrilinosWrappers::MPI::BlockVector &solution_delta_OUT);
2241 *  
2242 * @endcode
2243 *
2244 * Solve the linearized equations using a direct solver
2245 *
2246 * @code
2247 *   void solve_linear_system ( TrilinosWrappers::MPI::BlockVector &newton_update_OUT);
2248 *  
2249 * @endcode
2250 *
2251 * Retrieve the solution
2252 *
2253 * @code
2254 *   TrilinosWrappers::MPI::BlockVector
2255 *   get_total_solution(const TrilinosWrappers::MPI::BlockVector &solution_delta_IN) const;
2256 *  
2257 * @endcode
2258 *
2259 * Store the converged values of the internal variables at the end of each timestep
2260 *
2261 * @code
2262 *   void update_end_timestep();
2263 *  
2264 * @endcode
2265 *
2266 * Post-processing and writing data to files
2267 *
2268 * @code
2269 *   void output_results_to_vtu(const unsigned int timestep,
2270 *   const double current_time,
2271 *   TrilinosWrappers::MPI::BlockVector solution) const;
2272 *   void output_results_to_plot(const unsigned int timestep,
2273 *   const double current_time,
2274 *   TrilinosWrappers::MPI::BlockVector solution,
2275 *   std::vector<Point<dim> > &tracked_vertices,
2276 *   std::ofstream &pointfile) const;
2277 *  
2278 * @endcode
2279 *
2280 * Headers and footer for the output files
2281 *
2282 * @code
2283 *   void print_console_file_header( std::ofstream &outfile) const;
2284 *   void print_plot_file_header(std::vector<Point<dim> > &tracked_vertices,
2285 *   std::ofstream &pointfile) const;
2286 *   void print_console_file_footer(std::ofstream &outfile) const;
2287 *   void print_plot_file_footer( std::ofstream &pointfile) const;
2288 *  
2289 * @endcode
2290 *
2291 * For parallel communication
2292 *
2293 * @code
2294 *   MPI_Comm mpi_communicator;
2295 *   const unsigned int n_mpi_processes;
2296 *   const unsigned int this_mpi_process;
2297 *   mutable ConditionalOStream pcout;
2298 *  
2299 * @endcode
2300 *
2301 * A collection of the parameters used to describe the problem setup
2302 *
2303 * @code
2304 *   const Parameters::AllParameters &parameters;
2305 *  
2306 * @endcode
2307 *
2308 * Declare an instance of dealii Triangulation class (mesh)
2309 *
2310 * @code
2311 *   parallel::shared::Triangulation<dim> triangulation;
2312 *  
2313 * @endcode
2314 *
2315 * Keep track of the current time and the time spent evaluating certain functions
2316 *
2317 * @code
2318 *   Time time;
2319 *   TimerOutput timerconsole;
2320 *   TimerOutput timerfile;
2321 *  
2322 * @endcode
2323 *
2324 * A storage object for quadrature point information.
2325 *
2326 * @code
2327 *   CellDataStorage<typename Triangulation<dim>::cell_iterator, PointHistory<dim,ADNumberType> > quadrature_point_history;
2328 *  
2329 * @endcode
2330 *
2331 * Integers to store polynomial degree (needed for output)
2332 *
2333 * @code
2334 *   const unsigned int degree_displ;
2335 *   const unsigned int degree_pore;
2336 *  
2337 * @endcode
2338 *
2339 * Declare an instance of dealii FESystem class (finite element definition)
2340 *
2341 * @code
2342 *   const FESystem<dim> fe;
2343 *  
2344 * @endcode
2345 *
2346 * Declare an instance of dealii DoFHandler class (assign DoFs to mesh)
2347 *
2348 * @code
2349 *   DoFHandler<dim> dof_handler_ref;
2350 *  
2351 * @endcode
2352 *
2353 * Integer to store DoFs per element (this value will be used often)
2354 *
2355 * @code
2356 *   const unsigned int dofs_per_cell;
2357 *  
2358 * @endcode
2359 *
2360 * Declare an instance of dealii Extractor objects used to retrieve information from the solution vectors
2361 * We will use "u_fe" and "p_fluid_fe"as subscript in operator [] expressions on FEValues and FEFaceValues
2362 * objects to extract the components of the displacement vector and fluid pressure, respectively.
2363 *
2364 * @code
2365 *   const FEValuesExtractors::Vector u_fe;
2366 *   const FEValuesExtractors::Scalar p_fluid_fe;
2367 *  
2368 * @endcode
2369 *
2370 * Description of how the block-system is arranged. There are 3 blocks:
2371 * 0 - vector DOF displacements u
2372 * 1 - scalar DOF fluid pressure p_fluid
2373 *
2374 * @code
2375 *   static const unsigned int n_blocks = 2;
2376 *   static const unsigned int n_components = dim+1;
2377 *   static const unsigned int first_u_component = 0;
2378 *   static const unsigned int p_fluid_component = dim;
2379 *  
2380 *   enum
2381 *   {
2382 *   u_block = 0,
2383 *   p_fluid_block = 1
2384 *   };
2385 *  
2386 * @endcode
2387 *
2388 * Extractors
2389 *
2390 * @code
2391 *   const FEValuesExtractors::Scalar x_displacement;
2392 *   const FEValuesExtractors::Scalar y_displacement;
2393 *   const FEValuesExtractors::Scalar z_displacement;
2394 *   const FEValuesExtractors::Scalar pressure;
2395 *  
2396 * @endcode
2397 *
2398 * Block data
2399 *
2400 * @code
2401 *   std::vector<unsigned int> block_component;
2402 *  
2403 * @endcode
2404 *
2405 * DoF index data
2406 *
2407 * @code
2408 *   std::vector<IndexSet> all_locally_owned_dofs;
2409 *   IndexSet locally_owned_dofs;
2410 *   IndexSet locally_relevant_dofs;
2411 *   std::vector<IndexSet> locally_owned_partitioning;
2412 *   std::vector<IndexSet> locally_relevant_partitioning;
2413 *  
2414 *   std::vector<types::global_dof_index> dofs_per_block;
2415 *   std::vector<types::global_dof_index> element_indices_u;
2416 *   std::vector<types::global_dof_index> element_indices_p_fluid;
2417 *  
2418 * @endcode
2419 *
2420 * Declare an instance of dealii QGauss class (The Gauss-Legendre family of quadrature rules for numerical integration)
2421 * Gauss Points in element, with n quadrature points (in each of the dim space directions)
2422 *
2423 * @code
2424 *   const QGauss<dim> qf_cell;
2425 * @endcode
2426 *
2427 * Gauss Points on element faces (used for definition of BCs)
2428 *
2429 * @code
2430 *   const QGauss<dim - 1> qf_face;
2431 * @endcode
2432 *
2433 * Integer to store num GPs per element (this value will be used often)
2434 *
2435 * @code
2436 *   const unsigned int n_q_points;
2437 * @endcode
2438 *
2439 * Integer to store num GPs per face (this value will be used often)
2440 *
2441 * @code
2442 *   const unsigned int n_q_points_f;
2443 *  
2444 * @endcode
2445 *
2446 * Declare an instance of dealii AffineConstraints class (linear constraints on DoFs due to hanging nodes or BCs)
2447 *
2448 * @code
2449 *   AffineConstraints<double> constraints;
2450 *  
2451 * @endcode
2452 *
2453 * Declare an instance of dealii classes necessary for FE system set-up and assembly
2454 * Store elements of tangent matrix (indicated by SparsityPattern class) as sparse matrix (more efficient)
2455 *
2456 * @code
2457 *   TrilinosWrappers::BlockSparseMatrix tangent_matrix;
2458 *   TrilinosWrappers::BlockSparseMatrix tangent_matrix_preconditioner;
2459 * @endcode
2460 *
2461 * Right hand side vector of forces
2462 *
2463 * @code
2464 *   TrilinosWrappers::MPI::BlockVector system_rhs;
2465 * @endcode
2466 *
2467 * Total displacement values + pressure (accumulated solution to FE system)
2468 *
2469 * @code
2470 *   TrilinosWrappers::MPI::BlockVector solution_n;
2471 *  
2472 * @endcode
2473 *
2474 * Non-block system for the direct solver. We will copy the block system into these to solve the linearized system of equations.
2475 *
2476 * @code
2477 *   TrilinosWrappers::SparseMatrix tangent_matrix_nb;
2478 *   TrilinosWrappers::MPI::Vector system_rhs_nb;
2479 *  
2480 * @endcode
2481 *
2482 * We define variables to store norms and update norms and normalisation factors.
2483 *
2484 * @code
2485 *   struct Errors
2486 *   {
2487 *   Errors()
2488 *   :
2489 *   norm(1.0), u(1.0), p_fluid(1.0)
2490 *   {}
2491 *  
2492 *   void reset()
2493 *   {
2494 *   norm = 1.0;
2495 *   u = 1.0;
2496 *   p_fluid = 1.0;
2497 *   }
2498 *   void normalise(const Errors &rhs)
2499 *   {
2500 *   if (rhs.norm != 0.0)
2501 *   norm /= rhs.norm;
2502 *   if (rhs.u != 0.0)
2503 *   u /= rhs.u;
2504 *   if (rhs.p_fluid != 0.0)
2505 *   p_fluid /= rhs.p_fluid;
2506 *   }
2507 *  
2508 *   double norm, u, p_fluid;
2509 *   };
2510 *  
2511 * @endcode
2512 *
2513 * Declare several instances of the "Error" structure
2514 *
2515 * @code
2516 *   Errors error_residual, error_residual_0, error_residual_norm, error_update,
2517 *   error_update_0, error_update_norm;
2518 *  
2519 * @endcode
2520 *
2521 * Methods to calculate error measures
2522 *
2523 * @code
2524 *   void get_error_residual(Errors &error_residual_OUT);
2525 *   void get_error_update
2526 *   (const TrilinosWrappers::MPI::BlockVector &newton_update_IN,
2527 *   Errors &error_update_OUT);
2528 *  
2529 * @endcode
2530 *
2531 * Print information to screen
2532 *
2533 * @code
2534 *   void print_conv_header();
2535 *   void print_conv_footer();
2536 *  
2537 * @endcode
2538 *
2539 * NOTE: In all functions, we pass by reference (&), so these functions work on the original copy (not a clone copy),
2540 * modifying the input variables inside the functions will change them outside the function.
2541 *
2542 * @code
2543 *   };
2544 *  
2545 * @endcode
2546 *
2547 *
2548 * <a name="nonlinear-poro-viscoelasticity.cc-ImplementationofthecodeSolidcodeclass"></a>
2549 * <h3>Implementation of the <code>Solid</code> class</h3>
2550 *
2551 * <a name="nonlinear-poro-viscoelasticity.cc-Publicinterface"></a>
2552 * <h4>Public interface</h4>
2553 * We initialise the Solid class using data extracted from the parameter file.
2554 *
2555 * @code
2556 *   template <int dim>
2557 *   Solid<dim>::Solid(const Parameters::AllParameters &parameters)
2558 *   :
2559 *   mpi_communicator(MPI_COMM_WORLD),
2560 *   n_mpi_processes (Utilities::MPI::n_mpi_processes(mpi_communicator)),
2561 *   this_mpi_process (Utilities::MPI::this_mpi_process(mpi_communicator)),
2562 *   pcout(std::cout, this_mpi_process == 0),
2563 *   parameters(parameters),
2564 *   triangulation(mpi_communicator,Triangulation<dim>::maximum_smoothing),
2565 *   time(parameters.end_time, parameters.delta_t),
2566 *   timerconsole( mpi_communicator,
2567 *   pcout,
2570 *   timerfile( mpi_communicator,
2571 *   outfile,
2574 *   degree_displ(parameters.poly_degree_displ),
2575 *   degree_pore(parameters.poly_degree_pore),
2576 *   fe( FE_Q<dim>(parameters.poly_degree_displ), dim,
2577 *   FE_Q<dim>(parameters.poly_degree_pore), 1 ),
2578 *   dof_handler_ref(triangulation),
2579 *   dofs_per_cell (fe.dofs_per_cell),
2580 *   u_fe(first_u_component),
2581 *   p_fluid_fe(p_fluid_component),
2582 *   x_displacement(first_u_component),
2583 *   y_displacement(first_u_component+1),
2584 *   z_displacement(first_u_component+2),
2585 *   pressure(p_fluid_component),
2586 *   dofs_per_block(n_blocks),
2587 *   qf_cell(parameters.quad_order),
2588 *   qf_face(parameters.quad_order),
2589 *   n_q_points (qf_cell.size()),
2590 *   n_q_points_f (qf_face.size())
2591 *   {
2592 *   Assert(dim==3, ExcMessage("This problem only works in 3 space dimensions."));
2593 *   determine_component_extractors();
2594 *   }
2595 *  
2596 * @endcode
2597 *
2598 * The class destructor simply clears the data held by the DOFHandler
2599 *
2600 * @code
2601 *   template <int dim>
2602 *   Solid<dim>::~Solid()
2603 *   {
2604 *   dof_handler_ref.clear();
2605 *   }
2606 *  
2607 * @endcode
2608 *
2609 * Runs the 3D solid problem
2610 *
2611 * @code
2612 *   template <int dim>
2613 *   void Solid<dim>::run()
2614 *   {
2615 * @endcode
2616 *
2617 * The current solution increment is defined as a block vector to reflect the structure
2618 * of the PDE system, with multiple solution components
2619 *
2620 * @code
2621 *   TrilinosWrappers::MPI::BlockVector solution_delta;
2622 *  
2623 * @endcode
2624 *
2625 * Open file
2626 *
2627 * @code
2628 *   if (this_mpi_process == 0)
2629 *   {
2630 *   outfile.open("console-output.sol");
2631 *   print_console_file_header(outfile);
2632 *   }
2633 *  
2634 * @endcode
2635 *
2636 * Generate mesh
2637 *
2638 * @code
2639 *   make_grid();
2640 *  
2641 * @endcode
2642 *
2643 * Assign DOFs and create the stiffness and right-hand-side force vector
2644 *
2645 * @code
2646 *   system_setup(solution_delta);
2647 *  
2648 * @endcode
2649 *
2650 * Define points for post-processing
2651 *
2652 * @code
2653 *   std::vector<Point<dim> > tracked_vertices (2);
2654 *   define_tracked_vertices(tracked_vertices);
2655 *   std::vector<Point<dim>> reaction_force;
2656 *  
2657 *   if (this_mpi_process == 0)
2658 *   {
2659 *   pointfile.open("data-for-gnuplot.sol");
2660 *   print_plot_file_header(tracked_vertices, pointfile);
2661 *   }
2662 *  
2663 * @endcode
2664 *
2665 * Print results to output file
2666 *
2667 * @code
2668 *   if (parameters.outfiles_requested == "true")
2669 *   {
2670 *   output_results_to_vtu(time.get_timestep(),
2671 *   time.get_current(),
2672 *   solution_n );
2673 *   }
2674 *  
2675 *   output_results_to_plot(time.get_timestep(),
2676 *   time.get_current(),
2677 *   solution_n,
2678 *   tracked_vertices,
2679 *   pointfile);
2680 *  
2681 * @endcode
2682 *
2683 * Increment time step (=load step)
2684 * NOTE: In solving the quasi-static problem, the time becomes a loading parameter,
2685 * i.e. we increase the loading linearly with time, making the two concepts interchangeable.
2686 *
2687 * @code
2688 *   time.increment_time();
2689 *  
2690 * @endcode
2691 *
2692 * Print information on screen
2693 *
2694 * @code
2695 *   pcout << "\nSolver:";
2696 *   pcout << "\n CST = make constraints";
2697 *   pcout << "\n ASM_SYS = assemble system";
2698 *   pcout << "\n SLV = linear solver \n";
2699 *  
2700 * @endcode
2701 *
2702 * Print information on file
2703 *
2704 * @code
2705 *   outfile << "\nSolver:";
2706 *   outfile << "\n CST = make constraints";
2707 *   outfile << "\n ASM_SYS = assemble system";
2708 *   outfile << "\n SLV = linear solver \n";
2709 *  
2710 *   while ( (time.get_end() - time.get_current()) > -1.0*parameters.tol_u )
2711 *   {
2712 * @endcode
2713 *
2714 * Initialize the current solution increment to zero
2715 *
2716 * @code
2717 *   solution_delta = 0.0;
2718 *  
2719 * @endcode
2720 *
2721 * Solve the non-linear system using a Newton-Rapshon scheme
2722 *
2723 * @code
2724 *   solve_nonlinear_timestep(solution_delta);
2725 *  
2726 * @endcode
2727 *
2728 * Add the computed solution increment to total solution
2729 *
2730 * @code
2731 *   solution_n += solution_delta;
2732 *  
2733 * @endcode
2734 *
2735 * Store the converged values of the internal variables
2736 *
2737 * @code
2738 *   update_end_timestep();
2739 *  
2740 * @endcode
2741 *
2742 * Output results
2743 *
2744 * @code
2745 *   if (( (time.get_timestep()%parameters.timestep_output) == 0 )
2746 *   && (parameters.outfiles_requested == "true") )
2747 *   {
2748 *   output_results_to_vtu(time.get_timestep(),
2749 *   time.get_current(),
2750 *   solution_n );
2751 *   }
2752 *  
2753 *   output_results_to_plot(time.get_timestep(),
2754 *   time.get_current(),
2755 *   solution_n,
2756 *   tracked_vertices,
2757 *   pointfile);
2758 *  
2759 * @endcode
2760 *
2761 * Increment the time step (=load step)
2762 *
2763 * @code
2764 *   time.increment_time();
2765 *   }
2766 *  
2767 * @endcode
2768 *
2769 * Print the footers and close files
2770 *
2771 * @code
2772 *   if (this_mpi_process == 0)
2773 *   {
2774 *   print_plot_file_footer(pointfile);
2775 *   pointfile.close ();
2776 *   print_console_file_footer(outfile);
2777 *  
2778 * @endcode
2779 *
2780 * NOTE: ideally, we should close the outfile here [ >> outfile.close (); ]
2781 * But if we do, then the timer output will not be printed. That is why we leave it open.
2782 *
2783 * @code
2784 *   }
2785 *   }
2786 *  
2787 * @endcode
2788 *
2789 *
2790 * <a name="nonlinear-poro-viscoelasticity.cc-Privateinterface"></a>
2791 * <h4>Private interface</h4>
2792 * We define the structures needed for parallelization with Threading Building Blocks (TBB)
2793 * Tangent matrix and right-hand side force vector assembly structures.
2794 * PerTaskData_ASM stores local contributions
2795 *
2796 * @code
2797 *   template <int dim>
2798 *   struct Solid<dim>::PerTaskData_ASM
2799 *   {
2800 *   FullMatrix<double> cell_matrix;
2801 *   Vector<double> cell_rhs;
2802 *   std::vector<types::global_dof_index> local_dof_indices;
2803 *  
2804 *   PerTaskData_ASM(const unsigned int dofs_per_cell)
2805 *   :
2806 *   cell_matrix(dofs_per_cell, dofs_per_cell),
2807 *   cell_rhs(dofs_per_cell),
2808 *   local_dof_indices(dofs_per_cell)
2809 *   {}
2810 *  
2811 *   void reset()
2812 *   {
2813 *   cell_matrix = 0.0;
2814 *   cell_rhs = 0.0;
2815 *   }
2816 *   };
2817 *  
2818 * @endcode
2819 *
2820 * ScratchData_ASM stores larger objects used during the assembly
2821 *
2822 * @code
2823 *   template <int dim>
2824 *   template <typename NumberType>
2825 *   struct Solid<dim>::ScratchData_ASM
2826 *   {
2827 *   const TrilinosWrappers::MPI::BlockVector &solution_total;
2828 *  
2829 * @endcode
2830 *
2831 * Integration helper
2832 *
2833 * @code
2834 *   FEValues<dim> fe_values_ref;
2835 *   FEFaceValues<dim> fe_face_values_ref;
2836 *  
2837 * @endcode
2838 *
2839 * Quadrature point solution
2840 *
2841 * @code
2842 *   std::vector<NumberType> local_dof_values;
2843 *   std::vector<Tensor<2, dim, NumberType> > solution_grads_u_total;
2844 *   std::vector<NumberType> solution_values_p_fluid_total;
2845 *   std::vector<Tensor<1, dim, NumberType> > solution_grads_p_fluid_total;
2846 *   std::vector<Tensor<1, dim, NumberType> > solution_grads_face_p_fluid_total;
2847 *  
2848 * @endcode
2849 *
2850 * shape function values
2851 *
2852 * @code
2853 *   std::vector<std::vector<Tensor<1,dim>>> Nx;
2854 *   std::vector<std::vector<double>> Nx_p_fluid;
2855 * @endcode
2856 *
2857 * shape function gradients
2858 *
2859 * @code
2860 *   std::vector<std::vector<Tensor<2,dim, NumberType>>> grad_Nx;
2861 *   std::vector<std::vector<SymmetricTensor<2,dim, NumberType>>> symm_grad_Nx;
2862 *   std::vector<std::vector<Tensor<1,dim, NumberType>>> grad_Nx_p_fluid;
2863 *  
2864 *   ScratchData_ASM(const FiniteElement<dim> &fe_cell,
2865 *   const QGauss<dim> &qf_cell, const UpdateFlags uf_cell,
2866 *   const QGauss<dim - 1> & qf_face, const UpdateFlags uf_face,
2867 *   const TrilinosWrappers::MPI::BlockVector &solution_total )
2868 *   :
2869 *   solution_total (solution_total),
2870 *   fe_values_ref(fe_cell, qf_cell, uf_cell),
2871 *   fe_face_values_ref(fe_cell, qf_face, uf_face),
2872 *   local_dof_values(fe_cell.dofs_per_cell),
2873 *   solution_grads_u_total(qf_cell.size()),
2874 *   solution_values_p_fluid_total(qf_cell.size()),
2875 *   solution_grads_p_fluid_total(qf_cell.size()),
2876 *   solution_grads_face_p_fluid_total(qf_face.size()),
2877 *   Nx(qf_cell.size(), std::vector<Tensor<1,dim>>(fe_cell.dofs_per_cell)),
2878 *   Nx_p_fluid(qf_cell.size(), std::vector<double>(fe_cell.dofs_per_cell)),
2879 *   grad_Nx(qf_cell.size(), std::vector<Tensor<2, dim, NumberType>>(fe_cell.dofs_per_cell)),
2880 *   symm_grad_Nx(qf_cell.size(), std::vector<SymmetricTensor<2, dim, NumberType>> (fe_cell.dofs_per_cell)),
2881 *   grad_Nx_p_fluid(qf_cell.size(), std::vector<Tensor<1, dim, NumberType>>(fe_cell.dofs_per_cell))
2882 *   {}
2883 *  
2884 *   ScratchData_ASM(const ScratchData_ASM &rhs)
2885 *   :
2886 *   solution_total (rhs.solution_total),
2887 *   fe_values_ref(rhs.fe_values_ref.get_fe(),
2888 *   rhs.fe_values_ref.get_quadrature(),
2889 *   rhs.fe_values_ref.get_update_flags()),
2890 *   fe_face_values_ref(rhs.fe_face_values_ref.get_fe(),
2891 *   rhs.fe_face_values_ref.get_quadrature(),
2892 *   rhs.fe_face_values_ref.get_update_flags()),
2893 *   local_dof_values(rhs.local_dof_values),
2894 *   solution_grads_u_total(rhs.solution_grads_u_total),
2895 *   solution_values_p_fluid_total(rhs.solution_values_p_fluid_total),
2896 *   solution_grads_p_fluid_total(rhs.solution_grads_p_fluid_total),
2897 *   solution_grads_face_p_fluid_total(rhs.solution_grads_face_p_fluid_total),
2898 *   Nx(rhs.Nx),
2899 *   Nx_p_fluid(rhs.Nx_p_fluid),
2900 *   grad_Nx(rhs.grad_Nx),
2901 *   symm_grad_Nx(rhs.symm_grad_Nx),
2902 *   grad_Nx_p_fluid(rhs.grad_Nx_p_fluid)
2903 *   {}
2904 *  
2905 *   void reset()
2906 *   {
2907 *   const unsigned int n_q_points = Nx_p_fluid.size();
2908 *   const unsigned int n_dofs_per_cell = Nx_p_fluid[0].size();
2909 *  
2910 *   Assert(local_dof_values.size() == n_dofs_per_cell, ExcInternalError());
2911 *  
2912 *   for (unsigned int k = 0; k < n_dofs_per_cell; ++k)
2913 *   {
2914 *   local_dof_values[k] = 0.0;
2915 *   }
2916 *  
2917 *   Assert(solution_grads_u_total.size() == n_q_points, ExcInternalError());
2918 *   Assert(solution_values_p_fluid_total.size() == n_q_points, ExcInternalError());
2919 *   Assert(solution_grads_p_fluid_total.size() == n_q_points, ExcInternalError());
2920 *  
2921 *   Assert(Nx.size() == n_q_points, ExcInternalError());
2922 *   Assert(grad_Nx.size() == n_q_points, ExcInternalError());
2923 *   Assert(symm_grad_Nx.size() == n_q_points, ExcInternalError());
2924 *  
2925 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2926 *   {
2927 *   Assert( Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2928 *   Assert( grad_Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2929 *   Assert( symm_grad_Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2930 *  
2931 *   solution_grads_u_total[q_point] = 0.0;
2932 *   solution_values_p_fluid_total[q_point] = 0.0;
2933 *   solution_grads_p_fluid_total[q_point] = 0.0;
2934 *  
2935 *   for (unsigned int k = 0; k < n_dofs_per_cell; ++k)
2936 *   {
2937 *   Nx[q_point][k] = 0.0;
2938 *   Nx_p_fluid[q_point][k] = 0.0;
2939 *   grad_Nx[q_point][k] = 0.0;
2940 *   symm_grad_Nx[q_point][k] = 0.0;
2941 *   grad_Nx_p_fluid[q_point][k] = 0.0;
2942 *   }
2943 *   }
2944 *  
2945 *   const unsigned int n_f_q_points = solution_grads_face_p_fluid_total.size();
2946 *   Assert(solution_grads_face_p_fluid_total.size() == n_f_q_points, ExcInternalError());
2947 *  
2948 *   for (unsigned int f_q_point = 0; f_q_point < n_f_q_points; ++f_q_point)
2949 *   solution_grads_face_p_fluid_total[f_q_point] = 0.0;
2950 *   }
2951 *   };
2952 *  
2953 * @endcode
2954 *
2955 * Define the boundary conditions on the mesh
2956 *
2957 * @code
2958 *   template <int dim>
2959 *   void Solid<dim>::make_constraints(const int &it_nr_IN)
2960 *   {
2961 *   pcout << " CST " << std::flush;
2962 *   outfile << " CST " << std::flush;
2963 *  
2964 *   if (it_nr_IN > 1) return;
2965 *  
2966 *   const bool apply_dirichlet_bc = (it_nr_IN == 0);
2967 *  
2968 *   if (apply_dirichlet_bc)
2969 *   {
2970 *   constraints.clear();
2971 *   #if DEAL_II_VERSION_GTE(9, 6, 0)
2972 *   constraints.reinit(locally_owned_dofs, locally_relevant_dofs);
2973 *   #else
2974 *   constraints.reinit(locally_relevant_dofs);
2975 *   #endif
2976 *   make_dirichlet_constraints(constraints);
2977 *   }
2978 *   else
2979 *   {
2980 *   for (const auto i : locally_relevant_dofs)
2981 *   if (constraints.is_inhomogeneously_constrained(i) == true)
2982 *   constraints.set_inhomogeneity(i,0.0);
2983 *   }
2984 *   constraints.close();
2985 *   }
2986 *  
2987 * @endcode
2988 *
2989 * Set-up the FE system
2990 *
2991 * @code
2992 *   template <int dim>
2993 *   void Solid<dim>::system_setup(TrilinosWrappers::MPI::BlockVector &solution_delta_OUT)
2994 *   {
2995 *   timerconsole.enter_subsection("Setup system");
2996 *   timerfile.enter_subsection("Setup system");
2997 *  
2998 * @endcode
2999 *
3000 * Determine number of components per block
3001 *
3002 * @code
3003 *   std::vector<unsigned int> block_component(n_components, u_block);
3004 *   block_component[p_fluid_component] = p_fluid_block;
3005 *  
3006 * @endcode
3007 *
3008 * The DOF handler is initialised and we renumber the grid in an efficient manner.
3009 *
3010 * @code
3011 *   dof_handler_ref.distribute_dofs(fe);
3012 *   DoFRenumbering::Cuthill_McKee(dof_handler_ref);
3013 *   DoFRenumbering::component_wise(dof_handler_ref, block_component);
3014 *  
3015 * @endcode
3016 *
3017 * Count the number of DoFs in each block
3018 *
3019 * @code
3020 *   dofs_per_block = DoFTools::count_dofs_per_fe_block(dof_handler_ref, block_component);
3021 *  
3022 * @endcode
3023 *
3024 * Setup the sparsity pattern and tangent matrix
3025 *
3026 * @code
3027 *   all_locally_owned_dofs = DoFTools::locally_owned_dofs_per_subdomain (dof_handler_ref);
3028 *   std::vector<IndexSet> all_locally_relevant_dofs
3029 *   = DoFTools::locally_relevant_dofs_per_subdomain (dof_handler_ref);
3030 *  
3031 *   locally_owned_dofs.clear();
3032 *   locally_owned_partitioning.clear();
3033 *   Assert(all_locally_owned_dofs.size() > this_mpi_process, ExcInternalError());
3034 *   locally_owned_dofs = all_locally_owned_dofs[this_mpi_process];
3035 *  
3036 *   locally_relevant_dofs.clear();
3037 *   locally_relevant_partitioning.clear();
3038 *   Assert(all_locally_relevant_dofs.size() > this_mpi_process, ExcInternalError());
3039 *   locally_relevant_dofs = all_locally_relevant_dofs[this_mpi_process];
3040 *  
3041 *   locally_owned_partitioning.reserve(n_blocks);
3042 *   locally_relevant_partitioning.reserve(n_blocks);
3043 *  
3044 *   for (unsigned int b=0; b<n_blocks; ++b)
3045 *   {
3046 *   const types::global_dof_index idx_begin
3047 *   = std::accumulate(dofs_per_block.begin(),
3048 *   std::next(dofs_per_block.begin(),b), 0);
3049 *   const types::global_dof_index idx_end
3050 *   = std::accumulate(dofs_per_block.begin(),
3051 *   std::next(dofs_per_block.begin(),b+1), 0);
3052 *   locally_owned_partitioning.push_back(locally_owned_dofs.get_view(idx_begin, idx_end));
3053 *   locally_relevant_partitioning.push_back(locally_relevant_dofs.get_view(idx_begin, idx_end));
3054 *   }
3055 *  
3056 * @endcode
3057 *
3058 * Print information on screen
3059 *
3060 * @code
3061 *   pcout << "\nTriangulation:\n"
3062 *   << " Number of active cells: "
3063 *   << triangulation.n_active_cells()
3064 *   << " (by partition:";
3065 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
3066 *   pcout << (p==0 ? ' ' : '+')
3067 *   << (GridTools::count_cells_with_subdomain_association (triangulation,p));
3068 *   pcout << ")"
3069 *   << std::endl;
3070 *   pcout << " Number of degrees of freedom: "
3071 *   << dof_handler_ref.n_dofs()
3072 *   << " (by partition:";
3073 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
3074 *   pcout << (p==0 ? ' ' : '+')
3075 *   << (DoFTools::count_dofs_with_subdomain_association (dof_handler_ref,p));
3076 *   pcout << ")"
3077 *   << std::endl;
3078 *   pcout << " Number of degrees of freedom per block: "
3079 *   << "[n_u, n_p_fluid] = ["
3080 *   << dofs_per_block[u_block]
3081 *   << ", "
3082 *   << dofs_per_block[p_fluid_block]
3083 *   << "]"
3084 *   << std::endl;
3085 *  
3086 * @endcode
3087 *
3088 * Print information to file
3089 *
3090 * @code
3091 *   outfile << "\nTriangulation:\n"
3092 *   << " Number of active cells: "
3093 *   << triangulation.n_active_cells()
3094 *   << " (by partition:";
3095 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
3096 *   outfile << (p==0 ? ' ' : '+')
3097 *   << (GridTools::count_cells_with_subdomain_association (triangulation,p));
3098 *   outfile << ")"
3099 *   << std::endl;
3100 *   outfile << " Number of degrees of freedom: "
3101 *   << dof_handler_ref.n_dofs()
3102 *   << " (by partition:";
3103 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
3104 *   outfile << (p==0 ? ' ' : '+')
3105 *   << (DoFTools::count_dofs_with_subdomain_association (dof_handler_ref,p));
3106 *   outfile << ")"
3107 *   << std::endl;
3108 *   outfile << " Number of degrees of freedom per block: "
3109 *   << "[n_u, n_p_fluid] = ["
3110 *   << dofs_per_block[u_block]
3111 *   << ", "
3112 *   << dofs_per_block[p_fluid_block]
3113 *   << "]"
3114 *   << std::endl;
3115 *  
3116 * @endcode
3117 *
3118 * We optimise the sparsity pattern to reflect this structure and prevent
3119 * unnecessary data creation for the right-diagonal block components.
3120 *
3121 * @code
3122 *   Table<2, DoFTools::Coupling> coupling(n_components, n_components);
3123 *   for (unsigned int ii = 0; ii < n_components; ++ii)
3124 *   for (unsigned int jj = 0; jj < n_components; ++jj)
3125 *  
3126 * @endcode
3127 *
3128 * Identify "zero" matrix components of FE-system (The two components do not couple)
3129 *
3130 * @code
3131 *   if (((ii == p_fluid_component) && (jj < p_fluid_component))
3132 *   || ((ii < p_fluid_component) && (jj == p_fluid_component)) )
3133 *   coupling[ii][jj] = DoFTools::none;
3134 *  
3135 * @endcode
3136 *
3137 * The rest of components always couple
3138 *
3139 * @code
3140 *   else
3141 *   coupling[ii][jj] = DoFTools::always;
3142 *  
3143 *   TrilinosWrappers::BlockSparsityPattern bsp (locally_owned_partitioning,
3144 *   mpi_communicator);
3145 *  
3146 *   DoFTools::make_sparsity_pattern (dof_handler_ref, bsp, constraints,
3147 *   false, this_mpi_process);
3148 *   bsp.compress();
3149 *  
3150 * @endcode
3151 *
3152 * Reinitialize the (sparse) tangent matrix with the given sparsity pattern.
3153 *
3154 * @code
3155 *   tangent_matrix.reinit (bsp);
3156 *  
3157 * @endcode
3158 *
3159 * Initialize the right hand side and solution vectors with number of DoFs
3160 *
3161 * @code
3162 *   system_rhs.reinit(locally_owned_partitioning, mpi_communicator);
3163 *   solution_n.reinit(locally_owned_partitioning, mpi_communicator);
3164 *   solution_delta_OUT.reinit(locally_owned_partitioning, mpi_communicator);
3165 *  
3166 * @endcode
3167 *
3168 * Non-block system
3169 *
3170 * @code
3171 *   TrilinosWrappers::SparsityPattern sp (locally_owned_dofs,
3172 *   mpi_communicator);
3173 *   DoFTools::make_sparsity_pattern (dof_handler_ref, sp, constraints,
3174 *   false, this_mpi_process);
3175 *   sp.compress();
3176 *   tangent_matrix_nb.reinit (sp);
3177 *   system_rhs_nb.reinit(locally_owned_dofs, mpi_communicator);
3178 *  
3179 * @endcode
3180 *
3181 * Set up the quadrature point history
3182 *
3183 * @code
3184 *   setup_qph();
3185 *  
3186 *   timerconsole.leave_subsection();
3187 *   timerfile.leave_subsection();
3188 *   }
3189 *  
3190 * @endcode
3191 *
3192 * Component extractors: used to extract sub-blocks from the global matrix
3193 * Description of which local element DOFs are attached to which block component
3194 *
3195 * @code
3196 *   template <int dim>
3197 *   void Solid<dim>::determine_component_extractors()
3198 *   {
3199 *   element_indices_u.clear();
3200 *   element_indices_p_fluid.clear();
3201 *  
3202 *   for (unsigned int k = 0; k < fe.dofs_per_cell; ++k)
3203 *   {
3204 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
3205 *   if (k_group == u_block)
3206 *   element_indices_u.push_back(k);
3207 *   else if (k_group == p_fluid_block)
3208 *   element_indices_p_fluid.push_back(k);
3209 *   else
3210 *   {
3211 *   Assert(k_group <= p_fluid_block, ExcInternalError());
3212 *   }
3213 *   }
3214 *   }
3215 *  
3216 * @endcode
3217 *
3218 * Set-up quadrature point history (QPH) data objects
3219 *
3220 * @code
3221 *   template <int dim>
3222 *   void Solid<dim>::setup_qph()
3223 *   {
3224 *   pcout << "\nSetting up quadrature point data..." << std::endl;
3225 *   outfile << "\nSetting up quadrature point data..." << std::endl;
3226 *  
3227 * @endcode
3228 *
3229 * Create QPH data objects.
3230 *
3231 * @code
3232 *   quadrature_point_history.initialize(triangulation.begin_active(),
3233 *   triangulation.end(), n_q_points);
3234 *  
3235 * @endcode
3236 *
3237 * Setup the initial quadrature point data using the info stored in parameters
3238 *
3239 * @code
3242 *   dof_handler_ref.begin_active()),
3244 *   dof_handler_ref.end());
3245 *   for (; cell!=endc; ++cell)
3246 *   {
3247 *   Assert(cell->is_locally_owned(), ExcInternalError());
3248 *   Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
3249 *  
3250 *   const std::vector<std::shared_ptr<PointHistory<dim, ADNumberType> > >
3251 *   lqph = quadrature_point_history.get_data(cell);
3252 *   Assert(lqph.size() == n_q_points, ExcInternalError());
3253 *  
3254 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3255 *   lqph[q_point]->setup_lqp(parameters, time);
3256 *   }
3257 *   }
3258 *  
3259 * @endcode
3260 *
3261 * Solve the non-linear system using a Newton-Raphson scheme
3262 *
3263 * @code
3264 *   template <int dim>
3265 *   void Solid<dim>::solve_nonlinear_timestep(TrilinosWrappers::MPI::BlockVector &solution_delta_OUT)
3266 *   {
3267 * @endcode
3268 *
3269 * Print the load step
3270 *
3271 * @code
3272 *   pcout << std::endl
3273 *   << "\nTimestep "
3274 *   << time.get_timestep()
3275 *   << " @ "
3276 *   << time.get_current()
3277 *   << "s"
3278 *   << std::endl;
3279 *   outfile << std::endl
3280 *   << "\nTimestep "
3281 *   << time.get_timestep()
3282 *   << " @ "
3283 *   << time.get_current()
3284 *   << "s"
3285 *   << std::endl;
3286 *  
3287 * @endcode
3288 *
3289 * Declare newton_update vector (solution of a Newton iteration),
3290 * which must have as many positions as global DoFs.
3291 *
3292 * @code
3293 *   TrilinosWrappers::MPI::BlockVector newton_update
3294 *   (locally_owned_partitioning, mpi_communicator);
3295 *  
3296 * @endcode
3297 *
3298 * Reset the error storage objects
3299 *
3300 * @code
3301 *   error_residual.reset();
3302 *   error_residual_0.reset();
3303 *   error_residual_norm.reset();
3304 *   error_update.reset();
3305 *   error_update_0.reset();
3306 *   error_update_norm.reset();
3307 *  
3308 *   print_conv_header();
3309 *  
3310 * @endcode
3311 *
3312 * Declare and initialize iterator for the Newton-Raphson algorithm steps
3313 *
3314 * @code
3315 *   unsigned int newton_iteration = 0;
3316 *  
3317 * @endcode
3318 *
3319 * Iterate until error is below tolerance or max number iterations are reached
3320 *
3321 * @code
3322 *   while(newton_iteration < parameters.max_iterations_NR)
3323 *   {
3324 *   pcout << " " << std::setw(2) << newton_iteration << " " << std::flush;
3325 *   outfile << " " << std::setw(2) << newton_iteration << " " << std::flush;
3326 *  
3327 * @endcode
3328 *
3329 * Initialize global stiffness matrix and global force vector to zero
3330 *
3331 * @code
3332 *   tangent_matrix = 0.0;
3333 *   system_rhs = 0.0;
3334 *  
3335 *   tangent_matrix_nb = 0.0;
3336 *   system_rhs_nb = 0.0;
3337 *  
3338 * @endcode
3339 *
3340 * Apply boundary conditions
3341 *
3342 * @code
3343 *   make_constraints(newton_iteration);
3344 *   assemble_system(solution_delta_OUT);
3345 *  
3346 * @endcode
3347 *
3348 * Compute the rhs residual (error between external and internal forces in FE system)
3349 *
3350 * @code
3351 *   get_error_residual(error_residual);
3352 *  
3353 * @endcode
3354 *
3355 * error_residual in first iteration is stored to normalize posterior error measures
3356 *
3357 * @code
3358 *   if (newton_iteration == 0)
3359 *   error_residual_0 = error_residual;
3360 *  
3361 * @endcode
3362 *
3363 * Determine the normalised residual error
3364 *
3365 * @code
3366 *   error_residual_norm = error_residual;
3367 *   error_residual_norm.normalise(error_residual_0);
3368 *  
3369 * @endcode
3370 *
3371 * If both errors are below the tolerances, exit the loop.
3372 * We need to check the residual vector directly for convergence
3373 * in the load steps where no external forces or displacements are imposed.
3374 *
3375 * @code
3376 *   if ( ((newton_iteration > 0)
3377 *   && (error_update_norm.u <= parameters.tol_u)
3378 *   && (error_update_norm.p_fluid <= parameters.tol_p_fluid)
3379 *   && (error_residual_norm.u <= parameters.tol_f)
3380 *   && (error_residual_norm.p_fluid <= parameters.tol_f))
3381 *   || ( (newton_iteration > 0)
3382 *   && system_rhs.l2_norm() <= parameters.tol_f) )
3383 *   {
3384 *   pcout << "\n ***** CONVERGED! ***** "
3385 *   << system_rhs.l2_norm() << " "
3386 *   << " " << error_residual_norm.norm
3387 *   << " " << error_residual_norm.u
3388 *   << " " << error_residual_norm.p_fluid
3389 *   << " " << error_update_norm.norm
3390 *   << " " << error_update_norm.u
3391 *   << " " << error_update_norm.p_fluid
3392 *   << " " << std::endl;
3393 *   outfile << "\n ***** CONVERGED! ***** "
3394 *   << system_rhs.l2_norm() << " "
3395 *   << " " << error_residual_norm.norm
3396 *   << " " << error_residual_norm.u
3397 *   << " " << error_residual_norm.p_fluid
3398 *   << " " << error_update_norm.norm
3399 *   << " " << error_update_norm.u
3400 *   << " " << error_update_norm.p_fluid
3401 *   << " " << std::endl;
3402 *   print_conv_footer();
3403 *  
3404 *   break;
3405 *   }
3406 *  
3407 * @endcode
3408 *
3409 * Solve the linearized system
3410 *
3411 * @code
3412 *   solve_linear_system(newton_update);
3413 *   constraints.distribute(newton_update);
3414 *  
3415 * @endcode
3416 *
3417 * Compute the displacement error
3418 *
3419 * @code
3420 *   get_error_update(newton_update, error_update);
3421 *  
3422 * @endcode
3423 *
3424 * error_update in first iteration is stored to normalize posterior error measures
3425 *
3426 * @code
3427 *   if (newton_iteration == 0)
3428 *   error_update_0 = error_update;
3429 *  
3430 * @endcode
3431 *
3432 * Determine the normalised Newton update error
3433 *
3434 * @code
3435 *   error_update_norm = error_update;
3436 *   error_update_norm.normalise(error_update_0);
3437 *  
3438 * @endcode
3439 *
3440 * Determine the normalised residual error
3441 *
3442 * @code
3443 *   error_residual_norm = error_residual;
3444 *   error_residual_norm.normalise(error_residual_0);
3445 *  
3446 * @endcode
3447 *
3448 * Print error values
3449 *
3450 * @code
3451 *   pcout << " | " << std::fixed << std::setprecision(3)
3452 *   << std::setw(7) << std::scientific
3453 *   << system_rhs.l2_norm()
3454 *   << " " << error_residual_norm.norm
3455 *   << " " << error_residual_norm.u
3456 *   << " " << error_residual_norm.p_fluid
3457 *   << " " << error_update_norm.norm
3458 *   << " " << error_update_norm.u
3459 *   << " " << error_update_norm.p_fluid
3460 *   << " " << std::endl;
3461 *  
3462 *   outfile << " | " << std::fixed << std::setprecision(3)
3463 *   << std::setw(7) << std::scientific
3464 *   << system_rhs.l2_norm()
3465 *   << " " << error_residual_norm.norm
3466 *   << " " << error_residual_norm.u
3467 *   << " " << error_residual_norm.p_fluid
3468 *   << " " << error_update_norm.norm
3469 *   << " " << error_update_norm.u
3470 *   << " " << error_update_norm.p_fluid
3471 *   << " " << std::endl;
3472 *  
3473 * @endcode
3474 *
3475 * Update
3476 *
3477 * @code
3478 *   solution_delta_OUT += newton_update;
3479 *   newton_update = 0.0;
3480 *   newton_iteration++;
3481 *   }
3482 *  
3483 * @endcode
3484 *
3485 * If maximum allowed number of iterations for Newton algorithm are reached, print non-convergence message and abort program
3486 *
3487 * @code
3488 *   AssertThrow (newton_iteration < parameters.max_iterations_NR, ExcMessage("No convergence in nonlinear solver!"));
3489 *   }
3490 *  
3491 * @endcode
3492 *
3493 * Prints the header for convergence info on console
3494 *
3495 * @code
3496 *   template <int dim>
3497 *   void Solid<dim>::print_conv_header()
3498 *   {
3499 *   static const unsigned int l_width = 120;
3500 *  
3501 *   for (unsigned int i = 0; i < l_width; ++i)
3502 *   {
3503 *   pcout << "_";
3504 *   outfile << "_";
3505 *   }
3506 *  
3507 *   pcout << std::endl;
3508 *   outfile << std::endl;
3509 *  
3510 *   pcout << "\n SOLVER STEP | SYS_RES "
3511 *   << "RES_NORM RES_U RES_P "
3512 *   << "NU_NORM NU_U NU_P " << std::endl;
3513 *   outfile << "\n SOLVER STEP | SYS_RES "
3514 *   << "RES_NORM RES_U RES_P "
3515 *   << "NU_NORM NU_U NU_P " << std::endl;
3516 *  
3517 *   for (unsigned int i = 0; i < l_width; ++i)
3518 *   {
3519 *   pcout << "_";
3520 *   outfile << "_";
3521 *   }
3522 *   pcout << std::endl << std::endl;
3523 *   outfile << std::endl << std::endl;
3524 *   }
3525 *  
3526 * @endcode
3527 *
3528 * Prints the footer for convergence info on console
3529 *
3530 * @code
3531 *   template <int dim>
3532 *   void Solid<dim>::print_conv_footer()
3533 *   {
3534 *   static const unsigned int l_width = 120;
3535 *  
3536 *   for (unsigned int i = 0; i < l_width; ++i)
3537 *   {
3538 *   pcout << "_";
3539 *   outfile << "_";
3540 *   }
3541 *   pcout << std::endl << std::endl;
3542 *   outfile << std::endl << std::endl;
3543 *  
3544 *   pcout << "Relative errors:" << std::endl
3545 *   << "Displacement: "
3546 *   << error_update.u / error_update_0.u << std::endl
3547 *   << "Force (displ): "
3548 *   << error_residual.u / error_residual_0.u << std::endl
3549 *   << "Pore pressure: "
3550 *   << error_update.p_fluid / error_update_0.p_fluid << std::endl
3551 *   << "Force (pore): "
3552 *   << error_residual.p_fluid / error_residual_0.p_fluid << std::endl;
3553 *   outfile << "Relative errors:" << std::endl
3554 *   << "Displacement: "
3555 *   << error_update.u / error_update_0.u << std::endl
3556 *   << "Force (displ): "
3557 *   << error_residual.u / error_residual_0.u << std::endl
3558 *   << "Pore pressure: "
3559 *   << error_update.p_fluid / error_update_0.p_fluid << std::endl
3560 *   << "Force (pore): "
3561 *   << error_residual.p_fluid / error_residual_0.p_fluid << std::endl;
3562 *   }
3563 *  
3564 * @endcode
3565 *
3566 * Determine the true residual error for the problem
3567 *
3568 * @code
3569 *   template <int dim>
3570 *   void Solid<dim>::get_error_residual(Errors &error_residual_OUT)
3571 *   {
3572 *   TrilinosWrappers::MPI::BlockVector error_res(system_rhs);
3573 *   constraints.set_zero(error_res);
3574 *  
3575 *   error_residual_OUT.norm = error_res.l2_norm();
3576 *   error_residual_OUT.u = error_res.block(u_block).l2_norm();
3577 *   error_residual_OUT.p_fluid = error_res.block(p_fluid_block).l2_norm();
3578 *   }
3579 *  
3580 * @endcode
3581 *
3582 * Determine the true Newton update error for the problem
3583 *
3584 * @code
3585 *   template <int dim>
3586 *   void Solid<dim>::get_error_update
3587 *   (const TrilinosWrappers::MPI::BlockVector &newton_update_IN,
3588 *   Errors &error_update_OUT)
3589 *   {
3590 *   TrilinosWrappers::MPI::BlockVector error_ud(newton_update_IN);
3591 *   constraints.set_zero(error_ud);
3592 *  
3593 *   error_update_OUT.norm = error_ud.l2_norm();
3594 *   error_update_OUT.u = error_ud.block(u_block).l2_norm();
3595 *   error_update_OUT.p_fluid = error_ud.block(p_fluid_block).l2_norm();
3596 *   }
3597 *  
3598 * @endcode
3599 *
3600 * Compute the total solution, which is valid at any Newton step. This is required as, to reduce
3601 * computational error, the total solution is only updated at the end of the timestep.
3602 *
3603 * @code
3604 *   template <int dim>
3606 *   Solid<dim>::get_total_solution(const TrilinosWrappers::MPI::BlockVector &solution_delta_IN) const
3607 *   {
3608 * @endcode
3609 *
3610 * Cell interpolation -> Ghosted vector
3611 *
3612 * @code
3614 *   solution_total (locally_owned_partitioning,
3615 *   locally_relevant_partitioning,
3616 *   mpi_communicator,
3617 *   /*vector_writable = */ false);
3618 *   TrilinosWrappers::MPI::BlockVector tmp (solution_total);
3619 *   solution_total = solution_n;
3620 *   tmp = solution_delta_IN;
3621 *   solution_total += tmp;
3622 *   return solution_total;
3623 *   }
3624 *  
3625 * @endcode
3626 *
3627 * Compute elemental stiffness tensor and right-hand side force vector, and assemble into global ones
3628 *
3629 * @code
3630 *   template <int dim>
3631 *   void Solid<dim>::assemble_system( const TrilinosWrappers::MPI::BlockVector &solution_delta )
3632 *   {
3633 *   timerconsole.enter_subsection("Assemble system");
3634 *   timerfile.enter_subsection("Assemble system");
3635 *   pcout << " ASM_SYS " << std::flush;
3636 *   outfile << " ASM_SYS " << std::flush;
3637 *  
3638 *   const TrilinosWrappers::MPI::BlockVector solution_total(get_total_solution(solution_delta));
3639 *  
3640 * @endcode
3641 *
3642 * Info given to FEValues and FEFaceValues constructors, to indicate which data will be needed at each element.
3643 *
3644 * @code
3645 *   const UpdateFlags uf_cell(update_values |
3646 *   update_gradients |
3647 *   update_JxW_values);
3648 *   const UpdateFlags uf_face(update_values |
3649 *   update_gradients |
3652 *   update_JxW_values );
3653 *  
3654 * @endcode
3655 *
3656 * Setup a copy of the data structures required for the process and pass them, along with the
3657 * memory addresses of the assembly functions to the WorkStream object for processing
3658 *
3659 * @code
3660 *   PerTaskData_ASM per_task_data(dofs_per_cell);
3661 *   ScratchData_ASM<ADNumberType> scratch_data(fe, qf_cell, uf_cell,
3662 *   qf_face, uf_face,
3663 *   solution_total);
3664 *  
3667 *   dof_handler_ref.begin_active()),
3669 *   dof_handler_ref.end());
3670 *   for (; cell != endc; ++cell)
3671 *   {
3672 *   Assert(cell->is_locally_owned(), ExcInternalError());
3673 *   Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
3674 *  
3675 *   assemble_system_one_cell(cell, scratch_data, per_task_data);
3676 *   copy_local_to_global_system(per_task_data);
3677 *   }
3678 *   tangent_matrix.compress(VectorOperation::add);
3679 *   system_rhs.compress(VectorOperation::add);
3680 *  
3681 *   tangent_matrix_nb.compress(VectorOperation::add);
3682 *   system_rhs_nb.compress(VectorOperation::add);
3683 *  
3684 *   timerconsole.leave_subsection();
3685 *   timerfile.leave_subsection();
3686 *   }
3687 *  
3688 * @endcode
3689 *
3690 * Add the local elemental contribution to the global stiffness tensor
3691 * We do it twice, for the block and the non-block systems
3692 *
3693 * @code
3694 *   template <int dim>
3695 *   void Solid<dim>::copy_local_to_global_system (const PerTaskData_ASM &data)
3696 *   {
3697 *   constraints.distribute_local_to_global(data.cell_matrix,
3698 *   data.cell_rhs,
3699 *   data.local_dof_indices,
3700 *   tangent_matrix,
3701 *   system_rhs);
3702 *  
3703 *   constraints.distribute_local_to_global(data.cell_matrix,
3704 *   data.cell_rhs,
3705 *   data.local_dof_indices,
3706 *   tangent_matrix_nb,
3707 *   system_rhs_nb);
3708 *   }
3709 *  
3710 * @endcode
3711 *
3712 * Compute stiffness matrix and corresponding rhs for one element
3713 *
3714 * @code
3715 *   template <int dim>
3716 *   void Solid<dim>::assemble_system_one_cell
3717 *   (const typename DoFHandler<dim>::active_cell_iterator &cell,
3718 *   ScratchData_ASM<ADNumberType> &scratch,
3719 *   PerTaskData_ASM &data) const
3720 *   {
3721 *   Assert(cell->is_locally_owned(), ExcInternalError());
3722 *  
3723 *   data.reset();
3724 *   scratch.reset();
3725 *   scratch.fe_values_ref.reinit(cell);
3726 *   cell->get_dof_indices(data.local_dof_indices);
3727 *  
3728 * @endcode
3729 *
3730 * Setup automatic differentiation
3731 *
3732 * @code
3733 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
3734 *   {
3735 * @endcode
3736 *
3737 * Initialise the dofs for the cell using the current solution.
3738 *
3739 * @code
3740 *   scratch.local_dof_values[k] = scratch.solution_total[data.local_dof_indices[k]];
3741 * @endcode
3742 *
3743 * Mark this cell DoF as an independent variable
3744 *
3745 * @code
3746 *   scratch.local_dof_values[k].diff(k, dofs_per_cell);
3747 *   }
3748 *  
3749 * @endcode
3750 *
3751 * Update the quadrature point solution
3752 * Compute the values and gradients of the solution in terms of the AD variables
3753 *
3754 * @code
3755 *   for (unsigned int q = 0; q < n_q_points; ++q)
3756 *   {
3757 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
3758 *   {
3759 *   const unsigned int k_group = fe.system_to_base_index(k).first.first;
3760 *   if (k_group == u_block)
3761 *   {
3762 *   const Tensor<2, dim> Grad_Nx_u =
3763 *   scratch.fe_values_ref[u_fe].gradient(k, q);
3764 *   for (unsigned int dd = 0; dd < dim; ++dd)
3765 *   {
3766 *   for (unsigned int ee = 0; ee < dim; ++ee)
3767 *   {
3768 *   scratch.solution_grads_u_total[q][dd][ee]
3769 *   += scratch.local_dof_values[k] * Grad_Nx_u[dd][ee];
3770 *   }
3771 *   }
3772 *   }
3773 *   else if (k_group == p_fluid_block)
3774 *   {
3775 *   const double Nx_p = scratch.fe_values_ref[p_fluid_fe].value(k, q);
3776 *   const Tensor<1, dim> Grad_Nx_p =
3777 *   scratch.fe_values_ref[p_fluid_fe].gradient(k, q);
3778 *  
3779 *   scratch.solution_values_p_fluid_total[q]
3780 *   += scratch.local_dof_values[k] * Nx_p;
3781 *   for (unsigned int dd = 0; dd < dim; ++dd)
3782 *   {
3783 *   scratch.solution_grads_p_fluid_total[q][dd]
3784 *   += scratch.local_dof_values[k] * Grad_Nx_p[dd];
3785 *   }
3786 *   }
3787 *   else
3788 *   Assert(k_group <= p_fluid_block, ExcInternalError());
3789 *   }
3790 *   }
3791 *  
3792 * @endcode
3793 *
3794 * Set up pointer "lgph" to the PointHistory object of this element
3795 *
3796 * @code
3797 *   const std::vector<std::shared_ptr<const PointHistory<dim, ADNumberType> > >
3798 *   lqph = quadrature_point_history.get_data(cell);
3799 *   Assert(lqph.size() == n_q_points, ExcInternalError());
3800 *  
3801 *  
3802 * @endcode
3803 *
3804 * Precalculate the element shape function values and gradients
3805 *
3806 * @code
3807 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3808 *   {
3809 *   Tensor<2, dim, ADNumberType> F_AD = scratch.solution_grads_u_total[q_point];
3811 *   Assert(determinant(F_AD) > 0, ExcMessage("Invalid deformation map"));
3812 *   const Tensor<2, dim, ADNumberType> F_inv_AD = invert(F_AD);
3813 *  
3814 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3815 *   {
3816 *   const unsigned int i_group = fe.system_to_base_index(i).first.first;
3817 *  
3818 *   if (i_group == u_block)
3819 *   {
3820 *   scratch.Nx[q_point][i] =
3821 *   scratch.fe_values_ref[u_fe].value(i, q_point);
3822 *   scratch.grad_Nx[q_point][i] =
3823 *   scratch.fe_values_ref[u_fe].gradient(i, q_point)*F_inv_AD;
3824 *   scratch.symm_grad_Nx[q_point][i] =
3825 *   symmetrize(scratch.grad_Nx[q_point][i]);
3826 *   }
3827 *   else if (i_group == p_fluid_block)
3828 *   {
3829 *   scratch.Nx_p_fluid[q_point][i] =
3830 *   scratch.fe_values_ref[p_fluid_fe].value(i, q_point);
3831 *   scratch.grad_Nx_p_fluid[q_point][i] =
3832 *   scratch.fe_values_ref[p_fluid_fe].gradient(i, q_point)*F_inv_AD;
3833 *   }
3834 *   else
3835 *   Assert(i_group <= p_fluid_block, ExcInternalError());
3836 *   }
3837 *   }
3838 *  
3839 * @endcode
3840 *
3841 * Assemble the stiffness matrix and rhs vector
3842 *
3843 * @code
3844 *   std::vector<ADNumberType> residual_ad (dofs_per_cell, ADNumberType(0.0));
3845 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3846 *   {
3847 *   Tensor<2, dim, ADNumberType> F_AD = scratch.solution_grads_u_total[q_point];
3849 *   const ADNumberType det_F_AD = determinant(F_AD);
3850 *  
3851 *   Assert(det_F_AD > 0, ExcInternalError());
3852 *   const Tensor<2, dim, ADNumberType> F_inv_AD = invert(F_AD); //inverse of def. gradient tensor
3853 *  
3854 *   const ADNumberType p_fluid = scratch.solution_values_p_fluid_total[q_point];
3855 *  
3856 *   {
3857 *   PointHistory<dim, ADNumberType> *lqph_q_point_nc =
3858 *   const_cast<PointHistory<dim, ADNumberType>*>(lqph[q_point].get());
3859 *   lqph_q_point_nc->update_internal_equilibrium(F_AD);
3860 *   }
3861 *  
3862 * @endcode
3863 *
3864 * Get some info from constitutive model of solid
3865 *
3866 * @code
3867 *   static const SymmetricTensor< 2, dim, double>
3870 *   tau_E = lqph[q_point]->get_tau_E(F_AD);
3871 *   SymmetricTensor<2, dim, ADNumberType> tau_fluid_vol (I);
3872 *   tau_fluid_vol *= -1.0 * p_fluid * det_F_AD;
3873 *  
3874 * @endcode
3875 *
3876 * Get some info from constitutive model of fluid
3877 *
3878 * @code
3879 *   const ADNumberType det_F_aux = lqph[q_point]->get_converged_det_F();
3880 *   const double det_F_converged = Tensor<0,dim,double>(det_F_aux); //Needs to be double, not AD number
3881 *   const Tensor<1, dim, ADNumberType> overall_body_force
3882 *   = lqph[q_point]->get_overall_body_force(F_AD, parameters);
3883 *  
3884 * @endcode
3885 *
3886 * Define some aliases to make the assembly process easier to follow
3887 *
3888 * @code
3889 *   const std::vector<Tensor<1,dim>> &Nu = scratch.Nx[q_point];
3890 *   const std::vector<SymmetricTensor<2, dim, ADNumberType>>
3891 *   &symm_grad_Nu = scratch.symm_grad_Nx[q_point];
3892 *   const std::vector<double> &Np = scratch.Nx_p_fluid[q_point];
3893 *   const std::vector<Tensor<1, dim, ADNumberType> > &grad_Np
3894 *   = scratch.grad_Nx_p_fluid[q_point];
3895 *   const Tensor<1, dim, ADNumberType> grad_p
3896 *   = scratch.solution_grads_p_fluid_total[q_point]*F_inv_AD;
3897 *   const double JxW = scratch.fe_values_ref.JxW(q_point);
3898 *  
3899 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3900 *   {
3901 *   const unsigned int i_group = fe.system_to_base_index(i).first.first;
3902 *  
3903 *   if (i_group == u_block)
3904 *   {
3905 *   residual_ad[i] += symm_grad_Nu[i] * ( tau_E + tau_fluid_vol ) * JxW;
3906 *   residual_ad[i] -= Nu[i] * overall_body_force * JxW;
3907 *   }
3908 *   else if (i_group == p_fluid_block)
3909 *   {
3910 *   const Tensor<1, dim, ADNumberType> seepage_vel_current
3911 *   = lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p);
3912 *   residual_ad[i] += Np[i] * (det_F_AD - det_F_converged) * JxW;
3913 *   residual_ad[i] -= time.get_delta_t() * grad_Np[i]
3914 *   * seepage_vel_current * JxW;
3915 *   }
3916 *   else
3917 *   Assert(i_group <= p_fluid_block, ExcInternalError());
3918 *   }
3919 *   }
3920 *  
3921 * @endcode
3922 *
3923 * Assemble the Neumann contribution (external force contribution).
3924 *
3925 * @code
3926 *   for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face) //Loop over faces in element
3927 *   {
3928 *   if (cell->face(face)->at_boundary() == true)
3929 *   {
3930 *   scratch.fe_face_values_ref.reinit(cell, face);
3931 *  
3932 *   for (unsigned int f_q_point = 0; f_q_point < n_q_points_f; ++f_q_point)
3933 *   {
3934 *   const Tensor<1, dim> &N
3935 *   = scratch.fe_face_values_ref.normal_vector(f_q_point);
3936 *   const Point<dim> &pt
3937 *   = scratch.fe_face_values_ref.quadrature_point(f_q_point);
3938 *   const Tensor<1, dim> traction
3939 *   = get_neumann_traction(cell->face(face)->boundary_id(), pt, N);
3940 *   const double flow
3941 *   = get_prescribed_fluid_flow(cell->face(face)->boundary_id(), pt);
3942 *  
3943 *   if ( (traction.norm() < 1e-12) && (std::abs(flow) < 1e-12) ) continue;
3944 *  
3945 *   const double JxW_f = scratch.fe_face_values_ref.JxW(f_q_point);
3946 *  
3947 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3948 *   {
3949 *   const unsigned int i_group = fe.system_to_base_index(i).first.first;
3950 *  
3951 *   if ((i_group == u_block) && (traction.norm() > 1e-12))
3952 *   {
3953 *   const unsigned int component_i
3954 *   = fe.system_to_component_index(i).first;
3955 *   const double Nu_f
3956 *   = scratch.fe_face_values_ref.shape_value(i, f_q_point);
3957 *   residual_ad[i] -= (Nu_f * traction[component_i]) * JxW_f;
3958 *   }
3959 *   if ((i_group == p_fluid_block) && (std::abs(flow) > 1e-12))
3960 *   {
3961 *   const double Nu_p
3962 *   = scratch.fe_face_values_ref.shape_value(i, f_q_point);
3963 *   residual_ad[i] -= (Nu_p * flow) * JxW_f;
3964 *   }
3965 *   }
3966 *   }
3967 *   }
3968 *   }
3969 *  
3970 * @endcode
3971 *
3972 * Linearise the residual
3973 *
3974 * @code
3975 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
3976 *   {
3977 *   const ADNumberType &R_i = residual_ad[i];
3978 *  
3979 *   data.cell_rhs(i) -= R_i.val();
3980 *   for (unsigned int j=0; j<dofs_per_cell; ++j)
3981 *   data.cell_matrix(i,j) += R_i.fastAccessDx(j);
3982 *   }
3983 *   }
3984 *  
3985 * @endcode
3986 *
3987 * Store the converged values of the internal variables
3988 *
3989 * @code
3990 *   template <int dim>
3991 *   void Solid<dim>::update_end_timestep()
3992 *   {
3995 *   dof_handler_ref.begin_active()),
3997 *   dof_handler_ref.end());
3998 *   for (; cell!=endc; ++cell)
3999 *   {
4000 *   Assert(cell->is_locally_owned(), ExcInternalError());
4001 *   Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
4002 *  
4003 *   const std::vector<std::shared_ptr<PointHistory<dim, ADNumberType> > >
4004 *   lqph = quadrature_point_history.get_data(cell);
4005 *   Assert(lqph.size() == n_q_points, ExcInternalError());
4006 *   for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
4007 *   lqph[q_point]->update_end_timestep();
4008 *   }
4009 *   }
4010 *  
4011 *  
4012 * @endcode
4013 *
4014 * Solve the linearized equations
4015 *
4016 * @code
4017 *   template <int dim>
4018 *   void Solid<dim>::solve_linear_system( TrilinosWrappers::MPI::BlockVector &newton_update_OUT)
4019 *   {
4020 *  
4021 *   timerconsole.enter_subsection("Linear solver");
4022 *   timerfile.enter_subsection("Linear solver");
4023 *   pcout << " SLV " << std::flush;
4024 *   outfile << " SLV " << std::flush;
4025 *  
4026 *   TrilinosWrappers::MPI::Vector newton_update_nb;
4027 *   newton_update_nb.reinit(locally_owned_dofs, mpi_communicator);
4028 *  
4029 *   SolverControl solver_control (tangent_matrix_nb.m(),
4030 *   1.0e-6 * system_rhs_nb.l2_norm());
4031 *   TrilinosWrappers::SolverDirect solver (solver_control);
4032 *   solver.solve(tangent_matrix_nb, newton_update_nb, system_rhs_nb);
4033 *  
4034 * @endcode
4035 *
4036 * Copy the non-block solution back to block system
4037 *
4038 * @code
4039 *   for (unsigned int i=0; i<locally_owned_dofs.n_elements(); ++i)
4040 *   {
4041 *   const types::global_dof_index idx_i
4042 *   = locally_owned_dofs.nth_index_in_set(i);
4043 *   newton_update_OUT(idx_i) = newton_update_nb(idx_i);
4044 *   }
4045 *   newton_update_OUT.compress(VectorOperation::insert);
4046 *  
4047 *   timerconsole.leave_subsection();
4048 *   timerfile.leave_subsection();
4049 *   }
4050 *  
4051 * @endcode
4052 *
4053 * Class to compute gradient of the pressure
4054 *
4055 * @code
4056 *   template <int dim>
4057 *   class GradientPostprocessor : public DataPostprocessorVector<dim>
4058 *   {
4059 *   public:
4060 *   GradientPostprocessor (const unsigned int p_fluid_component)
4061 *   :
4062 *   DataPostprocessorVector<dim> ("grad_p",
4063 *   update_gradients),
4064 *   p_fluid_component (p_fluid_component)
4065 *   {}
4066 *  
4067 *   virtual ~GradientPostprocessor(){}
4068 *  
4069 *   virtual void
4070 *   evaluate_vector_field
4071 *   (const DataPostprocessorInputs::Vector<dim> &input_data,
4072 *   std::vector<Vector<double> > &computed_quantities) const override
4073 *   {
4074 *   AssertDimension (input_data.solution_gradients.size(),
4075 *   computed_quantities.size());
4076 *   for (unsigned int p=0; p<input_data.solution_gradients.size(); ++p)
4077 *   {
4078 *   AssertDimension (computed_quantities[p].size(), dim);
4079 *   for (unsigned int d=0; d<dim; ++d)
4080 *   computed_quantities[p][d]
4081 *   = input_data.solution_gradients[p][p_fluid_component][d];
4082 *   }
4083 *   }
4084 *  
4085 *   private:
4086 *   const unsigned int p_fluid_component;
4087 *   };
4088 *  
4089 *  
4090 * @endcode
4091 *
4092 * Print results to vtu file
4093 *
4094 * @code
4095 *   template <int dim> void Solid<dim>::output_results_to_vtu
4096 *   (const unsigned int timestep,
4097 *   const double current_time,
4098 *   TrilinosWrappers::MPI::BlockVector solution_IN) const
4099 *   {
4100 *   TrilinosWrappers::MPI::BlockVector solution_total(locally_owned_partitioning,
4101 *   locally_relevant_partitioning,
4102 *   mpi_communicator,
4103 *   false);
4104 *   solution_total = solution_IN;
4106 *   material_id.reinit(triangulation.n_active_cells());
4107 *   std::vector<types::subdomain_id> partition_int(triangulation.n_active_cells());
4108 *   GradientPostprocessor<dim> gradient_postprocessor(p_fluid_component);
4109 *  
4110 * @endcode
4111 *
4112 * Declare local variables with number of stress components
4113 * & assign value according to "dim" value
4114 *
4115 * @code
4116 *   unsigned int num_comp_symm_tensor = 6;
4117 *  
4118 * @endcode
4119 *
4120 * Declare local vectors to store values
4121 * OUTPUT AVERAGED ON ELEMENTS -------------------------------------------
4122 *
4123 * @code
4124 *   std::vector<Vector<double>>cauchy_stresses_total_elements
4125 *   (num_comp_symm_tensor,
4126 *   Vector<double> (triangulation.n_active_cells()));
4127 *   std::vector<Vector<double>>cauchy_stresses_E_elements
4128 *   (num_comp_symm_tensor,
4129 *   Vector<double> (triangulation.n_active_cells()));
4130 *   std::vector<Vector<double>>stretches_elements
4131 *   (dim,
4132 *   Vector<double> (triangulation.n_active_cells()));
4133 *   std::vector<Vector<double>>seepage_velocity_elements
4134 *   (dim,
4135 *   Vector<double> (triangulation.n_active_cells()));
4136 *   Vector<double> porous_dissipation_elements
4137 *   (triangulation.n_active_cells());
4138 *   Vector<double> viscous_dissipation_elements
4139 *   (triangulation.n_active_cells());
4140 *   Vector<double> solid_vol_fraction_elements
4141 *   (triangulation.n_active_cells());
4142 *  
4143 * @endcode
4144 *
4145 * OUTPUT AVERAGED ON NODES ----------------------------------------------
4146 * We need to create a new FE space with a single dof per node to avoid
4147 * duplication of the output on nodes for our problem with dim+1 dofs.
4148 *
4149 * @code
4150 *   FE_Q<dim> fe_vertex(1);
4151 *   DoFHandler<dim> vertex_handler_ref(triangulation);
4152 *   vertex_handler_ref.distribute_dofs(fe_vertex);
4153 *   AssertThrow(vertex_handler_ref.n_dofs() == triangulation.n_vertices(),
4154 *   ExcDimensionMismatch(vertex_handler_ref.n_dofs(),
4155 *   triangulation.n_vertices()));
4156 *  
4157 *   Vector<double> counter_on_vertices_mpi
4158 *   (vertex_handler_ref.n_dofs());
4159 *   Vector<double> sum_counter_on_vertices
4160 *   (vertex_handler_ref.n_dofs());
4161 *  
4162 *   std::vector<Vector<double>>cauchy_stresses_total_vertex_mpi
4163 *   (num_comp_symm_tensor,
4164 *   Vector<double>(vertex_handler_ref.n_dofs()));
4165 *   std::vector<Vector<double>>sum_cauchy_stresses_total_vertex
4166 *   (num_comp_symm_tensor,
4167 *   Vector<double>(vertex_handler_ref.n_dofs()));
4168 *   std::vector<Vector<double>>cauchy_stresses_E_vertex_mpi
4169 *   (num_comp_symm_tensor,
4170 *   Vector<double>(vertex_handler_ref.n_dofs()));
4171 *   std::vector<Vector<double>>sum_cauchy_stresses_E_vertex
4172 *   (num_comp_symm_tensor,
4173 *   Vector<double>(vertex_handler_ref.n_dofs()));
4174 *   std::vector<Vector<double>>stretches_vertex_mpi
4175 *   (dim,
4176 *   Vector<double>(vertex_handler_ref.n_dofs()));
4177 *   std::vector<Vector<double>>sum_stretches_vertex
4178 *   (dim,
4179 *   Vector<double>(vertex_handler_ref.n_dofs()));
4180 *   Vector<double> porous_dissipation_vertex_mpi(vertex_handler_ref.n_dofs());
4181 *   Vector<double> sum_porous_dissipation_vertex(vertex_handler_ref.n_dofs());
4182 *   Vector<double> viscous_dissipation_vertex_mpi(vertex_handler_ref.n_dofs());
4183 *   Vector<double> sum_viscous_dissipation_vertex(vertex_handler_ref.n_dofs());
4184 *   Vector<double> solid_vol_fraction_vertex_mpi(vertex_handler_ref.n_dofs());
4185 *   Vector<double> sum_solid_vol_fraction_vertex(vertex_handler_ref.n_dofs());
4186 *  
4187 * @endcode
4188 *
4189 * We need to create a new FE space with a dim dof per node to
4190 * be able to output data on nodes in vector form
4191 *
4192 * @code
4193 *   FESystem<dim> fe_vertex_vec(FE_Q<dim>(1),dim);
4194 *   DoFHandler<dim> vertex_vec_handler_ref(triangulation);
4195 *   vertex_vec_handler_ref.distribute_dofs(fe_vertex_vec);
4196 *   AssertThrow(vertex_vec_handler_ref.n_dofs() == (dim*triangulation.n_vertices()),
4197 *   ExcDimensionMismatch(vertex_vec_handler_ref.n_dofs(),
4198 *   (dim*triangulation.n_vertices())));
4199 *  
4200 *   Vector<double> seepage_velocity_vertex_vec_mpi(vertex_vec_handler_ref.n_dofs());
4201 *   Vector<double> sum_seepage_velocity_vertex_vec(vertex_vec_handler_ref.n_dofs());
4202 *   Vector<double> counter_on_vertices_vec_mpi(vertex_vec_handler_ref.n_dofs());
4203 *   Vector<double> sum_counter_on_vertices_vec(vertex_vec_handler_ref.n_dofs());
4204 * @endcode
4205 *
4206 * -----------------------------------------------------------------------
4207 *
4208
4209 *
4210 * Declare and initialize local unit vectors (to construct tensor basis)
4211 *
4212 * @code
4213 *   std::vector<Tensor<1,dim>> basis_vectors (dim, Tensor<1,dim>() );
4214 *   for (unsigned int i=0; i<dim; ++i)
4215 *   basis_vectors[i][i] = 1;
4216 *  
4217 * @endcode
4218 *
4219 * Declare an instance of the material class object
4220 *
4221 * @code
4222 *   if (parameters.mat_type == "Neo-Hooke")
4223 *   NeoHooke<dim,ADNumberType> material(parameters,time);
4224 *   else if (parameters.mat_type == "Ogden")
4225 *   Ogden<dim,ADNumberType> material(parameters,time);
4226 *   else if (parameters.mat_type == "visco-Ogden")
4227 *   visco_Ogden <dim,ADNumberType>material(parameters,time);
4228 *   else
4229 *   Assert (false, ExcMessage("Material type not implemented"));
4230 *  
4231 * @endcode
4232 *
4233 * Define a local instance of FEValues to compute updated values required
4234 * to calculate stresses
4235 *
4236 * @code
4237 *   const UpdateFlags uf_cell(update_values | update_gradients |
4238 *   update_JxW_values);
4239 *   FEValues<dim> fe_values_ref (fe, qf_cell, uf_cell);
4240 *  
4241 * @endcode
4242 *
4243 * Iterate through elements (cells) and Gauss Points
4244 *
4245 * @code
4248 *   dof_handler_ref.begin_active()),
4250 *   dof_handler_ref.end()),
4252 *   vertex_handler_ref.begin_active()),
4253 *   cell_v_vec(IteratorFilters::LocallyOwnedCell(),
4254 *   vertex_vec_handler_ref.begin_active());
4255 * @endcode
4256 *
4257 * start cell loop
4258 *
4259 * @code
4260 *   for (; cell!=endc; ++cell, ++cell_v, ++cell_v_vec)
4261 *   {
4262 *   Assert(cell->is_locally_owned(), ExcInternalError());
4263 *   Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
4264 *  
4265 *   material_id(cell->active_cell_index())=
4266 *   static_cast<int>(cell->material_id());
4267 *  
4268 *   fe_values_ref.reinit(cell);
4269 *  
4270 *   std::vector<Tensor<2,dim>> solution_grads_u(n_q_points);
4271 *   fe_values_ref[u_fe].get_function_gradients(solution_total,
4272 *   solution_grads_u);
4273 *  
4274 *   std::vector<double> solution_values_p_fluid_total(n_q_points);
4275 *   fe_values_ref[p_fluid_fe].get_function_values(solution_total,
4276 *   solution_values_p_fluid_total);
4277 *  
4278 *   std::vector<Tensor<1,dim>> solution_grads_p_fluid_AD (n_q_points);
4279 *   fe_values_ref[p_fluid_fe].get_function_gradients(solution_total,
4280 *   solution_grads_p_fluid_AD);
4281 *  
4282 * @endcode
4283 *
4284 * start gauss point loop
4285 *
4286 * @code
4287 *   for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
4288 *   {
4290 *   F_AD = Physics::Elasticity::Kinematics::F(solution_grads_u[q_point]);
4291 *   ADNumberType det_F_AD = determinant(F_AD);
4292 *   const double det_F = Tensor<0,dim,double>(det_F_AD);
4293 *  
4294 *   const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
4295 *   lqph = quadrature_point_history.get_data(cell);
4296 *   Assert(lqph.size() == n_q_points, ExcInternalError());
4297 *  
4298 *   const double p_fluid = solution_values_p_fluid_total[q_point];
4299 *  
4300 * @endcode
4301 *
4302 * Cauchy stress
4303 *
4304 * @code
4305 *   static const SymmetricTensor<2,dim,double>
4307 *   SymmetricTensor<2,dim> sigma_E;
4308 *   const SymmetricTensor<2,dim,ADNumberType> sigma_E_AD =
4309 *   lqph[q_point]->get_Cauchy_E(F_AD);
4310 *  
4311 *   for (unsigned int i=0; i<dim; ++i)
4312 *   for (unsigned int j=0; j<dim; ++j)
4313 *   sigma_E[i][j] = Tensor<0,dim,double>(sigma_E_AD[i][j]);
4314 *  
4315 *   SymmetricTensor<2,dim> sigma_fluid_vol (I);
4316 *   sigma_fluid_vol *= -p_fluid;
4317 *   const SymmetricTensor<2,dim> sigma = sigma_E + sigma_fluid_vol;
4318 *  
4319 * @endcode
4320 *
4321 * Volumes
4322 *
4323 * @code
4324 *   const double solid_vol_fraction = (parameters.solid_vol_frac)/det_F;
4325 *  
4326 * @endcode
4327 *
4328 * Green-Lagrange strain
4329 *
4330 * @code
4331 *   const Tensor<2,dim> E_strain = 0.5*(transpose(F_AD)*F_AD - I);
4332 *  
4333 * @endcode
4334 *
4335 * Seepage velocity
4336 *
4337 * @code
4338 *   const Tensor<2,dim,ADNumberType> F_inv = invert(F_AD);
4339 *   const Tensor<1,dim,ADNumberType> grad_p_fluid_AD =
4340 *   solution_grads_p_fluid_AD[q_point]*F_inv;
4341 *   const Tensor<1,dim,ADNumberType> seepage_vel_AD =
4342 *   lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p_fluid_AD);
4343 *  
4344 * @endcode
4345 *
4346 * Dissipations
4347 *
4348 * @code
4349 *   const double porous_dissipation =
4350 *   lqph[q_point]->get_porous_dissipation(F_AD, grad_p_fluid_AD);
4351 *   const double viscous_dissipation =
4352 *   lqph[q_point]->get_viscous_dissipation();
4353 *  
4354 * @endcode
4355 *
4356 * OUTPUT AVERAGED ON ELEMENTS -------------------------------------------
4357 * Both average on elements and on nodes is NOT weighted with the
4358 * integration point volume, i.e., we assume equal contribution of each
4359 * integration point to the average. Ideally, it should be weighted,
4360 * but I haven't invested time in getting it to work properly.
4361 *
4362 * @code
4363 *   if (parameters.outtype == "elements")
4364 *   {
4365 *   for (unsigned int j=0; j<dim; ++j)
4366 *   {
4367 *   cauchy_stresses_total_elements[j](cell->active_cell_index())
4368 *   += ((sigma*basis_vectors[j])*basis_vectors[j])/n_q_points;
4369 *   cauchy_stresses_E_elements[j](cell->active_cell_index())
4370 *   += ((sigma_E*basis_vectors[j])*basis_vectors[j])/n_q_points;
4371 *   stretches_elements[j](cell->active_cell_index())
4372 *   += std::sqrt(1.0+2.0*Tensor<0,dim,double>(E_strain[j][j]))
4373 *   /n_q_points;
4374 *   seepage_velocity_elements[j](cell->active_cell_index())
4375 *   += Tensor<0,dim,double>(seepage_vel_AD[j])/n_q_points;
4376 *   }
4377 *  
4378 *   porous_dissipation_elements(cell->active_cell_index())
4379 *   += porous_dissipation/n_q_points;
4380 *   viscous_dissipation_elements(cell->active_cell_index())
4381 *   += viscous_dissipation/n_q_points;
4382 *   solid_vol_fraction_elements(cell->active_cell_index())
4383 *   += solid_vol_fraction/n_q_points;
4384 *  
4385 *   cauchy_stresses_total_elements[3](cell->active_cell_index())
4386 *   += ((sigma*basis_vectors[0])*basis_vectors[1])/n_q_points; //sig_xy
4387 *   cauchy_stresses_total_elements[4](cell->active_cell_index())
4388 *   += ((sigma*basis_vectors[0])*basis_vectors[2])/n_q_points;//sig_xz
4389 *   cauchy_stresses_total_elements[5](cell->active_cell_index())
4390 *   += ((sigma*basis_vectors[1])*basis_vectors[2])/n_q_points;//sig_yz
4391 *  
4392 *   cauchy_stresses_E_elements[3](cell->active_cell_index())
4393 *   += ((sigma_E*basis_vectors[0])* basis_vectors[1])/n_q_points; //sig_xy
4394 *   cauchy_stresses_E_elements[4](cell->active_cell_index())
4395 *   += ((sigma_E*basis_vectors[0])* basis_vectors[2])/n_q_points;//sig_xz
4396 *   cauchy_stresses_E_elements[5](cell->active_cell_index())
4397 *   += ((sigma_E*basis_vectors[1])* basis_vectors[2])/n_q_points;//sig_yz
4398 *  
4399 *   }
4400 * @endcode
4401 *
4402 * OUTPUT AVERAGED ON NODES -------------------------------------------
4403 *
4404 * @code
4405 *   else if (parameters.outtype == "nodes")
4406 *   {
4407 *   for (unsigned int v=0; v<(GeometryInfo<dim>::vertices_per_cell); ++v)
4408 *   {
4409 *   types::global_dof_index local_vertex_indices =
4410 *   cell_v->vertex_dof_index(v, 0);
4411 *   counter_on_vertices_mpi(local_vertex_indices) += 1;
4412 *   for (unsigned int k=0; k<dim; ++k)
4413 *   {
4414 *   cauchy_stresses_total_vertex_mpi[k](local_vertex_indices)
4415 *   += (sigma*basis_vectors[k])*basis_vectors[k];
4416 *   cauchy_stresses_E_vertex_mpi[k](local_vertex_indices)
4417 *   += (sigma_E*basis_vectors[k])*basis_vectors[k];
4418 *   stretches_vertex_mpi[k](local_vertex_indices)
4419 *   += std::sqrt(1.0+2.0*Tensor<0,dim,double>(E_strain[k][k]));
4420 *  
4421 *   types::global_dof_index local_vertex_vec_indices =
4422 *   cell_v_vec->vertex_dof_index(v, k);
4423 *   counter_on_vertices_vec_mpi(local_vertex_vec_indices) += 1;
4424 *   seepage_velocity_vertex_vec_mpi(local_vertex_vec_indices)
4425 *   += Tensor<0,dim,double>(seepage_vel_AD[k]);
4426 *   }
4427 *  
4428 *   porous_dissipation_vertex_mpi(local_vertex_indices)
4429 *   += porous_dissipation;
4430 *   viscous_dissipation_vertex_mpi(local_vertex_indices)
4431 *   += viscous_dissipation;
4432 *   solid_vol_fraction_vertex_mpi(local_vertex_indices)
4433 *   += solid_vol_fraction;
4434 *  
4435 *   cauchy_stresses_total_vertex_mpi[3](local_vertex_indices)
4436 *   += (sigma*basis_vectors[0])*basis_vectors[1]; //sig_xy
4437 *   cauchy_stresses_total_vertex_mpi[4](local_vertex_indices)
4438 *   += (sigma*basis_vectors[0])*basis_vectors[2];//sig_xz
4439 *   cauchy_stresses_total_vertex_mpi[5](local_vertex_indices)
4440 *   += (sigma*basis_vectors[1])*basis_vectors[2]; //sig_yz
4441 *  
4442 *   cauchy_stresses_E_vertex_mpi[3](local_vertex_indices)
4443 *   += (sigma_E*basis_vectors[0])*basis_vectors[1]; //sig_xy
4444 *   cauchy_stresses_E_vertex_mpi[4](local_vertex_indices)
4445 *   += (sigma_E*basis_vectors[0])*basis_vectors[2];//sig_xz
4446 *   cauchy_stresses_E_vertex_mpi[5](local_vertex_indices)
4447 *   += (sigma_E*basis_vectors[1])*basis_vectors[2]; //sig_yz
4448 *   }
4449 *   }
4450 * @endcode
4451 *
4452 * ---------------------------------------------------------------
4453 *
4454 * @code
4455 *   } //end gauss point loop
4456 *   }//end cell loop
4457 *  
4458 * @endcode
4459 *
4460 * Different nodes might have different amount of contributions, e.g.,
4461 * corner nodes have less integration points contributing to the averaged.
4462 * This is why we need a counter and divide at the end, outside the cell loop.
4463 *
4464 * @code
4465 *   if (parameters.outtype == "nodes")
4466 *   {
4467 *   for (unsigned int d=0; d<(vertex_handler_ref.n_dofs()); ++d)
4468 *   {
4469 *   sum_counter_on_vertices[d] =
4470 *   Utilities::MPI::sum(counter_on_vertices_mpi[d],
4471 *   mpi_communicator);
4472 *   sum_porous_dissipation_vertex[d] =
4473 *   Utilities::MPI::sum(porous_dissipation_vertex_mpi[d],
4474 *   mpi_communicator);
4475 *   sum_viscous_dissipation_vertex[d] =
4476 *   Utilities::MPI::sum(viscous_dissipation_vertex_mpi[d],
4477 *   mpi_communicator);
4478 *   sum_solid_vol_fraction_vertex[d] =
4479 *   Utilities::MPI::sum(solid_vol_fraction_vertex_mpi[d],
4480 *   mpi_communicator);
4481 *  
4482 *   for (unsigned int k=0; k<num_comp_symm_tensor; ++k)
4483 *   {
4484 *   sum_cauchy_stresses_total_vertex[k][d] =
4485 *   Utilities::MPI::sum(cauchy_stresses_total_vertex_mpi[k][d],
4486 *   mpi_communicator);
4487 *   sum_cauchy_stresses_E_vertex[k][d] =
4488 *   Utilities::MPI::sum(cauchy_stresses_E_vertex_mpi[k][d],
4489 *   mpi_communicator);
4490 *   }
4491 *   for (unsigned int k=0; k<dim; ++k)
4492 *   {
4493 *   sum_stretches_vertex[k][d] =
4494 *   Utilities::MPI::sum(stretches_vertex_mpi[k][d],
4495 *   mpi_communicator);
4496 *   }
4497 *   }
4498 *  
4499 *   for (unsigned int d=0; d<(vertex_vec_handler_ref.n_dofs()); ++d)
4500 *   {
4501 *   sum_counter_on_vertices_vec[d] =
4502 *   Utilities::MPI::sum(counter_on_vertices_vec_mpi[d],
4503 *   mpi_communicator);
4504 *   sum_seepage_velocity_vertex_vec[d] =
4505 *   Utilities::MPI::sum(seepage_velocity_vertex_vec_mpi[d],
4506 *   mpi_communicator);
4507 *   }
4508 *  
4509 *   for (unsigned int d=0; d<(vertex_handler_ref.n_dofs()); ++d)
4510 *   {
4511 *   if (sum_counter_on_vertices[d]>0)
4512 *   {
4513 *   for (unsigned int i=0; i<num_comp_symm_tensor; ++i)
4514 *   {
4515 *   sum_cauchy_stresses_total_vertex[i][d] /= sum_counter_on_vertices[d];
4516 *   sum_cauchy_stresses_E_vertex[i][d] /= sum_counter_on_vertices[d];
4517 *   }
4518 *   for (unsigned int i=0; i<dim; ++i)
4519 *   {
4520 *   sum_stretches_vertex[i][d] /= sum_counter_on_vertices[d];
4521 *   }
4522 *   sum_porous_dissipation_vertex[d] /= sum_counter_on_vertices[d];
4523 *   sum_viscous_dissipation_vertex[d] /= sum_counter_on_vertices[d];
4524 *   sum_solid_vol_fraction_vertex[d] /= sum_counter_on_vertices[d];
4525 *   }
4526 *   }
4527 *  
4528 *   for (unsigned int d=0; d<(vertex_vec_handler_ref.n_dofs()); ++d)
4529 *   {
4530 *   if (sum_counter_on_vertices_vec[d]>0)
4531 *   {
4532 *   sum_seepage_velocity_vertex_vec[d] /= sum_counter_on_vertices_vec[d];
4533 *   }
4534 *   }
4535 *  
4536 *   }
4537 *  
4538 * @endcode
4539 *
4540 * Add the results to the solution to create the output file for Paraview
4541 *
4542 * @code
4543 *   DataOut<dim> data_out;
4544 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
4545 *   comp_type(dim,
4546 *   DataComponentInterpretation::component_is_part_of_vector);
4547 *   comp_type.push_back(DataComponentInterpretation::component_is_scalar);
4548 *  
4549 *   GridTools::get_subdomain_association(triangulation, partition_int);
4550 *  
4551 *   std::vector<std::string> solution_name(dim, "displacement");
4552 *   solution_name.push_back("pore_pressure");
4553 *  
4554 *   data_out.attach_dof_handler(dof_handler_ref);
4555 *   data_out.add_data_vector(solution_total,
4556 *   solution_name,
4557 *   DataOut<dim>::type_dof_data,
4558 *   comp_type);
4559 *  
4560 *   data_out.add_data_vector(solution_total,
4561 *   gradient_postprocessor);
4562 *  
4563 *   const Vector<double> partitioning(partition_int.begin(),
4564 *   partition_int.end());
4565 *  
4566 *   data_out.add_data_vector(partitioning, "partitioning");
4567 *   data_out.add_data_vector(material_id, "material_id");
4568 *  
4569 * @endcode
4570 *
4571 * Integration point results -----------------------------------------------------------
4572 *
4573 * @code
4574 *   if (parameters.outtype == "elements")
4575 *   {
4576 *   data_out.add_data_vector(cauchy_stresses_total_elements[0], "cauchy_xx");
4577 *   data_out.add_data_vector(cauchy_stresses_total_elements[1], "cauchy_yy");
4578 *   data_out.add_data_vector(cauchy_stresses_total_elements[2], "cauchy_zz");
4579 *   data_out.add_data_vector(cauchy_stresses_total_elements[3], "cauchy_xy");
4580 *   data_out.add_data_vector(cauchy_stresses_total_elements[4], "cauchy_xz");
4581 *   data_out.add_data_vector(cauchy_stresses_total_elements[5], "cauchy_yz");
4582 *  
4583 *   data_out.add_data_vector(cauchy_stresses_E_elements[0], "cauchy_E_xx");
4584 *   data_out.add_data_vector(cauchy_stresses_E_elements[1], "cauchy_E_yy");
4585 *   data_out.add_data_vector(cauchy_stresses_E_elements[2], "cauchy_E_zz");
4586 *   data_out.add_data_vector(cauchy_stresses_E_elements[3], "cauchy_E_xy");
4587 *   data_out.add_data_vector(cauchy_stresses_E_elements[4], "cauchy_E_xz");
4588 *   data_out.add_data_vector(cauchy_stresses_E_elements[5], "cauchy_E_yz");
4589 *  
4590 *   data_out.add_data_vector(stretches_elements[0], "stretch_xx");
4591 *   data_out.add_data_vector(stretches_elements[1], "stretch_yy");
4592 *   data_out.add_data_vector(stretches_elements[2], "stretch_zz");
4593 *  
4594 *   data_out.add_data_vector(seepage_velocity_elements[0], "seepage_vel_x");
4595 *   data_out.add_data_vector(seepage_velocity_elements[1], "seepage_vel_y");
4596 *   data_out.add_data_vector(seepage_velocity_elements[2], "seepage_vel_z");
4597 *  
4598 *   data_out.add_data_vector(porous_dissipation_elements, "dissipation_porous");
4599 *   data_out.add_data_vector(viscous_dissipation_elements, "dissipation_viscous");
4600 *   data_out.add_data_vector(solid_vol_fraction_elements, "solid_vol_fraction");
4601 *   }
4602 *   else if (parameters.outtype == "nodes")
4603 *   {
4604 *   data_out.add_data_vector(vertex_handler_ref,
4605 *   sum_cauchy_stresses_total_vertex[0],
4606 *   "cauchy_xx");
4607 *   data_out.add_data_vector(vertex_handler_ref,
4608 *   sum_cauchy_stresses_total_vertex[1],
4609 *   "cauchy_yy");
4610 *   data_out.add_data_vector(vertex_handler_ref,
4611 *   sum_cauchy_stresses_total_vertex[2],
4612 *   "cauchy_zz");
4613 *   data_out.add_data_vector(vertex_handler_ref,
4614 *   sum_cauchy_stresses_total_vertex[3],
4615 *   "cauchy_xy");
4616 *   data_out.add_data_vector(vertex_handler_ref,
4617 *   sum_cauchy_stresses_total_vertex[4],
4618 *   "cauchy_xz");
4619 *   data_out.add_data_vector(vertex_handler_ref,
4620 *   sum_cauchy_stresses_total_vertex[5],
4621 *   "cauchy_yz");
4622 *  
4623 *   data_out.add_data_vector(vertex_handler_ref,
4624 *   sum_cauchy_stresses_E_vertex[0],
4625 *   "cauchy_E_xx");
4626 *   data_out.add_data_vector(vertex_handler_ref,
4627 *   sum_cauchy_stresses_E_vertex[1],
4628 *   "cauchy_E_yy");
4629 *   data_out.add_data_vector(vertex_handler_ref,
4630 *   sum_cauchy_stresses_E_vertex[2],
4631 *   "cauchy_E_zz");
4632 *   data_out.add_data_vector(vertex_handler_ref,
4633 *   sum_cauchy_stresses_E_vertex[3],
4634 *   "cauchy_E_xy");
4635 *   data_out.add_data_vector(vertex_handler_ref,
4636 *   sum_cauchy_stresses_E_vertex[4],
4637 *   "cauchy_E_xz");
4638 *   data_out.add_data_vector(vertex_handler_ref,
4639 *   sum_cauchy_stresses_E_vertex[5],
4640 *   "cauchy_E_yz");
4641 *  
4642 *   data_out.add_data_vector(vertex_handler_ref,
4643 *   sum_stretches_vertex[0],
4644 *   "stretch_xx");
4645 *   data_out.add_data_vector(vertex_handler_ref,
4646 *   sum_stretches_vertex[1],
4647 *   "stretch_yy");
4648 *   data_out.add_data_vector(vertex_handler_ref,
4649 *   sum_stretches_vertex[2],
4650 *   "stretch_zz");
4651 *  
4652 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
4653 *   comp_type_vec(dim,
4654 *   DataComponentInterpretation::component_is_part_of_vector);
4655 *   std::vector<std::string> solution_name_vec(dim,"seepage_velocity");
4656 *  
4657 *   data_out.add_data_vector(vertex_vec_handler_ref,
4658 *   sum_seepage_velocity_vertex_vec,
4659 *   solution_name_vec,
4660 *   comp_type_vec);
4661 *  
4662 *   data_out.add_data_vector(vertex_handler_ref,
4663 *   sum_porous_dissipation_vertex,
4664 *   "dissipation_porous");
4665 *   data_out.add_data_vector(vertex_handler_ref,
4666 *   sum_viscous_dissipation_vertex,
4667 *   "dissipation_viscous");
4668 *   data_out.add_data_vector(vertex_handler_ref,
4669 *   sum_solid_vol_fraction_vertex,
4670 *   "solid_vol_fraction");
4671 *   }
4672 * @endcode
4673 *
4674 * ---------------------------------------------------------------------
4675 *
4676
4677 *
4678 *
4679 * @code
4680 *   data_out.build_patches(degree_displ);
4681 *  
4682 *   struct Filename
4683 *   {
4684 *   static std::string get_filename_vtu(unsigned int process,
4685 *   unsigned int timestep,
4686 *   const unsigned int n_digits = 5)
4687 *   {
4688 *   std::ostringstream filename_vtu;
4689 *   filename_vtu
4690 *   << "solution."
4691 *   << Utilities::int_to_string(process, n_digits)
4692 *   << "."
4693 *   << Utilities::int_to_string(timestep, n_digits)
4694 *   << ".vtu";
4695 *   return filename_vtu.str();
4696 *   }
4697 *  
4698 *   static std::string get_filename_pvtu(unsigned int timestep,
4699 *   const unsigned int n_digits = 5)
4700 *   {
4701 *   std::ostringstream filename_vtu;
4702 *   filename_vtu
4703 *   << "solution."
4704 *   << Utilities::int_to_string(timestep, n_digits)
4705 *   << ".pvtu";
4706 *   return filename_vtu.str();
4707 *   }
4708 *  
4709 *   static std::string get_filename_pvd (void)
4710 *   {
4711 *   std::ostringstream filename_vtu;
4712 *   filename_vtu
4713 *   << "solution.pvd";
4714 *   return filename_vtu.str();
4715 *   }
4716 *   };
4717 *  
4718 *   const std::string filename_vtu = Filename::get_filename_vtu(this_mpi_process,
4719 *   timestep);
4720 *   std::ofstream output(filename_vtu.c_str());
4721 *   data_out.write_vtu(output);
4722 *  
4723 * @endcode
4724 *
4725 * We have a collection of files written in parallel
4726 * This next set of steps should only be performed by master process
4727 *
4728 * @code
4729 *   if (this_mpi_process == 0)
4730 *   {
4731 * @endcode
4732 *
4733 * List of all files written out at this timestep by all processors
4734 *
4735 * @code
4736 *   std::vector<std::string> parallel_filenames_vtu;
4737 *   for (unsigned int p=0; p<n_mpi_processes; ++p)
4738 *   {
4739 *   parallel_filenames_vtu.push_back(Filename::get_filename_vtu(p, timestep));
4740 *   }
4741 *  
4742 *   const std::string filename_pvtu(Filename::get_filename_pvtu(timestep));
4743 *   std::ofstream pvtu_master(filename_pvtu.c_str());
4744 *   data_out.write_pvtu_record(pvtu_master,
4745 *   parallel_filenames_vtu);
4746 *  
4747 * @endcode
4748 *
4749 * Time dependent data master file
4750 *
4751 * @code
4752 *   static std::vector<std::pair<double,std::string>> time_and_name_history;
4753 *   time_and_name_history.push_back(std::make_pair(current_time,
4754 *   filename_pvtu));
4755 *   const std::string filename_pvd(Filename::get_filename_pvd());
4756 *   std::ofstream pvd_output(filename_pvd.c_str());
4757 *   DataOutBase::write_pvd_record(pvd_output, time_and_name_history);
4758 *   }
4759 *   }
4760 *  
4761 *  
4762 * @endcode
4763 *
4764 * Print results to plotting file
4765 *
4766 * @code
4767 *   template <int dim>
4768 *   void Solid<dim>::output_results_to_plot(
4769 *   const unsigned int timestep,
4770 *   const double current_time,
4771 *   TrilinosWrappers::MPI::BlockVector solution_IN,
4772 *   std::vector<Point<dim> > &tracked_vertices_IN,
4773 *   std::ofstream &plotpointfile) const
4774 *   {
4775 *   TrilinosWrappers::MPI::BlockVector solution_total(locally_owned_partitioning,
4776 *   locally_relevant_partitioning,
4777 *   mpi_communicator,
4778 *   false);
4779 *  
4780 *   (void) timestep;
4781 *   solution_total = solution_IN;
4782 *  
4783 * @endcode
4784 *
4785 * Variables needed to print the solution file for plotting
4786 *
4787 * @code
4788 *   Point<dim> reaction_force;
4789 *   Point<dim> reaction_force_pressure;
4790 *   Point<dim> reaction_force_extra;
4791 *   double total_fluid_flow = 0.0;
4792 *   double total_porous_dissipation = 0.0;
4793 *   double total_viscous_dissipation = 0.0;
4794 *   double total_solid_vol = 0.0;
4795 *   double total_vol_current = 0.0;
4796 *   double total_vol_reference = 0.0;
4797 *   std::vector<Point<dim+1>> solution_vertices(tracked_vertices_IN.size());
4798 *  
4799 * @endcode
4800 *
4801 * Auxiliary variables needed for mpi processing
4802 *
4803 * @code
4804 *   Tensor<1,dim> sum_reaction_mpi;
4805 *   Tensor<1,dim> sum_reaction_pressure_mpi;
4806 *   Tensor<1,dim> sum_reaction_extra_mpi;
4807 *   sum_reaction_mpi = 0.0;
4808 *   sum_reaction_pressure_mpi = 0.0;
4809 *   sum_reaction_extra_mpi = 0.0;
4810 *   double sum_total_flow_mpi = 0.0;
4811 *   double sum_porous_dissipation_mpi = 0.0;
4812 *   double sum_viscous_dissipation_mpi = 0.0;
4813 *   double sum_solid_vol_mpi = 0.0;
4814 *   double sum_vol_current_mpi = 0.0;
4815 *   double sum_vol_reference_mpi = 0.0;
4816 *  
4817 * @endcode
4818 *
4819 * Declare an instance of the material class object
4820 *
4821 * @code
4822 *   if (parameters.mat_type == "Neo-Hooke")
4823 *   NeoHooke<dim,ADNumberType> material(parameters,time);
4824 *   else if (parameters.mat_type == "Ogden")
4825 *   Ogden<dim,ADNumberType> material(parameters, time);
4826 *   else if (parameters.mat_type == "visco-Ogden")
4827 *   visco_Ogden <dim,ADNumberType>material(parameters,time);
4828 *   else
4829 *   Assert (false, ExcMessage("Material type not implemented"));
4830 *  
4831 * @endcode
4832 *
4833 * Define a local instance of FEValues to compute updated values required
4834 * to calculate stresses
4835 *
4836 * @code
4837 *   const UpdateFlags uf_cell(update_values | update_gradients |
4838 *   update_JxW_values);
4839 *   FEValues<dim> fe_values_ref (fe, qf_cell, uf_cell);
4840 *  
4841 * @endcode
4842 *
4843 * Iterate through elements (cells) and Gauss Points
4844 *
4845 * @code
4846 *   FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
4847 *   cell(IteratorFilters::LocallyOwnedCell(),
4848 *   dof_handler_ref.begin_active()),
4849 *   endc(IteratorFilters::LocallyOwnedCell(),
4850 *   dof_handler_ref.end());
4851 * @endcode
4852 *
4853 * start cell loop
4854 *
4855 * @code
4856 *   for (; cell!=endc; ++cell)
4857 *   {
4858 *   Assert(cell->is_locally_owned(), ExcInternalError());
4859 *   Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
4860 *  
4861 *   fe_values_ref.reinit(cell);
4862 *  
4863 *   std::vector<Tensor<2,dim>> solution_grads_u(n_q_points);
4864 *   fe_values_ref[u_fe].get_function_gradients(solution_total,
4865 *   solution_grads_u);
4866 *  
4867 *   std::vector<double> solution_values_p_fluid_total(n_q_points);
4868 *   fe_values_ref[p_fluid_fe].get_function_values(solution_total,
4869 *   solution_values_p_fluid_total);
4870 *  
4871 *   std::vector<Tensor<1,dim >> solution_grads_p_fluid_AD(n_q_points);
4872 *   fe_values_ref[p_fluid_fe].get_function_gradients(solution_total,
4873 *   solution_grads_p_fluid_AD);
4874 *  
4875 * @endcode
4876 *
4877 * start gauss point loop
4878 *
4879 * @code
4880 *   for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
4881 *   {
4882 *   const Tensor<2,dim,ADNumberType>
4883 *   F_AD = Physics::Elasticity::Kinematics::F(solution_grads_u[q_point]);
4884 *   ADNumberType det_F_AD = determinant(F_AD);
4885 *   const double det_F = Tensor<0,dim,double>(det_F_AD);
4886 *  
4887 *   const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
4888 *   lqph = quadrature_point_history.get_data(cell);
4889 *   Assert(lqph.size() == n_q_points, ExcInternalError());
4890 *  
4891 *   double JxW = fe_values_ref.JxW(q_point);
4892 *  
4893 * @endcode
4894 *
4895 * Volumes
4896 *
4897 * @code
4898 *   sum_vol_current_mpi += det_F * JxW;
4899 *   sum_vol_reference_mpi += JxW;
4900 *   sum_solid_vol_mpi += parameters.solid_vol_frac * JxW * det_F;
4901 *  
4902 * @endcode
4903 *
4904 * Seepage velocity
4905 *
4906 * @code
4907 *   const Tensor<2,dim,ADNumberType> F_inv = invert(F_AD);
4908 *   const Tensor<1,dim,ADNumberType>
4909 *   grad_p_fluid_AD = solution_grads_p_fluid_AD[q_point]*F_inv;
4910 *   const Tensor<1,dim,ADNumberType> seepage_vel_AD
4911 *   = lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p_fluid_AD);
4912 *  
4913 * @endcode
4914 *
4915 * Dissipations
4916 *
4917 * @code
4918 *   const double porous_dissipation =
4919 *   lqph[q_point]->get_porous_dissipation(F_AD, grad_p_fluid_AD);
4920 *   sum_porous_dissipation_mpi += porous_dissipation * det_F * JxW;
4921 *  
4922 *   const double viscous_dissipation = lqph[q_point]->get_viscous_dissipation();
4923 *   sum_viscous_dissipation_mpi += viscous_dissipation * det_F * JxW;
4924 *  
4925 * @endcode
4926 *
4927 * ---------------------------------------------------------------
4928 *
4929 * @code
4930 *   } //end gauss point loop
4931 *  
4932 * @endcode
4933 *
4934 * Compute reaction force on load boundary & total fluid flow across
4935 * drained boundary.
4936 * Define a local instance of FEFaceValues to compute values required
4937 * to calculate reaction force
4938 *
4939 * @code
4940 *   const UpdateFlags uf_face( update_values | update_gradients |
4941 *   update_normal_vectors | update_JxW_values );
4942 *   FEFaceValues<dim> fe_face_values_ref(fe, qf_face, uf_face);
4943 *  
4944 * @endcode
4945 *
4946 * start face loop
4947 *
4948 * @code
4949 *   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
4950 *   {
4951 * @endcode
4952 *
4953 * Reaction force
4954 *
4955 * @code
4956 *   if (cell->face(face)->at_boundary() == true &&
4957 *   cell->face(face)->boundary_id() == get_reaction_boundary_id_for_output() )
4958 *   {
4959 *   fe_face_values_ref.reinit(cell, face);
4960 *  
4961 * @endcode
4962 *
4963 * Get displacement gradients for current face
4964 *
4965 * @code
4966 *   std::vector<Tensor<2,dim> > solution_grads_u_f(n_q_points_f);
4967 *   fe_face_values_ref[u_fe].get_function_gradients
4968 *   (solution_total,
4969 *   solution_grads_u_f);
4970 *  
4971 * @endcode
4972 *
4973 * Get pressure for current element
4974 *
4975 * @code
4976 *   std::vector< double > solution_values_p_fluid_total_f(n_q_points_f);
4977 *   fe_face_values_ref[p_fluid_fe].get_function_values
4978 *   (solution_total,
4979 *   solution_values_p_fluid_total_f);
4980 *  
4981 * @endcode
4982 *
4983 * start gauss points on faces loop
4984 *
4985 * @code
4986 *   for (unsigned int f_q_point=0; f_q_point<n_q_points_f; ++f_q_point)
4987 *   {
4988 *   const Tensor<1,dim> &N = fe_face_values_ref.normal_vector(f_q_point);
4989 *   const double JxW_f = fe_face_values_ref.JxW(f_q_point);
4990 *  
4991 * @endcode
4992 *
4993 * Compute deformation gradient from displacements gradient
4994 * (present configuration)
4995 *
4996 * @code
4997 *   const Tensor<2,dim,ADNumberType> F_AD =
4998 *   Physics::Elasticity::Kinematics::F(solution_grads_u_f[f_q_point]);
4999 *  
5000 *   const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
5001 *   lqph = quadrature_point_history.get_data(cell);
5002 *   Assert(lqph.size() == n_q_points, ExcInternalError());
5003 *  
5004 *   const double p_fluid = solution_values_p_fluid_total[f_q_point];
5005 *  
5006 * @endcode
5007 *
5008 * Cauchy stress
5009 *
5010 * @code
5011 *   static const SymmetricTensor<2,dim,double>
5012 *   I (Physics::Elasticity::StandardTensors<dim>::I);
5013 *   SymmetricTensor<2,dim> sigma_E;
5014 *   const SymmetricTensor<2,dim,ADNumberType> sigma_E_AD =
5015 *   lqph[f_q_point]->get_Cauchy_E(F_AD);
5016 *  
5017 *   for (unsigned int i=0; i<dim; ++i)
5018 *   for (unsigned int j=0; j<dim; ++j)
5019 *   sigma_E[i][j] = Tensor<0,dim,double>(sigma_E_AD[i][j]);
5020 *  
5021 *   SymmetricTensor<2,dim> sigma_fluid_vol(I);
5022 *   sigma_fluid_vol *= -1.0*p_fluid;
5023 *   const SymmetricTensor<2,dim> sigma = sigma_E+sigma_fluid_vol;
5024 *   sum_reaction_mpi += sigma * N * JxW_f;
5025 *   sum_reaction_pressure_mpi += sigma_fluid_vol * N * JxW_f;
5026 *   sum_reaction_extra_mpi += sigma_E * N * JxW_f;
5027 *   }//end gauss points on faces loop
5028 *   }
5029 *  
5030 * @endcode
5031 *
5032 * Fluid flow
5033 *
5034 * @code
5035 *   if (cell->face(face)->at_boundary() == true &&
5036 *   (cell->face(face)->boundary_id() ==
5037 *   get_drained_boundary_id_for_output().first ||
5038 *   cell->face(face)->boundary_id() ==
5039 *   get_drained_boundary_id_for_output().second ) )
5040 *   {
5041 *   fe_face_values_ref.reinit(cell, face);
5042 *  
5043 * @endcode
5044 *
5045 * Get displacement gradients for current face
5046 *
5047 * @code
5048 *   std::vector<Tensor<2,dim>> solution_grads_u_f(n_q_points_f);
5049 *   fe_face_values_ref[u_fe].get_function_gradients
5050 *   (solution_total,
5051 *   solution_grads_u_f);
5052 *  
5053 * @endcode
5054 *
5055 * Get pressure gradients for current face
5056 *
5057 * @code
5058 *   std::vector<Tensor<1,dim>> solution_grads_p_f(n_q_points_f);
5059 *   fe_face_values_ref[p_fluid_fe].get_function_gradients
5060 *   (solution_total,
5061 *   solution_grads_p_f);
5062 *  
5063 * @endcode
5064 *
5065 * start gauss points on faces loop
5066 *
5067 * @code
5068 *   for (unsigned int f_q_point=0; f_q_point<n_q_points_f; ++f_q_point)
5069 *   {
5070 *   const Tensor<1,dim> &N =
5071 *   fe_face_values_ref.normal_vector(f_q_point);
5072 *   const double JxW_f = fe_face_values_ref.JxW(f_q_point);
5073 *  
5074 * @endcode
5075 *
5076 * Deformation gradient and inverse from displacements gradient
5077 * (present configuration)
5078 *
5079 * @code
5080 *   const Tensor<2,dim,ADNumberType> F_AD
5081 *   = Physics::Elasticity::Kinematics::F(solution_grads_u_f[f_q_point]);
5082 *  
5083 *   const Tensor<2,dim,ADNumberType> F_inv_AD = invert(F_AD);
5084 *   ADNumberType det_F_AD = determinant(F_AD);
5085 *  
5086 *   const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
5087 *   lqph = quadrature_point_history.get_data(cell);
5088 *   Assert(lqph.size() == n_q_points, ExcInternalError());
5089 *  
5090 * @endcode
5091 *
5092 * Seepage velocity
5093 *
5094 * @code
5095 *   Tensor<1,dim> seepage;
5096 *   double det_F = Tensor<0,dim,double>(det_F_AD);
5097 *   const Tensor<1,dim,ADNumberType> grad_p
5098 *   = solution_grads_p_f[f_q_point]*F_inv_AD;
5099 *   const Tensor<1,dim,ADNumberType> seepage_AD
5100 *   = lqph[f_q_point]->get_seepage_velocity_current(F_AD, grad_p);
5101 *  
5102 *   for (unsigned int i=0; i<dim; ++i)
5103 *   seepage[i] = Tensor<0,dim,double>(seepage_AD[i]);
5104 *  
5105 *   sum_total_flow_mpi += (seepage/det_F) * N * JxW_f;
5106 *   }//end gauss points on faces loop
5107 *   }
5108 *   }//end face loop
5109 *   }//end cell loop
5110 *  
5111 * @endcode
5112 *
5113 * Sum the results from different MPI process and then add to the reaction_force vector
5114 * In theory, the solution on each surface (each cell) only exists in one MPI process
5115 * so, we add all MPI process, one will have the solution and the others will be zero
5116 *
5117 * @code
5118 *   for (unsigned int d=0; d<dim; ++d)
5119 *   {
5120 *   reaction_force[d] = Utilities::MPI::sum(sum_reaction_mpi[d],
5121 *   mpi_communicator);
5122 *   reaction_force_pressure[d] = Utilities::MPI::sum(sum_reaction_pressure_mpi[d],
5123 *   mpi_communicator);
5124 *   reaction_force_extra[d] = Utilities::MPI::sum(sum_reaction_extra_mpi[d],
5125 *   mpi_communicator);
5126 *   }
5127 *  
5128 * @endcode
5129 *
5130 * Same for total fluid flow, and for porous and viscous dissipations
5131 *
5132 * @code
5133 *   total_fluid_flow = Utilities::MPI::sum(sum_total_flow_mpi,
5134 *   mpi_communicator);
5135 *   total_porous_dissipation = Utilities::MPI::sum(sum_porous_dissipation_mpi,
5136 *   mpi_communicator);
5137 *   total_viscous_dissipation = Utilities::MPI::sum(sum_viscous_dissipation_mpi,
5138 *   mpi_communicator);
5139 *   total_solid_vol = Utilities::MPI::sum(sum_solid_vol_mpi,
5140 *   mpi_communicator);
5141 *   total_vol_current = Utilities::MPI::sum(sum_vol_current_mpi,
5142 *   mpi_communicator);
5143 *   total_vol_reference = Utilities::MPI::sum(sum_vol_reference_mpi,
5144 *   mpi_communicator);
5145 *  
5146 * @endcode
5147 *
5148 * Extract solution for tracked vectors
5149 * Copying an MPI::BlockVector into MPI::Vector is not possible,
5150 * so we copy each block of MPI::BlockVector into an MPI::Vector
5151 * And then we copy the MPI::Vector into "normal" Vectors
5152 *
5153 * @code
5154 *   TrilinosWrappers::MPI::Vector solution_vector_u_MPI(solution_total.block(u_block));
5155 *   TrilinosWrappers::MPI::Vector solution_vector_p_MPI(solution_total.block(p_fluid_block));
5156 *   Vector<double> solution_u_vector(solution_vector_u_MPI);
5157 *   Vector<double> solution_p_vector(solution_vector_p_MPI);
5158 *  
5159 *   if (this_mpi_process == 0)
5160 *   {
5161 * @endcode
5162 *
5163 * Append the pressure solution vector to the displacement solution vector,
5164 * creating a single solution vector equivalent to the original BlockVector
5165 * so FEFieldFunction will work with the dof_handler_ref.
5166 *
5167 * @code
5168 *   Vector<double> solution_vector(solution_p_vector.size()
5169 *   +solution_u_vector.size());
5170 *  
5171 *   for (unsigned int d=0; d<(solution_u_vector.size()); ++d)
5172 *   solution_vector[d] = solution_u_vector[d];
5173 *  
5174 *   for (unsigned int d=0; d<(solution_p_vector.size()); ++d)
5175 *   solution_vector[solution_u_vector.size()+d] = solution_p_vector[d];
5176 *  
5177 *   Functions::FEFieldFunction<dim,Vector<double>>
5178 *   find_solution(dof_handler_ref, solution_vector);
5179 *  
5180 *   for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5181 *   {
5182 *   Vector<double> update(dim+1);
5183 *   Point<dim> pt_ref;
5184 *  
5185 *   pt_ref[0]= tracked_vertices_IN[p][0];
5186 *   pt_ref[1]= tracked_vertices_IN[p][1];
5187 *   pt_ref[2]= tracked_vertices_IN[p][2];
5188 *  
5189 *   find_solution.vector_value(pt_ref, update);
5190 *  
5191 *   for (unsigned int d=0; d<(dim+1); ++d)
5192 *   {
5193 * @endcode
5194 *
5195 * For values close to zero, set to 0.0
5196 *
5197 * @code
5198 *   if (abs(update[d])<1.5*parameters.tol_u)
5199 *   update[d] = 0.0;
5200 *   solution_vertices[p][d] = update[d];
5201 *   }
5202 *   }
5203 * @endcode
5204 *
5205 * Write the results to the plotting file.
5206 * Add two blank lines between cycles in the cyclic loading examples so GNUPLOT can detect each cycle as a different block
5207 *
5208 * @code
5209 *   if (( (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")||
5210 *   (parameters.geom_type == "Budday_cube_tension_compression")||
5211 *   (parameters.geom_type == "Budday_cube_shear_fully_fixed") ) &&
5212 *   ( (abs(current_time - parameters.end_time/3.) <0.9*parameters.delta_t)||
5213 *   (abs(current_time - 2.*parameters.end_time/3.)<0.9*parameters.delta_t) ) &&
5214 *   parameters.num_cycle_sets == 1 )
5215 *   {
5216 *   plotpointfile << std::endl<< std::endl;
5217 *   }
5218 *   if (( (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")||
5219 *   (parameters.geom_type == "Budday_cube_tension_compression")||
5220 *   (parameters.geom_type == "Budday_cube_shear_fully_fixed") ) &&
5221 *   ( (abs(current_time - parameters.end_time/9.) <0.9*parameters.delta_t)||
5222 *   (abs(current_time - 2.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5223 *   (abs(current_time - 3.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5224 *   (abs(current_time - 5.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5225 *   (abs(current_time - 7.*parameters.end_time/9.)<0.9*parameters.delta_t) ) &&
5226 *   parameters.num_cycle_sets == 2 )
5227 *   {
5228 *   plotpointfile << std::endl<< std::endl;
5229 *   }
5230 *  
5231 *   plotpointfile << std::setprecision(6) << std::scientific;
5232 *   plotpointfile << std::setw(16) << current_time << ","
5233 *   << std::setw(15) << total_vol_reference << ","
5234 *   << std::setw(15) << total_vol_current << ","
5235 *   << std::setw(15) << total_solid_vol << ",";
5236 *  
5237 *   if (current_time == 0.0)
5238 *   {
5239 *   for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5240 *   {
5241 *   for (unsigned int d=0; d<dim; ++d)
5242 *   plotpointfile << std::setw(15) << 0.0 << ",";
5243 *  
5244 *   plotpointfile << std::setw(15) << parameters.drained_pressure << ",";
5245 *   }
5246 *   for (unsigned int d=0; d<(3*dim+2); ++d)
5247 *   plotpointfile << std::setw(15) << 0.0 << ",";
5248 *  
5249 *   plotpointfile << std::setw(15) << 0.0;
5250 *   }
5251 *   else
5252 *   {
5253 *   for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5254 *   for (unsigned int d=0; d<(dim+1); ++d)
5255 *   plotpointfile << std::setw(15) << solution_vertices[p][d]<< ",";
5256 *  
5257 *   for (unsigned int d=0; d<dim; ++d)
5258 *   plotpointfile << std::setw(15) << reaction_force[d] << ",";
5259 *  
5260 *   for (unsigned int d=0; d<dim; ++d)
5261 *   plotpointfile << std::setw(15) << reaction_force_pressure[d] << ",";
5262 *  
5263 *   for (unsigned int d=0; d<dim; ++d)
5264 *   plotpointfile << std::setw(15) << reaction_force_extra[d] << ",";
5265 *  
5266 *   plotpointfile << std::setw(15) << total_fluid_flow << ","
5267 *   << std::setw(15) << total_porous_dissipation<< ","
5268 *   << std::setw(15) << total_viscous_dissipation;
5269 *   }
5270 *   plotpointfile << std::endl;
5271 *   }
5272 *   }
5273 *  
5274 * @endcode
5275 *
5276 * Header for console output file
5277 *
5278 * @code
5279 *   template <int dim>
5280 *   void Solid<dim>::print_console_file_header(std::ofstream &outputfile) const
5281 *   {
5282 *   outputfile << "/*-----------------------------------------------------------------------------------------";
5283 *   outputfile << "\n\n Poro-viscoelastic formulation to solve nonlinear solid mechanics problems using deal.ii";
5284 *   outputfile << "\n\n Problem setup by E Comellas and J-P Pelteret, University of Erlangen-Nuremberg, 2018";
5285 *   outputfile << "\n\n/*-----------------------------------------------------------------------------------------";
5286 *   outputfile << "\n\nCONSOLE OUTPUT: \n\n";
5287 *   }
5288 *  
5289 * @endcode
5290 *
5291 * Header for plotting output file
5292 *
5293 * @code
5294 *   template <int dim>
5295 *   void Solid<dim>::print_plot_file_header(std::vector<Point<dim> > &tracked_vertices,
5296 *   std::ofstream &plotpointfile) const
5297 *   {
5298 *   plotpointfile << "#\n# *** Solution history for tracked vertices -- DOF: 0 = Ux, 1 = Uy, 2 = Uz, 3 = P ***"
5299 *   << std::endl;
5300 *  
5301 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
5302 *   {
5303 *   plotpointfile << "# Point " << p << " coordinates: ";
5304 *   for (unsigned int d=0; d<dim; ++d)
5305 *   {
5306 *   plotpointfile << tracked_vertices[p][d];
5307 *   if (!( (p == tracked_vertices.size()-1) && (d == dim-1) ))
5308 *   plotpointfile << ", ";
5309 *   }
5310 *   plotpointfile << std::endl;
5311 *   }
5312 *   plotpointfile << "# The reaction force is the integral over the loaded surfaces in the "
5313 *   << "undeformed configuration of the Cauchy stress times the normal surface unit vector.\n"
5314 *   << "# reac(p) corresponds to the volumetric part of the Cauchy stress due to the pore fluid pressure"
5315 *   << " and reac(E) corresponds to the extra part of the Cauchy stress due to the solid contribution."
5316 *   << std::endl
5317 *   << "# The fluid flow is the integral over the drained surfaces in the "
5318 *   << "undeformed configuration of the seepage velocity times the normal surface unit vector."
5319 *   << std::endl
5320 *   << "# Column number:"
5321 *   << std::endl
5322 *   << "#";
5323 *  
5324 *   unsigned int columns = 24;
5325 *   for (unsigned int d=1; d<columns; ++d)
5326 *   plotpointfile << std::setw(15)<< d <<",";
5327 *  
5328 *   plotpointfile << std::setw(15)<< columns
5329 *   << std::endl
5330 *   << "#"
5331 *   << std::right << std::setw(16) << "Time,"
5332 *   << std::right << std::setw(16) << "ref vol,"
5333 *   << std::right << std::setw(16) << "def vol,"
5334 *   << std::right << std::setw(16) << "solid vol,";
5335 *   for (unsigned int p=0; p<tracked_vertices.size(); ++p)
5336 *   for (unsigned int d=0; d<(dim+1); ++d)
5337 *   plotpointfile << std::right<< std::setw(11)
5338 *   <<"P" << p << "[" << d << "],";
5339 *  
5340 *   for (unsigned int d=0; d<dim; ++d)
5341 *   plotpointfile << std::right<< std::setw(13)
5342 *   << "reaction [" << d << "],";
5343 *  
5344 *   for (unsigned int d=0; d<dim; ++d)
5345 *   plotpointfile << std::right<< std::setw(13)
5346 *   << "reac(p) [" << d << "],";
5347 *  
5348 *   for (unsigned int d=0; d<dim; ++d)
5349 *   plotpointfile << std::right<< std::setw(13)
5350 *   << "reac(E) [" << d << "],";
5351 *  
5352 *   plotpointfile << std::right<< std::setw(16)<< "fluid flow,"
5353 *   << std::right<< std::setw(16)<< "porous dissip,"
5354 *   << std::right<< std::setw(15)<< "viscous dissip"
5355 *   << std::endl;
5356 *   }
5357 *  
5358 * @endcode
5359 *
5360 * Footer for console output file
5361 *
5362 * @code
5363 *   template <int dim>
5364 *   void Solid<dim>::print_console_file_footer(std::ofstream &outputfile) const
5365 *   {
5366 * @endcode
5367 *
5368 * Copy "parameters" file at end of output file.
5369 *
5370 * @code
5371 *   std::ifstream infile("parameters.prm");
5372 *   std::string content = "";
5373 *   int i;
5374 *  
5375 *   for(i=0 ; infile.eof()!=true ; i++)
5376 *   {
5377 *   char aux = infile.get();
5378 *   content += aux;
5379 *   if(aux=='\n') content += '#';
5380 *   }
5381 *  
5382 *   i--;
5383 *   content.erase(content.end()-1);
5384 *   infile.close();
5385 *  
5386 *   outputfile << "\n\n\n\n PARAMETERS FILE USED IN THIS COMPUTATION: \n#"
5387 *   << std::endl
5388 *   << content;
5389 *   }
5390 *  
5391 * @endcode
5392 *
5393 * Footer for plotting output file
5394 *
5395 * @code
5396 *   template <int dim>
5397 *   void Solid<dim>::print_plot_file_footer(std::ofstream &plotpointfile) const
5398 *   {
5399 * @endcode
5400 *
5401 * Copy "parameters" file at end of output file.
5402 *
5403 * @code
5404 *   std::ifstream infile("parameters.prm");
5405 *   std::string content = "";
5406 *   int i;
5407 *  
5408 *   for(i=0 ; infile.eof()!=true ; i++)
5409 *   {
5410 *   char aux = infile.get();
5411 *   content += aux;
5412 *   if(aux=='\n') content += '#';
5413 *   }
5414 *  
5415 *   i--;
5416 *   content.erase(content.end()-1);
5417 *   infile.close();
5418 *  
5419 *   plotpointfile << "#"<< std::endl
5420 *   << "#"<< std::endl
5421 *   << "# PARAMETERS FILE USED IN THIS COMPUTATION:" << std::endl
5422 *   << "#"<< std::endl
5423 *   << content;
5424 *   }
5425 *  
5426 *  
5427 * @endcode
5428 *
5429 *
5430 * <a name="nonlinear-poro-viscoelasticity.cc-VerificationexamplesfromEhlersandEipper1999"></a>
5431 * <h3>Verification examples from Ehlers and Eipper 1999</h3>
5432 * We group the definition of the geometry, boundary and loading conditions specific to
5433 * the verification examples from Ehlers and Eipper 1999 into specific classes.
5434 *
5435
5436 *
5437 *
5438 * <a name="nonlinear-poro-viscoelasticity.cc-BaseclassTubegeometryandboundaryconditions"></a>
5439 * <h4>Base class: Tube geometry and boundary conditions</h4>
5440 *
5441 * @code
5442 *   template <int dim>
5443 *   class VerificationEhlers1999TubeBase
5444 *   : public Solid<dim>
5445 *   {
5446 *   public:
5447 *   VerificationEhlers1999TubeBase (const Parameters::AllParameters &parameters)
5448 *   : Solid<dim> (parameters)
5449 *   {}
5450 *  
5451 *   virtual ~VerificationEhlers1999TubeBase () {}
5452 *  
5453 *   private:
5454 *   virtual void make_grid() override
5455 *   {
5456 *   GridGenerator::cylinder( this->triangulation,
5457 *   0.1,
5458 *   0.5);
5459 *  
5460 *   const double rot_angle = 3.0*numbers::PI/2.0;
5461 *   GridTools::rotate( Point<3>::unit_vector(1), rot_angle, this->triangulation);
5462 *  
5463 *   this->triangulation.reset_manifold(0);
5464 *   static const CylindricalManifold<dim> manifold_description_3d(2);
5465 *   this->triangulation.set_manifold (0, manifold_description_3d);
5466 *   GridTools::scale(this->parameters.scale, this->triangulation);
5467 *   this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5468 *   this->triangulation.reset_manifold(0);
5469 *   }
5470 *  
5471 *   virtual void define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5472 *   {
5473 *   tracked_vertices[0][0] = 0.0*this->parameters.scale;
5474 *   tracked_vertices[0][1] = 0.0*this->parameters.scale;
5475 *   tracked_vertices[0][2] = 0.5*this->parameters.scale;
5476 *  
5477 *   tracked_vertices[1][0] = 0.0*this->parameters.scale;
5478 *   tracked_vertices[1][1] = 0.0*this->parameters.scale;
5479 *   tracked_vertices[1][2] = -0.5*this->parameters.scale;
5480 *   }
5481 *  
5482 *   virtual void make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5483 *   {
5484 *   if (this->time.get_timestep() < 2)
5485 *   {
5486 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5487 *   2,
5488 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5489 *   constraints,
5490 *   (this->fe.component_mask(this->pressure)));
5491 *   }
5492 *   else
5493 *   {
5494 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5495 *   2,
5496 *   Functions::ZeroFunction<dim>(this->n_components),
5497 *   constraints,
5498 *   (this->fe.component_mask(this->pressure)));
5499 *   }
5500 *  
5501 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5502 *   0,
5503 *   Functions::ZeroFunction<dim>(this->n_components),
5504 *   constraints,
5505 *   (this->fe.component_mask(this->x_displacement)|
5506 *   this->fe.component_mask(this->y_displacement) ) );
5507 *  
5508 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5509 *   1,
5510 *   Functions::ZeroFunction<dim>(this->n_components),
5511 *   constraints,
5512 *   (this->fe.component_mask(this->x_displacement) |
5513 *   this->fe.component_mask(this->y_displacement) |
5514 *   this->fe.component_mask(this->z_displacement) ));
5515 *   }
5516 *  
5517 *   virtual double
5518 *   get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5519 *   const Point<dim> &pt) const override
5520 *   {
5521 *   (void)pt;
5522 *   (void)boundary_id;
5523 *   return 0.0;
5524 *   }
5525 *  
5526 *   virtual types::boundary_id
5527 *   get_reaction_boundary_id_for_output() const override
5528 *   {
5529 *   return 2;
5530 *   }
5531 *  
5532 *   virtual std::pair<types::boundary_id,types::boundary_id>
5533 *   get_drained_boundary_id_for_output() const override
5534 *   {
5535 *   return std::make_pair(2,2);
5536 *   }
5537 *  
5538 *   virtual std::vector<double>
5539 *   get_dirichlet_load(const types::boundary_id &boundary_id,
5540 *   const int &direction) const override
5541 *   {
5542 *   std::vector<double> displ_incr(dim, 0.0);
5543 *   (void)boundary_id;
5544 *   (void)direction;
5545 *   AssertThrow(false, ExcMessage("Displacement loading not implemented for Ehlers verification examples."));
5546 *  
5547 *   return displ_incr;
5548 *   }
5549 *   };
5550 *  
5551 * @endcode
5552 *
5553 *
5554 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassSteploadexample"></a>
5555 * <h4>Derived class: Step load example</h4>
5556 *
5557 * @code
5558 *   template <int dim>
5559 *   class VerificationEhlers1999StepLoad
5560 *   : public VerificationEhlers1999TubeBase<dim>
5561 *   {
5562 *   public:
5563 *   VerificationEhlers1999StepLoad (const Parameters::AllParameters &parameters)
5564 *   : VerificationEhlers1999TubeBase<dim> (parameters)
5565 *   {}
5566 *  
5567 *   virtual ~VerificationEhlers1999StepLoad () {}
5568 *  
5569 *   private:
5570 *   virtual Tensor<1,dim>
5571 *   get_neumann_traction (const types::boundary_id &boundary_id,
5572 *   const Point<dim> &pt,
5573 *   const Tensor<1,dim> &N) const override
5574 *   {
5575 *   if (this->parameters.load_type == "pressure")
5576 *   {
5577 *   if (boundary_id == 2)
5578 *   {
5579 *   return this->parameters.load * N;
5580 *   }
5581 *   }
5582 *  
5583 *   (void)pt;
5584 *  
5585 *   return Tensor<1,dim>();
5586 *   }
5587 *   };
5588 *  
5589 * @endcode
5590 *
5591 *
5592 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassLoadincreasingexample"></a>
5593 * <h4>Derived class: Load increasing example</h4>
5594 *
5595 * @code
5596 *   template <int dim>
5597 *   class VerificationEhlers1999IncreaseLoad
5598 *   : public VerificationEhlers1999TubeBase<dim>
5599 *   {
5600 *   public:
5601 *   VerificationEhlers1999IncreaseLoad (const Parameters::AllParameters &parameters)
5602 *   : VerificationEhlers1999TubeBase<dim> (parameters)
5603 *   {}
5604 *  
5605 *   virtual ~VerificationEhlers1999IncreaseLoad () {}
5606 *  
5607 *   private:
5608 *   virtual Tensor<1,dim>
5609 *   get_neumann_traction (const types::boundary_id &boundary_id,
5610 *   const Point<dim> &pt,
5611 *   const Tensor<1,dim> &N) const override
5612 *   {
5613 *   if (this->parameters.load_type == "pressure")
5614 *   {
5615 *   if (boundary_id == 2)
5616 *   {
5617 *   const double initial_load = this->parameters.load;
5618 *   const double final_load = 20.0*initial_load;
5619 *   const double initial_time = this->time.get_delta_t();
5620 *   const double final_time = this->time.get_end();
5621 *   const double current_time = this->time.get_current();
5622 *   const double load = initial_load + (final_load-initial_load)*(current_time-initial_time)/(final_time-initial_time);
5623 *   return load * N;
5624 *   }
5625 *   }
5626 *  
5627 *   (void)pt;
5628 *  
5629 *   return Tensor<1,dim>();
5630 *   }
5631 *   };
5632 *  
5633 * @endcode
5634 *
5635 *
5636 * <a name="nonlinear-poro-viscoelasticity.cc-ClassConsolidationcube"></a>
5637 * <h4>Class: Consolidation cube</h4>
5638 *
5639 * @code
5640 *   template <int dim>
5641 *   class VerificationEhlers1999CubeConsolidation
5642 *   : public Solid<dim>
5643 *   {
5644 *   public:
5645 *   VerificationEhlers1999CubeConsolidation (const Parameters::AllParameters &parameters)
5646 *   : Solid<dim> (parameters)
5647 *   {}
5648 *  
5649 *   virtual ~VerificationEhlers1999CubeConsolidation () {}
5650 *  
5651 *   private:
5652 *   virtual void
5653 *   make_grid() override
5654 *   {
5655 *   GridGenerator::hyper_rectangle(this->triangulation,
5656 *   Point<dim>(0.0, 0.0, 0.0),
5657 *   Point<dim>(1.0, 1.0, 1.0),
5658 *   true);
5659 *  
5660 *   GridTools::scale(this->parameters.scale, this->triangulation);
5661 *   this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5662 *  
5663 *   typename Triangulation<dim>::active_cell_iterator cell =
5664 *   this->triangulation.begin_active(), endc = this->triangulation.end();
5665 *   for (; cell != endc; ++cell)
5666 *   {
5667 *   for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
5668 *   if (cell->face(face)->at_boundary() == true &&
5669 *   cell->face(face)->center()[2] == 1.0 * this->parameters.scale)
5670 *   {
5671 *   if (cell->face(face)->center()[0] < 0.5 * this->parameters.scale &&
5672 *   cell->face(face)->center()[1] < 0.5 * this->parameters.scale)
5673 *   cell->face(face)->set_boundary_id(100);
5674 *   else
5675 *   cell->face(face)->set_boundary_id(101);
5676 *   }
5677 *   }
5678 *   }
5679 *  
5680 *   virtual void
5681 *   define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5682 *   {
5683 *   tracked_vertices[0][0] = 0.0*this->parameters.scale;
5684 *   tracked_vertices[0][1] = 0.0*this->parameters.scale;
5685 *   tracked_vertices[0][2] = 1.0*this->parameters.scale;
5686 *  
5687 *   tracked_vertices[1][0] = 0.0*this->parameters.scale;
5688 *   tracked_vertices[1][1] = 0.0*this->parameters.scale;
5689 *   tracked_vertices[1][2] = 0.0*this->parameters.scale;
5690 *   }
5691 *  
5692 *   virtual void
5693 *   make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5694 *   {
5695 *   if (this->time.get_timestep() < 2)
5696 *   {
5697 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5698 *   101,
5699 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5700 *   constraints,
5701 *   (this->fe.component_mask(this->pressure)));
5702 *   }
5703 *   else
5704 *   {
5705 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5706 *   101,
5707 *   Functions::ZeroFunction<dim>(this->n_components),
5708 *   constraints,
5709 *   (this->fe.component_mask(this->pressure)));
5710 *   }
5711 *  
5712 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5713 *   0,
5714 *   Functions::ZeroFunction<dim>(this->n_components),
5715 *   constraints,
5716 *   this->fe.component_mask(this->x_displacement));
5717 *  
5718 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5719 *   1,
5720 *   Functions::ZeroFunction<dim>(this->n_components),
5721 *   constraints,
5722 *   this->fe.component_mask(this->x_displacement));
5723 *  
5724 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5725 *   2,
5726 *   Functions::ZeroFunction<dim>(this->n_components),
5727 *   constraints,
5728 *   this->fe.component_mask(this->y_displacement));
5729 *  
5730 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5731 *   3,
5732 *   Functions::ZeroFunction<dim>(this->n_components),
5733 *   constraints,
5734 *   this->fe.component_mask(this->y_displacement));
5735 *  
5736 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5737 *   4,
5738 *   Functions::ZeroFunction<dim>(this->n_components),
5739 *   constraints,
5740 *   ( this->fe.component_mask(this->x_displacement) |
5741 *   this->fe.component_mask(this->y_displacement) |
5742 *   this->fe.component_mask(this->z_displacement) ));
5743 *   }
5744 *  
5745 *   virtual Tensor<1,dim>
5746 *   get_neumann_traction (const types::boundary_id &boundary_id,
5747 *   const Point<dim> &pt,
5748 *   const Tensor<1,dim> &N) const override
5749 *   {
5750 *   if (this->parameters.load_type == "pressure")
5751 *   {
5752 *   if (boundary_id == 100)
5753 *   {
5754 *   return this->parameters.load * N;
5755 *   }
5756 *   }
5757 *  
5758 *   (void)pt;
5759 *  
5760 *   return Tensor<1,dim>();
5761 *   }
5762 *  
5763 *   virtual double
5764 *   get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5765 *   const Point<dim> &pt) const override
5766 *   {
5767 *   (void)pt;
5768 *   (void)boundary_id;
5769 *   return 0.0;
5770 *   }
5771 *  
5772 *   virtual types::boundary_id
5773 *   get_reaction_boundary_id_for_output() const override
5774 *   {
5775 *   return 100;
5776 *   }
5777 *  
5778 *   virtual std::pair<types::boundary_id,types::boundary_id>
5779 *   get_drained_boundary_id_for_output() const override
5780 *   {
5781 *   return std::make_pair(101,101);
5782 *   }
5783 *  
5784 *   virtual std::vector<double>
5785 *   get_dirichlet_load(const types::boundary_id &boundary_id,
5786 *   const int &direction) const override
5787 *   {
5788 *   std::vector<double> displ_incr(dim, 0.0);
5789 *   (void)boundary_id;
5790 *   (void)direction;
5791 *   AssertThrow(false, ExcMessage("Displacement loading not implemented for Ehlers verification examples."));
5792 *  
5793 *   return displ_incr;
5794 *   }
5795 *   };
5796 *  
5797 * @endcode
5798 *
5799 *
5800 * <a name="nonlinear-poro-viscoelasticity.cc-Franceschiniexperiments"></a>
5801 * <h4>Franceschini experiments</h4>
5802 *
5803 * @code
5804 *   template <int dim>
5805 *   class Franceschini2006Consolidation
5806 *   : public Solid<dim>
5807 *   {
5808 *   public:
5809 *   Franceschini2006Consolidation (const Parameters::AllParameters &parameters)
5810 *   : Solid<dim> (parameters)
5811 *   {}
5812 *  
5813 *   virtual ~Franceschini2006Consolidation () {}
5814 *  
5815 *   private:
5816 *   virtual void make_grid() override
5817 *   {
5818 *   const Point<dim-1> mesh_center(0.0, 0.0);
5819 *   const double radius = 0.5;
5820 * @endcode
5821 *
5822 * const double height = 0.27; //8.1 mm for 30 mm radius
5823 *
5824 * @code
5825 *   const double height = 0.23; //6.9 mm for 30 mm radius
5826 *   Triangulation<dim-1> triangulation_in;
5827 *   GridGenerator::hyper_ball( triangulation_in,
5828 *   mesh_center,
5829 *   radius);
5830 *  
5831 *   GridGenerator::extrude_triangulation(triangulation_in,
5832 *   2,
5833 *   height,
5834 *   this->triangulation);
5835 *  
5836 *   const CylindricalManifold<dim> cylinder_3d(2);
5837 *   const types::manifold_id cylinder_id = 0;
5838 *  
5839 *  
5840 *   this->triangulation.set_manifold(cylinder_id, cylinder_3d);
5841 *  
5842 *   for (auto cell : this->triangulation.active_cell_iterators())
5843 *   {
5844 *   for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
5845 *   {
5846 *   if (cell->face(face)->at_boundary() == true)
5847 *   {
5848 *   if (cell->face(face)->center()[2] == 0.0)
5849 *   cell->face(face)->set_boundary_id(1);
5850 *  
5851 *   else if (cell->face(face)->center()[2] == height)
5852 *   cell->face(face)->set_boundary_id(2);
5853 *  
5854 *   else
5855 *   {
5856 *   cell->face(face)->set_boundary_id(0);
5857 *   cell->face(face)->set_all_manifold_ids(cylinder_id);
5858 *   }
5859 *   }
5860 *   }
5861 *   }
5862 *  
5863 *   GridTools::scale(this->parameters.scale, this->triangulation);
5864 *   this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5865 *   }
5866 *  
5867 *   virtual void define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5868 *   {
5869 *   tracked_vertices[0][0] = 0.0*this->parameters.scale;
5870 *   tracked_vertices[0][1] = 0.0*this->parameters.scale;
5871 * @endcode
5872 *
5873 * tracked_vertices[0][2] = 0.27*this->parameters.scale;
5874 *
5875 * @code
5876 *   tracked_vertices[0][2] = 0.23*this->parameters.scale;
5877 *  
5878 *   tracked_vertices[1][0] = 0.0*this->parameters.scale;
5879 *   tracked_vertices[1][1] = 0.0*this->parameters.scale;
5880 *   tracked_vertices[1][2] = 0.0*this->parameters.scale;
5881 *   }
5882 *  
5883 *   virtual void make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5884 *   {
5885 *   if (this->time.get_timestep() < 2)
5886 *   {
5887 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5888 *   1,
5889 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5890 *   constraints,
5891 *   (this->fe.component_mask(this->pressure)));
5892 *  
5893 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5894 *   2,
5895 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5896 *   constraints,
5897 *   (this->fe.component_mask(this->pressure)));
5898 *   }
5899 *   else
5900 *   {
5901 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5902 *   1,
5903 *   Functions::ZeroFunction<dim>(this->n_components),
5904 *   constraints,
5905 *   (this->fe.component_mask(this->pressure)));
5906 *  
5907 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5908 *   2,
5909 *   Functions::ZeroFunction<dim>(this->n_components),
5910 *   constraints,
5911 *   (this->fe.component_mask(this->pressure)));
5912 *   }
5913 *  
5914 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5915 *   0,
5916 *   Functions::ZeroFunction<dim>(this->n_components),
5917 *   constraints,
5918 *   (this->fe.component_mask(this->x_displacement)|
5919 *   this->fe.component_mask(this->y_displacement) ) );
5920 *  
5921 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5922 *   1,
5923 *   Functions::ZeroFunction<dim>(this->n_components),
5924 *   constraints,
5925 *   (this->fe.component_mask(this->x_displacement) |
5926 *   this->fe.component_mask(this->y_displacement) |
5927 *   this->fe.component_mask(this->z_displacement) ));
5928 *  
5929 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5930 *   2,
5931 *   Functions::ZeroFunction<dim>(this->n_components),
5932 *   constraints,
5933 *   (this->fe.component_mask(this->x_displacement) |
5934 *   this->fe.component_mask(this->y_displacement) ));
5935 *   }
5936 *  
5937 *   virtual double
5938 *   get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5939 *   const Point<dim> &pt) const override
5940 *   {
5941 *   (void)pt;
5942 *   (void)boundary_id;
5943 *   return 0.0;
5944 *   }
5945 *  
5946 *   virtual types::boundary_id
5947 *   get_reaction_boundary_id_for_output() const override
5948 *   {
5949 *   return 2;
5950 *   }
5951 *  
5952 *   virtual std::pair<types::boundary_id,types::boundary_id>
5953 *   get_drained_boundary_id_for_output() const override
5954 *   {
5955 *   return std::make_pair(1,2);
5956 *   }
5957 *  
5958 *   virtual std::vector<double>
5959 *   get_dirichlet_load(const types::boundary_id &boundary_id,
5960 *   const int &direction) const override
5961 *   {
5962 *   std::vector<double> displ_incr(dim, 0.0);
5963 *   (void)boundary_id;
5964 *   (void)direction;
5965 *   AssertThrow(false, ExcMessage("Displacement loading not implemented for Franceschini examples."));
5966 *  
5967 *   return displ_incr;
5968 *   }
5969 *  
5970 *   virtual Tensor<1,dim>
5971 *   get_neumann_traction (const types::boundary_id &boundary_id,
5972 *   const Point<dim> &pt,
5973 *   const Tensor<1,dim> &N) const override
5974 *   {
5975 *   if (this->parameters.load_type == "pressure")
5976 *   {
5977 *   if (boundary_id == 2)
5978 *   {
5979 *   return (this->parameters.load * N);
5980 *   /*
5981 *   const double final_load = this->parameters.load;
5982 *   const double final_load_time = 10 * this->time.get_delta_t();
5983 *   const double current_time = this->time.get_current();
5984 *  
5985 *  
5986 *   const double c = final_load_time / 2.0;
5987 *   const double r = 200.0 * 0.03 / c;
5988 *  
5989 *   const double load = final_load * std::exp(r * current_time)
5990 *   / ( std::exp(c * current_time) + std::exp(r * current_time));
5991 *   return load * N;
5992 *   */
5993 *   }
5994 *   }
5995 *  
5996 *   (void)pt;
5997 *  
5998 *   return Tensor<1,dim>();
5999 *   }
6000 *   };
6001 *  
6002 * @endcode
6003 *
6004 *
6005 * <a name="nonlinear-poro-viscoelasticity.cc-ExamplestoreproduceexperimentsbyBuddayetal2017"></a>
6006 * <h3>Examples to reproduce experiments by Budday et al. 2017</h3>
6007 * We group the definition of the geometry, boundary and loading conditions specific to
6008 * the examples to reproduce experiments by Budday et al. 2017 into specific classes.
6009 *
6010
6011 *
6012 *
6013 * <a name="nonlinear-poro-viscoelasticity.cc-BaseclassCubegeometryandloadingpattern"></a>
6014 * <h4>Base class: Cube geometry and loading pattern</h4>
6015 *
6016 * @code
6017 *   template <int dim>
6018 *   class BrainBudday2017BaseCube
6019 *   : public Solid<dim>
6020 *   {
6021 *   public:
6022 *   BrainBudday2017BaseCube (const Parameters::AllParameters &parameters)
6023 *   : Solid<dim> (parameters)
6024 *   {}
6025 *  
6026 *   virtual ~BrainBudday2017BaseCube () {}
6027 *  
6028 *   private:
6029 *   virtual void
6030 *   make_grid() override
6031 *   {
6032 *   GridGenerator::hyper_cube(this->triangulation,
6033 *   0.0,
6034 *   1.0,
6035 *   true);
6036 *  
6037 *   typename Triangulation<dim>::active_cell_iterator cell =
6038 *   this->triangulation.begin_active(), endc = this->triangulation.end();
6039 *   for (; cell != endc; ++cell)
6040 *   {
6041 *   for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
6042 *   if (cell->face(face)->at_boundary() == true &&
6043 *   ( cell->face(face)->boundary_id() == 0 ||
6044 *   cell->face(face)->boundary_id() == 1 ||
6045 *   cell->face(face)->boundary_id() == 2 ||
6046 *   cell->face(face)->boundary_id() == 3 ) )
6047 *  
6048 *   cell->face(face)->set_boundary_id(100);
6049 *  
6050 *   }
6051 *  
6052 *   GridTools::scale(this->parameters.scale, this->triangulation);
6053 *   this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
6054 *   }
6055 *  
6056 *   virtual double
6057 *   get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
6058 *   const Point<dim> &pt) const override
6059 *   {
6060 *   (void)pt;
6061 *   (void)boundary_id;
6062 *   return 0.0;
6063 *   }
6064 *  
6065 *   virtual std::pair<types::boundary_id,types::boundary_id>
6066 *   get_drained_boundary_id_for_output() const override
6067 *   {
6068 *   return std::make_pair(100,100);
6069 *   }
6070 *   };
6071 *  
6072 * @endcode
6073 *
6074 *
6075 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassUniaxialboundaryconditions"></a>
6076 * <h4>Derived class: Uniaxial boundary conditions</h4>
6077 *
6078 * @code
6079 *   template <int dim>
6080 *   class BrainBudday2017CubeTensionCompression
6081 *   : public BrainBudday2017BaseCube<dim>
6082 *   {
6083 *   public:
6084 *   BrainBudday2017CubeTensionCompression (const Parameters::AllParameters &parameters)
6085 *   : BrainBudday2017BaseCube<dim> (parameters)
6086 *   {}
6087 *  
6088 *   virtual ~BrainBudday2017CubeTensionCompression () {}
6089 *  
6090 *   private:
6091 *   virtual void
6092 *   define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6093 *   {
6094 *   tracked_vertices[0][0] = 0.5*this->parameters.scale;
6095 *   tracked_vertices[0][1] = 0.5*this->parameters.scale;
6096 *   tracked_vertices[0][2] = 1.0*this->parameters.scale;
6097 *  
6098 *   tracked_vertices[1][0] = 0.5*this->parameters.scale;
6099 *   tracked_vertices[1][1] = 0.5*this->parameters.scale;
6100 *   tracked_vertices[1][2] = 0.5*this->parameters.scale;
6101 *   }
6102 *  
6103 *   virtual void
6104 *   make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6105 *   {
6106 *   if (this->time.get_timestep() < 2)
6107 *   {
6108 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6109 *   100,
6110 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6111 *   constraints,
6112 *   (this->fe.component_mask(this->pressure)));
6113 *   }
6114 *   else
6115 *   {
6116 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6117 *   100,
6118 *   Functions::ZeroFunction<dim>(this->n_components),
6119 *   constraints,
6120 *   (this->fe.component_mask(this->pressure)));
6121 *   }
6122 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6123 *   4,
6124 *   Functions::ZeroFunction<dim>(this->n_components),
6125 *   constraints,
6126 *   this->fe.component_mask(this->z_displacement) );
6127 *  
6128 *   Point<dim> fix_node(0.5*this->parameters.scale, 0.5*this->parameters.scale, 0.0);
6129 *   typename DoFHandler<dim>::active_cell_iterator
6130 *   cell = this->dof_handler_ref.begin_active(), endc = this->dof_handler_ref.end();
6131 *   for (; cell != endc; ++cell)
6132 *   for (unsigned int node = 0; node < GeometryInfo<dim>::vertices_per_cell; ++node)
6133 *   {
6134 *   if ( (abs(cell->vertex(node)[2]-fix_node[2]) < (1e-6 * this->parameters.scale))
6135 *   && (abs(cell->vertex(node)[0]-fix_node[0]) < (1e-6 * this->parameters.scale)))
6136 *   constraints.add_line(cell->vertex_dof_index(node, 0));
6137 *  
6138 *   if ( (abs(cell->vertex(node)[2]-fix_node[2]) < (1e-6 * this->parameters.scale))
6139 *   && (abs(cell->vertex(node)[1]-fix_node[1]) < (1e-6 * this->parameters.scale)))
6140 *   constraints.add_line(cell->vertex_dof_index(node, 1));
6141 *   }
6142 *  
6143 *   if (this->parameters.load_type == "displacement")
6144 *   {
6145 *   const std::vector<double> value = get_dirichlet_load(5,2);
6146 *   const FEValuesExtractors::Scalar direction(this->z_displacement);
6147 *  
6148 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6149 *   5,
6150 *   Functions::ConstantFunction<dim>(value[2],this->n_components),
6151 *   constraints,
6152 *   this->fe.component_mask(direction));
6153 *   }
6154 *   }
6155 *  
6156 *   virtual Tensor<1,dim>
6157 *   get_neumann_traction (const types::boundary_id &boundary_id,
6158 *   const Point<dim> &pt,
6159 *   const Tensor<1,dim> &N) const override
6160 *   {
6161 *   if (this->parameters.load_type == "pressure")
6162 *   {
6163 *   if (boundary_id == 5)
6164 *   {
6165 *   const double final_load = this->parameters.load;
6166 *   const double current_time = this->time.get_current();
6167 *   const double final_time = this->time.get_end();
6168 *   const double num_cycles = 3.0;
6169 *  
6170 *   return final_load/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5))) * N;
6171 *   }
6172 *   }
6173 *  
6174 *   (void)pt;
6175 *  
6176 *   return Tensor<1,dim>();
6177 *   }
6178 *  
6179 *   virtual types::boundary_id
6180 *   get_reaction_boundary_id_for_output() const override
6181 *   {
6182 *   return 5;
6183 *   }
6184 *  
6185 *   virtual std::vector<double>
6186 *   get_dirichlet_load(const types::boundary_id &boundary_id,
6187 *   const int &direction) const override
6188 *   {
6189 *   std::vector<double> displ_incr(dim,0.0);
6190 *  
6191 *   if ( (boundary_id == 5) && (direction == 2) )
6192 *   {
6193 *   const double final_displ = this->parameters.load;
6194 *   const double current_time = this->time.get_current();
6195 *   const double final_time = this->time.get_end();
6196 *   const double delta_time = this->time.get_delta_t();
6197 *   const double num_cycles = 3.0;
6198 *   double current_displ = 0.0;
6199 *   double previous_displ = 0.0;
6200 *  
6201 *   if (this->parameters.num_cycle_sets == 1)
6202 *   {
6203 *   current_displ = final_displ/2.0 * (1.0
6204 *   - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5)));
6205 *   previous_displ = final_displ/2.0 * (1.0
6206 *   - std::sin(numbers::PI * (2.0*num_cycles*(current_time-delta_time)/final_time + 0.5)));
6207 *   }
6208 *   else
6209 *   {
6210 *   if ( current_time <= (final_time*1.0/3.0) )
6211 *   {
6212 *   current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6213 *   (2.0*num_cycles*current_time/(final_time*1.0/3.0) + 0.5)));
6214 *   previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6215 *   (2.0*num_cycles*(current_time-delta_time)/(final_time*1.0/3.0) + 0.5)));
6216 *   }
6217 *   else
6218 *   {
6219 *   current_displ = final_displ * (1.0 - std::sin(numbers::PI *
6220 *   (2.0*num_cycles*current_time / (final_time*2.0/3.0)
6221 *   - (num_cycles - 0.5) )));
6222 *   previous_displ = final_displ * (1.0 - std::sin(numbers::PI *
6223 *   (2.0*num_cycles*(current_time-delta_time) / (final_time*2.0/3.0)
6224 *   - (num_cycles - 0.5))));
6225 *   }
6226 *   }
6227 *   displ_incr[2] = current_displ - previous_displ;
6228 *   }
6229 *   return displ_incr;
6230 *   }
6231 *   };
6232 *  
6233 * @endcode
6234 *
6235 *
6236 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassNolateraldisplacementinloadingsurfaces"></a>
6237 * <h4>Derived class: No lateral displacement in loading surfaces</h4>
6238 *
6239 * @code
6240 *   template <int dim>
6241 *   class BrainBudday2017CubeTensionCompressionFullyFixed
6242 *   : public BrainBudday2017BaseCube<dim>
6243 *   {
6244 *   public:
6245 *   BrainBudday2017CubeTensionCompressionFullyFixed (const Parameters::AllParameters &parameters)
6246 *   : BrainBudday2017BaseCube<dim> (parameters)
6247 *   {}
6248 *  
6249 *   virtual ~BrainBudday2017CubeTensionCompressionFullyFixed () {}
6250 *  
6251 *   private:
6252 *   virtual void
6253 *   define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6254 *   {
6255 *   tracked_vertices[0][0] = 0.5*this->parameters.scale;
6256 *   tracked_vertices[0][1] = 0.5*this->parameters.scale;
6257 *   tracked_vertices[0][2] = 1.0*this->parameters.scale;
6258 *  
6259 *   tracked_vertices[1][0] = 0.5*this->parameters.scale;
6260 *   tracked_vertices[1][1] = 0.5*this->parameters.scale;
6261 *   tracked_vertices[1][2] = 0.5*this->parameters.scale;
6262 *   }
6263 *  
6264 *   virtual void
6265 *   make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6266 *   {
6267 *   if (this->time.get_timestep() < 2)
6268 *   {
6269 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6270 *   100,
6271 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6272 *   constraints,
6273 *   (this->fe.component_mask(this->pressure)));
6274 *   }
6275 *   else
6276 *   {
6277 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6278 *   100,
6279 *   Functions::ZeroFunction<dim>(this->n_components),
6280 *   constraints,
6281 *   (this->fe.component_mask(this->pressure)));
6282 *   }
6283 *  
6284 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6285 *   4,
6286 *   Functions::ZeroFunction<dim>(this->n_components),
6287 *   constraints,
6288 *   (this->fe.component_mask(this->x_displacement) |
6289 *   this->fe.component_mask(this->y_displacement) |
6290 *   this->fe.component_mask(this->z_displacement) ));
6291 *  
6292 *  
6293 *   if (this->parameters.load_type == "displacement")
6294 *   {
6295 *   const std::vector<double> value = get_dirichlet_load(5,2);
6296 *   const FEValuesExtractors::Scalar direction(this->z_displacement);
6297 *  
6298 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6299 *   5,
6300 *   Functions::ConstantFunction<dim>(value[2],this->n_components),
6301 *   constraints,
6302 *   this->fe.component_mask(direction) );
6303 *  
6304 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6305 *   5,
6306 *   Functions::ZeroFunction<dim>(this->n_components),
6307 *   constraints,
6308 *   (this->fe.component_mask(this->x_displacement) |
6309 *   this->fe.component_mask(this->y_displacement) ));
6310 *   }
6311 *   }
6312 *  
6313 *   virtual Tensor<1,dim>
6314 *   get_neumann_traction (const types::boundary_id &boundary_id,
6315 *   const Point<dim> &pt,
6316 *   const Tensor<1,dim> &N) const override
6317 *   {
6318 *   if (this->parameters.load_type == "pressure")
6319 *   {
6320 *   if (boundary_id == 5)
6321 *   {
6322 *   const double final_load = this->parameters.load;
6323 *   const double current_time = this->time.get_current();
6324 *   const double final_time = this->time.get_end();
6325 *   const double num_cycles = 3.0;
6326 *  
6327 *   return final_load/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5))) * N;
6328 *   }
6329 *   }
6330 *  
6331 *   (void)pt;
6332 *  
6333 *   return Tensor<1,dim>();
6334 *   }
6335 *  
6336 *   virtual types::boundary_id
6337 *   get_reaction_boundary_id_for_output() const override
6338 *   {
6339 *   return 5;
6340 *   }
6341 *  
6342 *   virtual std::vector<double>
6343 *   get_dirichlet_load(const types::boundary_id &boundary_id,
6344 *   const int &direction) const override
6345 *   {
6346 *   std::vector<double> displ_incr(dim,0.0);
6347 *  
6348 *   if ( (boundary_id == 5) && (direction == 2) )
6349 *   {
6350 *   const double final_displ = this->parameters.load;
6351 *   const double current_time = this->time.get_current();
6352 *   const double final_time = this->time.get_end();
6353 *   const double delta_time = this->time.get_delta_t();
6354 *   const double num_cycles = 3.0;
6355 *   double current_displ = 0.0;
6356 *   double previous_displ = 0.0;
6357 *  
6358 *   if (this->parameters.num_cycle_sets == 1)
6359 *   {
6360 *   current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5)));
6361 *   previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*(current_time-delta_time)/final_time + 0.5)));
6362 *   }
6363 *   else
6364 *   {
6365 *   if ( current_time <= (final_time*1.0/3.0) )
6366 *   {
6367 *   current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6368 *   (2.0*num_cycles*current_time/(final_time*1.0/3.0) + 0.5)));
6369 *   previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6370 *   (2.0*num_cycles*(current_time-delta_time)/(final_time*1.0/3.0) + 0.5)));
6371 *   }
6372 *   else
6373 *   {
6374 *   current_displ = final_displ * (1.0 - std::sin(numbers::PI *
6375 *   (2.0*num_cycles*current_time / (final_time*2.0/3.0)
6376 *   - (num_cycles - 0.5) )));
6377 *   previous_displ = final_displ * (1.0 - std::sin(numbers::PI *
6378 *   (2.0*num_cycles*(current_time-delta_time) / (final_time*2.0/3.0)
6379 *   - (num_cycles - 0.5))));
6380 *   }
6381 *   }
6382 *   displ_incr[2] = current_displ - previous_displ;
6383 *   }
6384 *   return displ_incr;
6385 *   }
6386 *   };
6387 *  
6388 * @endcode
6389 *
6390 *
6391 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassNolateralorverticaldisplacementinloadingsurface"></a>
6392 * <h4>Derived class: No lateral or vertical displacement in loading surface</h4>
6393 *
6394 * @code
6395 *   template <int dim>
6396 *   class BrainBudday2017CubeShearFullyFixed
6397 *   : public BrainBudday2017BaseCube<dim>
6398 *   {
6399 *   public:
6400 *   BrainBudday2017CubeShearFullyFixed (const Parameters::AllParameters &parameters)
6401 *   : BrainBudday2017BaseCube<dim> (parameters)
6402 *   {}
6403 *  
6404 *   virtual ~BrainBudday2017CubeShearFullyFixed () {}
6405 *  
6406 *   private:
6407 *   virtual void
6408 *   define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6409 *   {
6410 *   tracked_vertices[0][0] = 0.75*this->parameters.scale;
6411 *   tracked_vertices[0][1] = 0.5*this->parameters.scale;
6412 *   tracked_vertices[0][2] = 0.0*this->parameters.scale;
6413 *  
6414 *   tracked_vertices[1][0] = 0.25*this->parameters.scale;
6415 *   tracked_vertices[1][1] = 0.5*this->parameters.scale;
6416 *   tracked_vertices[1][2] = 0.0*this->parameters.scale;
6417 *   }
6418 *  
6419 *   virtual void
6420 *   make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6421 *   {
6422 *   if (this->time.get_timestep() < 2)
6423 *   {
6424 *   VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6425 *   100,
6426 *   Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6427 *   constraints,
6428 *   (this->fe.component_mask(this->pressure)));
6429 *   }
6430 *   else
6431 *   {
6432 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6433 *   100,
6434 *   Functions::ZeroFunction<dim>(this->n_components),
6435 *   constraints,
6436 *   (this->fe.component_mask(this->pressure)));
6437 *   }
6438 *  
6439 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6440 *   5,
6441 *   Functions::ZeroFunction<dim>(this->n_components),
6442 *   constraints,
6443 *   (this->fe.component_mask(this->x_displacement) |
6444 *   this->fe.component_mask(this->y_displacement) |
6445 *   this->fe.component_mask(this->z_displacement) ));
6446 *  
6447 *  
6448 *   if (this->parameters.load_type == "displacement")
6449 *   {
6450 *   const std::vector<double> value = get_dirichlet_load(4,0);
6451 *   const FEValuesExtractors::Scalar direction(this->x_displacement);
6452 *  
6453 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6454 *   4,
6455 *   Functions::ConstantFunction<dim>(value[0],this->n_components),
6456 *   constraints,
6457 *   this->fe.component_mask(direction));
6458 *  
6459 *   VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6460 *   4,
6461 *   Functions::ZeroFunction<dim>(this->n_components),
6462 *   constraints,
6463 *   (this->fe.component_mask(this->y_displacement) |
6464 *   this->fe.component_mask(this->z_displacement) ));
6465 *   }
6466 *   }
6467 *  
6468 *   virtual Tensor<1,dim>
6469 *   get_neumann_traction (const types::boundary_id &boundary_id,
6470 *   const Point<dim> &pt,
6471 *   const Tensor<1,dim> &N) const override
6472 *   {
6473 *   if (this->parameters.load_type == "pressure")
6474 *   {
6475 *   if (boundary_id == 4)
6476 *   {
6477 *   const double final_load = this->parameters.load;
6478 *   const double current_time = this->time.get_current();
6479 *   const double final_time = this->time.get_end();
6480 *   const double num_cycles = 3.0;
6481 *   const Tensor<1,3> axis ({0.0,1.0,0.0});
6482 *   const double angle = numbers::PI;
6483 *   static const Tensor< 2, dim, double> R(Physics::Transformations::Rotations::rotation_matrix_3d(axis,angle));
6484 *  
6485 *   return (final_load * (std::sin(2.0*(numbers::PI)*num_cycles*current_time/final_time)) * (R * N));
6486 *   }
6487 *   }
6488 *  
6489 *   (void)pt;
6490 *  
6491 *   return Tensor<1,dim>();
6492 *   }
6493 *  
6494 *   virtual types::boundary_id
6495 *   get_reaction_boundary_id_for_output() const override
6496 *   {
6497 *   return 4;
6498 *   }
6499 *  
6500 *   virtual std::vector<double>
6501 *   get_dirichlet_load(const types::boundary_id &boundary_id,
6502 *   const int &direction) const override
6503 *   {
6504 *   std::vector<double> displ_incr (dim, 0.0);
6505 *  
6506 *   if ( (boundary_id == 4) && (direction == 0) )
6507 *   {
6508 *   const double final_displ = this->parameters.load;
6509 *   const double current_time = this->time.get_current();
6510 *   const double final_time = this->time.get_end();
6511 *   const double delta_time = this->time.get_delta_t();
6512 *   const double num_cycles = 3.0;
6513 *   double current_displ = 0.0;
6514 *   double previous_displ = 0.0;
6515 *  
6516 *   if (this->parameters.num_cycle_sets == 1)
6517 *   {
6518 *   current_displ = final_displ * (std::sin(2.0*(numbers::PI)*num_cycles*current_time/final_time));
6519 *   previous_displ = final_displ * (std::sin(2.0*(numbers::PI)*num_cycles*(current_time-delta_time)/final_time));
6520 *   }
6521 *   else
6522 *   {
6523 *   AssertThrow(false, ExcMessage("Problem type not defined. Budday shear experiments implemented only for one set of cycles."));
6524 *   }
6525 *   displ_incr[0] = current_displ - previous_displ;
6526 *   }
6527 *   return displ_incr;
6528 *   }
6529 *   };
6530 *  
6531 *   }
6532 *  
6533 * @endcode
6534 *
6535 *
6536 * <a name="nonlinear-poro-viscoelasticity.cc-Mainfunction"></a>
6537 * <h3>Main function</h3>
6538 * Lastly we provide the main driver function which is similar to the other tutorials.
6539 *
6540 * @code
6541 *   int main (int argc, char *argv[])
6542 *   {
6543 *   using namespace dealii;
6544 *   using namespace NonLinearPoroViscoElasticity;
6545 *  
6546 *   const unsigned int n_tbb_processes = 1;
6547 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, n_tbb_processes);
6548 *  
6549 *   try
6550 *   {
6551 *   Parameters::AllParameters parameters ("parameters.prm");
6552 *   if (parameters.geom_type == "Ehlers_tube_step_load")
6553 *   {
6554 *   VerificationEhlers1999StepLoad<3> solid_3d(parameters);
6555 *   solid_3d.run();
6556 *   }
6557 *   else if (parameters.geom_type == "Ehlers_tube_increase_load")
6558 *   {
6559 *   VerificationEhlers1999IncreaseLoad<3> solid_3d(parameters);
6560 *   solid_3d.run();
6561 *   }
6562 *   else if (parameters.geom_type == "Ehlers_cube_consolidation")
6563 *   {
6564 *   VerificationEhlers1999CubeConsolidation<3> solid_3d(parameters);
6565 *   solid_3d.run();
6566 *   }
6567 *   else if (parameters.geom_type == "Franceschini_consolidation")
6568 *   {
6569 *   Franceschini2006Consolidation<3> solid_3d(parameters);
6570 *   solid_3d.run();
6571 *   }
6572 *   else if (parameters.geom_type == "Budday_cube_tension_compression")
6573 *   {
6574 *   BrainBudday2017CubeTensionCompression<3> solid_3d(parameters);
6575 *   solid_3d.run();
6576 *   }
6577 *   else if (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")
6578 *   {
6579 *   BrainBudday2017CubeTensionCompressionFullyFixed<3> solid_3d(parameters);
6580 *   solid_3d.run();
6581 *   }
6582 *   else if (parameters.geom_type == "Budday_cube_shear_fully_fixed")
6583 *   {
6584 *   BrainBudday2017CubeShearFullyFixed<3> solid_3d(parameters);
6585 *   solid_3d.run();
6586 *   }
6587 *   else
6588 *   {
6589 *   AssertThrow(false, ExcMessage("Problem type not defined. Current setting: " + parameters.geom_type));
6590 *   }
6591 *  
6592 *   }
6593 *   catch (std::exception &exc)
6594 *   {
6595 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
6596 *   {
6597 *   std::cerr << std::endl << std::endl
6598 *   << "----------------------------------------------------"
6599 *   << std::endl;
6600 *   std::cerr << "Exception on processing: " << std::endl << exc.what()
6601 *   << std::endl << "Aborting!" << std::endl
6602 *   << "----------------------------------------------------"
6603 *   << std::endl;
6604 *  
6605 *   return 1;
6606 *   }
6607 *   }
6608 *   catch (...)
6609 *   {
6610 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
6611 *   {
6612 *   std::cerr << std::endl << std::endl
6613 *   << "----------------------------------------------------"
6614 *   << std::endl;
6615 *   std::cerr << "Unknown exception!" << std::endl << "Aborting!"
6616 *   << std::endl
6617 *   << "----------------------------------------------------"
6618 *   << std::endl;
6619 *   return 1;
6620 *   }
6621 *   }
6622 *   return 0;
6623 *   }
6624 * @endcode
6625
6626
6627*/
DataPostprocessorVector(const std::string &name, const UpdateFlags update_flags)
Definition fe_q.h:554
static constexpr const SymmetricTensor< 2, dim > I
Definition point.h:113
@ wall_times
Definition timer.h:651
@ maximum_smoothing
Definition tria.h:1539
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
#define DEAL_II_VERSION_GTE(major, minor, subminor)
Definition config.h:407
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
const bool IsBlockVector< VectorType >::value
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
UpdateFlags
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
const Event initial
Definition event.cc:68
void approximate(const SynchronousIterators< std::tuple< typename DoFHandler< dim, spacedim >::active_cell_iterator, Vector< float >::iterator > > &cell, const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component)
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >())
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
unsigned int count_cells_with_subdomain_association(const Triangulation< dim, spacedim > &triangulation, const types::subdomain_id subdomain)
double volume(const Triangulation< dim, spacedim > &tria)
@ valid
Iterator points to a valid object.
@ matrix
Contents is actually a matrix.
@ diagonal
Matrix is diagonal.
constexpr char N
constexpr types::blas_int zero
constexpr char A
constexpr types::blas_int one
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
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:193
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(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)
constexpr ReturnType< rank, T >::value_type & extract(T &t, const ArrayType &indices)
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:105
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:120
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
void abort(const ExceptionBase &exc) noexcept
bool check(const ConstraintKinds kind_in, const unsigned int dim)
long double gamma(const unsigned int n)
void copy(const T *begin, const T *end, U *dest)
int(& functions)(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition loop.h:70
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int material_id
Definition types.h:184
unsigned int boundary_id
Definition types.h:161
unsigned int global_dof_index
Definition types.h:94
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 > &)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
SymmetricTensorEigenvectorMethod