Flow123d  JS_before_hm-983-gccfd4b3
elasticity.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file transport_dg.cc
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #include "system/sys_profiler.hh"
20 #include "mechanics/elasticity.hh"
21 
22 #include "io/output_time.hh"
24 #include "fem/fe_values.hh"
25 #include "fem/fe_p.hh"
26 #include "fem/fe_rt.hh"
27 #include "fem/fe_system.hh"
28 #include "fields/field_fe.hh"
29 #include "la/linsys_PETSC.hh"
30 #include "coupling/balance.hh"
31 #include "mesh/neighbours.h"
32 
33 #include "fields/multi_field.hh"
34 #include "fields/generic_field.hh"
35 #include "input/factory.hh"
36 
37 
38 
39 
40 using namespace Input::Type;
41 
42 
44  std::string equation_name = std::string(Elasticity::EqData::name()) + "_FE";
45  return IT::Record(
46  std::string(equation_name),
47  "FEM for linear elasticity.")
49  .declare_key("balance", Balance::get_input_type(), Default("{}"),
50  "Settings for computing balance.")
52  "Parameters of output stream.")
54  "Linear solver for elasticity.")
55  .declare_key("input_fields", Array(
57  .make_field_descriptor_type(equation_name)),
59  "Input fields of the equation.")
60  .declare_key("output",
61  EqData().output_fields.make_output_type(equation_name, ""),
62  IT::Default("{ \"fields\": [ "+Elasticity::EqData::default_output_field()+" ] }"),
63  "Setting of the field output.")
64  .close();
65 }
66 
67 const int Elasticity::registrar =
68  Input::register_class< Elasticity, Mesh &, const Input::Record>(std::string(Elasticity::EqData::name()) + "_FE") +
70 
71 
72 
73 namespace Mechanics {
74 
75 FEObjects::FEObjects(Mesh *mesh_, unsigned int fe_order)
76 : q_(QGauss::make_array(2))
77 {
78  switch (fe_order)
79  {
80  case 1: {
81  MixedPtr<FE_P> fe_p(1);
82  fe_ = mixed_fe_system(fe_p, FEVector, 3);
83  break;
84  }
85  default:
86  xprintf(PrgErr, "Unsupported polynomial order %d for finite elements in Elasticity", fe_order);
87  break;
88  }
89 
90 
91  ds_ = std::make_shared<EqualOrderDiscreteSpace>(mesh_, fe_);
92  dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
93 
94  dh_->distribute_dofs(ds_);
95 
96 
97  MixedPtr<FE_P_disc> fe_p(0);
98  dh_scalar_ = make_shared<DOFHandlerMultiDim>(*mesh_);
99  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_p);
100  dh_scalar_->distribute_dofs(ds);
101 
102 
104  dh_tensor_ = make_shared<DOFHandlerMultiDim>(*mesh_);
105  std::shared_ptr<DiscreteSpace> dst = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_t);
106  dh_tensor_->distribute_dofs(dst);
107 }
108 
109 
111 {
112 }
113 
114 template<> std::shared_ptr<FiniteElement<0>> FEObjects::fe<0>() { return fe_[0_d]; }
115 template<> std::shared_ptr<FiniteElement<1>> FEObjects::fe<1>() { return fe_[1_d]; }
116 template<> std::shared_ptr<FiniteElement<2>> FEObjects::fe<2>() { return fe_[2_d]; }
117 template<> std::shared_ptr<FiniteElement<3>> FEObjects::fe<3>() { return fe_[3_d]; }
118 
119 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh() { return dh_; }
120 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh_scalar() { return dh_scalar_; }
121 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh_tensor() { return dh_tensor_; }
122 
123 
124 
125 } // namespace Mechanics
126 
127 
128 double lame_mu(double young, double poisson)
129 {
130  return young*0.5/(poisson+1.);
131 }
132 
133 
134 double lame_lambda(double young, double poisson)
135 {
136  return young*poisson/((poisson+1.)*(1.-2.*poisson));
137 }
138 
139 
140 
141 
142 
143 
145  return Selection("Elasticity_BC_Type", "Types of boundary conditions for mechanics.")
146  .add_value(bc_type_displacement, "displacement",
147  "Prescribed displacement.")
148  .add_value(bc_type_displacement_normal, "displacement_n",
149  "Prescribed displacement in the normal direction to the boundary.")
150  .add_value(bc_type_traction, "traction",
151  "Prescribed traction.")
152  .close();
153 }
154 
155 
157 {
158  *this+=bc_type
159  .name("bc_type")
160  .description(
161  "Type of boundary condition.")
162  .units( UnitSI::dimensionless() )
163  .input_default("\"traction\"")
164  .input_selection( get_bc_type_selection() )
166 
167  *this+=bc_displacement
168  .name("bc_displacement")
169  .description("Prescribed displacement on boundary.")
170  .units( UnitSI().m() )
171  .input_default("0.0")
172  .flags_add(in_rhs);
173 
174  *this+=bc_traction
175  .name("bc_traction")
176  .description("Prescribed traction on boundary.")
177  .units( UnitSI().Pa() )
178  .input_default("0.0")
179  .flags_add(in_rhs);
180 
181  *this+=load
182  .name("load")
183  .description("Prescribed bulk load.")
184  .units( UnitSI().N().m(-3) )
185  .input_default("0.0")
186  .flags_add(in_rhs);
187 
188  *this+=young_modulus
189  .name("young_modulus")
190  .description("Young's modulus.")
191  .units( UnitSI().Pa() )
192  .input_default("0.0")
193  .flags_add(in_main_matrix & in_rhs);
194 
195  *this+=poisson_ratio
196  .name("poisson_ratio")
197  .description("Poisson's ratio.")
198  .units( UnitSI().dimensionless() )
199  .input_default("0.0")
200  .flags_add(in_main_matrix & in_rhs);
201 
202  *this+=fracture_sigma
203  .name("fracture_sigma")
204  .description(
205  "Coefficient of transfer of forces through fractures.")
206  .units( UnitSI::dimensionless() )
207  .input_default("1.0")
208  .flags_add(in_main_matrix & in_rhs);
209 
210  *this += region_id.name("region_id")
211  .units( UnitSI::dimensionless())
213 
214  *this += subdomain.name("subdomain")
215  .units( UnitSI::dimensionless() )
217 
218  *this+=cross_section
219  .name("cross_section")
220  .units( UnitSI().m(3).md() )
221  .flags(input_copy & in_time_term & in_main_matrix & in_rhs);
222 
223  *this+=potential_load
224  .name("potential_load")
225  .units( UnitSI().m() )
226  .flags(input_copy & in_rhs);
227 
228  *this+=output_field
229  .name("displacement")
230  .description("Displacement vector field output.")
231  .units( UnitSI().m() )
232  .flags(equation_result);
233 
234  *this += output_stress
235  .name("stress")
236  .description("Stress tensor output.")
237  .units( UnitSI().Pa() )
238  .flags(equation_result);
239 
240  *this += output_von_mises_stress
241  .name("von_mises_stress")
242  .description("von Mises stress output.")
243  .units( UnitSI().Pa() )
244  .flags(equation_result);
245 
246  *this += output_cross_section
247  .name("cross_section_updated")
248  .description("Cross-section after deformation - output.")
249  .units( UnitSI().m() )
250  .flags(equation_result);
251 
252  *this += output_divergence
253  .name("displacement_divergence")
254  .description("Displacement divergence output.")
255  .units( UnitSI().dimensionless() )
256  .flags(equation_result);
257 
258  // add all input fields to the output list
259  output_fields += *this;
260 
261 }
262 
263 Elasticity::Elasticity(Mesh & init_mesh, const Input::Record in_rec, TimeGovernor *tm)
264  : EquationBase(init_mesh, in_rec),
265  input_rec(in_rec),
266  allocation_done(false)
267 {
268  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
269  // due to constexpr optimization.
271 
272  this->eq_data_ = &data_;
273 
274  auto time_rec = in_rec.val<Input::Record>("time");
275  if (tm == nullptr)
276  {
277  time_ = new TimeGovernor(time_rec);
278  }
279  else
280  {
281  TimeGovernor time_from_rec(time_rec);
282  ASSERT( time_from_rec.is_default() ).error("Duplicate key 'time', time in elasticity is already initialized from parent class!");
283  time_ = tm;
284  }
285 
286 
287  // Set up physical parameters.
288  data_.set_mesh(init_mesh);
291 
292  // create finite element structures and distribute DOFs
293  feo = new Mechanics::FEObjects(mesh_, 1);
294  DebugOut().fmt("Mechanics: solution size {}\n", feo->dh()->n_global_dofs());
295 
296 }
297 
298 
300 {
301  output_stream_ = OutputTime::create_output_stream("mechanics", input_rec.val<Input::Record>("output_stream"), time().get_unit_string());
302 
303  data_.set_components({"displacement"});
304  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), time() );
305 
306 // balance_ = std::make_shared<Balance>("mechanics", mesh_);
307 // balance_->init_from_input(input_rec.val<Input::Record>("balance"), *time_);
308  // initialization of balance object
309 // data_.balance_idx_ = {
310 // balance_->add_quantity("force_x"),
311 // balance_->add_quantity("force_y"),
312 // balance_->add_quantity("force_z")
313 // };
314 // balance_->units(UnitSI().kg().m().s(-2));
315 
316  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
317  data_.output_field_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(feo->dh());
319 
320  // setup output stress
321  data_.output_stress_ptr = create_field_fe<3, FieldValue<3>::TensorFixed>(feo->dh_tensor());
323 
324  // setup output von Mises stress
325  data_.output_von_mises_stress_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
327 
328  // setup output cross-section
329  data_.output_cross_section_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
331 
332  // setup output divergence
333  data_.output_div_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
335 
337 
338  // set time marks for writing the output
340 
341  // equation default PETSc solver options
342  std::string petsc_default_opts;
343  petsc_default_opts = "-ksp_type cg -pc_type hypre -pc_hypre_type boomeramg";
344 
345  // allocate matrix and vector structures
346  ls = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
348  ls->set_solution(data_.output_field_ptr->vec().petsc_vec());
349 
350  // initialization of balance object
351 // balance_->allocate(feo->dh()->distr()->lsize(),
352 // max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
353 
354 }
355 
356 
358 {
359 // delete time_;
360 
361  MatDestroy(&stiffness_matrix);
362  VecDestroy(&rhs);
363  delete ls;
364  delete feo;
365 
366 }
367 
368 
369 
371 {
372  // update ghost values of solution vector
373  data_.output_field_ptr->vec().local_to_ghost_begin();
374  data_.output_field_ptr->vec().local_to_ghost_end();
375 
376  // compute new output fields depending on solution (stress, divergence etc.)
377  data_.output_stress_ptr->vec().zero_entries();
378  data_.output_cross_section_ptr->vec().zero_entries();
379  data_.output_div_ptr->vec().zero_entries();
380  compute_output_fields<1>();
381  compute_output_fields<2>();
382  compute_output_fields<3>();
383 
384  // update ghost values of computed fields
385  data_.output_stress_ptr->vec().local_to_ghost_begin();
386  data_.output_stress_ptr->vec().local_to_ghost_end();
387  data_.output_von_mises_stress_ptr->vec().local_to_ghost_begin();
388  data_.output_von_mises_stress_ptr->vec().local_to_ghost_end();
389  data_.output_cross_section_ptr->vec().local_to_ghost_begin();
390  data_.output_cross_section_ptr->vec().local_to_ghost_end();
391  data_.output_div_ptr->vec().local_to_ghost_begin();
392  data_.output_div_ptr->vec().local_to_ghost_end();
393 }
394 
395 
396 
397 
399 {
403  std::stringstream ss; // print warning message with table of uninitialized fields
404  if ( FieldCommon::print_message_table(ss, "mechanics") ) {
405  WarningOut() << ss.str();
406  }
407 
408  ( (LinSys_PETSC *)ls )->set_initial_guess_nonzero();
409 
410  // check first time assembly - needs preallocation
412 
413 
414  // after preallocation we assemble the matrices and vectors required for balance of forces
415 // for (auto subst_idx : data_.balance_idx_)
416 // balance_->calculate_instant(subst_idx, ls->get_solution());
417 
418 // update_solution();
420  MatSetOption(*ls->get_matrix(), MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
421  ls->mat_zero_entries();
422  ls->rhs_zero_entries();
424  assemble_rhs();
425  ls->finish_assembly();
426  LinSys::SolveInfo si = ls->solve();
427  MessageOut().fmt("[mech solver] lin. it: {}, reason: {}, residual: {}\n",
429  output_data();
430 }
431 
432 
433 
435 {
436  // preallocate system matrix
437  ls->start_allocation();
438  stiffness_matrix = NULL;
439  rhs = NULL;
440 
442  assemble_rhs();
443 
444  allocation_done = true;
445 }
446 
447 
448 
449 
451 {
452  START_TIMER("DG-ONE STEP");
453 
454  next_time();
456 
458 
459  output_data();
460 
461  END_TIMER("DG-ONE STEP");
462 }
463 
465 {
466  time_->next_time();
467  time_->view("MECH");
468 
469 }
470 
471 
472 
474 {
475  START_TIMER("data reinit");
477  END_TIMER("data reinit");
478 
479  // assemble stiffness matrix
480  if (stiffness_matrix == NULL
482  {
483  DebugOut() << "Mechanics: Assembling matrix.\n";
485  ls->mat_zero_entries();
487  ls->finish_assembly();
488 
489  if (stiffness_matrix == NULL)
490  MatConvert(*( ls->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix);
491  else
492  MatCopy(*( ls->get_matrix() ), stiffness_matrix, DIFFERENT_NONZERO_PATTERN);
493  }
494 
495  // assemble right hand side (due to sources and boundary conditions)
496  if (rhs == NULL
497  || data_.subset(FieldFlag::in_rhs).changed())
498  {
499  DebugOut() << "Mechanics: Assembling right hand side.\n";
501  ls->rhs_zero_entries();
502  assemble_rhs();
503  ls->finish_assembly();
504 
505  if (rhs == nullptr) VecDuplicate(*( ls->get_rhs() ), &rhs);
506  VecCopy(*( ls->get_rhs() ), rhs);
507  }
508 
509  START_TIMER("solve");
510  LinSys::SolveInfo si = ls->solve();
511  MessageOut().fmt("[mech solver] lin. it: {}, reason: {}, residual: {}\n",
513  END_TIMER("solve");
514 }
515 
516 
517 
518 
519 
520 
522 {
523  START_TIMER("MECH-OUTPUT");
524 
525  // gather the solution from all processors
526  data_.output_fields.set_time( this->time().step(), LimitSide::left);
527  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
529  data_.output_fields.output(this->time().step());
530  output_stream_->write_time_frame();
531 
532 // START_TIMER("MECH-balance");
533 // balance_->calculate_instant(subst_idx, ls->get_solution());
534 // balance_->output();
535 // END_TIMER("MECH-balance");
536 
537  END_TIMER("MECH-OUTPUT");
538 }
539 
540 
541 template<unsigned int dim>
543 {
544  QGauss q(dim, 0), q_sub(dim-1, 0);
545  FEValues<3> fv(q, *feo->fe<dim>(),
547  FEValues<3> fsv(q_sub, *feo->fe<dim>(),
549  const unsigned int ndofs = feo->fe<dim>()->n_dofs();
550  auto vec = fv.vector_view(0);
551  auto vec_side = fsv.vector_view(0);
552  VectorMPI output_vec = data_.output_field_ptr->vec();
553  VectorMPI output_stress_vec = data_.output_stress_ptr->vec();
554  VectorMPI output_von_mises_stress_vec = data_.output_von_mises_stress_ptr->vec();
555  VectorMPI output_cross_sec_vec = data_.output_cross_section_ptr->vec();
556  VectorMPI output_div_vec = data_.output_div_ptr->vec();
557 
558  DHCellAccessor cell_tensor = *feo->dh_tensor()->own_range().begin();
559  DHCellAccessor cell_scalar = *feo->dh_scalar()->own_range().begin();
560  for (auto cell : feo->dh()->own_range())
561  {
562  if (cell.dim() == dim)
563  {
564  auto elm = cell.elm();
565 
566  double poisson = data_.poisson_ratio.value(elm.centre(), elm);
567  double young = data_.young_modulus.value(elm.centre(), elm);
568  double mu = lame_mu(young, poisson);
569  double lambda = lame_lambda(young, poisson);
570 
571  fv.reinit(elm);
572  LocDofVec dof_indices = cell.get_loc_dof_indices();
573  LocDofVec dof_indices_scalar = cell_scalar.get_loc_dof_indices();
574  LocDofVec dof_indices_tensor = cell_tensor.get_loc_dof_indices();
575 
576  arma::mat33 stress = arma::zeros(3,3);
577  double div = 0;
578  for (unsigned int i=0; i<ndofs; i++)
579  {
580  stress += (2*mu*vec.sym_grad(i,0) + lambda*vec.divergence(i,0)*arma::eye(3,3))*output_vec[dof_indices[i]];
581  div += vec.divergence(i,0)*output_vec[dof_indices[i]];
582  }
583 
584  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
585  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
586  output_div_vec[dof_indices_scalar[0]] += div;
587 
588  for (unsigned int i=0; i<3; i++)
589  for (unsigned int j=0; j<3; j++)
590  output_stress_vec[dof_indices_tensor[i*3+j]] += stress(i,j);
591  output_von_mises_stress_vec[dof_indices_scalar[0]] = von_mises_stress;
592 
593  output_cross_sec_vec[dof_indices_scalar[0]] += data_.cross_section.value(fv.point(0), elm);
594  }
595  else if (cell.dim() == dim-1)
596  {
597  auto elm = cell.elm();
598  double normal_displacement = 0;
599  double csection = data_.cross_section.value(fsv.point(0), elm);
600  arma::mat33 normal_stress = arma::zeros(3,3);
601 
602  double poisson = data_.poisson_ratio.value(elm.centre(), elm);
603  double young = data_.young_modulus.value(elm.centre(), elm);
604  double mu = lame_mu(young, poisson);
605  double lambda = lame_lambda(young, poisson);
606 
607  for (unsigned int inb=0; inb<elm->n_neighs_vb(); inb++)
608  {
609  auto side = elm->neigh_vb[inb]->side();
610  auto cell_side = side->element();
611  fsv.reinit(*side);
612  LocDofVec side_dof_indices =
613  feo->dh()->cell_accessor_from_element(cell_side.idx()).get_loc_dof_indices();
614 
615  for (unsigned int i=0; i<ndofs; i++)
616  {
617  normal_displacement -= arma::dot(vec_side.value(i,0)*output_vec[side_dof_indices[i]], fsv.normal_vector(0));
618  arma::mat33 grad = -arma::kron(vec_side.value(i,0)*output_vec[side_dof_indices[i]], fsv.normal_vector(0).t()) / csection;
619  normal_stress += mu*(grad+grad.t()) + lambda*arma::trace(grad)*arma::eye(3,3);
620  }
621  }
622  LocDofVec dof_indices_scalar = cell_scalar.get_loc_dof_indices();
623  LocDofVec dof_indices_tensor = cell_tensor.get_loc_dof_indices();
624  for (unsigned int i=0; i<3; i++)
625  for (unsigned int j=0; j<3; j++)
626  output_stress_vec[dof_indices_tensor[i*3+j]] += normal_stress(i,j);
627  output_cross_sec_vec[dof_indices_scalar[0]] += normal_displacement;
628  output_div_vec[dof_indices_scalar[0]] += normal_displacement / csection;
629  }
630  cell_scalar.inc();
631  cell_tensor.inc();
632  }
633 }
634 
635 
636 
637 
638 
640 {
641 // if (balance_->cumulative())
642 // {
643 // balance_->calculate_cumulative(subst_idx, ls->get_solution());
644 // }
645 }
646 
647 
648 
649 
650 
651 
652 
653 
655 {
656  START_TIMER("assemble_stiffness");
657  START_TIMER("assemble_volume_integrals");
658  assemble_volume_integrals<1>();
659  assemble_volume_integrals<2>();
660  assemble_volume_integrals<3>();
661  END_TIMER("assemble_volume_integrals");
662 
663  START_TIMER("assemble_fluxes_boundary");
664  assemble_fluxes_boundary<1>();
665  assemble_fluxes_boundary<2>();
666  assemble_fluxes_boundary<3>();
667  END_TIMER("assemble_fluxes_boundary");
668 
669  START_TIMER("assemble_matrix_elem_side");
670  assemble_matrix_element_side<1>();
671  assemble_matrix_element_side<2>();
672  assemble_matrix_element_side<3>();
673  END_TIMER("assemble_matrix_elem_side");
674  END_TIMER("assemble_stiffness");
675 }
676 
677 
678 
679 template<unsigned int dim>
681 {
682  FEValues<3> fe_values(*feo->q<dim>(), *feo->fe<dim>(),
684  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
685  vector<int> dof_indices(ndofs);
686  vector<arma::vec3> velocity(qsize);
687  vector<double> young(qsize), poisson(qsize), csection(qsize);
688  PetscScalar local_matrix[ndofs*ndofs];
689  auto vec = fe_values.vector_view(0);
690 
691  // assemble integral over elements
692  for (auto cell : feo->dh()->own_range())
693  {
694  if (cell.dim() != dim) continue;
695  ElementAccessor<3> elm_acc = cell.elm();
696 
697  fe_values.reinit(elm_acc);
698  cell.get_dof_indices(dof_indices);
699 
700  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
701  data_.young_modulus.value_list(fe_values.point_list(), elm_acc, young);
702  data_.poisson_ratio.value_list(fe_values.point_list(), elm_acc, poisson);
703 
704  // assemble the local stiffness matrix
705  for (unsigned int i=0; i<ndofs; i++)
706  for (unsigned int j=0; j<ndofs; j++)
707  local_matrix[i*ndofs+j] = 0;
708 
709  for (unsigned int k=0; k<qsize; k++)
710  {
711  double mu = lame_mu(young[k], poisson[k]);
712  double lambda = lame_lambda(young[k], poisson[k]);
713  for (unsigned int i=0; i<ndofs; i++)
714  {
715  for (unsigned int j=0; j<ndofs; j++)
716  local_matrix[i*ndofs+j] += csection[k]*(
717  2*mu*arma::dot(vec.sym_grad(j,k), vec.sym_grad(i,k))
718  + lambda*vec.divergence(j,k)*vec.divergence(i,k)
719  )*fe_values.JxW(k);
720  }
721  }
722  ls->mat_set_values(ndofs, dof_indices.data(), ndofs, dof_indices.data(), local_matrix);
723  }
724 }
725 
726 
727 
729 {
730  START_TIMER("assemble_rhs");
731 // balance_->start_source_assembly(subst_idx);
732 // balance_->start_flux_assembly(subst_idx);
733  assemble_sources<1>();
734  assemble_sources<2>();
735  assemble_sources<3>();
736  assemble_rhs_element_side<1>();
737  assemble_rhs_element_side<2>();
738  assemble_rhs_element_side<3>();
739  assemble_boundary_conditions<1>();
740  assemble_boundary_conditions<2>();
741  assemble_boundary_conditions<3>();
742 // balance_->finish_flux_assembly(subst_idx);
743 // balance_->finish_source_assembly(subst_idx);
744  END_TIMER("assemble_rhs");
745 }
746 
747 
748 template<unsigned int dim>
750 {
751  FEValues<3> fe_values(*feo->q<dim>(), *feo->fe<dim>(),
753  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
754  vector<arma::vec3> load(qsize);
755  vector<double> csection(qsize), potential(qsize);
756  vector<int> dof_indices(ndofs);
757  PetscScalar local_rhs[ndofs];
758  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
759  auto vec = fe_values.vector_view(0);
760 
761  // assemble integral over elements
762  for (auto cell : feo->dh()->own_range())
763  {
764  if (cell.dim() != dim) continue;
765  ElementAccessor<3> elm_acc = cell.elm();
766 
767  fe_values.reinit(elm_acc);
768  cell.get_dof_indices(dof_indices);
769 
770  // assemble the local stiffness matrix
771  fill_n(local_rhs, ndofs, 0);
772  local_source_balance_vector.assign(ndofs, 0);
773  local_source_balance_rhs.assign(ndofs, 0);
774 
775  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
776  data_.load.value_list(fe_values.point_list(), elm_acc, load);
777  data_.potential_load.value_list(fe_values.point_list(), elm_acc, potential);
778 
779  // compute sources
780  for (unsigned int k=0; k<qsize; k++)
781  {
782  for (unsigned int i=0; i<ndofs; i++)
783  local_rhs[i] += (
784  arma::dot(load[k], vec.value(i,k))
785  -potential[k]*vec.divergence(i,k)
786  )*csection[k]*fe_values.JxW(k);
787  }
788  ls->rhs_set_values(ndofs, dof_indices.data(), local_rhs);
789 
790 // for (unsigned int i=0; i<ndofs; i++)
791 // {
792 // for (unsigned int k=0; k<qsize; k++)
793 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_values[vec].value(i,k)*fe_values.JxW(k);
794 //
795 // local_source_balance_rhs[i] += local_rhs[i];
796 // }
797 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_vector);
798 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_rhs);
799  }
800 }
801 
802 
803 
805 {
806  double young = data_.young_modulus.value(side->centre(), side->element());
807  double poisson = data_.poisson_ratio.value(side->centre(), side->element());
808  return 1e3*(2*lame_mu(young, poisson) + lame_lambda(young, poisson)) / side->measure();
809 }
810 
811 
812 
813 template<unsigned int dim>
815 {
816  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
818  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
819  vector<int> side_dof_indices(ndofs);
820  PetscScalar local_matrix[ndofs*ndofs];
821  auto vec = fe_values_side.vector_view(0);
822 
823  // assemble boundary integral
824  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
825  {
826  Edge edg = mesh_->edge(feo->dh()->edge_index(iedg));
827  if (edg.n_sides() > 1) continue;
828  // check spatial dimension
829  if (edg.side(0)->dim() != dim-1) continue;
830  // skip edges lying not on the boundary
831  if (! edg.side(0)->is_boundary()) continue;
832 
833  SideIter side = edg.side(0);
834  ElementAccessor<3> cell = side->element();
835  DHCellAccessor dh_cell = feo->dh()->cell_accessor_from_element(cell.idx());
836  DHCellSide dh_side(dh_cell, side->side_idx());
837  dh_cell.get_dof_indices(side_dof_indices);
838  unsigned int bc_type = data_.bc_type.value(side->centre(), side->cond().element_accessor());
839  fe_values_side.reinit(*side);
840 
841  for (unsigned int i=0; i<ndofs; i++)
842  for (unsigned int j=0; j<ndofs; j++)
843  local_matrix[i*ndofs+j] = 0;
844 
845  if (bc_type == EqData::bc_type_displacement)
846  {
847  for (unsigned int k=0; k<qsize; k++)
848  for (unsigned int i=0; i<ndofs; i++)
849  for (unsigned int j=0; j<ndofs; j++)
850  local_matrix[i*ndofs+j] += dirichlet_penalty(side)*arma::dot(vec.value(i,k),vec.value(j,k))*fe_values_side.JxW(k);
851  }
852  else if (bc_type == EqData::bc_type_displacement_normal)
853  {
854  for (unsigned int k=0; k<qsize; k++)
855  for (unsigned int i=0; i<ndofs; i++)
856  for (unsigned int j=0; j<ndofs; j++)
857  local_matrix[i*ndofs+j] += dirichlet_penalty(side)*arma::dot(vec.value(i,k),fe_values_side.normal_vector(k))*arma::dot(vec.value(j,k),fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
858  }
859 
860  ls->mat_set_values(ndofs, side_dof_indices.data(), ndofs, side_dof_indices.data(), local_matrix);
861  }
862 }
863 
864 
866 {
867  arma::mat33 mt = m - m*arma::kron(n,n.t());
868  return mt;
869 }
870 
871 
872 template<unsigned int dim>
874 {
875  if (dim == 1) return;
876  FEValues<3> fe_values_sub(*feo->q<dim-1>(), *feo->fe<dim-1>(),
878  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
880 
881  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
882  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
883  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
884  vector<vector<int> > side_dof_indices(2);
885  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
886  vector<double> frac_sigma(qsize);
887  vector<double> csection_lower(qsize), csection_higher(qsize), young(qsize), poisson(qsize), alpha(qsize);
888  PetscScalar local_matrix[2][2][(ndofs_side)*(ndofs_side)];
889  auto vec_side = fe_values_side.vector_view(0);
890  auto vec_sub = fe_values_sub.vector_view(0);
891 
892  // index 0 = element with lower dimension,
893  // index 1 = side of element with higher dimension
894  side_dof_indices[0].resize(ndofs_sub);
895  side_dof_indices[1].resize(ndofs_side);
896 
897  // assemble integral over sides
898  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
899  {
900  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
901  // skip neighbours of different dimension
902  if (nb->element()->dim() != dim-1) continue;
903 
904  ElementAccessor<3> cell_sub = nb->element();
905  DHCellAccessor dh_cell_sub = feo->dh()->cell_accessor_from_element(cell_sub.idx());
906  dh_cell_sub.get_dof_indices(side_dof_indices[0]);
907  fe_values_sub.reinit(cell_sub);
908 
909  ElementAccessor<3> cell = nb->side()->element();
910  DHCellAccessor dh_cell = feo->dh()->cell_accessor_from_element(cell.idx());
911  DHCellSide dh_side(dh_cell, nb->side()->side_idx());
912  dh_cell.get_dof_indices(side_dof_indices[1]);
913  fe_values_side.reinit(dh_side.side());
914 
915  // Element id's for testing if they belong to local partition.
916  bool own_element_id[2];
917  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
918  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
919 
920  data_.cross_section.value_list(fe_values_sub.point_list(), cell_sub, csection_lower);
921  data_.cross_section.value_list(fe_values_sub.point_list(), cell, csection_higher);
922  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
923  data_.young_modulus.value_list(fe_values_sub.point_list(), cell_sub, young);
924  data_.poisson_ratio.value_list(fe_values_sub.point_list(), cell_sub, poisson);
925 
926  for (unsigned int n=0; n<2; ++n)
927  for (unsigned int i=0; i<ndofs_side; i++)
928  for (unsigned int m=0; m<2; ++m)
929  for (unsigned int j=0; j<ndofs_side; j++)
930  local_matrix[n][m][i*(ndofs_side)+j] = 0;
931 
932  // set transmission conditions
933  for (unsigned int k=0; k<qsize; k++)
934  {
935  arma::vec3 nv = fe_values_side.normal_vector(k);
936  double mu = lame_mu(young[k], poisson[k]);
937  double lambda = lame_lambda(young[k], poisson[k]);
938 
939  for (int n=0; n<2; n++)
940  {
941  if (!own_element_id[n]) continue;
942 
943  for (unsigned int i=0; i<n_dofs[n]; i++)
944  {
945  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
946  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
947  arma::mat33 gvft = (n==0)?vec_sub.grad(i,k):arma::zeros(3,3);
948 
949  for (int m=0; m<2; m++)
950  for (unsigned int j=0; j<n_dofs[m]; j++) {
951  arma::vec3 ui = (m==0)?arma::zeros(3):vec_side.value(j,k);
952  arma::vec3 uf = (m==1)?arma::zeros(3):vec_sub.value(j,k);
953  arma::mat33 guit = (m==1)?mat_t(vec_side.grad(j,k),nv):arma::zeros(3,3);
954  double divuit = (m==1)?arma::trace(guit):0;
955 
956  local_matrix[n][m][i*n_dofs[m] + j] +=
957  frac_sigma[k]*csection_higher[k]*(
958  arma::dot(vf-vi,
959  2*csection_higher[k]/csection_lower[k]*(mu*(uf-ui)+(mu+lambda)*(arma::dot(uf-ui,nv)*nv))
960  + mu*arma::trans(guit)*nv
961  + lambda*divuit*nv
962  )
963  - arma::dot(gvft, mu*arma::kron(nv,ui.t()) + lambda*arma::dot(ui,nv)*arma::eye(3,3))
964  )*fe_values_sub.JxW(k);
965  }
966 
967  }
968  }
969  }
970 
971  for (unsigned int n=0; n<2; ++n)
972  for (unsigned int m=0; m<2; ++m)
973  ls->mat_set_values(n_dofs[n], side_dof_indices[n].data(), n_dofs[m], side_dof_indices[m].data(), local_matrix[n][m]);
974  }
975 }
976 
977 
978 
979 template<unsigned int dim>
981 {
982  if (dim == 1) return;
983  FEValues<3> fe_values_sub(*feo->q<dim-1>(), *feo->fe<dim-1>(),
985  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
987 
988  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
989  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
990  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
991  vector<vector<int> > side_dof_indices(2);
992  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
993  vector<double> frac_sigma(qsize), potential(qsize), csection_higher(qsize);
994  PetscScalar local_rhs[2][ndofs_side];
995  auto vec_side = fe_values_side.vector_view(0);
996  auto vec_sub = fe_values_sub.vector_view(0);
997 
998  // index 0 = element with lower dimension,
999  // index 1 = side of element with higher dimension
1000  side_dof_indices[0].resize(ndofs_sub);
1001  side_dof_indices[1].resize(ndofs_side);
1002 
1003  // assemble integral over sides
1004  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
1005  {
1006  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
1007  // skip neighbours of different dimension
1008  if (nb->element()->dim() != dim-1) continue;
1009 
1010  ElementAccessor<3> cell_sub = nb->element();
1011  feo->dh()->cell_accessor_from_element(cell_sub.idx()).get_dof_indices(side_dof_indices[0]);
1012  fe_values_sub.reinit(cell_sub);
1013 
1014  ElementAccessor<3> cell = nb->side()->element();
1015  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices[1]);
1016  fe_values_side.reinit(*nb->side());
1017 
1018  // Element id's for testing if they belong to local partition.
1019  bool own_element_id[2];
1020  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
1021  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
1022 
1023  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
1024  data_.potential_load.value_list(fe_values_sub.point_list(), cell, potential);
1025  data_.cross_section.value_list(fe_values_sub.point_list(), cell, csection_higher);
1026 
1027  for (unsigned int n=0; n<2; ++n)
1028  for (unsigned int i=0; i<ndofs_side; i++)
1029  local_rhs[n][i] = 0;
1030 
1031  // set transmission conditions
1032  for (unsigned int k=0; k<qsize; k++)
1033  {
1034  arma::vec3 nv = fe_values_side.normal_vector(k);
1035 
1036  for (int n=0; n<2; n++)
1037  {
1038  if (!own_element_id[n]) continue;
1039 
1040  for (unsigned int i=0; i<n_dofs[n]; i++)
1041  {
1042  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
1043  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
1044 
1045  local_rhs[n][i] -= frac_sigma[k]*csection_higher[k]*arma::dot(vf-vi,potential[k]*nv)*fe_values_sub.JxW(k);
1046  }
1047  }
1048  }
1049 
1050  for (unsigned int n=0; n<2; ++n)
1051  ls->rhs_set_values(n_dofs[n], side_dof_indices[n].data(), local_rhs[n]);
1052  }
1053 }
1054 
1055 
1056 
1057 
1058 
1059 template<unsigned int dim>
1061 {
1062  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
1064  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1065  vector<int> side_dof_indices(ndofs);
1066  // unsigned int loc_b=0;
1067  double local_rhs[ndofs];
1068  vector<PetscScalar> local_flux_balance_vector(ndofs);
1069  // PetscScalar local_flux_balance_rhs;
1070  vector<arma::vec3> bc_values(qsize), bc_traction(qsize);
1071  vector<double> csection(qsize), bc_potential(qsize);
1072  auto vec = fe_values_side.vector_view(0);
1073 
1074  for (DHCellAccessor cell : feo->dh()->own_range())
1075  {
1076  ElementAccessor<3> elm = cell.elm();
1077  if (elm->boundary_idx_ == nullptr) continue;
1078 
1079  for (unsigned int si=0; si<elm->n_sides(); ++si)
1080  {
1081  Edge edg = elm.side(si)->edge();
1082  Side side = *cell.elm().side(si);
1083 
1084  if (edg.n_sides() > 1) continue;
1085  // skip edges lying not on the boundary
1086  if ( ! side.is_boundary()) continue;
1087 
1088  if (side.dim() != dim-1)
1089  {
1090  // if (edg.side(0)->cond() != nullptr) ++loc_b;
1091  continue;
1092  }
1093 
1094  ElementAccessor<3> bc_cell = side.cond().element_accessor();
1095  unsigned int bc_type = data_.bc_type.value(side.centre(), bc_cell);
1096 
1097  fe_values_side.reinit(side);
1098 
1099  data_.cross_section.value_list(fe_values_side.point_list(), elm, csection);
1100  // The b.c. data are fetched for all possible b.c. types since we allow
1101  // different bc_type for each substance.
1102  data_.bc_displacement.value_list(fe_values_side.point_list(), bc_cell, bc_values);
1103  data_.bc_traction.value_list(fe_values_side.point_list(), bc_cell, bc_traction);
1104  data_.potential_load.value_list(fe_values_side.point_list(), elm, bc_potential);
1105 
1106  cell.get_dof_indices(side_dof_indices);
1107 
1108  fill_n(local_rhs, ndofs, 0);
1109  local_flux_balance_vector.assign(ndofs, 0);
1110  // local_flux_balance_rhs = 0;
1111 
1112  if (bc_type == EqData::bc_type_displacement)
1113  {
1114  for (unsigned int k=0; k<qsize; k++)
1115  for (unsigned int i=0; i<ndofs; i++)
1116  local_rhs[i] += dirichlet_penalty(side)*arma::dot(bc_values[k],vec.value(i,k))*fe_values_side.JxW(k);
1117  }
1118  else if (bc_type == EqData::bc_type_displacement_normal)
1119  {
1120  for (unsigned int k=0; k<qsize; k++)
1121  for (unsigned int i=0; i<ndofs; i++)
1122  local_rhs[i] += dirichlet_penalty(side)*arma::dot(bc_values[k],fe_values_side.normal_vector(k))*arma::dot(vec.value(i,k),fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1123  }
1124  else if (bc_type == EqData::bc_type_traction)
1125  {
1126  for (unsigned int k=0; k<qsize; k++)
1127  {
1128  for (unsigned int i=0; i<ndofs; i++)
1129  local_rhs[i] += csection[k]*arma::dot(vec.value(i,k),bc_traction[k] + bc_potential[k]*fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1130  }
1131  }
1132  ls->rhs_set_values(ndofs, side_dof_indices.data(), local_rhs);
1133 
1134 
1135 
1136 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
1137 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
1138  // ++loc_b;
1139  }
1140  }
1141 }
1142 
1143 
1144 
1145 
TimeGovernor & time()
Definition: equation.hh:149
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:111
std::shared_ptr< DiscreteSpace > ds_
Definition: elasticity.hh:72
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:112
Field< 3, FieldValue< 3 >::TensorFixed > output_stress
Definition: elasticity.hh:126
Shape function values.
Definition: update_flags.hh:87
unsigned int n_sides() const
Returns number of sides aligned with the edge.
Definition: accessors.hh:300
Field< 3, FieldValue< 3 >::Scalar > young_modulus
Definition: elasticity.hh:115
FieldSet * eq_data_
Definition: equation.hh:227
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:50
static auto subdomain(Mesh &mesh) -> IndexField
static constexpr const char * name()
Definition: elasticity.hh:102
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
Definition: elasticity.hh:133
Field< 3, FieldValue< 3 >::Scalar > output_von_mises_stress
Definition: elasticity.hh:127
Transformed quadrature weight for cell sides.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:328
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:254
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
ArmaVec< double, N > vec
Definition: armor.hh:861
unsigned int * boundary_idx_
Definition: elements.h:80
void assemble_matrix_element_side()
Assembles the fluxes between elements of different dimensions depending on displacement.
Definition: elasticity.cc:873
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void preallocate()
Definition: elasticity.cc:434
void inc()
Iterates to next local element.
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: elasticity.hh:122
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:415
Field< 3, FieldValue< 3 >::Scalar > output_cross_section
Definition: elasticity.hh:128
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:37
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Edge edge() const
Returns pointer to the edge connected to the side.
void set_from_input(const Input::Record in_rec) override
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:255
virtual void start_add_assembly()
Definition: linsys.hh:341
void output(TimeStep step)
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:264
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
void next_time()
Proceed to the next time according to current estimated time step.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: elasticity.cc:43
void assemble_rhs()
Assembles the right hand side (forces, boundary conditions, tractions).
Definition: elasticity.cc:728
arma::vec3 centre() const
Centre of side.
Definition: accessors.cc:116
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
static std::shared_ptr< OutputTime > create_output_stream(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:184
LinSys * ls
Linear algebra system for the transport equation.
Definition: elasticity.hh:300
Definition: mesh.h:78
Fields computed from the mesh data.
virtual void start_allocation()
Definition: linsys.hh:333
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
Definition: elasticity.hh:132
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
virtual void finish_assembly()=0
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:120
Class FESystem for compound finite elements.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
~Elasticity()
Destructor.
Definition: elasticity.cc:357
SideIter side(const unsigned int loc_index)
double lame_mu(double young, double poisson)
Definition: elasticity.cc:128
const RegionDB & region_db() const
Definition: mesh.h:143
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:223
const TimeStep & step(int index=-1) const
void update_solution() override
Computes the solution in one time instant.
Definition: elasticity.cc:450
std::shared_ptr< DOFHandlerMultiDim > dh()
Definition: elasticity.cc:119
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
Definition: elasticity.cc:865
Basic time management functionality for unsteady (and steady) solvers (class Equation).
void assemble_rhs_element_side()
Assemble fluxes between different dimensions that are independent of displacement.
Definition: elasticity.cc:980
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:548
void compute_output_fields()
Definition: elasticity.cc:542
Input::Record input_rec
Record with input specification.
Definition: elasticity.hh:311
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
Definition: elasticity.hh:121
Class for declaration of inputs sequences.
Definition: type_base.hh:339
Vec rhs
Vector of right hand side.
Definition: elasticity.hh:294
void view(const char *name="") const
Symmetric Gauss-Legendre quadrature formulae on simplices.
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:35
FEM for linear elasticity.
bool allocation_done
Indicates whether matrices have been preallocated.
Definition: elasticity.hh:323
unsigned int dim() const
Definition: elements.h:121
ElementAccessor< 3 > element()
Definition: neighbours.h:161
Transformed quadrature points.
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:113
bool is_boundary() const
Returns true for side on the boundary.
void solve_linear_system()
Solve without updating time step and without output.
Definition: elasticity.cc:473
void assemble_sources()
Assembles the right hand side vector due to volume sources.
Definition: elasticity.cc:749
Definitions of basic Lagrangean finite elements with polynomial shape functions.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
SideIter side()
Definition: neighbours.h:145
void update_output_fields()
Definition: elasticity.cc:370
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:132
MixedPtr< FiniteElement > fe_
Finite elements for the solution of the mechanics equation.
Definition: elasticity.hh:68
static auto region_id(Mesh &mesh) -> IndexField
#define xprintf(...)
Definition: system.hh:93
#define START_TIMER(tag)
Starts a timer with specified tag.
int converged_reason
Definition: linsys.hh:109
Mesh * mesh_
Definition: equation.hh:218
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:434
Field< 3, FieldValue< 3 >::Scalar > output_divergence
Definition: elasticity.hh:129
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
void zero_time_step() override
Initialize solution in the zero time.
Definition: elasticity.cc:398
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
Shape function gradients.
Definition: update_flags.hh:95
void output_data()
Postprocesses the solution and writes to output file.
Definition: elasticity.cc:521
Mat stiffness_matrix
The stiffness matrix.
Definition: elasticity.hh:297
std::shared_ptr< FiniteElement< dim > > fe()
Quadrature * q()
Definition: elasticity.hh:58
virtual void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list) const
Definition: field.hh:448
double measure() const
Calculate metrics of the side.
Definition: accessors.cc:27
virtual const Vec * get_rhs()
Definition: linsys.hh:203
Normal vectors.
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:114
void set_solution(Vec sol_vec)
Definition: linsys.hh:290
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Definition: elasticity.hh:117
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
Definition: elasticity.cc:814
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:221
vector< Neighbour > vb_neighbours_
Definition: mesh.h:267
void assemble_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space ...
Definition: elasticity.cc:1060
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:193
Field< 3, FieldValue< 3 >::Scalar > poisson_ratio
Definition: elasticity.hh:116
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
const Selection & close() const
Close the Selection, no more values can be added.
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
Definition: elasticity.hh:77
void set_components(const std::vector< string > &names)
Definition: field_set.hh:180
EquationOutput output_fields
Definition: elasticity.hh:137
Definition: system.hh:65
static const int registrar
Registrar of class to factory.
Definition: elasticity.hh:208
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:227
void next_time()
Pass to next time and update equation data.
Definition: elasticity.cc:464
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
void initialize() override
Definition: elasticity.cc:299
EqData data_
Field data for model parameters.
Definition: elasticity.hh:275
double lame_lambda(double young, double poisson)
Definition: elasticity.cc:134
static const Input::Type::Selection & get_bc_type_selection()
Definition: elasticity.cc:144
MixedPtr< FESystem > mixed_fe_system(MixedPtr< FiniteElement > fe, Args &&...args)
Definition: fe_system.hh:171
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:258
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
static string default_output_field()
Definition: elasticity.hh:104
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
Definition: elasticity.hh:131
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
Edge edge(uint edge_idx) const
Definition: mesh.cc:254
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
Definition: elasticity.hh:75
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:186
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
Definition: elasticity.cc:654
std::shared_ptr< DOFHandlerMultiDim > dh_scalar()
Definition: elasticity.cc:120
Record type proxy class.
Definition: type_record.hh:182
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:200
Class for representation SI units of Fields.
Definition: unit_si.hh:40
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
Definition: elasticity.hh:135
Boundary cond() const
std::shared_ptr< DOFHandlerMultiDim > dh_tensor()
Definition: elasticity.cc:121
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
Definition: elasticity.cc:680
virtual const Mat * get_matrix()
Definition: linsys.hh:187
void calculate_cumulative_balance()
Definition: elasticity.cc:639
Elasticity(Mesh &init_mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
Constructor.
Definition: elasticity.cc:263
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
Definition: elasticity.hh:134
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.hh:213
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:264
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
virtual double compute_residual()=0
Field< 3, FieldValue< 3 >::VectorFixed > output_field
Definition: elasticity.hh:125
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
Definition: elasticity.hh:76
SideIter side(const unsigned int i) const
Gets side iterator of the i -th side.
bool changed() const
Definition: field_set.cc:165
Template for classes storing finite set of named values.
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: elasticity.hh:123
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
FieldCommon & output_type(OutputTime::DiscreteSpace rt)
ElementAccessor< 3 > element_accessor()
std::shared_ptr< OutputTime > output_stream_
Definition: elasticity.hh:308
TimeGovernor * time_
Definition: equation.hh:219
Mechanics::FEObjects * feo
Finite element objects.
Definition: elasticity.hh:284
double dirichlet_penalty(SideIter side)
Penalty to enforce boundary value in weak sense.
Definition: elasticity.cc:804
Transformed quadrature weights.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:234