Flow123d  DF_patch_fe_mechanics-5faa023
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"
26 #include "fem/op_factory.hh"
27 #include "fem/patch_op_impl.hh"
29 #include "coupling/balance.hh"
30 #include "fem/element_cache_map.hh"
31 
32 
33 /**
34  * Auxiliary container class for Finite element and related objects of given dimension.
35  */
36 template <unsigned int dim>
38 {
39 public:
40  typedef typename Elasticity::EqFields EqFields;
41  typedef typename Elasticity::EqData EqData;
42 
43  static constexpr const char * name() { return "StiffnessAssemblyElasticity"; }
44 
45  /// Constructor.
46  StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
47  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data), // quad_order = 1
48  JxW_( this->bulk_values().JxW() ),
49  JxW_side_( this->side_values().JxW() ),
50  normal_( this->side_values().normal_vector() ),
51  deform_side_( this->side_values().vector_shape() ),
52  grad_deform_( this->bulk_values().grad_vector_shape() ),
53  sym_grad_deform_( this->bulk_values().vector_sym_grad() ),
54  div_deform_( this->bulk_values().vector_divergence() ),
55  deform_join_( this->join_values().vector_join_shape() ),
56  deform_join_grad_( this->join_values().gradient_vector_join_shape() ) {
64  }
65 
66  /// Destructor.
68 
69  /// Initialize auxiliary vectors and other data members
70  void initialize(ElementCacheMap *element_cache_map) {
71  //this->balance_ = eq_data_->balance_;
72  this->element_cache_map_ = element_cache_map;
73 
74  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
75  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
76  this->fe_values_->template initialize<dim>(*this->quad_);
77  this->fe_values_->template initialize<dim>(*this->quad_low_);
78 
79  n_dofs_ = this->n_dofs();
80  n_dofs_sub_ = fe_low->n_dofs();
82  dof_indices_.resize(n_dofs_);
83  side_dof_indices_.resize(2*n_dofs_);
84  local_matrix_.resize(4*n_dofs_*n_dofs_);
85  }
86 
87 
88  /// Assemble integral over element
89  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
90  {
91  if (cell.dim() != dim) return;
92 
93  // Fracture stiffness is assembled in dimjoin_integral.
94  if (cell.elm()->n_neighs_vb() > 0) return;
95 
97 
98  // assemble the local stiffness matrix
99  for (unsigned int i=0; i<n_dofs_; i++)
100  for (unsigned int j=0; j<n_dofs_; j++)
101  local_matrix_[i*n_dofs_+j] = 0;
102 
103  for (auto p : this->bulk_points(element_patch_idx) )
104  {
105  for (unsigned int i=0; i<n_dofs_; i++)
106  {
107  for (unsigned int j=0; j<n_dofs_; j++)
109  arma::dot(eq_fields_->stress_tensor(p,sym_grad_deform_.shape(j)(p)), sym_grad_deform_.shape(i)(p))
110  )*JxW_(p);
111  }
112  }
114  }
115 
116  /// Assembles boundary integral.
117  inline void boundary_side_integral(DHCellSide cell_side)
118  {
119  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
120  if (!cell_side.cell().is_own()) return;
121 
122  Side side = cell_side.side();
123  const DHCellAccessor &dh_cell = cell_side.cell();
124  dh_cell.get_dof_indices(dof_indices_);
125 
126  for (unsigned int i=0; i<n_dofs_; i++)
127  for (unsigned int j=0; j<n_dofs_; j++)
128  local_matrix_[i*n_dofs_+j] = 0;
129 
130  auto p_side = *( this->boundary_points(cell_side).begin() );
131  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
132  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
133  double side_measure = cell_side.measure();
134  if (bc_type == EqFields::bc_type_displacement)
135  {
136  for (auto p : this->boundary_points(cell_side) ) {
137  for (unsigned int i=0; i<n_dofs_; i++)
138  for (unsigned int j=0; j<n_dofs_; j++)
139  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
140  arma::dot(deform_side_.shape(i)(p),deform_side_.shape(j)(p)) * JxW_side_(p);
141  }
142  }
143  else if (bc_type == EqFields::bc_type_displacement_normal)
144  {
145  for (auto p : this->boundary_points(cell_side) ) {
146  for (unsigned int i=0; i<n_dofs_; i++)
147  for (unsigned int j=0; j<n_dofs_; j++)
148  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
149  arma::dot(deform_side_.shape(i)(p), normal_(p)) *
150  arma::dot(deform_side_.shape(j)(p), normal_(p)) * JxW_side_(p);
151  }
152  }
153 
155  }
156 
157 
158  /// Assembles between elements of different dimensions.
159  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
160  if (dim == 1) return;
161  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
162 
163  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
164  for(unsigned int i=0; i<n_indices; ++i) {
166  }
167 
168  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
169  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
170  for(unsigned int i=0; i<n_indices; ++i) {
172  }
173 
174  // Element id's for testing if they belong to local partition.
175  bool own_element_id[2];
176  own_element_id[0] = cell_lower_dim.is_own();
177  own_element_id[1] = cell_higher_dim.is_own();
178 
179  unsigned int n_neighs = cell_lower_dim.elm()->n_neighs_vb();
180 
181  for (unsigned int i=0; i<n_dofs_ngh_[0]+n_dofs_ngh_[1]; i++)
182  for (unsigned int j=0; j<n_dofs_ngh_[0]+n_dofs_ngh_[1]; j++)
183  local_matrix_[i*(n_dofs_ngh_[0]+n_dofs_ngh_[1])+j] = 0;
184 
185  // set transmission conditions
186  for (auto p_high : this->coupling_points(neighb_side) )
187  {
188  auto p_low = p_high.lower_dim(cell_lower_dim);
189  arma::vec3 nv = normal_(p_high);
190 
191  for (uint i=0; i<deform_join_.n_dofs_both(); ++i) {
192  uint is_high_i = deform_join_.is_high_dim(i);
193  if (!own_element_id[is_high_i]) continue;
194  arma::vec3 diff_deform_i = deform_join_.shape(i)(p_low) - deform_join_.shape(i)(p_high);
195  arma::mat33 grad_deform_i = deform_join_grad_.shape(i)(p_low); // low dim element
196  arma::mat33 semi_grad_i = grad_deform_i + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_i,nv.t());
197  arma::mat33 semi_sym_grad_i = 0.5*(semi_grad_i + semi_grad_i.t());
198 
199  for (uint j=0; j<deform_join_.n_dofs_both(); ++j) {
200  arma::vec3 deform_j_high = deform_join_.shape(j)(p_high);
201  arma::vec3 diff_deform_j = deform_join_.shape(j)(p_low) - deform_j_high;
202  arma::mat33 grad_deform_j = deform_join_grad_.shape(j)(p_low); // low dim element
203  arma::mat33 semi_grad_j = grad_deform_j + n_neighs/eq_fields_->cross_section(p_low)*arma::kron(diff_deform_j,nv.t());
204  arma::mat33 semi_sym_grad_j = 0.5*(semi_grad_j + semi_grad_j.t());
205 
206  local_matrix_[i * (n_dofs_ngh_[0]+n_dofs_ngh_[1]) + j] +=
207  (
208  eq_fields_->fracture_sigma(p_low)*eq_fields_->cross_section(p_low) / n_neighs
209  * arma::dot(semi_sym_grad_i, eq_fields_->stress_tensor(p_low,semi_sym_grad_j))
210 
211  // This term corrects the tangential part of the above product so that there is no
212  // dependence on fracture_sigma.
213  // TODO: Fracture_sigma should be possibly removed and replaced by anisotropic elasticity.
214  + (1-eq_fields_->fracture_sigma(p_low))*eq_fields_->cross_section(p_low) / n_neighs
215  * 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())))
216  )*JxW_side_(p_high);
217  }
218  }
219 
220  }
221 
223  }
224 
225 
226 
227 private:
228  inline arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
229  {
230  arma::mat33 mt = m - m*arma::kron(n,n.t());
231  return mt;
232  }
233 
234 
235  /// Data objects shared with Elasticity
238 
239  /// Sub field set contains fields used in calculation.
241 
242  unsigned int n_dofs_; ///< Number of dofs
243  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
244  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
245 
246  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
247  vector<LongIdx> side_dof_indices_; ///< vector of DOF indices in neighbour calculation.
248  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
249 
250  /// Following data members represent Element quantities and FE quantities
260 
261  template < template<IntDim...> class DimAssembly>
262  friend class GenericAssembly;
263 
264 };
265 
266 
267 template <unsigned int dim>
269 {
270 public:
271  typedef typename Elasticity::EqFields EqFields;
272  typedef typename Elasticity::EqData EqData;
273 
274  static constexpr const char * name() { return "RhsAssemblyElasticity"; }
275 
276  /// Constructor.
277  RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
278  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
279  JxW_( this->bulk_values().JxW() ),
280  JxW_side_( this->side_values().JxW() ),
281  normal_( this->side_values().normal_vector() ),
282  deform_( this->bulk_values().vector_shape() ),
283  deform_side_( this->side_values().vector_shape() ),
284  grad_deform_( this->bulk_values().grad_vector_shape() ),
285  div_deform_( this->bulk_values().vector_divergence() ),
286  deform_join_( this->join_values().vector_join_shape() ) {
289  this->used_fields_ += eq_fields_->load;
294  this->used_fields_ += eq_fields_->bc_type;
299  }
300 
301  /// Destructor.
303 
304  /// Initialize auxiliary vectors and other data members
305  void initialize(ElementCacheMap *element_cache_map) {
306  //this->balance_ = eq_data_->balance_;
307  this->element_cache_map_ = element_cache_map;
308 
309  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
310  shared_ptr<FiniteElement<dim-1>> fe_low = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
311  this->fe_values_->template initialize<dim>(*this->quad_);
312  this->fe_values_->template initialize<dim>(*this->quad_low_);
313 
314  n_dofs_ = this->n_dofs();
315  n_dofs_sub_ = fe_low->n_dofs();
317  dof_indices_.resize(n_dofs_);
319  local_rhs_.resize(2*n_dofs_);
320  }
321 
322 
323  /// Assemble integral over element
324  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
325  {
326  if (cell.dim() != dim) return;
327  if (!cell.is_own()) return;
328 
330 
331  // assemble the local stiffness matrix
332  fill_n(&(local_rhs_[0]), n_dofs_, 0);
333  //local_source_balance_vector.assign(n_dofs_, 0);
334  //local_source_balance_rhs.assign(n_dofs_, 0);
335 
336  // compute sources
337  for (auto p : this->bulk_points(element_patch_idx) )
338  {
339  for (unsigned int i=0; i<n_dofs_; i++)
340  local_rhs_[i] += (
341  arma::dot(eq_fields_->load(p), deform_.shape(i)(p))
343  -arma::dot(eq_fields_->initial_stress(p), grad_deform_.shape(i)(p))
344  )*eq_fields_->cross_section(p)*JxW_(p);
345  }
347 
348 // for (unsigned int i=0; i<n_dofs_; i++)
349 // {
350 // for (unsigned int k=0; k<qsize_; k++) // point range
351 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_vals_[vec_view_].value(i,k)*JxW_(k);
352 //
353 // local_source_balance_rhs[i] += local_rhs_[i];
354 // }
355 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_vector);
356 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_rhs);
357  }
358 
359  /// Assembles boundary integral.
360  inline void boundary_side_integral(DHCellSide cell_side)
361  {
362  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
363  if (!cell_side.cell().is_own()) return;
364 
365  const DHCellAccessor &dh_cell = cell_side.cell();
366  dh_cell.get_dof_indices(dof_indices_);
367 
368  auto p_side = *( this->boundary_points(cell_side).begin() );
369  auto p_bdr = p_side.point_bdr( cell_side.cond().element_accessor() );
370  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
371 
372  fill_n(&(local_rhs_[0]), n_dofs_, 0);
373  // local_flux_balance_vector.assign(n_dofs_, 0);
374  // local_flux_balance_rhs = 0;
375 
376  // addtion from initial stress
377  for (auto p : this->boundary_points(cell_side) )
378  {
379  for (unsigned int i=0; i<n_dofs_; i++)
381  arma::dot(( eq_fields_->initial_stress(p) * normal_(p)),
382  deform_side_.shape(i)(p)) *
383  JxW_side_(p);
384  }
385 
386  if (bc_type == EqFields::bc_type_displacement)
387  {
388  double side_measure = cell_side.measure();
389  for (auto p : this->boundary_points(cell_side) )
390  {
391  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
392  for (unsigned int i=0; i<n_dofs_; i++)
393  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
394  arma::dot(eq_fields_->bc_displacement(p_bdr), deform_side_.shape(i)(p)) *
395  JxW_side_(p);
396  }
397  }
398  else if (bc_type == EqFields::bc_type_displacement_normal)
399  {
400  double side_measure = cell_side.measure();
401  for (auto p : this->boundary_points(cell_side) )
402  {
403  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
404  for (unsigned int i=0; i<n_dofs_; i++)
405  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
406  arma::dot(eq_fields_->bc_displacement(p_bdr), normal_(p)) *
407  arma::dot(deform_side_.shape(i)(p), normal_(p)) *
408  JxW_side_(p);
409  }
410  }
411  else if (bc_type == EqFields::bc_type_traction)
412  {
413  for (auto p : this->boundary_points(cell_side) )
414  {
415  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
416  for (unsigned int i=0; i<n_dofs_; i++)
418  arma::dot(deform_side_.shape(i)(p), eq_fields_->bc_traction(p_bdr) + eq_fields_->ref_potential_load(p) * normal_(p)) *
419  JxW_side_(p);
420  }
421  }
422  else if (bc_type == EqFields::bc_type_stress)
423  {
424  for (auto p : this->boundary_points(cell_side) )
425  {
426  auto p_bdr = p.point_bdr( cell_side.cond().element_accessor() );
427  for (unsigned int i=0; i<n_dofs_; i++)
428  // stress is multiplied by inward normal to obtain traction
430  arma::dot(deform_side_.shape(i)(p), -eq_fields_->bc_stress(p_bdr)*normal_(p)
432  * JxW_side_(p);
433  }
434  }
436 
437 
438 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
439 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
440  // ++loc_b;
441  }
442 
443 
444  /// Assembles between elements of different dimensions.
445  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
446  if (dim == 1) return;
447  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
448 
449  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
450  for(unsigned int i=0; i<n_indices; ++i) {
452  }
453 
454  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
455  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
456  for(unsigned int i=0; i<n_indices; ++i) {
458  }
459 
460  // Element id's for testing if they belong to local partition.
461  bool own_element_id[2];
462  own_element_id[0] = cell_lower_dim.is_own();
463  own_element_id[1] = cell_higher_dim.is_own();
464 
465  for (unsigned int i=0; i<2*n_dofs_; i++)
466  local_rhs_[i] = 0;
467 
468  // set transmission conditions
469  for (auto p_high : this->coupling_points(neighb_side) )
470  {
471  auto p_low = p_high.lower_dim(cell_lower_dim);
472  arma::vec3 nv = normal_(p_high);
473 
474  for (uint i=0; i<deform_join_.n_dofs_both(); ++i) {
475  uint is_high_i = deform_join_.is_high_dim(i);
476  if (!own_element_id[is_high_i]) continue;
477 
478  arma::vec3 vi = deform_join_.shape(i)(p_high);
479  arma::vec3 vf = deform_join_.shape(i)(p_low);
480 
481  local_rhs_[i] -= eq_fields_->fracture_sigma(p_low) * eq_fields_->cross_section(p_high) *
482  arma::dot(vf-vi, eq_fields_->potential_load(p_high) * nv) * JxW_side_(p_high);
483  }
484  }
485 
487  }
488 
489 
490 
491 private:
492  /// Data objects shared with Elasticity
495 
496  /// Sub field set contains fields used in calculation.
498 
499  unsigned int n_dofs_; ///< Number of dofs
500  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
501  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
502 
503  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
504  vector<LongIdx> side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
505  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for assemble methods
506 
507  /// Following data members represent Element quantities and FE quantities
516 
517 
518  template < template<IntDim...> class DimAssembly>
519  friend class GenericAssembly;
520 
521 };
522 
523 template <unsigned int dim>
525 {
526 public:
527  typedef typename Elasticity::EqFields EqFields;
529 
530  static constexpr const char * name() { return "OutpuFieldsAssemblyElasticity"; }
531 
532  /// Constructor.
534  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
535  normal_( this->side_values().normal_vector() ),
536  deform_side_( this->side_values().vector_shape() ),
537  grad_deform_( this->bulk_values().grad_vector_shape() ),
538  sym_grad_deform_( this->bulk_values().vector_sym_grad() ),
539  div_deform_( this->bulk_values().vector_divergence() ) {
542  this->used_fields_ += eq_fields_->lame_mu;
545  }
546 
547  /// Destructor.
549 
550  /// Initialize auxiliary vectors and other data members
551  void initialize(ElementCacheMap *element_cache_map) {
552  //this->balance_ = eq_data_->balance_;
553  this->element_cache_map_ = element_cache_map;
554 
555  this->fe_values_->template initialize<dim>(*this->quad_);
556  this->fe_values_->template initialize<dim>(*this->quad_low_);
557 
558  n_dofs_ = this->n_dofs();
559 
566  }
567 
568 
569  /// Assemble integral over element
570  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
571  {
572  if (cell.dim() != dim) return;
573  if (!cell.is_own()) return;
574  DHCellAccessor cell_tensor = cell.cell_with_other_dh(eq_data_->dh_tensor_.get());
575  DHCellAccessor cell_scalar = cell.cell_with_other_dh(eq_data_->dh_scalar_.get());
576 
580 
581  auto p = *( this->bulk_points(element_patch_idx).begin() );
582 
584  double div = 0;
585  for (unsigned int i=0; i<n_dofs_; i++)
586  {
587  stress += eq_fields_->stress_tensor(p,sym_grad_deform_.shape(i)(p))
589  div += div_deform_.shape(i)(p)*output_vec_.get(dof_indices_[i]);
590  }
591 
592  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
593  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
594  double mean_stress = arma::trace(stress) / 3;
596 
597  for (unsigned int i=0; i<3; i++)
598  for (unsigned int j=0; j<3; j++)
599  output_stress_vec_.add( dof_indices_tensor_[i*3+j], stress(i,j) );
600  output_von_mises_stress_vec_.set( dof_indices_scalar_[0], von_mises_stress );
602 
604  }
605 
606 
607  /// Assembles between elements of different dimensions.
608  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
609  if (dim == 1) return;
610  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
611 
613  normal_stress_.zeros();
614 
615  DHCellAccessor cell_higher_dim = neighb_side.cell();
616  DHCellAccessor cell_tensor = cell_lower_dim.cell_with_other_dh(eq_data_->dh_tensor_.get());
617  DHCellAccessor cell_scalar = cell_lower_dim.cell_with_other_dh(eq_data_->dh_scalar_.get());
618 
619  dof_indices_ = cell_higher_dim.get_loc_dof_indices();
620  auto p_high = *( this->coupling_points(neighb_side).begin() );
621  auto p_low = p_high.lower_dim(cell_lower_dim);
622 
623  for (unsigned int i=0; i<n_dofs_; i++)
624  {
625  normal_displacement_ -= arma::dot(deform_side_.shape(i)(p_high)*output_vec_.get(dof_indices_[i]), normal_(p_high));
626  arma::mat33 grad = -arma::kron(deform_side_.shape(i)(p_high)*output_vec_.get(dof_indices_[i]), normal_(p_high).t()) / eq_fields_->cross_section(p_low);
627  normal_stress_ += eq_fields_->stress_tensor(p_low, 0.5*(grad+grad.t()));
628  }
629 
632  for (unsigned int i=0; i<3; i++)
633  for (unsigned int j=0; j<3; j++)
637  }
638 
639 
640 
641 private:
642  /// Data objects shared with Elasticity
645 
646  /// Sub field set contains fields used in calculation.
648 
649  unsigned int n_dofs_; ///< Number of dofs
650  LocDofVec dof_indices_; ///< Vector of local DOF indices of vector fields
651  LocDofVec dof_indices_scalar_; ///< Vector of local DOF indices of scalar fields
652  LocDofVec dof_indices_tensor_; ///< Vector of local DOF indices of tensor fields
653 
654  double normal_displacement_; ///< Holds constributions of normal displacement.
655  arma::mat33 normal_stress_; ///< Holds constributions of normal stress.
656 
657  /// Following data members represent Element quantities and FE quantities
663 
664  /// Data vectors of output fields (FieldFE).
671 
672  template < template<IntDim...> class DimAssembly>
673  friend class GenericAssembly;
674 
675 };
676 
677 
678 
679 /**
680  * Container class for assembly of constraint matrix for contact condition.
681  */
682 template <unsigned int dim>
684 {
685 public:
686  typedef typename Elasticity::EqFields EqFields;
687  typedef typename Elasticity::EqData EqData;
688 
689  static constexpr const char * name() { return "ConstraintAssemblyElasticity"; }
690 
691  /// Constructor.
693  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
694  JxW_side_( this->side_values().JxW() ),
695  normal_( this->side_values().normal_vector() ),
696  deform_side_( this->side_values().vector_shape() ) {
700  }
701 
702  /// Destructor.
704 
705  /// Initialize auxiliary vectors and other data members
706  void initialize(ElementCacheMap *element_cache_map) {
707  this->element_cache_map_ = element_cache_map;
708 
709  this->fe_values_->template initialize<dim>(*this->quad_);
710  this->fe_values_->template initialize<dim>(*this->quad_low_);
711 
712  n_dofs_ = this->n_dofs();
713  dof_indices_.resize(n_dofs_);
714  local_matrix_.resize(n_dofs_*n_dofs_);
715  }
716 
717 
718  /// Assembles between elements of different dimensions.
719  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
720  if (dim == 1) return;
721  if (!cell_lower_dim.is_own()) return;
722 
723  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
724 
725  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
726  cell_higher_dim.get_dof_indices(dof_indices_);
727 
728  for (unsigned int i=0; i<n_dofs_; i++)
729  local_matrix_[i] = 0;
730 
731  // Assemble matrix and vector for contact conditions in the form B*x <= c,
732  // where B*x is the average jump of normal displacements and c is the average cross-section on element.
733  // Positive value means that the fracture closes.
734  double local_vector = 0;
735  for (auto p_high : this->coupling_points(neighb_side) )
736  {
737  auto p_low = p_high.lower_dim(cell_lower_dim);
738  arma::vec3 nv = normal_(p_high);
739 
740  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();
741 
742  for (unsigned int i=0; i<n_dofs_; i++)
743  {
744  local_matrix_[i] += eq_fields_->cross_section(p_high)*arma::dot(deform_side_.shape(i)(p_high), nv)*JxW_side_(p_high) / cell_lower_dim.elm().measure();
745  }
746  }
747 
748  int arow[1] = { eq_data_->constraint_idx[cell_lower_dim.elm_idx()] };
749  MatSetValues(eq_data_->constraint_matrix, 1, arow, n_dofs_, dof_indices_.data(), &(local_matrix_[0]), ADD_VALUES);
750  VecSetValue(eq_data_->constraint_vec, arow[0], local_vector, ADD_VALUES);
751  }
752 
753 
754 
755 private:
756 
757 
758  /// Data objects shared with Elasticity
761 
762  /// Sub field set contains fields used in calculation.
764 
765  unsigned int n_dofs_; ///< Number of dofs
766  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
767  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
768  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
769 
770  /// Following data members represent Element quantities and FE quantities
774 
775 
776  template < template<IntDim...> class DimAssembly>
777  friend class GenericAssembly;
778 
779 };
780 
781 
782 #endif /* ASSEMBLY_ELASTICITY_HH_ */
783 
#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.
FeQ< Scalar > JxW_side_
Following data members represent Element quantities and FE quantities.
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.
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
FeQ< ValueType > shape(unsigned int i_shape_fn_idx) const
Definition: op_accessors.hh:95
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
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.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
FeQArray< Vector > deform_
EqFields * eq_fields_
Data objects shared with Elasticity.
Elasticity::EqFields EqFields
FeQArray< Scalar > div_deform_
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
FeQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
Elasticity::EqData EqData
FeQArray< Tensor > grad_deform_
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.
FeQArray< Vector > deform_side_
FeQJoin< Vector > deform_join_
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()
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)
FeQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
unsigned int n_dofs_
Number of dofs.
StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
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
Declares top level factory classes of FE operations.
Store finite element reinit functions.
Definitions of particular quadrature rules on simplices.