Flow123d  JS_constraints-e651b99
assembly_elasticity.hh
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 assembly_elasticity.hh
15  * @brief
16  */
17 
18 #ifndef ASSEMBLY_ELASTICITY_HH_
19 #define ASSEMBLY_ELASTICITY_HH_
20 
23 #include "mechanics/elasticity.hh"
24 #include "fem/fe_p.hh"
25 #include "fem/patch_op_impl.hh"
27 #include "coupling/balance.hh"
28 #include "fem/element_cache_map.hh"
29 
30 
31 /**
32  * Auxiliary container class for Finite element and related objects of given dimension.
33  */
34 template <unsigned int dim, class TEqData>
36 {
37 public:
38  typedef typename TEqData::EqFields EqFields;
39  typedef TEqData EqData;
40 
41  static constexpr const char * name() { return "Elasticity_Stiffness_Assembly"; }
42 
43  /// Constructor.
45  : AssemblyBasePatch<dim>(eq_data->quad_order(), asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data), // quad_order = 1
49  JxW_( bulk_integral_->JxW() ),
50  JxW_side_( bdr_integral_->JxW() ),
51  JxW_join_( coupling_integral_->JxW() ),
52  normal_( bdr_integral_->normal_vector() ),
53  normal_join_( coupling_integral_->normal_vector() ),
54  deform_side_( bdr_integral_->vector_shape() ),
55  grad_deform_( bulk_integral_->grad_vector_shape() ),
56  sym_grad_deform_( bulk_integral_->vector_sym_grad() ),
57  deform_join_( coupling_integral_->vector_join_shape() ),
58  deform_join_grad_( coupling_integral_->gradient_vector_join_shape() ) {
59  this->used_fields_ += eq_fields_->cross_section;
60  this->used_fields_ += eq_fields_->lame_mu;
61  this->used_fields_ += eq_fields_->lame_lambda;
62  this->used_fields_ += eq_fields_->dirichlet_penalty;
63  this->used_fields_ += eq_fields_->bc_type;
64  this->used_fields_ += eq_fields_->fracture_sigma;
65  }
66 
67  /// Destructor.
69 
70  /// Initialize auxiliary vectors and other data members
71  void initialize() {
72  //this->balance_ = eq_data_->balance_;
73  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
74  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
75 
76  n_dofs_ = this->n_dofs();
77  n_dofs_high_ = this->n_dofs_high();
79  dof_indices_.resize(n_dofs_high_);
82  }
83 
84 
85  /// Assemble integral over element
86  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
87  {
88  if (cell.dim() != dim) return;
89 
90  // Fracture stiffness is assembled in dimjoin_integral.
91  if (cell.elm()->n_neighs_vb() > 0) return;
92 
94 
95  // assemble the local stiffness matrix
96  for (unsigned int i=0; i<n_dofs_; i++)
97  for (unsigned int j=0; j<n_dofs_; j++)
98  local_matrix_[i*n_dofs_+j] = 0;
99 
100  for (auto p : bulk_integral_->points(element_patch_idx) )
101  {
102  for (unsigned int i=0; i<n_dofs_; i++)
103  {
104  for (unsigned int j=0; j<n_dofs_; j++)
105  local_matrix_[i*n_dofs_+j] += eq_fields_->cross_section(p)*(
106  arma::dot(eq_fields_->stress_tensor(p,sym_grad_deform_.shape(j)(p)), sym_grad_deform_.shape(i)(p))
107  )*JxW_(p);
108  }
109  }
110  eq_data_->ls->mat_set_values(n_dofs_, dof_indices_.data(), n_dofs_, dof_indices_.data(), &(local_matrix_[0]));
111  }
112 
113  /// Assembles boundary integral.
114  inline void boundary_side_integral(DHCellSide cell_side)
115  {
116  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
117  if (!cell_side.cell().is_own()) return;
118 
119  Side side = cell_side.side();
120  const DHCellAccessor &dh_cell = cell_side.cell();
121  dh_cell.get_dof_indices(dof_indices_);
122 
123  for (unsigned int i=0; i<n_dofs_; i++)
124  for (unsigned int j=0; j<n_dofs_; j++)
125  local_matrix_[i*n_dofs_+j] = 0;
126 
127  auto p_side = *( bdr_integral_->points(cell_side).begin() );
128  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
129  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
130  double side_measure = cell_side.measure();
131  if (bc_type == EqFields::bc_type_displacement)
132  {
133  for (auto p : bdr_integral_->points(cell_side) ) {
134  for (unsigned int i=0; i<n_dofs_; i++)
135  for (unsigned int j=0; j<n_dofs_; j++)
136  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
137  arma::dot(deform_side_.shape(i)(p),deform_side_.shape(j)(p)) * JxW_side_(p);
138  }
139  }
140  else if (bc_type == EqFields::bc_type_displacement_normal)
141  {
142  for (auto p : bdr_integral_->points(cell_side) ) {
143  for (unsigned int i=0; i<n_dofs_; i++)
144  for (unsigned int j=0; j<n_dofs_; j++)
145  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
146  arma::dot(deform_side_.shape(i)(p), normal_(p)) *
147  arma::dot(deform_side_.shape(j)(p), normal_(p)) * JxW_side_(p);
148  }
149  }
150 
151  eq_data_->ls->mat_set_values(n_dofs_, dof_indices_.data(), n_dofs_, dof_indices_.data(), &(local_matrix_[0]));
152  }
153 
154 
155  /// Assembles between elements of different dimensions.
156  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
157  if (dim == 3) return;
158  ASSERT_EQ(cell_lower_dim.dim(), dim).error("Dimension of element mismatch!");
159 
160  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
161  for(unsigned int i=0; i<n_indices; ++i) {
163  }
164 
165  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
166  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
167  for(unsigned int i=0; i<n_indices; ++i) {
169  }
170 
171  // Element id's for testing if they belong to local partition.
172  bool own_element_id[2];
173  own_element_id[0] = cell_lower_dim.is_own();
174  own_element_id[1] = cell_higher_dim.is_own();
175 
176  unsigned int n_neighs = cell_lower_dim.elm()->n_neighs_vb();
177 
178  for (unsigned int i=0; i<n_dofs_ngh_[0]+n_dofs_ngh_[1]; i++)
179  for (unsigned int j=0; j<n_dofs_ngh_[0]+n_dofs_ngh_[1]; j++)
180  local_matrix_[i*(n_dofs_ngh_[0]+n_dofs_ngh_[1])+j] = 0;
181 
182  // set transmission conditions
183  for (auto p_high : coupling_integral_->points(neighb_side) )
184  {
185  auto p_low = p_high.lower_dim(cell_lower_dim);
186  arma::vec3 nv = normal_join_(p_high);
187 
188  for (uint i=0; i<deform_join_.n_dofs_both(); ++i) {
189  uint is_high_i = deform_join_.is_high_dim(i);
190  if (!own_element_id[is_high_i]) continue;
191  arma::vec3 diff_deform_i = deform_join_.shape(i)(p_low) - deform_join_.shape(i)(p_high);
192  arma::mat33 grad_deform_i = deform_join_grad_.shape(i)(p_low); // low dim element
193  arma::mat33 semi_grad_i = grad_deform_i + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_i,nv.t());
194  arma::mat33 semi_sym_grad_i = 0.5*(semi_grad_i + semi_grad_i.t());
195 
196  for (uint j=0; j<deform_join_.n_dofs_both(); ++j) {
197  arma::vec3 deform_j_high = deform_join_.shape(j)(p_high);
198  arma::vec3 diff_deform_j = deform_join_.shape(j)(p_low) - deform_j_high;
199  arma::mat33 grad_deform_j = deform_join_grad_.shape(j)(p_low); // low dim element
200  arma::mat33 semi_grad_j = grad_deform_j + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_j,nv.t());
201  arma::mat33 semi_sym_grad_j = 0.5*(semi_grad_j + semi_grad_j.t());
202 
203  local_matrix_[i * (n_dofs_ngh_[0]+n_dofs_ngh_[1]) + j] +=
204  (
205  eq_fields_->fracture_sigma(p_low)*eq_fields_->cross_section(p_low) / n_neighs
206  * arma::dot(semi_sym_grad_i, eq_fields_->stress_tensor(p_low,semi_sym_grad_j))
207 
208  // This term corrects the tangential part of the above product so that there is no
209  // dependence on fracture_sigma.
210  // TODO: Fracture_sigma should be possibly removed and replaced by anisotropic elasticity.
211  + (1-eq_fields_->fracture_sigma(p_low))*eq_fields_->cross_section(p_low) / n_neighs
212  * arma::dot(0.5*(grad_deform_i+grad_deform_i.t()), eq_fields_->stress_tensor(p_low,0.5*(grad_deform_j+grad_deform_j.t())))
213  )*JxW_join_(p_high);
214  }
215  }
216 
217  }
218 
219  eq_data_->ls->mat_set_values(n_dofs_ngh_[0]+n_dofs_ngh_[1], &(side_dof_indices_[0]), n_dofs_ngh_[0]+n_dofs_ngh_[1], &(side_dof_indices_[0]), &(local_matrix_[0]));
220  }
221 
222 
223 
224 private:
225  inline arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
226  {
227  arma::mat33 mt = m - m*arma::kron(n,n.t());
228  return mt;
229  }
230 
231 
232  /// Data objects shared with Elasticity
235 
236  /// Sub field set contains fields used in calculation.
238 
239  unsigned int n_dofs_; ///< Number of dofs
240  unsigned int n_dofs_high_; ///< Number of dofs (on lower dim element)
241  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
242 
243  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
244  vector<LongIdx> side_dof_indices_; ///< vector of DOF indices in neighbour calculation.
245  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
246 
247  std::shared_ptr<BulkIntegralAcc<dim>> bulk_integral_; ///< Bulk integral of assembly class
248  std::shared_ptr<BoundaryIntegralAcc<dim>> bdr_integral_; ///< Boundary integral of assembly class
249  std::shared_ptr<CouplingIntegralAcc<dim>> coupling_integral_; ///< Coupling integral of assembly class
250 
251  /// Following data members represent Element quantities and FE quantities
262 
263  template < template<IntDim...> class DimAssembly>
264  friend class GenericAssembly;
265 
266 };
267 
268 
269 template <unsigned int dim, class TEqData>
271 {
272 public:
273  typedef typename TEqData::EqFields EqFields;
274  typedef TEqData EqData;
275 
276  static constexpr const char * name() { return "Elasticity_Rhs_Assembly"; }
277 
278  /// Constructor.
280  : AssemblyBasePatch<dim>(eq_data->quad_order(), asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
284  JxW_( bulk_integral_->JxW() ),
285  JxW_side_( bdr_integral_->JxW() ),
286  JxW_join_( coupling_integral_->JxW() ),
287  normal_( bdr_integral_->normal_vector() ),
288  normal_join_( coupling_integral_->normal_vector() ),
289  deform_( bulk_integral_->vector_shape() ),
290  deform_side_( bdr_integral_->vector_shape() ),
291  grad_deform_( bulk_integral_->grad_vector_shape() ),
292  div_deform_( bulk_integral_->vector_divergence() ),
293  deform_join_( coupling_integral_->vector_join_shape() ) {
294  this->used_fields_ += eq_fields_->cross_section;
295  this->used_fields_ += eq_fields_->load;
296  this->used_fields_ += eq_fields_->potential_load;
297  this->used_fields_ += eq_fields_->ref_potential_load;
298  this->used_fields_ += eq_fields_->fracture_sigma;
299  this->used_fields_ += eq_fields_->dirichlet_penalty;
300  this->used_fields_ += eq_fields_->bc_type;
301  this->used_fields_ += eq_fields_->bc_displacement;
302  this->used_fields_ += eq_fields_->bc_traction;
303  this->used_fields_ += eq_fields_->bc_stress;
304  this->used_fields_ += eq_fields_->initial_stress;
305  }
306 
307  /// Destructor.
309 
310  /// Initialize auxiliary vectors and other data members
311  void initialize() {
312  //this->balance_ = eq_data_->balance_;
313  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
314  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
315 
316  n_dofs_ = this->n_dofs();
317  n_dofs_high_ = this->n_dofs_high();
319  dof_indices_.resize(n_dofs_high_);
320  side_dof_indices_.resize(n_dofs_ + n_dofs_high_);
321  local_rhs_.resize(2*n_dofs_high_);
322  }
323 
324 
325  /// Assemble integral over element
326  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
327  {
328  if (cell.dim() != dim) return;
329  if (!cell.is_own()) return;
330 
332 
333  // assemble the local stiffness matrix
334  fill_n(&(local_rhs_[0]), n_dofs_, 0);
335  //local_source_balance_vector.assign(n_dofs_, 0);
336  //local_source_balance_rhs.assign(n_dofs_, 0);
337 
338  // compute sources
339  for (auto p : bulk_integral_->points(element_patch_idx) )
340  {
341  for (unsigned int i=0; i<n_dofs_; i++)
342  local_rhs_[i] += (
343  arma::dot(eq_fields_->load(p), deform_.shape(i)(p))
344  -eq_fields_->potential_load(p)*div_deform_.shape(i)(p)
345  -arma::dot(eq_fields_->initial_stress(p), grad_deform_.shape(i)(p))
346  )*eq_fields_->cross_section(p)*JxW_(p);
347  }
348  eq_data_->ls->rhs_set_values(n_dofs_, dof_indices_.data(), &(local_rhs_[0]));
349 
350 // for (unsigned int i=0; i<n_dofs_; i++)
351 // {
352 // for (unsigned int k=0; k<qsize_; k++) // point range
353 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_vals_[vec_view_].value(i,k)*JxW_(k);
354 //
355 // local_source_balance_rhs[i] += local_rhs_[i];
356 // }
357 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_vector);
358 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_rhs);
359  }
360 
361  /// Assembles boundary integral.
362  inline void boundary_side_integral(DHCellSide cell_side)
363  {
364  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
365  if (!cell_side.cell().is_own()) return;
366 
367  const DHCellAccessor &dh_cell = cell_side.cell();
368  dh_cell.get_dof_indices(dof_indices_);
369 
370  auto p_side = *( bdr_integral_->points(cell_side).begin() );
371  auto p_bdr = p_side.point_bdr( cell_side.cond().element_accessor() );
372  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
373 
374  fill_n(&(local_rhs_[0]), n_dofs_, 0);
375  // local_flux_balance_vector.assign(n_dofs_, 0);
376  // local_flux_balance_rhs = 0;
377 
378  // addtion from initial stress
379  for (auto p : bdr_integral_->points(cell_side) )
380  {
381  for (unsigned int i=0; i<n_dofs_; i++)
382  local_rhs_[i] += eq_fields_->cross_section(p) *
383  arma::dot(( eq_fields_->initial_stress(p) * normal_(p)),
384  deform_side_.shape(i)(p)) *
385  JxW_side_(p);
386  }
387 
388  if (bc_type == EqFields::bc_type_displacement)
389  {
390  double side_measure = cell_side.measure();
391  for (auto p : bdr_integral_->points(cell_side) )
392  {
393  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
394  for (unsigned int i=0; i<n_dofs_; i++)
395  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
396  arma::dot(eq_fields_->bc_displacement(p_bdr), deform_side_.shape(i)(p)) *
397  JxW_side_(p);
398  }
399  }
400  else if (bc_type == EqFields::bc_type_displacement_normal)
401  {
402  double side_measure = cell_side.measure();
403  for (auto p : bdr_integral_->points(cell_side) )
404  {
405  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
406  for (unsigned int i=0; i<n_dofs_; i++)
407  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
408  arma::dot(eq_fields_->bc_displacement(p_bdr), normal_(p)) *
409  arma::dot(deform_side_.shape(i)(p), normal_(p)) *
410  JxW_side_(p);
411  }
412  }
413  else if (bc_type == EqFields::bc_type_traction)
414  {
415  for (auto p : bdr_integral_->points(cell_side) )
416  {
417  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
418  for (unsigned int i=0; i<n_dofs_; i++)
419  local_rhs_[i] += eq_fields_->cross_section(p) *
420  arma::dot(deform_side_.shape(i)(p), eq_fields_->bc_traction(p_bdr) + eq_fields_->ref_potential_load(p) * normal_(p)) *
421  JxW_side_(p);
422  }
423  }
424  else if (bc_type == EqFields::bc_type_stress)
425  {
426  for (auto p : bdr_integral_->points(cell_side) )
427  {
428  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
429  for (unsigned int i=0; i<n_dofs_; i++)
430  // stress is multiplied by inward normal to obtain traction
431  local_rhs_[i] += eq_fields_->cross_section(p) *
432  arma::dot(deform_side_.shape(i)(p), -eq_fields_->bc_stress(p_bdr)*normal_(p)
433  + eq_fields_->ref_potential_load(p) * normal_(p))
434  * JxW_side_(p);
435  }
436  }
437  eq_data_->ls->rhs_set_values(n_dofs_, dof_indices_.data(), &(local_rhs_[0]));
438 
439 
440 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
441 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
442  // ++loc_b;
443  }
444 
445 
446  /// Assembles between elements of different dimensions.
447  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
448  if (dim == 3) return;
449  ASSERT_EQ(cell_lower_dim.dim(), dim).error("Dimension of element mismatch!");
450 
451  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
452  for(unsigned int i=0; i<n_indices; ++i) {
454  }
455 
456  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
457  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
458  for(unsigned int i=0; i<n_indices; ++i) {
460  }
461 
462  // Element id's for testing if they belong to local partition.
463  bool own_element_id[2];
464  own_element_id[0] = cell_lower_dim.is_own();
465  own_element_id[1] = cell_higher_dim.is_own();
466 
467  for (unsigned int i=0; i<2*n_dofs_; i++)
468  local_rhs_[i] = 0;
469 
470  // set transmission conditions
471  for (auto p_high : coupling_integral_->points(neighb_side) )
472  {
473  auto p_low = p_high.lower_dim(cell_lower_dim);
474  arma::vec3 nv = normal_join_(p_high);
475 
476  for (uint i=0; i<deform_join_.n_dofs_both(); ++i) {
477  uint is_high_i = deform_join_.is_high_dim(i);
478  if (!own_element_id[is_high_i]) continue;
479 
480  arma::vec3 vi = deform_join_.shape(i)(p_high);
481  arma::vec3 vf = deform_join_.shape(i)(p_low);
482 
483  local_rhs_[i] -= eq_fields_->fracture_sigma(p_low) * eq_fields_->cross_section(p_high) *
484  arma::dot(vf-vi, eq_fields_->potential_load(p_high) * nv) * JxW_join_(p_high);
485  }
486  }
487 
488  eq_data_->ls->rhs_set_values(n_dofs_ngh_[0]+n_dofs_ngh_[1], side_dof_indices_.data(), &(local_rhs_[0]));
489  }
490 
491 
492 
493 private:
494  /// Data objects shared with Elasticity
497 
498  /// Sub field set contains fields used in calculation.
500 
501  unsigned int n_dofs_; ///< Number of dofs
502  unsigned int n_dofs_high_; ///< Number of dofs (on higher dim element)
503  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
504 
505  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
506  vector<LongIdx> side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
507  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for assemble methods
508 
509  std::shared_ptr<BulkIntegralAcc<dim>> bulk_integral_; ///< Bulk integral of assembly class
510  std::shared_ptr<BoundaryIntegralAcc<dim>> bdr_integral_; ///< Boundary integral of assembly class
511  std::shared_ptr<CouplingIntegralAcc<dim>> coupling_integral_; ///< Coupling integral of assembly class
512 
513  /// Following data members represent Element quantities and FE quantities
524 
525 
526  template < template<IntDim...> class DimAssembly>
527  friend class GenericAssembly;
528 
529 };
530 
531 template <unsigned int dim, class TEqData>
533 {
534 public:
535  typedef typename TEqData::EqFields EqFields;
536  typedef TEqData EqData;
537 
538  static constexpr const char * name() { return "Elasticity_OutpuFields_Assembly"; }
539 
540  /// Constructor.
542  : AssemblyBasePatch<dim>(eq_data->quad_order(), asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
545  normal_join_( coupling_integral_->normal_vector() ),
546  deform_join_( coupling_integral_->vector_shape() ),
547  grad_deform_( bulk_integral_->grad_vector_shape() ),
548  sym_grad_deform_( bulk_integral_->vector_sym_grad() ),
549  div_deform_( bulk_integral_->vector_divergence() ) {
550  this->used_fields_ += eq_fields_->cross_section;
551  this->used_fields_ += eq_fields_->lame_mu;
552  this->used_fields_ += eq_fields_->lame_lambda;
553  this->used_fields_ += eq_fields_->initial_stress;
554  }
555 
556  /// Destructor.
558 
559  /// Initialize auxiliary vectors and other data members
560  void initialize() {
561  //this->balance_ = eq_data_->balance_;
562 
563  n_dofs_ = this->n_dofs();
564  n_dofs_high_ = this->n_dofs_high();
565 
566  output_vec_ = eq_fields_->output_field_ptr->vec();
567  output_stress_vec_ = eq_fields_->output_stress_ptr->vec();
568  output_von_mises_stress_vec_ = eq_fields_->output_von_mises_stress_ptr->vec();
569  output_mean_stress_vec_ = eq_fields_->output_mean_stress_ptr->vec();
570  output_cross_sec_vec_ = eq_fields_->output_cross_section_ptr->vec();
571  output_div_vec_ = eq_fields_->output_div_ptr->vec();
572  }
573 
574 
575  /// Assemble integral over element
576  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
577  {
578  if (cell.dim() != dim) return;
579  if (!cell.is_own()) return;
580  DHCellAccessor cell_tensor = cell.cell_with_other_dh(eq_data_->dh_tensor_.get());
581  DHCellAccessor cell_scalar = cell.cell_with_other_dh(eq_data_->dh_scalar_.get());
582 
586 
587  auto p = *( bulk_integral_->points(element_patch_idx).begin() );
588 
589  arma::mat33 stress = eq_fields_->initial_stress(p);
590  double div = 0;
591  for (unsigned int i=0; i<n_dofs_; i++)
592  {
593  stress += eq_fields_->stress_tensor(p,sym_grad_deform_.shape(i)(p))
595  div += div_deform_.shape(i)(p)*output_vec_.get(dof_indices_[i]);
596  }
597 
598  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
599  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
600  double mean_stress = arma::trace(stress) / 3;
602 
603  for (unsigned int i=0; i<3; i++)
604  for (unsigned int j=0; j<3; j++)
605  output_stress_vec_.add( dof_indices_tensor_[i*3+j], stress(i,j) );
606  output_von_mises_stress_vec_.set( dof_indices_scalar_[0], von_mises_stress );
608 
609  output_cross_sec_vec_.add( dof_indices_scalar_[0], eq_fields_->cross_section(p) );
610  }
611 
612 
613  /// Assembles between elements of different dimensions.
614  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
615  if (dim == 3) return;
616  ASSERT_EQ(cell_lower_dim.dim(), dim).error("Dimension of element mismatch!");
617 
619  normal_stress_.zeros();
620 
621  DHCellAccessor cell_higher_dim = neighb_side.cell();
622  DHCellAccessor cell_tensor = cell_lower_dim.cell_with_other_dh(eq_data_->dh_tensor_.get());
623  DHCellAccessor cell_scalar = cell_lower_dim.cell_with_other_dh(eq_data_->dh_scalar_.get());
624 
625  dof_indices_ = cell_higher_dim.get_loc_dof_indices();
626  auto p_high = *( coupling_integral_->points(neighb_side).begin() );
627  auto p_low = p_high.lower_dim(cell_lower_dim);
628 
629  for (unsigned int i=0; i<n_dofs_high_; i++)
630  {
631  normal_displacement_ -= arma::dot(deform_join_.shape(i)(p_high)*output_vec_.get(dof_indices_[i]), normal_join_(p_high));
632  arma::mat33 grad = -arma::kron(deform_join_.shape(i)(p_high)*output_vec_.get(dof_indices_[i]), normal_join_(p_high).t()) / eq_fields_->cross_section(p_low);
633  normal_stress_ += eq_fields_->stress_tensor(p_low, 0.5*(grad+grad.t()));
634  }
635 
638  for (unsigned int i=0; i<3; i++)
639  for (unsigned int j=0; j<3; j++)
643  }
644 
645 
646 
647 private:
648  /// Data objects shared with Elasticity
651 
652  /// Sub field set contains fields used in calculation.
654 
655  unsigned int n_dofs_; ///< Number of dofs
656  unsigned int n_dofs_high_; ///< Number of dofs of higher dim element
657  LocDofVec dof_indices_; ///< Vector of local DOF indices of vector fields
658  LocDofVec dof_indices_scalar_; ///< Vector of local DOF indices of scalar fields
659  LocDofVec dof_indices_tensor_; ///< Vector of local DOF indices of tensor fields
660 
661  double normal_displacement_; ///< Holds constributions of normal displacement.
662  arma::mat33 normal_stress_; ///< Holds constributions of normal stress.
663 
664  std::shared_ptr<BulkIntegralAcc<dim>> bulk_integral_; ///< Bulk integral of assembly class
665  std::shared_ptr<CouplingIntegralAcc<dim>> coupling_integral_; ///< Coupling integral of assembly class
666 
667  /// Following data members represent Element quantities and FE quantities
673 
674  /// Data vectors of output fields (FieldFE).
681 
682  template < template<IntDim...> class DimAssembly>
683  friend class GenericAssembly;
684 
685 };
686 
687 
688 
689 /**
690  * Container class for assembly of constraint matrix for contact condition.
691  */
692 template <unsigned int dim, class TEqData>
694 {
695 public:
696  typedef typename TEqData::EqFields EqFields;
697  typedef TEqData EqData;
698 
699  static constexpr const char * name() { return "Elasticity_Constraint_Assembly"; }
700 
701  /// Constructor.
703  : AssemblyBasePatch<dim>(eq_data->quad_order(), asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
705  JxW_join_( coupling_integral_->JxW() ),
706  normal_join_( coupling_integral_->normal_vector() ),
707  deform_join_( coupling_integral_->vector_shape() ) {
708  this->used_fields_ += eq_fields_->cross_section;
709  this->used_fields_ += eq_fields_->cross_section_min;
710  }
711 
712  /// Destructor.
714 
715  /// Initialize auxiliary vectors and other data members
716  void initialize() {
717  n_dofs_ = this->n_dofs_high();
718  dof_indices_.resize(n_dofs_);
719  local_matrix_.resize(n_dofs_*n_dofs_);
720  }
721 
722 
723  /// Assembles between elements of different dimensions.
724  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
725  if (dim == 3) return;
726  if (!cell_lower_dim.is_own()) return;
727 
728  ASSERT_EQ(cell_lower_dim.dim(), dim).error("Dimension of element mismatch!");
729 
730  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
731  cell_higher_dim.get_dof_indices(dof_indices_);
732 
733  for (unsigned int i=0; i<n_dofs_; i++)
734  local_matrix_[i] = 0;
735 
736  // Assemble matrix and vector for contact conditions in the form B*x <= c,
737  // where B*x is the average jump of normal displacements and c is the average cross-section on element.
738  // Positive value means that the fracture closes.
739  double local_vector = 0;
740  for (auto p_high : coupling_integral_->points(neighb_side) )
741  {
742  auto p_low = p_high.lower_dim(cell_lower_dim);
743  arma::vec3 nv = normal_join_(p_high);
744 
745  local_vector += (eq_fields_->cross_section(p_low) - eq_fields_->cross_section_min(p_low))*JxW_join_(p_high) / cell_lower_dim.elm().measure() / cell_lower_dim.elm()->n_neighs_vb();
746 
747  for (unsigned int i=0; i<n_dofs_; i++)
748  {
749  local_matrix_[i] += eq_fields_->cross_section(p_high)*arma::dot(deform_join_.shape(i)(p_high), nv)*JxW_join_(p_high) / cell_lower_dim.elm().measure();
750  }
751  }
752 
753  int arow[1] = { eq_data_->constraint_idx[cell_lower_dim.elm_idx()] };
754  MatSetValues(eq_data_->constraint_matrix, 1, arow, n_dofs_, dof_indices_.data(), &(local_matrix_[0]), ADD_VALUES);
755  VecSetValue(eq_data_->constraint_vec, arow[0], local_vector, ADD_VALUES);
756  }
757 
758 
759 
760 private:
761 
762 
763  /// Data objects shared with Elasticity
766 
767  /// Sub field set contains fields used in calculation.
769 
770  unsigned int n_dofs_; ///< Number of dofs
771  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
772  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
773  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
774 
775  std::shared_ptr<CouplingIntegralAcc<dim>> coupling_integral_; ///< Coupling integral of assembly class
776 
777  /// Following data members represent Element quantities and FE quantities
781 
782 
783  template < template<IntDim...> class DimAssembly>
784  friend class GenericAssembly;
785 
786 };
787 
788 
789 #endif /* ASSEMBLY_ELASTICITY_HH_ */
790 
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
std::shared_ptr< CouplingIntegralAcc< dim > > create_coupling_integral(Quadrature *quad)
unsigned int n_dofs_high()
Return number of DOFs of higher dim element.
std::shared_ptr< BulkIntegralAcc< dim > > create_bulk_integral(Quadrature *quad)
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
Quadrature * quad_
Quadrature used in assembling methods.
unsigned int n_dofs()
Return number of DOFs.
std::shared_ptr< BoundaryIntegralAcc< dim > > create_boundary_integral(Quadrature *quad)
ElementAccessor< 3 > element_accessor()
static constexpr const char * name()
FieldSet used_fields_
Sub field set contains fields used in calculation.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
vector< vector< LongIdx > > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
unsigned int n_dofs_
Number of dofs.
EqFields * eq_fields_
Data objects shared with Elasticity.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
void initialize()
Initialize auxiliary vectors and other data members.
std::shared_ptr< CouplingIntegralAcc< dim > > coupling_integral_
Coupling integral of assembly class.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
FeQ< Scalar > JxW_join_
Following data members represent Element quantities and FE quantities.
ConstraintAssemblyElasticity(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
Cell accessor allow iterate over DOF handler cells.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
DHCellAccessor cell_with_other_dh(const DOFHandlerMultiDim *dh) const
Create new accessor with same local idx and given DOF handler. Actual and given DOF handler must be c...
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
Side side() const
Return Side of given cell and side_idx.
Boundary cond() const
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
double measure() const
ElementAccessor< 3 > element() const
double measure() const
Computes the measure of the element.
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:215
unsigned int n_neighs_vb() const
Return number of neighbours.
Definition: elements.h:65
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:102
Conforming Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:86
FeQ< ValueType > shape(unsigned int i_shape_fn_idx) const
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:160
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Generic class of assemblation.
LocDofVec dof_indices_tensor_
Vector of local DOF indices of tensor fields.
std::shared_ptr< BulkIntegralAcc< dim > > bulk_integral_
Bulk integral of assembly class.
static constexpr const char * name()
FieldSet used_fields_
Sub field set contains fields used in calculation.
double normal_displacement_
Holds constributions of normal displacement.
OutpuFieldsAssemblyElasticity(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
LocDofVec dof_indices_
Vector of local DOF indices of vector fields.
EqFields * eq_fields_
Data objects shared with Elasticity.
void initialize()
Initialize auxiliary vectors and other data members.
LocDofVec dof_indices_scalar_
Vector of local DOF indices of scalar fields.
std::shared_ptr< CouplingIntegralAcc< dim > > coupling_integral_
Coupling integral of assembly class.
arma::mat33 normal_stress_
Holds constributions of normal stress.
unsigned int n_dofs_
Number of dofs.
ElQ< Vector > normal_join_
Following data members represent Element quantities and FE quantities.
VectorMPI output_vec_
Data vectors of output fields (FieldFE).
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
unsigned int n_dofs_high_
Number of dofs of higher dim element.
TEqData::EqFields EqFields
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
FeQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
FeQArray< Tensor > grad_deform_
std::shared_ptr< BoundaryIntegralAcc< dim > > bdr_integral_
Boundary integral of assembly class.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
RhsAssemblyElasticity(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
FeQArray< Vector > deform_
EqFields * eq_fields_
Data objects shared with Elasticity.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
FeQArray< Scalar > div_deform_
void initialize()
Initialize auxiliary vectors and other data members.
std::shared_ptr< CouplingIntegralAcc< dim > > coupling_integral_
Coupling integral of assembly class.
std::shared_ptr< BulkIntegralAcc< dim > > bulk_integral_
Bulk integral of assembly class.
FieldSet used_fields_
Sub field set contains fields used in calculation.
vector< LongIdx > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
FeQJoin< Vector > deform_join_
static constexpr const char * name()
vector< PetscScalar > local_rhs_
Auxiliary vector for assemble methods.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
unsigned int n_dofs_high_
Number of dofs (on higher dim element)
FeQArray< Vector > deform_side_
unsigned int n_dofs_
Number of dofs.
Boundary cond() const
FieldSet used_fields_
Sub field set contains fields used in calculation.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
static constexpr const char * name()
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
vector< LongIdx > side_dof_indices_
vector of DOF indices in neighbour calculation.
void initialize()
Initialize auxiliary vectors and other data members.
std::shared_ptr< BulkIntegralAcc< dim > > bulk_integral_
Bulk integral of assembly class.
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
unsigned int n_dofs_
Number of dofs.
StiffnessAssemblyElasticity(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
vector< LongIdx > dof_indices_
Vector of global DOF indices.
unsigned int n_dofs_high_
Number of dofs (on lower dim element)
EqFields * eq_fields_
Data objects shared with Elasticity.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
std::shared_ptr< BoundaryIntegralAcc< dim > > bdr_integral_
Boundary integral of assembly class.
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
FeQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
std::shared_ptr< CouplingIntegralAcc< dim > > coupling_integral_
Coupling integral of assembly class.
double get(unsigned int pos) const
Return value on given position.
Definition: vector_mpi.hh:105
void set(unsigned int pos, double val)
Set value on given position.
Definition: vector_mpi.hh:111
void add(unsigned int pos, double val)
Add value to item on given position.
Definition: vector_mpi.hh:125
FEM for linear elasticity.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
@ FEVector
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
unsigned int uint
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Store finite element reinit functions.
Definitions of particular quadrature rules on simplices.
Holds common data shared between GenericAssemblz and Assembly<dim> classes.