deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
CeresFE.h
Go to the documentation of this file.
1
131 *  
132 *   /*
133 *   Summary of output files:
134 *  
135 *   One per run:
136 *  
137 *   - initial_mesh.eps : Visualization of initially imported mesh
138 *   - physical_times.txt : Columns are (1) step number corresponding to other files, (2) physical times at the time when each calculation is run in sec, (3) number of the final plasticity iteration in each timestep. Written in do_elastic_steps() for elastic steps and do_flow_step() for viscous steps
139 *  
140 *   One per time step:
141 *  
142 *   - timeXX_elastic_displacements.txt : Vtk-readable file with columns (1) x, (2) y, (3) u_x, (4) u_y, (5) P. Written in output_results() function, which is run immediately after solve().
143 *   - timeXX_baseviscosities.txt : Columns (1) cell x, (2) cell y, (3) base viscosity in Pa s. Written in solution_stresses().
144 *   - timeXX_surface.txt : Surface (defined as where P=0 boundary condition is applied) vertices at the beginning of timestep, except for the final timestep. Written in write_vertices() function, which is called immediately after setup_dofs() except for the final iteration, when it is called after move_mesh()
145 *  
146 *   One per plasticity step:
147 *  
148 *   - timeXX_flowYY.txt : Same as timeXX_elastic_displacements.txt above
149 *   - timeXX_principalstressesYY.txt : Columns with sigma1 and sigma3 at each cell. Same order as timeXX_baseviscosities.txt. Written in solution_stresses().
150 *   - timeXX_stresstensorYY.txt : Columns with components 11, 22, 33, and 13 of stress tensor at each cell. Written in solution_stresses().
151 *   - timeXX_failurelocationsYY.txt : Gives x,y coordinates of all cells where failure occurred. Written in solution_stresses().
152 *   - timeXX_viscositiesregYY.txt : Gives smoothed and regularized (i.e., floor and ceiling-filtered) effective viscosities. Written at end of solution_stresses().
153 *  
154 *   */
155 *  
156 *   #include <deal.II/base/quadrature_lib.h>
157 *   #include <deal.II/base/logstream.h>
158 *   #include <deal.II/base/function.h>
159 *   #include <deal.II/base/utilities.h>
160 *  
161 *   #include <deal.II/lac/block_vector.h>
162 *   #include <deal.II/lac/full_matrix.h>
163 *   #include <deal.II/lac/block_sparse_matrix.h>
164 *   #include <deal.II/lac/solver_cg.h>
165 *   #include <deal.II/lac/precondition.h>
166 *   #include <deal.II/lac/affine_constraints.h>
167 *  
168 *   #include <deal.II/grid/tria.h>
169 *   #include <deal.II/grid/grid_generator.h>
170 *   #include <deal.II/grid/tria_accessor.h>
171 *   #include <deal.II/grid/tria_iterator.h>
172 *   #include <deal.II/grid/grid_tools.h>
173 *   #include <deal.II/grid/manifold_lib.h>
174 *   #include <deal.II/grid/grid_refinement.h>
175 *   #include <deal.II/grid/grid_in.h>
176 *  
177 *   #include <deal.II/dofs/dof_handler.h>
178 *   #include <deal.II/dofs/dof_renumbering.h>
179 *   #include <deal.II/dofs/dof_accessor.h>
180 *   #include <deal.II/dofs/dof_tools.h>
181 *  
182 *   #include <deal.II/fe/fe_q.h>
183 *   #include <deal.II/fe/fe_system.h>
184 *   #include <deal.II/fe/fe_values.h>
185 *  
186 *   #include <deal.II/numerics/vector_tools.h>
187 *   #include <deal.II/numerics/matrix_tools.h>
188 *   #include <deal.II/numerics/data_out.h>
189 *   #include <deal.II/numerics/error_estimator.h>
190 *   #include <deal.II/numerics/derivative_approximation.h>
191 *   #include <deal.II/numerics/fe_field_function.h>
192 *  
193 *   #include <deal.II/lac/sparse_direct.h>
194 *   #include <deal.II/lac/sparse_ilu.h>
195 *  
196 *   #include <iostream>
197 *   #include <fstream>
198 *   #include <sstream>
199 *   #include <string>
200 *   #include <time.h>
201 *   #include <math.h>
202 *  
204 *   #include <armadillo>
206 *  
207 *   #include "../support_code/ellipsoid_grav.h"
208 *   #include "../support_code/ellipsoid_fit.h"
209 *   #include "../support_code/config_in.h"
210 *  
211 * @endcode
212 *
213 * As in all programs, the namespace dealii
214 * is included:
215 *
216 * @code
217 *   namespace Step22
218 *   {
219 *  
220 *   using namespace dealii;
221 *   using namespace arma;
222 *  
223 *   template<int dim>
224 *   struct InnerPreconditioner;
225 *  
226 *   template<>
227 *   struct InnerPreconditioner<2>
228 *   {
229 *   typedef SparseDirectUMFPACK type;
230 *   };
231 *  
232 *   template<>
233 *   struct InnerPreconditioner<3>
234 *   {
235 *   typedef SparseILU<double> type;
236 *   };
237 *  
238 * @endcode
239 *
240 * Auxiliary functions
241 *
242
243 *
244 *
245 * @code
246 *   template<int dim>
247 *   class AuxFunctions
248 *   {
249 *   public:
250 *   Tensor<2, 2> get_rotation_matrix(const std::vector<Tensor<1, 2> > &grad_u);
251 *   };
252 *  
253 *   template<int dim>
254 *   Tensor<2, 2> AuxFunctions<dim>::get_rotation_matrix(
255 *   const std::vector<Tensor<1, 2> > &grad_u)
256 *   {
257 *   const double curl = (grad_u[1][0] - grad_u[0][1]);
258 *  
259 *   const double angle = std::atan(curl);
260 *  
261 *   const double t[2][2] = { { cos(angle), sin(angle) }, {
262 *   -sin(angle), cos(
263 *   angle)
264 *   }
265 *   };
266 *   return Tensor<2, 2>(t);
267 *   }
268 *  
269 * @endcode
270 *
271 * Class for remembering material state/properties at each quadrature point
272 *
273
274 *
275 *
276 * @code
277 *   template<int dim>
278 *   struct PointHistory
279 *   {
280 *   SymmetricTensor<2, dim> old_stress;
281 *   double old_phiphi_stress;
282 *   double first_eta;
283 *   double new_eta;
284 *   double G;
285 *   };
286 *  
287 * @endcode
288 *
289 * Primary class of this problem
290 *
291
292 *
293 *
294 * @code
295 *   template<int dim>
296 *   class StokesProblem
297 *   {
298 *   public:
299 *   StokesProblem(const unsigned int degree);
300 *   void run();
301 *  
302 *   private:
303 *   void setup_dofs();
304 *   void assemble_system();
305 *   void solve();
306 *   void output_results() const;
307 *   void refine_mesh();
308 *   void solution_stesses();
309 *   void smooth_eta_field(std::vector<bool> failing_cells);
310 *  
311 *   void setup_initial_mesh();
312 *   void do_elastic_steps();
313 *   void do_flow_step();
314 *   void update_time_interval();
315 *   void initialize_eta_and_G();
316 *   void move_mesh();
317 *   void do_ellipse_fits();
318 *   void append_physical_times(int max_plastic);
319 *   void write_vertices(unsigned char);
320 *   void write_mesh();
321 *   void setup_quadrature_point_history();
322 *   void update_quadrature_point_history();
323 *  
324 *   const unsigned int degree;
325 *  
326 *   Triangulation<dim> triangulation;
327 *   const MappingQ1<dim> mapping;
328 *   FESystem<dim> fe;
329 *   DoFHandler<dim> dof_handler;
330 *   unsigned int n_u = 0, n_p = 0;
331 *   unsigned int plastic_iteration = 0;
332 *   unsigned int last_max_plasticity = 0;
333 *  
334 *   QGauss<dim> quadrature_formula;
335 *   std::vector< std::vector <Vector<double> > > quad_viscosities; // Indices for this object are [cell][q][q coords, eta]
336 *   std::vector<double> cell_viscosities; // This vector is only used for output, not FE computations
337 *   std::vector<PointHistory<dim> > quadrature_point_history;
338 *  
339 *   AffineConstraints<double> constraints;
340 *  
341 *   BlockSparsityPattern sparsity_pattern;
342 *   BlockSparseMatrix<double> system_matrix;
343 *  
344 *   BlockVector<double> solution;
345 *   BlockVector<double> system_rhs;
346 *  
347 *   std::shared_ptr<typename InnerPreconditioner<dim>::type> A_preconditioner;
348 *  
349 *   ellipsoid_fit<dim> ellipsoid;
350 *   };
351 *  
352 * @endcode
353 *
354 * Class for boundary conditions and rhs
355 *
356
357 *
358 *
359 * @code
360 *   template<int dim>
361 *   class BoundaryValuesP: public Function<dim>
362 *   {
363 *   public:
364 *   BoundaryValuesP() :
365 *   Function<dim>(dim + 1)
366 *   {
367 *   }
368 *  
369 *   virtual double value(const Point<dim> &p,
370 *   const unsigned int component = 0) const override;
371 *  
372 *   virtual void vector_value(const Point<dim> &p, Vector<double> &value) const override;
373 *   };
374 *  
375 *   template<int dim>
376 *   double BoundaryValuesP<dim>::value(const Point<dim> &p,
377 *   const unsigned int component) const
378 *   {
379 *   Assert(component < this->n_components,
380 *   ExcIndexRange (component, 0, this->n_components));
381 *   (void)p;
382 *   (void)component;
383 *  
384 *   Assert(p[0] >= -10.0, ExcLowerRangeType<double>(p[0], -10.0)); //value of -10 is to permit some small numerical error moving nodes left of x=0; a << value is in fact sufficient
385 *  
386 *   return 0;
387 *   }
388 *  
389 *   template<int dim>
390 *   void BoundaryValuesP<dim>::vector_value(const Point<dim> &p,
391 *   Vector<double> &values) const
392 *   {
393 *   for (unsigned int c = 0; c < this->n_components; ++c)
394 *   values(c) = BoundaryValuesP<dim>::value(p, c);
395 *   }
396 *  
397 *   template<int dim>
398 *   class RightHandSide: public Function<dim>
399 *   {
400 *   public:
401 *   RightHandSide () : Function<dim>(dim+1) {}
402 *  
403 *   using Function<dim>::value;
404 *   using Function<dim>::vector_value;
405 *   using Function<dim>::vector_value_list;
406 *  
407 *   virtual double value(const Point<dim> &p, const unsigned int component,
408 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const;
409 *  
410 *   virtual void vector_value(const Point<dim> &p, Vector<double> &value,
411 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const;
412 *  
413 *   virtual void vector_value_list(const std::vector<Point<dim> > &points,
414 *   std::vector<Vector<double> > &values,
415 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const;
416 *  
417 *   };
418 *  
419 *   template<int dim>
420 *   double RightHandSide<dim>::value(const Point<dim> &p,
421 *   const unsigned int component,
422 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const
423 *   {
424 *  
425 *   std::vector<double> temp_vector(2);
426 *   aGrav->get_gravity(p, temp_vector);
427 *  
428 *   if (component == 0)
429 *   {
430 *   return temp_vector[0] + system_parameters::omegasquared * p[0]; // * 1.2805;
431 *   }
432 *   else
433 *   {
434 *   if (component == 1)
435 *   return temp_vector[1];
436 *   else
437 *   return 0;
438 *   }
439 *   }
440 *  
441 *   template<int dim>
442 *   void RightHandSide<dim>::vector_value(const Point<dim> &p,
443 *   Vector<double> &values,
444 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const
445 *   {
446 *   for (unsigned int c = 0; c < this->n_components; ++c)
447 *   values(c) = RightHandSide<dim>::value(p, c, aGrav);
448 *   }
449 *  
450 *   template<int dim>
451 *   void RightHandSide<dim>::vector_value_list(
452 *   const std::vector<Point<dim> > &points,
453 *   std::vector<Vector<double> > &values,
454 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav) const
455 *   {
456 * @endcode
457 *
458 * check whether component is in
459 * the valid range is up to the
460 * derived class
461 *
462 * @code
463 *   Assert(values.size() == points.size(),
464 *   ExcDimensionMismatch(values.size(), points.size()));
465 *  
466 *   for (unsigned int i = 0; i < points.size(); ++i)
467 *   this->vector_value(points[i], values[i], aGrav);
468 *   }
469 *  
470 * @endcode
471 *
472 * Class for linear solvers and preconditioners
473 *
474
475 *
476 *
477 * @code
478 *   template<class Matrix, class Preconditioner>
479 *   class InverseMatrix: public Subscriptor
480 *   {
481 *   public:
482 *   InverseMatrix(const Matrix &m, const Preconditioner &preconditioner);
483 *  
484 *   void vmult(Vector<double> &dst, const Vector<double> &src) const;
485 *  
486 *   private:
488 *   const SmartPointer<const Preconditioner> preconditioner;
489 *   };
490 *  
491 *   template<class Matrix, class Preconditioner>
492 *   InverseMatrix<Matrix, Preconditioner>::InverseMatrix(const Matrix &m,
493 *   const Preconditioner &preconditioner) :
494 *   matrix(&m), preconditioner(&preconditioner)
495 *   {
496 *   }
497 *  
498 *   template<class Matrix, class Preconditioner>
499 *   void InverseMatrix<Matrix, Preconditioner>::vmult(Vector<double> &dst,
500 *   const Vector<double> &src) const
501 *   {
502 *   SolverControl solver_control(1000 * src.size(), 1e-9 * src.l2_norm());
503 *  
504 *   SolverCG<> cg(solver_control);
505 *  
506 *   dst = 0;
507 *  
508 *   cg.solve(*matrix, dst, src, *preconditioner);
509 *   }
510 *  
511 * @endcode
512 *
513 * Class for the SchurComplement
514 *
515
516 *
517 *
518 * @code
519 *   template<class Preconditioner>
520 *   class SchurComplement: public Subscriptor
521 *   {
522 *   public:
523 *   SchurComplement(const BlockSparseMatrix<double> &system_matrix,
524 *   const InverseMatrix<SparseMatrix<double>, Preconditioner> &A_inverse);
525 *  
526 *   void vmult(Vector<double> &dst, const Vector<double> &src) const;
527 *  
528 *   private:
529 *   const SmartPointer<const BlockSparseMatrix<double> > system_matrix;
530 *   const SmartPointer<const InverseMatrix<SparseMatrix<double>, Preconditioner> > A_inverse;
531 *  
532 *   mutable Vector<double> tmp1, tmp2;
533 *   };
534 *  
535 *   template<class Preconditioner>
536 *   SchurComplement<Preconditioner>::SchurComplement(
537 *   const BlockSparseMatrix<double> &system_matrix,
538 *   const InverseMatrix<SparseMatrix<double>, Preconditioner> &A_inverse) :
539 *   system_matrix(&system_matrix), A_inverse(&A_inverse), tmp1(
540 *   system_matrix.block(0, 0).m()), tmp2(
541 *   system_matrix.block(0, 0).m())
542 *   {
543 *   }
544 *  
545 *   template<class Preconditioner>
546 *   void SchurComplement<Preconditioner>::vmult(Vector<double> &dst,
547 *   const Vector<double> &src) const
548 *   {
549 *   system_matrix->block(0, 1).vmult(tmp1, src);
550 *   A_inverse->vmult(tmp2, tmp1);
551 *   system_matrix->block(1, 0).vmult(dst, tmp2);
552 *   }
553 *  
554 * @endcode
555 *
556 * StokesProblem::StokesProblem
557 *
558
559 *
560 *
561 * @code
562 *   template<int dim>
563 *   StokesProblem<dim>::StokesProblem(const unsigned int degree) :
564 *   degree(degree),
565 *   triangulation(Triangulation<dim>::maximum_smoothing),
566 *   mapping(),
567 *   fe(FE_Q<dim>(degree + 1), dim, FE_Q<dim>(degree), 1),
568 *   dof_handler(triangulation),
569 *   quadrature_formula(degree + 2),
570 *   ellipsoid(&triangulation)
571 *   {}
572 *  
573 * @endcode
574 *
575 * Set up dofs
576 *
577
578 *
579 *
580 * @code
581 *   template<int dim>
582 *   void StokesProblem<dim>::setup_dofs()
583 *   {
584 *   A_preconditioner.reset();
585 *   system_matrix.clear();
586 *  
587 *   dof_handler.distribute_dofs(fe);
588 *   DoFRenumbering::Cuthill_McKee(dof_handler);
589 *  
590 *   std::vector<unsigned int> block_component(dim + 1, 0);
591 *   block_component[dim] = 1;
592 *   DoFRenumbering::component_wise(dof_handler, block_component);
593 *  
594 * @endcode
595 *
596 * ========================================Apply Boundary Conditions=====================================
597 *
598 * @code
599 *   {
600 *   constraints.clear();
601 *   std::vector<bool> component_maskP(dim + 1, false);
602 *   component_maskP[dim] = true;
603 *   DoFTools::make_hanging_node_constraints(dof_handler, constraints);
605 *   BoundaryValuesP<dim>(), constraints, component_maskP);
606 *   }
607 *   {
608 *   std::set<types::boundary_id> no_normal_flux_boundaries;
609 *   no_normal_flux_boundaries.insert(99);
611 *   no_normal_flux_boundaries, constraints);
612 *   }
613 *  
614 *   constraints.close();
615 *  
616 *   std::vector<types::global_dof_index> dofs_per_block = DoFTools::count_dofs_per_fe_block(dof_handler, block_component);
617 *   n_u = dofs_per_block[0];
618 *   n_p = dofs_per_block[1];
619 *  
620 *   std::cout << " Number of active cells: " << triangulation.n_active_cells()
621 *   << std::endl << " Number of degrees of freedom: "
622 *   << dof_handler.n_dofs() << " (" << n_u << '+' << n_p << ')'
623 *   << std::endl;
624 *  
625 *   {
626 *   BlockDynamicSparsityPattern csp(2, 2);
627 *  
628 *   csp.block(0, 0).reinit(n_u, n_u);
629 *   csp.block(1, 0).reinit(n_p, n_u);
630 *   csp.block(0, 1).reinit(n_u, n_p);
631 *   csp.block(1, 1).reinit(n_p, n_p);
632 *  
633 *   csp.collect_sizes();
634 *  
635 *   DoFTools::make_sparsity_pattern(dof_handler, csp, constraints, false);
636 *   sparsity_pattern.copy_from(csp);
637 *   }
638 *  
639 *   system_matrix.reinit(sparsity_pattern);
640 *  
641 *   solution.reinit(2);
642 *   solution.block(0).reinit(n_u);
643 *   solution.block(1).reinit(n_p);
644 *   solution.collect_sizes();
645 *  
646 *   system_rhs.reinit(2);
647 *   system_rhs.block(0).reinit(n_u);
648 *   system_rhs.block(1).reinit(n_p);
649 *   system_rhs.collect_sizes();
650 *  
651 *   }
652 *  
653 * @endcode
654 *
655 * Viscosity and Shear modulus functions
656 *
657
658 *
659 *
660 * @code
661 *   template<int dim>
662 *   class Rheology
663 *   {
664 *   public:
665 *   double get_eta(double &r, double &z);
666 *   double get_G(unsigned int mat_id);
667 *  
668 *   private:
669 *   std::vector<double> get_manual_eta_profile();
670 *  
671 *   };
672 *  
673 *   template<int dim>
674 *   std::vector<double> Rheology<dim>::get_manual_eta_profile()
675 *   {
676 *   vector<double> etas;
677 *  
678 *   for (unsigned int i=0; i < system_parameters::sizeof_depths_eta; ++i)
679 *   {
680 *   etas.push_back(system_parameters::depths_eta[i]);
681 *   etas.push_back(system_parameters::eta_kinks[i]);
682 *   }
683 *   return etas;
684 *   }
685 *  
686 *   template<int dim>
687 *   double Rheology<dim>::get_eta(double &r, double &z)
688 *   {
689 * @endcode
690 *
691 * compute local depth
692 *
693 * @code
694 *   double ecc = system_parameters::q_axes[0] / system_parameters::p_axes[0];
695 *   double Rminusr = system_parameters::q_axes[0] - system_parameters::p_axes[0];
696 *   double approx_a = std::sqrt(r * r + z * z * ecc * ecc);
697 *   double group1 = r * r + z * z - Rminusr * Rminusr;
698 *  
699 *   double a0 = approx_a;
700 *   double error = 10000;
701 * @endcode
702 *
703 * While loop finds the a axis of the "isodepth" ellipse for which the input point is on the surface.
704 * An "isodepth" ellipse is defined as one whose axes a,b are related to the global axes A, B by: A-h = B-h
705 *
706
707 *
708 *
709 * @code
710 *   if ((r > system_parameters::q_axes[0] - system_parameters::depths_eta.back()) ||
711 *   (z > system_parameters::p_axes[0] - system_parameters::depths_eta.back()))
712 *   {
713 *  
714 *   double eps = 10.0;
715 *   while (error >= eps)
716 *   {
717 *   double a02 = a0 * a0;
718 *   double a03 = a0 * a02;
719 *   double a04 = a0 * a03;
720 *   double fofa = a04 - (2 * Rminusr * a03) - (group1 * a02)
721 *   + (2 * r * r * Rminusr * a0) - (r * r * Rminusr * Rminusr);
722 *   double fprimeofa = 4 * a03 - (6 * Rminusr * a02) - (2 * group1 * a0)
723 *   + (2 * r * r * Rminusr);
724 *   double deltaa = -fofa / fprimeofa;
725 *   a0 += deltaa;
726 *   error = std::abs(deltaa);
727 * @endcode
728 *
729 * cout << "error = " << error << endl;
730 *
731 * @code
732 *   }
733 *   }
734 *   else
735 *   {
736 *   a0 = 0.0;
737 *   }
738 *  
739 *   double local_depth = system_parameters::q_axes[0] - a0;
740 *   if (local_depth < 0)
741 *   local_depth = 0;
742 *  
743 *   if (local_depth > system_parameters::depths_eta.back())
744 *   {
745 *   if (system_parameters::eta_kinks.back() < system_parameters::eta_floor)
746 *   return system_parameters::eta_floor;
747 *   else if (system_parameters::eta_kinks.back() > system_parameters::eta_ceiling)
748 *   return system_parameters::eta_ceiling;
749 *   else
750 *   return system_parameters::eta_kinks.back();
751 *   }
752 *  
753 *   std::vector<double> viscosity_function = get_manual_eta_profile();
754 *  
755 *   unsigned int n_visc_kinks = viscosity_function.size() / 2;
756 *  
757 * @endcode
758 *
759 * find the correct interval to do the interpolation in
760 *
761 * @code
762 *   int n_minus_one = -1;
763 *   for (unsigned int n = 1; n <= n_visc_kinks; ++n)
764 *   {
765 *   unsigned int ndeep = 2 * n - 2;
766 *   unsigned int nshallow = 2 * n;
767 *   if (local_depth >= viscosity_function[ndeep] && local_depth <= viscosity_function[nshallow])
768 *   n_minus_one = ndeep;
769 *   }
770 *  
771 * @endcode
772 *
773 * find the viscosity interpolation
774 *
775 * @code
776 *   if (n_minus_one == -1)
777 *   return system_parameters::eta_ceiling;
778 *   else
779 *   {
780 *   double visc_exponent =
781 *   (viscosity_function[n_minus_one]
782 *   - local_depth)
783 *   / (viscosity_function[n_minus_one]
784 *   - viscosity_function[n_minus_one + 2]);
785 *   double visc_base = viscosity_function[n_minus_one + 3]
786 *   / viscosity_function[n_minus_one + 1];
787 * @endcode
788 *
789 * This is the true viscosity given the thermal profile
790 *
791 * @code
792 *   double true_eta = viscosity_function[n_minus_one + 1] * std::pow(visc_base, visc_exponent);
793 *  
794 * @endcode
795 *
796 * Implement latitude-dependence viscosity
797 *
798 * @code
799 *   if (system_parameters::lat_dependence)
800 *   {
801 *   double lat = 180 / PI * std::atan(z / r);
802 *   if (lat > 80)
803 *   lat = 80;
804 *   double T_eq = 155;
805 *   double T_surf = T_eq * std::sqrt( std::sqrt( std::cos( PI / 180 * lat ) ) );
806 *   double taper_depth = 40000;
807 *   double surface_taper = (taper_depth - local_depth) / taper_depth;
808 *   if (surface_taper < 0)
809 *   surface_taper = 0;
810 *   double log_eta_contrast = surface_taper * system_parameters::eta_Ea * 52.5365 * (T_eq - T_surf) / T_eq / T_surf;
811 *   true_eta *= std::pow(10, log_eta_contrast);
812 *   }
813 *  
814 *   if (true_eta > system_parameters::eta_ceiling)
815 *   return system_parameters::eta_ceiling;
816 *   else if (true_eta < system_parameters::eta_floor)
817 *   return system_parameters::eta_floor;
818 *   else
819 *   return true_eta;
820 *   }
821 *   }
822 *  
823 *  
824 *   template<int dim>
825 *   double Rheology<dim>::get_G(unsigned int mat_id)
826 *   {
827 *   return system_parameters::G[mat_id];
828 *   }
829 *  
830 *  
831 * @endcode
832 *
833 * Initialize the eta and G parts of the quadrature_point_history object
834 *
835
836 *
837 *
838 * @code
839 *   template<int dim>
840 *   void StokesProblem<dim>::initialize_eta_and_G()
841 *   {
842 *   FEValues<dim> fe_values(fe, quadrature_formula, update_quadrature_points);
843 *  
844 *   const unsigned int n_q_points = quadrature_formula.size();
845 *   Rheology<dim> rheology;
846 *  
847 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
848 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
849 *   {
850 *   PointHistory<dim> *local_quadrature_points_history =
851 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
852 *   Assert(
853 *   local_quadrature_points_history >= &quadrature_point_history.front(),
854 *   ExcInternalError());
855 *   Assert(
856 *   local_quadrature_points_history < &quadrature_point_history.back(),
857 *   ExcInternalError());
858 *   fe_values.reinit(cell);
859 *  
860 *   for (unsigned int q = 0; q < n_q_points; ++q)
861 *   {
862 *  
863 *   double r_value = fe_values.quadrature_point(q)[0];
864 *   double z_value = fe_values.quadrature_point(q)[1];
865 *  
866 * @endcode
867 *
868 * defines local viscosity
869 *
870 * @code
871 *   double local_viscosity = 0;
872 *  
873 *   local_viscosity = rheology.get_eta(r_value, z_value);
874 *  
875 *   local_quadrature_points_history[q].first_eta = local_viscosity;
876 *   local_quadrature_points_history[q].new_eta = local_viscosity;
877 *  
878 * @endcode
879 *
880 * defines local shear modulus
881 *
882 * @code
883 *   double local_G = 0;
884 *  
885 *   unsigned int mat_id = cell->material_id();
886 *  
887 *   local_G = rheology.get_G(mat_id);
888 *   local_quadrature_points_history[q].G = local_G;
889 *  
890 * @endcode
891 *
892 * initializes the phi-phi stress
893 *
894 * @code
895 *   local_quadrature_points_history[q].old_phiphi_stress = 0;
896 *   }
897 *   }
898 *   }
899 *  
900 * @endcode
901 *
902 * ====================== ASSEMBLE THE SYSTEM ======================
903 *
904
905 *
906 *
907 * @code
908 *   template<int dim>
909 *   void StokesProblem<dim>::assemble_system()
910 *   {
911 *   system_matrix = 0;
912 *   system_rhs = 0;
913 *  
914 *   FEValues<dim> fe_values(fe, quadrature_formula,
916 *   | update_gradients);
917 *  
918 *   const unsigned int dofs_per_cell = fe.dofs_per_cell;
919 *  
920 *   const unsigned int n_q_points = quadrature_formula.size();
921 *  
922 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
923 *   Vector<double> local_rhs(dofs_per_cell);
924 *  
925 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
926 *  
927 * @endcode
928 *
929 * runs the gravity script function
930 *
931 * @code
932 *   const RightHandSide<dim> right_hand_side;
933 *  
934 *   A_Grav_namespace::AnalyticGravity<dim> *aGrav =
935 *   new A_Grav_namespace::AnalyticGravity<dim>;
936 *   std::vector<double> grav_parameters;
937 *   grav_parameters.push_back(system_parameters::q_axes[system_parameters::present_timestep * 2 + 0]);
938 *   grav_parameters.push_back(system_parameters::p_axes[system_parameters::present_timestep * 2 + 0]);
939 *   grav_parameters.push_back(system_parameters::q_axes[system_parameters::present_timestep * 2 + 1]);
940 *   grav_parameters.push_back(system_parameters::p_axes[system_parameters::present_timestep * 2 + 1]);
941 *   grav_parameters.push_back(system_parameters::rho[0]);
942 *   grav_parameters.push_back(system_parameters::rho[1]);
943 *  
944 *   std::cout << "Body parameters are: " ;
945 *   for (int i=0; i<6; ++i)
946 *   std::cout << grav_parameters[i] << " ";
947 *   std::cout << endl;
948 *  
949 *   aGrav->setup_vars(grav_parameters);
950 *  
951 *   std::vector<Vector<double> > rhs_values(n_q_points,
952 *   Vector<double>(dim + 1));
953 *  
954 *   const FEValuesExtractors::Vector velocities(0);
955 *   const FEValuesExtractors::Scalar pressure(dim);
956 *  
957 *   std::vector<SymmetricTensor<2, dim> > phi_grads_u(dofs_per_cell);
958 *   std::vector<double> div_phi_u(dofs_per_cell);
959 *   std::vector<Tensor<1, dim> > phi_u(dofs_per_cell);
960 *   std::vector<double> phi_p(dofs_per_cell);
961 *  
962 *   typename DoFHandler<dim>::active_cell_iterator cell =
963 *   dof_handler.begin_active(), first_cell = dof_handler.begin_active(),
964 *   endc = dof_handler.end();
965 *  
966 *   for (; cell != endc; ++cell)
967 *   {
968 *   PointHistory<dim> *local_quadrature_points_history =
969 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
970 *   Assert(
971 *   local_quadrature_points_history >= &quadrature_point_history.front(),
972 *   ExcInternalError());
973 *   Assert(
974 *   local_quadrature_points_history < &quadrature_point_history.back(),
975 *   ExcInternalError());
976 *  
977 *   double cell_area = cell->measure();
978 *  
979 *   if (cell_area<0)
980 *   append_physical_times(-1); // This writes final line to physical_times.txt if step terminates prematurely
981 *   AssertThrow(cell_area > 0
982 *   ,
983 *   ExcInternalError());
984 *  
985 *  
986 *   unsigned int m_id = cell->material_id();
987 *  
988 * @endcode
989 *
990 * initializes the rhs vector to the correct g values
991 *
992 * @code
993 *   fe_values.reinit(cell);
994 *   right_hand_side.vector_value_list(fe_values.get_quadrature_points(),
995 *   rhs_values, aGrav);
996 *  
997 *   std::vector<Vector<double> > new_viscosities(quadrature_formula.size(), Vector<double>(dim + 1));
998 *  
999 * @endcode
1000 *
1001 * Finds vertices where the radius is zero DIM
1002 *
1003 * @code
1004 *   bool is_singular = false;
1005 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
1006 *   if (cell->face(f)->center()[0] == 0)
1007 *   is_singular = true;
1008 *  
1009 *   if (is_singular == false || system_parameters::cylindrical == false)
1010 *   {
1011 *   local_matrix = 0;
1012 *   local_rhs = 0;
1013 *  
1014 * @endcode
1015 *
1016 * ===== outputs the local gravity
1017 *
1018 * @code
1019 *   std::vector<Point<dim> > quad_points_list(n_q_points);
1020 *   quad_points_list = fe_values.get_quadrature_points();
1021 *  
1022 *   if (plastic_iteration
1023 *   == (system_parameters::max_plastic_iterations - 1))
1024 *   {
1025 *   if (cell != first_cell)
1026 *   {
1027 *   std::ofstream fout("gravity_field.txt", std::ios::app);
1028 *   fout << quad_points_list[0] << " " << rhs_values[0];
1029 *   fout.close();
1030 *   }
1031 *   else
1032 *   {
1033 *   std::ofstream fout("gravity_field.txt");
1034 *   fout << quad_points_list[0] << " " << rhs_values[0];
1035 *   fout.close();
1036 *   }
1037 *   }
1038 *  
1039 *   for (unsigned int q = 0; q < n_q_points; ++q)
1040 *   {
1041 *   SymmetricTensor<2, dim> &old_stress =
1042 *   local_quadrature_points_history[q].old_stress;
1043 *   double &local_old_phiphi_stress =
1044 *   local_quadrature_points_history[q].old_phiphi_stress;
1045 *   double r_value = fe_values.quadrature_point(q)[0];
1046 *  
1047 * @endcode
1048 *
1049 * get local density based on mat id
1050 *
1051 * @code
1052 *   double local_density = system_parameters::rho[m_id];
1053 *  
1054 * @endcode
1055 *
1056 * defines local viscosities
1057 *
1058 * @code
1059 *   double local_viscosity = 0;
1060 *   if (plastic_iteration == 0)
1061 *   local_viscosity = local_quadrature_points_history[q].first_eta;
1062 *   else
1063 *   local_viscosity = local_quadrature_points_history[q].new_eta;
1064 *  
1065 * @endcode
1066 *
1067 * Define the local viscoelastic constants
1068 *
1069 * @code
1070 *   double local_eta_ve = 2
1071 *   / ((1 / local_viscosity)
1072 *   + (1 / local_quadrature_points_history[q].G
1073 *   / system_parameters::current_time_interval));
1074 *   double local_chi_ve = 1
1075 *   / (1
1076 *   + (local_quadrature_points_history[q].G
1077 *   * system_parameters::current_time_interval
1078 *   / local_viscosity));
1079 *  
1080 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
1081 *   {
1082 *   phi_grads_u[k] = fe_values[velocities].symmetric_gradient(k,
1083 *   q);
1084 *   div_phi_u[k] = (fe_values[velocities].divergence(k, q));
1085 *   phi_u[k] = (fe_values[velocities].value(k, q));
1086 *   if (system_parameters::cylindrical == true)
1087 *   {
1088 *   div_phi_u[k] *= (r_value);
1089 *   div_phi_u[k] += (phi_u[k][0]);
1090 *   }
1091 *   phi_p[k] = fe_values[pressure].value(k, q);
1092 *   }
1093 *  
1094 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1095 *   {
1096 *   for (unsigned int j = 0; j <= i; ++j)
1097 *   {
1098 *   if (system_parameters::cylindrical == true)
1099 *   {
1100 *   local_matrix(i, j) += (phi_grads_u[i]
1101 *   * phi_grads_u[j] * 2 * local_eta_ve
1102 *   * r_value
1103 *   + 2 * phi_u[i][0] * phi_u[j][0]
1104 *   * local_eta_ve / r_value
1105 *   - div_phi_u[i] * phi_p[j]
1106 *   * system_parameters::pressure_scale
1107 *   - phi_p[i] * div_phi_u[j]
1108 *   * system_parameters::pressure_scale
1109 *   + phi_p[i] * phi_p[j] * r_value
1110 *   * system_parameters::pressure_scale)
1111 *   * fe_values.JxW(q);
1112 *   }
1113 *   else
1114 *   {
1115 *   local_matrix(i, j) += (phi_grads_u[i]
1116 *   * phi_grads_u[j] * 2 * local_eta_ve
1117 *   - div_phi_u[i] * phi_p[j]
1118 *   * system_parameters::pressure_scale
1119 *   - phi_p[i] * div_phi_u[j]
1120 *   * system_parameters::pressure_scale
1121 *   + phi_p[i] * phi_p[j]) * fe_values.JxW(q);
1122 *   }
1123 *   }
1124 *   if (system_parameters::cylindrical == true)
1125 *   {
1126 *   const unsigned int component_i =
1127 *   fe.system_to_component_index(i).first;
1128 *   local_rhs(i) += (fe_values.shape_value(i, q)
1129 *   * rhs_values[q](component_i) * r_value
1130 *   * local_density
1131 *   - local_chi_ve * phi_grads_u[i] * old_stress
1132 *   * r_value
1133 *   - local_chi_ve * phi_u[i][0]
1134 *   * local_old_phiphi_stress)
1135 *   * fe_values.JxW(q);
1136 *   }
1137 *   else
1138 *   {
1139 *   const unsigned int component_i =
1140 *   fe.system_to_component_index(i).first;
1141 *   local_rhs(i) += fe_values.shape_value(i, q)
1142 *   * rhs_values[q](component_i) * fe_values.JxW(q)
1143 *   * local_density;
1144 *   }
1145 *   }
1146 *   }
1147 *   } // end of non-singular
1148 *   else
1149 *   {
1150 *   local_matrix = 0;
1151 *   local_rhs = 0;
1152 *  
1153 * @endcode
1154 *
1155 * ===== outputs the local gravity
1156 *
1157 * @code
1158 *   std::vector<Point<dim> > quad_points_list(n_q_points);
1159 *   quad_points_list = fe_values.get_quadrature_points();
1160 *  
1161 *   for (unsigned int q = 0; q < n_q_points; ++q)
1162 *   {
1163 *   const SymmetricTensor<2, dim> &old_stress =
1164 *   local_quadrature_points_history[q].old_stress;
1165 *   double &local_old_phiphi_stress =
1166 *   local_quadrature_points_history[q].old_phiphi_stress;
1167 *   double r_value = fe_values.quadrature_point(q)[0];
1168 *  
1169 *   double local_density = system_parameters::rho[m_id];
1170 *  
1171 * @endcode
1172 *
1173 * defines local viscosities
1174 *
1175 * @code
1176 *   double local_viscosity = 0;
1177 *   if (plastic_iteration == 0)
1178 *   {
1179 *   local_viscosity = local_quadrature_points_history[q].first_eta;
1180 *   }
1181 *   else
1182 *   local_viscosity = local_quadrature_points_history[q].new_eta;
1183 *  
1184 * @endcode
1185 *
1186 * Define the local viscoelastic constants
1187 *
1188 * @code
1189 *   double local_eta_ve = 2
1190 *   / ((1 / local_viscosity)
1191 *   + (1 / local_quadrature_points_history[q].G
1192 *   / system_parameters::current_time_interval));
1193 *   double local_chi_ve = 1
1194 *   / (1
1195 *   + (local_quadrature_points_history[q].G
1196 *   * system_parameters::current_time_interval
1197 *   / local_viscosity));
1198 *  
1199 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
1200 *   {
1201 *   phi_grads_u[k] = fe_values[velocities].symmetric_gradient(k,
1202 *   q);
1203 *   div_phi_u[k] = (fe_values[velocities].divergence(k, q));
1204 *   phi_u[k] = (fe_values[velocities].value(k, q));
1205 *   if (system_parameters::cylindrical == true)
1206 *   {
1207 *   div_phi_u[k] *= (r_value);
1208 *   div_phi_u[k] += (phi_u[k][0]);
1209 *   }
1210 *   phi_p[k] = fe_values[pressure].value(k, q);
1211 *   }
1212 *  
1213 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1214 *   {
1215 *   for (unsigned int j = 0; j <= i; ++j)
1216 *   {
1217 *   if (system_parameters::cylindrical == true)
1218 *   {
1219 *   local_matrix(i, j) += (phi_grads_u[i]
1220 *   * phi_grads_u[j] * 2 * local_eta_ve
1221 *   * r_value
1222 *   + 2 * phi_u[i][0] * phi_u[j][0]
1223 *   * local_eta_ve / r_value
1224 *   - div_phi_u[i] * phi_p[j]
1225 *   * system_parameters::pressure_scale
1226 *   - phi_p[i] * div_phi_u[j]
1227 *   * system_parameters::pressure_scale
1228 *   + phi_p[i] * phi_p[j] * r_value
1229 *   * system_parameters::pressure_scale)
1230 *   * fe_values.JxW(q);
1231 *   }
1232 *   else
1233 *   {
1234 *   local_matrix(i, j) += (phi_grads_u[i]
1235 *   * phi_grads_u[j] * 2 * local_eta_ve
1236 *   - div_phi_u[i] * phi_p[j]
1237 *   * system_parameters::pressure_scale
1238 *   - phi_p[i] * div_phi_u[j]
1239 *   * system_parameters::pressure_scale
1240 *   + phi_p[i] * phi_p[j]) * fe_values.JxW(q);
1241 *   }
1242 *   }
1243 *   if (system_parameters::cylindrical == true)
1244 *   {
1245 *   const unsigned int component_i =
1246 *   fe.system_to_component_index(i).first;
1247 *   local_rhs(i) += (fe_values.shape_value(i, q)
1248 *   * rhs_values[q](component_i) * r_value
1249 *   * local_density
1250 *   - local_chi_ve * phi_grads_u[i] * old_stress
1251 *   * r_value
1252 *   - local_chi_ve * phi_u[i][0]
1253 *   * local_old_phiphi_stress)
1254 *   * fe_values.JxW(q);
1255 *   }
1256 *   else
1257 *   {
1258 *   const unsigned int component_i =
1259 *   fe.system_to_component_index(i).first;
1260 *   local_rhs(i) += fe_values.shape_value(i, q)
1261 *   * rhs_values[q](component_i) * fe_values.JxW(q)
1262 *   * local_density;
1263 *   }
1264 *   }
1265 *   }
1266 *   } // end of singular
1267 *  
1268 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1269 *   for (unsigned int j = i + 1; j < dofs_per_cell; ++j)
1270 *   local_matrix(i, j) = local_matrix(j, i);
1271 *  
1272 *   cell->get_dof_indices(local_dof_indices);
1273 *   constraints.distribute_local_to_global(local_matrix, local_rhs,
1274 *   local_dof_indices, system_matrix, system_rhs);
1275 *   }
1276 *  
1277 *   std::cout << " Computing preconditioner..." << std::endl << std::flush;
1278 *  
1279 *   A_preconditioner = std::shared_ptr<
1280 *   typename InnerPreconditioner<dim>::type>(
1281 *   new typename InnerPreconditioner<dim>::type());
1282 *   A_preconditioner->initialize(system_matrix.block(0, 0),
1283 *   typename InnerPreconditioner<dim>::type::AdditionalData());
1284 *  
1285 *   delete aGrav;
1286 *   }
1287 *  
1288 * @endcode
1289 *
1290 * ====================== SOLVER ======================
1291 *
1292
1293 *
1294 *
1295 * @code
1296 *   template<int dim>
1297 *   void StokesProblem<dim>::solve()
1298 *   {
1299 *   const InverseMatrix<SparseMatrix<double>,
1300 *   typename InnerPreconditioner<dim>::type> A_inverse(
1301 *   system_matrix.block(0, 0), *A_preconditioner);
1302 *   Vector<double> tmp(solution.block(0).size());
1303 *  
1304 *   {
1305 *   Vector<double> schur_rhs(solution.block(1).size());
1306 *   A_inverse.vmult(tmp, system_rhs.block(0));
1307 *   system_matrix.block(1, 0).vmult(schur_rhs, tmp);
1308 *   schur_rhs -= system_rhs.block(1);
1309 *  
1310 *   SchurComplement<typename InnerPreconditioner<dim>::type> schur_complement(
1311 *   system_matrix, A_inverse);
1312 *  
1313 *   int n_iterations = system_parameters::iteration_coefficient
1314 *   * solution.block(1).size();
1315 *   double tolerance_goal = system_parameters::tolerance_coefficient
1316 *   * schur_rhs.l2_norm();
1317 *  
1318 *   SolverControl solver_control(n_iterations, tolerance_goal);
1319 *   SolverCG<> cg(solver_control);
1320 *  
1321 *   std::cout << "\nMax iterations and tolerance are: " << n_iterations
1322 *   << " and " << tolerance_goal << std::endl;
1323 *  
1324 *   SparseILU<double> preconditioner;
1325 *   preconditioner.initialize(system_matrix.block(1, 1),
1327 *  
1328 *   InverseMatrix<SparseMatrix<double>, SparseILU<double> > m_inverse(
1329 *   system_matrix.block(1, 1), preconditioner);
1330 *  
1331 *   cg.solve(schur_complement, solution.block(1), schur_rhs, m_inverse);
1332 *  
1333 *   constraints.distribute(solution);
1334 *  
1335 *  
1336 *   std::cout << " " << solver_control.last_step()
1337 *   << " outer CG Schur complement iterations for pressure"
1338 *   << std::endl;
1339 *   }
1340 *  
1341 *   {
1342 *   system_matrix.block(0, 1).vmult(tmp, solution.block(1));
1343 *   tmp *= -1;
1344 *   tmp += system_rhs.block(0);
1345 *  
1346 *   A_inverse.vmult(solution.block(0), tmp);
1347 *   constraints.distribute(solution);
1348 *   solution.block(1) *= (system_parameters::pressure_scale);
1349 *   }
1350 *   }
1351 *  
1352 * @endcode
1353 *
1354 * ====================== OUTPUT RESULTS ======================
1355 *
1356 * @code
1357 *   template<int dim>
1358 *   void StokesProblem<dim>::output_results() const
1359 *   {
1360 *   std::vector < std::string > solution_names(dim, "velocity");
1361 *   solution_names.push_back("pressure");
1362 *  
1363 *   std::vector<DataComponentInterpretation::DataComponentInterpretation> data_component_interpretation(
1365 *   data_component_interpretation.push_back(
1367 *  
1368 *   DataOut<dim> data_out;
1369 *   data_out.attach_dof_handler(dof_handler);
1370 *   data_out.add_data_vector(solution, solution_names,
1371 *   DataOut<dim>::type_dof_data, data_component_interpretation);
1372 *   data_out.build_patches();
1373 *  
1374 *   std::ostringstream filename;
1375 *   if (system_parameters::present_timestep < system_parameters::initial_elastic_iterations)
1376 *   {
1377 *   filename << system_parameters::output_folder << "/time"
1378 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1379 *   << "_elastic_displacements" << ".txt";
1380 *   }
1381 *   else
1382 *   {
1383 *   filename << system_parameters::output_folder << "/time"
1384 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1385 *   << "_flow" << Utilities::int_to_string(plastic_iteration, 2) << ".txt";
1386 *   }
1387 *  
1388 *   std::ofstream output(filename.str().c_str());
1389 *   data_out.write_gnuplot(output);
1390 *   }
1391 *  
1392 * @endcode
1393 *
1394 * ====================== FIND AND WRITE TO FILE THE STRESS TENSOR; IMPLEMENT PLASTICITY ======================
1395 *
1396
1397 *
1398 *
1399 * @code
1400 *   template<int dim>
1401 *   void StokesProblem<dim>::solution_stesses()
1402 *   {
1403 * @endcode
1404 *
1405 * note most of this section only works with dim=2
1406 *
1407
1408 *
1409 * name the output text files
1410 *
1411 * @code
1412 *   std::ostringstream stress_output;
1413 *   stress_output << system_parameters::output_folder << "/time"
1414 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1415 *   << "_principalstresses" << Utilities::int_to_string(plastic_iteration, 2)
1416 *   << ".txt";
1417 *   std::ofstream fout_snew(stress_output.str().c_str());
1418 *   fout_snew.close();
1419 *  
1420 *   std::ostringstream stresstensor_output;
1421 *   stresstensor_output << system_parameters::output_folder << "/time"
1422 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1423 *   << "_stresstensor" << Utilities::int_to_string(plastic_iteration, 2)
1424 *   << ".txt";
1425 *   std::ofstream fout_sfull(stresstensor_output.str().c_str());
1426 *   fout_sfull.close();
1427 *  
1428 *   std::ostringstream failed_cells_output;
1429 *   failed_cells_output << system_parameters::output_folder << "/time"
1430 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1431 *   << "_failurelocations" << Utilities::int_to_string(plastic_iteration, 2)
1432 *   << ".txt";
1433 *   std::ofstream fout_failed_cells(failed_cells_output.str().c_str());
1434 *   fout_failed_cells.close();
1435 *  
1436 *   std::ostringstream plastic_eta_output;
1437 *   plastic_eta_output << system_parameters::output_folder << "/time"
1438 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1439 *   << "_viscositiesreg" << Utilities::int_to_string(plastic_iteration, 2)
1440 *   << ".txt";
1441 *   std::ofstream fout_vrnew(plastic_eta_output.str().c_str());
1442 *   fout_vrnew.close();
1443 *  
1444 *   std::ostringstream initial_eta_output;
1445 *   if (plastic_iteration == 0)
1446 *   {
1447 *   initial_eta_output << system_parameters::output_folder << "/time"
1448 *   << Utilities::int_to_string(system_parameters::present_timestep, 2)
1449 *   << "_baseviscosities.txt";
1450 *   std::ofstream fout_baseeta(initial_eta_output.str().c_str());
1451 *   fout_baseeta.close();
1452 *   }
1453 *  
1454 *   std::cout << "Running stress calculations for plasticity iteration "
1455 *   << plastic_iteration << "...\n";
1456 *  
1457 * @endcode
1458 *
1459 * This makes the set of points at which the stress tensor is calculated
1460 *
1461 * @code
1462 *   std::vector<Point<dim> > points_list(0);
1463 *   std::vector<unsigned int> material_list(0);
1464 *   typename DoFHandler<dim>::active_cell_iterator cell =
1465 *   dof_handler.begin_active(), endc = dof_handler.end();
1466 * @endcode
1467 *
1468 * This loop gets the gradients of the velocity field and saves it in the tensor_gradient_? objects DIM
1469 *
1470 * @code
1471 *   for (; cell != endc; ++cell)
1472 *   {
1473 *   points_list.push_back(cell->center());
1474 *   material_list.push_back(cell->material_id());
1475 *   }
1476 * @endcode
1477 *
1478 * Make the FEValues object to evaluate values and derivatives at quadrature points
1479 *
1480 * @code
1481 *   FEValues<dim> fe_values(fe, quadrature_formula,
1483 *  
1484 * @endcode
1485 *
1486 * Make the object that will hold the velocities and velocity gradients at the quadrature points
1487 *
1488 * @code
1489 *   std::vector < std::vector<Tensor<1, dim> >> velocity_grads(quadrature_formula.size(),
1490 *   std::vector < Tensor<1, dim> > (dim + 1));
1491 *   std::vector<Vector<double> > velocities(quadrature_formula.size(),
1492 *   Vector<double>(dim + 1));
1493 * @endcode
1494 *
1495 * Make the object to find rheology
1496 *
1497 * @code
1498 *   Rheology<dim> rheology;
1499 *  
1500 * @endcode
1501 *
1502 * Write the solution flow velocity and derivative for each cell
1503 *
1504 * @code
1505 *   std::vector<Vector<double> > vector_values(0);
1506 *   std::vector < std::vector<Tensor<1, dim> > > gradient_values(0);
1507 *   std::vector<bool> failing_cells;
1508 * @endcode
1509 *
1510 * Write the stresses from the previous step into vectors
1511 *
1512 * @code
1513 *   std::vector<SymmetricTensor<2, dim>> old_stress;
1514 *   std::vector<double> old_phiphi_stress;
1515 *   std::vector<double> cell_Gs;
1516 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
1517 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
1518 *   {
1519 * @endcode
1520 *
1521 * Makes pointer to data in quadrature_point_history
1522 *
1523 * @code
1524 *   PointHistory<dim> *local_quadrature_points_history =
1525 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
1526 *  
1527 *   fe_values.reinit(cell);
1528 *   fe_values.get_function_gradients(solution, velocity_grads);
1529 *   fe_values.get_function_values(solution, velocities);
1530 *   Vector<double> current_cell_velocity(dim+1);
1531 *   std::vector<Tensor<1, dim>> current_cell_grads(dim+1);
1532 *   SymmetricTensor<2, dim> current_cell_old_stress;
1533 *   current_cell_old_stress = 0;
1534 *   double current_cell_old_phiphi_stress = 0;
1535 *   double cell_area = 0;
1536 *  
1537 * @endcode
1538 *
1539 * Averages across each cell to find mean velocities, gradients, and old stresses
1540 *
1541 * @code
1542 *   for (unsigned int q = 0; q < quadrature_formula.size(); ++q)
1543 *   {
1544 *   cell_area += fe_values.JxW(q);
1545 *   velocities[q] *= fe_values.JxW(q);
1546 *   current_cell_velocity += velocities[q];
1547 *   for (unsigned int i = 0; i < (dim+1); ++i)
1548 *   {
1549 *   velocity_grads[q][i] *= fe_values.JxW(q);
1550 *   current_cell_grads[i] += velocity_grads[q][i];
1551 *   }
1552 *   current_cell_old_stress += local_quadrature_points_history[q].old_stress * fe_values.JxW(q);
1553 *   current_cell_old_phiphi_stress += local_quadrature_points_history[q].old_phiphi_stress * fe_values.JxW(q);
1554 *   }
1555 *   current_cell_velocity /= cell_area;
1556 *   for (unsigned int i = 0; i < (dim+1); ++i)
1557 *   current_cell_grads[i] /= cell_area;
1558 *   current_cell_old_stress /= cell_area;
1559 *   current_cell_old_phiphi_stress /= cell_area;
1560 *  
1561 *   vector_values.push_back(current_cell_velocity);
1562 *   gradient_values.push_back(current_cell_grads);
1563 *   old_stress.push_back(current_cell_old_stress);
1564 *   old_phiphi_stress.push_back(current_cell_old_phiphi_stress);
1565 *  
1566 * @endcode
1567 *
1568 * Get cell shear modulus: assumes it's constant for the cell
1569 *
1570 * @code
1571 *   unsigned int mat_id = cell->material_id();
1572 *   double local_G = rheology.get_G(mat_id);
1573 *   cell_Gs.push_back(local_G);
1574 *   }
1575 *  
1576 * @endcode
1577 *
1578 * tracks where failure occurred
1579 *
1580 * @code
1581 *   std::vector<double> reduction_factor;
1582 *   unsigned int total_fails = 0;
1583 *   if (plastic_iteration == 0)
1584 *   cell_viscosities.resize(0);
1585 * @endcode
1586 *
1587 * loop across all the cells to find and adjust eta of failing cells
1588 *
1589 * @code
1590 *   for (unsigned int i = 0; i < triangulation.n_active_cells(); ++i)
1591 *   {
1592 *   double current_cell_viscosity = 0;
1593 *  
1594 * @endcode
1595 *
1596 * Fill viscosities vector, analytically if plastic_iteration == 0 and from previous viscosities for later iteration
1597 *
1598 * @code
1599 *   if (plastic_iteration == 0)
1600 *   {
1601 *   double local_viscosity;
1602 *   local_viscosity = rheology.get_eta(points_list[i][0], points_list[i][1]);
1603 *   current_cell_viscosity = local_viscosity;
1604 *   cell_viscosities.push_back(current_cell_viscosity);
1605 *   }
1606 *   else
1607 *   {
1608 *   current_cell_viscosity = cell_viscosities[i];
1609 *   }
1610 *  
1611 *  
1612 *   double cell_eta_ve = 2
1613 *   / ((1 / current_cell_viscosity)
1614 *   + (1 / cell_Gs[i]
1615 *   / system_parameters::current_time_interval));
1616 *   double cell_chi_ve = 1
1617 *   / (1
1618 *   + (cell_Gs[i]
1619 *   * system_parameters::current_time_interval
1620 *   / current_cell_viscosity));
1621 *  
1622 * @endcode
1623 *
1624 * find local pressure
1625 *
1626 * @code
1627 *   double cell_p = vector_values[i].operator()(2);
1628 * @endcode
1629 *
1630 * find stresses tensor
1631 * makes non-diagonalized local matrix A
1632 *
1633 * @code
1634 *   double sigma13 = 0.5
1635 *   * (gradient_values[i][0][1] + gradient_values[i][1][0]);
1636 *   mat A;
1637 *   A << gradient_values[i][0][0] << 0 << sigma13 << endr
1638 *   << 0 << vector_values[i].operator()(0) / points_list[i].operator()(0)<< 0 << endr
1639 *   << sigma13 << 0 << gradient_values[i][1][1] << endr;
1640 *   mat olddevstress;
1641 *   olddevstress << old_stress[i][0][0] << 0 << old_stress[i][0][1] << endr
1642 *   << 0 << old_phiphi_stress[i] << 0 << endr
1643 *   << old_stress[i][0][1] << 0 << old_stress[i][1][1] << endr;
1644 *   vec P;
1645 *   P << cell_p << cell_p << cell_p;
1646 *   mat Pmat = diagmat(P);
1647 *   mat B;
1648 *   B = (cell_eta_ve * A + cell_chi_ve * olddevstress) - Pmat;
1649 *  
1650 * @endcode
1651 *
1652 * finds principal stresses
1653 *
1654 * @code
1655 *   vec eigval;
1656 *   mat eigvec;
1657 *   eig_sym(eigval, eigvec, B);
1658 *   double sigma1 = -min(eigval);
1659 *   double sigma3 = -max(eigval);
1660 *  
1661 * @endcode
1662 *
1663 * Writes text files for principal stresses, full stress tensor, base viscosities
1664 *
1665 * @code
1666 *   std::ofstream fout_snew(stress_output.str().c_str(), std::ios::app);
1667 *   fout_snew << " " << sigma1 << " " << sigma3 << "\n";
1668 *   fout_snew.close();
1669 *  
1670 *   std::ofstream fout_sfull(stresstensor_output.str().c_str(), std::ios::app);
1671 *   fout_sfull << A(0,0) << " " << A(1,1) << " " << A(2,2) << " " << A(0,2) << "\n";
1672 *   fout_sfull.close();
1673 *  
1674 *   if (plastic_iteration == 0)
1675 *   {
1676 *   std::ofstream fout_baseeta(initial_eta_output.str().c_str(), std::ios::app);
1677 *   fout_baseeta << points_list[i]<< " " << current_cell_viscosity << "\n";
1678 *   fout_baseeta.close();
1679 *   }
1680 *  
1681 * @endcode
1682 *
1683 * Finds adjusted effective viscosity
1684 *
1685 * @code
1686 *   if (system_parameters::plasticity_on)
1687 *   {
1688 *   if (system_parameters::failure_criterion == 0) //Apply Byerlee's rule
1689 *   {
1690 *   if (sigma1 >= 5 * sigma3) // this guarantees that viscosities only go down, never up
1691 *   {
1692 *   failing_cells.push_back(true);
1693 *   double temp_reductionfactor = 1;
1694 *   if (sigma3 < 0)
1695 *   temp_reductionfactor = 100;
1696 *   else
1697 *   temp_reductionfactor = 1.9 * sigma1 / 5 / sigma3;
1698 *  
1699 *   reduction_factor.push_back(temp_reductionfactor);
1700 *   total_fails++;
1701 *  
1702 * @endcode
1703 *
1704 * Text file of all failure locations
1705 *
1706 * @code
1707 *   std::ofstream fout_failed_cells(failed_cells_output.str().c_str(), std::ios::app);
1708 *   fout_failed_cells << points_list[i] << "\n";
1709 *   fout_failed_cells.close();
1710 *   }
1711 *   else
1712 *   {
1713 *   reduction_factor.push_back(1);
1714 *   failing_cells.push_back(false);
1715 *   }
1716 *   }
1717 *   else
1718 *   {
1719 *   if (system_parameters::failure_criterion == 1) //Apply Schultz criterion for frozen sand, RMR=45
1720 *   {
1721 *   double temp_reductionfactor = 1;
1722 *   if (sigma3 < -114037)
1723 *   {
1724 * @endcode
1725 *
1726 * std::cout << " ext ";
1727 *
1728 * @code
1729 *   failing_cells.push_back(true);
1730 *   temp_reductionfactor = 10;
1731 *   reduction_factor.push_back(temp_reductionfactor);
1732 *   total_fails++;
1733 *  
1734 * @endcode
1735 *
1736 * Text file of all failure locations
1737 *
1738 * @code
1739 *   std::ofstream fout_failed_cells(failed_cells_output.str().c_str(), std::ios::app);
1740 *   fout_failed_cells << points_list[i] << "\n";
1741 *   fout_failed_cells.close();
1742 *   }
1743 *   else
1744 *   {
1745 *   double sigma_c = 160e6; //Unconfined compressive strength
1746 *   double yield_sigma1 = sigma3 + std::sqrt( (3.086 * sigma_c * sigma3) + (0.002 * sigma3 * sigma3) );
1747 *   if (sigma1 >= yield_sigma1)
1748 *   {
1749 * @endcode
1750 *
1751 * std::cout << " comp ";
1752 *
1753 * @code
1754 *   failing_cells.push_back(true);
1755 *   temp_reductionfactor = 1.0 * sigma1 / 5 / sigma3;
1756 *  
1757 *   reduction_factor.push_back(temp_reductionfactor);
1758 *   total_fails++;
1759 *  
1760 * @endcode
1761 *
1762 * Text file of all failure locations
1763 *
1764 * @code
1765 *   std::ofstream fout_failed_cells(failed_cells_output.str().c_str(), std::ios::app);
1766 *   fout_failed_cells << points_list[i] << "\n";
1767 *   fout_failed_cells.close();
1768 *   }
1769 *   else
1770 *   {
1771 *   reduction_factor.push_back(1);
1772 *   failing_cells.push_back(false);
1773 *   }
1774 *   }
1775 *   }
1776 *   else
1777 *   {
1778 *   std::cout << "Specified failure criterion not found\n";
1779 *   }
1780 *   }
1781 *   }
1782 *   else
1783 *   reduction_factor.push_back(1);
1784 *   }
1785 *  
1786 * @endcode
1787 *
1788 * If there are enough failed cells, update eta at all quadrature points and perform smoothing
1789 *
1790 * @code
1791 *   std::cout << " Number of failing cells: " << total_fails << "\n";
1792 *   double last_max_plasticity_double = last_max_plasticity;
1793 *   double total_fails_double = total_fails;
1794 *   double decrease_in_plasticity = ((last_max_plasticity_double - total_fails_double) / last_max_plasticity_double);
1795 *   if (plastic_iteration == 0)
1796 *   decrease_in_plasticity = 1;
1797 *   last_max_plasticity = total_fails;
1798 *   if (total_fails <= 100 || decrease_in_plasticity <= 0.2)
1799 *   {
1800 *   system_parameters::continue_plastic_iterations = false;
1801 *   for (unsigned int j=0; j < triangulation.n_active_cells(); ++j)
1802 *   {
1803 * @endcode
1804 *
1805 * Writes to file the undisturbed cell viscosities
1806 *
1807 * @code
1808 *   std::ofstream fout_vrnew(plastic_eta_output.str().c_str(), std::ios::app);
1809 *   fout_vrnew << " " << cell_viscosities[j] << "\n";
1810 *   fout_vrnew.close();
1811 *   }
1812 *   }
1813 *   else
1814 *   {
1815 *   quad_viscosities.resize(triangulation.n_active_cells());
1816 * @endcode
1817 *
1818 * Decrease the eta at quadrature points in failing cells
1819 *
1820 * @code
1821 *   unsigned int cell_no = 0;
1822 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
1823 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
1824 *   {
1825 *   fe_values.reinit(cell);
1826 * @endcode
1827 *
1828 * Make local_quadrature_points_history pointer to the cell data
1829 *
1830 * @code
1831 *   PointHistory<dim> *local_quadrature_points_history =
1832 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
1833 *   Assert(
1834 *   local_quadrature_points_history >= &quadrature_point_history.front(),
1835 *   ExcInternalError());
1836 *   Assert(
1837 *   local_quadrature_points_history < &quadrature_point_history.back(),
1838 *   ExcInternalError());
1839 *  
1840 *   quad_viscosities[cell_no].resize(quadrature_formula.size());
1841 *  
1842 *   for (unsigned int q = 0; q < quadrature_formula.size(); ++q)
1843 *   {
1844 *   if (plastic_iteration == 0)
1845 *   local_quadrature_points_history[q].new_eta = local_quadrature_points_history[q].first_eta;
1846 *   local_quadrature_points_history[q].new_eta /= reduction_factor[cell_no];
1847 * @endcode
1848 *
1849 * Prevents viscosities from dropping below the floor necessary for numerical stability
1850 *
1851 * @code
1852 *   if (local_quadrature_points_history[q].new_eta < system_parameters::eta_floor)
1853 *   local_quadrature_points_history[q].new_eta = system_parameters::eta_floor;
1854 *  
1855 *   quad_viscosities[cell_no][q].reinit(dim+1);
1856 *   for (unsigned int ii=0; ii<dim; ++ii)
1857 *   quad_viscosities[cell_no][q](ii) = fe_values.quadrature_point(q)[ii];
1858 *   quad_viscosities[cell_no][q](dim) = local_quadrature_points_history[q].new_eta;
1859 *   }
1860 *   cell_no++;
1861 *   }
1862 *   smooth_eta_field(failing_cells);
1863 *  
1864 * @endcode
1865 *
1866 * Writes to file the smoothed eta field (which is defined at each quadrature point) for each cell
1867 *
1868 * @code
1869 *   cell_no = 0;
1870 * @endcode
1871 *
1872 * cell_viscosities.resize(triangulation.n_active_cells(), 0);
1873 *
1874 * @code
1875 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
1876 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
1877 *   {
1878 *   if (failing_cells[cell_no])
1879 *   {
1880 *   fe_values.reinit(cell);
1881 * @endcode
1882 *
1883 * Averages across each cell to find mean eta
1884 *
1885 * @code
1886 *   double cell_area = 0;
1887 *   for (unsigned int q = 0; q < quadrature_formula.size(); ++q)
1888 *   {
1889 *   cell_area += fe_values.JxW(q);
1890 *   cell_viscosities[cell_no] += quad_viscosities[cell_no][q][dim] * fe_values.JxW(q);
1891 *   }
1892 *   cell_viscosities[cell_no] /= cell_area;
1893 *  
1894 * @endcode
1895 *
1896 * Writes to file
1897 *
1898 * @code
1899 *   std::ofstream fout_vrnew(plastic_eta_output.str().c_str(), std::ios::app);
1900 *   fout_vrnew << " " << cell_viscosities[cell_no] << "\n";
1901 *   fout_vrnew.close();
1902 *   }
1903 *   else
1904 *   {
1905 *   std::ofstream fout_vrnew(plastic_eta_output.str().c_str(), std::ios::app);
1906 *   fout_vrnew << " " << cell_viscosities[cell_no] << "\n";
1907 *   fout_vrnew.close();
1908 *   }
1909 *   cell_no++;
1910 *   }
1911 *   }
1912 *   }
1913 *  
1914 * @endcode
1915 *
1916 * ====================== SMOOTHES THE VISCOSITY FIELD AT ALL QUADRATURE POINTS ======================
1917 *
1918
1919 *
1920 *
1921 * @code
1922 *   template<int dim>
1923 *   void StokesProblem<dim>::smooth_eta_field(std::vector<bool> failing_cells)
1924 *   {
1925 *   std::cout << " Smoothing viscosity field...\n";
1926 *   unsigned int cell_no = 0;
1927 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
1928 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
1929 *   {
1930 *   if (failing_cells[cell_no])
1931 *   {
1932 *   FEValues<dim> fe_values(fe, quadrature_formula, update_quadrature_points);
1933 *   PointHistory<dim> *local_quadrature_points_history =
1934 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
1935 *  
1936 * @endcode
1937 *
1938 * Currently this algorithm does not permit refinement. To permit refinement, daughter cells of neighbors must be identified
1939 * Find pointers and indices of all cells within certain radius
1940 *
1941 * @code
1942 *   bool find_more_cells = true;
1943 *   std::vector<bool> cell_touched(triangulation.n_active_cells(), false);
1944 *   std::vector< TriaIterator< CellAccessor<dim> > > neighbor_cells;
1945 *   std::vector<int> neighbor_indices;
1946 *   unsigned int start_cell = 0; // Which cell in the neighbor_cells vector to start from
1947 *   int new_cells_found = 0;
1948 *   neighbor_cells.push_back(cell);
1949 *   neighbor_indices.push_back(cell_no);
1950 *   cell_touched[cell_no] = true;
1951 *   while (find_more_cells)
1952 *   {
1953 *   new_cells_found = 0;
1954 *   for (unsigned int i = start_cell; i<neighbor_cells.size(); ++i)
1955 *   {
1956 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
1957 *   {
1958 *   if (!neighbor_cells[i]->face(f)->at_boundary())
1959 *   {
1960 *   int test_cell_no = neighbor_cells[i]->neighbor_index(f);
1961 *   if (!cell_touched[test_cell_no])
1962 *   if (cell->center().distance(neighbor_cells[i]->neighbor(f)->center()) < 2 * system_parameters::smoothing_radius)
1963 *   {
1964 * @endcode
1965 *
1966 * What to do if another nearby cell is found that hasn't been found before
1967 *
1968 * @code
1969 *   neighbor_cells.push_back(neighbor_cells[i]->neighbor(f));
1970 *   neighbor_indices.push_back(test_cell_no);
1971 *   cell_touched[test_cell_no] = true;
1972 *   start_cell++;
1973 *   new_cells_found++;
1974 *   }
1975 *   }
1976 *   }
1977 *   }
1978 *   if (new_cells_found == 0)
1979 *   {
1980 *   find_more_cells = false;
1981 *   }
1982 *   else
1983 *   start_cell = neighbor_cells.size() - new_cells_found;
1984 *   }
1985 *  
1986 *   fe_values.reinit(cell);
1987 * @endcode
1988 *
1989 * Collect the viscosities at nearby quadrature points
1990 *
1991 * @code
1992 *   for (unsigned int q = 0; q < quadrature_formula.size(); ++q)
1993 *   {
1994 *   std::vector<double> nearby_etas_q;
1995 *   for (unsigned int i = 0; i<neighbor_indices.size(); ++i)
1996 *   for (unsigned int j=0; j<quadrature_formula.size(); ++j)
1997 *   {
1998 *   Point<dim> test_q;
1999 *   for (unsigned int d=0; d<dim; ++d)
2000 *   test_q(d) = quad_viscosities[neighbor_indices[i]][j][d];
2001 *   double qq_distance = fe_values.quadrature_point(q).distance(test_q);
2002 *   if (qq_distance < system_parameters::smoothing_radius)
2003 *   nearby_etas_q.push_back(quad_viscosities[neighbor_indices[i]][j][dim]);
2004 *   }
2005 * @endcode
2006 *
2007 * Write smoothed viscosities to quadrature_points_history; simple boxcar function is the smoothing kernel
2008 *
2009 * @code
2010 *   double mean_eta = 0;
2011 *   for (unsigned int l = 0; l<nearby_etas_q.size(); ++l)
2012 *   {
2013 *   mean_eta += nearby_etas_q[l];
2014 *   }
2015 *   mean_eta /= nearby_etas_q.size();
2016 *   local_quadrature_points_history[q].new_eta = mean_eta;
2017 * @endcode
2018 *
2019 * std::cout << local_quadrature_points_history[q].new_eta << " ";
2020 *
2021 * @code
2022 *   }
2023 *   }
2024 *   cell_no++;
2025 *   }
2026 *   }
2027 *  
2028 * @endcode
2029 *
2030 * ====================== SAVE STRESS TENSOR AT QUADRATURE POINTS ======================
2031 *
2032
2033 *
2034 *
2035 * @code
2036 *   template<int dim>
2037 *   void StokesProblem<dim>::update_quadrature_point_history()
2038 *   {
2039 *   std::cout << " Updating stress field...";
2040 *  
2041 *   FEValues<dim> fe_values(fe, quadrature_formula,
2042 *   update_values | update_gradients | update_quadrature_points);
2043 *  
2044 * @endcode
2045 *
2046 * Make the object that will hold the velocity gradients
2047 *
2048 * @code
2049 *   std::vector < std::vector<Tensor<1, dim> >> velocity_grads(quadrature_formula.size(),
2050 *   std::vector < Tensor<1, dim> > (dim + 1));
2051 *   std::vector<Vector<double> > velocities(quadrature_formula.size(),
2052 *   Vector<double>(dim + 1));
2053 *  
2054 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
2055 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
2056 *   {
2057 *   PointHistory<dim> *local_quadrature_points_history =
2058 *   reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
2059 *   Assert(
2060 *   local_quadrature_points_history >= &quadrature_point_history.front(),
2061 *   ExcInternalError());
2062 *   Assert(
2063 *   local_quadrature_points_history < &quadrature_point_history.back(),
2064 *   ExcInternalError());
2065 *  
2066 *   fe_values.reinit(cell);
2067 *   fe_values.get_function_gradients(solution, velocity_grads);
2068 *   fe_values.get_function_values(solution, velocities);
2069 *  
2070 *   for (unsigned int q = 0; q < quadrature_formula.size(); ++q)
2071 *   {
2072 * @endcode
2073 *
2074 * Define the local viscoelastic constants
2075 *
2076 * @code
2077 *   double local_eta_ve = 2
2078 *   / ((1 / local_quadrature_points_history[q].new_eta)
2079 *   + (1 / local_quadrature_points_history[q].G
2080 *   / system_parameters::current_time_interval));
2081 *   double local_chi_ve =
2082 *   1
2083 *   / (1
2084 *   + (local_quadrature_points_history[q].G
2085 *   * system_parameters::current_time_interval
2086 *   / local_quadrature_points_history[q].new_eta));
2087 *  
2088 * @endcode
2089 *
2090 * Compute new stress at each quadrature point
2091 *
2092 * @code
2093 *   SymmetricTensor<2, dim> new_stress;
2094 *   for (unsigned int i = 0; i < dim; ++i)
2095 *   new_stress[i][i] =
2096 *   local_eta_ve * velocity_grads[q][i][i]
2097 *   + local_chi_ve
2098 *   * local_quadrature_points_history[q].old_stress[i][i];
2099 *  
2100 *   for (unsigned int i = 0; i < dim; ++i)
2101 *   for (unsigned int j = i + 1; j < dim; ++j)
2102 *   new_stress[i][j] =
2103 *   local_eta_ve
2104 *   * (velocity_grads[q][i][j]
2105 *   + velocity_grads[q][j][i]) / 2
2106 *   + local_chi_ve
2107 *   * local_quadrature_points_history[q].old_stress[i][j];
2108 *  
2109 * @endcode
2110 *
2111 * Rotate new stress
2112 *
2113 * @code
2114 *   AuxFunctions<dim> rotation_object;
2115 *   const Tensor<2, dim> rotation = rotation_object.get_rotation_matrix(
2116 *   velocity_grads[q]);
2117 *   const SymmetricTensor<2, dim> rotated_new_stress = symmetrize(
2118 *   transpose(rotation)
2119 *   * static_cast<Tensor<2, dim> >(new_stress)
2120 *   * rotation);
2121 *   local_quadrature_points_history[q].old_stress = rotated_new_stress;
2122 *  
2123 * @endcode
2124 *
2125 * For axisymmetric case, make the phi-phi element of stress tensor
2126 *
2127 * @code
2128 *   local_quadrature_points_history[q].old_phiphi_stress =
2129 *   (2 * local_eta_ve * velocities[q](0)
2130 *   / fe_values.quadrature_point(q)[0]
2131 *   + local_chi_ve
2132 *   * local_quadrature_points_history[q].old_phiphi_stress);
2133 *   }
2134 *   }
2135 *   }
2136 *  
2137 * @endcode
2138 *
2139 * ====================== REDEFINE THE TIME INTERVAL FOR THE VISCOUS STEPS ======================
2140 *
2141 * @code
2142 *   template<int dim>
2143 *   void StokesProblem<dim>::update_time_interval()
2144 *   {
2145 *   double move_goal_per_step = system_parameters::initial_disp_target;
2146 *   if (system_parameters::present_timestep > system_parameters::initial_elastic_iterations)
2147 *   {
2148 *   move_goal_per_step = system_parameters::initial_disp_target -
2149 *   ((system_parameters::initial_disp_target - system_parameters::final_disp_target) /
2150 *   system_parameters::total_viscous_steps *
2151 *   (system_parameters::present_timestep - system_parameters::initial_elastic_iterations));
2152 *   }
2153 *  
2154 *   double zero_tolerance = 1e-3;
2155 *   double max_velocity = 0;
2156 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
2157 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)// loop over all cells
2158 *   {
2159 *   if (cell->at_boundary())
2160 *   {
2161 *   int zero_faces = 0;
2162 *   for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
2163 *   for (unsigned int i=0; i<dim; ++i)
2164 *   if (fabs(cell->face(f)->center()[i]) < zero_tolerance)
2165 *   zero_faces++;
2166 *   if (zero_faces==0)
2167 *   {
2168 *   for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
2169 *   {
2170 *   Point<dim> vertex_velocity;
2171 *   Point<dim> vertex_position;
2172 *   for (unsigned int d = 0; d < dim; ++d)
2173 *   {
2174 *   vertex_velocity[d] = solution(cell->vertex_dof_index(v, d));
2175 *   vertex_position[d] = cell->vertex(v)[d];
2176 *   }
2177 * @endcode
2178 *
2179 * velocity to be evaluated is the radial component of a surface vertex
2180 *
2181 * @code
2182 *   double local_velocity = 0;
2183 *   for (unsigned int d = 0; d < dim; ++d)
2184 *   {
2185 *   local_velocity += vertex_velocity[d] * vertex_position [d];
2186 *   }
2187 *   local_velocity /= std::sqrt( vertex_position.square() );
2188 *   if (local_velocity < 0)
2189 *   local_velocity *= -1;
2190 *   if (local_velocity > max_velocity)
2191 *   {
2192 *   max_velocity = local_velocity;
2193 *   }
2194 *   }
2195 *   }
2196 *   }
2197 *   }
2198 * @endcode
2199 *
2200 * NOTE: It is possible for this time interval to be very different from that used in the viscoelasticity calculation.
2201 *
2202 * @code
2203 *   system_parameters::current_time_interval = move_goal_per_step / max_velocity;
2204 *   double step_time_yr = system_parameters::current_time_interval / SECSINYEAR;
2205 *   std::cout << "Timestep interval changed to: "
2206 *   << step_time_yr
2207 *   << " years\n";
2208 *   }
2209 *  
2210 * @endcode
2211 *
2212 * ====================== MOVE MESH ======================
2213 *
2214
2215 *
2216 *
2217 * @code
2218 *   template<int dim>
2219 *   void StokesProblem<dim>::move_mesh()
2220 *   {
2221 *  
2222 *   std::cout << "\n" << " Moving mesh...\n";
2223 *   std::vector<bool> vertex_touched(triangulation.n_vertices(), false);
2224 *   for (typename DoFHandler<dim>::active_cell_iterator cell =
2225 *   dof_handler.begin_active(); cell != dof_handler.end(); ++cell)
2226 *   for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
2227 *   if (vertex_touched[cell->vertex_index(v)] == false)
2228 *   {
2229 *   vertex_touched[cell->vertex_index(v)] = true;
2230 *  
2231 *   Point<dim> vertex_displacement;
2232 *   for (unsigned int d = 0; d < dim; ++d)
2233 *   vertex_displacement[d] = solution(
2234 *   cell->vertex_dof_index(v, d));
2235 *   cell->vertex(v) += vertex_displacement
2236 *   * system_parameters::current_time_interval;
2237 *   }
2238 *   }
2239 *  
2240 * @endcode
2241 *
2242 * ====================== WRITE MESH TO FILE ======================
2243 *
2244
2245 *
2246 *
2247 * @code
2248 *   template<int dim>
2249 *   void StokesProblem<dim>::write_mesh()
2250 *   {
2251 * @endcode
2252 *
2253 * output mesh in ucd
2254 *
2255 * @code
2256 *   std::ostringstream initial_mesh_file;
2257 *   initial_mesh_file << system_parameters::output_folder << "/time" <<
2258 *   Utilities::int_to_string(system_parameters::present_timestep, 2) <<
2259 *   "_mesh.inp";
2260 *   std::ofstream out_ucd (initial_mesh_file.str().c_str());
2261 *   GridOut grid_out;
2262 *   grid_out.write_ucd (triangulation, out_ucd);
2263 *   }
2264 *  
2265 * @endcode
2266 *
2267 * ====================== FIT ELLIPSE TO SURFACE AND WRITE RADII TO FILE ======================
2268 *
2269
2270 *
2271 *
2272 * @code
2273 *   template<int dim>
2274 *   void StokesProblem<dim>::do_ellipse_fits()
2275 *   {
2276 *   std::ostringstream ellipses_filename;
2277 *   ellipses_filename << system_parameters::output_folder << "/ellipse_fits.txt";
2278 * @endcode
2279 *
2280 * Find ellipsoidal axes for all layers
2281 *
2282 * @code
2283 *   std::vector<double> ellipse_axes(0);
2284 * @endcode
2285 *
2286 * compute fit to boundary 0, 1, 2 ...
2287 *
2288 * @code
2289 *   std::cout << endl;
2290 *   for (unsigned int i = 0; i<system_parameters::sizeof_material_id; ++i)
2291 *   {
2292 *   ellipse_axes.clear();
2293 *   ellipsoid.compute_fit(ellipse_axes, system_parameters::material_id[i]);
2294 *   system_parameters::q_axes.push_back(ellipse_axes[0]);
2295 *   system_parameters::p_axes.push_back(ellipse_axes[1]);
2296 *  
2297 *   std::cout << "a_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[0]
2298 *   << " " << " c_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[1] << std::endl;
2299 *  
2300 *   std::ofstream fout_ellipses(ellipses_filename.str().c_str(), std::ios::app);
2301 *   fout_ellipses << system_parameters::present_timestep << " a_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[0]
2302 *   << " " << " c_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[1] << endl;
2303 *   fout_ellipses.close();
2304 *   }
2305 *   }
2306 *  
2307 * @endcode
2308 *
2309 * ====================== APPEND LINE TO PHYSICAL_TIMES.TXT FILE WITH STEP NUMBER, PHYSICAL TIME, AND # PLASTIC ITERATIONS ======================
2310 *
2311
2312 *
2313 *
2314 * @code
2315 *   template<int dim>
2316 *   void StokesProblem<dim>::append_physical_times(int max_plastic)
2317 *   {
2318 *   std::ostringstream times_filename;
2319 *   times_filename << system_parameters::output_folder << "/physical_times.txt";
2320 *   std::ofstream fout_times(times_filename.str().c_str(), std::ios::app);
2321 *   fout_times << system_parameters::present_timestep << " "
2322 *   << system_parameters::present_time/SECSINYEAR << " "
2323 *   << max_plastic << "\n";
2324 * @endcode
2325 *
2326 * << system_parameters::q_axes[0] << " " << system_parameters::p_axes[0] << " "
2327 * << system_parameters::q_axes[1] << " " << system_parameters::p_axes[1] << "\n";
2328 *
2329 * @code
2330 *   fout_times.close();
2331 *   }
2332 *  
2333 * @endcode
2334 *
2335 * ====================== WRITE VERTICES TO FILE ======================
2336 *
2337
2338 *
2339 *
2340 * @code
2341 *   template<int dim>
2342 *   void StokesProblem<dim>::write_vertices(unsigned char boundary_that_we_need)
2343 *   {
2344 *   std::ostringstream vertices_output;
2345 *   vertices_output << system_parameters::output_folder << "/time" <<
2346 *   Utilities::int_to_string(system_parameters::present_timestep, 2) << "_" <<
2347 *   Utilities::int_to_string(boundary_that_we_need, 2) <<
2348 *   "_surface.txt";
2349 *   std::ofstream fout_final_vertices(vertices_output.str().c_str());
2350 *   fout_final_vertices.close();
2351 *  
2352 *   std::vector<bool> vertex_touched(triangulation.n_vertices(), false);
2353 *  
2354 *   if (boundary_that_we_need == 0)
2355 *   {
2356 * @endcode
2357 *
2358 * Figure out if the vertex is on the boundary of the domain
2359 *
2360 * @code
2361 *   for (typename Triangulation<dim>::active_cell_iterator cell =
2362 *   triangulation.begin_active(); cell != triangulation.end(); ++cell)
2363 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
2364 *   {
2365 *   unsigned char boundary_ids = cell->face(f)->boundary_id();
2366 *   if (boundary_ids == boundary_that_we_need)
2367 *   {
2368 *   for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
2369 *   if (vertex_touched[cell->face(f)->vertex_index(v)] == false)
2370 *   {
2371 *   vertex_touched[cell->face(f)->vertex_index(v)] = true;
2372 *   std::ofstream fout_final_vertices(vertices_output.str().c_str(), std::ios::app);
2373 *   fout_final_vertices << cell->face(f)->vertex(v) << "\n";
2374 *   fout_final_vertices.close();
2375 *   }
2376 *   }
2377 *   }
2378 *   }
2379 *   else
2380 *   {
2381 * @endcode
2382 *
2383 * Figure out if the vertex is on an internal boundary
2384 *
2385 * @code
2386 *   for (typename Triangulation<dim>::active_cell_iterator cell =
2387 *   triangulation.begin_active(); cell != triangulation.end(); ++cell)
2388 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
2389 *   {
2390 *   if (cell->neighbor(f) != triangulation.end())
2391 *   {
2392 *   if (cell->material_id() != cell->neighbor(f)->material_id()) //finds face is at internal boundary
2393 *   {
2394 *   int high_mat_id = std::max(cell->material_id(),
2395 *   cell->neighbor(f)->material_id());
2396 *   if (high_mat_id == boundary_that_we_need) //finds faces at the correct internal boundary
2397 *   {
2398 *   for (unsigned int v = 0;
2399 *   v < GeometryInfo<dim>::vertices_per_face;
2400 *   ++v)
2401 *   if (vertex_touched[cell->face(f)->vertex_index(
2402 *   v)] == false)
2403 *   {
2404 *   vertex_touched[cell->face(f)->vertex_index(
2405 *   v)] = true;
2406 *   std::ofstream fout_final_vertices(vertices_output.str().c_str(), std::ios::app);
2407 *   fout_final_vertices << cell->face(f)->vertex(v) << "\n";
2408 *   fout_final_vertices.close();
2409 *   }
2410 *   }
2411 *   }
2412 *   }
2413 *   }
2414 *   }
2415 *   }
2416 *  
2417 * @endcode
2418 *
2419 * ====================== SETUP INITIAL MESH ======================
2420 *
2421
2422 *
2423 *
2424 * @code
2425 *   template<int dim>
2426 *   void StokesProblem<dim>::setup_initial_mesh()
2427 *   {
2428 *   GridIn<dim> grid_in;
2429 *   grid_in.attach_triangulation(triangulation);
2430 *   std::ifstream mesh_stream(system_parameters::mesh_filename,
2431 *   std::ifstream::in);
2432 *   grid_in.read_ucd(mesh_stream);
2433 *  
2434 * @endcode
2435 *
2436 * output initial mesh in eps
2437 *
2438 * @code
2439 *   std::ostringstream initial_mesh_file;
2440 *   initial_mesh_file << system_parameters::output_folder << "/initial_mesh.eps";
2441 *   std::ofstream out_eps (initial_mesh_file.str().c_str());
2442 *   GridOut grid_out;
2443 *   grid_out.write_eps (triangulation, out_eps);
2444 *   out_eps.close();
2445 *  
2446 * @endcode
2447 *
2448 * set boundary ids
2449 * boundary indicator 0 is outer free surface; 1, 2, 3 ... is boundary between layers, 99 is flat boundaries
2450 *
2451 * @code
2452 *   typename Triangulation<dim>::active_cell_iterator
2453 *   cell=triangulation.begin_active(), endc=triangulation.end();
2454 *  
2455 *   unsigned int how_many; // how many components away from cardinal planes
2456 *  
2457 *   std::ostringstream boundaries_file;
2458 *   boundaries_file << system_parameters::output_folder << "/boundaries.txt";
2459 *   std::ofstream fout_boundaries(boundaries_file.str().c_str());
2460 *   fout_boundaries.close();
2461 *  
2462 *   double zero_tolerance = 1e-3;
2463 *   for (; cell != endc; ++cell) // loop over all cells
2464 *   {
2465 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f) // loop over all vertices
2466 *   {
2467 *   if (cell->face(f)->at_boundary())
2468 *   {
2469 * @endcode
2470 *
2471 * print boundary
2472 *
2473 * @code
2474 *   std::ofstream fout_boundaries(boundaries_file.str().c_str(), std::ios::app);
2475 *   fout_boundaries << cell->face(f)->center()[0] << " " << cell->face(f)->center()[1]<< "\n";
2476 *   fout_boundaries.close();
2477 *  
2478 *   how_many = 0;
2479 *   for (unsigned int i=0; i<dim; ++i)
2480 *   if (fabs(cell->face(f)->center()[i]) > zero_tolerance)
2481 *   how_many++;
2482 *   if (how_many==dim)
2483 *   cell->face(f)->set_all_boundary_ids(0); // if face center coordinates > zero_tol, set bnry indicators to 0
2484 *   else
2485 *   cell->face(f)->set_all_boundary_ids(99);
2486 *   }
2487 *   }
2488 *   }
2489 *  
2490 *   std::ostringstream ellipses_filename;
2491 *   ellipses_filename << system_parameters::output_folder << "/ellipse_fits.txt";
2492 *   std::ofstream fout_ellipses(ellipses_filename.str().c_str());
2493 *   fout_ellipses.close();
2494 *  
2495 * @endcode
2496 *
2497 * Find ellipsoidal axes for all layers
2498 *
2499 * @code
2500 *   std::vector<double> ellipse_axes(0);
2501 * @endcode
2502 *
2503 * compute fit to boundary 0, 1, 2 ...
2504 *
2505 * @code
2506 *   std::cout << endl;
2507 *   for (unsigned int i = 0; i<system_parameters::sizeof_material_id; ++i)
2508 *   {
2509 *   ellipse_axes.clear();
2510 *   ellipsoid.compute_fit(ellipse_axes, system_parameters::material_id[i]);
2511 *   system_parameters::q_axes.push_back(ellipse_axes[0]);
2512 *   system_parameters::p_axes.push_back(ellipse_axes[1]);
2513 *  
2514 *   std::cout << "a_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[0]
2515 *   << " " << " c_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[1] << std::endl;
2516 *  
2517 *   std::ofstream fout_ellipses(ellipses_filename.str().c_str(), std::ios::app);
2518 *   fout_ellipses << system_parameters::present_timestep << " a_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[0]
2519 *   << " " << " c_"<< system_parameters::material_id[i] <<" = " << ellipse_axes[1] << endl;
2520 *   fout_ellipses.close();
2521 *   }
2522 *  
2523 *   triangulation.refine_global(system_parameters::global_refinement);
2524 *  
2525 *  
2526 * @endcode
2527 *
2528 * refines crustal region
2529 *
2530 * @code
2531 *   if (system_parameters::crustal_refinement != 0)
2532 *   {
2533 *   double a = system_parameters::q_axes[0] - system_parameters::crust_refine_region;
2534 *   double b = system_parameters::p_axes[0] - system_parameters::crust_refine_region;
2535 *  
2536 *  
2537 *   for (unsigned int step = 0;
2538 *   step < system_parameters::crustal_refinement; ++step)
2539 *   {
2540 *   typename ::Triangulation<dim>::active_cell_iterator cell =
2541 *   triangulation.begin_active(), endc = triangulation.end();
2542 *   for (; cell != endc; ++cell)
2543 *   for (unsigned int v = 0;
2544 *   v < GeometryInfo<dim>::vertices_per_cell; ++v)
2545 *   {
2546 *   Point<dim> current_vertex = cell->vertex(v);
2547 *  
2548 *   const double x_coord = current_vertex.operator()(0);
2549 *   const double y_coord = current_vertex.operator()(1);
2550 *   double expected_z = -1;
2551 *  
2552 *   if ((x_coord - a) < -1e-10)
2553 *   expected_z = b
2554 *   * std::sqrt(1 - (x_coord * x_coord / a / a));
2555 *  
2556 *   if (y_coord >= expected_z)
2557 *   {
2558 *   cell->set_refine_flag();
2559 *   break;
2560 *   }
2561 *   }
2562 *   triangulation.execute_coarsening_and_refinement();
2563 *   }
2564 *   }
2565 *  
2566 *  
2567 * @endcode
2568 *
2569 * output initial mesh in eps
2570 *
2571 * @code
2572 *   std::ostringstream refined_mesh_file;
2573 *   refined_mesh_file << system_parameters::output_folder << "/refined_mesh.eps";
2574 *   std::ofstream out_eps_refined (refined_mesh_file.str().c_str());
2575 *   GridOut grid_out_refined;
2576 *   grid_out_refined.write_eps (triangulation, out_eps_refined);
2577 *   out_eps_refined.close();
2578 *   write_vertices(0);
2579 *   write_vertices(1);
2580 *   write_mesh();
2581 *   }
2582 *  
2583 * @endcode
2584 *
2585 * ====================== REFINE MESH ======================
2586 *
2587
2588 *
2589 *
2590 * @code
2591 *   template<int dim>
2592 *   void StokesProblem<dim>::refine_mesh()
2593 *   {
2594 *   using FunctionMap = std::map<types::boundary_id, const Function<dim> *>;
2595 *   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
2596 *  
2597 *   std::vector<bool> component_mask(dim + 1, false);
2598 *   component_mask[dim] = true;
2599 *   KellyErrorEstimator<dim>::estimate(dof_handler, QGauss<dim - 1>(degree + 1),
2600 *   FunctionMap(), solution,
2601 *   estimated_error_per_cell, component_mask);
2602 *  
2603 *   GridRefinement::refine_and_coarsen_fixed_number(triangulation,
2604 *   estimated_error_per_cell, 0.3, 0.0);
2605 *   triangulation.execute_coarsening_and_refinement();
2606 *   }
2607 *  
2608 * @endcode
2609 *
2610 * ====================== SET UP THE DATA STRUCTURES TO REMEMBER STRESS FIELD ======================
2611 *
2612 * @code
2613 *   template<int dim>
2614 *   void StokesProblem<dim>::setup_quadrature_point_history()
2615 *   {
2616 *   unsigned int our_cells = 0;
2617 *   for (typename Triangulation<dim>::active_cell_iterator cell =
2618 *   triangulation.begin_active(); cell != triangulation.end(); ++cell)
2619 *   ++our_cells;
2620 *  
2621 *   triangulation.clear_user_data();
2622 *  
2623 *   quadrature_point_history.resize(our_cells * quadrature_formula.size());
2624 *  
2625 *   unsigned int history_index = 0;
2626 *   for (typename Triangulation<dim>::active_cell_iterator cell =
2627 *   triangulation.begin_active(); cell != triangulation.end(); ++cell)
2628 *   {
2629 *   cell->set_user_pointer(&quadrature_point_history[history_index]);
2630 *   history_index += quadrature_formula.size();
2631 *   }
2632 *  
2633 *   Assert(history_index == quadrature_point_history.size(), ExcInternalError());
2634 *   }
2635 *  
2636 * @endcode
2637 *
2638 * ====================== DOES ELASTIC STEPS ======================
2639 *
2640 * @code
2641 *   template<int dim>
2642 *   void StokesProblem<dim>::do_elastic_steps()
2643 *   {
2644 *   unsigned int elastic_iteration = 0;
2645 *  
2646 *   while (elastic_iteration < system_parameters::initial_elastic_iterations)
2647 *   {
2648 *  
2649 *   std::cout << "\n\nElastic iteration " << elastic_iteration
2650 *   << "\n";
2651 *   setup_dofs();
2652 *  
2653 *   if (system_parameters::present_timestep == 0)
2654 *   initialize_eta_and_G();
2655 *  
2656 *   if (elastic_iteration == 0)
2657 *   system_parameters::current_time_interval =
2658 *   system_parameters::viscous_time; //This is the time interval needed in assembling the problem
2659 *  
2660 *   std::cout << " Assembling..." << std::endl << std::flush;
2661 *   assemble_system();
2662 *  
2663 *   std::cout << " Solving..." << std::flush;
2664 *   solve();
2665 *  
2666 *   output_results();
2667 *   update_quadrature_point_history();
2668 *  
2669 *   append_physical_times(0);
2670 *   elastic_iteration++;
2671 *   system_parameters::present_timestep++;
2672 *   do_ellipse_fits();
2673 *   write_vertices(0);
2674 *   write_vertices(1);
2675 *   write_mesh();
2676 *   update_time_interval();
2677 *   }
2678 *   }
2679 *  
2680 * @endcode
2681 *
2682 * ====================== DO A SINGLE VISCOELASTOPLASTIC TIMESTEP ======================
2683 *
2684 * @code
2685 *   template<int dim>
2686 *   void StokesProblem<dim>::do_flow_step()
2687 *   {
2688 *   plastic_iteration = 0;
2689 *   while (plastic_iteration < system_parameters::max_plastic_iterations)
2690 *   {
2691 *   if (system_parameters::continue_plastic_iterations == true)
2692 *   {
2693 *   std::cout << "Plasticity iteration " << plastic_iteration << "\n";
2694 *   setup_dofs();
2695 *  
2696 *   std::cout << " Assembling..." << std::endl << std::flush;
2697 *   assemble_system();
2698 *  
2699 *   std::cout << " Solving..." << std::flush;
2700 *   solve();
2701 *  
2702 *   output_results();
2703 *   solution_stesses();
2704 *  
2705 *   if (system_parameters::continue_plastic_iterations == false)
2706 *   break;
2707 *  
2708 *   plastic_iteration++;
2709 *   }
2710 *   }
2711 *   }
2712 *  
2713 * @endcode
2714 *
2715 * ====================== RUN ======================
2716 *
2717
2718 *
2719 *
2720 * @code
2721 *   template<int dim>
2722 *   void StokesProblem<dim>::run()
2723 *   {
2724 * @endcode
2725 *
2726 * Sets up mesh and data structure for viscosity and stress at quadrature points
2727 *
2728 * @code
2729 *   setup_initial_mesh();
2730 *   setup_quadrature_point_history();
2731 *  
2732 * @endcode
2733 *
2734 * Makes the physical_times.txt file
2735 *
2736 * @code
2737 *   std::ostringstream times_filename;
2738 *   times_filename << system_parameters::output_folder << "/physical_times.txt";
2739 *   std::ofstream fout_times(times_filename.str().c_str());
2740 *   fout_times.close();
2741 *  
2742 * @endcode
2743 *
2744 * Computes elastic timesteps
2745 *
2746 * @code
2747 *   do_elastic_steps();
2748 * @endcode
2749 *
2750 * Computes viscous timesteps
2751 *
2752 * @code
2753 *   unsigned int VEPstep = 0;
2754 *   while (system_parameters::present_timestep
2755 *   < (system_parameters::initial_elastic_iterations
2756 *   + system_parameters::total_viscous_steps))
2757 *   {
2758 *  
2759 *   if (system_parameters::continue_plastic_iterations == false)
2760 *   system_parameters::continue_plastic_iterations = true;
2761 *   std::cout << "\n\nViscoelastoplastic iteration " << VEPstep << "\n\n";
2762 * @endcode
2763 *
2764 * Computes plasticity
2765 *
2766 * @code
2767 *   do_flow_step();
2768 *   update_quadrature_point_history();
2769 *   move_mesh();
2770 *   append_physical_times(plastic_iteration);
2771 *   system_parameters::present_timestep++;
2772 *   system_parameters::present_time = system_parameters::present_time + system_parameters::current_time_interval;
2773 *   do_ellipse_fits();
2774 *   write_vertices(0);
2775 *   write_vertices(1);
2776 *   write_mesh();
2777 *   VEPstep++;
2778 *   }
2779 *   append_physical_times(-1);
2780 *   }
2781 *  
2782 *   }
2783 *  
2784 * @endcode
2785 *
2786 * ====================== MAIN ======================
2787 *
2788
2789 *
2790 *
2791 * @code
2792 *   int main(int argc, char *argv[])
2793 *   {
2794 *  
2795 * @endcode
2796 *
2797 * output program name
2798 *
2799 * @code
2800 *   std::cout << "Running: " << argv[0] << std::endl;
2801 *  
2802 *   char *cfg_filename = new char[120];
2803 *  
2804 *   if (argc == 1) // if no input parameters (as if launched from eclipse)
2805 *   {
2806 *   std::strcpy(cfg_filename,"config/sample_CeresFE_config.cfg");
2807 *   }
2808 *   else
2809 *   std::strcpy(cfg_filename,argv[1]);
2810 *  
2811 *   try
2812 *   {
2813 *   using namespace dealii;
2814 *   using namespace Step22;
2815 *   config_in cfg(cfg_filename);
2816 *  
2817 *   std::clock_t t1;
2818 *   std::clock_t t2;
2819 *   t1 = std::clock();
2820 *  
2821 *   deallog.depth_console(0);
2822 *  
2823 *   StokesProblem<2> flow_problem(1);
2824 *   flow_problem.run();
2825 *  
2826 *   std::cout << std::endl << "\a";
2827 *  
2828 *   t2 = std::clock();
2829 *   float diff (((float)t2 - (float)t1) / (float)CLOCKS_PER_SEC);
2830 *   std::cout << "\n Program run in: " << diff << " seconds" << endl;
2831 *   }
2832 *   catch (std::exception &exc)
2833 *   {
2834 *   std::cerr << std::endl << std::endl
2835 *   << "----------------------------------------------------"
2836 *   << std::endl;
2837 *   std::cerr << "Exception on processing: " << std::endl << exc.what()
2838 *   << std::endl << "Aborting!" << std::endl
2839 *   << "----------------------------------------------------"
2840 *   << std::endl;
2841 *  
2842 *   return 1;
2843 *   }
2844 *   catch (...)
2845 *   {
2846 *   std::cerr << std::endl << std::endl
2847 *   << "----------------------------------------------------"
2848 *   << std::endl;
2849 *   std::cerr << "Unknown exception!" << std::endl << "Aborting!"
2850 *   << std::endl
2851 *   << "----------------------------------------------------"
2852 *   << std::endl;
2853 *   return 1;
2854 *   }
2855 *  
2856 *   return 0;
2857 *   }
2858 * @endcode
2859
2860
2861<a name="ann-support_code/config_in.h"></a>
2862<h1>Annotated version of support_code/config_in.h</h1>
2863 *
2864 *
2865 *
2866 *
2867 * @code
2868 *   /* -----------------------------------------------------------------------------
2869 *   *
2870 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
2871 *   * Copyright (C) 2015 by Anton Ermakov
2872 *   *
2873 *   * This file is part of the deal.II code gallery.
2874 *   *
2875 *   * -----------------------------------------------------------------------------
2876 *   *
2877 *   * Author: antonermakov
2878 *   */
2879 *  
2880 *  
2881 *   #include <iostream>
2882 *   #include <iomanip>
2883 *   #include <cstdlib>
2884 *   #include <string>
2885 *   #include <libconfig.h++>
2886 *  
2887 *   #include "local_math.h"
2888 *  
2889 *   using namespace std;
2890 *   using namespace libconfig;
2891 *  
2892 *   namespace Step22
2893 *   {
2894 *   namespace system_parameters
2895 *   {
2896 *  
2897 * @endcode
2898 *
2899 * Mesh file name
2900 *
2901 * @code
2902 *   string mesh_filename;
2903 *   string output_folder;
2904 *  
2905 * @endcode
2906 *
2907 * Body parameters
2908 *
2909 * @code
2910 *   double r_mean;
2911 *   double period;
2912 *   double omegasquared;
2913 *   double beta;
2914 *   double intercept;
2915 *  
2916 * @endcode
2917 *
2918 * Rheology parameters
2919 *
2920 * @code
2921 *   vector<double> depths_eta;
2922 *   vector<double> eta_kinks;
2923 *   vector<double> depths_rho;
2924 *   vector<double> rho;
2925 *   vector<int> material_id;
2926 *   vector<double> G;
2927 *  
2928 *   double eta_ceiling;
2929 *   double eta_floor;
2930 *   double eta_Ea;
2931 *   bool lat_dependence;
2932 *  
2933 *   unsigned int sizeof_depths_eta;
2934 *   unsigned int sizeof_depths_rho;
2935 *   unsigned int sizeof_rho;
2936 *   unsigned int sizeof_eta_kinks;
2937 *   unsigned int sizeof_material_id;
2938 *   unsigned int sizeof_G;
2939 *  
2940 *   double pressure_scale;
2941 *   double q;
2942 *   bool cylindrical;
2943 *   bool continue_plastic_iterations;
2944 *  
2945 * @endcode
2946 *
2947 * plasticity variables
2948 *
2949 * @code
2950 *   bool plasticity_on;
2951 *   unsigned int failure_criterion;
2952 *   unsigned int max_plastic_iterations;
2953 *   double smoothing_radius;
2954 *  
2955 * @endcode
2956 *
2957 * viscoelasticity variables
2958 *
2959 * @code
2960 *   unsigned int initial_elastic_iterations;
2961 *   double elastic_time;
2962 *   double viscous_time;
2963 *   double initial_disp_target;
2964 *   double final_disp_target;
2965 *   double current_time_interval;
2966 *  
2967 * @endcode
2968 *
2969 * mesh refinement variables
2970 *
2971 * @code
2972 *   unsigned int global_refinement;
2973 *   unsigned int small_r_refinement;
2974 *   unsigned int crustal_refinement;
2975 *   double crust_refine_region;
2976 *   unsigned int surface_refinement;
2977 *  
2978 * @endcode
2979 *
2980 * solver variables
2981 *
2982 * @code
2983 *   int iteration_coefficient;
2984 *   double tolerance_coefficient;
2985 *  
2986 * @endcode
2987 *
2988 * time step variables
2989 *
2990 * @code
2991 *   double present_time;
2992 *   unsigned int present_timestep;
2993 *   unsigned int total_viscous_steps;
2994 *  
2995 *  
2996 * @endcode
2997 *
2998 * ellipse axes
2999 *
3000 * @code
3001 *   vector<double> q_axes;
3002 *   vector<double> p_axes;
3003 *  
3004 *   }
3005 *  
3006 *   class config_in
3007 *   {
3008 *   public:
3009 *   config_in(char *);
3010 *  
3011 *   private:
3012 *   void write_config();
3013 *   };
3014 *  
3015 *   void config_in::write_config()
3016 *   {
3017 *   std::ostringstream config_parameters;
3018 *   config_parameters << system_parameters::output_folder << "/run_parameters.txt";
3019 *   std::ofstream fout_config(config_parameters.str().c_str());
3020 *  
3021 * @endcode
3022 *
3023 * mesh filename
3024 *
3025 * @code
3026 *   fout_config << "mesh filename: " << system_parameters::mesh_filename << endl << endl;
3027 *  
3028 * @endcode
3029 *
3030 * body parameters
3031 *
3032 * @code
3033 *   fout_config << "r_mean = " << system_parameters::r_mean << endl;
3034 *   fout_config << "period = " << system_parameters::period << endl;
3035 *   fout_config << "omegasquared = " << system_parameters::omegasquared << endl;
3036 *   fout_config << "beta = " << system_parameters::beta << endl;
3037 *   fout_config << "intercept = " << system_parameters::intercept << endl;
3038 *  
3039 * @endcode
3040 *
3041 * rheology parameters
3042 *
3043
3044 *
3045 *
3046 * @code
3047 *   for (unsigned int i=0; i<system_parameters::sizeof_depths_eta; ++i)
3048 *   fout_config << "depths_eta[" << i << "] = " << system_parameters::depths_eta[i] << endl;
3049 *  
3050 *   for (unsigned int i=0; i<system_parameters::sizeof_eta_kinks; ++i)
3051 *   fout_config << "eta_kinks[" << i << "] = " << system_parameters::eta_kinks[i] << endl;
3052 *  
3053 *   for (unsigned int i=0; i<system_parameters::sizeof_depths_rho; ++i)
3054 *   fout_config << "depths_rho[" << i << "] = " << system_parameters::depths_rho[i] << endl;
3055 *  
3056 *   for (unsigned int i=0; i<system_parameters::sizeof_rho; ++i)
3057 *   fout_config << "rho[" << i << "] = " << system_parameters::rho[i] << endl;
3058 *  
3059 *   for (unsigned int i=0; i<system_parameters::sizeof_material_id; ++i)
3060 *   fout_config << "material_id[" << i << "] = " << system_parameters::material_id[i] << endl;
3061 *  
3062 *   for (unsigned int i=0; i<system_parameters::sizeof_G; ++i)
3063 *   fout_config << "G[" << i << "] = " << system_parameters::G[i] << endl;
3064 *  
3065 *   fout_config << "eta_ceiling = " << system_parameters::eta_ceiling << endl;
3066 *   fout_config << "eta_floor = " << system_parameters::eta_floor << endl;
3067 *   fout_config << "eta_Ea = " << system_parameters::eta_Ea << endl;
3068 *   fout_config << "lat_dependence = " << system_parameters::lat_dependence << endl;
3069 *   fout_config << "pressure_scale = " << system_parameters::pressure_scale << endl;
3070 *   fout_config << "q = " << system_parameters::q << endl;
3071 *   fout_config << "cylindrical = " << system_parameters::cylindrical << endl;
3072 *   fout_config << "continue_plastic_iterations = " << system_parameters::continue_plastic_iterations << endl;
3073 *  
3074 * @endcode
3075 *
3076 * Plasticity parameters
3077 *
3078 * @code
3079 *   fout_config << "plasticity_on = " << system_parameters::plasticity_on << endl;
3080 *   fout_config << "failure_criterion = " << system_parameters::failure_criterion << endl;
3081 *   fout_config << "max_plastic_iterations = " << system_parameters::max_plastic_iterations << endl;
3082 *   fout_config << "smoothing_radius = " << system_parameters::smoothing_radius << endl;
3083 *  
3084 * @endcode
3085 *
3086 * Viscoelasticity parameters
3087 *
3088 * @code
3089 *   fout_config << "initial_elastic_iterations = " << system_parameters::initial_elastic_iterations << endl;
3090 *   fout_config << "elastic_time = " << system_parameters::elastic_time << endl;
3091 *   fout_config << "viscous_time = " << system_parameters::viscous_time << endl;
3092 *   fout_config << "initial_disp_target = " << system_parameters::initial_disp_target << endl;
3093 *   fout_config << "final_disp_target = " << system_parameters::final_disp_target << endl;
3094 *   fout_config << "current_time_interval = " << system_parameters::current_time_interval << endl;
3095 *  
3096 * @endcode
3097 *
3098 * Mesh refinement parameters
3099 *
3100 * @code
3101 *   fout_config << "global_refinement = " << system_parameters::global_refinement << endl;
3102 *   fout_config << "small_r_refinement = " << system_parameters::small_r_refinement << endl;
3103 *   fout_config << "crustal_refinement = " << system_parameters::crustal_refinement << endl;
3104 *   fout_config << "crust_refine_region = " << system_parameters::crust_refine_region << endl;
3105 *   fout_config << "surface_refinement = " << system_parameters::surface_refinement << endl;
3106 *  
3107 * @endcode
3108 *
3109 * Solver parameters
3110 *
3111 * @code
3112 *   fout_config << "iteration_coefficient = " << system_parameters::iteration_coefficient << endl;
3113 *   fout_config << "tolerance_coefficient = " << system_parameters::tolerance_coefficient << endl;
3114 *  
3115 * @endcode
3116 *
3117 * Time step parameters
3118 *
3119 * @code
3120 *   fout_config << "present_time = " << system_parameters::present_time << endl;
3121 *   fout_config << "present_timestep = " << system_parameters::present_timestep << endl;
3122 *   fout_config << "total_viscous_steps = " << system_parameters::total_viscous_steps << endl;
3123 *  
3124 *   fout_config.close();
3125 *   }
3126 *  
3127 *   config_in::config_in(char *filename)
3128 *   {
3129 *  
3130 * @endcode
3131 *
3132 * This example reads the configuration file 'example.cfg' and displays
3133 * some of its contents.
3134 *
3135
3136 *
3137 *
3138 * @code
3139 *   Config cfg;
3140 *  
3141 * @endcode
3142 *
3143 * Read the file. If there is an error, report it and exit.
3144 *
3145 * @code
3146 *   try
3147 *   {
3148 *   cfg.readFile(filename);
3149 *   }
3150 *   catch (const FileIOException &fioex)
3151 *   {
3152 *   std::cerr << "I/O error while reading file:" << filename << std::endl;
3153 *   }
3154 *   catch (const ParseException &pex)
3155 *   {
3156 *   std::cerr << "Parse error at " << pex.getFile() << ":" << pex.getLine()
3157 *   << " - " << pex.getError() << std::endl;
3158 *   }
3159 *  
3160 * @endcode
3161 *
3162 * Get mesh name.
3163 *
3164 * @code
3165 *   try
3166 *   {
3167 *   string msh = cfg.lookup("mesh_filename");
3168 *   system_parameters::mesh_filename = msh;
3169 *   }
3170 *   catch (const SettingNotFoundException &nfex)
3171 *   {
3172 *   cerr << "No 'mesh_filename' setting in configuration file." << endl;
3173 *   }
3174 *  
3175 * @endcode
3176 *
3177 * get output folder
3178 *
3179
3180 *
3181 *
3182 * @code
3183 *   try
3184 *   {
3185 *   string output = cfg.lookup("output_folder");
3186 *   system_parameters::output_folder = output;
3187 *  
3188 *   std::cout << "Writing to folder: " << output << endl;
3189 *   }
3190 *   catch (const SettingNotFoundException &nfex)
3191 *   {
3192 *   cerr << "No 'output_folder' setting in configuration file." << endl;
3193 *   }
3194 *  
3195 * @endcode
3196 *
3197 * get radii
3198 *
3199
3200 *
3201 *
3202 * @code
3203 *   const Setting &root = cfg.getRoot();
3204 *  
3205 *  
3206 * @endcode
3207 *
3208 * get body parameters
3209 *
3210 * @code
3211 *   try
3212 *   {
3213 *   const Setting &body_parameters = root["body_parameters"];
3214 *  
3215 *   body_parameters.lookupValue("period", system_parameters::period);
3216 *   system_parameters::omegasquared = pow(TWOPI / 3600.0 / system_parameters::period, 2.0);
3217 *   body_parameters.lookupValue("r_mean", system_parameters::r_mean);
3218 *   body_parameters.lookupValue("beta", system_parameters::beta);
3219 *   body_parameters.lookupValue("intercept", system_parameters::intercept);
3220 *   }
3221 *   catch (const SettingNotFoundException &nfex)
3222 *   {
3223 *   cerr << "We've got a problem in the body parameters block" << endl;
3224 *  
3225 *   }
3226 *  
3227 * @endcode
3228 *
3229 * Rheology parameters
3230 *
3231
3232 *
3233 *
3234 * @code
3235 *   try
3236 *   {
3237 * @endcode
3238 *
3239 * get depths_eta ---------------------
3240 *
3241 * @code
3242 *   const Setting &set_depths_eta = cfg.lookup("rheology_parameters.depths_eta");
3243 *  
3244 *   unsigned int ndepths_eta = set_depths_eta.getLength();
3245 *   system_parameters::sizeof_depths_eta = ndepths_eta;
3246 *  
3247 *   for (unsigned int i=0; i<ndepths_eta; ++i)
3248 *   {
3249 *   system_parameters::depths_eta.push_back(set_depths_eta[i]);
3250 *   cout << "depth_eta[" << i << "] = " << system_parameters::depths_eta[i] << endl;
3251 *   }
3252 *  
3253 * @endcode
3254 *
3255 * get eta_kinks -------------------------
3256 *
3257 * @code
3258 *   const Setting &set_eta_kinks = cfg.lookup("rheology_parameters.eta_kinks");
3259 *  
3260 *   unsigned int neta_kinks = set_eta_kinks.getLength();
3261 *   system_parameters::sizeof_eta_kinks = neta_kinks;
3262 *  
3263 * @endcode
3264 *
3265 * cout << "Number of depth = " << ndepths << endl;
3266 *
3267
3268 *
3269 *
3270 * @code
3271 *   for (unsigned int i=0; i<neta_kinks; ++i)
3272 *   {
3273 *   system_parameters::eta_kinks.push_back(set_eta_kinks[i]);
3274 *   cout << "eta_kinks[" << i << "] = " << system_parameters::eta_kinks[i] << endl;
3275 *   }
3276 *  
3277 * @endcode
3278 *
3279 * get depths_rho -------------------------
3280 *
3281 * @code
3282 *   const Setting &set_depths_rho = cfg.lookup("rheology_parameters.depths_rho");
3283 *  
3284 *   unsigned int ndepths_rho = set_depths_rho.getLength();
3285 *   system_parameters::sizeof_depths_rho = ndepths_rho;
3286 *  
3287 * @endcode
3288 *
3289 * cout << "Number of depth = " << ndepths << endl;
3290 *
3291
3292 *
3293 *
3294 * @code
3295 *   for (unsigned int i=0; i<ndepths_rho; ++i)
3296 *   {
3297 *   system_parameters::depths_rho.push_back(set_depths_rho[i]);
3298 *   cout << "depths_rho[" << i << "] = " << system_parameters::depths_rho[i] << endl;
3299 *   }
3300 *  
3301 * @endcode
3302 *
3303 * get rho -------------------------
3304 *
3305 * @code
3306 *   const Setting &set_rho = cfg.lookup("rheology_parameters.rho");
3307 *  
3308 *   unsigned int nrho = set_rho.getLength();
3309 *   system_parameters::sizeof_rho = nrho;
3310 *  
3311 * @endcode
3312 *
3313 * cout << "Number of depth = " << ndepths << endl;
3314 *
3315
3316 *
3317 *
3318 * @code
3319 *   for (unsigned int i=0; i<nrho; ++i)
3320 *   {
3321 *   system_parameters::rho.push_back(set_rho[i]);
3322 *   cout << "rho[" << i << "] = " << system_parameters::rho[i] << endl;
3323 *   }
3324 *  
3325 * @endcode
3326 *
3327 * get material_id -------------------------
3328 *
3329 * @code
3330 *   const Setting &set_material_id = cfg.lookup("rheology_parameters.material_id");
3331 *  
3332 *   unsigned int nmaterial_id = set_material_id.getLength();
3333 *   system_parameters::sizeof_material_id = nmaterial_id;
3334 *  
3335 * @endcode
3336 *
3337 * cout << "Number of depth = " << ndepths << endl;
3338 *
3339
3340 *
3341 *
3342 * @code
3343 *   for (unsigned int i=0; i<nmaterial_id; ++i)
3344 *   {
3345 *   system_parameters::material_id.push_back(set_material_id[i]);
3346 *   cout << "material_id[" << i << "] = " << system_parameters::material_id[i] << endl;
3347 *   }
3348 *  
3349 * @endcode
3350 *
3351 * get G -------------------------
3352 *
3353 * @code
3354 *   const Setting &set_G = cfg.lookup("rheology_parameters.G");
3355 *  
3356 *   unsigned int nG = set_G.getLength();
3357 *   system_parameters::sizeof_G = nG;
3358 *  
3359 * @endcode
3360 *
3361 * cout << "Number of depth = " << ndepths << endl;
3362 *
3363
3364 *
3365 *
3366 * @code
3367 *   for (unsigned int i=0; i<nG; ++i)
3368 *   {
3369 *   system_parameters::G.push_back(set_G[i]);
3370 *   cout << "G[" << i << "] = " << system_parameters::G[i] << endl;
3371 *   }
3372 *  
3373 *   const Setting &rheology_parameters = root["rheology_parameters"];
3374 *   rheology_parameters.lookupValue("eta_ceiling", system_parameters::eta_ceiling);
3375 *   rheology_parameters.lookupValue("eta_floor", system_parameters::eta_floor);
3376 *   rheology_parameters.lookupValue("eta_Ea", system_parameters::eta_Ea);
3377 *   rheology_parameters.lookupValue("lat_dependence", system_parameters::lat_dependence);
3378 *   rheology_parameters.lookupValue("pressure_scale", system_parameters::pressure_scale);
3379 *   rheology_parameters.lookupValue("q", system_parameters::q);
3380 *   rheology_parameters.lookupValue("cylindrical", system_parameters::cylindrical);
3381 *   rheology_parameters.lookupValue("continue_plastic_iterations", system_parameters::continue_plastic_iterations);
3382 *   }
3383 *   catch (const SettingNotFoundException &nfex)
3384 *   {
3385 *   cerr << "We've got a problem in the rheology parameters block" << endl;
3386 *   }
3387 *  
3388 * @endcode
3389 *
3390 * Plasticity parameters
3391 *
3392 * @code
3393 *   try
3394 *   {
3395 *  
3396 *   const Setting &plasticity_parameters = root["plasticity_parameters"];
3397 *   plasticity_parameters.lookupValue("plasticity_on", system_parameters::plasticity_on);
3398 *   plasticity_parameters.lookupValue("failure_criterion", system_parameters::failure_criterion);
3399 *   plasticity_parameters.lookupValue("max_plastic_iterations", system_parameters::max_plastic_iterations);
3400 *   plasticity_parameters.lookupValue("smoothing_radius", system_parameters::smoothing_radius);
3401 *   }
3402 *   catch (const SettingNotFoundException &nfex)
3403 *   {
3404 *   cerr << "We've got a problem in the plasticity parameters block" << endl;
3405 *   }
3406 *  
3407 * @endcode
3408 *
3409 * Viscoelasticity parameters
3410 *
3411
3412 *
3413 *
3414 * @code
3415 *   try
3416 *   {
3417 *  
3418 *   const Setting &viscoelasticity_parameters = root["viscoelasticity_parameters"];
3419 *   viscoelasticity_parameters.lookupValue("initial_elastic_iterations", system_parameters::initial_elastic_iterations);
3420 *   viscoelasticity_parameters.lookupValue("elastic_time", system_parameters::elastic_time);
3421 *   viscoelasticity_parameters.lookupValue("viscous_time", system_parameters::viscous_time);
3422 *   viscoelasticity_parameters.lookupValue("initial_disp_target", system_parameters::initial_disp_target);
3423 *   viscoelasticity_parameters.lookupValue("final_disp_target", system_parameters::final_disp_target);
3424 *   viscoelasticity_parameters.lookupValue("current_time_interval", system_parameters::current_time_interval);
3425 *  
3426 *   system_parameters::viscous_time *= SECSINYEAR;
3427 *   }
3428 *   catch (const SettingNotFoundException &nfex)
3429 *   {
3430 *   cerr << "We've got a problem in the viscoelasticity parameters block" << endl;
3431 *   }
3432 *  
3433 * @endcode
3434 *
3435 * Mesh refinement parameters
3436 *
3437 * @code
3438 *   try
3439 *   {
3440 *  
3441 *   const Setting &mesh_refinement_parameters = root["mesh_refinement_parameters"];
3442 *   mesh_refinement_parameters.lookupValue("global_refinement", system_parameters::global_refinement);
3443 *   mesh_refinement_parameters.lookupValue("small_r_refinement", system_parameters::small_r_refinement);
3444 *   mesh_refinement_parameters.lookupValue("crustal_refinement", system_parameters::crustal_refinement);
3445 *   mesh_refinement_parameters.lookupValue("crust_refine_region", system_parameters::crust_refine_region);
3446 *   mesh_refinement_parameters.lookupValue("surface_refinement", system_parameters::surface_refinement);
3447 *   }
3448 *   catch (const SettingNotFoundException &nfex)
3449 *   {
3450 *   cerr << "We've got a problem in the mesh refinement parameters block" << endl;
3451 *   }
3452 *  
3453 * @endcode
3454 *
3455 * Solver parameters
3456 *
3457 * @code
3458 *   try
3459 *   {
3460 *   const Setting &solve_parameters = root["solve_parameters"];
3461 *   solve_parameters.lookupValue("iteration_coefficient", system_parameters::iteration_coefficient);
3462 *   solve_parameters.lookupValue("tolerance_coefficient", system_parameters::tolerance_coefficient);
3463 *  
3464 *  
3465 *   }
3466 *   catch (const SettingNotFoundException &nfex)
3467 *   {
3468 *   cerr << "We've got a problem in the solver parameters block" << endl;
3469 *   }
3470 *  
3471 * @endcode
3472 *
3473 * Time step parameters
3474 *
3475 * @code
3476 *   try
3477 *   {
3478 *   const Setting &time_step_parameters = root["time_step_parameters"];
3479 *   time_step_parameters.lookupValue("present_time", system_parameters::present_time);
3480 *   time_step_parameters.lookupValue("present_timestep", system_parameters::present_timestep);
3481 *   time_step_parameters.lookupValue("total_viscous_steps", system_parameters::total_viscous_steps);
3482 *   }
3483 *   catch (const SettingNotFoundException &nfex)
3484 *   {
3485 *   cerr << "We've got a problem in the time step parameters block" << endl;
3486 *   }
3487 *  
3488 *   write_config();
3489 *   }
3490 *   }
3491 *  
3492 *  
3493 *  
3494 *  
3495 *  
3496 * @endcode
3497
3498
3499<a name="ann-support_code/ellipsoid_fit.h"></a>
3500<h1>Annotated version of support_code/ellipsoid_fit.h</h1>
3501 *
3502 *
3503 *
3504 *
3505 * @code
3506 *   /* -----------------------------------------------------------------------------
3507 *   *
3508 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
3509 *   * Copyright (C) 2015 by Anton Ermakov
3510 *   *
3511 *   * This file is part of the deal.II code gallery.
3512 *   *
3513 *   * -----------------------------------------------------------------------------
3514 *   *
3515 *   * Author: antonermakov
3516 *   */
3517 *  
3518 *   #include <deal.II/base/quadrature_lib.h>
3519 *   #include <deal.II/base/function.h>
3520 *   #include <deal.II/base/logstream.h>
3521 *   #include <deal.II/lac/vector.h>
3522 *   #include <deal.II/lac/full_matrix.h>
3523 *   #include <deal.II/lac/sparse_matrix.h>
3524 *   #include <deal.II/lac/solver_cg.h>
3525 *   #include <deal.II/lac/precondition.h>
3526 *   #include <deal.II/grid/tria.h>
3527 *   #include <deal.II/dofs/dof_handler.h>
3528 *   #include <deal.II/grid/tria_accessor.h>
3529 *   #include <deal.II/grid/tria_iterator.h>
3530 *   #include <deal.II/dofs/dof_accessor.h>
3531 *   #include <deal.II/dofs/dof_tools.h>
3532 *   #include <deal.II/fe/fe_q.h>
3533 *   #include <deal.II/fe/fe_values.h>
3534 *   #include <deal.II/numerics/vector_tools.h>
3535 *   #include <deal.II/numerics/matrix_tools.h>
3536 *   #include <deal.II/numerics/data_out.h>
3537 *   #include <deal.II/grid/grid_in.h>
3538 *   #include <deal.II/grid/grid_out.h>
3539 *   #include <deal.II/base/point.h>
3540 *   #include <deal.II/grid/grid_generator.h>
3541 *  
3542 *   #include <fstream>
3543 *   #include <sstream>
3544 *   #include <iostream>
3545 *   #include <iomanip>
3546 *   #include <cstdlib>
3547 *  
3548 *  
3549 *   #include "local_math.h"
3550 *  
3551 *   using namespace dealii;
3552 *  
3553 *   template <int dim>
3554 *   class ellipsoid_fit
3555 *   {
3556 *   public:
3557 *   inline ellipsoid_fit (Triangulation<dim,dim> *pi)
3558 *   {
3559 *   p_triangulation = pi;
3560 *   };
3561 *   void compute_fit(std::vector<double> &ell, unsigned char bndry);
3562 *  
3563 *  
3564 *   private:
3565 *   Triangulation<dim,dim> *p_triangulation;
3566 *  
3567 *   };
3568 *  
3569 *  
3570 * @endcode
3571 *
3572 * This function computes ellipsoid fit to a set of vertices that lie on the
3573 * boundary_that_we_need
3574 *
3575 * @code
3576 *   template <int dim>
3577 *   void ellipsoid_fit<dim>::compute_fit(std::vector<double> &ell, unsigned char boundary_that_we_need)
3578 *   {
3579 *   typename Triangulation<dim>::active_cell_iterator cell = p_triangulation->begin_active();
3580 *   typename Triangulation<dim>::active_cell_iterator endc = p_triangulation->end();
3581 *  
3582 *   FullMatrix<double> A(p_triangulation->n_vertices(),dim);
3583 *   Vector<double> x(dim);
3584 *   Vector<double> b(p_triangulation->n_vertices());
3585 *  
3586 *   std::vector<bool> vertex_touched (p_triangulation->n_vertices(),
3587 *   false);
3588 *  
3589 *   unsigned int j = 0;
3590 *   unsigned char boundary_ids;
3591 *   std::vector<unsigned int> ind_bnry_row;
3592 *   std::vector<unsigned int> ind_bnry_col;
3593 *  
3594 * @endcode
3595 *
3596 * assemble the sensitivity matrix and r.h.s.
3597 *
3598 * @code
3599 *   for (; cell != endc; ++cell)
3600 *   {
3601 *   if (boundary_that_we_need != 0)
3602 *   cell->set_manifold_id(cell->material_id());
3603 *   for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
3604 *   {
3605 *   if (boundary_that_we_need == 0) //if this is the outer surface, then look for boundary ID 0; otherwise look for material ID change.
3606 *   {
3607 *   boundary_ids = cell->face(f)->boundary_id();
3608 *   if (boundary_ids == boundary_that_we_need)
3609 *   {
3610 *   for (unsigned int v = 0;
3611 *   v < GeometryInfo<dim>::vertices_per_face; ++v)
3612 *   if (vertex_touched[cell->face(f)->vertex_index(v)]
3613 *   == false)
3614 *   {
3615 *   vertex_touched[cell->face(f)->vertex_index(v)] =
3616 *   true;
3617 *   for (unsigned int i = 0; i < dim; ++i)
3618 *   {
3619 * @endcode
3620 *
3621 * stiffness matrix entry
3622 *
3623 * @code
3624 *   A(j, i) = pow(cell->face(f)->vertex(v)[i], 2);
3625 * @endcode
3626 *
3627 * r.h.s. entry
3628 *
3629 * @code
3630 *   b[j] = 1.0;
3631 * @endcode
3632 *
3633 * if mesh if not full: set the indicator
3634 *
3635 * @code
3636 *   }
3637 *   ind_bnry_row.push_back(j);
3638 *   j++;
3639 *   }
3640 *   }
3641 *   }
3642 *   else //find the faces that are at the boundary between materials, get the vertices, and write them into the stiffness matrix
3643 *   {
3644 *   if (cell->neighbor(f) != endc)
3645 *   {
3646 *   if (cell->material_id() != cell->neighbor(f)->material_id()) //finds face is at internal boundary
3647 *   {
3648 *   int high_mat_id = std::max(cell->material_id(),
3649 *   cell->neighbor(f)->material_id());
3650 *   if (high_mat_id == boundary_that_we_need) //finds faces at the correct internal boundary
3651 *   {
3652 *   for (unsigned int v = 0;
3653 *   v < GeometryInfo<dim>::vertices_per_face;
3654 *   ++v)
3655 *   if (vertex_touched[cell->face(f)->vertex_index(
3656 *   v)] == false)
3657 *   {
3658 *   vertex_touched[cell->face(f)->vertex_index(
3659 *   v)] = true;
3660 *   for (unsigned int i = 0; i < dim; ++i)
3661 *   {
3662 * @endcode
3663 *
3664 * stiffness matrix entry
3665 *
3666 * @code
3667 *   A(j, i) = pow(
3668 *   cell->face(f)->vertex(v)[i], 2);
3669 * @endcode
3670 *
3671 * r.h.s. entry
3672 *
3673 * @code
3674 *   b[j] = 1.0;
3675 * @endcode
3676 *
3677 * if mesh if not full: set the indicator
3678 *
3679 * @code
3680 *   }
3681 *   ind_bnry_row.push_back(j);
3682 *   j++;
3683 *   }
3684 *   }
3685 *   }
3686 *   }
3687 *   }
3688 *   }
3689 *   }
3690 *   if (ind_bnry_row.size()>0)
3691 *   {
3692 *  
3693 * @endcode
3694 *
3695 * maxtrix A'*A and vector A'*b; A'*A*x = A'*b -- normal system of equations
3696 *
3697 * @code
3698 *   FullMatrix<double> AtA(dim,dim);
3699 *   Vector<double> Atb(dim);
3700 *  
3701 *   FullMatrix<double> A_out(ind_bnry_row.size(),dim);
3702 *   Vector<double> b_out(ind_bnry_row.size());
3703 *  
3704 *   for (unsigned int i=0; i<dim; ++i)
3705 *   ind_bnry_col.push_back(i);
3706 *  
3707 *   for (unsigned int i=0; i<ind_bnry_row.size(); ++i)
3708 *   b_out(i) = 1;
3709 *  
3710 *   A_out.extract_submatrix_from(A, ind_bnry_row, ind_bnry_col);
3711 *   A_out.Tmmult(AtA,A_out,true);
3712 *   A_out.Tvmult(Atb,b_out,true);
3713 *  
3714 * @endcode
3715 *
3716 * solve normal system of equations
3717 *
3718 * @code
3719 *   SolverControl solver_control (1000, 1e-12);
3720 *   SolverCG<> solver (solver_control);
3721 *   solver.solve (AtA, x, Atb, PreconditionIdentity());
3722 *  
3723 * @endcode
3724 *
3725 * find ellipsoidal axes
3726 *
3727 * @code
3728 *   for (unsigned int i=0; i<dim; ++i)
3729 *   ell.push_back(sqrt(1.0/x[i]));
3730 *   }
3731 *   else
3732 *   std::cerr << "fit_ellipsoid: no points to fit" << std::endl;
3733 *  
3734 *   }
3735 *  
3736 * @endcode
3737
3738
3739<a name="ann-support_code/ellipsoid_grav.h"></a>
3740<h1>Annotated version of support_code/ellipsoid_grav.h</h1>
3741 *
3742 *
3743 *
3744 *
3745 * @code
3746 *   /* -----------------------------------------------------------------------------
3747 *   *
3748 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
3749 *   * Copyright (C) 2012 by Roger R. Fu
3750 *   *
3751 *   * This file is part of the deal.II code gallery.
3752 *   *
3753 *   * -----------------------------------------------------------------------------
3754 *   *
3755 *   * Author: Roger R. Fu
3756 *   * Reference Pohanka 2011, Contrib. Geophys. Geodes.
3757 *   */
3758 *  
3759 *   #include <math.h>
3760 *   #include <deal.II/base/point.h>
3761 *   #include <fstream>
3762 *   #include <iostream>
3763 *  
3764 *   namespace A_Grav_namespace
3765 *   {
3766 *   namespace system_parameters
3767 *   {
3768 *   double mantle_rho;
3769 *   double core_rho;
3770 *   double excess_rho;
3771 *   double r_eq;
3772 *   double r_polar;
3773 *  
3774 *   double r_core_eq;
3775 *   double r_core_polar;
3776 *   }
3777 *  
3778 *   template <int dim>
3779 *   class AnalyticGravity
3780 *   {
3781 *   public:
3782 *   void setup_vars (std::vector<double> v);
3783 *   void get_gravity (const ::Point<dim> &p, std::vector<double> &g);
3784 *  
3785 *   private:
3786 *   double ecc;
3787 *   double eV;
3788 *   double ke;
3789 *   double r00;
3790 *   double r01;
3791 *   double r11;
3792 *   double ecc_c;
3793 *   double eV_c;
3794 *   double ke_c;
3795 *   double r00_c;
3796 *   double r01_c;
3797 *   double r11_c;
3798 *   double g_coeff;
3799 *   double g_coeff_c;
3800 *   };
3801 *  
3802 *   template <int dim>
3803 *   void AnalyticGravity<dim>::get_gravity (const ::Point<dim> &p, std::vector<double> &g)
3804 *   {
3805 *   double rsph = std::sqrt(p[0] * p[0] + p[1] * p[1]);
3806 *   double thetasph = std::atan2(p[0], p[1]);
3807 *   double costhetasph = std::cos(thetasph);
3808 *  
3809 * @endcode
3810 *
3811 * convert to elliptical coordinates for silicates
3812 *
3813 * @code
3814 *   double stemp = std::sqrt((rsph * rsph - eV * eV + std::sqrt((rsph * rsph - eV * eV) * (rsph * rsph - eV * eV)
3815 *   + 4 * eV * eV * rsph * rsph * costhetasph *costhetasph)) / 2);
3816 *   double vout = stemp / system_parameters::r_eq / std::sqrt(1 - ecc * ecc);
3817 *   double eout = std::acos(rsph * costhetasph / stemp);
3818 *  
3819 * @endcode
3820 *
3821 * convert to elliptical coordinates for core correction
3822 *
3823 * @code
3824 *   double stemp_c = std::sqrt((rsph * rsph - eV_c * eV_c + std::sqrt((rsph * rsph - eV_c * eV_c) * (rsph * rsph - eV_c * eV_c)
3825 *   + 4 * eV_c * eV_c * rsph * rsph * costhetasph *costhetasph)) / 2);
3826 *   double vout_c = stemp_c / system_parameters::r_core_eq / std::sqrt(1 - ecc_c * ecc_c);
3827 *   double eout_c = std::acos(rsph * costhetasph / stemp_c);
3828 *  
3829 * @endcode
3830 *
3831 * shell contribution
3832 *
3833 * @code
3834 *   g[0] = g_coeff * r11 * std::sqrt((1 - ecc * ecc) * vout * vout + ecc * ecc) * std::sin(eout);
3835 *   g[1] = g_coeff * r01 * vout * std::cos(eout) / std::sqrt(1 - ecc * ecc);
3836 *  
3837 * @endcode
3838 *
3839 * core contribution
3840 *
3841 * @code
3842 *   double expected_y = system_parameters::r_core_polar * std::sqrt(1 -
3843 *   (p[0] * p[0] / system_parameters::r_core_eq / system_parameters::r_core_eq));
3844 *  
3845 *  
3846 *   if (p[1] <= expected_y)
3847 *   {
3848 *   g[0] += g_coeff_c * r11_c * std::sqrt((1 - ecc_c * ecc_c) * vout_c * vout_c + ecc_c * ecc_c) * std::sin(eout_c);
3849 *   g[1] += g_coeff_c * r01_c * vout_c * std::cos(eout_c) / std::sqrt(1 - ecc_c * ecc_c);
3850 *   }
3851 *   else
3852 *   {
3853 *   double g_coeff_co = - 2.795007963255562e-10 * system_parameters::excess_rho * system_parameters::r_core_eq
3854 *   / vout_c / vout_c;
3855 *   double r00_co = 0;
3856 *   double r01_co = 0;
3857 *   double r11_co = 0;
3858 *  
3859 *   if (system_parameters::r_core_polar == system_parameters::r_core_eq)
3860 *   {
3861 *   r00_co = 1;
3862 *   r01_co = 1;
3863 *   r11_co = 1;
3864 *   }
3865 *   else
3866 *   {
3867 *   r00_co = ke_c * vout_c * std::atan2(1, ke_c * vout_c);
3868 *   double ke_co2 = ke_c * ke_c * vout_c * vout_c;
3869 *   r01_co = 3 * ke_co2 * (1 - r00_co);
3870 *   r11_co = 3 * ((ke_co2 + 1) * r00_co - ke_co2) / 2;
3871 *   }
3872 *   g[0] += g_coeff_co * vout_c * r11_co / std::sqrt((1 - ecc_c* ecc_c) * vout_c * vout_c + ecc_c * ecc_c) * std::sin(eout_c);
3873 *   g[1] += g_coeff_co * r01_co * std::cos(eout_c) / std::sqrt(1 - ecc_c * ecc_c);
3874 *   }
3875 *   }
3876 *  
3877 *   template <int dim>
3878 *   void AnalyticGravity<dim>::setup_vars (std::vector<double> v)
3879 *   {
3880 *   system_parameters::r_eq = v[0];
3881 *   system_parameters::r_polar = v[1];
3882 *   system_parameters::r_core_eq = v[2];
3883 *   system_parameters::r_core_polar = v[3];
3884 *   system_parameters::mantle_rho = v[4];
3885 *   system_parameters::core_rho = v[5];
3886 *   system_parameters::excess_rho = system_parameters::core_rho - system_parameters::mantle_rho;
3887 *  
3888 *  
3889 * @endcode
3890 *
3891 * Shell
3892 *
3893 * @code
3894 *   if (system_parameters::r_polar > system_parameters::r_eq)
3895 *   {
3896 * @endcode
3897 *
3898 * This makes the gravity field nearly that of a sphere in case the body becomes prolate
3899 *
3900 * @code
3901 *   std::cout << "\nWarning: The model body has become prolate. \n";
3902 *   ecc = 0.001;
3903 *   }
3904 *   else
3905 *   {
3906 *   ecc = std::sqrt(1 - (system_parameters::r_polar * system_parameters::r_polar / system_parameters::r_eq / system_parameters::r_eq));
3907 *   }
3908 *  
3909 *   eV = ecc * system_parameters::r_eq;
3910 *   ke = std::sqrt(1 - (ecc * ecc)) / ecc;
3911 *   r00 = ke * std::atan2(1, ke);
3912 *   double ke2 = ke * ke;
3913 *   r01 = 3 * ke2 * (1 - r00);
3914 *   r11 = 3 * ((ke2 + 1) * r00 - ke2) / 2;
3915 *   g_coeff = - 2.795007963255562e-10 * system_parameters::mantle_rho * system_parameters::r_eq;
3916 *  
3917 * @endcode
3918 *
3919 * Core
3920 *
3921 * @code
3922 *   if (system_parameters::r_core_polar > system_parameters::r_core_eq)
3923 *   {
3924 *   std::cout << "\nWarning: The model core has become prolate. \n";
3925 *   ecc_c = 0.001;
3926 *   }
3927 *   else
3928 *   {
3929 *   ecc_c = std::sqrt(1 - (system_parameters::r_core_polar * system_parameters::r_core_polar / system_parameters::r_core_eq / system_parameters::r_core_eq));
3930 *   }
3931 *   eV_c = ecc_c * system_parameters::r_core_eq;
3932 *   if (system_parameters::r_core_polar == system_parameters::r_core_eq)
3933 *   {
3934 *   ke_c = 1;
3935 *   r00_c = 1;
3936 *   r01_c = 1;
3937 *   r11_c = 1;
3938 *   g_coeff_c = - 2.795007963255562e-10 * system_parameters::excess_rho * system_parameters::r_core_eq;
3939 *   }
3940 *   else
3941 *   {
3942 *   ke_c = std::sqrt(1 - (ecc_c * ecc_c)) / ecc_c;
3943 *   r00_c = ke_c * std::atan2(1, ke_c);
3944 *   double ke2_c = ke_c * ke_c;
3945 *   r01_c = 3 * ke2_c * (1 - r00_c);
3946 *   r11_c = 3 * ((ke2_c + 1) * r00_c - ke2_c) / 2;
3947 *   g_coeff_c = - 2.795007963255562e-10 * system_parameters::excess_rho * system_parameters::r_core_eq;
3948 *   }
3949 * @endcode
3950 *
3951 * std::cout << "Loaded variables: ecc = " << ecc_c << " ke = " << ke_c << " r00 = " << r00_c << " r01 = " << r01_c << " r11 = " << r11_c << "\n";
3952 *
3953 * @code
3954 *   }
3955 *   }
3956 * @endcode
3957
3958
3959<a name="ann-support_code/local_math.h"></a>
3960<h1>Annotated version of support_code/local_math.h</h1>
3961 *
3962 *
3963 *
3964 *
3965 * @code
3966 *   /* -----------------------------------------------------------------------------
3967 *   *
3968 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
3969 *   * Copyright (C) 2013 by Anton Ermakov
3970 *   *
3971 *   * This file is part of the deal.II code gallery.
3972 *   *
3973 *   * -----------------------------------------------------------------------------
3974 *   *
3975 *   * Author: antonermakov
3976 *   */
3977 *  
3978 *   #ifndef LOCAL_MATH_
3979 *   #define LOCAL_MATH_
3980 *  
3981 *   #define PI 3.14159265358979323846
3982 *   #define TWOPI 6.283185307179586476925287
3983 *   #define SECSINYEAR 3.155692608e+07
3984 *  
3985 * @endcode
3986 *
3987 * double factorial(int n)
3988 * {
3989 * if(n == 0) {
3990 * return(1.);
3991 * } else if(n == 1) {
3992 * return(1.);
3993 * } else if(n == 2) {
3994 * return(2.);
3995 * } else if(n == 3) {
3996 * return(6.);
3997 * } else if(n == 4) {
3998 * return(24.);
3999 * } else {
4000 * exit(-1);
4001 * }
4002 * }
4003 *
4004
4005 *
4006 * double fudge(int m)
4007 * {
4008 * if(m == 0) {
4009 * return(1.0);
4010 * } else {
4011 * return(2.0);
4012 * }
4013 * }
4014 *
4015
4016 *
4017 *
4018
4019 *
4020 * double sign(double x)
4021 * {
4022 * if(x > 0) {
4023 * return(1.0);
4024 * } else if(x < 0.0) {
4025 * return(-1.0);
4026 * } else {
4027 * return(0.0);
4028 * }
4029 * }
4030 *
4031
4032 *
4033 * double pv0(double x)
4034 * {
4035 * double ans;
4036 *
4037
4038 *
4039 * ans = x - TWOPI*floor(x/TWOPI);
4040 * if(ans > TWOPI/2.0) {
4041 * ans = ans - TWOPI;
4042 * }
4043 *
4044
4045 *
4046 * return(ans);
4047 * }
4048 *
4049
4050 *
4051 *
4052 * @code
4053 *   /* assumes l=2 */
4054 * @endcode
4055 *
4056 * double System::Plm(int m, double x)
4057 * {
4058 * if(m == 0) {
4059 * return(1.5*x*x - 0.5);
4060 * } else if(m == 1) {
4061 * return(3.0*x*sqrt(1.0 - x*x));
4062 * } else if(m == 2) {
4063 * return(3.0 - 3.0*x*x);
4064 * } else {
4065 * exit(-1);
4066 * }
4067 * }
4068 *
4069
4070 *
4071 * double System::DP(int m, double x)
4072 * {
4073 * if(m == 0) {
4074 * return(3.0*x);
4075 * } else if(m == 1) {
4076 * return((3.0 - 6.0*x*x)/sqrt(1.0 - x*x));
4077 * } else if(m == 2) {
4078 * return(- 6.0*x);
4079 * } else {
4080 * exit(-1);
4081 * }
4082 * }
4083 *
4084
4085 *
4086 *
4087
4088 *
4089 *
4090 * @code
4091 *   #endif /* LOCALMATH_H */
4092 *  
4093 * @endcode
4094
4095
4096*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
Definition fe_q.h:554
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< RangeNumberType > > &values) const
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
Definition point.h:113
typename SparseLUDecomposition< number >::AdditionalData AdditionalData
Definition sparse_ilu.h:79
void initialize(const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
EnableObserverPointer Subscriptor
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
LinearOperator< Range_2, Domain_2, Payload > schur_complement(const LinearOperator< Domain_1, Range_1, Payload > &A_inv, const LinearOperator< Range_1, Domain_2, Payload > &B, const LinearOperator< Range_2, Domain_1, Payload > &C, const LinearOperator< Range_2, Domain_2, Payload > &D)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
const bool IsBlockVector< VectorType >::value
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void compute_no_normal_flux_constraints(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, AffineConstraints< number > &constraints, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const bool use_manifold_for_normal=true)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
std::vector< value_type > l2_norm(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
Expression floor(const Expression &x)
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >())
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
@ valid
Iterator points to a valid object.
@ matrix
Contents is actually a matrix.
constexpr types::blas_int zero
constexpr char A
constexpr types::blas_int one
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:193
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Number angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:466
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
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)
bool check(const ConstraintKinds kind_in, const unsigned int dim)
int(& functions)(const void *v1, const void *v2)
constexpr double PI
Definition numbers.h:239
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
inline ::VectorizedArray< Number, width > atan(const ::VectorizedArray< Number, width > &x)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int material_id
Definition types.h:184
ObserverPointer< T, P > SmartPointer