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