Flow123d  DF_patch_fe_mechanics-a6ba684
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 
92 
93  // assemble the local stiffness matrix
94  for (unsigned int i=0; i<n_dofs_; i++)
95  for (unsigned int j=0; j<n_dofs_; j++)
96  local_matrix_[i*n_dofs_+j] = 0;
97 
98  for (auto p : this->bulk_points(element_patch_idx) )
99  {
100  for (unsigned int i=0; i<n_dofs_; i++)
101  {
102  for (unsigned int j=0; j<n_dofs_; j++)
104  2*eq_fields_->lame_mu(p)*arma::dot(sym_grad_deform_(j,p), sym_grad_deform_(i,p))
106  )*JxW_(p);
107  }
108  }
110  }
111 
112  /// Assembles boundary integral.
113  inline void boundary_side_integral(DHCellSide cell_side)
114  {
115  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
116  if (!cell_side.cell().is_own()) return;
117 
118  Side side = cell_side.side();
119  const DHCellAccessor &dh_cell = cell_side.cell();
120  dh_cell.get_dof_indices(dof_indices_);
121 
122  for (unsigned int i=0; i<n_dofs_; i++)
123  for (unsigned int j=0; j<n_dofs_; j++)
124  local_matrix_[i*n_dofs_+j] = 0;
125 
126  auto p_side = *( this->boundary_points(cell_side).begin() );
127  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
128  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
129  double side_measure = cell_side.measure();
130  if (bc_type == EqFields::bc_type_displacement)
131  {
132  for (auto p : this->boundary_points(cell_side) ) {
133  for (unsigned int i=0; i<n_dofs_; i++)
134  for (unsigned int j=0; j<n_dofs_; j++)
135  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
136  arma::dot(deform_side_(i,p),deform_side_(j,p)) * JxW_side_(p);
137  }
138  }
139  else if (bc_type == EqFields::bc_type_displacement_normal)
140  {
141  for (auto p : this->boundary_points(cell_side) ) {
142  for (unsigned int i=0; i<n_dofs_; i++)
143  for (unsigned int j=0; j<n_dofs_; j++)
144  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
145  arma::dot(deform_side_(i,p), normal_(p)) *
146  arma::dot(deform_side_(j,p), normal_(p)) * JxW_side_(p);
147  }
148  }
149 
151  }
152 
153 
154  /// Assembles between elements of different dimensions.
155  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
156  if (dim == 1) return;
157  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
158 
159  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
160  for(unsigned int i=0; i<n_indices; ++i) {
162  }
163 
164  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
165  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
166  for(unsigned int i=0; i<n_indices; ++i) {
168  }
169 
170  // Element id's for testing if they belong to local partition.
171  bool own_element_id[2];
172  own_element_id[0] = cell_lower_dim.is_own();
173  own_element_id[1] = cell_higher_dim.is_own();
174 
175  for (unsigned int i=0; i<n_dofs_ngh_[0]+n_dofs_ngh_[1]; i++)
176  for (unsigned int j=0; j<n_dofs_ngh_[0]+n_dofs_ngh_[1]; j++)
177  local_matrix_[i*(n_dofs_ngh_[0]+n_dofs_ngh_[1])+j] = 0;
178 
179  // set transmission conditions
180  for (auto p_high : this->coupling_points(neighb_side) )
181  {
182  auto p_low = p_high.lower_dim(cell_lower_dim);
183  arma::vec3 nv = normal_(p_high);
184 
185  auto deform_shape_i = deform_join_.begin();
186  auto deform_grad_i = deform_join_grad_.begin();
187  for( ; deform_shape_i != deform_join_.end() && deform_grad_i != deform_join_grad_.end(); ++deform_shape_i, ++deform_grad_i) {
188  uint is_high_i = deform_shape_i->is_high_dim();
189  if (!own_element_id[is_high_i]) continue;
190  uint i_mat_idx = deform_shape_i->join_idx();
191  arma::vec3 diff_deform_i = (*deform_shape_i)(p_low) - (*deform_shape_i)(p_high);
192  arma::mat33 grad_deform_i = (*deform_grad_i)(p_low); // low dim element
193 
194  auto deform_shape_j = deform_join_.begin();
195  auto deform_grad_j = deform_join_grad_.begin();
196  for( ; deform_shape_j != deform_join_.end() && deform_grad_j != deform_join_grad_.end(); ++deform_shape_j, ++deform_grad_j) {
197  uint j_mat_idx = deform_shape_j->join_idx();
198  arma::vec3 deform_j_high = (*deform_shape_j)(p_high);
199  arma::vec3 diff_deform_j = (*deform_shape_j)(p_low) - (*deform_shape_j)(p_high);
200  arma::mat33 grad_deform_j = mat_t( arma::trans((*deform_grad_j)(p_low)), nv); // low dim element
201  double div_deform_j = arma::trace(grad_deform_j);
202 
203  local_matrix_[i_mat_idx * (n_dofs_ngh_[0]+n_dofs_ngh_[1]) + j_mat_idx] +=
204  eq_fields_->fracture_sigma(p_low)*(
205  arma::dot(diff_deform_i,
206  2/eq_fields_->cross_section(p_low)*(eq_fields_->lame_mu(p_low)*(diff_deform_j)+(eq_fields_->lame_mu(p_low)+eq_fields_->lame_lambda(p_low))*(arma::dot(diff_deform_j,nv)*nv))
207  + eq_fields_->lame_mu(p_low)*arma::trans(grad_deform_j)*nv
208  + eq_fields_->lame_lambda(p_low)*div_deform_j*nv
209  )
210  - arma::dot(arma::trans(grad_deform_i), eq_fields_->lame_mu(p_low)*arma::kron(nv,deform_j_high.t()) + eq_fields_->lame_lambda(p_low)*arma::dot(deform_j_high,nv)*arma::eye(3,3))
211  )*JxW_side_(p_high);
212  }
213  }
214  }
215 
217  }
218 
219 
220 
221 private:
222  inline arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
223  {
224  arma::mat33 mt = m - m*arma::kron(n,n.t());
225  return mt;
226  }
227 
228 
229  /// Data objects shared with Elasticity
232 
233  /// Sub field set contains fields used in calculation.
235 
236  unsigned int n_dofs_; ///< Number of dofs
237  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
238  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
239 
240  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
241  vector<LongIdx> side_dof_indices_; ///< vector of DOF indices in neighbour calculation.
242  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
243 
244  /// Following data members represent Element quantities and FE quantities
254 
255  template < template<IntDim...> class DimAssembly>
256  friend class GenericAssembly;
257 
258 };
259 
260 
261 template <unsigned int dim>
263 {
264 public:
265  typedef typename Elasticity::EqFields EqFields;
266  typedef typename Elasticity::EqData EqData;
267 
268  static constexpr const char * name() { return "RhsAssemblyElasticity"; }
269 
270  /// Constructor.
271  RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
272  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
273  JxW_( this->bulk_values().JxW() ),
274  JxW_side_( this->side_values().JxW() ),
275  normal_( this->side_values().normal_vector() ),
276  deform_( this->bulk_values().vector_shape() ),
277  deform_side_( this->side_values().vector_shape() ),
278  gras_deform_( this->bulk_values().grad_vector_shape() ),
279  div_deform_( this->bulk_values().vector_divergence() ),
280  deform_join_( Range< JoinShapeAccessor<Vector> >( this->join_values().vector_join_shape() ) ) {
283  this->used_fields_ += eq_fields_->load;
288  this->used_fields_ += eq_fields_->bc_type;
293  }
294 
295  /// Destructor.
297 
298  /// Initialize auxiliary vectors and other data members
299  void initialize(ElementCacheMap *element_cache_map) {
300  //this->balance_ = eq_data_->balance_;
301  this->element_cache_map_ = element_cache_map;
302 
303  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
304  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
305  this->fe_values_->template initialize<dim>(*this->quad_);
306  this->fe_values_->template initialize<dim>(*this->quad_low_);
307 
308  n_dofs_ = this->n_dofs();
309  n_dofs_sub_ = fe_low->n_dofs();
311  dof_indices_.resize(n_dofs_);
313  local_rhs_.resize(2*n_dofs_);
314  }
315 
316 
317  /// Assemble integral over element
318  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
319  {
320  if (cell.dim() != dim) return;
321  if (!cell.is_own()) return;
322 
324 
325  // assemble the local stiffness matrix
326  fill_n(&(local_rhs_[0]), n_dofs_, 0);
327  //local_source_balance_vector.assign(n_dofs_, 0);
328  //local_source_balance_rhs.assign(n_dofs_, 0);
329 
330  // compute sources
331  for (auto p : this->bulk_points(element_patch_idx) )
332  {
333  for (unsigned int i=0; i<n_dofs_; i++)
334  local_rhs_[i] += (
335  arma::dot(eq_fields_->load(p), deform_(i,p))
337  -arma::dot(eq_fields_->initial_stress(p), gras_deform_(i,p))
338  )*eq_fields_->cross_section(p)*JxW_(p);
339  }
341 
342 // for (unsigned int i=0; i<n_dofs_; i++)
343 // {
344 // for (unsigned int k=0; k<qsize_; k++) // point range
345 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_vals_[vec_view_].value(i,k)*JxW_(k);
346 //
347 // local_source_balance_rhs[i] += local_rhs_[i];
348 // }
349 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_vector);
350 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_rhs);
351  }
352 
353  /// Assembles boundary integral.
354  inline void boundary_side_integral(DHCellSide cell_side)
355  {
356  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
357  if (!cell_side.cell().is_own()) return;
358 
359  const DHCellAccessor &dh_cell = cell_side.cell();
360  dh_cell.get_dof_indices(dof_indices_);
361 
362  auto p_side = *( this->boundary_points(cell_side).begin() );
363  auto p_bdr = p_side.point_bdr( cell_side.cond().element_accessor() );
364  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
365 
366  fill_n(&(local_rhs_[0]), n_dofs_, 0);
367  // local_flux_balance_vector.assign(n_dofs_, 0);
368  // local_flux_balance_rhs = 0;
369 
370  // addtion from initial stress
371  for (auto p : this->boundary_points(cell_side) )
372  {
373  for (unsigned int i=0; i<n_dofs_; i++)
375  arma::dot(( eq_fields_->initial_stress(p) * normal_(p)),
376  deform_side_(i,p)) *
377  JxW_side_(p);
378  }
379 
380  if (bc_type == EqFields::bc_type_displacement)
381  {
382  double side_measure = cell_side.measure();
383  for (auto p : this->boundary_points(cell_side) )
384  {
385  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
386  for (unsigned int i=0; i<n_dofs_; i++)
387  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
388  arma::dot(eq_fields_->bc_displacement(p_bdr), deform_side_(i,p)) *
389  JxW_side_(p);
390  }
391  }
392  else if (bc_type == EqFields::bc_type_displacement_normal)
393  {
394  double side_measure = cell_side.measure();
395  for (auto p : this->boundary_points(cell_side) )
396  {
397  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
398  for (unsigned int i=0; i<n_dofs_; i++)
399  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
400  arma::dot(eq_fields_->bc_displacement(p_bdr), normal_(p)) *
401  arma::dot(deform_side_(i,p), normal_(p)) *
402  JxW_side_(p);
403  }
404  }
405  else if (bc_type == EqFields::bc_type_traction)
406  {
407  for (auto p : this->boundary_points(cell_side) )
408  {
409  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
410  for (unsigned int i=0; i<n_dofs_; i++)
412  arma::dot(deform_side_(i,p), eq_fields_->bc_traction(p_bdr) + eq_fields_->ref_potential_load(p) * normal_(p)) *
413  JxW_side_(p);
414  }
415  }
416  else if (bc_type == EqFields::bc_type_stress)
417  {
418  for (auto p : this->boundary_points(cell_side) )
419  {
420  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
421  for (unsigned int i=0; i<n_dofs_; i++)
422  // stress is multiplied by inward normal to obtain traction
424  arma::dot(deform_side_(i,p), -eq_fields_->bc_stress(p_bdr)*normal_(p)
426  * JxW_side_(p);
427  }
428  }
430 
431 
432 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
433 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
434  // ++loc_b;
435  }
436 
437 
438  /// Assembles between elements of different dimensions.
439  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
440  if (dim == 1) return;
441  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
442 
443  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
444  for(unsigned int i=0; i<n_indices; ++i) {
446  }
447 
448  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
449  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
450  for(unsigned int i=0; i<n_indices; ++i) {
452  }
453 
454  // Element id's for testing if they belong to local partition.
455  bool own_element_id[2];
456  own_element_id[0] = cell_lower_dim.is_own();
457  own_element_id[1] = cell_higher_dim.is_own();
458 
459  for (unsigned int i=0; i<2*n_dofs_; i++)
460  local_rhs_[i] = 0;
461 
462  // set transmission conditions
463  for (auto p_high : this->coupling_points(neighb_side) )
464  {
465  auto p_low = p_high.lower_dim(cell_lower_dim);
466  arma::vec3 nv = normal_(p_high);
467 
468  for( auto join_shape_i : deform_join_) {
469  uint is_high_i = join_shape_i.is_high_dim();
470  if (!own_element_id[is_high_i]) continue;
471 
472  arma::vec3 vi = join_shape_i(p_high);
473  arma::vec3 vf = join_shape_i(p_low);
474 
475  local_rhs_[join_shape_i.join_idx()] -= eq_fields_->fracture_sigma(p_low) * eq_fields_->cross_section(p_high) *
476  arma::dot(vf-vi, eq_fields_->potential_load(p_high) * nv) * JxW_side_(p_high);
477  }
478  }
479 
481  }
482 
483 
484 
485 private:
486  /// Data objects shared with Elasticity
489 
490  /// Sub field set contains fields used in calculation.
492 
493  unsigned int n_dofs_; ///< Number of dofs
494  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
495  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
496 
497  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
498  vector<LongIdx> side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
499  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for assemble methods
500 
501  /// Following data members represent Element quantities and FE quantities
510 
511 
512  template < template<IntDim...> class DimAssembly>
513  friend class GenericAssembly;
514 
515 };
516 
517 template <unsigned int dim>
519 {
520 public:
521  typedef typename Elasticity::EqFields EqFields;
523 
524  static constexpr const char * name() { return "OutpuFieldsAssemblyElasticity"; }
525 
526  /// Constructor.
528  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
529  normal_( this->side_values().normal_vector() ),
530  deform_side_( this->side_values().vector_shape() ),
531  gras_deform_( this->bulk_values().grad_vector_shape() ),
532  sym_grad_deform_( this->bulk_values().vector_sym_grad() ),
533  div_deform_( this->bulk_values().vector_divergence() ) {
536  this->used_fields_ += eq_fields_->lame_mu;
539  }
540 
541  /// Destructor.
543 
544  /// Initialize auxiliary vectors and other data members
545  void initialize(ElementCacheMap *element_cache_map) {
546  //this->balance_ = eq_data_->balance_;
547  this->element_cache_map_ = element_cache_map;
548 
549  this->fe_values_->template initialize<dim>(*this->quad_);
550  this->fe_values_->template initialize<dim>(*this->quad_low_);
551 
552  n_dofs_ = this->n_dofs();
553 
560  }
561 
562 
563  /// Assemble integral over element
564  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
565  {
566  if (cell.dim() != dim) return;
567  if (!cell.is_own()) return;
568  DHCellAccessor cell_tensor = cell.cell_with_other_dh(eq_data_->dh_tensor_.get());
569  DHCellAccessor cell_scalar = cell.cell_with_other_dh(eq_data_->dh_scalar_.get());
570 
574 
575  auto p = *( this->bulk_points(element_patch_idx).begin() );
576 
578  double div = 0;
579  for (unsigned int i=0; i<n_dofs_; i++)
580  {
581  stress += (2*eq_fields_->lame_mu(p)*sym_grad_deform_(i,p) + eq_fields_->lame_lambda(p)*div_deform_(i,p)*arma::eye(3,3))
583  div += div_deform_(i,p)*output_vec_.get(dof_indices_[i]);
584  }
585 
586  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
587  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
588  double mean_stress = arma::trace(stress) / 3;
590 
591  for (unsigned int i=0; i<3; i++)
592  for (unsigned int j=0; j<3; j++)
593  output_stress_vec_.add( dof_indices_tensor_[i*3+j], stress(i,j) );
594  output_von_mises_stress_vec_.set( dof_indices_scalar_[0], von_mises_stress );
596 
598  }
599 
600 
601  /// Assembles between elements of different dimensions.
602  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
603  if (dim == 1) return;
604  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
605 
607  normal_stress_.zeros();
608 
609  DHCellAccessor cell_higher_dim = neighb_side.cell();
610  DHCellAccessor cell_tensor = cell_lower_dim.cell_with_other_dh(eq_data_->dh_tensor_.get());
611  DHCellAccessor cell_scalar = cell_lower_dim.cell_with_other_dh(eq_data_->dh_scalar_.get());
612 
613  dof_indices_ = cell_higher_dim.get_loc_dof_indices();
614  auto p_high = *( this->coupling_points(neighb_side).begin() );
615  auto p_low = p_high.lower_dim(cell_lower_dim);
616 
617  for (unsigned int i=0; i<n_dofs_; i++)
618  {
619  normal_displacement_ -= arma::dot(deform_side_(i,p_high)*output_vec_.get(dof_indices_[i]), normal_(p_high));
620  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);
621  normal_stress_ += eq_fields_->lame_mu(p_low)*(grad+grad.t()) + eq_fields_->lame_lambda(p_low)*arma::trace(grad)*arma::eye(3,3);
622  }
623 
626  for (unsigned int i=0; i<3; i++)
627  for (unsigned int j=0; j<3; j++)
631  }
632 
633 
634 
635 private:
636  /// Data objects shared with Elasticity
639 
640  /// Sub field set contains fields used in calculation.
642 
643  unsigned int n_dofs_; ///< Number of dofs
644  LocDofVec dof_indices_; ///< Vector of local DOF indices of vector fields
645  LocDofVec dof_indices_scalar_; ///< Vector of local DOF indices of scalar fields
646  LocDofVec dof_indices_tensor_; ///< Vector of local DOF indices of tensor fields
647 
648  double normal_displacement_; ///< Holds constributions of normal displacement.
649  arma::mat33 normal_stress_; ///< Holds constributions of normal stress.
650 
651  /// Following data members represent Element quantities and FE quantities
657 
658  /// Data vectors of output fields (FieldFE).
665 
666  template < template<IntDim...> class DimAssembly>
667  friend class GenericAssembly;
668 
669 };
670 
671 
672 
673 /**
674  * Container class for assembly of constraint matrix for contact condition.
675  */
676 template <unsigned int dim>
678 {
679 public:
680  typedef typename Elasticity::EqFields EqFields;
681  typedef typename Elasticity::EqData EqData;
682 
683  static constexpr const char * name() { return "ConstraintAssemblyElasticity"; }
684 
685  /// Constructor.
687  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
688  JxW_side_( this->side_values().JxW() ),
689  normal_( this->side_values().normal_vector() ),
690  deform_side_( this->side_values().vector_shape() ) {
694  }
695 
696  /// Destructor.
698 
699  /// Initialize auxiliary vectors and other data members
700  void initialize(ElementCacheMap *element_cache_map) {
701  this->element_cache_map_ = element_cache_map;
702 
703  this->fe_values_->template initialize<dim>(*this->quad_);
704  this->fe_values_->template initialize<dim>(*this->quad_low_);
705 
706  n_dofs_ = this->n_dofs();
707  dof_indices_.resize(n_dofs_);
708  local_matrix_.resize(n_dofs_*n_dofs_);
709  }
710 
711 
712  /// Assembles between elements of different dimensions.
713  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
714  if (dim == 1) return;
715  if (!cell_lower_dim.is_own()) return;
716 
717  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
718 
719  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
720  cell_higher_dim.get_dof_indices(dof_indices_);
721 
722  for (unsigned int i=0; i<n_dofs_; i++)
723  local_matrix_[i] = 0;
724 
725  // Assemble matrix and vector for contact conditions in the form B*x <= c,
726  // where B*x is the average jump of normal displacements and c is the average cross-section on element.
727  // Positive value means that the fracture closes.
728  double local_vector = 0;
729  for (auto p_high : this->coupling_points(neighb_side) )
730  {
731  auto p_low = p_high.lower_dim(cell_lower_dim);
732  arma::vec3 nv = normal_(p_high);
733 
734  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();
735 
736  for (unsigned int i=0; i<n_dofs_; i++)
737  {
738  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();
739  }
740  }
741 
742  int arow[1] = { eq_data_->constraint_idx[cell_lower_dim.elm_idx()] };
743  MatSetValues(eq_data_->constraint_matrix, 1, arow, n_dofs_, dof_indices_.data(), &(local_matrix_[0]), ADD_VALUES);
744  VecSetValue(eq_data_->constraint_vec, arow[0], local_vector, ADD_VALUES);
745  }
746 
747 
748 
749 private:
750 
751 
752  /// Data objects shared with Elasticity
755 
756  /// Sub field set contains fields used in calculation.
758 
759  unsigned int n_dofs_; ///< Number of dofs
760  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
761  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
762  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
763 
764  /// Following data members represent Element quantities and FE quantities
768 
769 
770  template < template<IntDim...> class DimAssembly>
771  friend class GenericAssembly;
772 
773 };
774 
775 
776 #endif /* ASSEMBLY_ELASTICITY_HH_ */
777 
#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:112
LinSys * ls
Linear algebraic system.
Definition: elasticity.hh:139
std::shared_ptr< DOFHandlerMultiDim > dh_
Objects for distribution of dofs.
Definition: elasticity.hh:133
std::map< LongIdx, LongIdx > constraint_idx
Definition: elasticity.hh:145
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
Definition: elasticity.hh:79
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
Definition: elasticity.hh:101
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
Definition: elasticity.hh:102
Field< 3, FieldValue< 3 >::Scalar > cross_section_min
Definition: elasticity.hh:78
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
Definition: elasticity.hh:104
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
Definition: elasticity.hh:100
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:66
Field< 3, FieldValue< 3 >::Scalar > lame_lambda
Definition: elasticity.hh:95
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
Definition: elasticity.hh:105
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:67
Field< 3, FieldValue< 3 >::TensorFixed > initial_stress
Definition: elasticity.hh:74
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:68
Field< 3, FieldValue< 3 >::Scalar > ref_potential_load
Potential of reference external load on boundary. TODO: Switch to BCField when possible.
Definition: elasticity.hh:80
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:70
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_mean_stress_ptr
Definition: elasticity.hh:103
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Definition: elasticity.hh:73
Field< 3, FieldValue< 3 >::Scalar > dirichlet_penalty
Definition: elasticity.hh:96
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:77
BCField< 3, FieldValue< 3 >::TensorFixed > bc_stress
Definition: elasticity.hh:69
Field< 3, FieldValue< 3 >::Scalar > lame_mu
Definition: elasticity.hh:94
Data of output parameters.
Definition: elasticity.hh:156
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
Objects for distribution of dofs.
Definition: elasticity.hh:173
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
Definition: elasticity.hh:174
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.