Flow123d  release_3.0.0-1264-g45bfb2a
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 nullptr; }
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<0,3> *FEObjects::mapping<0>() { return nullptr; }
127 template<> MappingP1<1,3> *FEObjects::mapping<1>() { return map1_; }
128 template<> MappingP1<2,3> *FEObjects::mapping<2>() { return map2_; }
129 template<> MappingP1<3,3> *FEObjects::mapping<3>() { return map3_; }
130 
131 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh() { return dh_; }
132 
133 
134 
135 } // namespace Mechanics
136 
137 
138 
139 
140 
141 
142 
144  return Selection("Elasticity_BC_Type", "Types of boundary conditions for heat transfer model.")
145  .add_value(bc_type_displacement, "displacement",
146  "Prescribed displacement.")
147  .add_value(bc_type_traction, "traction",
148  "Prescribed traction.")
149  .close();
150 }
151 
152 
154 {
155  *this+=bc_type
156  .name("bc_type")
157  .description(
158  "Type of boundary condition.")
159  .units( UnitSI::dimensionless() )
160  .input_default("\"traction\"")
161  .input_selection( get_bc_type_selection() )
163 
164  *this+=bc_displacement
165  .name("bc_displacement")
166  .description("Prescribed displacement.")
167  .units( UnitSI().m() )
168  .input_default("0.0")
169  .flags_add(in_rhs);
170 
171  *this+=bc_traction
172  .name("bc_traction")
173  .description("Prescribed traction.")
174  .units( UnitSI().Pa() )
175  .input_default("0.0")
176  .flags_add(in_rhs);
177 
178  *this+=load
179  .name("load")
180  .description("Prescribed load.")
181  .units( UnitSI().N().m(-3) )
182  .input_default("0.0")
183  .flags_add(in_rhs);
184 
185  *this+=young_modulus
186  .name("young_modulus")
187  .description("Young modulus.")
188  .units( UnitSI().Pa() )
189  .input_default("0.0")
190  .flags_add(in_main_matrix);
191 
192  *this+=poisson_ratio
193  .name("poisson_ratio")
194  .description("Poisson ratio.")
195  .units( UnitSI().dimensionless() )
196  .input_default("0.0")
197  .flags_add(in_main_matrix);
198 
199  *this+=fracture_sigma
200  .name("fracture_sigma")
201  .description(
202  "Coefficient of diffusive transfer through fractures (for each substance).")
203  .units( UnitSI::dimensionless() )
204  .input_default("1.0")
205  .flags_add(FieldFlag::in_main_matrix);
206 
207  *this += region_id.name("region_id")
208  .units( UnitSI::dimensionless())
210 
211  *this += subdomain.name("subdomain")
212  .units( UnitSI::dimensionless() )
214 
215  *this+=cross_section
216  .name("cross_section")
217  .units( UnitSI().m(3).md() )
218  .flags(input_copy & in_time_term & in_main_matrix);
219 
220  *this+=output_field
221  .name("displacement")
222  .units( UnitSI().m() )
223  .flags(equation_result);
224 
225 
226  // add all input fields to the output list
227  output_fields += *this;
228 
229 }
230 
231 Elasticity::Elasticity(Mesh & init_mesh, const Input::Record in_rec)
232  : EquationBase(init_mesh, in_rec),
233  input_rec(in_rec),
234  allocation_done(false)
235 {
236  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
237  // due to constexpr optimization.
239 
240  this->eq_data_ = &data_;
241  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
242 
243 
244  // Set up physical parameters.
245  data_.set_mesh(init_mesh);
248 
249  // create finite element structures and distribute DOFs
250  feo = new Mechanics::FEObjects(mesh_, 1);
251  DebugOut().fmt("Mechanics: solution size {}\n", feo->dh()->n_global_dofs());
252 
253 }
254 
255 
257 {
258  output_stream_ = OutputTime::create_output_stream("mechanics", input_rec.val<Input::Record>("output_stream"), time().get_unit_string());
259 
260  data_.set_components({"displacement"});
261  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), time() );
262 
263 // balance_ = std::make_shared<Balance>("mechanics", mesh_);
264 // balance_->init_from_input(input_rec.val<Input::Record>("balance"), *time_);
265  // initialization of balance object
266 // data_.balance_idx_ = {
267 // balance_->add_quantity("force_x"),
268 // balance_->add_quantity("force_y"),
269 // balance_->add_quantity("force_z")
270 // };
271 // balance_->units(UnitSI().kg().m().s(-2));
272 
273  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
274  std::shared_ptr<FieldFE<3, FieldValue<3>::VectorFixed> > output_field_ptr(new FieldFE<3, FieldValue<3>::VectorFixed>);
275  output_vec = output_field_ptr->set_fe_data(feo->dh());
276  data_.output_field.set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0.);
278 
279  // set time marks for writing the output
281 
282  // equation default PETSc solver options
283  std::string petsc_default_opts;
284  petsc_default_opts = "-ksp_type cg -pc_type hypre -pc_hypre_type boomeramg";
285 
286  // allocate matrix and vector structures
287  ls = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
290 
291  // initialization of balance object
292 // balance_->allocate(feo->dh()->distr()->lsize(),
293 // max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
294 
295 }
296 
297 
299 {
300 // delete time_;
301 
302  MatDestroy(&stiffness_matrix);
303  VecDestroy(&rhs);
304  delete ls;
305  delete feo;
306 
307 }
308 
309 
310 
312 {
313  VecScatter output_scatter;
314  VecScatterCreateToZero(ls->get_solution(), &output_scatter, PETSC_NULL);
315  // gather solution to output_vec
316  VecScatterBegin(output_scatter, ls->get_solution(), output_vec.petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
317  VecScatterEnd(output_scatter, ls->get_solution(), output_vec.petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
318  VecScatterDestroy(&(output_scatter));
319 }
320 
321 
322 
323 
325 {
329  std::stringstream ss; // print warning message with table of uninitialized fields
330  if ( FieldCommon::print_message_table(ss, "mechanics") ) {
331  WarningOut() << ss.str();
332  }
333 
334  ( (LinSys_PETSC *)ls )->set_initial_guess_nonzero();
335 
336  // check first time assembly - needs preallocation
338 
339 
340  // after preallocation we assemble the matrices and vectors required for balance of forces
341 // for (auto subst_idx : data_.balance_idx_)
342 // balance_->calculate_instant(subst_idx, ls->get_solution());
343 
344 // update_solution();
346  MatSetOption(*ls->get_matrix(), MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
347  ls->mat_zero_entries();
349  set_sources();
351  ls->finish_assembly();
352  ls->apply_constrains(1.0);
353  ls->solve();
354  output_data();
355 }
356 
357 
358 
360 {
361  // preallocate system matrix
362  ls->start_allocation();
363  stiffness_matrix = NULL;
364  rhs = NULL;
365 
367  set_sources();
369 
370  allocation_done = true;
371 }
372 
373 
374 
375 
377 {
378  START_TIMER("DG-ONE STEP");
379 
380  time_->next_time();
381  time_->view("MECH");
382 
383  START_TIMER("data reinit");
385  END_TIMER("data reinit");
386 
387  // assemble stiffness matrix
388  if (stiffness_matrix == NULL
390  || flux_changed)
391  {
393  ls->mat_zero_entries();
395  ls->finish_assembly();
396 
397  if (stiffness_matrix == NULL)
398  MatConvert(*( ls->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix);
399  else
400  MatCopy(*( ls->get_matrix() ), stiffness_matrix, DIFFERENT_NONZERO_PATTERN);
401  }
402 
403  // assemble right hand side (due to sources and boundary conditions)
404  if (rhs == NULL
405  || data_.subset(FieldFlag::in_rhs).changed()
406  || flux_changed)
407  {
409  ls->rhs_zero_entries();
410  set_sources();
412  ls->finish_assembly();
413  ls->apply_constrains(1.0);
414 
415  if (rhs == nullptr) VecDuplicate(*( ls->get_rhs() ), &rhs);
416  VecCopy(*( ls->get_rhs() ), rhs);
417  }
418 
419  flux_changed = false;
420 
421 
422  START_TIMER("solve");
423  ls->solve();
424  END_TIMER("solve");
425 
427 
428  output_data();
429 
430  END_TIMER("DG-ONE STEP");
431 }
432 
433 
434 
435 
437 {
438  START_TIMER("MECH-OUTPUT");
439 
440  // gather the solution from all processors
441  data_.output_fields.set_time( this->time().step(), LimitSide::left);
442  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
444  data_.output_fields.output(this->time().step());
445  output_stream_->write_time_frame();
446 
447 // START_TIMER("MECH-balance");
448 // balance_->calculate_instant(subst_idx, ls->get_solution());
449 // balance_->output();
450 // END_TIMER("MECH-balance");
451 
452  END_TIMER("MECH-OUTPUT");
453 }
454 
455 
456 
458 {
459 // if (balance_->cumulative())
460 // {
461 // balance_->calculate_cumulative(subst_idx, ls->get_solution());
462 // }
463 }
464 
465 
466 
467 
468 
469 
470 
471 
473 {
474  START_TIMER("assemble_stiffness");
475  START_TIMER("assemble_volume_integrals");
476  assemble_volume_integrals<1>();
477  assemble_volume_integrals<2>();
478  assemble_volume_integrals<3>();
479  END_TIMER("assemble_volume_integrals");
480 
481  START_TIMER("assemble_fluxes_boundary");
482  assemble_fluxes_boundary<1>();
483  assemble_fluxes_boundary<2>();
484  assemble_fluxes_boundary<3>();
485  END_TIMER("assemble_fluxes_boundary");
486 
487  START_TIMER("assemble_fluxes_elem_side");
488  assemble_fluxes_element_side<1>();
489  assemble_fluxes_element_side<2>();
490  assemble_fluxes_element_side<3>();
491  END_TIMER("assemble_fluxes_elem_side");
492  END_TIMER("assemble_stiffness");
493 }
494 
495 
496 double lame_mu(double young, double poisson)
497 {
498  return young*0.5/(poisson+1.);
499 }
500 
501 
502 double lame_lambda(double young, double poisson)
503 {
504  return young*poisson/((poisson+1.)*(1.-2.*poisson));
505 }
506 
507 
508 template<unsigned int dim>
510 {
511  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
513  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
514  vector<int> dof_indices(ndofs);
515  vector<arma::vec3> velocity(qsize);
516  vector<double> young(qsize), poisson(qsize), csection(qsize);
517  PetscScalar local_matrix[ndofs*ndofs];
518  auto vec = fe_values.vector_view(0);
519 
520  // assemble integral over elements
521  for (auto cell : feo->dh()->own_range())
522  {
523  if (cell.dim() != dim) continue;
524  ElementAccessor<3> elm_acc = cell.elm();
525 
526  fe_values.reinit(elm_acc);
527  cell.get_dof_indices(dof_indices);
528 
529  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
530  data_.young_modulus.value_list(fe_values.point_list(), elm_acc, young);
531  data_.poisson_ratio.value_list(fe_values.point_list(), elm_acc, poisson);
532 
533  // assemble the local stiffness matrix
534  for (unsigned int i=0; i<ndofs; i++)
535  for (unsigned int j=0; j<ndofs; j++)
536  local_matrix[i*ndofs+j] = 0;
537 
538  for (unsigned int k=0; k<qsize; k++)
539  {
540  double mu = lame_mu(young[k], poisson[k]);
541  double lambda = lame_lambda(young[k], poisson[k]);
542  for (unsigned int i=0; i<ndofs; i++)
543  {
544  for (unsigned int j=0; j<ndofs; j++)
545  local_matrix[i*ndofs+j] += csection[k]*(
546  2*mu*arma::dot(vec.sym_grad(j,k), vec.sym_grad(i,k))
547  + lambda*vec.divergence(j,k)*vec.divergence(i,k)
548  )*fe_values.JxW(k);
549  }
550  }
551  ls->mat_set_values(ndofs, dof_indices.data(), ndofs, dof_indices.data(), local_matrix);
552  }
553 }
554 
555 
556 
558 {
559  START_TIMER("assemble_sources");
560 // balance_->start_source_assembly(subst_idx);
561  set_sources<1>();
562  set_sources<2>();
563  set_sources<3>();
564 // balance_->finish_source_assembly(subst_idx);
565  END_TIMER("assemble_sources");
566 }
567 
568 
569 template<unsigned int dim>
571 {
572  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
574  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
575  vector<arma::vec3> load(qsize);
576  vector<double> csection(qsize);
577  vector<int> dof_indices(ndofs);
578  PetscScalar local_rhs[ndofs];
579  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
580  auto vec = fe_values.vector_view(0);
581 
582  // assemble integral over elements
583  for (auto cell : feo->dh()->own_range())
584  {
585  if (cell.dim() != dim) continue;
586  ElementAccessor<3> elm_acc = cell.elm();
587 
588  fe_values.reinit(elm_acc);
589  cell.get_dof_indices(dof_indices);
590 
591  // assemble the local stiffness matrix
592  fill_n(local_rhs, ndofs, 0);
593  local_source_balance_vector.assign(ndofs, 0);
594  local_source_balance_rhs.assign(ndofs, 0);
595 
596  data_.cross_section.value_list(fe_values.point_list(), elm_acc, csection);
597  data_.load.value_list(fe_values.point_list(), elm_acc, load);
598 
599  // compute sources
600  for (unsigned int k=0; k<qsize; k++)
601  {
602  for (unsigned int i=0; i<ndofs; i++)
603  local_rhs[i] += arma::dot(load[k], vec.value(i,k))*csection[k]*fe_values.JxW(k);
604  }
605  ls->rhs_set_values(ndofs, &(dof_indices[0]), local_rhs);
606 
607 // for (unsigned int i=0; i<ndofs; i++)
608 // {
609 // for (unsigned int k=0; k<qsize; k++)
610 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_values[vec].value(i,k)*fe_values.JxW(k);
611 //
612 // local_source_balance_rhs[i] += local_rhs[i];
613 // }
614 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_vector);
615 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices, local_source_balance_rhs);
616  }
617 }
618 
619 
620 
621 
622 
623 
624 
625 template<unsigned int dim>
627 {
628  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
630  const unsigned int ndofs = feo->fe<dim>()->n_dofs();
631  vector<int> side_dof_indices(ndofs);
632  // PetscScalar local_matrix[ndofs*ndofs];
633 
634  // assemble boundary integral
635  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
636  {
637  Edge edg = mesh_->edge(feo->dh()->edge_index(iedg));
638  if (edg.n_sides() > 1) continue;
639  // check spatial dimension
640  if (edg.side(0)->dim() != dim-1) continue;
641  // skip edges lying not on the boundary
642  if (! edg.side(0)->is_boundary()) continue;
643 
644  SideIter side = edg.side(0);
645  ElementAccessor<3> cell = side->element();
646  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices);
647  fe_values_side.reinit(cell, side->side_idx());
648 // unsigned int bc_type = data_.bc_type.value(side->centre(), side->cond()->element_accessor());
649 
650 // for (unsigned int i=0; i<ndofs; i++)
651 // for (unsigned int j=0; j<ndofs; j++)
652 // local_matrix[i*ndofs+j] = 0;
653 //
654 // ls->mat_set_values(ndofs, side_dof_indices.data(), ndofs, side_dof_indices.data(), local_matrix);
655  }
656 }
657 
658 
660 {
661  arma::mat33 mt = m - m*arma::kron(n,n.t());
662  return mt;
663 }
664 
665 
666 template<unsigned int dim>
668 {
669  if (dim == 1) return;
670  FEValues<dim-1,3> fe_values_sub(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
672  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
674 
675  vector<FEValuesSpaceBase<3>*> fv_sb(2);
676  const unsigned int ndofs_side = feo->fe<dim>()->n_dofs(); // number of local dofs
677  const unsigned int ndofs_sub = feo->fe<dim-1>()->n_dofs();
678  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
679  vector<vector<int> > side_dof_indices(2);
680  vector<unsigned int> n_dofs = { ndofs_sub, ndofs_side };
681  vector<double> frac_sigma(qsize);
682  vector<double> csection_lower(qsize), csection_higher(qsize), young(qsize), poisson(qsize), alpha(qsize);
683  PetscScalar local_matrix[2][2][(ndofs_side)*(ndofs_side)];
684  PetscScalar local_rhs[2][ndofs_side];
685  auto vec_side = fe_values_side.vector_view(0);
686  auto vec_sub = fe_values_sub.vector_view(0);
687 
688  // index 0 = element with lower dimension,
689  // index 1 = side of element with higher dimension
690  side_dof_indices[0].resize(ndofs_sub);
691  side_dof_indices[1].resize(ndofs_side);
692  fv_sb[0] = &fe_values_sub;
693  fv_sb[1] = &fe_values_side;
694 
695  // assemble integral over sides
696  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
697  {
698  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
699  // skip neighbours of different dimension
700  if (nb->element()->dim() != dim-1) continue;
701 
702  ElementAccessor<3> cell_sub = nb->element();
703  feo->dh()->cell_accessor_from_element(cell_sub.idx()).get_dof_indices(side_dof_indices[0]);
704  fe_values_sub.reinit(cell_sub);
705 
706  ElementAccessor<3> cell = nb->side()->element();
707  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices[1]);
708  fe_values_side.reinit(cell, nb->side()->side_idx());
709 
710  // Element id's for testing if they belong to local partition.
711  bool own_element_id[2];
712  own_element_id[0] = feo->dh()->cell_accessor_from_element(cell_sub.idx()).is_own();
713  own_element_id[1] = feo->dh()->cell_accessor_from_element(cell.idx()).is_own();
714 
715  data_.cross_section.value_list(fe_values_sub.point_list(), cell_sub, csection_lower);
716  data_.cross_section.value_list(fe_values_sub.point_list(), cell, csection_higher);
717  data_.fracture_sigma.value_list(fe_values_sub.point_list(), cell_sub, frac_sigma);
718  data_.young_modulus.value_list(fe_values_sub.point_list(), cell_sub, young);
719  data_.poisson_ratio.value_list(fe_values_sub.point_list(), cell_sub, poisson);
720 
721  for (unsigned int n=0; n<2; ++n)
722  {
723  for (unsigned int i=0; i<ndofs_side; i++)
724  {
725  for (unsigned int m=0; m<2; ++m)
726  {
727  for (unsigned int j=0; j<ndofs_side; j++)
728  local_matrix[n][m][i*(ndofs_side)+j] = 0;
729  }
730  local_rhs[n][i] = 0;
731  }
732  }
733 
734  // set transmission conditions
735  for (unsigned int k=0; k<qsize; k++)
736  {
737  arma::vec3 nv = fe_values_side.normal_vector(k);
738  double mu = lame_mu(young[k], poisson[k]);
739  double lambda = lame_lambda(young[k], poisson[k]);
740 
741  for (int n=0; n<2; n++)
742  {
743  if (!own_element_id[n]) continue;
744 
745  for (unsigned int i=0; i<n_dofs[n]; i++)
746  {
747  arma::vec3 vi = (n==0)?arma::zeros(3):vec_side.value(i,k);
748  arma::vec3 vf = (n==1)?arma::zeros(3):vec_sub.value(i,k);
749  arma::mat33 gvft = (n==0)?vec_sub.grad(i,k):arma::zeros(3,3);
750 
751  for (int m=0; m<2; m++)
752  for (unsigned int j=0; j<n_dofs[m]; j++) {
753  arma::vec3 ui = (m==0)?arma::zeros(3):vec_side.value(j,k);
754  arma::vec3 uf = (m==1)?arma::zeros(3):vec_sub.value(j,k);
755  arma::mat33 guit = (m==1)?mat_t(vec_side.grad(j,k),nv):arma::zeros(3,3);
756  double divuit = (m==1)?arma::trace(guit):0;
757 
758  local_matrix[n][m][i*n_dofs[m] + j] +=
759  (
760  frac_sigma[k]*arma::dot(vf-vi,
761  2/csection_lower[k]*(mu*(uf-ui)+(mu+lambda)*(arma::dot(uf-ui,nv)*nv))
762  + mu*arma::trans(guit)*nv
763  + lambda*divuit*nv
764  )
765  - arma::dot(gvft, mu*arma::kron(nv,ui.t()) + lambda*arma::dot(ui,nv)*arma::eye(3,3))
766  )*fv_sb[0]->JxW(k);
767  }
768  }
769  }
770  }
771 
772  for (unsigned int n=0; n<2; ++n)
773  {
774  for (unsigned int m=0; m<2; ++m)
775  ls->mat_set_values(n_dofs[n], side_dof_indices[n].data(), n_dofs[m], side_dof_indices[m].data(), local_matrix[n][m]);
776 
777  ls->rhs_set_values(n_dofs[n], side_dof_indices[n].data(), local_rhs[n]);
778  }
779  }
780 }
781 
782 
783 
784 
785 
786 
787 
789 {
790  START_TIMER("assemble_bc");
791 // balance_->start_flux_assembly(subst_idx);
792  set_boundary_conditions<1>();
793  set_boundary_conditions<2>();
794  set_boundary_conditions<3>();
795 // balance_->finish_flux_assembly(subst_idx);
796  END_TIMER("assemble_bc");
797 }
798 
799 
800 
801 template<unsigned int dim>
803 {
804  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
806  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
807  vector<int> side_dof_indices(ndofs);
808  // unsigned int loc_b=0;
809  double local_rhs[ndofs];
810  vector<PetscScalar> local_flux_balance_vector(ndofs);
811  // PetscScalar local_flux_balance_rhs;
812  vector<arma::vec3> bc_values(qsize), bc_traction(qsize);
813  vector<double> csection(qsize);
814  auto vec = fe_values_side.vector_view(0);
815 
816  for (auto cell : feo->dh()->own_range())
817  {
818  ElementAccessor<3> elm = cell.elm();
819  if (elm->boundary_idx_ == nullptr) continue;
820 
821  for (unsigned int si=0; si<elm->n_sides(); ++si)
822  {
823  Edge edg = elm.side(si)->edge();
824  if (edg.n_sides() > 1) continue;
825  // skip edges lying not on the boundary
826  if (!edg.side(0)->is_boundary()) continue;
827 
828  if (edg.side(0)->dim() != dim-1)
829  {
830  // if (edg.side(0)->cond() != nullptr) ++loc_b;
831  continue;
832  }
833 
834  SideIter side = edg.side(0);
835  ElementAccessor<3> cell = side->element();
836  ElementAccessor<3> bc_cell = side->cond().element_accessor();
837 
838  unsigned int bc_type = data_.bc_type.value(side->centre(), bc_cell);
839 
840  fe_values_side.reinit(cell, side->side_idx());
841 
842  data_.cross_section.value_list(fe_values_side.point_list(), cell, csection);
843  // The b.c. data are fetched for all possible b.c. types since we allow
844  // different bc_type for each substance.
845  data_.bc_displacement.value_list(fe_values_side.point_list(), bc_cell, bc_values);
846  data_.bc_traction.value_list(fe_values_side.point_list(), bc_cell, bc_traction);
847 
848  feo->dh()->cell_accessor_from_element(cell.idx()).get_dof_indices(side_dof_indices);
849 
850  fill_n(local_rhs, ndofs, 0);
851  local_flux_balance_vector.assign(ndofs, 0);
852  // local_flux_balance_rhs = 0;
853 
854  if (bc_type == EqData::bc_type_displacement)
855  {
856  // We solve a local problem to determine which local dof
857  // to assign the boundary value to.
858  // TODO: Should be possibly optimized by inverting the matrix only once within reference element.
859  arma::mat d_mat(ndofs,ndofs);
860  arma::vec d_vec(ndofs);
861  for (unsigned int i=0; i<ndofs; i++)
862  {
863  d_vec(i) = 0;
864  for (unsigned int k=0; k<qsize; k++)
865  d_vec(i) += arma::dot(vec.value(i,k),bc_values[k])*fe_values_side.JxW(k);
866  for (unsigned int j=i; j<ndofs; j++)
867  {
868  d_mat(i,j) = 0;
869  for (unsigned int k=0; k<qsize; k++)
870  d_mat(i,j) += arma::dot(vec.value(i,k),vec.value(j,k))*fe_values_side.JxW(k);
871  d_mat(j,i) = d_mat(i,j);
872  }
873  }
874  arma::vec d_val = pinv(d_mat)*d_vec;
875 
876  for (unsigned int i=0; i<ndofs; i++)
877  if (norm(d_mat.row(i)) > 0)
878  ls->add_constraint(side_dof_indices[i], d_val(i));
879  }
880  else if (bc_type == EqData::bc_type_traction)
881  {
882  for (unsigned int k=0; k<qsize; k++)
883  {
884  for (unsigned int i=0; i<ndofs; i++)
885  local_rhs[i] += csection[k]*arma::dot(vec.value(i,k),bc_traction[k])*fe_values_side.JxW(k);
886  }
887  }
888  ls->rhs_set_values(ndofs, &(side_dof_indices[0]), local_rhs);
889 
890 
891 
892 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
893 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
894  // ++loc_b;
895  }
896  }
897 }
898 
899 
900 
901 
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
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
unsigned int n_sides() const
Returns number of sides aligned with the edge.
Definition: accessors.hh:298
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:48
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:82
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:359
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: elasticity.hh:131
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:413
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:37
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
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: accessors.cc:118
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:788
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: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.
~Elasticity()
Destructor.
Definition: elasticity.cc:298
SideIter side(const unsigned int loc_index)
double lame_mu(double young, double poisson)
Definition: elasticity.cc:496
const RegionDB & region_db() const
Definition: mesh.h:141
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:231
void update_solution() override
Computes the solution in one time instant.
Definition: elasticity.cc:376
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
Definition: elasticity.cc:659
Basic time management functionality for unsteady (and steady) solvers (class Equation).
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
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:123
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:123
bool is_boundary() const
Returns true for side on the boundary.
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:145
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:134
static auto region_id(Mesh &mesh) -> IndexField
#define xprintf(...)
Definition: system.hh:93
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:131
void set_sources()
Assembles the right hand side due to volume sources.
Definition: elasticity.cc:557
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:324
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:436
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:626
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:218
vector< Neighbour > vb_neighbours_
Definition: mesh.h:268
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:667
EquationOutput output_fields
Definition: elasticity.hh:136
Definition: system.hh:65
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
void initialize() override
Definition: elasticity.cc:256
EqData data_
Field data for model parameters.
Definition: elasticity.hh:262
double lame_lambda(double young, double poisson)
Definition: elasticity.cc:502
static const Input::Type::Selection & get_bc_type_selection()
Definition: elasticity.cc:143
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
Edge edge(uint edge_idx) const
Definition: mesh.cc:248
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
Definition: elasticity.cc:472
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
Boundary cond() const
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
Definition: elasticity.cc:509
virtual const Mat * get_matrix()
Definition: linsys.hh:187
void calculate_cumulative_balance()
Definition: elasticity.cc:457
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:180
Field< 3, FieldValue< 3 >::VectorFixed > output_field
Definition: elasticity.hh:134
SideIter side(const unsigned int i) const
Gets side iterator of the i -th side.
void output_vector_gather()
Definition: elasticity.cc:311
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()
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