Flow123d  release_3.0.0-1165-gffc6a5f
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/mapping_p1.hh"
25 #include "fem/fe_values.hh"
26 #include "fem/fe_p.hh"
27 #include "fem/fe_rt.hh"
28 #include "fem/fe_system.hh"
29 #include "fields/field_fe.hh"
30 #include "la/linsys_PETSC.hh"
31 #include "coupling/balance.hh"
32 #include "mesh/neighbours.h"
33 
34 #include "fields/multi_field.hh"
35 #include "fields/generic_field.hh"
36 #include "input/factory.hh"
37 
38 
39 
40 
41 using namespace Input::Type;
42 
43 
45  std::string equation_name = std::string(Elasticity::EqData::name()) + "_FE";
46  return IT::Record(
47  std::string(equation_name),
48  "FEM for linear elasticity.")
50  "Time governor setting for the secondary equation.")
51  .declare_key("balance", Balance::get_input_type(), Default("{}"),
52  "Settings for computing balance.")
54  "Parameters of output stream.")
56  "Linear solver for elasticity.")
57  .declare_key("input_fields", Array(
59  .make_field_descriptor_type(equation_name)),
61  "Input fields of the equation.")
62  .declare_key("output",
63  EqData().output_fields.make_output_type(equation_name, ""),
64  IT::Default("{ \"fields\": [ "+Elasticity::EqData::default_output_field()+" ] }"),
65  "Setting of the field output.")
66  .close();
67 }
68 
69 const int Elasticity::registrar =
70  Input::register_class< Elasticity, Mesh &, const Input::Record>(std::string(Elasticity::EqData::name()) + "_FE") +
72 
73 
74 
75 namespace Mechanics {
76 
77 FEObjects::FEObjects(Mesh *mesh_, unsigned int fe_order)
78 {
79  unsigned int q_order;
80 
81  switch (fe_order)
82  {
83  case 1:
84  q_order = 2;
85  fe0_ = new FESystem<0>(std::make_shared<FE_P<0> >(1), FEVector, 3);
86  fe1_ = new FESystem<1>(std::make_shared<FE_P<1> >(1), FEVector, 3);
87  fe2_ = new FESystem<2>(std::make_shared<FE_P<2> >(1), FEVector, 3);
88  fe3_ = new FESystem<3>(std::make_shared<FE_P<3> >(1), FEVector, 3);
89  break;
90 
91  default:
92  q_order=0;
93  xprintf(PrgErr, "Unsupported polynomial order %d for finite elements in Elasticity", fe_order);
94  break;
95  }
96 
97  for (unsigned int dim = 0; dim < 4; dim++) q_[dim] = new QGauss(dim, q_order);
98 
99  map1_ = new MappingP1<1,3>;
100  map2_ = new MappingP1<2,3>;
101  map3_ = new MappingP1<3,3>;
102 
103  ds_ = std::make_shared<EqualOrderDiscreteSpace>(mesh_, fe0_, fe1_, fe2_, fe3_);
104  dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
105 
106  dh_->distribute_dofs(ds_);
107 }
108 
109 
111 {
112  delete fe1_;
113  delete fe2_;
114  delete fe3_;
115  for (unsigned int dim=0; dim < 4; dim++) delete q_[dim];
116  delete map1_;
117  delete map2_;
118  delete map3_;
119 }
120 
121 template<> FiniteElement<0> *FEObjects::fe<0>() { return 0; }
122 template<> FiniteElement<1> *FEObjects::fe<1>() { return fe1_; }
123 template<> FiniteElement<2> *FEObjects::fe<2>() { return fe2_; }
124 template<> FiniteElement<3> *FEObjects::fe<3>() { return fe3_; }
125 
126 template<> MappingP1<1,3> *FEObjects::mapping<1>() { return map1_; }
127 template<> MappingP1<2,3> *FEObjects::mapping<2>() { return map2_; }
128 template<> MappingP1<3,3> *FEObjects::mapping<3>() { return map3_; }
129 
130 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh() { return dh_; }
131 
132 
133 
134 } // namespace Mechanics
135 
136 
137 
138 
139 
140 
141 
143  return Selection("Elasticity_BC_Type", "Types of boundary conditions for heat transfer model.")
144  .add_value(bc_type_displacement, "displacement",
145  "Prescribed displacement.")
146  .add_value(bc_type_traction, "traction",
147  "Prescribed traction.")
148  .close();
149 }
150 
151 
153 {
154  *this+=bc_type
155  .name("bc_type")
156  .description(
157  "Type of boundary condition.")
158  .units( UnitSI::dimensionless() )
159  .input_default("\"traction\"")
160  .input_selection( get_bc_type_selection() )
162 
163  *this+=bc_displacement
164  .name("bc_displacement")
165  .description("Prescribed displacement.")
166  .units( UnitSI().m() )
167  .input_default("0.0")
168  .flags_add(in_rhs);
169 
170  *this+=bc_traction
171  .name("bc_traction")
172  .description("Prescribed traction.")
173  .units( UnitSI().Pa() )
174  .input_default("0.0")
175  .flags_add(in_rhs);
176 
177  *this+=load
178  .name("load")
179  .description("Prescribed load.")
180  .units( UnitSI().N().m(-3) )
181  .input_default("0.0")
182  .flags_add(in_rhs);
183 
184  *this+=young_modulus
185  .name("young_modulus")
186  .description("Young modulus.")
187  .units( UnitSI().Pa() )
188  .input_default("0.0")
189  .flags_add(in_main_matrix);
190 
191  *this+=poisson_ratio
192  .name("poisson_ratio")
193  .description("Poisson ratio.")
194  .units( UnitSI().dimensionless() )
195  .input_default("0.0")
196  .flags_add(in_main_matrix);
197 
198  *this+=fracture_sigma
199  .name("fracture_sigma")
200  .description(
201  "Coefficient of diffusive transfer through fractures (for each substance).")
202  .units( UnitSI::dimensionless() )
203  .input_default("1.0")
204  .flags_add(FieldFlag::in_main_matrix);
205 
206  *this += region_id.name("region_id")
207  .units( UnitSI::dimensionless())
209 
210  *this += subdomain.name("subdomain")
211  .units( UnitSI::dimensionless() )
213 
214  *this+=cross_section
215  .name("cross_section")
216  .units( UnitSI().m(3).md() )
217  .flags(input_copy & in_time_term & in_main_matrix);
218 
219  *this+=output_field
220  .name("displacement")
221  .units( UnitSI().m() )
222  .flags(equation_result);
223 
224 
225  // add all input fields to the output list
226  output_fields += *this;
227 
228 }
229 
230 Elasticity::Elasticity(Mesh & init_mesh, const Input::Record in_rec)
231  : EquationBase(init_mesh, in_rec),
232  input_rec(in_rec),
233  allocation_done(false)
234 {
235  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
236  // due to constexpr optimization.
238 
239  this->eq_data_ = &data_;
240  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
241 
242 
243  // Set up physical parameters.
244  data_.set_mesh(init_mesh);
247 
248  // create finite element structures and distribute DOFs
249  feo = new Mechanics::FEObjects(mesh_, 1);
250  DebugOut().fmt("Mechanics: solution size {}\n", feo->dh()->n_global_dofs());
251 
252 }
253 
254 
256 {
257  output_stream_ = OutputTime::create_output_stream("mechanics", input_rec.val<Input::Record>("output_stream"), time().get_unit_string());
258 
259  data_.set_components({"displacement"});
260  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), time() );
261 
262 // balance_ = std::make_shared<Balance>("mechanics", mesh_);
263 // balance_->init_from_input(input_rec.val<Input::Record>("balance"), *time_);
264  // initialization of balance object
265 // data_.balance_idx_ = {
266 // balance_->add_quantity("force_x"),
267 // balance_->add_quantity("force_y"),
268 // balance_->add_quantity("force_z")
269 // };
270 // balance_->units(UnitSI().kg().m().s(-2));
271 
272  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
273  std::shared_ptr<FieldFE<3, FieldValue<3>::VectorFixed> > output_field_ptr(new FieldFE<3, FieldValue<3>::VectorFixed>);
274  output_vec = output_field_ptr->set_fe_data(feo->dh());
275  data_.output_field.set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0.);
277 
278  // set time marks for writing the output
280 
281  // equation default PETSc solver options
282  std::string petsc_default_opts;
283  petsc_default_opts = "-ksp_type cg -pc_type hypre -pc_hypre_type boomeramg";
284 
285  // allocate matrix and vector structures
286  ls = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
289 
290  // initialization of balance object
291 // balance_->allocate(feo->dh()->distr()->lsize(),
292 // max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
293 
294 }
295 
296 
298 {
299 // delete time_;
300 
301  MatDestroy(&stiffness_matrix);
302  VecDestroy(&rhs);
303  delete ls;
304  delete feo;
305 
306 }
307 
308 
309 
311 {
312  VecScatter output_scatter;
313  VecScatterCreateToZero(ls->get_solution(), &output_scatter, PETSC_NULL);
314  // gather solution to output_vec
315  VecScatterBegin(output_scatter, ls->get_solution(), output_vec.petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
316  VecScatterEnd(output_scatter, ls->get_solution(), output_vec.petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
317  VecScatterDestroy(&(output_scatter));
318 }
319 
320 
321 
322 
324 {
328  std::stringstream ss; // print warning message with table of uninitialized fields
329  if ( FieldCommon::print_message_table(ss, "mechanics") ) {
330  WarningOut() << ss.str();
331  }
332 
333  ( (LinSys_PETSC *)ls )->set_initial_guess_nonzero();
334 
335  // check first time assembly - needs preallocation
337 
338 
339  // after preallocation we assemble the matrices and vectors required for balance of forces
340 // for (auto subst_idx : data_.balance_idx_)
341 // balance_->calculate_instant(subst_idx, ls->get_solution());
342 
343 // update_solution();
345  MatSetOption(*ls->get_matrix(), MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
346  ls->mat_zero_entries();
348  set_sources();
350  ls->finish_assembly();
351  ls->apply_constrains(1.0);
352  ls->solve();
353  output_data();
354 }
355 
356 
357 
359 {
360  // preallocate system matrix
361  ls->start_allocation();
362  stiffness_matrix = NULL;
363  rhs = NULL;
364 
366  set_sources();
368 
369  allocation_done = true;
370 }
371 
372 
373 
374 
376 {
377  START_TIMER("DG-ONE STEP");
378 
379  time_->next_time();
380  time_->view("MECH");
381 
382  START_TIMER("data reinit");
384  END_TIMER("data reinit");
385 
386  // assemble stiffness matrix
387  if (stiffness_matrix == NULL
389  || flux_changed)
390  {
392  ls->mat_zero_entries();
394  ls->finish_assembly();
395 
396  if (stiffness_matrix == NULL)
397  MatConvert(*( ls->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix);
398  else
399  MatCopy(*( ls->get_matrix() ), stiffness_matrix, DIFFERENT_NONZERO_PATTERN);
400  }
401 
402  // assemble right hand side (due to sources and boundary conditions)
403  if (rhs == NULL
404  || data_.subset(FieldFlag::in_rhs).changed()
405  || flux_changed)
406  {
408  ls->rhs_zero_entries();
409  set_sources();
411  ls->finish_assembly();
412  ls->apply_constrains(1.0);
413 
414  if (rhs == nullptr) VecDuplicate(*( ls->get_rhs() ), &rhs);
415  VecCopy(*( ls->get_rhs() ), rhs);
416  }
417 
418  flux_changed = false;
419 
420 
421  START_TIMER("solve");
422  ls->solve();
423  END_TIMER("solve");
424 
426 
427  output_data();
428 
429  END_TIMER("DG-ONE STEP");
430 }
431 
432 
433 
434 
436 {
437  START_TIMER("MECH-OUTPUT");
438 
439  // gather the solution from all processors
440  data_.output_fields.set_time( this->time().step(), LimitSide::left);
441  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
443  data_.output_fields.output(this->time().step());
444  output_stream_->write_time_frame();
445 
446 // START_TIMER("MECH-balance");
447 // balance_->calculate_instant(subst_idx, ls->get_solution());
448 // balance_->output();
449 // END_TIMER("MECH-balance");
450 
451  END_TIMER("MECH-OUTPUT");
452 }
453 
454 
455 
457 {
458 // if (balance_->cumulative())
459 // {
460 // balance_->calculate_cumulative(subst_idx, ls->get_solution());
461 // }
462 }
463 
464 
465 
466 
467 
468 
469 
470 
472 {
473  START_TIMER("assemble_stiffness");
474  START_TIMER("assemble_volume_integrals");
475  assemble_volume_integrals<1>();
476  assemble_volume_integrals<2>();
477  assemble_volume_integrals<3>();
478  END_TIMER("assemble_volume_integrals");
479 
480  START_TIMER("assemble_fluxes_boundary");
481  assemble_fluxes_boundary<1>();
482  assemble_fluxes_boundary<2>();
483  assemble_fluxes_boundary<3>();
484  END_TIMER("assemble_fluxes_boundary");
485 
486  START_TIMER("assemble_fluxes_elem_side");
487  assemble_fluxes_element_side<1>();
488  assemble_fluxes_element_side<2>();
489  assemble_fluxes_element_side<3>();
490  END_TIMER("assemble_fluxes_elem_side");
491  END_TIMER("assemble_stiffness");
492 }
493 
494 
495 double lame_mu(double young, double poisson)
496 {
497  return young*0.5/(poisson+1.);
498 }
499 
500 
501 double lame_lambda(double young, double poisson)
502 {
503  return young*poisson/((poisson+1.)*(1.-2.*poisson));
504 }
505 
506 
507 template<unsigned int dim>
509 {
510  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
512  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
513  vector<int> dof_indices(ndofs);
514  vector<arma::vec3> velocity(qsize);
515  vector<double> young(qsize), poisson(qsize), csection(qsize);
516  PetscScalar local_matrix[ndofs*ndofs];
517  auto vec = fe_values.vector_view(0);
518 
519  // assemble integral over elements
520  for (auto cell : feo->dh()->own_range())
521  {
522  if (cell.dim() != dim) continue;
523  ElementAccessor<3> elm_acc = cell.elm();
524 
525  fe_values.reinit(elm_acc);
526  cell.get_dof_indices(dof_indices);
527 
528  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
529  data_.young_modulus.value_list(fe_values.point_list(), elm_acc, young);
530  data_.poisson_ratio.value_list(fe_values.point_list(), elm_acc, poisson);
531 
532  // assemble the local stiffness matrix
533  for (unsigned int i=0; i<ndofs; i++)
534  for (unsigned int j=0; j<ndofs; j++)
535  local_matrix[i*ndofs+j] = 0;
536 
537  for (unsigned int k=0; k<qsize; k++)
538  {
539  double mu = lame_mu(young[k], poisson[k]);
540  double lambda = lame_lambda(young[k], poisson[k]);
541  for (unsigned int i=0; i<ndofs; i++)
542  {
543  for (unsigned int j=0; j<ndofs; j++)
544  local_matrix[i*ndofs+j] += csection[k]*(
545  2*mu*arma::dot(vec.sym_grad(j,k), vec.sym_grad(i,k))
546  + lambda*vec.divergence(j,k)*vec.divergence(i,k)
547  )*fe_values.JxW(k);
548  }
549  }
550  ls->mat_set_values(ndofs, dof_indices.data(), ndofs, dof_indices.data(), local_matrix);
551  }
552 }
553 
554 
555 
557 {
558  START_TIMER("assemble_sources");
559 // balance_->start_source_assembly(subst_idx);
560  set_sources<1>();
561  set_sources<2>();
562  set_sources<3>();
563 // balance_->finish_source_assembly(subst_idx);
564  END_TIMER("assemble_sources");
565 }
566 
567 
568 template<unsigned int dim>
570 {
571  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
573  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
574  vector<arma::vec3> load(qsize);
575  vector<double> csection(qsize);
576  vector<int> dof_indices(ndofs);
577  PetscScalar local_rhs[ndofs];
578  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
579  auto vec = fe_values.vector_view(0);
580 
581  // assemble integral over elements
582  for (auto cell : feo->dh()->own_range())
583  {
584  if (cell.dim() != dim) continue;
585  ElementAccessor<3> elm_acc = cell.elm();
586 
587  fe_values.reinit(elm_acc);
588  cell.get_dof_indices(dof_indices);
589 
590  // assemble the local stiffness matrix
591  fill_n(local_rhs, ndofs, 0);
592  local_source_balance_vector.assign(ndofs, 0);
593  local_source_balance_rhs.assign(ndofs, 0);
594 
595  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
596  data_.load.value_list(fe_values.point_list(), elm_acc, load);
597 
598  // compute sources
599  for (unsigned int k=0; k<qsize; k++)
600  {
601  for (unsigned int i=0; i<ndofs; i++)
602  local_rhs[i] += arma::dot(load[k], vec.value(i,k))*csection[k]*fe_values.JxW(k);
603  }
604  ls->rhs_set_values(ndofs, &(dof_indices[0]), local_rhs);
605 
606 // for (unsigned int i=0; i<ndofs; i++)
607 // {
608 // for (unsigned int k=0; k<qsize; k++)
609 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_values[vec].value(i,k)*fe_values.JxW(k);
610 //
611 // local_source_balance_rhs[i] += local_rhs[i];
612 // }
613 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_vector);
614 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_rhs);
615  }
616 }
617 
618 
619 
620 
621 
622 
623 
624 template<unsigned int dim>
626 {
627  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
629  const unsigned int ndofs = feo->fe<dim>()->n_dofs();
630  vector<int> side_dof_indices(ndofs);
631  PetscScalar local_matrix[ndofs*ndofs];
632 
633  // assemble boundary integral
634  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
635  {
636  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
637  if (edg->n_sides > 1) continue;
638  // check spatial dimension
639  if (edg->side(0)->dim() != dim-1) continue;
640  // skip edges lying not on the boundary
641  if (edg->side(0)->cond() == NULL) continue;
642 
643  SideIter side = edg->side(0);
644  ElementAccessor<3> cell = side->element();
645  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices);
646  fe_values_side.reinit(cell, side->side_idx());
647 // unsigned int bc_type = data_.bc_type.value(side->centre(), side->cond()->element_accessor());
648 
649 // for (unsigned int i=0; i<ndofs; i++)
650 // for (unsigned int j=0; j<ndofs; j++)
651 // local_matrix[i*ndofs+j] = 0;
652 //
653 // ls->mat_set_values(ndofs, side_dof_indices.data(), ndofs, side_dof_indices.data(), local_matrix);
654  }
655 }
656 
657 
659 {
660  arma::mat33 mt = m - m*arma::kron(n,n.t());
661  return mt;
662 }
663 
664 
665 template<unsigned int dim>
667 {
668  if (dim == 1) return;
669  FEValues<dim-1,3> fe_values_sub(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
671  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
673 
674  vector<FEValuesSpaceBase<3>*> fv_sb(2);
675  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
676  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
677  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
678  vector<vector<int> > side_dof_indices(2);
679  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
680  vector<double> frac_sigma(qsize);
681  vector<double> csection_lower(qsize), csection_higher(qsize), young(qsize), poisson(qsize), alpha(qsize);
682  PetscScalar local_matrix[2][2][(ndofs_side)*(ndofs_side)];
683  PetscScalar local_rhs[2][ndofs_side];
684  auto vec_side = fe_values_side.vector_view(0);
685  auto vec_sub = fe_values_sub.vector_view(0);
686 
687  // index 0 = element with lower dimension,
688  // index 1 = side of element with higher dimension
689  side_dof_indices[0].resize(ndofs_sub);
690  side_dof_indices[1].resize(ndofs_side);
691  fv_sb[0] = &fe_values_sub;
692  fv_sb[1] = &fe_values_side;
693 
694  // assemble integral over sides
695  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
696  {
697  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
698  // skip neighbours of different dimension
699  if (nb->element()->dim() != dim-1) continue;
700 
701  ElementAccessor<3> cell_sub = nb->element();
702  feo->dh()->cell_accessor_from_element(cell_sub.idx()).get_dof_indices(side_dof_indices[0]);
703  fe_values_sub.reinit(cell_sub);
704 
705  ElementAccessor<3> cell = nb->side()->element();
706  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices[1]);
707  fe_values_side.reinit(cell, nb->side()->side_idx());
708 
709  // Element id's for testing if they belong to local partition.
710  bool own_element_id[2];
711  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
712  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
713 
714  data_.cross_section.value_list(fe_values_sub.point_list(), cell_sub, csection_lower);
715  data_.cross_section.value_list(fe_values_sub.point_list(), cell, csection_higher);
716  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
717  data_.young_modulus.value_list(fe_values_sub.point_list(), cell_sub, young);
718  data_.poisson_ratio.value_list(fe_values_sub.point_list(), cell_sub, poisson);
719 
720  for (unsigned int n=0; n<2; ++n)
721  {
722  for (unsigned int i=0; i<ndofs_side; i++)
723  {
724  for (unsigned int m=0; m<2; ++m)
725  {
726  for (unsigned int j=0; j<ndofs_side; j++)
727  local_matrix[n][m][i*(ndofs_side)+j] = 0;
728  }
729  local_rhs[n][i] = 0;
730  }
731  }
732 
733  // set transmission conditions
734  for (unsigned int k=0; k<qsize; k++)
735  {
736  arma::vec3 nv = fe_values_side.normal_vector(k);
737  double mu = lame_mu(young[k], poisson[k]);
738  double lambda = lame_lambda(young[k], poisson[k]);
739 
740  for (int n=0; n<2; n++)
741  {
742  if (!own_element_id[n]) continue;
743 
744  for (unsigned int i=0; i<n_dofs[n]; i++)
745  {
746  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
747  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
748  arma::mat33 gvft = (n==0)?vec_sub.grad(i,k):arma::zeros(3,3);
749 
750  for (int m=0; m<2; m++)
751  for (unsigned int j=0; j<n_dofs[m]; j++) {
752  arma::vec3 ui = (m==0)?arma::zeros(3):vec_side.value(j,k);
753  arma::vec3 uf = (m==1)?arma::zeros(3):vec_sub.value(j,k);
754  arma::mat33 guit = (m==1)?mat_t(vec_side.grad(j,k),nv):arma::zeros(3,3);
755  double divuit = (m==1)?arma::trace(guit):0;
756 
757  local_matrix[n][m][i*n_dofs[m] + j] +=
758  (
759  frac_sigma[k]*arma::dot(vf-vi,
760  2/csection_lower[k]*(mu*(uf-ui)+(mu+lambda)*(arma::dot(uf-ui,nv)*nv))
761  + mu*arma::trans(guit)*nv
762  + lambda*divuit*nv
763  )
764  - arma::dot(gvft, mu*arma::kron(nv,ui.t()) + lambda*arma::dot(ui,nv)*arma::eye(3,3))
765  )*fv_sb[0]->JxW(k);
766  }
767  }
768  }
769  }
770 
771  for (unsigned int n=0; n<2; ++n)
772  {
773  for (unsigned int m=0; m<2; ++m)
774  ls->mat_set_values(n_dofs[n], side_dof_indices[n].data(), n_dofs[m], side_dof_indices[m].data(), local_matrix[n][m]);
775 
776  ls->rhs_set_values(n_dofs[n], side_dof_indices[n].data(), local_rhs[n]);
777  }
778  }
779 }
780 
781 
782 
783 
784 
785 
786 
788 {
789  START_TIMER("assemble_bc");
790 // balance_->start_flux_assembly(subst_idx);
791  set_boundary_conditions<1>();
792  set_boundary_conditions<2>();
793  set_boundary_conditions<3>();
794 // balance_->finish_flux_assembly(subst_idx);
795  END_TIMER("assemble_bc");
796 }
797 
798 
799 
800 template<unsigned int dim>
802 {
803  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
805  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
806  vector<int> side_dof_indices(ndofs);
807  unsigned int loc_b=0;
808  double local_rhs[ndofs];
809  vector<PetscScalar> local_flux_balance_vector(ndofs);
810  PetscScalar local_flux_balance_rhs;
811  vector<arma::vec3> bc_values(qsize), bc_traction(qsize);
812  vector<double> csection(qsize);
813  auto vec = fe_values_side.vector_view(0);
814 
815  for (auto cell : feo->dh()->own_range())
816  {
817  ElementAccessor<3> elm = cell.elm();
818  if (elm->boundary_idx_ == nullptr) continue;
819 
820  for (unsigned int si=0; si<elm->n_sides(); ++si)
821  {
822  const Edge *edg = elm.side(si)->edge();
823  if (edg->n_sides > 1) continue;
824  // skip edges lying not on the boundary
825  if (edg->side(0)->cond() == NULL) continue;
826 
827  if (edg->side(0)->dim() != dim-1)
828  {
829  if (edg->side(0)->cond() != nullptr) ++loc_b;
830  continue;
831  }
832 
833  SideIter side = edg->side(0);
834  ElementAccessor<3> cell = side->element();
835  ElementAccessor<3> bc_cell = side->cond()->element_accessor();
836 
837  unsigned int bc_type = data_.bc_type.value(side->centre(), bc_cell);
838 
839  fe_values_side.reinit(cell, side->side_idx());
840 
841  data_.cross_section.value_list(fe_values_side.point_list(), cell, csection);
842  // The b.c. data are fetched for all possible b.c. types since we allow
843  // different bc_type for each substance.
844  data_.bc_displacement.value_list(fe_values_side.point_list(), bc_cell, bc_values);
845  data_.bc_traction.value_list(fe_values_side.point_list(), bc_cell, bc_traction);
846 
847  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices);
848 
849  fill_n(local_rhs, ndofs, 0);
850  local_flux_balance_vector.assign(ndofs, 0);
851  local_flux_balance_rhs = 0;
852 
853  if (bc_type == EqData::bc_type_displacement)
854  {
855  // We solve a local problem to determine which local dof
856  // to assign the boundary value to.
857  // TODO: Should be possibly optimized by inverting the matrix only once within reference element.
858  arma::mat d_mat(ndofs,ndofs);
859  arma::vec d_vec(ndofs);
860  for (unsigned int i=0; i<ndofs; i++)
861  {
862  d_vec(i) = 0;
863  for (unsigned int k=0; k<qsize; k++)
864  d_vec(i) += arma::dot(vec.value(i,k),bc_values[k])*fe_values_side.JxW(k);
865  for (unsigned int j=i; j<ndofs; j++)
866  {
867  d_mat(i,j) = 0;
868  for (unsigned int k=0; k<qsize; k++)
869  d_mat(i,j) += arma::dot(vec.value(i,k),vec.value(j,k))*fe_values_side.JxW(k);
870  d_mat(j,i) = d_mat(i,j);
871  }
872  }
873  arma::vec d_val = pinv(d_mat)*d_vec;
874 
875  for (unsigned int i=0; i<ndofs; i++)
876  if (norm(d_mat.row(i)) > 0)
877  ls->add_constraint(side_dof_indices[i], d_val(i));
878  }
879  else if (bc_type == EqData::bc_type_traction)
880  {
881  for (unsigned int k=0; k<qsize; k++)
882  {
883  for (unsigned int i=0; i<ndofs; i++)
884  local_rhs[i] += csection[k]*arma::dot(vec.value(i,k),bc_traction[k])*fe_values_side.JxW(k);
885  }
886  }
887  ls->rhs_set_values(ndofs, &(side_dof_indices[0]), local_rhs);
888 
889 
890 
891 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
892 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
893  ++loc_b;
894  }
895  }
896 }
897 
898 
899 
900 
TimeGovernor & time()
Definition: equation.hh:148
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:121
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:211
const Edge * edge() const
Definition: side_impl.hh:66
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:122
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
Field< 3, FieldValue< 3 >::Scalar > young_modulus
Definition: elasticity.hh:125
FieldSet * eq_data_
Definition: equation.hh:232
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:49
static auto subdomain(Mesh &mesh) -> IndexField
static constexpr const char * name()
Definition: elasticity.hh:112
Transformed quadrature weight for cell sides.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
Mat< double, N, M > mat
Definition: armor.hh:214
unsigned int * boundary_idx_
Definition: elements.h:83
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:615
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void preallocate()
Definition: elasticity.cc:358
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: elasticity.hh:131
unsigned int side_idx() const
Definition: side_impl.hh:82
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:38
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
void set_from_input(const Input::Record in_rec) override
virtual void start_add_assembly()
Definition: linsys.hh:341
void output(TimeStep step)
Boundary * cond() const
Definition: side_impl.hh:71
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:264
VectorMPI output_vec
Vector of solution data.
Definition: elasticity.hh:296
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:44
arma::vec3 centre() const
Centre of side.
Definition: sides.cc:121
FiniteElement< dim > * fe()
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
Definition: elasticity.cc:787
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:185
LinSys * ls
Linear algebra system for the transport equation.
Definition: elasticity.hh:287
Definition: mesh.h:76
Fields computed from the mesh data.
virtual void start_allocation()
Definition: linsys.hh:333
Class FEValues calculates finite element data on the actual cells such as shape function values...
Mat< double, N, 1 > vec
Definition: armor.hh:211
virtual void finish_assembly()=0
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:130
Class FESystem for compound finite elements.
int n_sides
Definition: edges.h:36
~Elasticity()
Destructor.
Definition: elasticity.cc:297
Definition: edges.h:26
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
double lame_mu(double young, double poisson)
Definition: elasticity.cc:495
const RegionDB & region_db() const
Definition: mesh.h:143
std::shared_ptr< DOFHandlerMultiDim > dh()
const TimeStep & step(int index=-1) const
void add_constraint(int row, double value)
Definition: linsys.hh:494
Elasticity(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
Definition: elasticity.cc:230
void update_solution() override
Computes the solution in one time instant.
Definition: elasticity.cc:375
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
Definition: elasticity.cc:658
Basic time management functionality for unsteady (and steady) solvers (class Equation).
ElementAccessor< 3 > element() const
Definition: side_impl.hh:53
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list) const
Definition: field.hh:403
bool flux_changed
Indicator of change in advection vector field.
Definition: elasticity.hh:316
Input::Record input_rec
Record with input specification.
Definition: elasticity.hh:301
Class for declaration of inputs sequences.
Definition: type_base.hh:346
virtual void apply_constrains(double scalar)=0
Vec rhs
Vector of right hand side.
Definition: elasticity.hh:281
void view(const char *name="") const
Symmetric Gauss-Legendre quadrature formulae on simplices.
FEM for linear elasticity.
bool allocation_done
Indicates whether matrices have been preallocated.
Definition: elasticity.hh:313
unsigned int dim() const
Definition: elements.h:124
ElementAccessor< 3 > element()
Definition: neighbours.h:163
Transformed quadrature points.
unsigned int dim() const
Definition: side_impl.hh:38
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:123
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
const Vec & get_solution()
Definition: linsys.hh:282
double JxW(const unsigned int point_no) override
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:336
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:147
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:135
static auto region_id(Mesh &mesh) -> IndexField
#define xprintf(...)
Definition: system.hh:92
Type dot(const Mat< Type, nRows, nCols > &a, const Mat< Type, nRows, nCols > &b)
Definition: armor.hh:176
#define START_TIMER(tag)
Starts a timer with specified tag.
std::shared_ptr< DOFHandlerMultiDim > dh()
Definition: elasticity.cc:130
void set_sources()
Assembles the right hand side due to volume sources.
Definition: elasticity.cc:556
Mesh * mesh_
Definition: equation.hh:223
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:389
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:323
Vec & petsc_vec()
Getter for PETSC vector of output data (e.g. can be used by scatters).
Definition: vector_mpi.hh:120
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:435
Mat stiffness_matrix
The stiffness matrix.
Definition: elasticity.hh:284
Quadrature * q()
Definition: elasticity.hh:60
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:388
EqData & data()
Definition: elasticity.hh:184
virtual const Vec * get_rhs()
Definition: linsys.hh:203
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no) override
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:368
Normal vectors.
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:124
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:127
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
Definition: elasticity.cc:625
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:218
vector< Neighbour > vb_neighbours_
Definition: mesh.h:273
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:190
Field< 3, FieldValue< 3 >::Scalar > poisson_ratio
Definition: elasticity.hh:126
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.
void set_components(const std::vector< string > &names)
Definition: field_set.hh:177
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
Definition: elasticity.cc:666
EquationOutput output_fields
Definition: elasticity.hh:136
Definition: system.hh:64
static const int registrar
Registrar of class to factory.
Definition: elasticity.hh:194
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:198
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:252
void initialize() override
Definition: elasticity.cc:255
EqData data_
Field data for model parameters.
Definition: elasticity.hh:262
double lame_lambda(double young, double poisson)
Definition: elasticity.cc:501
static const Input::Type::Selection & get_bc_type_selection()
Definition: elasticity.cc:142
static const Input::Type::Record & get_input_type()
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
static string default_output_field()
Definition: elasticity.hh:114
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:357
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
Definition: elasticity.cc:471
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
Class for representation SI units of Fields.
Definition: unit_si.hh:40
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
Definition: elasticity.cc:508
virtual const Mat * get_matrix()
Definition: linsys.hh:187
void calculate_cumulative_balance()
Definition: elasticity.cc:456
MappingP1< dim, 3 > * mapping()
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:252
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:111
Field< 3, FieldValue< 3 >::VectorFixed > output_field
Definition: elasticity.hh:134
SideIter side(const unsigned int i) const
Definition: edges.h:31
void output_vector_gather()
Definition: elasticity.cc:310
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:132
Definitions of Raviart-Thomas finite elements.
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:48
std::shared_ptr< OutputTime > output_stream_
Definition: elasticity.hh:298
FEObjects(Mesh *mesh_, unsigned int fe_order)
TimeGovernor * time_
Definition: equation.hh:224
Definition: field.hh:56
Mechanics::FEObjects * feo
Finite element objects.
Definition: elasticity.hh:271
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:541
Transformed quadrature weights.
Calculates finite element data on a side.
Definition: fe_values.hh:582