Flow123d  JS_before_hm-919-g5f1bbbf
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 heat transfer model.")
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.")
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.")
177  .units( UnitSI().Pa() )
178  .input_default("0.0")
179  .flags_add(in_rhs);
180 
181  *this+=load
182  .name("load")
183  .description("Prescribed 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 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 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 diffusive transfer through fractures (for each substance).")
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  .units( UnitSI().m() )
231  .flags(equation_result);
232 
233  *this += output_stress
234  .name("stress")
235  .units( UnitSI().Pa() )
236  .flags(equation_result);
237 
238  *this += output_von_mises_stress
239  .name("von_mises_stress")
240  .units( UnitSI().Pa() )
241  .flags(equation_result);
242 
243  *this += output_cross_section
244  .name("cross_section_updated")
245  .units( UnitSI().m() )
246  .flags(equation_result);
247 
248  *this += output_divergence
249  .name("displacement_divergence")
250  .units( UnitSI().dimensionless() )
251  .flags(equation_result);
252 
253  // add all input fields to the output list
254  output_fields += *this;
255 
256 }
257 
258 Elasticity::Elasticity(Mesh & init_mesh, const Input::Record in_rec, TimeGovernor *tm)
259  : EquationBase(init_mesh, in_rec),
260  input_rec(in_rec),
261  allocation_done(false)
262 {
263  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
264  // due to constexpr optimization.
266 
267  this->eq_data_ = &data_;
268 
269  auto time_rec = in_rec.val<Input::Record>("time");
270  if (tm == nullptr)
271  {
272  time_ = new TimeGovernor(time_rec);
273  }
274  else
275  {
276  TimeGovernor time_from_rec(time_rec);
277  ASSERT( time_from_rec.is_default() ).error("Duplicate key 'time', time in elasticity is already initialized from parent class!");
278  time_ = tm;
279  }
280 
281 
282  // Set up physical parameters.
283  data_.set_mesh(init_mesh);
286 
287  // create finite element structures and distribute DOFs
288  feo = new Mechanics::FEObjects(mesh_, 1);
289  DebugOut().fmt("Mechanics: solution size {}\n", feo->dh()->n_global_dofs());
290 
291 }
292 
293 
295 {
296  output_stream_ = OutputTime::create_output_stream("mechanics", input_rec.val<Input::Record>("output_stream"), time().get_unit_string());
297 
298  data_.set_components({"displacement"});
299  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), time() );
300 
301 // balance_ = std::make_shared<Balance>("mechanics", mesh_);
302 // balance_->init_from_input(input_rec.val<Input::Record>("balance"), *time_);
303  // initialization of balance object
304 // data_.balance_idx_ = {
305 // balance_->add_quantity("force_x"),
306 // balance_->add_quantity("force_y"),
307 // balance_->add_quantity("force_z")
308 // };
309 // balance_->units(UnitSI().kg().m().s(-2));
310 
311  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
312  data_.output_field_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(feo->dh());
314 
315  // setup output stress
316  data_.output_stress_ptr = create_field_fe<3, FieldValue<3>::TensorFixed>(feo->dh_tensor());
318 
319  // setup output von Mises stress
320  data_.output_von_mises_stress_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
322 
323  // setup output cross-section
324  data_.output_cross_section_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
326 
327  // setup output divergence
328  data_.output_div_ptr = create_field_fe<3, FieldValue<3>::Scalar>(feo->dh_scalar());
330 
332 
333  // set time marks for writing the output
335 
336  // equation default PETSc solver options
337  std::string petsc_default_opts;
338  petsc_default_opts = "-ksp_type cg -pc_type hypre -pc_hypre_type boomeramg";
339 
340  // allocate matrix and vector structures
341  ls = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
343  ls->set_solution(data_.output_field_ptr->get_data_vec().petsc_vec());
344 
345  // initialization of balance object
346 // balance_->allocate(feo->dh()->distr()->lsize(),
347 // max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
348 
349 }
350 
351 
353 {
354 // delete time_;
355 
356  MatDestroy(&stiffness_matrix);
357  VecDestroy(&rhs);
358  delete ls;
359  delete feo;
360 
361 }
362 
363 
364 
366 {
367  // update ghost values of solution vector
368  data_.output_field_ptr->get_data_vec().local_to_ghost_begin();
369  data_.output_field_ptr->get_data_vec().local_to_ghost_end();
370 
371  // compute new output fields depending on solution (stress, divergence etc.)
372  data_.output_stress_ptr->get_data_vec().zero_entries();
373  data_.output_cross_section_ptr->get_data_vec().zero_entries();
374  data_.output_div_ptr->get_data_vec().zero_entries();
375  compute_output_fields<1>();
376  compute_output_fields<2>();
377  compute_output_fields<3>();
378 
379  // update ghost values of computed fields
380  data_.output_stress_ptr->get_data_vec().local_to_ghost_begin();
381  data_.output_stress_ptr->get_data_vec().local_to_ghost_end();
382  data_.output_von_mises_stress_ptr->get_data_vec().local_to_ghost_begin();
383  data_.output_von_mises_stress_ptr->get_data_vec().local_to_ghost_end();
384  data_.output_cross_section_ptr->get_data_vec().local_to_ghost_begin();
385  data_.output_cross_section_ptr->get_data_vec().local_to_ghost_end();
386  data_.output_div_ptr->get_data_vec().local_to_ghost_begin();
387  data_.output_div_ptr->get_data_vec().local_to_ghost_end();
388 }
389 
390 
391 
392 
394 {
398  std::stringstream ss; // print warning message with table of uninitialized fields
399  if ( FieldCommon::print_message_table(ss, "mechanics") ) {
400  WarningOut() << ss.str();
401  }
402 
403  ( (LinSys_PETSC *)ls )->set_initial_guess_nonzero();
404 
405  // check first time assembly - needs preallocation
407 
408 
409  // after preallocation we assemble the matrices and vectors required for balance of forces
410 // for (auto subst_idx : data_.balance_idx_)
411 // balance_->calculate_instant(subst_idx, ls->get_solution());
412 
413 // update_solution();
415  MatSetOption(*ls->get_matrix(), MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
416  ls->mat_zero_entries();
417  ls->rhs_zero_entries();
419  assemble_rhs();
420  ls->finish_assembly();
421  ls->solve();
422  output_data();
423 }
424 
425 
426 
428 {
429  // preallocate system matrix
430  ls->start_allocation();
431  stiffness_matrix = NULL;
432  rhs = NULL;
433 
435  assemble_rhs();
436 
437  allocation_done = true;
438 }
439 
440 
441 
442 
444 {
445  START_TIMER("DG-ONE STEP");
446 
447  next_time();
449 
451 
452  output_data();
453 
454  END_TIMER("DG-ONE STEP");
455 }
456 
458 {
459  time_->next_time();
460  time_->view("MECH");
461 
462 }
463 
464 
465 
467 {
468  START_TIMER("data reinit");
470  END_TIMER("data reinit");
471 
472  // assemble stiffness matrix
473  if (stiffness_matrix == NULL
475  {
476  DebugOut() << "Mechanics: Assembling matrix.\n";
478  ls->mat_zero_entries();
480  ls->finish_assembly();
481 
482  if (stiffness_matrix == NULL)
483  MatConvert(*( ls->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix);
484  else
485  MatCopy(*( ls->get_matrix() ), stiffness_matrix, DIFFERENT_NONZERO_PATTERN);
486  }
487 
488  // assemble right hand side (due to sources and boundary conditions)
489  if (rhs == NULL
490  || data_.subset(FieldFlag::in_rhs).changed())
491  {
492  DebugOut() << "Mechanics: Assembling right hand side.\n";
494  ls->rhs_zero_entries();
495  assemble_rhs();
496  ls->finish_assembly();
497 
498  if (rhs == nullptr) VecDuplicate(*( ls->get_rhs() ), &rhs);
499  VecCopy(*( ls->get_rhs() ), rhs);
500  }
501 
502  START_TIMER("solve");
503  ls->solve();
504  END_TIMER("solve");
505 }
506 
507 
508 
509 
510 
511 
513 {
514  START_TIMER("MECH-OUTPUT");
515 
516  // gather the solution from all processors
517  data_.output_fields.set_time( this->time().step(), LimitSide::left);
518  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
520  data_.output_fields.output(this->time().step());
521  output_stream_->write_time_frame();
522 
523 // START_TIMER("MECH-balance");
524 // balance_->calculate_instant(subst_idx, ls->get_solution());
525 // balance_->output();
526 // END_TIMER("MECH-balance");
527 
528  END_TIMER("MECH-OUTPUT");
529 }
530 
531 
532 template<unsigned int dim>
534 {
535  QGauss q(dim, 0), q_sub(dim-1, 0);
536  FEValues<3> fv(q, *feo->fe<dim>(),
538  FEValues<3> fsv(q_sub, *feo->fe<dim>(),
540  const unsigned int ndofs = feo->fe<dim>()->n_dofs();
541  auto vec = fv.vector_view(0);
542  auto vec_side = fsv.vector_view(0);
543  VectorMPI output_vec = data_.output_field_ptr->get_data_vec();
544  VectorMPI output_stress_vec = data_.output_stress_ptr->get_data_vec();
545  VectorMPI output_von_mises_stress_vec = data_.output_von_mises_stress_ptr->get_data_vec();
546  VectorMPI output_cross_sec_vec = data_.output_cross_section_ptr->get_data_vec();
547  VectorMPI output_div_vec = data_.output_div_ptr->get_data_vec();
548 
549  DHCellAccessor cell_tensor = *feo->dh_tensor()->own_range().begin();
550  DHCellAccessor cell_scalar = *feo->dh_scalar()->own_range().begin();
551  for (auto cell : feo->dh()->own_range())
552  {
553  if (cell.dim() == dim)
554  {
555  auto elm = cell.elm();
556 
557  double poisson = data_.poisson_ratio.value(elm.centre(), elm);
558  double young = data_.young_modulus.value(elm.centre(), elm);
559  double mu = lame_mu(young, poisson);
560  double lambda = lame_lambda(young, poisson);
561 
562  fv.reinit(elm);
563  LocDofVec dof_indices = cell.get_loc_dof_indices();
564  LocDofVec dof_indices_scalar = cell_scalar.get_loc_dof_indices();
565  LocDofVec dof_indices_tensor = cell_tensor.get_loc_dof_indices();
566 
567  arma::mat33 stress = arma::zeros(3,3);
568  double div = 0;
569  for (unsigned int i=0; i<ndofs; i++)
570  {
571  stress += (2*mu*vec.sym_grad(i,0) + lambda*vec.divergence(i,0)*arma::eye(3,3))*output_vec[dof_indices[i]];
572  div += vec.divergence(i,0)*output_vec[dof_indices[i]];
573  }
574 
575  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
576  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
577  output_div_vec[dof_indices_scalar[0]] += div;
578 
579  for (unsigned int i=0; i<3; i++)
580  for (unsigned int j=0; j<3; j++)
581  output_stress_vec[dof_indices_tensor[i*3+j]] += stress(i,j);
582  output_von_mises_stress_vec[dof_indices_scalar[0]] = von_mises_stress;
583 
584  output_cross_sec_vec[dof_indices_scalar[0]] += data_.cross_section.value(fv.point(0), elm);
585  }
586  else if (cell.dim() == dim-1)
587  {
588  auto elm = cell.elm();
589  double normal_displacement = 0;
590  double csection = data_.cross_section.value(fsv.point(0), elm);
591  arma::mat33 normal_stress = arma::zeros(3,3);
592 
593  double poisson = data_.poisson_ratio.value(elm.centre(), elm);
594  double young = data_.young_modulus.value(elm.centre(), elm);
595  double mu = lame_mu(young, poisson);
596  double lambda = lame_lambda(young, poisson);
597 
598  for (unsigned int inb=0; inb<elm->n_neighs_vb(); inb++)
599  {
600  auto side = elm->neigh_vb[inb]->side();
601  auto cell_side = side->element();
602  fsv.reinit(*side);
603  LocDofVec side_dof_indices =
604  feo->dh()->cell_accessor_from_element(cell_side.idx()).get_loc_dof_indices();
605 
606  for (unsigned int i=0; i<ndofs; i++)
607  {
608  normal_displacement -= arma::dot(vec_side.value(i,0)*output_vec[side_dof_indices[i]], fsv.normal_vector(0));
609  arma::mat33 grad = -arma::kron(vec_side.value(i,0)*output_vec[side_dof_indices[i]], fsv.normal_vector(0).t()) / csection;
610  normal_stress += mu*(grad+grad.t()) + lambda*arma::trace(grad)*arma::eye(3,3);
611  }
612  }
613  LocDofVec dof_indices_scalar = cell_scalar.get_loc_dof_indices();
614  LocDofVec dof_indices_tensor = cell_tensor.get_loc_dof_indices();
615  for (unsigned int i=0; i<3; i++)
616  for (unsigned int j=0; j<3; j++)
617  output_stress_vec[dof_indices_tensor[i*3+j]] += normal_stress(i,j);
618  output_cross_sec_vec[dof_indices_scalar[0]] += normal_displacement;
619  output_div_vec[dof_indices_scalar[0]] += normal_displacement / csection;
620  }
621  cell_scalar.inc();
622  cell_tensor.inc();
623  }
624 }
625 
626 
627 
628 
629 
631 {
632 // if (balance_->cumulative())
633 // {
634 // balance_->calculate_cumulative(subst_idx, ls->get_solution());
635 // }
636 }
637 
638 
639 
640 
641 
642 
643 
644 
646 {
647  START_TIMER("assemble_stiffness");
648  START_TIMER("assemble_volume_integrals");
649  assemble_volume_integrals<1>();
650  assemble_volume_integrals<2>();
651  assemble_volume_integrals<3>();
652  END_TIMER("assemble_volume_integrals");
653 
654  START_TIMER("assemble_fluxes_boundary");
655  assemble_fluxes_boundary<1>();
656  assemble_fluxes_boundary<2>();
657  assemble_fluxes_boundary<3>();
658  END_TIMER("assemble_fluxes_boundary");
659 
660  START_TIMER("assemble_matrix_elem_side");
661  assemble_matrix_element_side<1>();
662  assemble_matrix_element_side<2>();
663  assemble_matrix_element_side<3>();
664  END_TIMER("assemble_matrix_elem_side");
665  END_TIMER("assemble_stiffness");
666 }
667 
668 
669 
670 template<unsigned int dim>
672 {
673  FEValues<3> fe_values(*feo->q<dim>(), *feo->fe<dim>(),
675  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
676  vector<int> dof_indices(ndofs);
677  vector<arma::vec3> velocity(qsize);
678  vector<double> young(qsize), poisson(qsize), csection(qsize);
679  PetscScalar local_matrix[ndofs*ndofs];
680  auto vec = fe_values.vector_view(0);
681 
682  // assemble integral over elements
683  for (auto cell : feo->dh()->own_range())
684  {
685  if (cell.dim() != dim) continue;
686  ElementAccessor<3> elm_acc = cell.elm();
687 
688  fe_values.reinit(elm_acc);
689  cell.get_dof_indices(dof_indices);
690 
691  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
692  data_.young_modulus.value_list(fe_values.point_list(), elm_acc, young);
693  data_.poisson_ratio.value_list(fe_values.point_list(), elm_acc, poisson);
694 
695  // assemble the local stiffness matrix
696  for (unsigned int i=0; i<ndofs; i++)
697  for (unsigned int j=0; j<ndofs; j++)
698  local_matrix[i*ndofs+j] = 0;
699 
700  for (unsigned int k=0; k<qsize; k++)
701  {
702  double mu = lame_mu(young[k], poisson[k]);
703  double lambda = lame_lambda(young[k], poisson[k]);
704  for (unsigned int i=0; i<ndofs; i++)
705  {
706  for (unsigned int j=0; j<ndofs; j++)
707  local_matrix[i*ndofs+j] += csection[k]*(
708  2*mu*arma::dot(vec.sym_grad(j,k), vec.sym_grad(i,k))
709  + lambda*vec.divergence(j,k)*vec.divergence(i,k)
710  )*fe_values.JxW(k);
711  }
712  }
713  ls->mat_set_values(ndofs, dof_indices.data(), ndofs, dof_indices.data(), local_matrix);
714  }
715 }
716 
717 
718 
720 {
721  START_TIMER("assemble_rhs");
722 // balance_->start_source_assembly(subst_idx);
723 // balance_->start_flux_assembly(subst_idx);
724  assemble_sources<1>();
725  assemble_sources<2>();
726  assemble_sources<3>();
727  assemble_rhs_element_side<1>();
728  assemble_rhs_element_side<2>();
729  assemble_rhs_element_side<3>();
730  assemble_boundary_conditions<1>();
731  assemble_boundary_conditions<2>();
732  assemble_boundary_conditions<3>();
733 // balance_->finish_flux_assembly(subst_idx);
734 // balance_->finish_source_assembly(subst_idx);
735  END_TIMER("assemble_rhs");
736 }
737 
738 
739 template<unsigned int dim>
741 {
742  FEValues<3> fe_values(*feo->q<dim>(), *feo->fe<dim>(),
744  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
745  vector<arma::vec3> load(qsize);
746  vector<double> csection(qsize), potential(qsize);
747  vector<int> dof_indices(ndofs);
748  PetscScalar local_rhs[ndofs];
749  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
750  auto vec = fe_values.vector_view(0);
751 
752  // assemble integral over elements
753  for (auto cell : feo->dh()->own_range())
754  {
755  if (cell.dim() != dim) continue;
756  ElementAccessor<3> elm_acc = cell.elm();
757 
758  fe_values.reinit(elm_acc);
759  cell.get_dof_indices(dof_indices);
760 
761  // assemble the local stiffness matrix
762  fill_n(local_rhs, ndofs, 0);
763  local_source_balance_vector.assign(ndofs, 0);
764  local_source_balance_rhs.assign(ndofs, 0);
765 
766  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
767  data_.load.value_list(fe_values.point_list(), elm_acc, load);
768  data_.potential_load.value_list(fe_values.point_list(), elm_acc, potential);
769 
770  // compute sources
771  for (unsigned int k=0; k<qsize; k++)
772  {
773  for (unsigned int i=0; i<ndofs; i++)
774  local_rhs[i] += (
775  arma::dot(load[k], vec.value(i,k))
776  -potential[k]*vec.divergence(i,k)
777  )*csection[k]*fe_values.JxW(k);
778  }
779  ls->rhs_set_values(ndofs, dof_indices.data(), local_rhs);
780 
781 // for (unsigned int i=0; i<ndofs; i++)
782 // {
783 // for (unsigned int k=0; k<qsize; k++)
784 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_values[vec].value(i,k)*fe_values.JxW(k);
785 //
786 // local_source_balance_rhs[i] += local_rhs[i];
787 // }
788 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_vector);
789 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_rhs);
790  }
791 }
792 
793 
794 
796 {
797  double young = data_.young_modulus.value(side->centre(), side->element());
798  double poisson = data_.poisson_ratio.value(side->centre(), side->element());
799  return 1e3*(2*lame_mu(young, poisson) + lame_lambda(young, poisson)) / side->measure();
800 }
801 
802 
803 
804 template<unsigned int dim>
806 {
807  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
809  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
810  vector<int> side_dof_indices(ndofs);
811  PetscScalar local_matrix[ndofs*ndofs];
812  auto vec = fe_values_side.vector_view(0);
813 
814  // assemble boundary integral
815  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
816  {
817  Edge edg = mesh_->edge(feo->dh()->edge_index(iedg));
818  if (edg.n_sides() > 1) continue;
819  // check spatial dimension
820  if (edg.side(0)->dim() != dim-1) continue;
821  // skip edges lying not on the boundary
822  if (! edg.side(0)->is_boundary()) continue;
823 
824  SideIter side = edg.side(0);
825  ElementAccessor<3> cell = side->element();
826  DHCellAccessor dh_cell = feo->dh()->cell_accessor_from_element(cell.idx());
827  DHCellSide dh_side(dh_cell, side->side_idx());
828  dh_cell.get_dof_indices(side_dof_indices);
829  unsigned int bc_type = data_.bc_type.value(side->centre(), side->cond().element_accessor());
830  fe_values_side.reinit(*side);
831 
832  for (unsigned int i=0; i<ndofs; i++)
833  for (unsigned int j=0; j<ndofs; j++)
834  local_matrix[i*ndofs+j] = 0;
835 
836  if (bc_type == EqData::bc_type_displacement)
837  {
838  for (unsigned int k=0; k<qsize; k++)
839  for (unsigned int i=0; i<ndofs; i++)
840  for (unsigned int j=0; j<ndofs; j++)
841  local_matrix[i*ndofs+j] += dirichlet_penalty(side)*arma::dot(vec.value(i,k),vec.value(j,k))*fe_values_side.JxW(k);
842  }
843  else if (bc_type == EqData::bc_type_displacement_normal)
844  {
845  for (unsigned int k=0; k<qsize; k++)
846  for (unsigned int i=0; i<ndofs; i++)
847  for (unsigned int j=0; j<ndofs; j++)
848  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);
849  }
850 
851  ls->mat_set_values(ndofs, side_dof_indices.data(), ndofs, side_dof_indices.data(), local_matrix);
852  }
853 }
854 
855 
857 {
858  arma::mat33 mt = m - m*arma::kron(n,n.t());
859  return mt;
860 }
861 
862 
863 template<unsigned int dim>
865 {
866  if (dim == 1) return;
867  FEValues<3> fe_values_sub(*feo->q<dim-1>(), *feo->fe<dim-1>(),
869  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
871 
872  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
873  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
874  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
875  vector<vector<int> > side_dof_indices(2);
876  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
877  vector<double> frac_sigma(qsize);
878  vector<double> csection_lower(qsize), csection_higher(qsize), young(qsize), poisson(qsize), alpha(qsize);
879  PetscScalar local_matrix[2][2][(ndofs_side)*(ndofs_side)];
880  auto vec_side = fe_values_side.vector_view(0);
881  auto vec_sub = fe_values_sub.vector_view(0);
882 
883  // index 0 = element with lower dimension,
884  // index 1 = side of element with higher dimension
885  side_dof_indices[0].resize(ndofs_sub);
886  side_dof_indices[1].resize(ndofs_side);
887 
888  // assemble integral over sides
889  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
890  {
891  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
892  // skip neighbours of different dimension
893  if (nb->element()->dim() != dim-1) continue;
894 
895  ElementAccessor<3> cell_sub = nb->element();
896  DHCellAccessor dh_cell_sub = feo->dh()->cell_accessor_from_element(cell_sub.idx());
897  dh_cell_sub.get_dof_indices(side_dof_indices[0]);
898  fe_values_sub.reinit(cell_sub);
899 
900  ElementAccessor<3> cell = nb->side()->element();
901  DHCellAccessor dh_cell = feo->dh()->cell_accessor_from_element(cell.idx());
902  DHCellSide dh_side(dh_cell, nb->side()->side_idx());
903  dh_cell.get_dof_indices(side_dof_indices[1]);
904  fe_values_side.reinit(dh_side.side());
905 
906  // Element id's for testing if they belong to local partition.
907  bool own_element_id[2];
908  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
909  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
910 
911  data_.cross_section.value_list(fe_values_sub.point_list(), cell_sub, csection_lower);
912  data_.cross_section.value_list(fe_values_sub.point_list(), cell, csection_higher);
913  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
914  data_.young_modulus.value_list(fe_values_sub.point_list(), cell_sub, young);
915  data_.poisson_ratio.value_list(fe_values_sub.point_list(), cell_sub, poisson);
916 
917  for (unsigned int n=0; n<2; ++n)
918  for (unsigned int i=0; i<ndofs_side; i++)
919  for (unsigned int m=0; m<2; ++m)
920  for (unsigned int j=0; j<ndofs_side; j++)
921  local_matrix[n][m][i*(ndofs_side)+j] = 0;
922 
923  // set transmission conditions
924  for (unsigned int k=0; k<qsize; k++)
925  {
926  arma::vec3 nv = fe_values_side.normal_vector(k);
927  double mu = lame_mu(young[k], poisson[k]);
928  double lambda = lame_lambda(young[k], poisson[k]);
929 
930  for (int n=0; n<2; n++)
931  {
932  if (!own_element_id[n]) continue;
933 
934  for (unsigned int i=0; i<n_dofs[n]; i++)
935  {
936  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
937  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
938  arma::mat33 gvft = (n==0)?vec_sub.grad(i,k):arma::zeros(3,3);
939 
940  for (int m=0; m<2; m++)
941  for (unsigned int j=0; j<n_dofs[m]; j++) {
942  arma::vec3 ui = (m==0)?arma::zeros(3):vec_side.value(j,k);
943  arma::vec3 uf = (m==1)?arma::zeros(3):vec_sub.value(j,k);
944  arma::mat33 guit = (m==1)?mat_t(vec_side.grad(j,k),nv):arma::zeros(3,3);
945  double divuit = (m==1)?arma::trace(guit):0;
946 
947  local_matrix[n][m][i*n_dofs[m] + j] +=
948  frac_sigma[k]*(
949  arma::dot(vf-vi,
950  2/csection_lower[k]*(mu*(uf-ui)+(mu+lambda)*(arma::dot(uf-ui,nv)*nv))
951  + mu*arma::trans(guit)*nv
952  + lambda*divuit*nv
953  )
954  - arma::dot(gvft, mu*arma::kron(nv,ui.t()) + lambda*arma::dot(ui,nv)*arma::eye(3,3))
955  )*fe_values_sub.JxW(k);
956  }
957 
958  }
959  }
960  }
961 
962  for (unsigned int n=0; n<2; ++n)
963  for (unsigned int m=0; m<2; ++m)
964  ls->mat_set_values(n_dofs[n], side_dof_indices[n].data(), n_dofs[m], side_dof_indices[m].data(), local_matrix[n][m]);
965  }
966 }
967 
968 
969 
970 template<unsigned int dim>
972 {
973  if (dim == 1) return;
974  FEValues<3> fe_values_sub(*feo->q<dim-1>(), *feo->fe<dim-1>(),
976  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
978 
979  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
980  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
981  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
982  vector<vector<int> > side_dof_indices(2);
983  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
984  vector<double> frac_sigma(qsize), potential(qsize);
985  PetscScalar local_rhs[2][ndofs_side];
986  auto vec_side = fe_values_side.vector_view(0);
987  auto vec_sub = fe_values_sub.vector_view(0);
988 
989  // index 0 = element with lower dimension,
990  // index 1 = side of element with higher dimension
991  side_dof_indices[0].resize(ndofs_sub);
992  side_dof_indices[1].resize(ndofs_side);
993 
994  // assemble integral over sides
995  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
996  {
997  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
998  // skip neighbours of different dimension
999  if (nb->element()->dim() != dim-1) continue;
1000 
1001  ElementAccessor<3> cell_sub = nb->element();
1002  feo->dh()->cell_accessor_from_element(cell_sub.idx()).get_dof_indices(side_dof_indices[0]);
1003  fe_values_sub.reinit(cell_sub);
1004 
1005  ElementAccessor<3> cell = nb->side()->element();
1006  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices[1]);
1007  fe_values_side.reinit(*nb->side());
1008 
1009  // Element id's for testing if they belong to local partition.
1010  bool own_element_id[2];
1011  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
1012  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
1013 
1014  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
1015  data_.potential_load.value_list(fe_values_sub.point_list(), cell, potential);
1016 
1017  for (unsigned int n=0; n<2; ++n)
1018  for (unsigned int i=0; i<ndofs_side; i++)
1019  local_rhs[n][i] = 0;
1020 
1021  // set transmission conditions
1022  for (unsigned int k=0; k<qsize; k++)
1023  {
1024  arma::vec3 nv = fe_values_side.normal_vector(k);
1025 
1026  for (int n=0; n<2; n++)
1027  {
1028  if (!own_element_id[n]) continue;
1029 
1030  for (unsigned int i=0; i<n_dofs[n]; i++)
1031  {
1032  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
1033  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
1034 
1035  local_rhs[n][i] -= frac_sigma[k]*arma::dot(vf-vi,potential[k]*nv)*fe_values_sub.JxW(k);
1036  }
1037  }
1038  }
1039 
1040  for (unsigned int n=0; n<2; ++n)
1041  ls->rhs_set_values(n_dofs[n], side_dof_indices[n].data(), local_rhs[n]);
1042  }
1043 }
1044 
1045 
1046 
1047 
1048 
1049 template<unsigned int dim>
1051 {
1052  FEValues<3> fe_values_side(*feo->q<dim-1>(), *feo->fe<dim>(),
1054  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1055  vector<int> side_dof_indices(ndofs);
1056  // unsigned int loc_b=0;
1057  double local_rhs[ndofs];
1058  vector<PetscScalar> local_flux_balance_vector(ndofs);
1059  // PetscScalar local_flux_balance_rhs;
1060  vector<arma::vec3> bc_values(qsize), bc_traction(qsize);
1061  vector<double> csection(qsize), bc_potential(qsize);
1062  auto vec = fe_values_side.vector_view(0);
1063 
1064  for (DHCellAccessor cell : feo->dh()->own_range())
1065  {
1066  ElementAccessor<3> elm = cell.elm();
1067  if (elm->boundary_idx_ == nullptr) continue;
1068 
1069  for (unsigned int si=0; si<elm->n_sides(); ++si)
1070  {
1071  Edge edg = elm.side(si)->edge();
1072  Side side = *cell.elm().side(si);
1073 
1074  if (edg.n_sides() > 1) continue;
1075  // skip edges lying not on the boundary
1076  if ( ! side.is_boundary()) continue;
1077 
1078  if (side.dim() != dim-1)
1079  {
1080  // if (edg.side(0)->cond() != nullptr) ++loc_b;
1081  continue;
1082  }
1083 
1084  ElementAccessor<3> bc_cell = side.cond().element_accessor();
1085  unsigned int bc_type = data_.bc_type.value(side.centre(), bc_cell);
1086 
1087  fe_values_side.reinit(side);
1088 
1089  data_.cross_section.value_list(fe_values_side.point_list(), elm, csection);
1090  // The b.c. data are fetched for all possible b.c. types since we allow
1091  // different bc_type for each substance.
1092  data_.bc_displacement.value_list(fe_values_side.point_list(), bc_cell, bc_values);
1093  data_.bc_traction.value_list(fe_values_side.point_list(), bc_cell, bc_traction);
1094  data_.potential_load.value_list(fe_values_side.point_list(), elm, bc_potential);
1095 
1096  cell.get_dof_indices(side_dof_indices);
1097 
1098  fill_n(local_rhs, ndofs, 0);
1099  local_flux_balance_vector.assign(ndofs, 0);
1100  // local_flux_balance_rhs = 0;
1101 
1102  if (bc_type == EqData::bc_type_displacement)
1103  {
1104  for (unsigned int k=0; k<qsize; k++)
1105  for (unsigned int i=0; i<ndofs; i++)
1106  local_rhs[i] += dirichlet_penalty(side)*arma::dot(bc_values[k],vec.value(i,k))*fe_values_side.JxW(k);
1107  }
1108  else if (bc_type == EqData::bc_type_displacement_normal)
1109  {
1110  for (unsigned int k=0; k<qsize; k++)
1111  for (unsigned int i=0; i<ndofs; i++)
1112  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);
1113  }
1114  else if (bc_type == EqData::bc_type_traction)
1115  {
1116  for (unsigned int k=0; k<qsize; k++)
1117  {
1118  for (unsigned int i=0; i<ndofs; i++)
1119  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);
1120  }
1121  }
1122  ls->rhs_set_values(ndofs, side_dof_indices.data(), local_rhs);
1123 
1124 
1125 
1126 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
1127 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
1128  // ++loc_b;
1129  }
1130  }
1131 }
1132 
1133 
1134 
1135 
TimeGovernor & time()
Definition: equation.hh:151
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:229
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:329
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:567
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:864
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:427
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
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:719
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:352
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:443
std::shared_ptr< DOFHandlerMultiDim > dh()
Definition: elasticity.cc:119
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
Definition: elasticity.cc:856
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:971
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:533
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:346
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:36
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:466
void assemble_sources()
Assembles the right hand side vector due to volume sources.
Definition: elasticity.cc:740
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:365
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
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.
Mesh * mesh_
Definition: equation.hh:220
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:431
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:393
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:512
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:445
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:805
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:260
void assemble_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space ...
Definition: elasticity.cc:1050
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:229
void next_time()
Pass to next time and update equation data.
Definition: elasticity.cc:457
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
void initialize() override
Definition: elasticity.cc:294
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:645
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:671
virtual const Mat * get_matrix()
Definition: linsys.hh:187
void calculate_cumulative_balance()
Definition: elasticity.cc:630
Elasticity(Mesh &init_mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
Constructor.
Definition: elasticity.cc:258
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
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:221
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:795
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