Flow123d  DF_patch_fe_mechanics-8d0e131
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/fe_values.hh"
27 #include "coupling/balance.hh"
29 
30 
31 /**
32  * Auxiliary container class for Finite element and related objects of given dimension.
33  */
34 template <unsigned int dim>
36 {
37 public:
38  typedef typename Elasticity::EqFields EqFields;
39  typedef typename Elasticity::EqData EqData;
40 
41  static constexpr const char * name() { return "StiffnessAssemblyElasticity"; }
42 
43  /// Constructor.
44  StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
45  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data), // quad_order = 1
46  JxW_( this->bulk_values().JxW() ),
47  JxW_side_( this->side_values().JxW() ),
48  normal_( this->side_values().normal_vector() ),
49  deform_side_( this->side_values().vector_shape() ),
50  gras_deform_( this->bulk_values().grad_vector_shape() ),
51  sym_grad_deform_( this->bulk_values().vector_sym_grad() ),
52  div_deform_( this->bulk_values().vector_divergence() ),
53  deform_join_( Range< JoinShapeAccessor<Vector> >( this->join_values().vector_join_shape() ) ),
54  deform_join_grad_( Range< JoinShapeAccessor<Tensor> >( this->join_values().gradient_vector_join_shape() ) ) {
62  }
63 
64  /// Destructor.
66 
67  /// Initialize auxiliary vectors and other data members
68  void initialize(ElementCacheMap *element_cache_map) {
69  //this->balance_ = eq_data_->balance_;
70  this->element_cache_map_ = element_cache_map;
71 
72  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
73  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
74  this->fe_values_->template initialize<dim>(*this->quad_);
75  this->fe_values_->template initialize<dim>(*this->quad_low_);
76 
77  n_dofs_ = this->n_dofs();
78  n_dofs_sub_ = fe_low->n_dofs();
80  dof_indices_.resize(n_dofs_);
81  side_dof_indices_.resize(2*n_dofs_);
82  local_matrix_.resize(4*n_dofs_*n_dofs_);
83  }
84 
85 
86  /// Assemble integral over element
87  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
88  {
89  if (cell.dim() != dim) return;
90 
91  // Fracture stiffness is assembled in dimjoin_integral.
92  if (cell.elm()->n_neighs_vb() > 0) return;
93 
95 
96  // assemble the local stiffness matrix
97  for (unsigned int i=0; i<n_dofs_; i++)
98  for (unsigned int j=0; j<n_dofs_; j++)
99  local_matrix_[i*n_dofs_+j] = 0;
100 
101  for (auto p : this->bulk_points(element_patch_idx) )
102  {
103  for (unsigned int i=0; i<n_dofs_; i++)
104  {
105  for (unsigned int j=0; j<n_dofs_; j++)
108  )*JxW_(p);
109  }
110  }
112  }
113 
114  /// Assembles boundary integral.
115  inline void boundary_side_integral(DHCellSide cell_side)
116  {
117  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
118  if (!cell_side.cell().is_own()) return;
119 
120  Side side = cell_side.side();
121  const DHCellAccessor &dh_cell = cell_side.cell();
122  dh_cell.get_dof_indices(dof_indices_);
123 
124  for (unsigned int i=0; i<n_dofs_; i++)
125  for (unsigned int j=0; j<n_dofs_; j++)
126  local_matrix_[i*n_dofs_+j] = 0;
127 
128  auto p_side = *( this->boundary_points(cell_side).begin() );
129  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
130  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
131  double side_measure = cell_side.measure();
132  if (bc_type == EqFields::bc_type_displacement)
133  {
134  for (auto p : this->boundary_points(cell_side) ) {
135  for (unsigned int i=0; i<n_dofs_; i++)
136  for (unsigned int j=0; j<n_dofs_; j++)
137  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
138  arma::dot(deform_side_(i,p),deform_side_(j,p)) * JxW_side_(p);
139  }
140  }
141  else if (bc_type == EqFields::bc_type_displacement_normal)
142  {
143  for (auto p : this->boundary_points(cell_side) ) {
144  for (unsigned int i=0; i<n_dofs_; i++)
145  for (unsigned int j=0; j<n_dofs_; j++)
146  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
147  arma::dot(deform_side_(i,p), normal_(p)) *
148  arma::dot(deform_side_(j,p), normal_(p)) * JxW_side_(p);
149  }
150  }
151 
153  }
154 
155 
156  /// Assembles between elements of different dimensions.
157  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
158  if (dim == 1) return;
159  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
160 
161  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
162  for(unsigned int i=0; i<n_indices; ++i) {
164  }
165 
166  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
167  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
168  for(unsigned int i=0; i<n_indices; ++i) {
170  }
171 
172  // Element id's for testing if they belong to local partition.
173  bool own_element_id[2];
174  own_element_id[0] = cell_lower_dim.is_own();
175  own_element_id[1] = cell_higher_dim.is_own();
176 
177  unsigned int n_neighs = cell_lower_dim.elm()->n_neighs_vb();
178 
179  for (unsigned int i=0; i<n_dofs_ngh_[0]+n_dofs_ngh_[1]; i++)
180  for (unsigned int j=0; j<n_dofs_ngh_[0]+n_dofs_ngh_[1]; j++)
181  local_matrix_[i*(n_dofs_ngh_[0]+n_dofs_ngh_[1])+j] = 0;
182 
183  // set transmission conditions
184  for (auto p_high : this->coupling_points(neighb_side) )
185  {
186  auto p_low = p_high.lower_dim(cell_lower_dim);
187  arma::vec3 nv = normal_(p_high);
188 
189  auto deform_shape_i = deform_join_.begin();
190  auto deform_grad_i = deform_join_grad_.begin();
191  for( ; deform_shape_i != deform_join_.end() && deform_grad_i != deform_join_grad_.end(); ++deform_shape_i, ++deform_grad_i) {
192  uint is_high_i = deform_shape_i->is_high_dim();
193  if (!own_element_id[is_high_i]) continue;
194  uint i_mat_idx = deform_shape_i->join_idx();
195  arma::vec3 diff_deform_i = (*deform_shape_i)(p_low) - (*deform_shape_i)(p_high);
196  arma::mat33 grad_deform_i = (*deform_grad_i)(p_low); // low dim element
197  arma::mat33 semi_grad_i = grad_deform_i + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_i,nv.t());
198  arma::mat33 semi_sym_grad_i = 0.5*(semi_grad_i + semi_grad_i.t());
199 
200  auto deform_shape_j = deform_join_.begin();
201  auto deform_grad_j = deform_join_grad_.begin();
202  for( ; deform_shape_j != deform_join_.end() && deform_grad_j != deform_join_grad_.end(); ++deform_shape_j, ++deform_grad_j) {
203  uint j_mat_idx = deform_shape_j->join_idx();
204  arma::vec3 deform_j_high = (*deform_shape_j)(p_high);
205  arma::vec3 diff_deform_j = (*deform_shape_j)(p_low) - (*deform_shape_j)(p_high);
206  arma::mat33 grad_deform_j = (*deform_grad_j)(p_low); // low dim element
207  arma::mat33 semi_grad_j = grad_deform_j + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_j,nv.t());
208  arma::mat33 semi_sym_grad_j = 0.5*(semi_grad_j + semi_grad_j.t());
209 
210  local_matrix_[i_mat_idx * (n_dofs_ngh_[0]+n_dofs_ngh_[1]) + j_mat_idx] +=
211  (
212  eq_fields_->fracture_sigma(p_low)*eq_fields_->cross_section(p_low) / n_neighs
213  * arma::dot(semi_sym_grad_i, eq_fields_->stress_tensor(p_low,semi_sym_grad_j))
214 
215  // This term corrects the tangential part of the above product so that there is no
216  // dependence on fracture_sigma.
217  // TODO: Fracture_sigma should be possibly removed and replaced by anisotropic elasticity.
218  + (1-eq_fields_->fracture_sigma(p_low))*eq_fields_->cross_section(p_low) / n_neighs
219  * 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())))
220  )*JxW_side_(p_high);
221  }
222  }
223  }
224 
226  }
227 
228 
229 
230 private:
231  inline arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
232  {
233  arma::mat33 mt = m - m*arma::kron(n,n.t());
234  return mt;
235  }
236 
237 
238  /// Data objects shared with Elasticity
241 
242  /// Sub field set contains fields used in calculation.
244 
245  unsigned int n_dofs_; ///< Number of dofs
246  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
247  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
248 
249  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
250  vector<LongIdx> side_dof_indices_; ///< vector of DOF indices in neighbour calculation.
251  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
252 
253  /// Following data members represent Element quantities and FE quantities
263 
264  template < template<IntDim...> class DimAssembly>
265  friend class GenericAssembly;
266 
267 };
268 
269 
270 template <unsigned int dim>
272 {
273 public:
274  typedef typename Elasticity::EqFields EqFields;
275  typedef typename Elasticity::EqData EqData;
276 
277  static constexpr const char * name() { return "RhsAssemblyElasticity"; }
278 
279  /// Constructor.
280  RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
281  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
282  JxW_( this->bulk_values().JxW() ),
283  JxW_side_( this->side_values().JxW() ),
284  normal_( this->side_values().normal_vector() ),
285  deform_( this->bulk_values().vector_shape() ),
286  deform_side_( this->side_values().vector_shape() ),
287  gras_deform_( this->bulk_values().grad_vector_shape() ),
288  div_deform_( this->bulk_values().vector_divergence() ),
289  deform_join_( Range< JoinShapeAccessor<Vector> >( this->join_values().vector_join_shape() ) ) {
292  this->used_fields_ += eq_fields_->load;
297  this->used_fields_ += eq_fields_->bc_type;
302  }
303 
304  /// Destructor.
306 
307  /// Initialize auxiliary vectors and other data members
308  void initialize(ElementCacheMap *element_cache_map) {
309  //this->balance_ = eq_data_->balance_;
310  this->element_cache_map_ = element_cache_map;
311 
312  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
313  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
314  this->fe_values_->template initialize<dim>(*this->quad_);
315  this->fe_values_->template initialize<dim>(*this->quad_low_);
316 
317  n_dofs_ = this->n_dofs();
318  n_dofs_sub_ = fe_low->n_dofs();
320  dof_indices_.resize(n_dofs_);
322  local_rhs_.resize(2*n_dofs_);
323  }
324 
325 
326  /// Assemble integral over element
327  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
328  {
329  if (cell.dim() != dim) return;
330  if (!cell.is_own()) return;
331 
333 
334  // assemble the local stiffness matrix
335  fill_n(&(local_rhs_[0]), n_dofs_, 0);
336  //local_source_balance_vector.assign(n_dofs_, 0);
337  //local_source_balance_rhs.assign(n_dofs_, 0);
338 
339  // compute sources
340  for (auto p : this->bulk_points(element_patch_idx) )
341  {
342  for (unsigned int i=0; i<n_dofs_; i++)
343  local_rhs_[i] += (
344  arma::dot(eq_fields_->load(p), deform_(i,p))
346  -arma::dot(eq_fields_->initial_stress(p), gras_deform_(i,p))
347  )*eq_fields_->cross_section(p)*JxW_(p);
348  }
350 
351 // for (unsigned int i=0; i<n_dofs_; i++)
352 // {
353 // for (unsigned int k=0; k<qsize_; k++) // point range
354 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_vals_[vec_view_].value(i,k)*JxW_(k);
355 //
356 // local_source_balance_rhs[i] += local_rhs_[i];
357 // }
358 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_vector);
359 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_rhs);
360  }
361 
362  /// Assembles boundary integral.
363  inline void boundary_side_integral(DHCellSide cell_side)
364  {
365  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
366  if (!cell_side.cell().is_own()) return;
367 
368  const DHCellAccessor &dh_cell = cell_side.cell();
369  dh_cell.get_dof_indices(dof_indices_);
370 
371  auto p_side = *( this->boundary_points(cell_side).begin() );
372  auto p_bdr = p_side.point_bdr( cell_side.cond().element_accessor() );
373  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
374 
375  fill_n(&(local_rhs_[0]), n_dofs_, 0);
376  // local_flux_balance_vector.assign(n_dofs_, 0);
377  // local_flux_balance_rhs = 0;
378 
379  // addtion from initial stress
380  for (auto p : this->boundary_points(cell_side) )
381  {
382  for (unsigned int i=0; i<n_dofs_; i++)
384  arma::dot(( eq_fields_->initial_stress(p) * normal_(p)),
385  deform_side_(i,p)) *
386  JxW_side_(p);
387  }
388 
389  if (bc_type == EqFields::bc_type_displacement)
390  {
391  double side_measure = cell_side.measure();
392  for (auto p : this->boundary_points(cell_side) )
393  {
394  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
395  for (unsigned int i=0; i<n_dofs_; i++)
396  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
397  arma::dot(eq_fields_->bc_displacement(p_bdr), deform_side_(i,p)) *
398  JxW_side_(p);
399  }
400  }
401  else if (bc_type == EqFields::bc_type_displacement_normal)
402  {
403  double side_measure = cell_side.measure();
404  for (auto p : this->boundary_points(cell_side) )
405  {
406  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
407  for (unsigned int i=0; i<n_dofs_; i++)
408  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
409  arma::dot(eq_fields_->bc_displacement(p_bdr), normal_(p)) *
410  arma::dot(deform_side_(i,p), normal_(p)) *
411  JxW_side_(p);
412  }
413  }
414  else if (bc_type == EqFields::bc_type_traction)
415  {
416  for (auto p : this->boundary_points(cell_side) )
417  {
418  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
419  for (unsigned int i=0; i<n_dofs_; i++)
421  arma::dot(deform_side_(i,p), eq_fields_->bc_traction(p_bdr) + eq_fields_->ref_potential_load(p) * normal_(p)) *
422  JxW_side_(p);
423  }
424  }
425  else if (bc_type == EqFields::bc_type_stress)
426  {
427  for (auto p : this->boundary_points(cell_side) )
428  {
429  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
430  for (unsigned int i=0; i<n_dofs_; i++)
431  // stress is multiplied by inward normal to obtain traction
433  arma::dot(deform_side_(i,p), -eq_fields_->bc_stress(p_bdr)*normal_(p)
435  * JxW_side_(p);
436  }
437  }
439 
440 
441 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
442 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
443  // ++loc_b;
444  }
445 
446 
447  /// Assembles between elements of different dimensions.
448  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
449  if (dim == 1) return;
450  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
451 
452  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
453  for(unsigned int i=0; i<n_indices; ++i) {
455  }
456 
457  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
458  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
459  for(unsigned int i=0; i<n_indices; ++i) {
461  }
462 
463  // Element id's for testing if they belong to local partition.
464  bool own_element_id[2];
465  own_element_id[0] = cell_lower_dim.is_own();
466  own_element_id[1] = cell_higher_dim.is_own();
467 
468  for (unsigned int i=0; i<2*n_dofs_; i++)
469  local_rhs_[i] = 0;
470 
471  // set transmission conditions
472  for (auto p_high : this->coupling_points(neighb_side) )
473  {
474  auto p_low = p_high.lower_dim(cell_lower_dim);
475  arma::vec3 nv = normal_(p_high);
476 
477  for( auto join_shape_i : deform_join_) {
478  uint is_high_i = join_shape_i.is_high_dim();
479  if (!own_element_id[is_high_i]) continue;
480 
481  arma::vec3 vi = join_shape_i(p_high);
482  arma::vec3 vf = join_shape_i(p_low);
483 
484  local_rhs_[join_shape_i.join_idx()] -= eq_fields_->fracture_sigma(p_low) * eq_fields_->cross_section(p_high) *
485  arma::dot(vf-vi, eq_fields_->potential_load(p_high) * nv) * JxW_side_(p_high);
486  }
487  }
488 
490  }
491 
492 
493 
494 private:
495  /// Data objects shared with Elasticity
498 
499  /// Sub field set contains fields used in calculation.
501 
502  unsigned int n_dofs_; ///< Number of dofs
503  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
504  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
505 
506  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
507  vector<LongIdx> side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
508  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for assemble methods
509 
510  /// Following data members represent Element quantities and FE quantities
519 
520 
521  template < template<IntDim...> class DimAssembly>
522  friend class GenericAssembly;
523 
524 };
525 
526 template <unsigned int dim>
528 {
529 public:
530  typedef typename Elasticity::EqFields EqFields;
532 
533  static constexpr const char * name() { return "OutpuFieldsAssemblyElasticity"; }
534 
535  /// Constructor.
537  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
538  normal_( this->side_values().normal_vector() ),
539  deform_side_( this->side_values().vector_shape() ),
540  gras_deform_( this->bulk_values().grad_vector_shape() ),
541  sym_grad_deform_( this->bulk_values().vector_sym_grad() ),
542  div_deform_( this->bulk_values().vector_divergence() ) {
545  this->used_fields_ += eq_fields_->lame_mu;
548  }
549 
550  /// Destructor.
552 
553  /// Initialize auxiliary vectors and other data members
554  void initialize(ElementCacheMap *element_cache_map) {
555  //this->balance_ = eq_data_->balance_;
556  this->element_cache_map_ = element_cache_map;
557 
558  this->fe_values_->template initialize<dim>(*this->quad_);
559  this->fe_values_->template initialize<dim>(*this->quad_low_);
560 
561  n_dofs_ = this->n_dofs();
562 
569  }
570 
571 
572  /// Assemble integral over element
573  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
574  {
575  if (cell.dim() != dim) return;
576  if (!cell.is_own()) return;
577  DHCellAccessor cell_tensor = cell.cell_with_other_dh(eq_data_->dh_tensor_.get());
578  DHCellAccessor cell_scalar = cell.cell_with_other_dh(eq_data_->dh_scalar_.get());
579 
583 
584  auto p = *( this->bulk_points(element_patch_idx).begin() );
585 
587  double div = 0;
588  for (unsigned int i=0; i<n_dofs_; i++)
589  {
590  stress += eq_fields_->stress_tensor(p,sym_grad_deform_(i,p))
592  div += div_deform_(i,p)*output_vec_.get(dof_indices_[i]);
593  }
594 
595  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
596  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
597  double mean_stress = arma::trace(stress) / 3;
599 
600  for (unsigned int i=0; i<3; i++)
601  for (unsigned int j=0; j<3; j++)
602  output_stress_vec_.add( dof_indices_tensor_[i*3+j], stress(i,j) );
603  output_von_mises_stress_vec_.set( dof_indices_scalar_[0], von_mises_stress );
605 
607  }
608 
609 
610  /// Assembles between elements of different dimensions.
611  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
612  if (dim == 1) return;
613  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
614 
616  normal_stress_.zeros();
617 
618  DHCellAccessor cell_higher_dim = neighb_side.cell();
619  DHCellAccessor cell_tensor = cell_lower_dim.cell_with_other_dh(eq_data_->dh_tensor_.get());
620  DHCellAccessor cell_scalar = cell_lower_dim.cell_with_other_dh(eq_data_->dh_scalar_.get());
621 
622  dof_indices_ = cell_higher_dim.get_loc_dof_indices();
623  auto p_high = *( this->coupling_points(neighb_side).begin() );
624  auto p_low = p_high.lower_dim(cell_lower_dim);
625 
626  for (unsigned int i=0; i<n_dofs_; i++)
627  {
628  normal_displacement_ -= arma::dot(deform_side_(i,p_high)*output_vec_.get(dof_indices_[i]), normal_(p_high));
629  arma::mat33 grad = -arma::kron(deform_side_(i,p_high)*output_vec_.get(dof_indices_[i]), normal_(p_high).t()) / eq_fields_->cross_section(p_low);
630  normal_stress_ += eq_fields_->stress_tensor(p_low, 0.5*(grad+grad.t()));
631  }
632 
635  for (unsigned int i=0; i<3; i++)
636  for (unsigned int j=0; j<3; j++)
640  }
641 
642 
643 
644 private:
645  /// Data objects shared with Elasticity
648 
649  /// Sub field set contains fields used in calculation.
651 
652  unsigned int n_dofs_; ///< Number of dofs
653  LocDofVec dof_indices_; ///< Vector of local DOF indices of vector fields
654  LocDofVec dof_indices_scalar_; ///< Vector of local DOF indices of scalar fields
655  LocDofVec dof_indices_tensor_; ///< Vector of local DOF indices of tensor fields
656 
657  double normal_displacement_; ///< Holds constributions of normal displacement.
658  arma::mat33 normal_stress_; ///< Holds constributions of normal stress.
659 
660  /// Following data members represent Element quantities and FE quantities
666 
667  /// Data vectors of output fields (FieldFE).
674 
675  template < template<IntDim...> class DimAssembly>
676  friend class GenericAssembly;
677 
678 };
679 
680 
681 
682 /**
683  * Container class for assembly of constraint matrix for contact condition.
684  */
685 template <unsigned int dim>
687 {
688 public:
689  typedef typename Elasticity::EqFields EqFields;
690  typedef typename Elasticity::EqData EqData;
691 
692  static constexpr const char * name() { return "ConstraintAssemblyElasticity"; }
693 
694  /// Constructor.
696  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
697  JxW_side_( this->side_values().JxW() ),
698  normal_( this->side_values().normal_vector() ),
699  deform_side_( this->side_values().vector_shape() ) {
703  }
704 
705  /// Destructor.
707 
708  /// Initialize auxiliary vectors and other data members
709  void initialize(ElementCacheMap *element_cache_map) {
710  this->element_cache_map_ = element_cache_map;
711 
712  this->fe_values_->template initialize<dim>(*this->quad_);
713  this->fe_values_->template initialize<dim>(*this->quad_low_);
714 
715  n_dofs_ = this->n_dofs();
716  dof_indices_.resize(n_dofs_);
717  local_matrix_.resize(n_dofs_*n_dofs_);
718  }
719 
720 
721  /// Assembles between elements of different dimensions.
722  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
723  if (dim == 1) return;
724  if (!cell_lower_dim.is_own()) return;
725 
726  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
727 
728  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
729  cell_higher_dim.get_dof_indices(dof_indices_);
730 
731  for (unsigned int i=0; i<n_dofs_; i++)
732  local_matrix_[i] = 0;
733 
734  // Assemble matrix and vector for contact conditions in the form B*x <= c,
735  // where B*x is the average jump of normal displacements and c is the average cross-section on element.
736  // Positive value means that the fracture closes.
737  double local_vector = 0;
738  for (auto p_high : this->coupling_points(neighb_side) )
739  {
740  auto p_low = p_high.lower_dim(cell_lower_dim);
741  arma::vec3 nv = normal_(p_high);
742 
743  local_vector += (eq_fields_->cross_section(p_low) - eq_fields_->cross_section_min(p_low))*JxW_side_(p_high) / cell_lower_dim.elm().measure() / cell_lower_dim.elm()->n_neighs_vb();
744 
745  for (unsigned int i=0; i<n_dofs_; i++)
746  {
747  local_matrix_[i] += eq_fields_->cross_section(p_high)*arma::dot(deform_side_(i,p_high), nv)*JxW_side_(p_high) / cell_lower_dim.elm().measure();
748  }
749  }
750 
751  int arow[1] = { eq_data_->constraint_idx[cell_lower_dim.elm_idx()] };
752  MatSetValues(eq_data_->constraint_matrix, 1, arow, n_dofs_, dof_indices_.data(), &(local_matrix_[0]), ADD_VALUES);
753  VecSetValue(eq_data_->constraint_vec, arow[0], local_vector, ADD_VALUES);
754  }
755 
756 
757 
758 private:
759 
760 
761  /// Data objects shared with Elasticity
764 
765  /// Sub field set contains fields used in calculation.
767 
768  unsigned int n_dofs_; ///< Number of dofs
769  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
770  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
771  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
772 
773  /// Following data members represent Element quantities and FE quantities
777 
778 
779  template < template<IntDim...> class DimAssembly>
780  friend class GenericAssembly;
781 
782 };
783 
784 
785 #endif /* ASSEMBLY_ELASTICITY_HH_ */
786 
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
JoinValues< dim > join_values()
Return JoinValues object.
PatchFEValues< 3 > * fe_values_
Common FEValues object over all dimensions.
BulkValues< dim > bulk_values()
Return BulkValues object.
unsigned int n_dofs()
Return BulkValues object.
SideValues< dim > side_values()
Return SideValues object.
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
Quadrature * quad_
Quadrature used in assembling methods.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
int active_integrals_
Holds mask of active integrals.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
ElementAccessor< 3 > element_accessor()
EqFields * eq_fields_
Data objects shared with Elasticity.
ConstraintAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
static constexpr const char * name()
vector< vector< LongIdx > > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
unsigned int n_dofs_
Number of dofs.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
FieldSet used_fields_
Sub field set contains fields used in calculation.
ElQ< Scalar > JxW_side_
Following data members represent Element quantities and FE quantities.
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
Data of equation parameters.
Definition: elasticity.hh:114
LinSys * ls
Linear algebraic system.
Definition: elasticity.hh:141
std::shared_ptr< DOFHandlerMultiDim > dh_
Objects for distribution of dofs.
Definition: elasticity.hh:135
std::map< LongIdx, LongIdx > constraint_idx
Definition: elasticity.hh:147
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
Definition: elasticity.hh:81
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
Definition: elasticity.hh:103
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
Definition: elasticity.hh:104
Field< 3, FieldValue< 3 >::Scalar > cross_section_min
Definition: elasticity.hh:80
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
Definition: elasticity.hh:106
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
Definition: elasticity.hh:102
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:68
Field< 3, FieldValue< 3 >::Scalar > lame_lambda
Definition: elasticity.hh:97
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
Definition: elasticity.hh:107
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:69
Field< 3, FieldValue< 3 >::TensorFixed > initial_stress
Definition: elasticity.hh:76
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:70
Field< 3, FieldValue< 3 >::Scalar > ref_potential_load
Potential of reference external load on boundary. TODO: Switch to BCField when possible.
Definition: elasticity.hh:82
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:72
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_mean_stress_ptr
Definition: elasticity.hh:105
arma::mat33 stress_tensor(BulkPoint &p, const arma::mat33 &strain_tensor)
Definition: elasticity.cc:113
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Definition: elasticity.hh:75
Field< 3, FieldValue< 3 >::Scalar > dirichlet_penalty
Definition: elasticity.hh:98
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:79
BCField< 3, FieldValue< 3 >::TensorFixed > bc_stress
Definition: elasticity.hh:71
Field< 3, FieldValue< 3 >::Scalar > lame_mu
Definition: elasticity.hh:96
Data of output parameters.
Definition: elasticity.hh:158
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
Objects for distribution of dofs.
Definition: elasticity.hh:175
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
Definition: elasticity.hh:176
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
Directing class of FieldValueCache.
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
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Generic class of assemblation.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
LocDofVec dof_indices_scalar_
Vector of local DOF indices of scalar fields.
VectorMPI output_vec_
Data vectors of output fields (FieldFE).
unsigned int n_dofs_
Number of dofs.
FieldSet used_fields_
Sub field set contains fields used in calculation.
static constexpr const char * name()
LocDofVec dof_indices_
Vector of local DOF indices of vector fields.
LocDofVec dof_indices_tensor_
Vector of local DOF indices of tensor fields.
EqFields * eq_fields_
Data objects shared with Elasticity.
arma::mat33 normal_stress_
Holds constributions of normal stress.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
double normal_displacement_
Holds constributions of normal displacement.
ElQ< Vector > normal_
Following data members represent Element quantities and FE quantities.
OutpuFieldsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Elasticity::OutputEqData EqData
Range helper class.
vector< LongIdx > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
ElQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Range< JoinShapeAccessor< Vector > > deform_join_
EqFields * eq_fields_
Data objects shared with Elasticity.
Elasticity::EqFields EqFields
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
Elasticity::EqData EqData
RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
vector< PetscScalar > local_rhs_
Auxiliary vector for assemble methods.
unsigned int n_dofs_
Number of dofs.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
static constexpr const char * name()
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
FieldSet used_fields_
Sub field set contains fields used in calculation.
Boundary cond() const
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
FieldSet used_fields_
Sub field set contains fields used in calculation.
vector< LongIdx > side_dof_indices_
vector of DOF indices in neighbour calculation.
static constexpr const char * name()
ElQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
EqFields * eq_fields_
Data objects shared with Elasticity.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
Range< JoinShapeAccessor< Vector > > deform_join_
unsigned int n_dofs_
Number of dofs.
StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Range< JoinShapeAccessor< Tensor > > deform_join_grad_
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
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.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ FEVector
@ coupling
@ boundary
@ bulk
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
unsigned int uint
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Definitions of particular quadrature rules on simplices.