224 * #include <deal.II/base/mpi.h>
225 * #include <deal.II/base/function.h>
226 * #include <deal.II/base/parameter_handler.h>
227 * #include <deal.II/base/
point.h>
228 * #include <deal.II/base/quadrature_lib.h>
229 * #include <deal.II/base/symmetric_tensor.h>
230 * #include <deal.II/base/tensor.h>
231 * #include <deal.II/base/timer.h>
232 * #include <deal.II/base/work_stream.h>
233 * #include <deal.II/dofs/dof_renumbering.h>
234 * #include <deal.II/dofs/dof_tools.h>
235 * #include <deal.II/base/quadrature_point_data.h>
236 * #include <deal.II/grid/filtered_iterator.h>
237 * #include <deal.II/grid/grid_generator.h>
238 * #include <deal.II/grid/grid_tools.h>
239 * #include <deal.II/grid/grid_in.h>
240 * #include <deal.II/grid/manifold_lib.h>
241 * #include <deal.II/grid/tria.h>
242 * #include <deal.II/fe/fe_dgp_monomial.h>
243 * #include <deal.II/fe/fe_q.h>
244 * #include <deal.II/fe/fe_system.h>
245 * #include <deal.II/fe/fe_tools.h>
246 * #include <deal.II/fe/fe_values.h>
247 * #include <deal.II/fe/mapping_q_eulerian.h>
248 * #include <deal.II/lac/block_sparsity_pattern.h>
249 * #include <deal.II/lac/dynamic_sparsity_pattern.h>
250 * #include <deal.II/lac/affine_constraints.h>
251 * #include <deal.II/lac/full_matrix.h>
252 * #include <deal.II/lac/solver_selector.h>
253 * #include <deal.II/lac/trilinos_block_sparse_matrix.h>
254 * #include <deal.II/lac/trilinos_precondition.h>
255 * #include <deal.II/lac/trilinos_sparsity_pattern.h>
256 * #include <deal.II/lac/trilinos_sparse_matrix.h>
257 * #include <deal.II/lac/trilinos_vector.h>
259 * #include <deal.II/lac/packaged_operation.h>
260 * #include <deal.II/lac/trilinos_linear_operator.h>
261 * #include <deal.II/numerics/data_out.h>
262 * #include <deal.II/numerics/vector_tools.h>
263 * #include <deal.II/physics/transformations.h>
264 * #include <deal.II/physics/elasticity/kinematics.h>
265 * #include <deal.II/physics/elasticity/standard_tensors.h>
266 * #include <iostream>
270 * #include <deal.II/grid/grid_out.h>
272 *
namespace ViscoElasStripHole
276 *
namespace Parameters
278 *
struct BoundaryConditions
280 * BoundaryConditions();
282 * std::string driver;
301 * BoundaryConditions::BoundaryConditions()
303 * driver (
"Neumann"),
307 * boundary_id_minus_X (1),
308 * boundary_id_plus_X (2),
309 * boundary_id_minus_Y (3),
310 * boundary_id_plus_Y (4),
311 * boundary_id_minus_Z (5),
312 * boundary_id_plus_Z (6),
313 * boundary_id_hole (10),
314 * manifold_id_hole (10)
318 * prm.enter_subsection(
"Boundary conditions");
320 * prm.declare_entry(
"Driver",
"Dirichlet",
322 *
"Driver boundary condition for the problem");
323 * prm.declare_entry(
"Final stretch",
"2.0",
325 *
"Positive stretch applied length-ways to the strip");
326 * prm.declare_entry(
"Applied pressure",
"0.0",
328 *
"Hydrostatic pressure applied (in the referential configuration) to the interior surface of the hole");
329 * prm.declare_entry(
"Load time",
"2.5",
331 *
"Total time over which the stretch/pressure is ramped up");
333 * prm.leave_subsection();
337 * prm.enter_subsection(
"Boundary conditions");
339 * driver = prm.get(
"Driver");
340 * stretch = prm.get_double(
"Final stretch");
341 * pressure = prm.get_double(
"Applied pressure");
342 * load_time = prm.get_double(
"Load time");
344 * prm.leave_subsection();
348 *
unsigned int poly_degree;
349 *
unsigned int quad_order;
357 * prm.enter_subsection(
"Finite element system");
359 * prm.declare_entry(
"Polynomial degree",
"2",
361 *
"Displacement system polynomial order");
362 * prm.declare_entry(
"Quadrature order",
"3",
364 *
"Gauss quadrature order");
366 * prm.leave_subsection();
370 * prm.enter_subsection(
"Finite element system");
372 * poly_degree = prm.get_integer(
"Polynomial degree");
373 * quad_order = prm.get_integer(
"Quadrature order");
375 * prm.leave_subsection();
382 *
double hole_diameter;
383 *
double hole_division_fraction;
384 *
unsigned int n_repetitions_xy;
385 *
unsigned int n_repetitions_z;
386 *
unsigned int global_refinement;
395 * prm.enter_subsection(
"Geometry");
397 * prm.declare_entry(
"Length",
"100.0",
399 *
"Total sample length");
400 * prm.declare_entry(
"Width",
"50.0",
402 *
"Total sample width");
403 * prm.declare_entry(
"Thickness",
"5.0",
405 *
"Total sample thickness");
406 * prm.declare_entry(
"Hole diameter",
"20.0",
409 * prm.declare_entry(
"Hole division fraction",
"0.5",
411 *
"A geometric factor affecting the discretisation near the hole");
412 * prm.declare_entry(
"Number of subdivisions in cross-section",
"2",
414 *
"A factor defining the number of initial grid subdivisions in the cross-section");
415 * prm.declare_entry(
"Number of subdivisions thickness",
"6",
417 *
"A factor defining the number of initial grid subdivisions through the thickness");
418 * prm.declare_entry(
"Global refinement",
"2",
420 *
"Global refinement level");
421 * prm.declare_entry(
"Grid scale",
"1e-3",
423 *
"Global grid scaling factor");
425 * prm.leave_subsection();
429 * prm.enter_subsection(
"Geometry");
431 * length = prm.get_double(
"Length");
432 * width = prm.get_double(
"Width");
433 * thickness = prm.get_double(
"Thickness");
434 * hole_diameter = prm.get_double(
"Hole diameter");
435 * hole_division_fraction = prm.get_double(
"Hole division fraction");
436 * n_repetitions_xy = prm.get_integer(
"Number of subdivisions in cross-section");
437 * n_repetitions_z = prm.get_integer(
"Number of subdivisions thickness");
438 * global_refinement = prm.get_integer(
"Global refinement");
439 *
scale = prm.get_double(
"Grid scale");
441 * prm.leave_subsection();
456 * prm.enter_subsection(
"Material properties");
458 * prm.declare_entry(
"Poisson's ratio",
"0.4999",
460 *
"Poisson's ratio");
461 * prm.declare_entry(
"Elastic shear modulus",
"80.194e6",
463 *
"Elastic shear modulus");
464 * prm.declare_entry(
"Viscous shear modulus",
"80.194e6",
466 *
"Viscous shear modulus");
467 * prm.declare_entry(
"Viscous relaxation time",
"2.0",
469 *
"Viscous relaxation time");
471 * prm.leave_subsection();
475 * prm.enter_subsection(
"Material properties");
477 * nu_e = prm.get_double(
"Poisson's ratio");
478 * mu_e = prm.get_double(
"Elastic shear modulus");
479 * mu_v = prm.get_double(
"Viscous shear modulus");
480 * tau_v = prm.get_double(
"Viscous relaxation time");
482 * prm.leave_subsection();
484 *
struct LinearSolver
486 * std::string type_lin;
488 *
double max_iterations_lin;
496 * prm.enter_subsection(
"Linear solver");
498 * prm.declare_entry(
"Solver type",
"cg",
500 *
"Type of solver used to solve the linear system");
501 * prm.declare_entry(
"Residual",
"1e-6",
503 *
"Linear solver residual (scaled by residual norm)");
504 * prm.declare_entry(
"Max iteration multiplier",
"1",
506 *
"Linear solver iterations (multiples of the system matrix size)");
508 * prm.leave_subsection();
512 * prm.enter_subsection(
"Linear solver");
514 * type_lin = prm.get(
"Solver type");
515 * tol_lin = prm.get_double(
"Residual");
516 * max_iterations_lin = prm.get_double(
"Max iteration multiplier");
518 * prm.leave_subsection();
520 *
struct NonlinearSolver
522 *
unsigned int max_iterations_NR;
532 * prm.enter_subsection(
"Nonlinear solver");
534 * prm.declare_entry(
"Max iterations Newton-Raphson",
"10",
536 *
"Number of Newton-Raphson iterations allowed");
537 * prm.declare_entry(
"Tolerance displacement",
"1.0e-6",
539 *
"Displacement error tolerance");
540 * prm.declare_entry(
"Tolerance force",
"1.0e-9",
542 *
"Force residual tolerance");
544 * prm.leave_subsection();
548 * prm.enter_subsection(
"Nonlinear solver");
550 * max_iterations_NR = prm.get_integer(
"Max iterations Newton-Raphson");
551 * tol_f = prm.get_double(
"Tolerance force");
552 * tol_u = prm.get_double(
"Tolerance displacement");
554 * prm.leave_subsection();
567 * prm.enter_subsection(
"Time");
569 * prm.declare_entry(
"End time",
"1",
572 * prm.declare_entry(
"Time step size",
"0.1",
576 * prm.leave_subsection();
580 * prm.enter_subsection(
"Time");
582 * end_time = prm.get_double(
"End time");
583 * delta_t = prm.get_double(
"Time step size");
585 * prm.leave_subsection();
587 *
struct AllParameters
588 * :
public BoundaryConditions,
592 *
public LinearSolver,
593 *
public NonlinearSolver,
596 * AllParameters(
const std::string &input_file);
602 * AllParameters::AllParameters(
const std::string &input_file)
605 * declare_parameters(prm);
606 * prm.parse_input(input_file);
607 * parse_parameters(prm);
611 * BoundaryConditions::declare_parameters(prm);
612 * FESystem::declare_parameters(prm);
613 * Geometry::declare_parameters(prm);
614 * Materials::declare_parameters(prm);
615 * LinearSolver::declare_parameters(prm);
616 * NonlinearSolver::declare_parameters(prm);
617 * Time::declare_parameters(prm);
621 * BoundaryConditions::parse_parameters(prm);
622 * FESystem::parse_parameters(prm);
623 * Geometry::parse_parameters(prm);
624 * Materials::parse_parameters(prm);
625 * LinearSolver::parse_parameters(prm);
626 * NonlinearSolver::parse_parameters(prm);
627 * Time::parse_parameters(prm);
633 * Time (
const double time_end,
634 *
const double delta_t)
638 * time_end(time_end),
643 *
double current() const
645 *
return time_current;
651 *
double get_delta_t() const
655 *
unsigned int get_timestep() const
661 * time_current += delta_t;
665 *
unsigned int timestep;
666 *
double time_current;
667 *
const double time_end;
668 *
const double delta_t;
671 *
class Material_Compressible_Three_Field_Linear_Viscoelastic
674 * Material_Compressible_Three_Field_Linear_Viscoelastic(
const double mu_e,
677 *
const double tau_v,
680 * kappa((2.0 * mu_e * (1.0 + nu_e)) / (3.0 * (1.0 - 2.0 * nu_e))),
688 *
Assert(kappa > 0, ExcInternalError());
690 * ~Material_Compressible_Three_Field_Linear_Viscoelastic()
695 *
const double &p_tilde)
const
697 *
return get_tau_iso(F) + get_tau_vol(F,p_tilde);
700 *
const double &p_tilde)
const
702 *
return get_Jc_iso(F) + get_Jc_vol(F,p_tilde);
705 * get_dPsi_vol_dJ(
const double &J_tilde)
const
707 *
return (kappa / 2.0) * (J_tilde - 1.0 / J_tilde);
710 * get_d2Psi_vol_dJ2(
const double &J_tilde)
const
712 *
return ( (kappa / 2.0) * (1.0 + 1.0 / (J_tilde * J_tilde)));
724 * Assumes first-oder backward Euler time discretisation
727 * Q_n_t = (1.0/(1.0 + time.get_delta_t()/tau_v))*(Q_t1 + (time.get_delta_t()/tau_v)*
invert(C_bar));
730 * update_end_timestep()
736 *
const double kappa;
739 *
const double tau_v;
746 *
const double &p_tilde)
const
765 * Elastic Neo-Hookean + Linder2011 eq 47
768 *
return mu_e * b_bar
772 *
const double &p_tilde)
const
775 *
return p_tilde * det_F
790 *
return (2.0 / dim) *
trace(tau_bar)
792 * - (2.0 / dim) * (tau_iso_x_I + I_x_tau_iso)
800 * Elastic Neo-Hookean + Linder2011 eq 56
812 *
virtual ~PointHistory()
815 * setup_lqp (
const Parameters::AllParameters ¶meters,
818 * material.reset(
new Material_Compressible_Three_Field_Linear_Viscoelastic<dim>(
819 * parameters.mu_e, parameters.nu_e,
820 * parameters.mu_v, parameters.tau_v,
826 *
const double &p_tilde)
const
828 *
return material->get_tau(F, p_tilde);
832 *
const double &p_tilde)
const
834 *
return material->get_Jc(F, p_tilde);
837 * get_dPsi_vol_dJ(
const double &J_tilde)
const
839 *
return material->get_dPsi_vol_dJ(J_tilde);
842 * get_d2Psi_vol_dJ2(
const double &J_tilde)
const
844 *
return material->get_d2Psi_vol_dJ2(J_tilde);
848 *
const double &p_tilde,
849 *
const double &J_tilde)
851 * material->update_internal_equilibrium(F,p_tilde,J_tilde);
854 * update_end_timestep()
856 * material->update_end_timestep();
859 * std::shared_ptr< Material_Compressible_Three_Field_Linear_Viscoelastic<dim> > material;
865 * Solid(
const std::string &input_file);
871 *
struct PerTaskData_ASM;
872 *
struct ScratchData_ASM;
877 *
const double half_length,
878 *
const double half_width,
879 *
const double hole_radius,
880 *
const unsigned int n_repetitions_xy = 1,
881 *
const double hole_division_fraction = 0.25);
883 * setup_system(LA::MPI::BlockVector &solution_delta);
885 * determine_component_extractors();
887 * assemble_system(
const LA::MPI::BlockVector &solution_delta);
890 * ScratchData_ASM &scratch,
891 * PerTaskData_ASM &data)
const;
893 * copy_local_to_global_system(
const PerTaskData_ASM &data);
895 * make_constraints(
const int &it_nr);
899 * solve_nonlinear_timestep(LA::MPI::BlockVector &solution_delta);
900 * std::pair<unsigned int, double>
901 * solve_linear_system(LA::MPI::BlockVector &newton_update);
902 * LA::MPI::BlockVector
903 * get_solution_total(
const LA::MPI::BlockVector &solution_delta)
const;
905 * update_end_timestep();
907 * output_results(
const unsigned int timestep,
908 *
const double current_time)
const;
910 * compute_vertex_positions(std::vector<double> &real_time,
911 * std::vector<std::vector<
Point<dim> > > &tracked_vertices,
912 *
const LA::MPI::BlockVector &solution_total)
const;
916 * Parallel communication
924 * Parameters::AllParameters parameters;
929 * PointHistory<dim> > quadrature_point_history;
930 *
const unsigned int degree;
933 *
const unsigned int dofs_per_cell;
937 *
static const unsigned int n_blocks = 3;
938 *
static const unsigned int n_components = dim + 2;
939 *
static const unsigned int first_u_component = 0;
940 *
static const unsigned int p_component = dim;
941 *
static const unsigned int J_component = dim + 1;
953 * std::vector<unsigned int> block_component;
960 * std::vector<IndexSet> all_locally_owned_dofs;
963 * std::vector<IndexSet> locally_owned_partitioning;
964 * std::vector<IndexSet> locally_relevant_partitioning;
965 * std::vector<types::global_dof_index> dofs_per_block;
966 * std::vector<types::global_dof_index> element_indices_u;
967 * std::vector<types::global_dof_index> element_indices_p;
968 * std::vector<types::global_dof_index> element_indices_J;
970 *
const QGauss<dim - 1> qf_face;
971 *
const unsigned int n_q_points;
972 *
const unsigned int n_q_points_f;
974 * LA::BlockSparseMatrix tangent_matrix;
975 * LA::MPI::BlockVector system_rhs;
976 * LA::MPI::BlockVector solution_n;
981 *
norm(1.0), u(1.0), p(1.0), J(1.0)
990 *
void normalise(
const Errors &rhs)
992 *
if (rhs.norm != 0.0)
1001 *
double norm, u, p, J;
1003 * Errors error_residual, error_residual_0, error_residual_norm, error_update,
1004 * error_update_0, error_update_norm;
1006 * get_error_residual(Errors &error_residual);
1008 * get_error_update(
const LA::MPI::BlockVector &newton_update,
1009 * Errors &error_update);
1010 * std::pair<double, std::pair<double,double> >
1011 * get_error_dilation(
const LA::MPI::BlockVector &solution_total)
const;
1013 * print_conv_header();
1015 * print_conv_footer(
const LA::MPI::BlockVector &solution_delta);
1017 *
template <
int dim>
1018 * Solid<dim>::Solid(
const std::string &input_file)
1020 * mpi_communicator(MPI_COMM_WORLD),
1023 * pcout(std::cout, this_mpi_process == 0),
1024 * parameters(input_file),
1026 * time(parameters.end_time, parameters.delta_t),
1027 * timer(mpi_communicator,
1031 * degree(parameters.poly_degree),
1032 * fe(
FE_Q<dim>(parameters.poly_degree), dim,
1035 * dof_handler(triangulation),
1036 * dofs_per_cell (fe.dofs_per_cell),
1037 * u_fe(first_u_component),
1038 * p_fe(p_component),
1039 * J_fe(J_component),
1040 * dofs_per_block(n_blocks),
1041 * qf_cell(parameters.quad_order),
1042 * qf_face(parameters.quad_order),
1043 * n_q_points (qf_cell.size()),
1044 * n_q_points_f (qf_face.size())
1046 *
Assert(dim==2 || dim==3, ExcMessage(
"This problem only works in 2 or 3 space dimensions."));
1047 * determine_component_extractors();
1049 *
template <
int dim>
1050 * Solid<dim>::~Solid()
1052 * dof_handler.clear();
1054 *
template <
int dim>
1055 *
void Solid<dim>::run()
1057 * LA::MPI::BlockVector solution_delta;
1060 * setup_system(solution_delta);
1064 * Set the J component to
one, and all others to
zero:
1068 * solution_n.block(J_block) = 1.;
1069 * output_results(time.get_timestep(), time.current());
1074 * Some points
for post-processing
1077 * std::vector<double> real_time;
1078 * real_time.push_back(0);
1079 * std::vector<std::vector<Point<dim> > > tracked_vertices (4);
1082 * p[1] = parameters.length/2.0;
1083 * tracked_vertices[0].push_back(p*parameters.scale);
1087 * p[1] = parameters.hole_diameter/2.0;
1088 * tracked_vertices[1].push_back(p*parameters.scale);
1092 * p[0] = parameters.hole_diameter/2.0;
1093 * tracked_vertices[2].push_back(p*parameters.scale);
1097 * p[0] = parameters.width/2.0;
1098 * tracked_vertices[3].push_back(p*parameters.scale);
1101 *
while (time.current() < time.end()+0.01*time.get_delta_t())
1103 * solve_nonlinear_timestep(solution_delta);
1104 * solution_n += solution_delta;
1105 * solution_delta = 0.0;
1106 * output_results(time.get_timestep(), time.current());
1107 * compute_vertex_positions(real_time,
1109 * get_solution_total(solution_delta));
1110 * update_end_timestep();
1114 * pcout <<
"\n\n*** Spatial position history for tracked vertices ***" << std::endl;
1115 *
for (
unsigned int t=0; t<real_time.size(); ++t)
1120 *
for (
unsigned int p=0; p<tracked_vertices.size(); ++p)
1122 *
for (
unsigned int d=0;
d<dim; ++
d)
1124 * pcout <<
"Point " << p <<
" [" <<
d <<
"]";
1125 *
if (!(p == tracked_vertices.size()-1 && d == dim-1))
1129 * pcout << std::endl;
1132 * pcout << std::setprecision(6);
1133 * pcout << real_time[t] <<
",";
1134 *
for (
unsigned int p=0; p<tracked_vertices.size(); ++p)
1136 *
Assert(tracked_vertices[p].size() == real_time.size(),
1137 * ExcMessage(
"Vertex not tracked at each timestep"));
1138 *
for (
unsigned int d=0;
d<dim; ++
d)
1140 * pcout << tracked_vertices[p][t][
d];
1141 *
if (!(p == tracked_vertices.size()-1 && d == dim-1))
1145 * pcout << std::endl;
1148 *
template <
int dim>
1149 *
struct Solid<dim>::PerTaskData_ASM
1152 * Vector<double> cell_rhs;
1153 * std::vector<types::global_dof_index> local_dof_indices;
1154 * PerTaskData_ASM(
const unsigned int dofs_per_cell)
1157 * cell_rhs(dofs_per_cell),
1158 * local_dof_indices(dofs_per_cell)
1166 *
template <
int dim>
1167 *
struct Solid<dim>::ScratchData_ASM
1169 *
const LA::MPI::BlockVector &solution_total;
1173 * Integration helper
1176 * FEValues<dim> fe_values_ref;
1177 * FEFaceValues<dim> fe_face_values_ref;
1181 * Quadrature
point solution
1184 * std::vector<Tensor<2, dim> > solution_grads_u_total;
1185 * std::vector<double> solution_values_p_total;
1186 * std::vector<double> solution_values_J_total;
1193 * std::vector<std::vector<double> > Nx;
1194 * std::vector<std::vector<Tensor<2, dim> > > grad_Nx;
1195 * std::vector<std::vector<SymmetricTensor<2, dim> > > symm_grad_Nx;
1197 * ScratchData_ASM(
const FiniteElement<dim> &fe_cell,
1198 *
const QGauss<dim> &qf_cell,
const UpdateFlags uf_cell,
1199 *
const QGauss<dim - 1> & qf_face,
const UpdateFlags uf_face,
1200 *
const LA::MPI::BlockVector &solution_total)
1202 * solution_total (solution_total),
1203 * fe_values_ref(fe_cell, qf_cell, uf_cell),
1204 * fe_face_values_ref(fe_cell, qf_face, uf_face),
1205 * solution_grads_u_total(qf_cell.size()),
1206 * solution_values_p_total(qf_cell.size()),
1207 * solution_values_J_total(qf_cell.size()),
1208 * Nx(qf_cell.size(),
1209 * std::vector<double>(fe_cell.dofs_per_cell)),
1210 * grad_Nx(qf_cell.size(),
1211 * std::vector<Tensor<2, dim> >(fe_cell.dofs_per_cell)),
1212 * symm_grad_Nx(qf_cell.size(),
1213 * std::vector<SymmetricTensor<2, dim> >
1214 * (fe_cell.dofs_per_cell))
1216 * ScratchData_ASM(
const ScratchData_ASM &rhs)
1218 * solution_total (rhs.solution_total),
1219 * fe_values_ref(rhs.fe_values_ref.get_fe(),
1220 * rhs.fe_values_ref.get_quadrature(),
1221 * rhs.fe_values_ref.get_update_flags()),
1222 * fe_face_values_ref(rhs.fe_face_values_ref.get_fe(),
1223 * rhs.fe_face_values_ref.get_quadrature(),
1224 * rhs.fe_face_values_ref.get_update_flags()),
1225 * solution_grads_u_total(rhs.solution_grads_u_total),
1226 * solution_values_p_total(rhs.solution_values_p_total),
1227 * solution_values_J_total(rhs.solution_values_J_total),
1229 * grad_Nx(rhs.grad_Nx),
1230 * symm_grad_Nx(rhs.symm_grad_Nx)
1234 *
const unsigned int n_q_points = solution_grads_u_total.size();
1235 *
const unsigned int n_dofs_per_cell = Nx[0].size();
1237 *
Assert(solution_grads_u_total.size() == n_q_points,
1238 * ExcInternalError());
1239 *
Assert(solution_values_p_total.size() == n_q_points,
1240 * ExcInternalError());
1241 *
Assert(solution_values_J_total.size() == n_q_points,
1242 * ExcInternalError());
1243 *
Assert(Nx.size() == n_q_points,
1244 * ExcInternalError());
1245 *
Assert(grad_Nx.size() == n_q_points,
1246 * ExcInternalError());
1247 *
Assert(symm_grad_Nx.size() == n_q_points,
1248 * ExcInternalError());
1250 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1252 *
Assert( Nx[q_point].size() == n_dofs_per_cell, ExcInternalError());
1253 *
Assert( grad_Nx[q_point].size() == n_dofs_per_cell,
1254 * ExcInternalError());
1255 *
Assert( symm_grad_Nx[q_point].size() == n_dofs_per_cell,
1256 * ExcInternalError());
1258 * solution_grads_u_total[q_point] = 0.0;
1259 * solution_values_p_total[q_point] = 0.0;
1260 * solution_values_J_total[q_point] = 0.0;
1261 *
for (
unsigned int k = 0; k < n_dofs_per_cell; ++k)
1263 * Nx[q_point][k] = 0.0;
1264 * grad_Nx[q_point][k] = 0.0;
1265 * symm_grad_Nx[q_point][k] = 0.0;
1271 *
void Solid<2>::make_grid()
1273 *
const int dim = 2;
1274 *
const double tol = 1
e-12;
1275 * make_2d_quarter_plate_with_hole(triangulation,
1276 * parameters.length/2.0,
1277 * parameters.width/2.0,
1278 * parameters.hole_diameter/2.0,
1279 * parameters.n_repetitions_xy,
1280 * parameters.hole_division_fraction);
1284 * Clear boundary ID
's
1287 * for (typename Triangulation<dim>::active_cell_iterator
1288 * cell = triangulation.begin_active();
1289 * cell != triangulation.end(); ++cell)
1291 * for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1292 * if (cell->face(face)->at_boundary())
1294 * cell->face(face)->set_all_boundary_ids(0);
1300 * Set boundary IDs and and manifolds
1303 * const Point<dim> centre (0,0);
1304 * for (typename Triangulation<dim>::active_cell_iterator
1305 * cell = triangulation.begin_active();
1306 * cell != triangulation.end(); ++cell)
1308 * for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1309 * if (cell->face(face)->at_boundary())
1316 * if (std::abs(cell->face(face)->center()[0] - 0.0) < tol)
1318 * cell->face(face)->set_boundary_id(parameters.boundary_id_minus_X);
1320 * else if (std::abs(cell->face(face)->center()[0] - parameters.width/2.0) < tol)
1322 * cell->face(face)->set_boundary_id(parameters.boundary_id_plus_X);
1324 * else if (std::abs(cell->face(face)->center()[1] - 0.0) < tol)
1326 * cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Y);
1328 * else if (std::abs(cell->face(face)->center()[1] - parameters.length/2.0) < tol)
1330 * cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Y);
1334 * for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1335 * if (std::abs(cell->vertex(vertex).distance(centre) - parameters.hole_diameter/2.0) < tol)
1337 * cell->face(face)->set_boundary_id(parameters.boundary_id_hole);
1347 * for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1348 * if (std::abs(cell->vertex(vertex).distance(centre) - parameters.hole_diameter/2.0) < tol)
1350 * cell->face(face)->set_manifold_id(parameters.manifold_id_hole);
1355 * static SphericalManifold<dim> spherical_manifold (centre);
1356 * triangulation.set_manifold(parameters.manifold_id_hole,spherical_manifold);
1357 * triangulation.refine_global(parameters.global_refinement);
1358 * GridTools::scale(parameters.scale,triangulation);
1361 * void Solid<3>::make_grid()
1363 * const int dim = 3;
1364 * const double tol = 1e-12;
1365 * Triangulation<2> tria_2d;
1366 * make_2d_quarter_plate_with_hole(tria_2d,
1367 * parameters.length/2.0,
1368 * parameters.width/2.0,
1369 * parameters.hole_diameter/2.0,
1370 * parameters.n_repetitions_xy,
1371 * parameters.hole_division_fraction);
1372 * GridGenerator::extrude_triangulation(tria_2d,
1373 * parameters.n_repetitions_z+1,
1374 * parameters.thickness/2.0,
1379 * Clear boundary ID's
1383 * cell = triangulation.begin_active();
1384 * cell != triangulation.end(); ++cell)
1386 *
for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1387 *
if (cell->face(face)->at_boundary())
1389 * cell->face(face)->set_all_boundary_ids(0);
1395 * Set boundary IDs and and manifolds
1401 * cell = triangulation.begin_active();
1402 * cell != triangulation.end(); ++cell)
1404 *
for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
1405 *
if (cell->face(face)->at_boundary())
1412 *
if (
std::abs(cell->face(face)->center()[0] - 0.0) < tol)
1414 * cell->face(face)->set_boundary_id(parameters.boundary_id_minus_X);
1416 *
else if (
std::abs(cell->face(face)->center()[0] - parameters.width/2.0) < tol)
1418 * cell->face(face)->set_boundary_id(parameters.boundary_id_plus_X);
1420 *
else if (
std::abs(cell->face(face)->center()[1] - 0.0) < tol)
1422 * cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Y);
1424 *
else if (
std::abs(cell->face(face)->center()[1] - parameters.length/2.0) < tol)
1426 * cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Y);
1428 *
else if (
std::abs(cell->face(face)->center()[2] - 0.0) < tol)
1430 * cell->face(face)->set_boundary_id(parameters.boundary_id_minus_Z);
1432 *
else if (
std::abs(cell->face(face)->center()[2] - parameters.thickness/2.0) < tol)
1434 * cell->face(face)->set_boundary_id(parameters.boundary_id_plus_Z);
1438 *
for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1442 * Project the cell vertex to the XY plane and
1443 * test the distance from the
cylinder axis
1446 *
Point<dim> vertex_proj = cell->vertex(vertex);
1447 * vertex_proj[2] = 0.0;
1448 *
if (
std::abs(vertex_proj.distance(centre) - parameters.hole_diameter/2.0) < tol)
1450 * cell->face(face)->set_boundary_id(parameters.boundary_id_hole);
1461 *
for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
1465 * Project the cell vertex to the XY plane and
1466 * test the distance from the
cylinder axis
1469 *
Point<dim> vertex_proj = cell->vertex(vertex);
1470 * vertex_proj[2] = 0.0;
1471 *
if (
std::abs(vertex_proj.distance(centre) - parameters.hole_diameter/2.0) < 1e-12)
1475 * Set manifold ID on face and edges
1478 * cell->face(face)->set_all_manifold_ids(parameters.manifold_id_hole);
1485 * triangulation.set_manifold(parameters.manifold_id_hole,cylindrical_manifold);
1486 * triangulation.refine_global(parameters.global_refinement);
1489 *
template <
int dim>
1490 *
void Solid<dim>::make_2d_quarter_plate_with_hole(
Triangulation<2> &tria_2d,
1491 *
const double half_length,
1492 *
const double half_width,
1493 *
const double hole_radius,
1494 *
const unsigned int n_repetitions_xy,
1495 *
const double hole_division_fraction)
1497 *
const double length = 2.0*half_length;
1498 *
const double width = 2.0*half_width;
1499 *
const double hole_diameter = 2.0*hole_radius;
1501 *
const double internal_width = hole_diameter + hole_division_fraction*(width - hole_diameter);
1506 * hole_diameter/2.0,
1507 * internal_width/2.0);
1509 * std::set<typename Triangulation<2>::active_cell_iterator > cells_to_remove;
1511 * cell = tria_plate_hole.begin_active();
1512 * cell != tria_plate_hole.end(); ++cell)
1516 * Remove all cells that are not in the first quadrant
1519 *
if (cell->center()[0] < 0.0 || cell->center()[1] < 0.0)
1520 * cells_to_remove.insert(cell);
1522 *
Assert(cells_to_remove.size() > 0, ExcInternalError());
1523 *
Assert(cells_to_remove.size() != tria_plate_hole.n_active_cells(), ExcInternalError());
1532 * Subdivide the plate so that we
're left one
1533 * cell to remove (we'll replace
this with the
1534 * plate with the hole) and then make the
1535 * rest of the subdivisions so that we
're left
1536 * with cells with a decent aspect ratio
1539 * std::vector<std::vector<double> > step_sizes;
1541 * std::vector<double> subdivision_width;
1542 * subdivision_width.push_back(internal_width/2.0);
1543 * const double width_remaining = (width - internal_width)/2.0;
1544 * const unsigned int n_subs = static_cast<unsigned int>(std::max(1.0,std::ceil(width_remaining/(internal_width/2.0))));
1545 * Assert(n_subs>0, ExcInternalError());
1546 * for (unsigned int s=0; s<n_subs; ++s)
1547 * subdivision_width.push_back(width_remaining/n_subs);
1548 * step_sizes.push_back(subdivision_width);
1550 * const double sum_half_width = std::accumulate(subdivision_width.begin(), subdivision_width.end(), 0.0);
1551 * (void)sum_half_width;
1552 * Assert(std::abs(sum_half_width-width/2.0) < 1e-12, ExcInternalError());
1555 * std::vector<double> subdivision_length;
1556 * subdivision_length.push_back(internal_width/2.0);
1557 * const double length_remaining = (length - internal_width)/2.0;
1558 * const unsigned int n_subs = static_cast<unsigned int>(std::max(1.0,std::ceil(length_remaining/(internal_width/2.0))));
1559 * Assert(n_subs>0, ExcInternalError());
1560 * for (unsigned int s=0; s<n_subs; ++s)
1561 * subdivision_length.push_back(length_remaining/n_subs);
1562 * step_sizes.push_back(subdivision_length);
1564 * const double sum_half_length = std::accumulate(subdivision_length.begin(), subdivision_length.end(), 0.0);
1565 * (void)sum_half_length;
1566 * Assert(std::abs(sum_half_length-length/2.0) < 1e-12, ExcInternalError());
1569 * GridGenerator::subdivided_hyper_rectangle(tria_plate,
1571 * Point<2>(0.0, 0.0),
1572 * Point<2>(width/2.0, length/2.0));
1574 * std::set<typename Triangulation<2>::active_cell_iterator > cells_to_remove;
1575 * for (typename Triangulation<2>::active_cell_iterator
1576 * cell = tria_plate.begin_active();
1577 * cell != tria_plate.end(); ++cell)
1581 * Remove all cells that are in the first quadrant
1584 * if (cell->center()[0] < internal_width/2.0 && cell->center()[1] < internal_width/2.0)
1585 * cells_to_remove.insert(cell);
1587 * Assert(cells_to_remove.size() > 0, ExcInternalError());
1588 * Assert(cells_to_remove.size() != tria_plate.n_active_cells(), ExcInternalError());
1589 * GridGenerator::create_triangulation_with_removed_cells(tria_plate,cells_to_remove,tria_cut_plate);
1592 * Triangulation<2> tria_2d_not_flat;
1593 * GridGenerator::merge_triangulations(tria_quarter_plate_hole,
1595 * tria_2d_not_flat);
1599 * Attach a manifold to the curved boundary and refine
1600 * Note: We can only guarantee that the vertices sit on
1601 * the curve, so we must test with their position instead
1602 * of the cell centre.
1605 * const Point<2> centre_2d (0,0);
1606 * for (typename Triangulation<2>::active_cell_iterator
1607 * cell = tria_2d_not_flat.begin_active();
1608 * cell != tria_2d_not_flat.end(); ++cell)
1610 * for (unsigned int face=0; face<GeometryInfo<2>::faces_per_cell; ++face)
1611 * if (cell->face(face)->at_boundary())
1612 * for (unsigned int vertex=0; vertex<GeometryInfo<2>::vertices_per_face; ++vertex)
1613 * if (std::abs(cell->vertex(vertex).distance(centre_2d) - hole_diameter/2.0) < 1e-12)
1615 * cell->face(face)->set_manifold_id(10);
1619 * SphericalManifold<2> spherical_manifold_2d (centre_2d);
1620 * tria_2d_not_flat.set_manifold(10,spherical_manifold_2d);
1621 * tria_2d_not_flat.refine_global(std::max (1U, n_repetitions_xy));
1622 * tria_2d_not_flat.reset_manifold(10); // Clear manifold
1624 * GridGenerator::flatten_triangulation(tria_2d_not_flat,tria_2d);
1626 * template <int dim>
1627 * void Solid<dim>::setup_system(LA::MPI::BlockVector &solution_delta)
1629 * timer.enter_subsection("Setup system");
1630 * pcout << "Setting up linear system..." << std::endl;
1632 * block_component = std::vector<unsigned int> (n_components, u_block); // Displacement
1633 * block_component[p_component] = p_block; // Pressure
1634 * block_component[J_component] = J_block; // Dilatation
1635 * dof_handler.distribute_dofs(fe);
1636 * DoFRenumbering::Cuthill_McKee(dof_handler);
1637 * DoFRenumbering::component_wise(dof_handler, block_component);
1641 * Count DoFs in each block
1644 * dofs_per_block = DoFTools::count_dofs_per_fe_block(dof_handler, block_component);
1646 * all_locally_owned_dofs = DoFTools::locally_owned_dofs_per_subdomain (dof_handler);
1647 * std::vector<IndexSet> all_locally_relevant_dofs
1648 * = DoFTools::locally_relevant_dofs_per_subdomain (dof_handler);
1650 * locally_owned_dofs.clear();
1651 * locally_owned_partitioning.clear();
1652 * Assert(all_locally_owned_dofs.size() > this_mpi_process, ExcInternalError());
1653 * locally_owned_dofs = all_locally_owned_dofs[this_mpi_process];
1655 * locally_relevant_dofs.clear();
1656 * locally_relevant_partitioning.clear();
1657 * Assert(all_locally_relevant_dofs.size() > this_mpi_process, ExcInternalError());
1658 * locally_relevant_dofs = all_locally_relevant_dofs[this_mpi_process];
1660 * locally_owned_partitioning.reserve(n_blocks);
1661 * locally_relevant_partitioning.reserve(n_blocks);
1662 * for (unsigned int b=0; b<n_blocks; ++b)
1664 * const types::global_dof_index idx_begin
1665 * = std::accumulate(dofs_per_block.begin(),
1666 * std::next(dofs_per_block.begin(),b), 0);
1667 * const types::global_dof_index idx_end
1668 * = std::accumulate(dofs_per_block.begin(),
1669 * std::next(dofs_per_block.begin(),b+1), 0);
1670 * locally_owned_partitioning.push_back(locally_owned_dofs.get_view(idx_begin, idx_end));
1671 * locally_relevant_partitioning.push_back(locally_relevant_dofs.get_view(idx_begin, idx_end));
1675 * << " Number of active cells: " << triangulation.n_active_cells()
1676 * << " (by partition:";
1677 * for (unsigned int p=0; p<n_mpi_processes; ++p)
1679 * << (p==0 ? ' ' : '+
')
1680 * << (GridTools::count_cells_with_subdomain_association (triangulation,p));
1681 * pcout << ")" << std::endl;
1684 * << " Number of degrees of freedom: " << dof_handler.n_dofs()
1685 * << " (by partition:";
1686 * for (unsigned int p=0; p<n_mpi_processes; ++p)
1688 * << (p==0 ? ' ' : '+
')
1689 * << (DoFTools::count_dofs_with_subdomain_association (dof_handler,p));
1690 * pcout << ")" << std::endl;
1692 * << " Number of degrees of freedom per block: "
1693 * << "[n_u, n_p, n_J] = ["
1694 * << dofs_per_block[u_block] << ", "
1695 * << dofs_per_block[p_block] << ", "
1696 * << dofs_per_block[J_block] << "]"
1700 * Table<2, DoFTools::Coupling> coupling(n_components, n_components);
1701 * for (unsigned int ii = 0; ii < n_components; ++ii)
1702 * for (unsigned int jj = 0; jj < n_components; ++jj)
1703 * if (((ii < p_component) && (jj == J_component))
1704 * || ((ii == J_component) && (jj < p_component))
1705 * || ((ii == p_component) && (jj == p_component)))
1706 * coupling[ii][jj] = DoFTools::none;
1708 * coupling[ii][jj] = DoFTools::always;
1710 * TrilinosWrappers::BlockSparsityPattern bsp (locally_owned_partitioning,
1711 * locally_owned_partitioning,
1712 * locally_relevant_partitioning,
1713 * mpi_communicator);
1714 * DoFTools::make_sparsity_pattern (dof_handler, bsp,
1715 * constraints, false,
1716 * this_mpi_process);
1718 * tangent_matrix.reinit (bsp);
1722 * We then set up storage vectors
1725 * system_rhs.reinit(locally_owned_partitioning,
1726 * mpi_communicator);
1727 * solution_n.reinit(locally_owned_partitioning,
1728 * mpi_communicator);
1729 * solution_delta.reinit(locally_owned_partitioning,
1730 * mpi_communicator);
1732 * timer.leave_subsection();
1734 * template <int dim>
1736 * Solid<dim>::determine_component_extractors()
1738 * element_indices_u.clear();
1739 * element_indices_p.clear();
1740 * element_indices_J.clear();
1741 * for (unsigned int k = 0; k < fe.dofs_per_cell; ++k)
1743 * const unsigned int k_group = fe.system_to_base_index(k).first.first;
1744 * if (k_group == u_block)
1745 * element_indices_u.push_back(k);
1746 * else if (k_group == p_block)
1747 * element_indices_p.push_back(k);
1748 * else if (k_group == J_block)
1749 * element_indices_J.push_back(k);
1752 * Assert(k_group <= J_block, ExcInternalError());
1756 * template <int dim>
1757 * void Solid<dim>::setup_qph()
1759 * pcout << "Setting up quadrature point data..." << std::endl;
1760 * quadrature_point_history.initialize(triangulation.begin_active(),
1761 * triangulation.end(),
1763 * FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1764 * cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1765 * dof_handler.begin_active()),
1766 * endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1767 * dof_handler.end());
1768 * for (; cell!=endc; ++cell)
1770 * Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1771 * const std::vector<std::shared_ptr<PointHistory<dim> > > lqph =
1772 * quadrature_point_history.get_data(cell);
1773 * Assert(lqph.size() == n_q_points, ExcInternalError());
1774 * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1775 * lqph[q_point]->setup_lqp(parameters, time);
1778 * template <int dim>
1780 * Solid<dim>::solve_nonlinear_timestep(LA::MPI::BlockVector &solution_delta)
1782 * pcout << std::endl
1783 * << "Timestep " << time.get_timestep() << " @ "
1784 * << time.current() << "s of "
1785 * << time.end() << "s" << std::endl;
1786 * LA::MPI::BlockVector newton_update(locally_owned_partitioning,
1787 * mpi_communicator);
1788 * error_residual.reset();
1789 * error_residual_0.reset();
1790 * error_residual_norm.reset();
1791 * error_update.reset();
1792 * error_update_0.reset();
1793 * error_update_norm.reset();
1794 * print_conv_header();
1795 * unsigned int newton_iteration = 0;
1796 * for (; newton_iteration < parameters.max_iterations_NR;
1797 * ++newton_iteration)
1799 * pcout << " " << std::setw(2) << newton_iteration << " " << std::flush;
1800 * make_constraints(newton_iteration);
1801 * assemble_system(solution_delta);
1802 * get_error_residual(error_residual);
1803 * if (newton_iteration == 0)
1804 * error_residual_0 = error_residual;
1805 * error_residual_norm = error_residual;
1806 * error_residual_norm.normalise(error_residual_0);
1807 * if (newton_iteration > 0 &&
1808 * (error_update_norm.u <= parameters.tol_u &&
1809 * error_residual_norm.u <= parameters.tol_f) )
1811 * pcout << " CONVERGED! " << std::endl;
1812 * print_conv_footer(solution_delta);
1815 * const std::pair<unsigned int, double>
1816 * lin_solver_output = solve_linear_system(newton_update);
1817 * get_error_update(newton_update, error_update);
1818 * if (newton_iteration == 0)
1819 * error_update_0 = error_update;
1820 * error_update_norm = error_update;
1821 * error_update_norm.normalise(error_update_0);
1822 * solution_delta += newton_update;
1823 * newton_update = 0.0;
1824 * pcout << " | " << std::fixed << std::setprecision(3) << std::setw(7)
1825 * << std::scientific << lin_solver_output.first << " "
1826 * << lin_solver_output.second << " " << error_residual_norm.norm
1827 * << " " << error_residual_norm.u << " "
1828 * << error_residual_norm.p << " " << error_residual_norm.J
1829 * << " " << error_update_norm.norm << " " << error_update_norm.u
1830 * << " " << error_update_norm.p << " " << error_update_norm.J
1831 * << " " << std::endl;
1833 * AssertThrow (newton_iteration <= parameters.max_iterations_NR,
1834 * ExcMessage("No convergence in nonlinear solver!"));
1836 * template <int dim>
1837 * void Solid<dim>::print_conv_header()
1839 * pcout << std::string(132,'_
') << std::endl;
1840 * pcout << " SOLVER STEP "
1841 * << " | LIN_IT LIN_RES RES_NORM "
1842 * << " RES_U RES_P RES_J NU_NORM "
1843 * << " NU_U NU_P NU_J " << std::endl;
1844 * pcout << std::string(132,'_
') << std::endl;
1846 * template <int dim>
1847 * void Solid<dim>::print_conv_footer(const LA::MPI::BlockVector &solution_delta)
1849 * pcout << std::string(132,'_
') << std::endl;
1850 * const std::pair<double,std::pair<double,double> > error_dil = get_error_dilation(get_solution_total(solution_delta));
1851 * pcout << "Relative errors:" << std::endl
1852 * << "Displacement:\t" << error_update.u / error_update_0.u << std::endl
1853 * << "Force: \t\t" << error_residual.u / error_residual_0.u << std::endl
1854 * << "Dilatation:\t" << error_dil.first << std::endl
1855 * << "v / V_0:\t" << error_dil.second.second << " / " << error_dil.second.first
1856 * << " = " << (error_dil.second.second/error_dil.second.first) << std::endl;
1858 * template <int dim>
1859 * std::pair<double,std::pair<double,double> >
1860 * Solid<dim>::get_error_dilation(const LA::MPI::BlockVector &solution_total) const
1862 * double vol_reference = 0.0;
1863 * double vol_current = 0.0;
1864 * double dil_L2_error = 0.0;
1865 * FEValues<dim> fe_values_ref(fe, qf_cell,
1866 * update_values | update_gradients | update_JxW_values);
1867 * std::vector<Tensor<2, dim> > solution_grads_u_total (qf_cell.size());
1868 * std::vector<double> solution_values_J_total (qf_cell.size());
1869 * FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1870 * cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1871 * dof_handler.begin_active()),
1872 * endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1873 * dof_handler.end());
1874 * for (; cell != endc; ++cell)
1876 * Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1877 * fe_values_ref.reinit(cell);
1878 * fe_values_ref[u_fe].get_function_gradients(solution_total,
1879 * solution_grads_u_total);
1880 * fe_values_ref[J_fe].get_function_values(solution_total,
1881 * solution_values_J_total);
1882 * const std::vector<std::shared_ptr<const PointHistory<dim> > > lqph =
1883 * quadrature_point_history.get_data(cell);
1884 * Assert(lqph.size() == n_q_points, ExcInternalError());
1885 * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
1887 * const double det_F_qp = determinant(Physics::Elasticity::Kinematics::F(solution_grads_u_total[q_point]));
1888 * const double J_tilde_qp = solution_values_J_total[q_point];
1889 * const double the_error_qp_squared = std::pow((det_F_qp - J_tilde_qp),
1891 * const double JxW = fe_values_ref.JxW(q_point);
1892 * dil_L2_error += the_error_qp_squared * JxW;
1893 * vol_reference += JxW;
1894 * vol_current += det_F_qp * JxW;
1897 * Assert(vol_current > 0.0, ExcInternalError());
1900 * Sum across all processors
1903 * dil_L2_error = Utilities::MPI::sum(dil_L2_error,mpi_communicator);
1904 * vol_reference = Utilities::MPI::sum(vol_reference,mpi_communicator);
1905 * vol_current = Utilities::MPI::sum(vol_current,mpi_communicator);
1907 * return std::make_pair(std::sqrt(dil_L2_error),
1908 * std::make_pair(vol_reference,vol_current));
1910 * template <int dim>
1911 * void Solid<dim>::get_error_residual(Errors &error_residual)
1915 * Construct a residual vector that has the values for all of its
1916 * constrained DoFs set to zero.
1919 * LA::MPI::BlockVector error_res (system_rhs);
1920 * constraints.set_zero(error_res);
1921 * error_residual.norm = error_res.l2_norm();
1922 * error_residual.u = error_res.block(u_block).l2_norm();
1923 * error_residual.p = error_res.block(p_block).l2_norm();
1924 * error_residual.J = error_res.block(J_block).l2_norm();
1926 * template <int dim>
1927 * void Solid<dim>::get_error_update(const LA::MPI::BlockVector &newton_update,
1928 * Errors &error_update)
1932 * Construct a update vector that has the values for all of its
1933 * constrained DoFs set to zero.
1936 * LA::MPI::BlockVector error_ud (newton_update);
1937 * constraints.set_zero(error_ud);
1938 * error_update.norm = error_ud.l2_norm();
1939 * error_update.u = error_ud.block(u_block).l2_norm();
1940 * error_update.p = error_ud.block(p_block).l2_norm();
1941 * error_update.J = error_ud.block(J_block).l2_norm();
1943 * template <int dim>
1944 * LA::MPI::BlockVector
1945 * Solid<dim>::get_solution_total(const LA::MPI::BlockVector &solution_delta) const
1949 * Cell interpolation -> Ghosted vector
1952 * LA::MPI::BlockVector solution_total (locally_owned_partitioning,
1953 * locally_relevant_partitioning,
1955 * /*vector_writable = */ false);
1956 * LA::MPI::BlockVector tmp (solution_total);
1957 * solution_total = solution_n;
1958 * tmp = solution_delta;
1959 * solution_total += tmp;
1960 * return solution_total;
1962 * template <int dim>
1963 * void Solid<dim>::assemble_system(const LA::MPI::BlockVector &solution_delta)
1965 * timer.enter_subsection("Assemble system");
1966 * pcout << " ASM_SYS " << std::flush;
1967 * tangent_matrix = 0.0;
1969 * const LA::MPI::BlockVector solution_total(get_solution_total(solution_delta));
1970 * const UpdateFlags uf_cell(update_values |
1971 * update_gradients |
1972 * update_JxW_values);
1973 * const UpdateFlags uf_face(update_values |
1974 * update_normal_vectors |
1975 * update_JxW_values);
1976 * PerTaskData_ASM per_task_data(dofs_per_cell);
1977 * ScratchData_ASM scratch_data(fe, qf_cell, uf_cell, qf_face, uf_face, solution_total);
1979 * FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
1980 * cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1981 * dof_handler.begin_active()),
1982 * endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
1983 * dof_handler.end());
1984 * for (; cell != endc; ++cell)
1986 * Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
1987 * assemble_system_one_cell(cell, scratch_data, per_task_data);
1988 * copy_local_to_global_system(per_task_data);
1990 * tangent_matrix.compress(VectorOperation::add);
1991 * system_rhs.compress(VectorOperation::add);
1992 * timer.leave_subsection();
1994 * template <int dim>
1995 * void Solid<dim>::copy_local_to_global_system(const PerTaskData_ASM &data)
1997 * constraints.distribute_local_to_global(data.cell_matrix, data.cell_rhs,
1998 * data.local_dof_indices,
1999 * tangent_matrix, system_rhs);
2001 * template <int dim>
2003 * Solid<dim>::assemble_system_one_cell(const typename DoFHandler<dim>::active_cell_iterator &cell,
2004 * ScratchData_ASM &scratch,
2005 * PerTaskData_ASM &data) const
2007 * Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2011 * scratch.fe_values_ref.reinit(cell);
2012 * cell->get_dof_indices(data.local_dof_indices);
2013 * const std::vector<std::shared_ptr<const PointHistory<dim> > > lqph =
2014 * quadrature_point_history.get_data(cell);
2015 * Assert(lqph.size() == n_q_points, ExcInternalError());
2019 * Update quadrature point solution
2022 * scratch.fe_values_ref[u_fe].get_function_gradients(scratch.solution_total,
2023 * scratch.solution_grads_u_total);
2024 * scratch.fe_values_ref[p_fe].get_function_values(scratch.solution_total,
2025 * scratch.solution_values_p_total);
2026 * scratch.fe_values_ref[J_fe].get_function_values(scratch.solution_total,
2027 * scratch.solution_values_J_total);
2031 * Update shape functions and their gradients (push-forward)
2034 * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2036 * const Tensor<2, dim> F = Physics::Elasticity::Kinematics::F(scratch.solution_grads_u_total[q_point]);
2037 * const Tensor<2, dim> F_inv = invert(F);
2039 * for (unsigned int k = 0; k < dofs_per_cell; ++k)
2041 * const unsigned int k_group = fe.system_to_base_index(k).first.first;
2042 * if (k_group == u_block)
2044 * scratch.grad_Nx[q_point][k] = scratch.fe_values_ref[u_fe].gradient(k, q_point)
2046 * scratch.symm_grad_Nx[q_point][k] = symmetrize(scratch.grad_Nx[q_point][k]);
2048 * else if (k_group == p_block)
2049 * scratch.Nx[q_point][k] = scratch.fe_values_ref[p_fe].value(k,
2051 * else if (k_group == J_block)
2052 * scratch.Nx[q_point][k] = scratch.fe_values_ref[J_fe].value(k,
2055 * Assert(k_group <= J_block, ExcInternalError());
2058 * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2060 * const SymmetricTensor<2, dim> &I = Physics::Elasticity::StandardTensors<dim>::I;
2061 * const Tensor<2, dim> F = Physics::Elasticity::Kinematics::F(scratch.solution_grads_u_total[q_point]);
2062 * const double det_F = determinant(F);
2063 * const double &p_tilde = scratch.solution_values_p_total[q_point];
2064 * const double &J_tilde = scratch.solution_values_J_total[q_point];
2065 * Assert(det_F > 0, ExcInternalError());
2068 * PointHistory<dim> *lqph_q_point_nc = const_cast<PointHistory<dim>*>(lqph[q_point].get());
2069 * lqph_q_point_nc->update_internal_equilibrium(F,p_tilde,J_tilde);
2072 * const SymmetricTensor<2, dim> tau = lqph[q_point]->get_tau(F,p_tilde);
2073 * const Tensor<2, dim> tau_ns (tau);
2074 * const SymmetricTensor<4, dim> Jc = lqph[q_point]->get_Jc(F,p_tilde);
2075 * const double dPsi_vol_dJ = lqph[q_point]->get_dPsi_vol_dJ(J_tilde);
2076 * const double d2Psi_vol_dJ2 = lqph[q_point]->get_d2Psi_vol_dJ2(J_tilde);
2078 * const std::vector<double> &Nx = scratch.Nx[q_point];
2079 * const std::vector<Tensor<2, dim> > &grad_Nx = scratch.grad_Nx[q_point];
2080 * const std::vector<SymmetricTensor<2, dim> > &symm_grad_Nx = scratch.symm_grad_Nx[q_point];
2081 * const double JxW = scratch.fe_values_ref.JxW(q_point);
2083 * for (unsigned int i = 0; i < dofs_per_cell; ++i)
2085 * const unsigned int component_i = fe.system_to_component_index(i).first;
2086 * const unsigned int i_group = fe.system_to_base_index(i).first.first;
2087 * if (i_group == u_block)
2088 * data.cell_rhs(i) -= (symm_grad_Nx[i] * tau) * JxW;
2089 * else if (i_group == p_block)
2090 * data.cell_rhs(i) -= Nx[i] * (det_F - J_tilde) * JxW;
2091 * else if (i_group == J_block)
2092 * data.cell_rhs(i) -= Nx[i] * (dPsi_vol_dJ - p_tilde) * JxW;
2094 * Assert(i_group <= J_block, ExcInternalError());
2096 * for (unsigned int j = 0; j <= i; ++j)
2098 * const unsigned int component_j = fe.system_to_component_index(j).first;
2099 * const unsigned int j_group = fe.system_to_base_index(j).first.first;
2100 * if ((i_group == u_block) && (j_group == u_block))
2102 * data.cell_matrix(i, j) += symm_grad_Nx[i] * Jc // The material contribution:
2103 * * symm_grad_Nx[j] * JxW;
2104 * if (component_i == component_j) // geometrical stress contribution
2105 * data.cell_matrix(i, j) += grad_Nx[i][component_i] * tau_ns
2106 * * grad_Nx[j][component_j] * JxW;
2108 * else if ((i_group == u_block) && (j_group == p_block))
2110 * data.cell_matrix(i, j) += (symm_grad_Nx[i] * I)
2114 * else if ((i_group == p_block) && (j_group == u_block))
2116 * data.cell_matrix(i, j) += Nx[i] * det_F
2117 * * (symm_grad_Nx[j] * I)
2120 * else if ((i_group == p_block) && (j_group == J_block))
2121 * data.cell_matrix(i, j) -= Nx[i] * Nx[j] * JxW;
2122 * else if ((i_group == J_block) && (j_group == p_block))
2123 * data.cell_matrix(i, j) -= Nx[i] * Nx[j] * JxW;
2124 * else if ((i_group == J_block) && (j_group == J_block))
2125 * data.cell_matrix(i, j) += Nx[i] * d2Psi_vol_dJ2 * Nx[j] * JxW;
2127 * Assert((i_group <= J_block) && (j_group <= J_block),
2128 * ExcInternalError());
2133 * for (unsigned int i = 0; i < dofs_per_cell; ++i)
2134 * for (unsigned int j = i + 1; j < dofs_per_cell; ++j)
2135 * data.cell_matrix(i, j) = data.cell_matrix(j, i);
2137 * if (parameters.driver == "Neumann")
2138 * for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
2140 * if (cell->face(face)->at_boundary() == true
2141 * && cell->face(face)->boundary_id() == parameters.boundary_id_plus_Y)
2143 * scratch.fe_face_values_ref.reinit(cell, face);
2144 * for (unsigned int f_q_point = 0; f_q_point < n_q_points_f;
2147 * const Tensor<1, dim> &N =
2148 * scratch.fe_face_values_ref.normal_vector(f_q_point);
2149 * static const double pressure_nom = parameters.pressure
2150 * / (parameters.scale * parameters.scale);
2151 * const double time_ramp = (time.current() < parameters.load_time ?
2152 * time.current() / parameters.load_time : 1.0);
2153 * const double pressure = -pressure_nom * time_ramp;
2154 * const Tensor<1, dim> traction = pressure * N;
2155 * for (unsigned int i = 0; i < dofs_per_cell; ++i)
2157 * const unsigned int i_group =
2158 * fe.system_to_base_index(i).first.first;
2159 * if (i_group == u_block)
2161 * const unsigned int component_i =
2162 * fe.system_to_component_index(i).first;
2164 * scratch.fe_face_values_ref.shape_value(i,
2166 * const double JxW = scratch.fe_face_values_ref.JxW(
2168 * data.cell_rhs(i) += (Ni * traction[component_i])
2175 * template <int dim>
2176 * void Solid<dim>::make_constraints(const int &it_nr)
2178 * pcout << " CST " << std::flush;
2181 * constraints.clear();
2182 * #if DEAL_II_VERSION_GTE(9, 6, 0)
2183 * constraints.reinit (locally_owned_dofs,
2184 * locally_relevant_dofs);
2186 * constraints.reinit (locally_relevant_dofs);
2189 * const bool apply_dirichlet_bc = (it_nr == 0);
2190 * const FEValuesExtractors::Scalar x_displacement(0);
2191 * const FEValuesExtractors::Scalar y_displacement(1);
2193 * const int boundary_id = parameters.boundary_id_minus_X;
2194 * if (apply_dirichlet_bc == true)
2195 * VectorTools::interpolate_boundary_values(dof_handler,
2197 * Functions::ZeroFunction<dim>(n_components),
2199 * fe.component_mask(x_displacement));
2201 * VectorTools::interpolate_boundary_values(dof_handler,
2203 * Functions::ZeroFunction<dim>(n_components),
2205 * fe.component_mask(x_displacement));
2208 * const int boundary_id = parameters.boundary_id_minus_Y;
2209 * if (apply_dirichlet_bc == true)
2210 * VectorTools::interpolate_boundary_values(dof_handler,
2212 * Functions::ZeroFunction<dim>(n_components),
2214 * fe.component_mask(y_displacement));
2216 * VectorTools::interpolate_boundary_values(dof_handler,
2218 * Functions::ZeroFunction<dim>(n_components),
2220 * fe.component_mask(y_displacement));
2224 * const FEValuesExtractors::Scalar z_displacement(2);
2226 * const int boundary_id = parameters.boundary_id_minus_Z;
2227 * if (apply_dirichlet_bc == true)
2228 * VectorTools::interpolate_boundary_values(dof_handler,
2230 * Functions::ZeroFunction<dim>(n_components),
2232 * fe.component_mask(z_displacement));
2234 * VectorTools::interpolate_boundary_values(dof_handler,
2236 * Functions::ZeroFunction<dim>(n_components),
2238 * fe.component_mask(z_displacement));
2241 * const int boundary_id = parameters.boundary_id_plus_Z;
2242 * if (apply_dirichlet_bc == true)
2243 * VectorTools::interpolate_boundary_values(dof_handler,
2245 * Functions::ZeroFunction<dim>(n_components),
2247 * fe.component_mask(z_displacement));
2249 * VectorTools::interpolate_boundary_values(dof_handler,
2251 * Functions::ZeroFunction<dim>(n_components),
2253 * fe.component_mask(z_displacement));
2256 * if (parameters.driver == "Dirichlet")
2258 * const int boundary_id = parameters.boundary_id_plus_Y;
2259 * if (apply_dirichlet_bc == true)
2262 * if (time.current() < parameters.load_time+0.01*time.get_delta_t())
2264 * const double delta_length = parameters.length*(parameters.stretch - 1.0)*parameters.scale;
2265 * const unsigned int n_stretch_steps = static_cast<unsigned int>(parameters.load_time/time.get_delta_t());
2266 * const double delta_u_y = delta_length/2.0/n_stretch_steps;
2267 * VectorTools::interpolate_boundary_values(dof_handler,
2269 * Functions::ConstantFunction<dim>(delta_u_y,n_components),
2271 * fe.component_mask(y_displacement));
2274 * VectorTools::interpolate_boundary_values(dof_handler,
2276 * Functions::ZeroFunction<dim>(n_components),
2278 * fe.component_mask(y_displacement));
2281 * VectorTools::interpolate_boundary_values(dof_handler,
2283 * Functions::ZeroFunction<dim>(n_components),
2285 * fe.component_mask(y_displacement));
2287 * constraints.close();
2289 * template <int dim>
2290 * std::pair<unsigned int, double>
2291 * Solid<dim>::solve_linear_system(LA::MPI::BlockVector &newton_update)
2293 * unsigned int lin_it = 0;
2294 * double lin_res = 0.0;
2296 * timer.enter_subsection("Linear solver");
2297 * pcout << " SLV " << std::flush;
2299 * const LA::MPI::Vector &f_u = system_rhs.block(u_block);
2300 * const LA::MPI::Vector &f_p = system_rhs.block(p_block);
2301 * const LA::MPI::Vector &f_J = system_rhs.block(J_block);
2302 * LA::MPI::Vector &d_u = newton_update.block(u_block);
2303 * LA::MPI::Vector &d_p = newton_update.block(p_block);
2304 * LA::MPI::Vector &d_J = newton_update.block(J_block);
2305 * const auto K_uu = linear_operator<LA::MPI::Vector>(tangent_matrix.block(u_block, u_block));
2306 * const auto K_up = linear_operator<LA::MPI::Vector>(tangent_matrix.block(u_block, p_block));
2307 * const auto K_pu = linear_operator<LA::MPI::Vector>(tangent_matrix.block(p_block, u_block));
2308 * const auto K_Jp = linear_operator<LA::MPI::Vector>(tangent_matrix.block(J_block, p_block));
2309 * const auto K_JJ = linear_operator<LA::MPI::Vector>(tangent_matrix.block(J_block, J_block));
2311 * LA::PreconditionJacobi preconditioner_K_Jp_inv;
2312 * preconditioner_K_Jp_inv.initialize(
2313 * tangent_matrix.block(J_block, p_block),
2314 * LA::PreconditionJacobi::AdditionalData());
2315 * ReductionControl solver_control_K_Jp_inv (
2316 * static_cast<unsigned int>(tangent_matrix.block(J_block, p_block).m()
2317 * * parameters.max_iterations_lin),
2319 * ::SolverCG<LA::MPI::Vector> solver_K_Jp_inv (solver_control_K_Jp_inv);
2321 * const auto K_Jp_inv = inverse_operator(K_Jp,
2323 * preconditioner_K_Jp_inv);
2324 * const auto K_pJ_inv = transpose_operator(K_Jp_inv);
2325 * const auto K_pp_bar = K_Jp_inv * K_JJ * K_pJ_inv;
2326 * const auto K_uu_bar_bar = K_up * K_pp_bar * K_pu;
2327 * const auto K_uu_con = K_uu + K_uu_bar_bar;
2329 * LA::PreconditionAMG preconditioner_K_con_inv;
2330 * preconditioner_K_con_inv.initialize(
2331 * tangent_matrix.block(u_block, u_block),
2332 * LA::PreconditionAMG::AdditionalData(
2333 * true /*elliptic*/,
2334 * (parameters.poly_degree > 1 /*higher_order_elements*/)) );
2335 * ReductionControl solver_control_K_con_inv (
2336 * static_cast<unsigned int>(tangent_matrix.block(u_block, u_block).m()
2337 * * parameters.max_iterations_lin),
2338 * 1.0e-30, parameters.tol_lin);
2339 * ::SolverSelector<LA::MPI::Vector> solver_K_con_inv;
2340 * solver_K_con_inv.select(parameters.type_lin);
2341 * solver_K_con_inv.set_control(solver_control_K_con_inv);
2342 * const auto K_uu_con_inv = inverse_operator(K_uu_con,
2344 * preconditioner_K_con_inv);
2346 * d_u = K_uu_con_inv*(f_u - K_up*(K_Jp_inv*f_J - K_pp_bar*f_p));
2347 * lin_it = solver_control_K_con_inv.last_step();
2348 * lin_res = solver_control_K_con_inv.last_value();
2349 * timer.leave_subsection();
2351 * timer.enter_subsection("Linear solver postprocessing");
2352 * d_J = K_pJ_inv*(f_p - K_pu*d_u);
2353 * d_p = K_Jp_inv*(f_J - K_JJ*d_J);
2354 * timer.leave_subsection();
2356 * constraints.distribute(newton_update);
2357 * return std::make_pair(lin_it, lin_res);
2359 * template <int dim>
2361 * Solid<dim>::update_end_timestep ()
2363 * FilteredIterator<typename DoFHandler<dim>::active_cell_iterator>
2364 * cell (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2365 * dof_handler.begin_active()),
2366 * endc (IteratorFilters::SubdomainEqualTo(this_mpi_process),
2367 * dof_handler.end());
2368 * for (; cell != endc; ++cell)
2370 * Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2371 * const std::vector<std::shared_ptr<PointHistory<dim> > > lqph =
2372 * quadrature_point_history.get_data(cell);
2373 * Assert(lqph.size() == n_q_points, ExcInternalError());
2374 * for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
2375 * lqph[q_point]->update_end_timestep();
2380 * class FilteredDataOut : public DataOut<dim>
2383 * FilteredDataOut (const unsigned int subdomain_id)
2385 * subdomain_id (subdomain_id)
2388 * virtual ~FilteredDataOut() {}
2390 * virtual typename DataOut<dim>::cell_iterator
2393 * auto cell = this->dofs->begin_active();
2394 * while ((cell != this->dofs->end()) &&
2395 * (cell->subdomain_id() != subdomain_id))
2400 * virtual typename DataOut<dim>::cell_iterator
2401 * next_cell (const typename DataOut<dim>::cell_iterator &old_cell)
2403 * if (old_cell != this->dofs->end())
2405 * const IteratorFilters::SubdomainEqualTo predicate(subdomain_id);
2407 * ++(FilteredIterator
2408 * <typename DataOut<dim>::cell_iterator>
2409 * (predicate,old_cell));
2416 * const unsigned int subdomain_id;
2419 * template <int dim>
2420 * void Solid<dim>::output_results(const unsigned int timestep,
2421 * const double current_time) const
2425 * Output -> Ghosted vector
2428 * LA::MPI::BlockVector solution_total (locally_owned_partitioning,
2429 * locally_relevant_partitioning,
2431 * /*vector_writable = */ false);
2432 * LA::MPI::BlockVector residual (locally_owned_partitioning,
2433 * locally_relevant_partitioning,
2435 * /*vector_writable = */ false);
2436 * solution_total = solution_n;
2437 * residual = system_rhs;
2442 * --- Additional data ---
2445 * Vector<double> material_id;
2446 * Vector<double> polynomial_order;
2447 * material_id.reinit(triangulation.n_active_cells());
2448 * polynomial_order.reinit(triangulation.n_active_cells());
2449 * std::vector<types::subdomain_id> partition_int (triangulation.n_active_cells());
2451 * FilteredDataOut<dim> data_out(this_mpi_process);
2452 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
2453 * data_component_interpretation(dim,
2454 * DataComponentInterpretation::component_is_part_of_vector);
2455 * data_component_interpretation.push_back(DataComponentInterpretation::component_is_scalar);
2456 * data_component_interpretation.push_back(DataComponentInterpretation::component_is_scalar);
2458 * GridTools::get_subdomain_association (triangulation, partition_int);
2462 * Can't use filtered iterators here because the cell
2463 * count
"c" is incorrect
for the
parallel case
2466 *
unsigned int c = 0;
2468 * cell = dof_handler.begin_active(),
2469 * endc = dof_handler.end();
2470 *
for (; cell!=endc; ++cell, ++c)
2472 *
if (cell->subdomain_id() != this_mpi_process)
continue;
2474 *
material_id(c) =
static_cast<int>(cell->material_id());
2477 * std::vector<std::string> solution_name(n_components,
"solution_");
2478 * std::vector<std::string> residual_name(n_components,
"residual_");
2479 *
for (
unsigned int c=0; c<n_components; ++c)
2481 *
if (block_component[c] == u_block)
2483 * solution_name[c] +=
"u";
2484 * residual_name[c] +=
"u";
2486 *
else if (block_component[c] == p_block)
2488 * solution_name[c] +=
"p";
2489 * residual_name[c] +=
"p";
2491 *
else if (block_component[c] == J_block)
2493 * solution_name[c] +=
"J";
2494 * residual_name[c] +=
"J";
2498 *
Assert(c <= J_block, ExcInternalError());
2502 * data_out.attach_dof_handler(dof_handler);
2503 * data_out.add_data_vector(solution_total,
2506 * data_component_interpretation);
2507 * data_out.add_data_vector(residual,
2510 * data_component_interpretation);
2512 * partition_int.end());
2513 * data_out.add_data_vector (material_id,
"material_id");
2514 * data_out.add_data_vector (partitioning,
"partitioning");
2515 * data_out.build_patches(degree);
2519 *
static std::string get_filename_vtu (
unsigned int process,
2520 *
unsigned int timestep,
2521 *
const unsigned int n_digits = 4)
2523 * std::ostringstream filename_vtu;
2526 * << (std::to_string(dim) +
"d")
2532 *
return filename_vtu.str();
2535 *
static std::string get_filename_pvtu (
unsigned int timestep,
2536 *
const unsigned int n_digits = 4)
2538 * std::ostringstream filename_vtu;
2541 * << (std::to_string(dim) +
"d")
2545 *
return filename_vtu.str();
2548 *
static std::string get_filename_pvd (
void)
2550 * std::ostringstream filename_vtu;
2553 * << (std::to_string(dim) +
"d")
2555 *
return filename_vtu.str();
2561 * Write out main data file
2564 *
const std::string filename_vtu = Filename::get_filename_vtu(this_mpi_process, timestep);
2565 * std::ofstream output(filename_vtu.c_str());
2566 * data_out.write_vtu(output);
2570 * Collection of files written in
parallel
2571 * This next set of steps should only be performed
2575 *
if (this_mpi_process == 0)
2579 * List of all files written out at
this timestep by all processors
2582 * std::vector<std::string> parallel_filenames_vtu;
2585 * parallel_filenames_vtu.push_back(Filename::get_filename_vtu(p, timestep));
2588 *
const std::string filename_pvtu (Filename::get_filename_pvtu(timestep));
2589 * std::ofstream pvtu_master(filename_pvtu.c_str());
2590 * data_out.write_pvtu_record(pvtu_master,
2591 * parallel_filenames_vtu);
2595 * Time dependent data master file
2598 *
static std::vector<std::pair<double,std::string> > time_and_name_history;
2599 * time_and_name_history.push_back (std::make_pair (current_time,
2601 *
const std::string filename_pvd (Filename::get_filename_pvd());
2602 * std::ofstream pvd_output (filename_pvd.c_str());
2606 *
template <
int dim>
2607 *
void Solid<dim>::compute_vertex_positions(std::vector<double> &real_time,
2608 * std::vector<std::vector<
Point<dim> > > &tracked_vertices,
2609 *
const LA::MPI::BlockVector &solution_total)
const
2611 * real_time.push_back(time.current());
2613 * std::vector<bool> vertex_found (tracked_vertices.size(),
false);
2614 * std::vector<Tensor<1,dim> > vertex_update (tracked_vertices.size());
2618 * dof_handler.begin_active()),
2620 * dof_handler.end());
2621 *
for (; cell != endc; ++cell)
2623 *
Assert(cell->subdomain_id()==this_mpi_process, ExcInternalError());
2624 *
for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2626 *
for (
unsigned int p=0; p<tracked_vertices.size(); ++p)
2628 *
if (vertex_found[p] ==
true)
continue;
2630 *
const Point<dim> pt_ref = tracked_vertices[p][0];
2631 *
if (cell->vertex(v).distance(pt_ref) < 1e-6*parameters.scale)
2633 *
for (
unsigned int d=0;
d<dim; ++
d)
2634 * vertex_update[p][d] = solution_total(cell->vertex_dof_index(v,u_block+d));
2636 * vertex_found[p] =
true;
2642 *
for (
unsigned int p=0; p<tracked_vertices.size(); ++p)
2644 *
const int found_on_n_processes =
Utilities::MPI::sum(
int(vertex_found[p]), mpi_communicator);
2645 *
Assert(found_on_n_processes>0, ExcMessage(
"Vertex not found on any processor"));
2647 *
for (
unsigned int d=0;
d<dim; ++
d)
2649 * update /= found_on_n_processes;
2650 * tracked_vertices[p].push_back(tracked_vertices[p][0] + update);
2655 *
int main (
int argc,
char *argv[])
2657 *
using namespace dealii;
2658 *
using namespace ViscoElasStripHole;
2664 *
const unsigned int dim = 2;
2665 * Solid<dim> solid(
"parameters.prm");
2668 *
catch (std::exception &exc)
2672 * std::cerr << std::endl << std::endl
2673 * <<
"----------------------------------------------------"
2675 * std::cerr <<
"Exception on processing: " << std::endl << exc.what()
2676 * << std::endl <<
"Aborting!" << std::endl
2677 * <<
"----------------------------------------------------"
2686 * std::cerr << std::endl << std::endl
2687 * <<
"----------------------------------------------------"
2689 * std::cerr <<
"Unknown exception!" << std::endl <<
"Aborting!"
2691 * <<
"----------------------------------------------------"
static constexpr const SymmetricTensor< 4, dim > dev_P
static constexpr const SymmetricTensor< 2, dim > I
static constexpr const SymmetricTensor< 4, dim > S
static constexpr const SymmetricTensor< 4, dim > IxI
static std::string get_solver_names()
#define Assert(cond, exc)
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
void write_pvd_record(std::ostream &out, const std::vector< std::pair< double, std::string > > ×_and_names)
void hyper_cube_with_cylindrical_hole(Triangulation< dim, spacedim > &triangulation, const double inner_radius=.25, const double outer_radius=.5, const double L=.5, const unsigned int repetitions=1, const bool colorize=false)
void cylinder(Triangulation< dim > &tria, const double radius=1., const double half_length=1.)
void create_triangulation_with_removed_cells(const Triangulation< dim, spacedim > &input_triangulation, const std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > &cells_to_remove, Triangulation< dim, spacedim > &result)
constexpr types::blas_int zero
constexpr types::blas_int one
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
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 > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * end(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)