598 * We start by including all the necessary deal.II header files and some
C++
599 * related ones. They have been discussed in detail in previous tutorial
600 * programs, so you need only refer to past tutorials
for details.
606 * #include <deal.II/base/function.h>
607 * #include <deal.II/base/parameter_handler.h>
608 * #include <deal.II/base/
point.h>
609 * #include <deal.II/base/quadrature_lib.h>
610 * #include <deal.II/base/symmetric_tensor.h>
611 * #include <deal.II/base/tensor.h>
612 * #include <deal.II/base/timer.h>
613 * #include <deal.II/base/work_stream.h>
614 * #include <deal.II/base/mpi.h>
615 * #include <deal.II/base/quadrature_point_data.h>
617 * #include <deal.II/differentiation/ad.h>
619 * #include <deal.II/distributed/shared_tria.h>
621 * #include <deal.II/dofs/dof_renumbering.h>
622 * #include <deal.II/dofs/dof_tools.h>
623 * #include <deal.II/dofs/dof_accessor.h>
625 * #include <deal.II/grid/filtered_iterator.h>
626 * #include <deal.II/grid/grid_generator.h>
627 * #include <deal.II/grid/grid_tools.h>
628 * #include <deal.II/grid/grid_in.h>
629 * #include <deal.II/grid/grid_out.h>
630 * #include <deal.II/grid/manifold_lib.h>
631 * #include <deal.II/grid/tria_accessor.h>
632 * #include <deal.II/grid/tria_iterator.h>
634 * #include <deal.II/fe/fe_dgp_monomial.h>
635 * #include <deal.II/fe/fe_q.h>
636 * #include <deal.II/fe/fe_system.h>
637 * #include <deal.II/fe/fe_tools.h>
638 * #include <deal.II/fe/fe_values.h>
640 * #include <deal.II/lac/block_sparsity_pattern.h>
641 * #include <deal.II/lac/affine_constraints.h>
642 * #include <deal.II/lac/dynamic_sparsity_pattern.h>
643 * #include <deal.II/lac/full_matrix.h>
645 * #include <deal.II/lac/packaged_operation.h>
647 * #include <deal.II/lac/trilinos_block_sparse_matrix.h>
648 * #include <deal.II/lac/trilinos_linear_operator.h>
649 * #include <deal.II/lac/trilinos_parallel_block_vector.h>
650 * #include <deal.II/lac/trilinos_precondition.h>
651 * #include <deal.II/lac/trilinos_sparse_matrix.h>
652 * #include <deal.II/lac/trilinos_sparsity_pattern.h>
653 * #include <deal.II/lac/trilinos_solver.h>
654 * #include <deal.II/lac/trilinos_vector.h>
656 * #include <deal.II/lac/block_vector.h>
657 * #include <deal.II/lac/vector.h>
659 * #include <deal.II/numerics/data_postprocessor.h>
660 * #include <deal.II/numerics/data_out.h>
661 * #include <deal.II/numerics/data_out_faces.h>
662 * #include <deal.II/numerics/fe_field_function.h>
663 * #include <deal.II/numerics/vector_tools.h>
665 * #include <deal.II/physics/transformations.h>
666 * #include <deal.II/physics/elasticity/kinematics.h>
667 * #include <deal.II/physics/elasticity/standard_tensors.h>
669 * #include <iostream>
677 * We create a
namespace for everything that relates to
678 * the nonlinear poro-viscoelastic formulation,
679 * and
import all the deal.II function and class names into it:
682 * namespace NonLinearPoroViscoElasticity
684 * using namespace dealii;
689 * <a name=
"nonlinear-poro-viscoelasticity.cc-Runtimeparameters"></a>
690 * <h3>Run-time parameters</h3>
695 * introduced by the user through the file
"parameters.prm"
698 *
namespace Parameters
703 * <a name=
"nonlinear-poro-viscoelasticity.cc-FiniteElementsystem"></a>
704 * <h4>Finite Element system</h4>
705 * Here we specify the polynomial order used to
approximate the solution,
706 * both
for the displacements and pressure unknowns.
707 * The quadrature order should be adjusted accordingly.
712 *
unsigned int poly_degree_displ;
713 *
unsigned int poly_degree_pore;
714 *
unsigned int quad_order;
725 * prm.enter_subsection(
"Finite element system");
727 * prm.declare_entry(
"Polynomial degree displ",
"2",
729 *
"Displacement system polynomial order");
731 * prm.declare_entry(
"Polynomial degree pore",
"1",
733 *
"Pore pressure system polynomial order");
735 * prm.declare_entry(
"Quadrature order",
"3",
737 *
"Gauss quadrature order");
739 * prm.leave_subsection();
744 * prm.enter_subsection(
"Finite element system");
746 * poly_degree_displ = prm.get_integer(
"Polynomial degree displ");
747 * poly_degree_pore = prm.get_integer(
"Polynomial degree pore");
748 * quad_order = prm.get_integer(
"Quadrature order");
750 * prm.leave_subsection();
756 * <a name=
"nonlinear-poro-viscoelasticity.cc-Geometry"></a>
758 * These parameters are related to the geometry definition and mesh generation.
759 * We select the type of problem to solve and introduce the desired load values.
764 * std::string geom_type;
765 *
unsigned int global_refinement;
767 * std::string load_type;
769 *
unsigned int num_cycle_sets;
771 *
double drained_pressure;
782 * prm.enter_subsection(
"Geometry");
784 * prm.declare_entry(
"Geometry type",
"Ehlers_tube_step_load",
786 *
"|Ehlers_tube_increase_load"
787 *
"|Ehlers_cube_consolidation"
788 *
"|Franceschini_consolidation"
789 *
"|Budday_cube_tension_compression"
790 *
"|Budday_cube_tension_compression_fully_fixed"
791 *
"|Budday_cube_shear_fully_fixed"),
792 *
"Type of geometry used. "
793 *
"For Ehlers verification examples see Ehlers and Eipper (1999). "
794 *
"For Franceschini brain consolidation see Franceschini et al. (2006)"
795 *
"For Budday brain examples see Budday et al. (2017)");
797 * prm.declare_entry(
"Global refinement",
"1",
799 *
"Global refinement level");
801 * prm.declare_entry(
"Grid scale",
"1.0",
803 *
"Global grid scaling factor");
805 * prm.declare_entry(
"Load type",
"pressure",
807 *
"Type of loading");
809 * prm.declare_entry(
"Load value",
"-7.5e+6",
813 * prm.declare_entry(
"Number of cycle sets",
"1",
815 *
"Number of times each set of 3 cycles is repeated, only for "
816 *
"Budday_cube_tension_compression and Budday_cube_tension_compression_fully_fixed. "
817 *
"Load value is doubled in second set, load rate is kept constant."
818 *
"Final time indicates end of second cycle set.");
820 * prm.declare_entry(
"Fluid flow value",
"0.0",
822 *
"Prescribed fluid flow. Not implemented in any example yet.");
824 * prm.declare_entry(
"Drained pressure",
"0.0",
826 *
"Increase of pressure value at drained boundary w.r.t the atmospheric pressure.");
828 * prm.leave_subsection();
833 * prm.enter_subsection(
"Geometry");
835 * geom_type = prm.get(
"Geometry type");
836 * global_refinement = prm.get_integer(
"Global refinement");
837 *
scale = prm.get_double(
"Grid scale");
838 * load_type = prm.get(
"Load type");
839 * load = prm.get_double(
"Load value");
840 * num_cycle_sets = prm.get_integer(
"Number of cycle sets");
841 * fluid_flow = prm.get_double(
"Fluid flow value");
842 * drained_pressure = prm.get_double(
"Drained pressure");
844 * prm.leave_subsection();
850 * <a name=
"nonlinear-poro-viscoelasticity.cc-Materials"></a>
855 * Here we select the type of material
for the solid component
856 * and define the corresponding material parameters.
857 * Then we define he fluid data, including the type of
858 * seepage velocity definition to use.
863 * std::string mat_type;
869 *
double alpha1_infty;
870 *
double alpha2_infty;
871 *
double alpha3_infty;
875 *
double alpha1_mode_1;
876 *
double alpha2_mode_1;
877 *
double alpha3_mode_1;
878 *
double viscosity_mode_1;
879 * std::string fluid_type;
880 *
double solid_vol_frac;
881 *
double kappa_darcy;
882 *
double init_intrinsic_perm;
883 *
double viscosity_FR;
884 *
double init_darcy_coef;
887 *
int gravity_direction;
888 *
double gravity_value;
902 * prm.enter_subsection(
"Material properties");
904 * prm.declare_entry(
"material",
"Neo-Hooke",
906 *
"Type of material used in the problem");
908 * prm.declare_entry(
"lambda",
"8.375e6",
910 *
"First Lamé parameter for extension function related to compactation point in solid material [Pa].");
912 * prm.declare_entry(
"shear modulus",
"5.583e6",
914 *
"shear modulus for Neo-Hooke materials [Pa].");
916 * prm.declare_entry(
"eigen solver",
"QL Implicit Shifts",
918 *
"The type of eigen solver to be used for Ogden and visco-Ogden models.");
920 * prm.declare_entry(
"mu1",
"0.0",
922 *
"Shear material parameter 'mu1' for Ogden material [Pa].");
924 * prm.declare_entry(
"mu2",
"0.0",
926 *
"Shear material parameter 'mu2' for Ogden material [Pa].");
928 * prm.declare_entry(
"mu3",
"0.0",
930 *
"Shear material parameter 'mu1' for Ogden material [Pa].");
932 * prm.declare_entry(
"alpha1",
"1.0",
934 *
"Stiffness material parameter 'alpha1' for Ogden material [-].");
936 * prm.declare_entry(
"alpha2",
"1.0",
938 *
"Stiffness material parameter 'alpha2' for Ogden material [-].");
940 * prm.declare_entry(
"alpha3",
"1.0",
942 *
"Stiffness material parameter 'alpha3' for Ogden material [-].");
944 * prm.declare_entry(
"mu1_1",
"0.0",
946 *
"Shear material parameter 'mu1' for first viscous mode in Ogden material [Pa].");
948 * prm.declare_entry(
"mu2_1",
"0.0",
950 *
"Shear material parameter 'mu2' for first viscous mode in Ogden material [Pa].");
952 * prm.declare_entry(
"mu3_1",
"0.0",
954 *
"Shear material parameter 'mu1' for first viscous mode in Ogden material [Pa].");
956 * prm.declare_entry(
"alpha1_1",
"1.0",
958 *
"Stiffness material parameter 'alpha1' for first viscous mode in Ogden material [-].");
960 * prm.declare_entry(
"alpha2_1",
"1.0",
962 *
"Stiffness material parameter 'alpha2' for first viscous mode in Ogden material [-].");
964 * prm.declare_entry(
"alpha3_1",
"1.0",
966 *
"Stiffness material parameter 'alpha3' for first viscous mode in Ogden material [-].");
968 * prm.declare_entry(
"viscosity_1",
"1e-10",
970 *
"Deformation-independent viscosity parameter 'eta_1' for first viscous mode in Ogden material [-].");
972 * prm.declare_entry(
"seepage definition",
"Ehlers",
974 *
"Type of formulation used to define the seepage velocity in the problem. "
975 *
"Choose between Markert formulation of deformation-dependent intrinsic permeability "
976 *
"and Ehlers formulation of deformation-dependent Darcy flow coefficient.");
978 * prm.declare_entry(
"initial solid volume fraction",
"0.67",
980 *
"Initial porosity (solid volume fraction, 0 < n_0s < 1)");
982 * prm.declare_entry(
"kappa",
"0.0",
984 *
"Deformation-dependency control parameter for specific permeability (kappa >= 0)");
986 * prm.declare_entry(
"initial intrinsic permeability",
"0.0",
988 *
"Initial intrinsic permeability parameter [m^2] (isotropic permeability). To be used with Markert formulation.");
990 * prm.declare_entry(
"fluid viscosity",
"0.0",
992 *
"Effective shear viscosity parameter of the fluid [Pa·s, (N·s)/m^2]. To be used with Markert formulation.");
994 * prm.declare_entry(
"initial Darcy coefficient",
"1.0e-4",
996 *
"Initial Darcy flow coefficient [m/s] (isotropic permeability). To be used with Ehlers formulation.");
998 * prm.declare_entry(
"fluid weight",
"1.0e4",
1000 *
"Effective weight of the fluid [N/m^3]. To be used with Ehlers formulation.");
1002 * prm.declare_entry(
"gravity term",
"false",
1004 *
"Gravity term considered (true) or neglected (false)");
1006 * prm.declare_entry(
"fluid density",
"1.0",
1008 *
"Real (or effective) density of the fluid");
1010 * prm.declare_entry(
"solid density",
"1.0",
1012 *
"Real (or effective) density of the solid");
1014 * prm.declare_entry(
"gravity direction",
"2",
1016 *
"Direction of gravity (unit vector 0 for x, 1 for y, 2 for z)");
1018 * prm.declare_entry(
"gravity value",
"-9.81",
1020 *
"Value of gravity (be careful to have consistent units!)");
1022 * prm.leave_subsection();
1027 * prm.enter_subsection(
"Material properties");
1034 * mat_type = prm.get(
"material");
1035 *
lambda = prm.get_double(
"lambda");
1036 * mu = prm.get_double(
"shear modulus");
1037 * mu1_infty = prm.get_double(
"mu1");
1038 * mu2_infty = prm.get_double(
"mu2");
1039 * mu3_infty = prm.get_double(
"mu3");
1040 * alpha1_infty = prm.get_double(
"alpha1");
1041 * alpha2_infty = prm.get_double(
"alpha2");
1042 * alpha3_infty = prm.get_double(
"alpha3");
1043 * mu1_mode_1 = prm.get_double(
"mu1_1");
1044 * mu2_mode_1 = prm.get_double(
"mu2_1");
1045 * mu3_mode_1 = prm.get_double(
"mu3_1");
1046 * alpha1_mode_1 = prm.get_double(
"alpha1_1");
1047 * alpha2_mode_1 = prm.get_double(
"alpha2_1");
1048 * alpha3_mode_1 = prm.get_double(
"alpha3_1");
1049 * viscosity_mode_1 = prm.get_double(
"viscosity_1");
1055 * fluid_type = prm.get(
"seepage definition");
1056 * solid_vol_frac = prm.get_double(
"initial solid volume fraction");
1057 * kappa_darcy = prm.get_double(
"kappa");
1058 * init_intrinsic_perm = prm.get_double(
"initial intrinsic permeability");
1059 * viscosity_FR = prm.get_double(
"fluid viscosity");
1060 * init_darcy_coef = prm.get_double(
"initial Darcy coefficient");
1061 * weight_FR = prm.get_double(
"fluid weight");
1067 * gravity_term = prm.get_bool(
"gravity term");
1068 * density_FR = prm.get_double(
"fluid density");
1069 * density_SR = prm.get_double(
"solid density");
1070 * gravity_direction = prm.get_integer(
"gravity direction");
1071 * gravity_value = prm.get_double(
"gravity value");
1073 *
if ( (fluid_type ==
"Markert") && ((init_intrinsic_perm == 0.0) || (viscosity_FR == 0.0)) )
1074 *
AssertThrow(
false, ExcMessage(
"Markert seepage velocity formulation requires the definition of "
1075 *
"'initial intrinsic permeability' and 'fluid viscosity' greater than 0.0."));
1077 *
if ( (fluid_type ==
"Ehlers") && ((init_darcy_coef == 0.0) || (weight_FR == 0.0)) )
1078 *
AssertThrow(
false, ExcMessage(
"Ehler seepage velocity formulation requires the definition of "
1079 *
"'initial Darcy coefficient' and 'fluid weight' greater than 0.0."));
1081 *
const std::string eigen_solver_type = prm.get(
"eigen solver");
1082 *
if (eigen_solver_type ==
"QL Implicit Shifts")
1084 *
else if (eigen_solver_type ==
"Jacobi")
1088 *
AssertThrow(
false, ExcMessage(
"Unknown eigen solver selected."));
1091 * prm.leave_subsection();
1097 * <a name=
"nonlinear-poro-viscoelasticity.cc-Nonlinearsolver"></a>
1098 * <h4>Nonlinear solver</h4>
1102 * We now define the tolerances and the maximum number of iterations
for the
1103 * Newton-Raphson scheme used to solve the nonlinear system of governing equations.
1106 *
struct NonlinearSolver
1108 *
unsigned int max_iterations_NR;
1111 *
double tol_p_fluid;
1122 * prm.enter_subsection(
"Nonlinear solver");
1124 * prm.declare_entry(
"Max iterations Newton-Raphson",
"15",
1126 *
"Number of Newton-Raphson iterations allowed");
1128 * prm.declare_entry(
"Tolerance force",
"1.0e-8",
1130 *
"Force residual tolerance");
1132 * prm.declare_entry(
"Tolerance displacement",
"1.0e-6",
1134 *
"Displacement error tolerance");
1136 * prm.declare_entry(
"Tolerance pore pressure",
"1.0e-6",
1138 *
"Pore pressure error tolerance");
1140 * prm.leave_subsection();
1145 * prm.enter_subsection(
"Nonlinear solver");
1147 * max_iterations_NR = prm.get_integer(
"Max iterations Newton-Raphson");
1148 * tol_f = prm.get_double(
"Tolerance force");
1149 * tol_u = prm.get_double(
"Tolerance displacement");
1150 * tol_p_fluid = prm.get_double(
"Tolerance pore pressure");
1152 * prm.leave_subsection();
1158 * <a name=
"nonlinear-poro-viscoelasticity.cc-Time"></a>
1160 * Here we
set the timestep size @f$ \varDelta t @f$ and the simulation
end-time.
1176 * prm.enter_subsection(
"Time");
1178 * prm.declare_entry(
"End time",
"10.0",
1182 * prm.declare_entry(
"Time step size",
"0.002",
1184 *
"Time step size. The value must be larger than the displacement error tolerance defined.");
1186 * prm.leave_subsection();
1191 * prm.enter_subsection(
"Time");
1193 * end_time = prm.get_double(
"End time");
1194 * delta_t = prm.get_double(
"Time step size");
1196 * prm.leave_subsection();
1203 * <a name=
"nonlinear-poro-viscoelasticity.cc-Output"></a>
1205 * We can choose the frequency of the data
for the output files.
1208 *
struct OutputParam
1211 * std::string outfiles_requested;
1212 *
unsigned int timestep_output;
1213 * std::string outtype;
1224 * prm.enter_subsection(
"Output parameters");
1226 * prm.declare_entry(
"Output files",
"true",
1228 *
"Paraview output files to generate.");
1229 * prm.declare_entry(
"Time step number output",
"1",
1231 *
"Output data for time steps multiple of the given "
1232 *
"integer value.");
1233 * prm.declare_entry(
"Averaged results",
"nodes",
1235 *
"Output data associated with integration point values"
1236 *
" averaged on elements or on nodes.");
1238 * prm.leave_subsection();
1243 * prm.enter_subsection(
"Output parameters");
1245 * outfiles_requested = prm.get(
"Output files");
1246 * timestep_output = prm.get_integer(
"Time step number output");
1247 * outtype = prm.get(
"Averaged results");
1249 * prm.leave_subsection();
1255 * <a name=
"nonlinear-poro-viscoelasticity.cc-Allparameters"></a>
1256 * <h4>All parameters</h4>
1257 * We
finally consolidate all of the above structures into a single container that holds all the
run-time selections.
1260 *
struct AllParameters :
public FESystem,
1263 *
public NonlinearSolver,
1265 *
public OutputParam
1267 * AllParameters(
const std::string &input_file);
1276 * AllParameters::AllParameters(
const std::string &input_file)
1279 * declare_parameters(prm);
1280 * prm.parse_input(input_file);
1281 * parse_parameters(prm);
1286 * FESystem::declare_parameters(prm);
1287 * Geometry::declare_parameters(prm);
1288 * Materials::declare_parameters(prm);
1289 * NonlinearSolver::declare_parameters(prm);
1290 * Time::declare_parameters(prm);
1291 * OutputParam::declare_parameters(prm);
1296 * FESystem::parse_parameters(prm);
1297 * Geometry::parse_parameters(prm);
1298 * Materials::parse_parameters(prm);
1299 * NonlinearSolver::parse_parameters(prm);
1300 * Time::parse_parameters(prm);
1301 * OutputParam::parse_parameters(prm);
1308 * <a name=
"nonlinear-poro-viscoelasticity.cc-Timeclass"></a>
1309 * <h3>Time
class</h3>
1310 * A simple
class to store time data.
1311 * For simplicity we assume a
constant time step size.
1317 * Time (
const double time_end,
1318 *
const double delta_t)
1321 * time_current(0.0),
1322 * time_end(time_end),
1329 *
double get_current() const
1331 *
return time_current;
1333 *
double get_end() const
1337 *
double get_delta_t() const
1341 *
unsigned int get_timestep() const
1345 *
void increment_time ()
1347 * time_current += delta_t;
1352 *
unsigned int timestep;
1353 *
double time_current;
1355 *
const double delta_t;
1361 * <a name=
"nonlinear-poro-viscoelasticity.cc-Constitutiveequationforthesolidcomponentofthebiphasicmaterial"></a>
1362 * <h3>Constitutive equation
for the solid component of the biphasic material</h3>
1367 * <a name=
"nonlinear-poro-viscoelasticity.cc-Baseclassgenerichyperelasticmaterial"></a>
1368 * <h4>Base
class:
generic hyperelastic material</h4>
1369 * The
"extra" Kirchhoff stress in the solid component is the
sum of isochoric
1370 * and a volumetric part:
1371 * @f$\mathbf{\tau} = \mathbf{\tau}_E^{(\bullet)} + \mathbf{\tau}^{\textrm{vol}}@f$.
1372 * The deviatoric part changes depending on the type of material model selected:
1373 * Neo-Hooken hyperelasticity, Ogden hyperelasticiy,
1374 * or a single-mode finite viscoelasticity based on the Ogden hyperelastic model.
1375 * In
this base
class we declare it as a virtual function,
1376 * and it will be defined
for each model type in the corresponding derived
class.
1377 * We define here the volumetric component, which depends on the
1378 * extension function @f$U(J_S)@f$ selected, and in
this case is the same
for all models.
1379 * We use the function proposed by
1380 * Ehlers & Eipper 1999 doi:10.1023/A:1006565509095.
1381 * We also define some
public functions to access and update the
internal variables.
1384 *
template <
int dim,
typename NumberType = Sacado::Fad::DFad<
double> >
1385 *
class Material_Hyperelastic
1388 * Material_Hyperelastic(
const Parameters::AllParameters ¶meters,
1391 * n_OS (parameters.solid_vol_frac),
1392 *
lambda (parameters.lambda),
1395 * det_F_converged (1.0),
1396 * eigen_solver (parameters.eigen_solver)
1398 * ~Material_Hyperelastic()
1404 *
return ( get_tau_E_base(F) + get_tau_E_ext_func(F) );
1411 *
Assert(det_F > 0, ExcInternalError());
1412 *
return get_tau_E(F)*NumberType(1/det_F);
1416 * get_converged_det_F() const
1418 *
return det_F_converged;
1422 * update_end_timestep()
1424 * det_F_converged = det_F;
1434 * get_viscous_dissipation( )
const = 0;
1436 *
const double n_OS;
1440 *
double det_F_converged;
1448 *
Assert(det_F > 0, ExcInternalError());
1452 *
return ( NumberType(lambda * (1.0-n_OS)*(1.0-n_OS)
1453 * * (det_F/(1.0-n_OS) - det_F/(det_F-n_OS))) * I );
1463 * <a name=
"nonlinear-poro-viscoelasticity.cc-DerivedclassNeoHookeanhyperelasticmaterial"></a>
1464 * <h4>Derived
class: Neo-Hookean hyperelastic material</h4>
1467 *
template <
int dim,
typename NumberType = Sacado::Fad::DFad<
double> >
1468 *
class NeoHooke :
public Material_Hyperelastic < dim, NumberType >
1471 * NeoHooke(
const Parameters::AllParameters ¶meters,
1474 * Material_Hyperelastic< dim, NumberType > (parameters,time),
1477 *
virtual ~NeoHooke()
1481 * get_viscous_dissipation() const override
1495 *
const bool use_standard_model =
true;
1497 *
if (use_standard_model)
1501 * Standard Neo-Hooke
1510 * Neo-Hooke in terms of principal stretches
1515 *
const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1519 *
for (
unsigned int d=0;
d<dim; ++
d)
1522 *
return ( mu*(B_ev-I) );
1530 * <a name=
"nonlinear-poro-viscoelasticity.cc-DerivedclassOgdenhyperelasticmaterial"></a>
1531 * <h4>Derived
class: Ogden hyperelastic material</h4>
1534 *
template <
int dim,
typename NumberType = Sacado::Fad::DFad<
double> >
1535 *
class Ogden :
public Material_Hyperelastic < dim, NumberType >
1538 * Ogden(
const Parameters::AllParameters ¶meters,
1541 * Material_Hyperelastic< dim, NumberType > (parameters,time),
1542 * mu({parameters.mu1_infty,
1543 * parameters.mu2_infty,
1544 * parameters.mu3_infty}),
1545 * alpha({parameters.alpha1_infty,
1546 * parameters.alpha2_infty,
1547 * parameters.alpha3_infty})
1553 * get_viscous_dissipation() const override
1559 * std::vector<double> mu;
1560 * std::vector<double> alpha;
1568 *
const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1575 *
for (
unsigned int i = 0; i < 3; ++i)
1577 *
for (
unsigned int A = 0; A < dim; ++A)
1581 * tau_aux1 *= mu[i]*
std::pow(eigen_B[A].
first, (alpha[i]/2.) );
1585 * tau_aux2 *= mu[i];
1595 * <a name=
"nonlinear-poro-viscoelasticity.cc-DerivedclassSinglemodeOgdenviscoelasticmaterial"></a>
1596 * <h4>Derived
class: Single-mode Ogden viscoelastic material</h4>
1597 * We use the finite viscoelastic model described in
1598 * Reese & Govindjee (1998) doi:10.1016/S0020-7683(97)00217-5
1599 * The algorithm for the implicit exponential time integration is given in
1600 * Budday et al. (2017) doi: 10.1016/j.actbio.2017.06.024
1603 * template <
int dim, typename NumberType = Sacado::Fad::DFad<
double> >
1604 * class visco_Ogden : public Material_Hyperelastic < dim, NumberType >
1607 * visco_Ogden(
const Parameters::AllParameters ¶meters,
1610 * Material_Hyperelastic< dim, NumberType > (parameters,time),
1611 * mu_infty({parameters.mu1_infty,
1612 * parameters.mu2_infty,
1613 * parameters.mu3_infty}),
1614 * alpha_infty({parameters.alpha1_infty,
1615 * parameters.alpha2_infty,
1616 * parameters.alpha3_infty}),
1617 * mu_mode_1({parameters.mu1_mode_1,
1618 * parameters.mu2_mode_1,
1619 * parameters.mu3_mode_1}),
1620 * alpha_mode_1({parameters.alpha1_mode_1,
1621 * parameters.alpha2_mode_1,
1622 * parameters.alpha3_mode_1}),
1623 * viscosity_mode_1(parameters.viscosity_mode_1),
1627 *
virtual ~visco_Ogden()
1633 * Material_Hyperelastic < dim, NumberType >::update_internal_equilibrium(
F);
1635 * this->Cinv_v_1 = this->Cinv_v_1_converged;
1638 *
const std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim >
1639 * eigen_B_e_1_tr =
eigenvectors(B_e_1_tr, this->eigen_solver);
1643 *
for (
int a = 0; a < dim; ++a)
1646 * epsilon_e_1_tr[a] =
std::log(lambdas_e_1_tr[a]);
1649 *
const double tolerance = 1
e-8;
1650 *
double residual_check = tolerance*10.0;
1656 * std::vector<NumberType> lambdas_e_1_iso(dim);
1658 *
int iteration = 0;
1662 * epsilon_e_1 = epsilon_e_1_tr;
1664 *
while(residual_check > tolerance)
1666 * NumberType aux_J_e_1 = 1.0;
1667 *
for (
unsigned int a = 0; a < dim; ++a)
1669 * lambdas_e_1[a] =
std::exp(epsilon_e_1[a]);
1670 * aux_J_e_1 *= lambdas_e_1[a];
1673 * J_e_1 = aux_J_e_1;
1675 *
for (
unsigned int a = 0; a < dim; ++a)
1676 * lambdas_e_1_iso[a] = lambdas_e_1[a]*
std::pow(J_e_1,-1.0/dim);
1678 *
for (
unsigned int a = 0; a < dim; ++a)
1680 * residual[a] = get_beta_mode_1(lambdas_e_1_iso, a);
1681 * residual[a] *= this->time.get_delta_t()/(2.0*viscosity_mode_1);
1682 * residual[a] += epsilon_e_1[a];
1683 * residual[a] -= epsilon_e_1_tr[a];
1685 *
for (
unsigned int b = 0;
b < dim; ++
b)
1687 * tangent[a][
b] = get_gamma_mode_1(lambdas_e_1_iso, a, b);
1688 * tangent[a][
b] *= this->time.get_delta_t()/(2.0*viscosity_mode_1);
1689 * tangent[a][
b] += I[a][
b];
1693 * epsilon_e_1 -=
invert(tangent)*residual;
1695 * residual_check = 0.0;
1696 *
for (
unsigned int a = 0; a < dim; ++a)
1698 *
if (
std::abs(residual[a]) > residual_check)
1702 *
if (iteration > 15 )
1703 *
AssertThrow(
false, ExcMessage(
"No convergence in local Newton iteration for the "
1704 *
"viscoelastic exponential time integration algorithm."));
1707 * NumberType aux_J_e_1 = 1.0;
1708 *
for (
unsigned int a = 0; a < dim; ++a)
1710 * lambdas_e_1[a] =
std::exp(epsilon_e_1[a]);
1711 * aux_J_e_1 *= lambdas_e_1[a];
1713 * J_e_1 = aux_J_e_1;
1715 *
for (
unsigned int a = 0; a < dim; ++a)
1716 * lambdas_e_1_iso[a] = lambdas_e_1[a]*
std::pow(J_e_1,-1.0/dim);
1718 *
for (
unsigned int a = 0; a < dim; ++a)
1722 * B_e_1_aux *= lambdas_e_1[a] * lambdas_e_1[a];
1723 * B_e_1 += B_e_1_aux;
1728 * this->tau_neq_1 = 0;
1729 *
for (
unsigned int a = 0; a < dim; ++a)
1733 * tau_neq_1_aux *= get_beta_mode_1(lambdas_e_1_iso, a);
1734 * this->tau_neq_1 += tau_neq_1_aux;
1742 *
for (
unsigned int a = 0; a < dim; ++a)
1743 *
for (
unsigned int b = 0;
b < dim; ++
b)
1747 *
void update_end_timestep() override
1749 * Material_Hyperelastic < dim, NumberType >::update_end_timestep();
1750 * this->Cinv_v_1_converged = this->Cinv_v_1;
1753 *
double get_viscous_dissipation() const override
1755 * NumberType dissipation_term = get_tau_E_neq() * get_tau_E_neq();
1756 * dissipation_term /= (2*viscosity_mode_1);
1758 *
return dissipation_term.val();
1762 * std::vector<double> mu_infty;
1763 * std::vector<double> alpha_infty;
1764 * std::vector<double> mu_mode_1;
1765 * std::vector<double> alpha_mode_1;
1766 *
double viscosity_mode_1;
1774 *
return ( get_tau_E_neq() + get_tau_E_eq(F) );
1782 * std::array< std::pair< NumberType, Tensor< 1, dim, NumberType > >, dim > eigen_B;
1789 *
for (
unsigned int i = 0; i < 3; ++i)
1791 *
for (
unsigned int A = 0; A < dim; ++A)
1795 * tau_aux1 *= mu_infty[i]*
std::pow(eigen_B[A].
first, (alpha_infty[i]/2.) );
1799 * tau_aux2 *= mu_infty[i];
1806 * get_tau_E_neq() const
1812 * get_beta_mode_1(std::vector< NumberType > &lambda,
const int &A)
const
1814 * NumberType beta = 0.0;
1816 *
for (
unsigned int i = 0; i < 3; ++i)
1819 * NumberType aux = 0.0;
1820 *
for (
int p = 0; p < dim; ++p)
1821 * aux +=
std::pow(lambda[p],alpha_mode_1[i]);
1824 * aux +=
std::pow(lambda[A], alpha_mode_1[i]);
1825 * aux *= mu_mode_1[i];
1833 * get_gamma_mode_1(std::vector< NumberType > &lambda,
1835 *
const int &B )
const
1837 * NumberType
gamma = 0.0;
1841 *
for (
unsigned int i = 0; i < 3; ++i)
1843 * NumberType aux = 0.0;
1844 *
for (
int p = 0; p < dim; ++p)
1845 * aux +=
std::pow(lambda[p],alpha_mode_1[i]);
1847 * aux *= 1.0/(dim*dim);
1848 * aux += 1.0/dim *
std::pow(lambda[A], alpha_mode_1[i]);
1849 * aux *= mu_mode_1[i]*alpha_mode_1[i];
1856 *
for (
unsigned int i = 0; i < 3; ++i)
1858 * NumberType aux = 0.0;
1859 *
for (
int p = 0; p < dim; ++p)
1860 * aux +=
std::pow(lambda[p],alpha_mode_1[i]);
1862 * aux *= 1.0/(dim*dim);
1863 * aux -= 1.0/dim *
std::pow(lambda[A], alpha_mode_1[i]);
1864 * aux -= 1.0/dim *
std::pow(lambda[B], alpha_mode_1[i]);
1865 * aux *= mu_mode_1[i]*alpha_mode_1[i];
1879 * <a name=
"nonlinear-poro-viscoelasticity.cc-Constitutiveequationforthefluidcomponentofthebiphasicmaterial"></a>
1880 * <h3>Constitutive equation
for the fluid component of the biphasic material</h3>
1881 * We consider two slightly different definitions to define the seepage velocity with a Darcy-like law.
1882 * Ehlers & Eipper 1999, doi:10.1023/A:1006565509095
1883 * Markert 2007, doi:10.1007/s11242-007-9107-6
1884 * The selection of one or another is made by the user via the parameters file.
1887 *
template <
int dim,
typename NumberType = Sacado::Fad::DFad<
double> >
1888 *
class Material_Darcy_Fluid
1891 * Material_Darcy_Fluid(
const Parameters::AllParameters ¶meters)
1893 * fluid_type(parameters.fluid_type),
1894 * n_OS(parameters.solid_vol_frac),
1895 * initial_intrinsic_permeability(parameters.init_intrinsic_perm),
1896 * viscosity_FR(parameters.viscosity_FR),
1897 * initial_darcy_coefficient(parameters.init_darcy_coef),
1898 * weight_FR(parameters.weight_FR),
1899 * kappa_darcy(parameters.kappa_darcy),
1900 * gravity_term(parameters.gravity_term),
1901 * density_FR(parameters.density_FR),
1902 * gravity_direction(parameters.gravity_direction),
1903 * gravity_value(parameters.gravity_value)
1905 *
Assert(kappa_darcy >= 0, ExcInternalError());
1907 * ~Material_Darcy_Fluid()
1915 *
Assert(det_F > 0.0, ExcInternalError());
1919 *
if (fluid_type ==
"Markert")
1920 * permeability_term = get_instrinsic_permeability_current(F) / viscosity_FR;
1922 *
else if (fluid_type ==
"Ehlers")
1923 * permeability_term = get_darcy_flow_current(F) / weight_FR;
1927 *
"Material_Darcy_Fluid --> Only Markert "
1928 *
"and Ehlers formulations have been implemented."));
1930 *
return ( -1.0 * permeability_term * det_F
1931 * * (grad_p_fluid - get_body_force_FR_current()) );
1937 * NumberType dissipation_term;
1942 *
Assert(det_F > 0.0, ExcInternalError());
1944 *
if (fluid_type ==
"Markert")
1946 * permeability_term = get_instrinsic_permeability_current(F) / viscosity_FR;
1947 * seepage_velocity = get_seepage_velocity_current(F,grad_p_fluid);
1949 *
else if (fluid_type ==
"Ehlers")
1951 * permeability_term = get_darcy_flow_current(F) / weight_FR;
1952 * seepage_velocity = get_seepage_velocity_current(F,grad_p_fluid);
1956 *
"Material_Darcy_Fluid --> Only Markert and Ehlers "
1957 *
"formulations have been implemented."));
1959 * dissipation_term = (
invert(permeability_term) * seepage_velocity ) * seepage_velocity;
1960 * dissipation_term *= 1.0/(det_F*det_F);
1965 *
const std::string fluid_type;
1966 *
const double n_OS;
1967 *
const double initial_intrinsic_permeability;
1968 *
const double viscosity_FR;
1969 *
const double initial_darcy_coefficient;
1970 *
const double weight_FR;
1971 *
const double kappa_darcy;
1972 *
const bool gravity_term;
1973 *
const double density_FR;
1974 *
const int gravity_direction;
1975 *
const double gravity_value;
1986 *
Assert(det_F > 0.0, ExcInternalError());
1988 *
const NumberType fraction = (det_F - n_OS)/(1 - n_OS);
1989 *
return ( NumberType (
std::pow(fraction, kappa_darcy))
1990 * * initial_instrinsic_permeability_tensor );
2002 *
Assert(det_F > 0.0, ExcInternalError());
2004 *
const NumberType fraction = (1.0 - (n_OS / det_F) )/(1.0 - n_OS);
2005 *
return ( NumberType (
std::pow(fraction, kappa_darcy))
2006 * * initial_darcy_flow_tensor);
2010 * get_body_force_FR_current() const
2014 *
if (gravity_term ==
true)
2017 * gravity_vector[gravity_direction] = gravity_value;
2018 * body_force_FR_current = density_FR * gravity_vector;
2020 *
return body_force_FR_current;
2027 * <a name=
"nonlinear-poro-viscoelasticity.cc-Quadraturepointhistory"></a>
2029 * As seen in @ref step_18
"step-18", the <code> PointHistory </code>
class offers a method
2030 *
for storing data at the quadrature points. Here each quadrature
point
2031 * holds a pointer to a material description. Thus, different material models
2032 * can be used in different regions of the domain. Among other data, we
2033 * choose to store the
"extra" Kirchhoff stress @f$\boldsymbol{\tau}_E@f$ and
2034 * the dissipation values @f$\mathcal{D}_p@f$ and @f$\mathcal{D}_v@f$.
2037 *
template <
int dim,
typename NumberType = Sacado::Fad::DFad<
double> >
2038 *
class PointHistory
2044 *
virtual ~PointHistory()
2047 *
void setup_lqp (
const Parameters::AllParameters ¶meters,
2050 *
if (parameters.mat_type ==
"Neo-Hooke")
2051 * solid_material.reset(
new NeoHooke<dim,NumberType>(parameters,time));
2052 *
else if (parameters.mat_type ==
"Ogden")
2053 * solid_material.reset(
new Ogden<dim,NumberType>(parameters,time));
2054 *
else if (parameters.mat_type ==
"visco-Ogden")
2055 * solid_material.reset(
new visco_Ogden<dim,NumberType>(parameters,time));
2057 *
Assert (
false, ExcMessage(
"Material type not implemented"));
2059 * fluid_material.reset(
new Material_Darcy_Fluid<dim,NumberType>(parameters));
2065 *
return solid_material->get_tau_E(F);
2071 *
return solid_material->get_Cauchy_E(F);
2075 * get_converged_det_F() const
2077 *
return solid_material->get_converged_det_F();
2081 * update_end_timestep()
2083 * solid_material->update_end_timestep();
2089 * solid_material->update_internal_equilibrium(F);
2093 * get_viscous_dissipation() const
2095 *
return solid_material->get_viscous_dissipation();
2102 *
return fluid_material->get_seepage_velocity_current(F, grad_p_fluid);
2109 *
return fluid_material->get_porous_dissipation(F, grad_p_fluid);
2114 *
const Parameters::AllParameters ¶meters)
const
2118 *
if (parameters.gravity_term ==
true)
2121 *
Assert(det_F_AD > 0.0, ExcInternalError());
2123 *
const NumberType overall_density_ref
2124 * = parameters.density_SR * parameters.solid_vol_frac
2125 * + parameters.density_FR
2126 * * (det_F_AD - parameters.solid_vol_frac);
2129 * gravity_vector[parameters.gravity_direction] = parameters.gravity_value;
2130 * body_force = overall_density_ref * gravity_vector;
2133 *
return body_force;
2136 * std::shared_ptr< Material_Hyperelastic<dim, NumberType> > solid_material;
2137 * std::shared_ptr< Material_Darcy_Fluid<dim, NumberType> > fluid_material;
2143 * <a name=
"nonlinear-poro-viscoelasticity.cc-Nonlinearporoviscoelasticsolid"></a>
2144 * <h3>Nonlinear poro-viscoelastic solid</h3>
2145 * The Solid
class is the central class as it represents the problem at hand:
2146 * the nonlinear poro-viscoelastic solid
2149 * template <int dim>
2153 * Solid(
const Parameters::AllParameters ¶meters);
2158 *
using ADNumberType = Sacado::Fad::DFad<double>;
2160 * std::ofstream outfile;
2161 * std::ofstream pointfile;
2163 *
struct PerTaskData_ASM;
2164 *
template<
typename NumberType =
double>
struct ScratchData_ASM;
2171 *
virtual void make_grid() = 0;
2175 * Define points
for post-processing
2178 *
virtual void define_tracked_vertices(std::vector<
Point<dim> > &tracked_vertices) = 0;
2182 * Set up the finite element system to be solved:
2189 * Extract sub-blocks from the global
matrix
2192 *
void determine_component_extractors();
2196 * Several
functions to
assemble the system and right hand side matrices
using multithreading.
2199 *
void assemble_system
2201 *
void assemble_system_one_cell
2203 * ScratchData_ASM<ADNumberType> &scratch,
2204 * PerTaskData_ASM &data)
const;
2205 *
void copy_local_to_global_system(
const PerTaskData_ASM &data);
2209 * Define boundary conditions
2212 *
virtual void make_constraints(
const int &it_nr);
2218 *
virtual double get_prescribed_fluid_flow
2222 * get_reaction_boundary_id_for_output ()
const = 0;
2223 *
virtual std::pair<types::boundary_id,types::boundary_id>
2224 * get_drained_boundary_id_for_output ()
const = 0;
2225 *
virtual std::vector<double> get_dirichlet_load
2227 *
const int &direction)
const = 0;
2231 * Create and update the quadrature points.
2238 * Solve non-linear system
using a Newton-Raphson scheme
2245 * Solve the linearized equations
using a direct solver
2252 * Retrieve the solution
2260 * Store the converged values of the
internal variables at the
end of each timestep
2263 *
void update_end_timestep();
2267 * Post-processing and writing data to files
2270 *
void output_results_to_vtu(
const unsigned int timestep,
2271 *
const double current_time,
2273 *
void output_results_to_plot(
const unsigned int timestep,
2274 *
const double current_time,
2276 * std::vector<
Point<dim> > &tracked_vertices,
2277 * std::ofstream &pointfile)
const;
2281 * Headers and footer
for the output files
2284 *
void print_console_file_header( std::ofstream &outfile)
const;
2285 *
void print_plot_file_header(std::vector<
Point<dim> > &tracked_vertices,
2286 * std::ofstream &pointfile)
const;
2287 *
void print_console_file_footer(std::ofstream &outfile)
const;
2288 *
void print_plot_file_footer( std::ofstream &pointfile)
const;
2302 * A collection of the parameters used to describe the problem setup
2305 *
const Parameters::AllParameters ¶meters;
2316 * Keep track of the current time and the time spent evaluating certain
functions
2325 * A storage
object for quadrature
point information.
2332 * Integers to store polynomial degree (needed
for output)
2335 *
const unsigned int degree_displ;
2336 *
const unsigned int degree_pore;
2340 * Declare an instance of
dealii FESystem class (finite element definition)
2354 * Integer to store DoFs per element (
this value will be used often)
2357 *
const unsigned int dofs_per_cell;
2361 * Declare an instance of
dealii Extractor objects used to retrieve information from the solution vectors
2362 * We will use
"u_fe" and
"p_fluid_fe"as subscript in
operator [] expressions on
FEValues and
FEFaceValues
2363 * objects to
extract the components of the displacement vector and fluid pressure, respectively.
2371 * Description of how the block-system is arranged. There are 3 blocks:
2372 * 0 - vector DOF displacements u
2373 * 1 -
scalar DOF fluid pressure p_fluid
2376 *
static const unsigned int n_blocks = 2;
2377 *
static const unsigned int n_components = dim+1;
2378 *
static const unsigned int first_u_component = 0;
2379 *
static const unsigned int p_fluid_component = dim;
2402 * std::vector<unsigned int> block_component;
2409 * std::vector<IndexSet> all_locally_owned_dofs;
2412 * std::vector<IndexSet> locally_owned_partitioning;
2413 * std::vector<IndexSet> locally_relevant_partitioning;
2415 * std::vector<types::global_dof_index> dofs_per_block;
2416 * std::vector<types::global_dof_index> element_indices_u;
2417 * std::vector<types::global_dof_index> element_indices_p_fluid;
2421 * Declare an instance of
dealii QGauss class (The Gauss-Legendre family of quadrature rules
for numerical integration)
2422 * Gauss Points in element, with n quadrature points (in each space direction <dim> )
2428 * Gauss Points on element faces (used
for definition of BCs)
2431 *
const QGauss<dim - 1> qf_face;
2434 * Integer to store num GPs per element (
this value will be used often)
2437 *
const unsigned int n_q_points;
2440 * Integer to store num GPs per face (
this value will be used often)
2443 *
const unsigned int n_q_points_f;
2447 * Declare an instance of
dealii AffineConstraints class (linear constraints on DoFs due to hanging nodes or BCs)
2454 * Declare an instance of
dealii classes necessary
for FE system
set-up and assembly
2462 * Right hand side vector of forces
2468 * Total displacement values + pressure (accumulated solution to FE system)
2475 * Non-block system
for the direct solver. We will
copy the block system into these to solve the linearized system of equations.
2483 * We define variables to store norms and update norms and normalisation factors.
2490 *
norm(1.0), u(1.0), p_fluid(1.0)
2499 *
void normalise(
const Errors &rhs)
2501 *
if (rhs.norm != 0.0)
2505 *
if (rhs.p_fluid != 0.0)
2506 * p_fluid /= rhs.p_fluid;
2509 *
double norm, u, p_fluid;
2514 * Declare several instances of the
"Error" structure
2517 * Errors error_residual, error_residual_0, error_residual_norm, error_update,
2518 * error_update_0, error_update_norm;
2522 * Methods to calculate error measures
2525 *
void get_error_residual(Errors &error_residual_OUT);
2526 *
void get_error_update
2528 * Errors &error_update_OUT);
2532 * Print information to screen
2535 *
void print_conv_header();
2536 *
void print_conv_footer();
2540 * NOTE: In all
functions, we pass by reference (&), so these
functions work on the original
copy (not a clone copy),
2549 * <a name=
"nonlinear-poro-viscoelasticity.cc-ImplementationofthecodeSolidcodeclass"></a>
2550 * <h3>Implementation of the <code>Solid</code>
class</h3>
2552 * <a name=
"nonlinear-poro-viscoelasticity.cc-Publicinterface"></a>
2553 * <h4>Public interface</h4>
2554 * We initialise the Solid
class using data extracted from the parameter file.
2557 *
template <
int dim>
2558 * Solid<dim>::Solid(
const Parameters::AllParameters ¶meters)
2560 * mpi_communicator(MPI_COMM_WORLD),
2563 * pcout(std::cout, this_mpi_process == 0),
2564 * parameters(parameters),
2566 * time(parameters.end_time, parameters.delta_t),
2567 * timerconsole( mpi_communicator,
2571 * timerfile( mpi_communicator,
2575 * degree_displ(parameters.poly_degree_displ),
2576 * degree_pore(parameters.poly_degree_pore),
2577 * fe(
FE_Q<dim>(parameters.poly_degree_displ), dim,
2578 *
FE_Q<dim>(parameters.poly_degree_pore), 1 ),
2580 * dofs_per_cell (fe.dofs_per_cell),
2581 * u_fe(first_u_component),
2582 * p_fluid_fe(p_fluid_component),
2583 * x_displacement(first_u_component),
2584 * y_displacement(first_u_component+1),
2585 * z_displacement(first_u_component+2),
2586 * pressure(p_fluid_component),
2587 * dofs_per_block(n_blocks),
2588 * qf_cell(parameters.quad_order),
2589 * qf_face(parameters.quad_order),
2590 * n_q_points (qf_cell.size()),
2591 * n_q_points_f (qf_face.size())
2593 *
Assert(dim==3, ExcMessage(
"This problem only works in 3 space dimensions."));
2594 * determine_component_extractors();
2599 * The
class destructor simply clears the data held by the DOFHandler
2602 *
template <
int dim>
2603 * Solid<dim>::~Solid()
2605 * dof_handler_ref.clear();
2610 * Runs the 3D solid problem
2613 *
template <
int dim>
2614 *
void Solid<dim>::run()
2618 * The current solution increment is defined as a block vector to reflect the structure
2619 * of the PDE system, with multiple solution components
2629 *
if (this_mpi_process == 0)
2631 * outfile.open(
"console-output.sol");
2632 * print_console_file_header(outfile);
2644 * Assign DOFs and create the stiffness and right-hand-side force vector
2647 * system_setup(solution_delta);
2651 * Define points
for post-processing
2654 * std::vector<Point<dim> > tracked_vertices (2);
2655 * define_tracked_vertices(tracked_vertices);
2656 * std::vector<Point<dim>> reaction_force;
2658 *
if (this_mpi_process == 0)
2660 * pointfile.open(
"data-for-gnuplot.sol");
2661 * print_plot_file_header(tracked_vertices, pointfile);
2666 * Print results to output file
2669 *
if (parameters.outfiles_requested ==
"true")
2671 * output_results_to_vtu(time.get_timestep(),
2672 * time.get_current(),
2676 * output_results_to_plot(time.get_timestep(),
2677 * time.get_current(),
2684 * Increment time step (=load step)
2685 * NOTE: In solving the quasi-
static problem, the time becomes a loading parameter,
2686 * i.e. we increase the loading linearly with time, making the two
concepts interchangeable.
2689 * time.increment_time();
2693 * Print information on screen
2696 * pcout <<
"\nSolver:";
2697 * pcout <<
"\n CST = make constraints";
2698 * pcout <<
"\n ASM_SYS = assemble system";
2699 * pcout <<
"\n SLV = linear solver \n";
2703 * Print information on file
2706 * outfile <<
"\nSolver:";
2707 * outfile <<
"\n CST = make constraints";
2708 * outfile <<
"\n ASM_SYS = assemble system";
2709 * outfile <<
"\n SLV = linear solver \n";
2711 *
while ( (time.get_end() - time.get_current()) > -1.0*parameters.tol_u )
2715 * Initialize the current solution increment to zero
2718 * solution_delta = 0.0;
2722 * Solve the non-linear system
using a Newton-Rapshon scheme
2725 * solve_nonlinear_timestep(solution_delta);
2729 * Add the computed solution increment to total solution
2732 * solution_n += solution_delta;
2736 * Store the converged values of the
internal variables
2739 * update_end_timestep();
2746 *
if (( (time.get_timestep()%parameters.timestep_output) == 0 )
2747 * && (parameters.outfiles_requested ==
"true") )
2749 * output_results_to_vtu(time.get_timestep(),
2750 * time.get_current(),
2754 * output_results_to_plot(time.get_timestep(),
2755 * time.get_current(),
2762 * Increment the time step (=load step)
2765 * time.increment_time();
2770 * Print the footers and close files
2773 *
if (this_mpi_process == 0)
2775 * print_plot_file_footer(pointfile);
2776 * pointfile.close ();
2777 * print_console_file_footer(outfile);
2781 * NOTE: ideally, we should close the outfile here [ >> outfile.close (); ]
2782 * But
if we
do, then the timer output will not be printed. That is why we leave it open.
2791 * <a name=
"nonlinear-poro-viscoelasticity.cc-Privateinterface"></a>
2792 * <h4>Private interface</h4>
2793 * We define the structures needed
for parallelization with Threading Building Blocks (TBB)
2794 * Tangent
matrix and right-hand side force vector assembly structures.
2795 * PerTaskData_ASM stores local contributions
2798 *
template <
int dim>
2799 *
struct Solid<dim>::PerTaskData_ASM
2803 * std::vector<types::global_dof_index> local_dof_indices;
2805 * PerTaskData_ASM(
const unsigned int dofs_per_cell)
2808 * cell_rhs(dofs_per_cell),
2809 * local_dof_indices(dofs_per_cell)
2821 * ScratchData_ASM stores larger objects used during the assembly
2824 *
template <
int dim>
2825 *
template <
typename NumberType>
2826 *
struct Solid<dim>::ScratchData_ASM
2832 * Integration helper
2843 * std::vector<NumberType> local_dof_values;
2844 * std::vector<Tensor<2, dim, NumberType> > solution_grads_u_total;
2845 * std::vector<NumberType> solution_values_p_fluid_total;
2846 * std::vector<Tensor<1, dim, NumberType> > solution_grads_p_fluid_total;
2847 * std::vector<Tensor<1, dim, NumberType> > solution_grads_face_p_fluid_total;
2851 * shape function values
2854 * std::vector<std::vector<Tensor<1,dim>>> Nx;
2855 * std::vector<std::vector<double>> Nx_p_fluid;
2861 * std::vector<std::vector<Tensor<2,dim, NumberType>>> grad_Nx;
2862 * std::vector<std::vector<SymmetricTensor<2,dim, NumberType>>> symm_grad_Nx;
2863 * std::vector<std::vector<Tensor<1,dim, NumberType>>> grad_Nx_p_fluid;
2870 * solution_total (solution_total),
2871 * fe_values_ref(fe_cell, qf_cell, uf_cell),
2872 * fe_face_values_ref(fe_cell, qf_face, uf_face),
2873 * local_dof_values(fe_cell.dofs_per_cell),
2874 * solution_grads_u_total(qf_cell.size()),
2875 * solution_values_p_fluid_total(qf_cell.size()),
2876 * solution_grads_p_fluid_total(qf_cell.size()),
2877 * solution_grads_face_p_fluid_total(qf_face.size()),
2878 * Nx(qf_cell.size(), std::vector<
Tensor<1,dim>>(fe_cell.dofs_per_cell)),
2879 * Nx_p_fluid(qf_cell.size(), std::vector<double>(fe_cell.dofs_per_cell)),
2885 * ScratchData_ASM(
const ScratchData_ASM &rhs)
2887 * solution_total (rhs.solution_total),
2888 * fe_values_ref(rhs.fe_values_ref.get_fe(),
2889 * rhs.fe_values_ref.get_quadrature(),
2890 * rhs.fe_values_ref.get_update_flags()),
2891 * fe_face_values_ref(rhs.fe_face_values_ref.get_fe(),
2892 * rhs.fe_face_values_ref.get_quadrature(),
2893 * rhs.fe_face_values_ref.get_update_flags()),
2894 * local_dof_values(rhs.local_dof_values),
2895 * solution_grads_u_total(rhs.solution_grads_u_total),
2896 * solution_values_p_fluid_total(rhs.solution_values_p_fluid_total),
2897 * solution_grads_p_fluid_total(rhs.solution_grads_p_fluid_total),
2898 * solution_grads_face_p_fluid_total(rhs.solution_grads_face_p_fluid_total),
2900 * Nx_p_fluid(rhs.Nx_p_fluid),
2901 * grad_Nx(rhs.grad_Nx),
2902 * symm_grad_Nx(rhs.symm_grad_Nx),
2903 * grad_Nx_p_fluid(rhs.grad_Nx_p_fluid)
2908 *
const unsigned int n_q_points = Nx_p_fluid.size();
2909 *
const unsigned int n_dofs_per_cell = Nx_p_fluid[0].size();
2911 *
Assert(local_dof_values.size() == n_dofs_per_cell, ExcInternalError());
2913 *
for (
unsigned int k = 0; k < n_dofs_per_cell; ++k)
2915 * local_dof_values[k] = 0.0;
2918 *
Assert(solution_grads_u_total.size() == n_q_points, ExcInternalError());
2919 *
Assert(solution_values_p_fluid_total.size() == n_q_points, ExcInternalError());
2920 *
Assert(solution_grads_p_fluid_total.size() == n_q_points, ExcInternalError());
2922 *
Assert(Nx.size() == n_q_points, ExcInternalError());
2923 *
Assert(grad_Nx.size() == n_q_points, ExcInternalError());
2924 *
Assert(symm_grad_Nx.size() == n_q_points, ExcInternalError());
2926 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2928 *
Assert( Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2929 *
Assert( grad_Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2930 *
Assert( symm_grad_Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
2932 * solution_grads_u_total[q_point] = 0.0;
2933 * solution_values_p_fluid_total[q_point] = 0.0;
2934 * solution_grads_p_fluid_total[q_point] = 0.0;
2936 *
for (
unsigned int k = 0; k < n_dofs_per_cell; ++k)
2938 * Nx[q_point][k] = 0.0;
2939 * Nx_p_fluid[q_point][k] = 0.0;
2940 * grad_Nx[q_point][k] = 0.0;
2941 * symm_grad_Nx[q_point][k] = 0.0;
2942 * grad_Nx_p_fluid[q_point][k] = 0.0;
2946 *
const unsigned int n_f_q_points = solution_grads_face_p_fluid_total.size();
2947 *
Assert(solution_grads_face_p_fluid_total.size() == n_f_q_points, ExcInternalError());
2949 *
for (
unsigned int f_q_point = 0; f_q_point < n_f_q_points; ++f_q_point)
2950 * solution_grads_face_p_fluid_total[f_q_point] = 0.0;
2956 * Define the boundary conditions on the mesh
2959 *
template <
int dim>
2960 *
void Solid<dim>::make_constraints(
const int &it_nr_IN)
2962 * pcout <<
" CST " << std::flush;
2963 * outfile <<
" CST " << std::flush;
2965 *
if (it_nr_IN > 1)
return;
2967 *
const bool apply_dirichlet_bc = (it_nr_IN == 0);
2969 *
if (apply_dirichlet_bc)
2971 * constraints.clear();
2972 * make_dirichlet_constraints(constraints);
2976 *
for (
unsigned int i=0; i<dof_handler_ref.n_dofs(); ++i)
2977 *
if (constraints.is_inhomogeneously_constrained(i) ==
true)
2978 * constraints.set_inhomogeneity(i,0.0);
2980 * constraints.close();
2985 * Set-up the FE system
2988 *
template <
int dim>
2991 * timerconsole.enter_subsection(
"Setup system");
2992 * timerfile.enter_subsection(
"Setup system");
2996 * Determine number of components per block
2999 * std::vector<unsigned int> block_component(n_components, u_block);
3000 * block_component[p_fluid_component] = p_fluid_block;
3004 * The DOF handler is initialised and we renumber the grid in an efficient manner.
3007 * dof_handler_ref.distribute_dofs(fe);
3013 * Count the number of DoFs in each block
3020 * Setup the sparsity pattern and tangent
matrix
3024 * std::vector<IndexSet> all_locally_relevant_dofs
3027 * locally_owned_dofs.clear();
3028 * locally_owned_partitioning.clear();
3029 *
Assert(all_locally_owned_dofs.size() > this_mpi_process, ExcInternalError());
3032 * locally_relevant_dofs.clear();
3033 * locally_relevant_partitioning.clear();
3034 *
Assert(all_locally_relevant_dofs.size() > this_mpi_process, ExcInternalError());
3037 * locally_owned_partitioning.reserve(n_blocks);
3038 * locally_relevant_partitioning.reserve(n_blocks);
3043 * = std::accumulate(dofs_per_block.begin(),
3044 * std::next(dofs_per_block.begin(),b), 0);
3046 * = std::accumulate(dofs_per_block.begin(),
3047 * std::next(dofs_per_block.begin(),b+1), 0);
3048 * locally_owned_partitioning.push_back(locally_owned_dofs.get_view(idx_begin, idx_end));
3049 * locally_relevant_partitioning.push_back(locally_relevant_dofs.get_view(idx_begin, idx_end));
3054 * Print information on screen
3057 * pcout <<
"\nTriangulation:\n"
3058 * <<
" Number of active cells: "
3060 * <<
" (by partition:";
3062 * pcout << (p==0 ?
' ' :
'+')
3066 * pcout <<
" Number of degrees of freedom: "
3067 * << dof_handler_ref.n_dofs()
3068 * <<
" (by partition:";
3070 * pcout << (p==0 ?
' ' :
'+')
3074 * pcout <<
" Number of degrees of freedom per block: "
3075 * <<
"[n_u, n_p_fluid] = ["
3076 * << dofs_per_block[u_block]
3078 * << dofs_per_block[p_fluid_block]
3084 * Print information to file
3087 * outfile <<
"\nTriangulation:\n"
3088 * <<
" Number of active cells: "
3090 * <<
" (by partition:";
3092 * outfile << (p==0 ?
' ' :
'+')
3096 * outfile <<
" Number of degrees of freedom: "
3097 * << dof_handler_ref.n_dofs()
3098 * <<
" (by partition:";
3100 * outfile << (p==0 ?
' ' :
'+')
3104 * outfile <<
" Number of degrees of freedom per block: "
3105 * <<
"[n_u, n_p_fluid] = ["
3106 * << dofs_per_block[u_block]
3108 * << dofs_per_block[p_fluid_block]
3114 * We optimise the sparsity pattern to reflect
this structure and prevent
3115 * unnecessary data creation
for the right-
diagonal block components.
3119 *
for (
unsigned int ii = 0; ii < n_components; ++ii)
3120 *
for (
unsigned int jj = 0; jj < n_components; ++jj)
3124 * Identify
"zero" matrix components of FE-system (The two components
do not couple)
3127 *
if (((ii == p_fluid_component) && (jj < p_fluid_component))
3128 * || ((ii < p_fluid_component) && (jj == p_fluid_component)) )
3133 * The rest of components
always couple
3140 * mpi_communicator);
3143 *
false, this_mpi_process);
3148 * Reinitialize the (sparse) tangent
matrix with the given sparsity pattern.
3151 * tangent_matrix.reinit (bsp);
3155 * Initialize the right hand side and solution vectors with number of DoFs
3158 * system_rhs.reinit(locally_owned_partitioning, mpi_communicator);
3159 * solution_n.reinit(locally_owned_partitioning, mpi_communicator);
3160 * solution_delta_OUT.reinit(locally_owned_partitioning, mpi_communicator);
3168 * mpi_communicator);
3170 *
false, this_mpi_process);
3172 * tangent_matrix_nb.reinit (sp);
3173 * system_rhs_nb.reinit(locally_owned_dofs, mpi_communicator);
3177 * Set up the quadrature
point history
3182 * timerconsole.leave_subsection();
3183 * timerfile.leave_subsection();
3188 * Component extractors: used to
extract sub-blocks from the global
matrix
3189 * Description of which local element DOFs are attached to which block component
3192 *
template <
int dim>
3193 *
void Solid<dim>::determine_component_extractors()
3195 * element_indices_u.clear();
3196 * element_indices_p_fluid.clear();
3198 *
for (
unsigned int k = 0; k < fe.dofs_per_cell; ++k)
3200 *
const unsigned int k_group = fe.system_to_base_index(k).first.first;
3201 *
if (k_group == u_block)
3202 * element_indices_u.push_back(k);
3203 *
else if (k_group == p_fluid_block)
3204 * element_indices_p_fluid.push_back(k);
3207 *
Assert(k_group <= p_fluid_block, ExcInternalError());
3214 * Set-up quadrature
point history (QPH) data objects
3217 *
template <
int dim>
3218 *
void Solid<dim>::setup_qph()
3220 * pcout <<
"\nSetting up quadrature point data..." << std::endl;
3221 * outfile <<
"\nSetting up quadrature point data..." << std::endl;
3225 * Create QPH data objects.
3233 * Setup the
initial quadrature
point data
using the info stored in parameters
3238 * dof_handler_ref.begin_active()),
3240 * dof_handler_ref.end());
3241 *
for (; cell!=endc; ++cell)
3243 *
Assert(cell->is_locally_owned(), ExcInternalError());
3244 *
Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
3246 *
const std::vector<std::shared_ptr<PointHistory<dim, ADNumberType> > >
3247 * lqph = quadrature_point_history.get_data(cell);
3248 *
Assert(lqph.size() == n_q_points, ExcInternalError());
3250 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3251 * lqph[q_point]->setup_lqp(parameters, time);
3257 * Solve the non-linear system
using a Newton-Raphson scheme
3260 *
template <
int dim>
3265 * Print the load step
3268 * pcout << std::endl
3270 * << time.get_timestep()
3272 * << time.get_current()
3275 * outfile << std::endl
3277 * << time.get_timestep()
3279 * << time.get_current()
3285 * Declare newton_update vector (solution of a Newton iteration),
3286 * which must have as many positions as global DoFs.
3290 * (locally_owned_partitioning, mpi_communicator);
3294 * Reset the error storage objects
3297 * error_residual.reset();
3298 * error_residual_0.reset();
3299 * error_residual_norm.reset();
3300 * error_update.reset();
3301 * error_update_0.reset();
3302 * error_update_norm.reset();
3304 * print_conv_header();
3308 * Declare and initialize iterator
for the Newton-Raphson algorithm steps
3311 *
unsigned int newton_iteration = 0;
3315 * Iterate until error is below tolerance or
max number iterations are reached
3318 *
while(newton_iteration < parameters.max_iterations_NR)
3320 * pcout <<
" " << std::setw(2) << newton_iteration <<
" " << std::flush;
3321 * outfile <<
" " << std::setw(2) << newton_iteration <<
" " << std::flush;
3325 * Initialize global stiffness
matrix and global force vector to zero
3328 * tangent_matrix = 0.0;
3331 * tangent_matrix_nb = 0.0;
3332 * system_rhs_nb = 0.0;
3336 * Apply boundary conditions
3339 * make_constraints(newton_iteration);
3340 * assemble_system(solution_delta_OUT);
3344 * Compute the rhs residual (error between external and
internal forces in FE system)
3347 * get_error_residual(error_residual);
3351 * error_residual in
first iteration is stored to normalize posterior error measures
3354 *
if (newton_iteration == 0)
3355 * error_residual_0 = error_residual;
3359 * Determine the normalised residual error
3362 * error_residual_norm = error_residual;
3363 * error_residual_norm.normalise(error_residual_0);
3367 * If both errors are below the tolerances, exit the
loop.
3368 * We need to
check the residual vector directly
for convergence
3369 * in the load steps where no external forces or displacements are imposed.
3372 *
if ( ((newton_iteration > 0)
3373 * && (error_update_norm.u <= parameters.tol_u)
3374 * && (error_update_norm.p_fluid <= parameters.tol_p_fluid)
3375 * && (error_residual_norm.u <= parameters.tol_f)
3376 * && (error_residual_norm.p_fluid <= parameters.tol_f))
3377 * || ( (newton_iteration > 0)
3378 * && system_rhs.l2_norm() <= parameters.tol_f) )
3380 * pcout <<
"\n ***** CONVERGED! ***** "
3381 * << system_rhs.l2_norm() <<
" "
3382 * <<
" " << error_residual_norm.norm
3383 * <<
" " << error_residual_norm.u
3384 * <<
" " << error_residual_norm.p_fluid
3385 * <<
" " << error_update_norm.norm
3386 * <<
" " << error_update_norm.u
3387 * <<
" " << error_update_norm.p_fluid
3388 * <<
" " << std::endl;
3389 * outfile <<
"\n ***** CONVERGED! ***** "
3390 * << system_rhs.l2_norm() <<
" "
3391 * <<
" " << error_residual_norm.norm
3392 * <<
" " << error_residual_norm.u
3393 * <<
" " << error_residual_norm.p_fluid
3394 * <<
" " << error_update_norm.norm
3395 * <<
" " << error_update_norm.u
3396 * <<
" " << error_update_norm.p_fluid
3397 * <<
" " << std::endl;
3398 * print_conv_footer();
3405 * Solve the linearized system
3408 * solve_linear_system(newton_update);
3409 * constraints.distribute(newton_update);
3413 * Compute the displacement error
3416 * get_error_update(newton_update, error_update);
3420 * error_update in
first iteration is stored to normalize posterior error measures
3423 *
if (newton_iteration == 0)
3424 * error_update_0 = error_update;
3428 * Determine the normalised Newton update error
3431 * error_update_norm = error_update;
3432 * error_update_norm.normalise(error_update_0);
3436 * Determine the normalised residual error
3439 * error_residual_norm = error_residual;
3440 * error_residual_norm.normalise(error_residual_0);
3444 * Print error values
3447 * pcout <<
" | " << std::fixed << std::setprecision(3)
3448 * << std::setw(7) << std::scientific
3449 * << system_rhs.l2_norm()
3450 * <<
" " << error_residual_norm.
norm
3451 * <<
" " << error_residual_norm.u
3452 * <<
" " << error_residual_norm.p_fluid
3453 * <<
" " << error_update_norm.norm
3454 * <<
" " << error_update_norm.u
3455 * <<
" " << error_update_norm.p_fluid
3456 * <<
" " << std::endl;
3458 * outfile <<
" | " << std::fixed << std::setprecision(3)
3459 * << std::setw(7) << std::scientific
3460 * << system_rhs.l2_norm()
3461 * <<
" " << error_residual_norm.norm
3462 * <<
" " << error_residual_norm.u
3463 * <<
" " << error_residual_norm.p_fluid
3464 * <<
" " << error_update_norm.norm
3465 * <<
" " << error_update_norm.u
3466 * <<
" " << error_update_norm.p_fluid
3467 * <<
" " << std::endl;
3474 * solution_delta_OUT += newton_update;
3475 * newton_update = 0.0;
3476 * newton_iteration++;
3481 * If maximum allowed number of iterations
for Newton algorithm are reached, print non-convergence message and
abort program
3484 *
AssertThrow (newton_iteration < parameters.max_iterations_NR, ExcMessage(
"No convergence in nonlinear solver!"));
3489 * Prints the header
for convergence info on console
3492 *
template <
int dim>
3493 *
void Solid<dim>::print_conv_header()
3495 *
static const unsigned int l_width = 120;
3497 *
for (
unsigned int i = 0; i < l_width; ++i)
3503 * pcout << std::endl;
3504 * outfile << std::endl;
3506 * pcout <<
"\n SOLVER STEP | SYS_RES "
3507 * <<
"RES_NORM RES_U RES_P "
3508 * <<
"NU_NORM NU_U NU_P " << std::endl;
3509 * outfile <<
"\n SOLVER STEP | SYS_RES "
3510 * <<
"RES_NORM RES_U RES_P "
3511 * <<
"NU_NORM NU_U NU_P " << std::endl;
3513 *
for (
unsigned int i = 0; i < l_width; ++i)
3518 * pcout << std::endl << std::endl;
3519 * outfile << std::endl << std::endl;
3524 * Prints the footer
for convergence info on console
3527 *
template <
int dim>
3528 *
void Solid<dim>::print_conv_footer()
3530 *
static const unsigned int l_width = 120;
3532 *
for (
unsigned int i = 0; i < l_width; ++i)
3537 * pcout << std::endl << std::endl;
3538 * outfile << std::endl << std::endl;
3540 * pcout <<
"Relative errors:" << std::endl
3541 * <<
"Displacement: "
3542 * << error_update.u / error_update_0.u << std::endl
3543 * <<
"Force (displ): "
3544 * << error_residual.u / error_residual_0.u << std::endl
3545 * <<
"Pore pressure: "
3546 * << error_update.p_fluid / error_update_0.p_fluid << std::endl
3547 * <<
"Force (pore): "
3548 * << error_residual.p_fluid / error_residual_0.p_fluid << std::endl;
3549 * outfile <<
"Relative errors:" << std::endl
3550 * <<
"Displacement: "
3551 * << error_update.u / error_update_0.u << std::endl
3552 * <<
"Force (displ): "
3553 * << error_residual.u / error_residual_0.u << std::endl
3554 * <<
"Pore pressure: "
3555 * << error_update.p_fluid / error_update_0.p_fluid << std::endl
3556 * <<
"Force (pore): "
3557 * << error_residual.p_fluid / error_residual_0.p_fluid << std::endl;
3562 * Determine the
true residual error
for the problem
3565 *
template <
int dim>
3566 *
void Solid<dim>::get_error_residual(Errors &error_residual_OUT)
3569 * constraints.set_zero(error_res);
3571 * error_residual_OUT.norm = error_res.l2_norm();
3572 * error_residual_OUT.u = error_res.block(u_block).l2_norm();
3573 * error_residual_OUT.p_fluid = error_res.block(p_fluid_block).l2_norm();
3578 * Determine the
true Newton update error
for the problem
3581 *
template <
int dim>
3582 *
void Solid<dim>::get_error_update
3584 * Errors &error_update_OUT)
3587 * constraints.set_zero(error_ud);
3589 * error_update_OUT.norm = error_ud.l2_norm();
3590 * error_update_OUT.u = error_ud.block(u_block).l2_norm();
3591 * error_update_OUT.p_fluid = error_ud.block(p_fluid_block).l2_norm();
3596 * Compute the total solution, which is
valid at any Newton step. This is required as, to
reduce
3597 * computational error, the total solution is only updated at the
end of the timestep.
3600 *
template <
int dim>
3606 * Cell interpolation -> Ghosted vector
3610 * solution_total (locally_owned_partitioning,
3611 * locally_relevant_partitioning,
3615 * solution_total = solution_n;
3616 * tmp = solution_delta_IN;
3617 * solution_total += tmp;
3618 *
return solution_total;
3623 * Compute elemental stiffness tensor and right-hand side force vector, and
assemble into global ones
3626 *
template <
int dim>
3629 * timerconsole.enter_subsection(
"Assemble system");
3630 * timerfile.enter_subsection(
"Assemble system");
3631 * pcout <<
" ASM_SYS " << std::flush;
3632 * outfile <<
" ASM_SYS " << std::flush;
3638 * Info given to
FEValues and
FEFaceValues constructors, to indicate which data will be needed at each element.
3652 * Setup a
copy of the data structures required
for the process and pass them, along with the
3656 * PerTaskData_ASM per_task_data(dofs_per_cell);
3657 * ScratchData_ASM<ADNumberType> scratch_data(fe, qf_cell, uf_cell,
3663 * dof_handler_ref.begin_active()),
3665 * dof_handler_ref.end());
3666 *
for (; cell != endc; ++cell)
3668 *
Assert(cell->is_locally_owned(), ExcInternalError());
3669 *
Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
3671 * assemble_system_one_cell(cell, scratch_data, per_task_data);
3672 * copy_local_to_global_system(per_task_data);
3680 * timerconsole.leave_subsection();
3681 * timerfile.leave_subsection();
3686 * Add the local elemental contribution to the global stiffness tensor
3687 * We
do it twice,
for the block and the non-block systems
3690 *
template <
int dim>
3691 *
void Solid<dim>::copy_local_to_global_system (
const PerTaskData_ASM &data)
3693 * constraints.distribute_local_to_global(data.cell_matrix,
3695 * data.local_dof_indices,
3699 * constraints.distribute_local_to_global(data.cell_matrix,
3701 * data.local_dof_indices,
3702 * tangent_matrix_nb,
3708 * Compute stiffness
matrix and corresponding rhs
for one element
3711 *
template <
int dim>
3712 *
void Solid<dim>::assemble_system_one_cell
3714 * ScratchData_ASM<ADNumberType> &scratch,
3715 * PerTaskData_ASM &data)
const
3717 *
Assert(cell->is_locally_owned(), ExcInternalError());
3721 * scratch.fe_values_ref.reinit(cell);
3722 * cell->get_dof_indices(data.local_dof_indices);
3726 * Setup automatic differentiation
3729 *
for (
unsigned int k = 0; k < dofs_per_cell; ++k)
3733 * Initialise the dofs
for the cell
using the current solution.
3736 * scratch.local_dof_values[k] = scratch.solution_total[data.local_dof_indices[k]];
3742 * scratch.local_dof_values[k].diff(k, dofs_per_cell);
3747 * Update the quadrature
point solution
3748 * Compute the values and
gradients of the solution in terms of the AD variables
3751 *
for (
unsigned int q = 0; q < n_q_points; ++q)
3753 *
for (
unsigned int k = 0; k < dofs_per_cell; ++k)
3755 *
const unsigned int k_group = fe.system_to_base_index(k).first.first;
3756 *
if (k_group == u_block)
3759 * scratch.fe_values_ref[u_fe].gradient(k, q);
3760 *
for (
unsigned int dd = 0; dd < dim; ++dd)
3762 *
for (
unsigned int ee = 0; ee < dim; ++ee)
3764 * scratch.solution_grads_u_total[q][dd][ee]
3765 * += scratch.local_dof_values[k] * Grad_Nx_u[dd][ee];
3769 *
else if (k_group == p_fluid_block)
3771 *
const double Nx_p = scratch.fe_values_ref[p_fluid_fe].value(k, q);
3773 * scratch.fe_values_ref[p_fluid_fe].gradient(k, q);
3775 * scratch.solution_values_p_fluid_total[q]
3776 * += scratch.local_dof_values[k] * Nx_p;
3777 *
for (
unsigned int dd = 0; dd < dim; ++dd)
3779 * scratch.solution_grads_p_fluid_total[q][dd]
3780 * += scratch.local_dof_values[k] * Grad_Nx_p[dd];
3784 *
Assert(k_group <= p_fluid_block, ExcInternalError());
3790 * Set up pointer
"lgph" to the PointHistory
object of
this element
3793 *
const std::vector<std::shared_ptr<const PointHistory<dim, ADNumberType> > >
3794 * lqph = quadrature_point_history.get_data(cell);
3795 *
Assert(lqph.size() == n_q_points, ExcInternalError());
3800 * Precalculate the element shape function values and
gradients
3803 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3810 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3812 *
const unsigned int i_group = fe.system_to_base_index(i).first.first;
3814 *
if (i_group == u_block)
3816 * scratch.Nx[q_point][i] =
3817 * scratch.fe_values_ref[u_fe].value(i, q_point);
3818 * scratch.grad_Nx[q_point][i] =
3819 * scratch.fe_values_ref[u_fe].gradient(i, q_point)*F_inv_AD;
3820 * scratch.symm_grad_Nx[q_point][i] =
3823 *
else if (i_group == p_fluid_block)
3825 * scratch.Nx_p_fluid[q_point][i] =
3826 * scratch.fe_values_ref[p_fluid_fe].value(i, q_point);
3827 * scratch.grad_Nx_p_fluid[q_point][i] =
3828 * scratch.fe_values_ref[p_fluid_fe].gradient(i, q_point)*F_inv_AD;
3831 *
Assert(i_group <= p_fluid_block, ExcInternalError());
3837 * Assemble the stiffness
matrix and rhs vector
3840 * std::vector<ADNumberType> residual_ad (dofs_per_cell, ADNumberType(0.0));
3841 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
3845 *
const ADNumberType det_F_AD =
determinant(F_AD);
3847 *
Assert(det_F_AD > 0, ExcInternalError());
3850 *
const ADNumberType p_fluid = scratch.solution_values_p_fluid_total[q_point];
3853 * PointHistory<dim, ADNumberType> *lqph_q_point_nc =
3854 *
const_cast<PointHistory<dim, ADNumberType>*
>(lqph[q_point].get());
3855 * lqph_q_point_nc->update_internal_equilibrium(F_AD);
3860 * Get some info from constitutive model of solid
3866 * tau_E = lqph[q_point]->get_tau_E(F_AD);
3868 * tau_fluid_vol *= -1.0 * p_fluid * det_F_AD;
3872 * Get some info from constitutive model of fluid
3875 *
const ADNumberType det_F_aux = lqph[q_point]->get_converged_det_F();
3878 * = lqph[q_point]->get_overall_body_force(F_AD, parameters);
3882 * Define some aliases to make the assembly process easier to follow
3885 *
const std::vector<Tensor<1,dim>> &Nu = scratch.Nx[q_point];
3886 *
const std::vector<SymmetricTensor<2, dim, ADNumberType>>
3887 * &symm_grad_Nu = scratch.symm_grad_Nx[q_point];
3888 *
const std::vector<double> &Np = scratch.Nx_p_fluid[q_point];
3889 *
const std::vector<Tensor<1, dim, ADNumberType> > &grad_Np
3890 * = scratch.grad_Nx_p_fluid[q_point];
3892 * = scratch.solution_grads_p_fluid_total[q_point]*F_inv_AD;
3893 *
const double JxW = scratch.fe_values_ref.JxW(q_point);
3895 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3897 *
const unsigned int i_group = fe.system_to_base_index(i).first.first;
3899 *
if (i_group == u_block)
3901 * residual_ad[i] += symm_grad_Nu[i] * ( tau_E + tau_fluid_vol ) * JxW;
3902 * residual_ad[i] -= Nu[i] * overall_body_force * JxW;
3904 *
else if (i_group == p_fluid_block)
3907 * = lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p);
3908 * residual_ad[i] += Np[i] * (det_F_AD - det_F_converged) * JxW;
3909 * residual_ad[i] -= time.get_delta_t() * grad_Np[i]
3910 * * seepage_vel_current * JxW;
3913 *
Assert(i_group <= p_fluid_block, ExcInternalError());
3919 * Assemble the Neumann contribution (external force contribution).
3922 *
for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
3924 *
if (cell->face(face)->at_boundary() ==
true)
3926 * scratch.fe_face_values_ref.reinit(cell, face);
3928 *
for (
unsigned int f_q_point = 0; f_q_point < n_q_points_f; ++f_q_point)
3931 * = scratch.fe_face_values_ref.normal_vector(f_q_point);
3933 * = scratch.fe_face_values_ref.quadrature_point(f_q_point);
3935 * = get_neumann_traction(cell->face(face)->boundary_id(), pt, N);
3937 * = get_prescribed_fluid_flow(cell->face(face)->boundary_id(), pt);
3939 *
if ( (traction.norm() < 1e-12) && (
std::abs(flow) < 1e-12) )
continue;
3941 *
const double JxW_f = scratch.fe_face_values_ref.JxW(f_q_point);
3943 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3945 *
const unsigned int i_group = fe.system_to_base_index(i).first.first;
3947 *
if ((i_group == u_block) && (traction.norm() > 1e-12))
3949 *
const unsigned int component_i
3950 * = fe.system_to_component_index(i).first;
3952 * = scratch.fe_face_values_ref.shape_value(i, f_q_point);
3953 * residual_ad[i] -= (Nu_f * traction[component_i]) * JxW_f;
3955 *
if ((i_group == p_fluid_block) && (
std::abs(flow) > 1
e-12))
3958 * = scratch.fe_face_values_ref.shape_value(i, f_q_point);
3959 * residual_ad[i] -= (Nu_p * flow) * JxW_f;
3968 * Linearise the residual
3971 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3973 *
const ADNumberType &R_i = residual_ad[i];
3975 * data.cell_rhs(i) -= R_i.val();
3976 *
for (
unsigned int j=0; j<dofs_per_cell; ++j)
3977 * data.cell_matrix(i,j) += R_i.fastAccessDx(j);
3983 * Store the converged values of the
internal variables
3986 *
template <
int dim>
3987 *
void Solid<dim>::update_end_timestep()
3991 * dof_handler_ref.begin_active()),
3993 * dof_handler_ref.end());
3994 *
for (; cell!=endc; ++cell)
3996 *
Assert(cell->is_locally_owned(), ExcInternalError());
3997 *
Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
3999 *
const std::vector<std::shared_ptr<PointHistory<dim, ADNumberType> > >
4000 * lqph = quadrature_point_history.get_data(cell);
4001 *
Assert(lqph.size() == n_q_points, ExcInternalError());
4002 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
4003 * lqph[q_point]->update_end_timestep();
4010 * Solve the linearized equations
4013 *
template <
int dim>
4017 * timerconsole.enter_subsection(
"Linear solver");
4018 * timerfile.enter_subsection(
"Linear solver");
4019 * pcout <<
" SLV " << std::flush;
4020 * outfile <<
" SLV " << std::flush;
4023 * newton_update_nb.
reinit(locally_owned_dofs, mpi_communicator);
4026 * 1.0e-6 * system_rhs_nb.l2_norm());
4028 * solver.solve(tangent_matrix_nb, newton_update_nb, system_rhs_nb);
4032 * Copy the non-block solution back to block system
4035 *
for (
unsigned int i=0; i<locally_owned_dofs.n_elements(); ++i)
4038 * = locally_owned_dofs.nth_index_in_set(i);
4039 * newton_update_OUT(idx_i) = newton_update_nb(idx_i);
4043 * timerconsole.leave_subsection();
4044 * timerfile.leave_subsection();
4049 * Class to compute
gradient of the pressure
4052 *
template <
int dim>
4056 * GradientPostprocessor (
const unsigned int p_fluid_component)
4060 * p_fluid_component (p_fluid_component)
4063 *
virtual ~GradientPostprocessor(){}
4066 * evaluate_vector_field
4068 * std::vector<
Vector<double> > &computed_quantities)
const override
4071 * computed_quantities.size());
4072 *
for (
unsigned int p=0; p<input_data.solution_gradients.size(); ++p)
4075 *
for (
unsigned int d=0;
d<dim; ++
d)
4076 * computed_quantities[p][d]
4077 * = input_data.solution_gradients[p][p_fluid_component][d];
4082 *
const unsigned int p_fluid_component;
4088 * Print results to
vtu file
4091 *
template <
int dim>
void Solid<dim>::output_results_to_vtu
4092 * (
const unsigned int timestep,
4093 *
const double current_time,
4097 * locally_relevant_partitioning,
4100 * solution_total = solution_IN;
4104 * GradientPostprocessor<dim> gradient_postprocessor(p_fluid_component);
4108 * Declare local variables with number of stress components
4112 *
unsigned int num_comp_symm_tensor = 6;
4116 * Declare local vectors to store values
4117 * OUTPUT AVERAGED ON ELEMENTS -------------------------------------------
4120 * std::vector<Vector<double>>cauchy_stresses_total_elements
4121 * (num_comp_symm_tensor,
4123 * std::vector<Vector<double>>cauchy_stresses_E_elements
4124 * (num_comp_symm_tensor,
4126 * std::vector<Vector<double>>stretches_elements
4129 * std::vector<Vector<double>>seepage_velocity_elements
4141 * OUTPUT AVERAGED ON NODES ----------------------------------------------
4142 * We need to create a
new FE space with a single dof per node to avoid
4143 * duplication of the output on nodes
for our problem with dim+1 dofs.
4148 * vertex_handler_ref.distribute_dofs(fe_vertex);
4150 * ExcDimensionMismatch(vertex_handler_ref.n_dofs(),
4154 * (vertex_handler_ref.n_dofs());
4156 * (vertex_handler_ref.n_dofs());
4158 * std::vector<Vector<double>>cauchy_stresses_total_vertex_mpi
4159 * (num_comp_symm_tensor,
4161 * std::vector<Vector<double>>sum_cauchy_stresses_total_vertex
4162 * (num_comp_symm_tensor,
4164 * std::vector<Vector<double>>cauchy_stresses_E_vertex_mpi
4165 * (num_comp_symm_tensor,
4167 * std::vector<Vector<double>>sum_cauchy_stresses_E_vertex
4168 * (num_comp_symm_tensor,
4170 * std::vector<Vector<double>>stretches_vertex_mpi
4173 * std::vector<Vector<double>>sum_stretches_vertex
4176 *
Vector<double> porous_dissipation_vertex_mpi(vertex_handler_ref.n_dofs());
4177 *
Vector<double> sum_porous_dissipation_vertex(vertex_handler_ref.n_dofs());
4178 *
Vector<double> viscous_dissipation_vertex_mpi(vertex_handler_ref.n_dofs());
4179 *
Vector<double> sum_viscous_dissipation_vertex(vertex_handler_ref.n_dofs());
4180 *
Vector<double> solid_vol_fraction_vertex_mpi(vertex_handler_ref.n_dofs());
4181 *
Vector<double> sum_solid_vol_fraction_vertex(vertex_handler_ref.n_dofs());
4185 * We need to create a
new FE space with a dim dof per node to
4186 * be able to output data on nodes in vector form
4191 * vertex_vec_handler_ref.distribute_dofs(fe_vertex_vec);
4193 * ExcDimensionMismatch(vertex_vec_handler_ref.n_dofs(),
4196 *
Vector<double> seepage_velocity_vertex_vec_mpi(vertex_vec_handler_ref.n_dofs());
4197 *
Vector<double> sum_seepage_velocity_vertex_vec(vertex_vec_handler_ref.n_dofs());
4198 *
Vector<double> counter_on_vertices_vec_mpi(vertex_vec_handler_ref.n_dofs());
4199 *
Vector<double> sum_counter_on_vertices_vec(vertex_vec_handler_ref.n_dofs());
4202 * -----------------------------------------------------------------------
4206 * Declare and initialize local unit vectors (to construct tensor basis)
4209 * std::vector<Tensor<1,dim>> basis_vectors (dim,
Tensor<1,dim>() );
4210 *
for (
unsigned int i=0; i<dim; ++i)
4211 * basis_vectors[i][i] = 1;
4215 * Declare an instance of the material
class object
4218 *
if (parameters.mat_type ==
"Neo-Hooke")
4219 * NeoHooke<dim,ADNumberType> material(parameters,time);
4220 *
else if (parameters.mat_type ==
"Ogden")
4221 * Ogden<dim,ADNumberType> material(parameters,time);
4222 *
else if (parameters.mat_type ==
"visco-Ogden")
4223 * visco_Ogden <dim,ADNumberType>material(parameters,time);
4225 *
Assert (
false, ExcMessage(
"Material type not implemented"));
4229 * Define a local instance of
FEValues to compute updated values required
4230 * to calculate stresses
4239 * Iterate through elements (cells) and Gauss Points
4244 * dof_handler_ref.begin_active()),
4246 * dof_handler_ref.end()),
4248 * vertex_handler_ref.begin_active()),
4250 * vertex_vec_handler_ref.begin_active());
4256 *
for (; cell!=endc; ++cell, ++cell_v, ++cell_v_vec)
4258 *
Assert(cell->is_locally_owned(), ExcInternalError());
4259 *
Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
4262 *
static_cast<int>(cell->material_id());
4264 * fe_values_ref.reinit(cell);
4266 * std::vector<Tensor<2,dim>> solution_grads_u(n_q_points);
4267 * fe_values_ref[u_fe].get_function_gradients(solution_total,
4268 * solution_grads_u);
4270 * std::vector<double> solution_values_p_fluid_total(n_q_points);
4271 * fe_values_ref[p_fluid_fe].get_function_values(solution_total,
4272 * solution_values_p_fluid_total);
4274 * std::vector<Tensor<1,dim>> solution_grads_p_fluid_AD (n_q_points);
4275 * fe_values_ref[p_fluid_fe].get_function_gradients(solution_total,
4276 * solution_grads_p_fluid_AD);
4283 *
for (
unsigned int q_point=0; q_point<n_q_points; ++q_point)
4290 *
const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
4291 * lqph = quadrature_point_history.get_data(cell);
4292 *
Assert(lqph.size() == n_q_points, ExcInternalError());
4294 *
const double p_fluid = solution_values_p_fluid_total[q_point];
4305 * lqph[q_point]->get_Cauchy_E(F_AD);
4307 *
for (
unsigned int i=0; i<dim; ++i)
4308 *
for (
unsigned int j=0; j<dim; ++j)
4312 * sigma_fluid_vol *= -p_fluid;
4320 *
const double solid_vol_fraction = (parameters.solid_vol_frac)/det_F;
4324 * Green-Lagrange strain
4336 * solution_grads_p_fluid_AD[q_point]*F_inv;
4338 * lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p_fluid_AD);
4345 *
const double porous_dissipation =
4346 * lqph[q_point]->get_porous_dissipation(F_AD, grad_p_fluid_AD);
4347 *
const double viscous_dissipation =
4348 * lqph[q_point]->get_viscous_dissipation();
4352 * OUTPUT AVERAGED ON ELEMENTS -------------------------------------------
4353 * Both average on elements and on nodes is NOT weighted with the
4355 * integration
point to the average. Ideally, it should be weighted,
4356 * but I haven
't invested time in getting it to work properly.
4359 * if (parameters.outtype == "elements")
4361 * for (unsigned int j=0; j<dim; ++j)
4363 * cauchy_stresses_total_elements[j](cell->active_cell_index())
4364 * += ((sigma*basis_vectors[j])*basis_vectors[j])/n_q_points;
4365 * cauchy_stresses_E_elements[j](cell->active_cell_index())
4366 * += ((sigma_E*basis_vectors[j])*basis_vectors[j])/n_q_points;
4367 * stretches_elements[j](cell->active_cell_index())
4368 * += std::sqrt(1.0+2.0*Tensor<0,dim,double>(E_strain[j][j]))
4370 * seepage_velocity_elements[j](cell->active_cell_index())
4371 * += Tensor<0,dim,double>(seepage_vel_AD[j])/n_q_points;
4374 * porous_dissipation_elements(cell->active_cell_index())
4375 * += porous_dissipation/n_q_points;
4376 * viscous_dissipation_elements(cell->active_cell_index())
4377 * += viscous_dissipation/n_q_points;
4378 * solid_vol_fraction_elements(cell->active_cell_index())
4379 * += solid_vol_fraction/n_q_points;
4381 * cauchy_stresses_total_elements[3](cell->active_cell_index())
4382 * += ((sigma*basis_vectors[0])*basis_vectors[1])/n_q_points; //sig_xy
4383 * cauchy_stresses_total_elements[4](cell->active_cell_index())
4384 * += ((sigma*basis_vectors[0])*basis_vectors[2])/n_q_points;//sig_xz
4385 * cauchy_stresses_total_elements[5](cell->active_cell_index())
4386 * += ((sigma*basis_vectors[1])*basis_vectors[2])/n_q_points;//sig_yz
4388 * cauchy_stresses_E_elements[3](cell->active_cell_index())
4389 * += ((sigma_E*basis_vectors[0])* basis_vectors[1])/n_q_points; //sig_xy
4390 * cauchy_stresses_E_elements[4](cell->active_cell_index())
4391 * += ((sigma_E*basis_vectors[0])* basis_vectors[2])/n_q_points;//sig_xz
4392 * cauchy_stresses_E_elements[5](cell->active_cell_index())
4393 * += ((sigma_E*basis_vectors[1])* basis_vectors[2])/n_q_points;//sig_yz
4398 * OUTPUT AVERAGED ON NODES -------------------------------------------
4401 * else if (parameters.outtype == "nodes")
4403 * for (unsigned int v=0; v<(GeometryInfo<dim>::vertices_per_cell); ++v)
4405 * types::global_dof_index local_vertex_indices =
4406 * cell_v->vertex_dof_index(v, 0);
4407 * counter_on_vertices_mpi(local_vertex_indices) += 1;
4408 * for (unsigned int k=0; k<dim; ++k)
4410 * cauchy_stresses_total_vertex_mpi[k](local_vertex_indices)
4411 * += (sigma*basis_vectors[k])*basis_vectors[k];
4412 * cauchy_stresses_E_vertex_mpi[k](local_vertex_indices)
4413 * += (sigma_E*basis_vectors[k])*basis_vectors[k];
4414 * stretches_vertex_mpi[k](local_vertex_indices)
4415 * += std::sqrt(1.0+2.0*Tensor<0,dim,double>(E_strain[k][k]));
4417 * types::global_dof_index local_vertex_vec_indices =
4418 * cell_v_vec->vertex_dof_index(v, k);
4419 * counter_on_vertices_vec_mpi(local_vertex_vec_indices) += 1;
4420 * seepage_velocity_vertex_vec_mpi(local_vertex_vec_indices)
4421 * += Tensor<0,dim,double>(seepage_vel_AD[k]);
4424 * porous_dissipation_vertex_mpi(local_vertex_indices)
4425 * += porous_dissipation;
4426 * viscous_dissipation_vertex_mpi(local_vertex_indices)
4427 * += viscous_dissipation;
4428 * solid_vol_fraction_vertex_mpi(local_vertex_indices)
4429 * += solid_vol_fraction;
4431 * cauchy_stresses_total_vertex_mpi[3](local_vertex_indices)
4432 * += (sigma*basis_vectors[0])*basis_vectors[1]; //sig_xy
4433 * cauchy_stresses_total_vertex_mpi[4](local_vertex_indices)
4434 * += (sigma*basis_vectors[0])*basis_vectors[2];//sig_xz
4435 * cauchy_stresses_total_vertex_mpi[5](local_vertex_indices)
4436 * += (sigma*basis_vectors[1])*basis_vectors[2]; //sig_yz
4438 * cauchy_stresses_E_vertex_mpi[3](local_vertex_indices)
4439 * += (sigma_E*basis_vectors[0])*basis_vectors[1]; //sig_xy
4440 * cauchy_stresses_E_vertex_mpi[4](local_vertex_indices)
4441 * += (sigma_E*basis_vectors[0])*basis_vectors[2];//sig_xz
4442 * cauchy_stresses_E_vertex_mpi[5](local_vertex_indices)
4443 * += (sigma_E*basis_vectors[1])*basis_vectors[2]; //sig_yz
4448 * ---------------------------------------------------------------
4451 * } //end gauss point loop
4456 * Different nodes might have different amount of contributions, e.g.,
4457 * corner nodes have less integration points contributing to the averaged.
4458 * This is why we need a counter and divide at the end, outside the cell loop.
4461 * if (parameters.outtype == "nodes")
4463 * for (unsigned int d=0; d<(vertex_handler_ref.n_dofs()); ++d)
4465 * sum_counter_on_vertices[d] =
4466 * Utilities::MPI::sum(counter_on_vertices_mpi[d],
4467 * mpi_communicator);
4468 * sum_porous_dissipation_vertex[d] =
4469 * Utilities::MPI::sum(porous_dissipation_vertex_mpi[d],
4470 * mpi_communicator);
4471 * sum_viscous_dissipation_vertex[d] =
4472 * Utilities::MPI::sum(viscous_dissipation_vertex_mpi[d],
4473 * mpi_communicator);
4474 * sum_solid_vol_fraction_vertex[d] =
4475 * Utilities::MPI::sum(solid_vol_fraction_vertex_mpi[d],
4476 * mpi_communicator);
4478 * for (unsigned int k=0; k<num_comp_symm_tensor; ++k)
4480 * sum_cauchy_stresses_total_vertex[k][d] =
4481 * Utilities::MPI::sum(cauchy_stresses_total_vertex_mpi[k][d],
4482 * mpi_communicator);
4483 * sum_cauchy_stresses_E_vertex[k][d] =
4484 * Utilities::MPI::sum(cauchy_stresses_E_vertex_mpi[k][d],
4485 * mpi_communicator);
4487 * for (unsigned int k=0; k<dim; ++k)
4489 * sum_stretches_vertex[k][d] =
4490 * Utilities::MPI::sum(stretches_vertex_mpi[k][d],
4491 * mpi_communicator);
4495 * for (unsigned int d=0; d<(vertex_vec_handler_ref.n_dofs()); ++d)
4497 * sum_counter_on_vertices_vec[d] =
4498 * Utilities::MPI::sum(counter_on_vertices_vec_mpi[d],
4499 * mpi_communicator);
4500 * sum_seepage_velocity_vertex_vec[d] =
4501 * Utilities::MPI::sum(seepage_velocity_vertex_vec_mpi[d],
4502 * mpi_communicator);
4505 * for (unsigned int d=0; d<(vertex_handler_ref.n_dofs()); ++d)
4507 * if (sum_counter_on_vertices[d]>0)
4509 * for (unsigned int i=0; i<num_comp_symm_tensor; ++i)
4511 * sum_cauchy_stresses_total_vertex[i][d] /= sum_counter_on_vertices[d];
4512 * sum_cauchy_stresses_E_vertex[i][d] /= sum_counter_on_vertices[d];
4514 * for (unsigned int i=0; i<dim; ++i)
4516 * sum_stretches_vertex[i][d] /= sum_counter_on_vertices[d];
4518 * sum_porous_dissipation_vertex[d] /= sum_counter_on_vertices[d];
4519 * sum_viscous_dissipation_vertex[d] /= sum_counter_on_vertices[d];
4520 * sum_solid_vol_fraction_vertex[d] /= sum_counter_on_vertices[d];
4524 * for (unsigned int d=0; d<(vertex_vec_handler_ref.n_dofs()); ++d)
4526 * if (sum_counter_on_vertices_vec[d]>0)
4528 * sum_seepage_velocity_vertex_vec[d] /= sum_counter_on_vertices_vec[d];
4536 * Add the results to the solution to create the output file for Paraview
4539 * DataOut<dim> data_out;
4540 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
4542 * DataComponentInterpretation::component_is_part_of_vector);
4543 * comp_type.push_back(DataComponentInterpretation::component_is_scalar);
4545 * GridTools::get_subdomain_association(triangulation, partition_int);
4547 * std::vector<std::string> solution_name(dim, "displacement");
4548 * solution_name.push_back("pore_pressure");
4550 * data_out.attach_dof_handler(dof_handler_ref);
4551 * data_out.add_data_vector(solution_total,
4553 * DataOut<dim>::type_dof_data,
4556 * data_out.add_data_vector(solution_total,
4557 * gradient_postprocessor);
4559 * const Vector<double> partitioning(partition_int.begin(),
4560 * partition_int.end());
4562 * data_out.add_data_vector(partitioning, "partitioning");
4563 * data_out.add_data_vector(material_id, "material_id");
4567 * Integration point results -----------------------------------------------------------
4570 * if (parameters.outtype == "elements")
4572 * data_out.add_data_vector(cauchy_stresses_total_elements[0], "cauchy_xx");
4573 * data_out.add_data_vector(cauchy_stresses_total_elements[1], "cauchy_yy");
4574 * data_out.add_data_vector(cauchy_stresses_total_elements[2], "cauchy_zz");
4575 * data_out.add_data_vector(cauchy_stresses_total_elements[3], "cauchy_xy");
4576 * data_out.add_data_vector(cauchy_stresses_total_elements[4], "cauchy_xz");
4577 * data_out.add_data_vector(cauchy_stresses_total_elements[5], "cauchy_yz");
4579 * data_out.add_data_vector(cauchy_stresses_E_elements[0], "cauchy_E_xx");
4580 * data_out.add_data_vector(cauchy_stresses_E_elements[1], "cauchy_E_yy");
4581 * data_out.add_data_vector(cauchy_stresses_E_elements[2], "cauchy_E_zz");
4582 * data_out.add_data_vector(cauchy_stresses_E_elements[3], "cauchy_E_xy");
4583 * data_out.add_data_vector(cauchy_stresses_E_elements[4], "cauchy_E_xz");
4584 * data_out.add_data_vector(cauchy_stresses_E_elements[5], "cauchy_E_yz");
4586 * data_out.add_data_vector(stretches_elements[0], "stretch_xx");
4587 * data_out.add_data_vector(stretches_elements[1], "stretch_yy");
4588 * data_out.add_data_vector(stretches_elements[2], "stretch_zz");
4590 * data_out.add_data_vector(seepage_velocity_elements[0], "seepage_vel_x");
4591 * data_out.add_data_vector(seepage_velocity_elements[1], "seepage_vel_y");
4592 * data_out.add_data_vector(seepage_velocity_elements[2], "seepage_vel_z");
4594 * data_out.add_data_vector(porous_dissipation_elements, "dissipation_porous");
4595 * data_out.add_data_vector(viscous_dissipation_elements, "dissipation_viscous");
4596 * data_out.add_data_vector(solid_vol_fraction_elements, "solid_vol_fraction");
4598 * else if (parameters.outtype == "nodes")
4600 * data_out.add_data_vector(vertex_handler_ref,
4601 * sum_cauchy_stresses_total_vertex[0],
4603 * data_out.add_data_vector(vertex_handler_ref,
4604 * sum_cauchy_stresses_total_vertex[1],
4606 * data_out.add_data_vector(vertex_handler_ref,
4607 * sum_cauchy_stresses_total_vertex[2],
4609 * data_out.add_data_vector(vertex_handler_ref,
4610 * sum_cauchy_stresses_total_vertex[3],
4612 * data_out.add_data_vector(vertex_handler_ref,
4613 * sum_cauchy_stresses_total_vertex[4],
4615 * data_out.add_data_vector(vertex_handler_ref,
4616 * sum_cauchy_stresses_total_vertex[5],
4619 * data_out.add_data_vector(vertex_handler_ref,
4620 * sum_cauchy_stresses_E_vertex[0],
4622 * data_out.add_data_vector(vertex_handler_ref,
4623 * sum_cauchy_stresses_E_vertex[1],
4625 * data_out.add_data_vector(vertex_handler_ref,
4626 * sum_cauchy_stresses_E_vertex[2],
4628 * data_out.add_data_vector(vertex_handler_ref,
4629 * sum_cauchy_stresses_E_vertex[3],
4631 * data_out.add_data_vector(vertex_handler_ref,
4632 * sum_cauchy_stresses_E_vertex[4],
4634 * data_out.add_data_vector(vertex_handler_ref,
4635 * sum_cauchy_stresses_E_vertex[5],
4638 * data_out.add_data_vector(vertex_handler_ref,
4639 * sum_stretches_vertex[0],
4641 * data_out.add_data_vector(vertex_handler_ref,
4642 * sum_stretches_vertex[1],
4644 * data_out.add_data_vector(vertex_handler_ref,
4645 * sum_stretches_vertex[2],
4648 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
4649 * comp_type_vec(dim,
4650 * DataComponentInterpretation::component_is_part_of_vector);
4651 * std::vector<std::string> solution_name_vec(dim,"seepage_velocity");
4653 * data_out.add_data_vector(vertex_vec_handler_ref,
4654 * sum_seepage_velocity_vertex_vec,
4655 * solution_name_vec,
4658 * data_out.add_data_vector(vertex_handler_ref,
4659 * sum_porous_dissipation_vertex,
4660 * "dissipation_porous");
4661 * data_out.add_data_vector(vertex_handler_ref,
4662 * sum_viscous_dissipation_vertex,
4663 * "dissipation_viscous");
4664 * data_out.add_data_vector(vertex_handler_ref,
4665 * sum_solid_vol_fraction_vertex,
4666 * "solid_vol_fraction");
4670 * ---------------------------------------------------------------------
4676 * data_out.build_patches(degree_displ);
4680 * static std::string get_filename_vtu(unsigned int process,
4681 * unsigned int timestep,
4682 * const unsigned int n_digits = 5)
4684 * std::ostringstream filename_vtu;
4687 * << Utilities::int_to_string(process, n_digits)
4689 * << Utilities::int_to_string(timestep, n_digits)
4691 * return filename_vtu.str();
4694 * static std::string get_filename_pvtu(unsigned int timestep,
4695 * const unsigned int n_digits = 5)
4697 * std::ostringstream filename_vtu;
4700 * << Utilities::int_to_string(timestep, n_digits)
4702 * return filename_vtu.str();
4705 * static std::string get_filename_pvd (void)
4707 * std::ostringstream filename_vtu;
4709 * << "solution.pvd";
4710 * return filename_vtu.str();
4714 * const std::string filename_vtu = Filename::get_filename_vtu(this_mpi_process,
4716 * std::ofstream output(filename_vtu.c_str());
4717 * data_out.write_vtu(output);
4721 * We have a collection of files written in parallel
4722 * This next set of steps should only be performed by master process
4725 * if (this_mpi_process == 0)
4729 * List of all files written out at this timestep by all processors
4732 * std::vector<std::string> parallel_filenames_vtu;
4733 * for (unsigned int p=0; p<n_mpi_processes; ++p)
4735 * parallel_filenames_vtu.push_back(Filename::get_filename_vtu(p, timestep));
4738 * const std::string filename_pvtu(Filename::get_filename_pvtu(timestep));
4739 * std::ofstream pvtu_master(filename_pvtu.c_str());
4740 * data_out.write_pvtu_record(pvtu_master,
4741 * parallel_filenames_vtu);
4745 * Time dependent data master file
4748 * static std::vector<std::pair<double,std::string>> time_and_name_history;
4749 * time_and_name_history.push_back(std::make_pair(current_time,
4751 * const std::string filename_pvd(Filename::get_filename_pvd());
4752 * std::ofstream pvd_output(filename_pvd.c_str());
4753 * DataOutBase::write_pvd_record(pvd_output, time_and_name_history);
4760 * Print results to plotting file
4763 * template <int dim>
4764 * void Solid<dim>::output_results_to_plot(
4765 * const unsigned int timestep,
4766 * const double current_time,
4767 * TrilinosWrappers::MPI::BlockVector solution_IN,
4768 * std::vector<Point<dim> > &tracked_vertices_IN,
4769 * std::ofstream &plotpointfile) const
4771 * TrilinosWrappers::MPI::BlockVector solution_total(locally_owned_partitioning,
4772 * locally_relevant_partitioning,
4777 * solution_total = solution_IN;
4781 * Variables needed to print the solution file for plotting
4784 * Point<dim> reaction_force;
4785 * Point<dim> reaction_force_pressure;
4786 * Point<dim> reaction_force_extra;
4787 * double total_fluid_flow = 0.0;
4788 * double total_porous_dissipation = 0.0;
4789 * double total_viscous_dissipation = 0.0;
4790 * double total_solid_vol = 0.0;
4791 * double total_vol_current = 0.0;
4792 * double total_vol_reference = 0.0;
4793 * std::vector<Point<dim+1>> solution_vertices(tracked_vertices_IN.size());
4797 * Auxiliary variables needed for mpi processing
4800 * Tensor<1,dim> sum_reaction_mpi;
4801 * Tensor<1,dim> sum_reaction_pressure_mpi;
4802 * Tensor<1,dim> sum_reaction_extra_mpi;
4803 * sum_reaction_mpi = 0.0;
4804 * sum_reaction_pressure_mpi = 0.0;
4805 * sum_reaction_extra_mpi = 0.0;
4806 * double sum_total_flow_mpi = 0.0;
4807 * double sum_porous_dissipation_mpi = 0.0;
4808 * double sum_viscous_dissipation_mpi = 0.0;
4809 * double sum_solid_vol_mpi = 0.0;
4810 * double sum_vol_current_mpi = 0.0;
4811 * double sum_vol_reference_mpi = 0.0;
4815 * Declare an instance of the material class object
4818 * if (parameters.mat_type == "Neo-Hooke")
4819 * NeoHooke<dim,ADNumberType> material(parameters,time);
4820 * else if (parameters.mat_type == "Ogden")
4821 * Ogden<dim,ADNumberType> material(parameters, time);
4822 * else if (parameters.mat_type == "visco-Ogden")
4823 * visco_Ogden <dim,ADNumberType>material(parameters,time);
4825 * Assert (false, ExcMessage("Material type not implemented"));
4829 * Define a local instance of FEValues to compute updated values required
4830 * to calculate stresses
4833 * const UpdateFlags uf_cell(update_values | update_gradients |
4834 * update_JxW_values);
4835 * FEValues<dim> fe_values_ref (fe, qf_cell, uf_cell);
4839 * Iterate through elements (cells) and Gauss Points
4842 * FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
4843 * cell(IteratorFilters::LocallyOwnedCell(),
4844 * dof_handler_ref.begin_active()),
4845 * endc(IteratorFilters::LocallyOwnedCell(),
4846 * dof_handler_ref.end());
4852 * for (; cell!=endc; ++cell)
4854 * Assert(cell->is_locally_owned(), ExcInternalError());
4855 * Assert(cell->subdomain_id() == this_mpi_process, ExcInternalError());
4857 * fe_values_ref.reinit(cell);
4859 * std::vector<Tensor<2,dim>> solution_grads_u(n_q_points);
4860 * fe_values_ref[u_fe].get_function_gradients(solution_total,
4861 * solution_grads_u);
4863 * std::vector<double> solution_values_p_fluid_total(n_q_points);
4864 * fe_values_ref[p_fluid_fe].get_function_values(solution_total,
4865 * solution_values_p_fluid_total);
4867 * std::vector<Tensor<1,dim >> solution_grads_p_fluid_AD(n_q_points);
4868 * fe_values_ref[p_fluid_fe].get_function_gradients(solution_total,
4869 * solution_grads_p_fluid_AD);
4873 * start gauss point loop
4876 * for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
4878 * const Tensor<2,dim,ADNumberType>
4879 * F_AD = Physics::Elasticity::Kinematics::F(solution_grads_u[q_point]);
4880 * ADNumberType det_F_AD = determinant(F_AD);
4881 * const double det_F = Tensor<0,dim,double>(det_F_AD);
4883 * const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
4884 * lqph = quadrature_point_history.get_data(cell);
4885 * Assert(lqph.size() == n_q_points, ExcInternalError());
4887 * double JxW = fe_values_ref.JxW(q_point);
4894 * sum_vol_current_mpi += det_F * JxW;
4895 * sum_vol_reference_mpi += JxW;
4896 * sum_solid_vol_mpi += parameters.solid_vol_frac * JxW * det_F;
4903 * const Tensor<2,dim,ADNumberType> F_inv = invert(F_AD);
4904 * const Tensor<1,dim,ADNumberType>
4905 * grad_p_fluid_AD = solution_grads_p_fluid_AD[q_point]*F_inv;
4906 * const Tensor<1,dim,ADNumberType> seepage_vel_AD
4907 * = lqph[q_point]->get_seepage_velocity_current(F_AD, grad_p_fluid_AD);
4914 * const double porous_dissipation =
4915 * lqph[q_point]->get_porous_dissipation(F_AD, grad_p_fluid_AD);
4916 * sum_porous_dissipation_mpi += porous_dissipation * det_F * JxW;
4918 * const double viscous_dissipation = lqph[q_point]->get_viscous_dissipation();
4919 * sum_viscous_dissipation_mpi += viscous_dissipation * det_F * JxW;
4923 * ---------------------------------------------------------------
4926 * } //end gauss point loop
4930 * Compute reaction force on load boundary & total fluid flow across
4932 * Define a local instance of FEFaceValues to compute values required
4933 * to calculate reaction force
4936 * const UpdateFlags uf_face( update_values | update_gradients |
4937 * update_normal_vectors | update_JxW_values );
4938 * FEFaceValues<dim> fe_face_values_ref(fe, qf_face, uf_face);
4945 * for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
4952 * if (cell->face(face)->at_boundary() == true &&
4953 * cell->face(face)->boundary_id() == get_reaction_boundary_id_for_output() )
4955 * fe_face_values_ref.reinit(cell, face);
4959 * Get displacement gradients for current face
4962 * std::vector<Tensor<2,dim> > solution_grads_u_f(n_q_points_f);
4963 * fe_face_values_ref[u_fe].get_function_gradients
4965 * solution_grads_u_f);
4969 * Get pressure for current element
4972 * std::vector< double > solution_values_p_fluid_total_f(n_q_points_f);
4973 * fe_face_values_ref[p_fluid_fe].get_function_values
4975 * solution_values_p_fluid_total_f);
4979 * start gauss points on faces loop
4982 * for (unsigned int f_q_point=0; f_q_point<n_q_points_f; ++f_q_point)
4984 * const Tensor<1,dim> &N = fe_face_values_ref.normal_vector(f_q_point);
4985 * const double JxW_f = fe_face_values_ref.JxW(f_q_point);
4989 * Compute deformation gradient from displacements gradient
4990 * (present configuration)
4993 * const Tensor<2,dim,ADNumberType> F_AD =
4994 * Physics::Elasticity::Kinematics::F(solution_grads_u_f[f_q_point]);
4996 * const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
4997 * lqph = quadrature_point_history.get_data(cell);
4998 * Assert(lqph.size() == n_q_points, ExcInternalError());
5000 * const double p_fluid = solution_values_p_fluid_total[f_q_point];
5007 * static const SymmetricTensor<2,dim,double>
5008 * I (Physics::Elasticity::StandardTensors<dim>::I);
5009 * SymmetricTensor<2,dim> sigma_E;
5010 * const SymmetricTensor<2,dim,ADNumberType> sigma_E_AD =
5011 * lqph[f_q_point]->get_Cauchy_E(F_AD);
5013 * for (unsigned int i=0; i<dim; ++i)
5014 * for (unsigned int j=0; j<dim; ++j)
5015 * sigma_E[i][j] = Tensor<0,dim,double>(sigma_E_AD[i][j]);
5017 * SymmetricTensor<2,dim> sigma_fluid_vol(I);
5018 * sigma_fluid_vol *= -1.0*p_fluid;
5019 * const SymmetricTensor<2,dim> sigma = sigma_E+sigma_fluid_vol;
5020 * sum_reaction_mpi += sigma * N * JxW_f;
5021 * sum_reaction_pressure_mpi += sigma_fluid_vol * N * JxW_f;
5022 * sum_reaction_extra_mpi += sigma_E * N * JxW_f;
5023 * }//end gauss points on faces loop
5031 * if (cell->face(face)->at_boundary() == true &&
5032 * (cell->face(face)->boundary_id() ==
5033 * get_drained_boundary_id_for_output().first ||
5034 * cell->face(face)->boundary_id() ==
5035 * get_drained_boundary_id_for_output().second ) )
5037 * fe_face_values_ref.reinit(cell, face);
5041 * Get displacement gradients for current face
5044 * std::vector<Tensor<2,dim>> solution_grads_u_f(n_q_points_f);
5045 * fe_face_values_ref[u_fe].get_function_gradients
5047 * solution_grads_u_f);
5051 * Get pressure gradients for current face
5054 * std::vector<Tensor<1,dim>> solution_grads_p_f(n_q_points_f);
5055 * fe_face_values_ref[p_fluid_fe].get_function_gradients
5057 * solution_grads_p_f);
5061 * start gauss points on faces loop
5064 * for (unsigned int f_q_point=0; f_q_point<n_q_points_f; ++f_q_point)
5066 * const Tensor<1,dim> &N =
5067 * fe_face_values_ref.normal_vector(f_q_point);
5068 * const double JxW_f = fe_face_values_ref.JxW(f_q_point);
5072 * Deformation gradient and inverse from displacements gradient
5073 * (present configuration)
5076 * const Tensor<2,dim,ADNumberType> F_AD
5077 * = Physics::Elasticity::Kinematics::F(solution_grads_u_f[f_q_point]);
5079 * const Tensor<2,dim,ADNumberType> F_inv_AD = invert(F_AD);
5080 * ADNumberType det_F_AD = determinant(F_AD);
5082 * const std::vector<std::shared_ptr<const PointHistory<dim,ADNumberType>>>
5083 * lqph = quadrature_point_history.get_data(cell);
5084 * Assert(lqph.size() == n_q_points, ExcInternalError());
5091 * Tensor<1,dim> seepage;
5092 * double det_F = Tensor<0,dim,double>(det_F_AD);
5093 * const Tensor<1,dim,ADNumberType> grad_p
5094 * = solution_grads_p_f[f_q_point]*F_inv_AD;
5095 * const Tensor<1,dim,ADNumberType> seepage_AD
5096 * = lqph[f_q_point]->get_seepage_velocity_current(F_AD, grad_p);
5098 * for (unsigned int i=0; i<dim; ++i)
5099 * seepage[i] = Tensor<0,dim,double>(seepage_AD[i]);
5101 * sum_total_flow_mpi += (seepage/det_F) * N * JxW_f;
5102 * }//end gauss points on faces loop
5109 * Sum the results from different MPI process and then add to the reaction_force vector
5110 * In theory, the solution on each surface (each cell) only exists in one MPI process
5111 * so, we add all MPI process, one will have the solution and the others will be zero
5114 * for (unsigned int d=0; d<dim; ++d)
5116 * reaction_force[d] = Utilities::MPI::sum(sum_reaction_mpi[d],
5117 * mpi_communicator);
5118 * reaction_force_pressure[d] = Utilities::MPI::sum(sum_reaction_pressure_mpi[d],
5119 * mpi_communicator);
5120 * reaction_force_extra[d] = Utilities::MPI::sum(sum_reaction_extra_mpi[d],
5121 * mpi_communicator);
5126 * Same for total fluid flow, and for porous and viscous dissipations
5129 * total_fluid_flow = Utilities::MPI::sum(sum_total_flow_mpi,
5130 * mpi_communicator);
5131 * total_porous_dissipation = Utilities::MPI::sum(sum_porous_dissipation_mpi,
5132 * mpi_communicator);
5133 * total_viscous_dissipation = Utilities::MPI::sum(sum_viscous_dissipation_mpi,
5134 * mpi_communicator);
5135 * total_solid_vol = Utilities::MPI::sum(sum_solid_vol_mpi,
5136 * mpi_communicator);
5137 * total_vol_current = Utilities::MPI::sum(sum_vol_current_mpi,
5138 * mpi_communicator);
5139 * total_vol_reference = Utilities::MPI::sum(sum_vol_reference_mpi,
5140 * mpi_communicator);
5144 * Extract solution for tracked vectors
5145 * Copying an MPI::BlockVector into MPI::Vector is not possible,
5146 * so we copy each block of MPI::BlockVector into an MPI::Vector
5147 * And then we copy the MPI::Vector into "normal" Vectors
5150 * TrilinosWrappers::MPI::Vector solution_vector_u_MPI(solution_total.block(u_block));
5151 * TrilinosWrappers::MPI::Vector solution_vector_p_MPI(solution_total.block(p_fluid_block));
5152 * Vector<double> solution_u_vector(solution_vector_u_MPI);
5153 * Vector<double> solution_p_vector(solution_vector_p_MPI);
5155 * if (this_mpi_process == 0)
5159 * Append the pressure solution vector to the displacement solution vector,
5160 * creating a single solution vector equivalent to the original BlockVector
5161 * so FEFieldFunction will work with the dof_handler_ref.
5164 * Vector<double> solution_vector(solution_p_vector.size()
5165 * +solution_u_vector.size());
5167 * for (unsigned int d=0; d<(solution_u_vector.size()); ++d)
5168 * solution_vector[d] = solution_u_vector[d];
5170 * for (unsigned int d=0; d<(solution_p_vector.size()); ++d)
5171 * solution_vector[solution_u_vector.size()+d] = solution_p_vector[d];
5173 * Functions::FEFieldFunction<dim,Vector<double>>
5174 * find_solution(dof_handler_ref, solution_vector);
5176 * for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5178 * Vector<double> update(dim+1);
5179 * Point<dim> pt_ref;
5181 * pt_ref[0]= tracked_vertices_IN[p][0];
5182 * pt_ref[1]= tracked_vertices_IN[p][1];
5183 * pt_ref[2]= tracked_vertices_IN[p][2];
5185 * find_solution.vector_value(pt_ref, update);
5187 * for (unsigned int d=0; d<(dim+1); ++d)
5191 * For values close to zero, set to 0.0
5194 * if (abs(update[d])<1.5*parameters.tol_u)
5196 * solution_vertices[p][d] = update[d];
5201 * Write the results to the plotting file.
5202 * Add two blank lines between cycles in the cyclic loading examples so GNUPLOT can detect each cycle as a different block
5205 * if (( (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")||
5206 * (parameters.geom_type == "Budday_cube_tension_compression")||
5207 * (parameters.geom_type == "Budday_cube_shear_fully_fixed") ) &&
5208 * ( (abs(current_time - parameters.end_time/3.) <0.9*parameters.delta_t)||
5209 * (abs(current_time - 2.*parameters.end_time/3.)<0.9*parameters.delta_t) ) &&
5210 * parameters.num_cycle_sets == 1 )
5212 * plotpointfile << std::endl<< std::endl;
5214 * if (( (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")||
5215 * (parameters.geom_type == "Budday_cube_tension_compression")||
5216 * (parameters.geom_type == "Budday_cube_shear_fully_fixed") ) &&
5217 * ( (abs(current_time - parameters.end_time/9.) <0.9*parameters.delta_t)||
5218 * (abs(current_time - 2.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5219 * (abs(current_time - 3.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5220 * (abs(current_time - 5.*parameters.end_time/9.)<0.9*parameters.delta_t)||
5221 * (abs(current_time - 7.*parameters.end_time/9.)<0.9*parameters.delta_t) ) &&
5222 * parameters.num_cycle_sets == 2 )
5224 * plotpointfile << std::endl<< std::endl;
5227 * plotpointfile << std::setprecision(6) << std::scientific;
5228 * plotpointfile << std::setw(16) << current_time << ","
5229 * << std::setw(15) << total_vol_reference << ","
5230 * << std::setw(15) << total_vol_current << ","
5231 * << std::setw(15) << total_solid_vol << ",";
5233 * if (current_time == 0.0)
5235 * for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5237 * for (unsigned int d=0; d<dim; ++d)
5238 * plotpointfile << std::setw(15) << 0.0 << ",";
5240 * plotpointfile << std::setw(15) << parameters.drained_pressure << ",";
5242 * for (unsigned int d=0; d<(3*dim+2); ++d)
5243 * plotpointfile << std::setw(15) << 0.0 << ",";
5245 * plotpointfile << std::setw(15) << 0.0;
5249 * for (unsigned int p=0; p<tracked_vertices_IN.size(); ++p)
5250 * for (unsigned int d=0; d<(dim+1); ++d)
5251 * plotpointfile << std::setw(15) << solution_vertices[p][d]<< ",";
5253 * for (unsigned int d=0; d<dim; ++d)
5254 * plotpointfile << std::setw(15) << reaction_force[d] << ",";
5256 * for (unsigned int d=0; d<dim; ++d)
5257 * plotpointfile << std::setw(15) << reaction_force_pressure[d] << ",";
5259 * for (unsigned int d=0; d<dim; ++d)
5260 * plotpointfile << std::setw(15) << reaction_force_extra[d] << ",";
5262 * plotpointfile << std::setw(15) << total_fluid_flow << ","
5263 * << std::setw(15) << total_porous_dissipation<< ","
5264 * << std::setw(15) << total_viscous_dissipation;
5266 * plotpointfile << std::endl;
5272 * Header for console output file
5275 * template <int dim>
5276 * void Solid<dim>::print_console_file_header(std::ofstream &outputfile) const
5278 * outputfile << "/*-----------------------------------------------------------------------------------------";
5279 * outputfile << "\n\n Poro-viscoelastic formulation to solve nonlinear solid mechanics problems using deal.ii";
5280 * outputfile << "\n\n Problem setup by E Comellas and J-P Pelteret, University of Erlangen-Nuremberg, 2018";
5281 * outputfile << "\n\n/*-----------------------------------------------------------------------------------------";
5282 * outputfile << "\n\nCONSOLE OUTPUT: \n\n";
5287 * Header for plotting output file
5290 * template <int dim>
5291 * void Solid<dim>::print_plot_file_header(std::vector<Point<dim> > &tracked_vertices,
5292 * std::ofstream &plotpointfile) const
5294 * plotpointfile << "#\n# *** Solution history for tracked vertices -- DOF: 0 = Ux, 1 = Uy, 2 = Uz, 3 = P ***"
5297 * for (unsigned int p=0; p<tracked_vertices.size(); ++p)
5299 * plotpointfile << "# Point " << p << " coordinates: ";
5300 * for (unsigned int d=0; d<dim; ++d)
5302 * plotpointfile << tracked_vertices[p][d];
5303 * if (!( (p == tracked_vertices.size()-1) && (d == dim-1) ))
5304 * plotpointfile << ", ";
5306 * plotpointfile << std::endl;
5308 * plotpointfile << "# The reaction force is the integral over the loaded surfaces in the "
5309 * << "undeformed configuration of the Cauchy stress times the normal surface unit vector.\n"
5310 * << "# reac(p) corresponds to the volumetric part of the Cauchy stress due to the pore fluid pressure"
5311 * << " and reac(E) corresponds to the extra part of the Cauchy stress due to the solid contribution."
5313 * << "# The fluid flow is the integral over the drained surfaces in the "
5314 * << "undeformed configuration of the seepage velocity times the normal surface unit vector."
5316 * << "# Column number:"
5320 * unsigned int columns = 24;
5321 * for (unsigned int d=1; d<columns; ++d)
5322 * plotpointfile << std::setw(15)<< d <<",";
5324 * plotpointfile << std::setw(15)<< columns
5327 * << std::right << std::setw(16) << "Time,"
5328 * << std::right << std::setw(16) << "ref vol,"
5329 * << std::right << std::setw(16) << "def vol,"
5330 * << std::right << std::setw(16) << "solid vol,";
5331 * for (unsigned int p=0; p<tracked_vertices.size(); ++p)
5332 * for (unsigned int d=0; d<(dim+1); ++d)
5333 * plotpointfile << std::right<< std::setw(11)
5334 * <<"P" << p << "[" << d << "],";
5336 * for (unsigned int d=0; d<dim; ++d)
5337 * plotpointfile << std::right<< std::setw(13)
5338 * << "reaction [" << d << "],";
5340 * for (unsigned int d=0; d<dim; ++d)
5341 * plotpointfile << std::right<< std::setw(13)
5342 * << "reac(p) [" << d << "],";
5344 * for (unsigned int d=0; d<dim; ++d)
5345 * plotpointfile << std::right<< std::setw(13)
5346 * << "reac(E) [" << d << "],";
5348 * plotpointfile << std::right<< std::setw(16)<< "fluid flow,"
5349 * << std::right<< std::setw(16)<< "porous dissip,"
5350 * << std::right<< std::setw(15)<< "viscous dissip"
5356 * Footer for console output file
5359 * template <int dim>
5360 * void Solid<dim>::print_console_file_footer(std::ofstream &outputfile) const
5364 * Copy "parameters" file at end of output file.
5367 * std::ifstream infile("parameters.prm");
5368 * std::string content = "";
5371 * for(i=0 ; infile.eof()!=true ; i++)
5373 * char aux = infile.get();
5375 * if(aux=='\n
') content += '#
';
5379 * content.erase(content.end()-1);
5382 * outputfile << "\n\n\n\n PARAMETERS FILE USED IN THIS COMPUTATION: \n#"
5389 * Footer for plotting output file
5392 * template <int dim>
5393 * void Solid<dim>::print_plot_file_footer(std::ofstream &plotpointfile) const
5397 * Copy "parameters" file at end of output file.
5400 * std::ifstream infile("parameters.prm");
5401 * std::string content = "";
5404 * for(i=0 ; infile.eof()!=true ; i++)
5406 * char aux = infile.get();
5408 * if(aux=='\n
') content += '#
';
5412 * content.erase(content.end()-1);
5415 * plotpointfile << "#"<< std::endl
5416 * << "#"<< std::endl
5417 * << "# PARAMETERS FILE USED IN THIS COMPUTATION:" << std::endl
5418 * << "#"<< std::endl
5426 * <a name="nonlinear-poro-viscoelasticity.cc-VerificationexamplesfromEhlersandEipper1999"></a>
5427 * <h3>Verification examples from Ehlers and Eipper 1999</h3>
5428 * We group the definition of the geometry, boundary and loading conditions specific to
5429 * the verification examples from Ehlers and Eipper 1999 into specific classes.
5434 * <a name="nonlinear-poro-viscoelasticity.cc-BaseclassTubegeometryandboundaryconditions"></a>
5435 * <h4>Base class: Tube geometry and boundary conditions</h4>
5438 * template <int dim>
5439 * class VerificationEhlers1999TubeBase
5440 * : public Solid<dim>
5443 * VerificationEhlers1999TubeBase (const Parameters::AllParameters ¶meters)
5444 * : Solid<dim> (parameters)
5447 * virtual ~VerificationEhlers1999TubeBase () {}
5450 * virtual void make_grid() override
5452 * GridGenerator::cylinder( this->triangulation,
5456 * const double rot_angle = 3.0*numbers::PI/2.0;
5457 * GridTools::rotate( Point<3>::unit_vector(1), rot_angle, this->triangulation);
5459 * this->triangulation.reset_manifold(0);
5460 * static const CylindricalManifold<dim> manifold_description_3d(2);
5461 * this->triangulation.set_manifold (0, manifold_description_3d);
5462 * GridTools::scale(this->parameters.scale, this->triangulation);
5463 * this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5464 * this->triangulation.reset_manifold(0);
5467 * virtual void define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5469 * tracked_vertices[0][0] = 0.0*this->parameters.scale;
5470 * tracked_vertices[0][1] = 0.0*this->parameters.scale;
5471 * tracked_vertices[0][2] = 0.5*this->parameters.scale;
5473 * tracked_vertices[1][0] = 0.0*this->parameters.scale;
5474 * tracked_vertices[1][1] = 0.0*this->parameters.scale;
5475 * tracked_vertices[1][2] = -0.5*this->parameters.scale;
5478 * virtual void make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5480 * if (this->time.get_timestep() < 2)
5482 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5484 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5486 * (this->fe.component_mask(this->pressure)));
5490 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5492 * Functions::ZeroFunction<dim>(this->n_components),
5494 * (this->fe.component_mask(this->pressure)));
5497 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5499 * Functions::ZeroFunction<dim>(this->n_components),
5501 * (this->fe.component_mask(this->x_displacement)|
5502 * this->fe.component_mask(this->y_displacement) ) );
5504 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5506 * Functions::ZeroFunction<dim>(this->n_components),
5508 * (this->fe.component_mask(this->x_displacement) |
5509 * this->fe.component_mask(this->y_displacement) |
5510 * this->fe.component_mask(this->z_displacement) ));
5514 * get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5515 * const Point<dim> &pt) const override
5518 * (void)boundary_id;
5522 * virtual types::boundary_id
5523 * get_reaction_boundary_id_for_output() const override
5528 * virtual std::pair<types::boundary_id,types::boundary_id>
5529 * get_drained_boundary_id_for_output() const override
5531 * return std::make_pair(2,2);
5534 * virtual std::vector<double>
5535 * get_dirichlet_load(const types::boundary_id &boundary_id,
5536 * const int &direction) const override
5538 * std::vector<double> displ_incr(dim, 0.0);
5539 * (void)boundary_id;
5541 * AssertThrow(false, ExcMessage("Displacement loading not implemented for Ehlers verification examples."));
5543 * return displ_incr;
5550 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassSteploadexample"></a>
5551 * <h4>Derived class: Step load example</h4>
5554 * template <int dim>
5555 * class VerificationEhlers1999StepLoad
5556 * : public VerificationEhlers1999TubeBase<dim>
5559 * VerificationEhlers1999StepLoad (const Parameters::AllParameters ¶meters)
5560 * : VerificationEhlers1999TubeBase<dim> (parameters)
5563 * virtual ~VerificationEhlers1999StepLoad () {}
5566 * virtual Tensor<1,dim>
5567 * get_neumann_traction (const types::boundary_id &boundary_id,
5568 * const Point<dim> &pt,
5569 * const Tensor<1,dim> &N) const override
5571 * if (this->parameters.load_type == "pressure")
5573 * if (boundary_id == 2)
5575 * return this->parameters.load * N;
5581 * return Tensor<1,dim>();
5588 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassLoadincreasingexample"></a>
5589 * <h4>Derived class: Load increasing example</h4>
5592 * template <int dim>
5593 * class VerificationEhlers1999IncreaseLoad
5594 * : public VerificationEhlers1999TubeBase<dim>
5597 * VerificationEhlers1999IncreaseLoad (const Parameters::AllParameters ¶meters)
5598 * : VerificationEhlers1999TubeBase<dim> (parameters)
5601 * virtual ~VerificationEhlers1999IncreaseLoad () {}
5604 * virtual Tensor<1,dim>
5605 * get_neumann_traction (const types::boundary_id &boundary_id,
5606 * const Point<dim> &pt,
5607 * const Tensor<1,dim> &N) const override
5609 * if (this->parameters.load_type == "pressure")
5611 * if (boundary_id == 2)
5613 * const double initial_load = this->parameters.load;
5614 * const double final_load = 20.0*initial_load;
5615 * const double initial_time = this->time.get_delta_t();
5616 * const double final_time = this->time.get_end();
5617 * const double current_time = this->time.get_current();
5618 * const double load = initial_load + (final_load-initial_load)*(current_time-initial_time)/(final_time-initial_time);
5625 * return Tensor<1,dim>();
5632 * <a name="nonlinear-poro-viscoelasticity.cc-ClassConsolidationcube"></a>
5633 * <h4>Class: Consolidation cube</h4>
5636 * template <int dim>
5637 * class VerificationEhlers1999CubeConsolidation
5638 * : public Solid<dim>
5641 * VerificationEhlers1999CubeConsolidation (const Parameters::AllParameters ¶meters)
5642 * : Solid<dim> (parameters)
5645 * virtual ~VerificationEhlers1999CubeConsolidation () {}
5649 * make_grid() override
5651 * GridGenerator::hyper_rectangle(this->triangulation,
5652 * Point<dim>(0.0, 0.0, 0.0),
5653 * Point<dim>(1.0, 1.0, 1.0),
5656 * GridTools::scale(this->parameters.scale, this->triangulation);
5657 * this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5659 * typename Triangulation<dim>::active_cell_iterator cell =
5660 * this->triangulation.begin_active(), endc = this->triangulation.end();
5661 * for (; cell != endc; ++cell)
5663 * for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
5664 * if (cell->face(face)->at_boundary() == true &&
5665 * cell->face(face)->center()[2] == 1.0 * this->parameters.scale)
5667 * if (cell->face(face)->center()[0] < 0.5 * this->parameters.scale &&
5668 * cell->face(face)->center()[1] < 0.5 * this->parameters.scale)
5669 * cell->face(face)->set_boundary_id(100);
5671 * cell->face(face)->set_boundary_id(101);
5677 * define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5679 * tracked_vertices[0][0] = 0.0*this->parameters.scale;
5680 * tracked_vertices[0][1] = 0.0*this->parameters.scale;
5681 * tracked_vertices[0][2] = 1.0*this->parameters.scale;
5683 * tracked_vertices[1][0] = 0.0*this->parameters.scale;
5684 * tracked_vertices[1][1] = 0.0*this->parameters.scale;
5685 * tracked_vertices[1][2] = 0.0*this->parameters.scale;
5689 * make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5691 * if (this->time.get_timestep() < 2)
5693 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5695 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5697 * (this->fe.component_mask(this->pressure)));
5701 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5703 * Functions::ZeroFunction<dim>(this->n_components),
5705 * (this->fe.component_mask(this->pressure)));
5708 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5710 * Functions::ZeroFunction<dim>(this->n_components),
5712 * this->fe.component_mask(this->x_displacement));
5714 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5716 * Functions::ZeroFunction<dim>(this->n_components),
5718 * this->fe.component_mask(this->x_displacement));
5720 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5722 * Functions::ZeroFunction<dim>(this->n_components),
5724 * this->fe.component_mask(this->y_displacement));
5726 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5728 * Functions::ZeroFunction<dim>(this->n_components),
5730 * this->fe.component_mask(this->y_displacement));
5732 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5734 * Functions::ZeroFunction<dim>(this->n_components),
5736 * ( this->fe.component_mask(this->x_displacement) |
5737 * this->fe.component_mask(this->y_displacement) |
5738 * this->fe.component_mask(this->z_displacement) ));
5741 * virtual Tensor<1,dim>
5742 * get_neumann_traction (const types::boundary_id &boundary_id,
5743 * const Point<dim> &pt,
5744 * const Tensor<1,dim> &N) const override
5746 * if (this->parameters.load_type == "pressure")
5748 * if (boundary_id == 100)
5750 * return this->parameters.load * N;
5756 * return Tensor<1,dim>();
5760 * get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5761 * const Point<dim> &pt) const override
5764 * (void)boundary_id;
5768 * virtual types::boundary_id
5769 * get_reaction_boundary_id_for_output() const override
5774 * virtual std::pair<types::boundary_id,types::boundary_id>
5775 * get_drained_boundary_id_for_output() const override
5777 * return std::make_pair(101,101);
5780 * virtual std::vector<double>
5781 * get_dirichlet_load(const types::boundary_id &boundary_id,
5782 * const int &direction) const override
5784 * std::vector<double> displ_incr(dim, 0.0);
5785 * (void)boundary_id;
5787 * AssertThrow(false, ExcMessage("Displacement loading not implemented for Ehlers verification examples."));
5789 * return displ_incr;
5796 * <a name="nonlinear-poro-viscoelasticity.cc-Franceschiniexperiments"></a>
5797 * <h4>Franceschini experiments</h4>
5800 * template <int dim>
5801 * class Franceschini2006Consolidation
5802 * : public Solid<dim>
5805 * Franceschini2006Consolidation (const Parameters::AllParameters ¶meters)
5806 * : Solid<dim> (parameters)
5809 * virtual ~Franceschini2006Consolidation () {}
5812 * virtual void make_grid() override
5814 * const Point<dim-1> mesh_center(0.0, 0.0);
5815 * const double radius = 0.5;
5818 * const double height = 0.27; //8.1 mm for 30 mm radius
5821 * const double height = 0.23; //6.9 mm for 30 mm radius
5822 * Triangulation<dim-1> triangulation_in;
5823 * GridGenerator::hyper_ball( triangulation_in,
5827 * GridGenerator::extrude_triangulation(triangulation_in,
5830 * this->triangulation);
5832 * const CylindricalManifold<dim> cylinder_3d(2);
5833 * const types::manifold_id cylinder_id = 0;
5836 * this->triangulation.set_manifold(cylinder_id, cylinder_3d);
5838 * for (auto cell : this->triangulation.active_cell_iterators())
5840 * for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
5842 * if (cell->face(face)->at_boundary() == true)
5844 * if (cell->face(face)->center()[2] == 0.0)
5845 * cell->face(face)->set_boundary_id(1);
5847 * else if (cell->face(face)->center()[2] == height)
5848 * cell->face(face)->set_boundary_id(2);
5852 * cell->face(face)->set_boundary_id(0);
5853 * cell->face(face)->set_all_manifold_ids(cylinder_id);
5859 * GridTools::scale(this->parameters.scale, this->triangulation);
5860 * this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
5863 * virtual void define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
5865 * tracked_vertices[0][0] = 0.0*this->parameters.scale;
5866 * tracked_vertices[0][1] = 0.0*this->parameters.scale;
5869 * tracked_vertices[0][2] = 0.27*this->parameters.scale;
5872 * tracked_vertices[0][2] = 0.23*this->parameters.scale;
5874 * tracked_vertices[1][0] = 0.0*this->parameters.scale;
5875 * tracked_vertices[1][1] = 0.0*this->parameters.scale;
5876 * tracked_vertices[1][2] = 0.0*this->parameters.scale;
5879 * virtual void make_dirichlet_constraints(AffineConstraints<double> &constraints) override
5881 * if (this->time.get_timestep() < 2)
5883 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5885 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5887 * (this->fe.component_mask(this->pressure)));
5889 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5891 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
5893 * (this->fe.component_mask(this->pressure)));
5897 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5899 * Functions::ZeroFunction<dim>(this->n_components),
5901 * (this->fe.component_mask(this->pressure)));
5903 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
5905 * Functions::ZeroFunction<dim>(this->n_components),
5907 * (this->fe.component_mask(this->pressure)));
5910 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5912 * Functions::ZeroFunction<dim>(this->n_components),
5914 * (this->fe.component_mask(this->x_displacement)|
5915 * this->fe.component_mask(this->y_displacement) ) );
5917 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5919 * Functions::ZeroFunction<dim>(this->n_components),
5921 * (this->fe.component_mask(this->x_displacement) |
5922 * this->fe.component_mask(this->y_displacement) |
5923 * this->fe.component_mask(this->z_displacement) ));
5925 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
5927 * Functions::ZeroFunction<dim>(this->n_components),
5929 * (this->fe.component_mask(this->x_displacement) |
5930 * this->fe.component_mask(this->y_displacement) ));
5934 * get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
5935 * const Point<dim> &pt) const override
5938 * (void)boundary_id;
5942 * virtual types::boundary_id
5943 * get_reaction_boundary_id_for_output() const override
5948 * virtual std::pair<types::boundary_id,types::boundary_id>
5949 * get_drained_boundary_id_for_output() const override
5951 * return std::make_pair(1,2);
5954 * virtual std::vector<double>
5955 * get_dirichlet_load(const types::boundary_id &boundary_id,
5956 * const int &direction) const override
5958 * std::vector<double> displ_incr(dim, 0.0);
5959 * (void)boundary_id;
5961 * AssertThrow(false, ExcMessage("Displacement loading not implemented for Franceschini examples."));
5963 * return displ_incr;
5966 * virtual Tensor<1,dim>
5967 * get_neumann_traction (const types::boundary_id &boundary_id,
5968 * const Point<dim> &pt,
5969 * const Tensor<1,dim> &N) const override
5971 * if (this->parameters.load_type == "pressure")
5973 * if (boundary_id == 2)
5975 * return (this->parameters.load * N);
5977 * const double final_load = this->parameters.load;
5978 * const double final_load_time = 10 * this->time.get_delta_t();
5979 * const double current_time = this->time.get_current();
5982 * const double c = final_load_time / 2.0;
5983 * const double r = 200.0 * 0.03 / c;
5985 * const double load = final_load * std::exp(r * current_time)
5986 * / ( std::exp(c * current_time) + std::exp(r * current_time));
5994 * return Tensor<1,dim>();
6001 * <a name="nonlinear-poro-viscoelasticity.cc-ExamplestoreproduceexperimentsbyBuddayetal2017"></a>
6002 * <h3>Examples to reproduce experiments by Budday et al. 2017</h3>
6003 * We group the definition of the geometry, boundary and loading conditions specific to
6004 * the examples to reproduce experiments by Budday et al. 2017 into specific classes.
6009 * <a name="nonlinear-poro-viscoelasticity.cc-BaseclassCubegeometryandloadingpattern"></a>
6010 * <h4>Base class: Cube geometry and loading pattern</h4>
6013 * template <int dim>
6014 * class BrainBudday2017BaseCube
6015 * : public Solid<dim>
6018 * BrainBudday2017BaseCube (const Parameters::AllParameters ¶meters)
6019 * : Solid<dim> (parameters)
6022 * virtual ~BrainBudday2017BaseCube () {}
6026 * make_grid() override
6028 * GridGenerator::hyper_cube(this->triangulation,
6033 * typename Triangulation<dim>::active_cell_iterator cell =
6034 * this->triangulation.begin_active(), endc = this->triangulation.end();
6035 * for (; cell != endc; ++cell)
6037 * for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
6038 * if (cell->face(face)->at_boundary() == true &&
6039 * ( cell->face(face)->boundary_id() == 0 ||
6040 * cell->face(face)->boundary_id() == 1 ||
6041 * cell->face(face)->boundary_id() == 2 ||
6042 * cell->face(face)->boundary_id() == 3 ) )
6044 * cell->face(face)->set_boundary_id(100);
6048 * GridTools::scale(this->parameters.scale, this->triangulation);
6049 * this->triangulation.refine_global(std::max (1U, this->parameters.global_refinement));
6053 * get_prescribed_fluid_flow (const types::boundary_id &boundary_id,
6054 * const Point<dim> &pt) const override
6057 * (void)boundary_id;
6061 * virtual std::pair<types::boundary_id,types::boundary_id>
6062 * get_drained_boundary_id_for_output() const override
6064 * return std::make_pair(100,100);
6071 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassUniaxialboundaryconditions"></a>
6072 * <h4>Derived class: Uniaxial boundary conditions</h4>
6075 * template <int dim>
6076 * class BrainBudday2017CubeTensionCompression
6077 * : public BrainBudday2017BaseCube<dim>
6080 * BrainBudday2017CubeTensionCompression (const Parameters::AllParameters ¶meters)
6081 * : BrainBudday2017BaseCube<dim> (parameters)
6084 * virtual ~BrainBudday2017CubeTensionCompression () {}
6088 * define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6090 * tracked_vertices[0][0] = 0.5*this->parameters.scale;
6091 * tracked_vertices[0][1] = 0.5*this->parameters.scale;
6092 * tracked_vertices[0][2] = 1.0*this->parameters.scale;
6094 * tracked_vertices[1][0] = 0.5*this->parameters.scale;
6095 * tracked_vertices[1][1] = 0.5*this->parameters.scale;
6096 * tracked_vertices[1][2] = 0.5*this->parameters.scale;
6100 * make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6102 * if (this->time.get_timestep() < 2)
6104 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6106 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6108 * (this->fe.component_mask(this->pressure)));
6112 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6114 * Functions::ZeroFunction<dim>(this->n_components),
6116 * (this->fe.component_mask(this->pressure)));
6118 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6120 * Functions::ZeroFunction<dim>(this->n_components),
6122 * this->fe.component_mask(this->z_displacement) );
6124 * Point<dim> fix_node(0.5*this->parameters.scale, 0.5*this->parameters.scale, 0.0);
6125 * typename DoFHandler<dim>::active_cell_iterator
6126 * cell = this->dof_handler_ref.begin_active(), endc = this->dof_handler_ref.end();
6127 * for (; cell != endc; ++cell)
6128 * for (unsigned int node = 0; node < GeometryInfo<dim>::vertices_per_cell; ++node)
6130 * if ( (abs(cell->vertex(node)[2]-fix_node[2]) < (1e-6 * this->parameters.scale))
6131 * && (abs(cell->vertex(node)[0]-fix_node[0]) < (1e-6 * this->parameters.scale)))
6132 * constraints.add_line(cell->vertex_dof_index(node, 0));
6134 * if ( (abs(cell->vertex(node)[2]-fix_node[2]) < (1e-6 * this->parameters.scale))
6135 * && (abs(cell->vertex(node)[1]-fix_node[1]) < (1e-6 * this->parameters.scale)))
6136 * constraints.add_line(cell->vertex_dof_index(node, 1));
6139 * if (this->parameters.load_type == "displacement")
6141 * const std::vector<double> value = get_dirichlet_load(5,2);
6142 * const FEValuesExtractors::Scalar direction(this->z_displacement);
6144 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6146 * Functions::ConstantFunction<dim>(value[2],this->n_components),
6148 * this->fe.component_mask(direction));
6152 * virtual Tensor<1,dim>
6153 * get_neumann_traction (const types::boundary_id &boundary_id,
6154 * const Point<dim> &pt,
6155 * const Tensor<1,dim> &N) const override
6157 * if (this->parameters.load_type == "pressure")
6159 * if (boundary_id == 5)
6161 * const double final_load = this->parameters.load;
6162 * const double current_time = this->time.get_current();
6163 * const double final_time = this->time.get_end();
6164 * const double num_cycles = 3.0;
6166 * return final_load/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5))) * N;
6172 * return Tensor<1,dim>();
6175 * virtual types::boundary_id
6176 * get_reaction_boundary_id_for_output() const override
6181 * virtual std::vector<double>
6182 * get_dirichlet_load(const types::boundary_id &boundary_id,
6183 * const int &direction) const override
6185 * std::vector<double> displ_incr(dim,0.0);
6187 * if ( (boundary_id == 5) && (direction == 2) )
6189 * const double final_displ = this->parameters.load;
6190 * const double current_time = this->time.get_current();
6191 * const double final_time = this->time.get_end();
6192 * const double delta_time = this->time.get_delta_t();
6193 * const double num_cycles = 3.0;
6194 * double current_displ = 0.0;
6195 * double previous_displ = 0.0;
6197 * if (this->parameters.num_cycle_sets == 1)
6199 * current_displ = final_displ/2.0 * (1.0
6200 * - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5)));
6201 * previous_displ = final_displ/2.0 * (1.0
6202 * - std::sin(numbers::PI * (2.0*num_cycles*(current_time-delta_time)/final_time + 0.5)));
6206 * if ( current_time <= (final_time*1.0/3.0) )
6208 * current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6209 * (2.0*num_cycles*current_time/(final_time*1.0/3.0) + 0.5)));
6210 * previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6211 * (2.0*num_cycles*(current_time-delta_time)/(final_time*1.0/3.0) + 0.5)));
6215 * current_displ = final_displ * (1.0 - std::sin(numbers::PI *
6216 * (2.0*num_cycles*current_time / (final_time*2.0/3.0)
6217 * - (num_cycles - 0.5) )));
6218 * previous_displ = final_displ * (1.0 - std::sin(numbers::PI *
6219 * (2.0*num_cycles*(current_time-delta_time) / (final_time*2.0/3.0)
6220 * - (num_cycles - 0.5))));
6223 * displ_incr[2] = current_displ - previous_displ;
6225 * return displ_incr;
6232 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassNolateraldisplacementinloadingsurfaces"></a>
6233 * <h4>Derived class: No lateral displacement in loading surfaces</h4>
6236 * template <int dim>
6237 * class BrainBudday2017CubeTensionCompressionFullyFixed
6238 * : public BrainBudday2017BaseCube<dim>
6241 * BrainBudday2017CubeTensionCompressionFullyFixed (const Parameters::AllParameters ¶meters)
6242 * : BrainBudday2017BaseCube<dim> (parameters)
6245 * virtual ~BrainBudday2017CubeTensionCompressionFullyFixed () {}
6249 * define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6251 * tracked_vertices[0][0] = 0.5*this->parameters.scale;
6252 * tracked_vertices[0][1] = 0.5*this->parameters.scale;
6253 * tracked_vertices[0][2] = 1.0*this->parameters.scale;
6255 * tracked_vertices[1][0] = 0.5*this->parameters.scale;
6256 * tracked_vertices[1][1] = 0.5*this->parameters.scale;
6257 * tracked_vertices[1][2] = 0.5*this->parameters.scale;
6261 * make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6263 * if (this->time.get_timestep() < 2)
6265 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6267 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6269 * (this->fe.component_mask(this->pressure)));
6273 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6275 * Functions::ZeroFunction<dim>(this->n_components),
6277 * (this->fe.component_mask(this->pressure)));
6280 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6282 * Functions::ZeroFunction<dim>(this->n_components),
6284 * (this->fe.component_mask(this->x_displacement) |
6285 * this->fe.component_mask(this->y_displacement) |
6286 * this->fe.component_mask(this->z_displacement) ));
6289 * if (this->parameters.load_type == "displacement")
6291 * const std::vector<double> value = get_dirichlet_load(5,2);
6292 * const FEValuesExtractors::Scalar direction(this->z_displacement);
6294 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6296 * Functions::ConstantFunction<dim>(value[2],this->n_components),
6298 * this->fe.component_mask(direction) );
6300 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6302 * Functions::ZeroFunction<dim>(this->n_components),
6304 * (this->fe.component_mask(this->x_displacement) |
6305 * this->fe.component_mask(this->y_displacement) ));
6309 * virtual Tensor<1,dim>
6310 * get_neumann_traction (const types::boundary_id &boundary_id,
6311 * const Point<dim> &pt,
6312 * const Tensor<1,dim> &N) const override
6314 * if (this->parameters.load_type == "pressure")
6316 * if (boundary_id == 5)
6318 * const double final_load = this->parameters.load;
6319 * const double current_time = this->time.get_current();
6320 * const double final_time = this->time.get_end();
6321 * const double num_cycles = 3.0;
6323 * return final_load/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5))) * N;
6329 * return Tensor<1,dim>();
6332 * virtual types::boundary_id
6333 * get_reaction_boundary_id_for_output() const override
6338 * virtual std::vector<double>
6339 * get_dirichlet_load(const types::boundary_id &boundary_id,
6340 * const int &direction) const override
6342 * std::vector<double> displ_incr(dim,0.0);
6344 * if ( (boundary_id == 5) && (direction == 2) )
6346 * const double final_displ = this->parameters.load;
6347 * const double current_time = this->time.get_current();
6348 * const double final_time = this->time.get_end();
6349 * const double delta_time = this->time.get_delta_t();
6350 * const double num_cycles = 3.0;
6351 * double current_displ = 0.0;
6352 * double previous_displ = 0.0;
6354 * if (this->parameters.num_cycle_sets == 1)
6356 * current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*current_time/final_time + 0.5)));
6357 * previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI * (2.0*num_cycles*(current_time-delta_time)/final_time + 0.5)));
6361 * if ( current_time <= (final_time*1.0/3.0) )
6363 * current_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6364 * (2.0*num_cycles*current_time/(final_time*1.0/3.0) + 0.5)));
6365 * previous_displ = final_displ/2.0 * (1.0 - std::sin(numbers::PI *
6366 * (2.0*num_cycles*(current_time-delta_time)/(final_time*1.0/3.0) + 0.5)));
6370 * current_displ = final_displ * (1.0 - std::sin(numbers::PI *
6371 * (2.0*num_cycles*current_time / (final_time*2.0/3.0)
6372 * - (num_cycles - 0.5) )));
6373 * previous_displ = final_displ * (1.0 - std::sin(numbers::PI *
6374 * (2.0*num_cycles*(current_time-delta_time) / (final_time*2.0/3.0)
6375 * - (num_cycles - 0.5))));
6378 * displ_incr[2] = current_displ - previous_displ;
6380 * return displ_incr;
6387 * <a name="nonlinear-poro-viscoelasticity.cc-DerivedclassNolateralorverticaldisplacementinloadingsurface"></a>
6388 * <h4>Derived class: No lateral or vertical displacement in loading surface</h4>
6391 * template <int dim>
6392 * class BrainBudday2017CubeShearFullyFixed
6393 * : public BrainBudday2017BaseCube<dim>
6396 * BrainBudday2017CubeShearFullyFixed (const Parameters::AllParameters ¶meters)
6397 * : BrainBudday2017BaseCube<dim> (parameters)
6400 * virtual ~BrainBudday2017CubeShearFullyFixed () {}
6404 * define_tracked_vertices(std::vector<Point<dim> > &tracked_vertices) override
6406 * tracked_vertices[0][0] = 0.75*this->parameters.scale;
6407 * tracked_vertices[0][1] = 0.5*this->parameters.scale;
6408 * tracked_vertices[0][2] = 0.0*this->parameters.scale;
6410 * tracked_vertices[1][0] = 0.25*this->parameters.scale;
6411 * tracked_vertices[1][1] = 0.5*this->parameters.scale;
6412 * tracked_vertices[1][2] = 0.0*this->parameters.scale;
6416 * make_dirichlet_constraints(AffineConstraints<double> &constraints) override
6418 * if (this->time.get_timestep() < 2)
6420 * VectorTools::interpolate_boundary_values(this->dof_handler_ref,
6422 * Functions::ConstantFunction<dim>(this->parameters.drained_pressure,this->n_components),
6424 * (this->fe.component_mask(this->pressure)));
6428 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6430 * Functions::ZeroFunction<dim>(this->n_components),
6432 * (this->fe.component_mask(this->pressure)));
6435 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6437 * Functions::ZeroFunction<dim>(this->n_components),
6439 * (this->fe.component_mask(this->x_displacement) |
6440 * this->fe.component_mask(this->y_displacement) |
6441 * this->fe.component_mask(this->z_displacement) ));
6444 * if (this->parameters.load_type == "displacement")
6446 * const std::vector<double> value = get_dirichlet_load(4,0);
6447 * const FEValuesExtractors::Scalar direction(this->x_displacement);
6449 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6451 * Functions::ConstantFunction<dim>(value[0],this->n_components),
6453 * this->fe.component_mask(direction));
6455 * VectorTools::interpolate_boundary_values( this->dof_handler_ref,
6457 * Functions::ZeroFunction<dim>(this->n_components),
6459 * (this->fe.component_mask(this->y_displacement) |
6460 * this->fe.component_mask(this->z_displacement) ));
6464 * virtual Tensor<1,dim>
6465 * get_neumann_traction (const types::boundary_id &boundary_id,
6466 * const Point<dim> &pt,
6467 * const Tensor<1,dim> &N) const override
6469 * if (this->parameters.load_type == "pressure")
6471 * if (boundary_id == 4)
6473 * const double final_load = this->parameters.load;
6474 * const double current_time = this->time.get_current();
6475 * const double final_time = this->time.get_end();
6476 * const double num_cycles = 3.0;
6477 * const Tensor<1,3> axis ({0.0,1.0,0.0});
6478 * const double angle = numbers::PI;
6479 * static const Tensor< 2, dim, double> R(Physics::Transformations::Rotations::rotation_matrix_3d(axis,angle));
6481 * return (final_load * (std::sin(2.0*(numbers::PI)*num_cycles*current_time/final_time)) * (R * N));
6487 * return Tensor<1,dim>();
6490 * virtual types::boundary_id
6491 * get_reaction_boundary_id_for_output() const override
6496 * virtual std::vector<double>
6497 * get_dirichlet_load(const types::boundary_id &boundary_id,
6498 * const int &direction) const override
6500 * std::vector<double> displ_incr (dim, 0.0);
6502 * if ( (boundary_id == 4) && (direction == 0) )
6504 * const double final_displ = this->parameters.load;
6505 * const double current_time = this->time.get_current();
6506 * const double final_time = this->time.get_end();
6507 * const double delta_time = this->time.get_delta_t();
6508 * const double num_cycles = 3.0;
6509 * double current_displ = 0.0;
6510 * double previous_displ = 0.0;
6512 * if (this->parameters.num_cycle_sets == 1)
6514 * current_displ = final_displ * (std::sin(2.0*(numbers::PI)*num_cycles*current_time/final_time));
6515 * previous_displ = final_displ * (std::sin(2.0*(numbers::PI)*num_cycles*(current_time-delta_time)/final_time));
6519 * AssertThrow(false, ExcMessage("Problem type not defined. Budday shear experiments implemented only for one set of cycles."));
6521 * displ_incr[0] = current_displ - previous_displ;
6523 * return displ_incr;
6532 * <a name="nonlinear-poro-viscoelasticity.cc-Mainfunction"></a>
6533 * <h3>Main function</h3>
6534 * Lastly we provide the main driver function which is similar to the other tutorials.
6537 * int main (int argc, char *argv[])
6539 * using namespace dealii;
6540 * using namespace NonLinearPoroViscoElasticity;
6542 * const unsigned int n_tbb_processes = 1;
6543 * Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, n_tbb_processes);
6547 * Parameters::AllParameters parameters ("parameters.prm");
6548 * if (parameters.geom_type == "Ehlers_tube_step_load")
6550 * VerificationEhlers1999StepLoad<3> solid_3d(parameters);
6553 * else if (parameters.geom_type == "Ehlers_tube_increase_load")
6555 * VerificationEhlers1999IncreaseLoad<3> solid_3d(parameters);
6558 * else if (parameters.geom_type == "Ehlers_cube_consolidation")
6560 * VerificationEhlers1999CubeConsolidation<3> solid_3d(parameters);
6563 * else if (parameters.geom_type == "Franceschini_consolidation")
6565 * Franceschini2006Consolidation<3> solid_3d(parameters);
6568 * else if (parameters.geom_type == "Budday_cube_tension_compression")
6570 * BrainBudday2017CubeTensionCompression<3> solid_3d(parameters);
6573 * else if (parameters.geom_type == "Budday_cube_tension_compression_fully_fixed")
6575 * BrainBudday2017CubeTensionCompressionFullyFixed<3> solid_3d(parameters);
6578 * else if (parameters.geom_type == "Budday_cube_shear_fully_fixed")
6580 * BrainBudday2017CubeShearFullyFixed<3> solid_3d(parameters);
6585 * AssertThrow(false, ExcMessage("Problem type not defined. Current setting: " + parameters.geom_type));
6589 * catch (std::exception &exc)
6591 * if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
6593 * std::cerr << std::endl << std::endl
6594 * << "----------------------------------------------------"
6596 * std::cerr << "Exception on processing: " << std::endl << exc.what()
6597 * << std::endl << "Aborting!" << std::endl
6598 * << "----------------------------------------------------"
6606 * if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
6608 * std::cerr << std::endl << std::endl
6609 * << "----------------------------------------------------"
6611 * std::cerr << "Unknown exception!" << std::endl << "Aborting!"
6613 * << "----------------------------------------------------"
DataPostprocessorVector(const std::string &name, const UpdateFlags update_flags)
numbers::NumberTraits< Number >::real_type norm() const
unsigned int n_active_cells() const
cell_iterator end() const
unsigned int n_vertices() const
active_cell_iterator begin_active(const unsigned int level=0) const
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
__global__ void set(Number *val, const Number s, const size_type N)
#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())
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
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 >())
@ valid
Iterator points to a valid object.
@ matrix
Contents is actually a matrix.
@ diagonal
Matrix is diagonal.
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.)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
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)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
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)
::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 > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)
SymmetricTensorEigenvectorMethod