Flow123d  JS_before_hm-1602-g5680f2c
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)
45  : AssemblyBase<dim>(1), eq_fields_(eq_fields), eq_data_(eq_data) {
53  }
54 
55  /// Destructor.
57 
58  /// Initialize auxiliary vectors and other data members
59  void initialize(ElementCacheMap *element_cache_map) {
60  //this->balance_ = eq_data_->balance_;
61  this->element_cache_map_ = element_cache_map;
62 
63  shared_ptr<FE_P<dim>> fe_p = std::make_shared< FE_P<dim> >(1);
64  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
65  fe_ = std::make_shared<FESystem<dim>>(fe_p, FEVector, 3);
66  fe_low_ = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
67  fe_values_.initialize(*this->quad_, *fe_,
71  fe_values_sub_.initialize(*this->quad_low_, *fe_low_,
73 
74  n_dofs_ = fe_->n_dofs();
75  n_dofs_sub_ = fe_low_->n_dofs();
77  dof_indices_.resize(n_dofs_);
78  side_dof_indices_.resize(2);
79  side_dof_indices_[0].resize(n_dofs_sub_); // index 0 = element with lower dimension,
80  side_dof_indices_[1].resize(n_dofs_); // index 1 = side of element with higher dimension
81  local_matrix_.resize(n_dofs_*n_dofs_);
82  local_matrix_ngh_.resize(2);
83  for (uint m=0; m<2; ++m) {
84  local_matrix_ngh_[m].resize(2);
85  for (uint n=0; n<2; ++n)
86  local_matrix_ngh_[m][n].resize(n_dofs_*n_dofs_);
87  }
90  if (dim>1) vec_view_sub_ = &fe_values_sub_.vector_view(0);
91  }
92 
93 
94  /// Assemble integral over element
95  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
96  {
97  if (cell.dim() != dim) return;
98 
99  ElementAccessor<3> elm_acc = cell.elm();
100 
101  fe_values_.reinit(elm_acc);
103 
104  // assemble the local stiffness matrix
105  for (unsigned int i=0; i<n_dofs_; i++)
106  for (unsigned int j=0; j<n_dofs_; j++)
107  local_matrix_[i*n_dofs_+j] = 0;
108 
109  unsigned int k=0;
110  for (auto p : this->bulk_points(element_patch_idx) )
111  {
112  for (unsigned int i=0; i<n_dofs_; i++)
113  {
114  for (unsigned int j=0; j<n_dofs_; j++)
115  local_matrix_[i*n_dofs_+j] += eq_fields_->cross_section(p)*(
116  2*eq_fields_->lame_mu(p)*arma::dot(vec_view_->sym_grad(j,k), vec_view_->sym_grad(i,k))
118  )*fe_values_.JxW(k);
119  }
120  k++;
121  }
123  }
124 
125  /// Assembles boundary integral.
126  inline void boundary_side_integral(DHCellSide cell_side)
127  {
128  ASSERT_EQ_DBG(cell_side.dim(), dim).error("Dimension of element mismatch!");
129  if (!cell_side.cell().is_own()) return;
130 
131  Side side = cell_side.side();
132  const DHCellAccessor &dh_cell = cell_side.cell();
133  dh_cell.get_dof_indices(dof_indices_);
134  fe_values_side_.reinit(side);
135 
136  for (unsigned int i=0; i<n_dofs_; i++)
137  for (unsigned int j=0; j<n_dofs_; j++)
138  local_matrix_[i*n_dofs_+j] = 0;
139 
140  auto p_side = *( this->boundary_points(cell_side).begin() );
141  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
142  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
143  double side_measure = cell_side.measure();
144  if (bc_type == EqFields::bc_type_displacement)
145  {
146  unsigned int k=0;
147  for (auto p : this->boundary_points(cell_side) ) {
148  for (unsigned int i=0; i<n_dofs_; i++)
149  for (unsigned int j=0; j<n_dofs_; j++)
150  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
151  arma::dot(vec_view_side_->value(i,k),vec_view_side_->value(j,k)) * fe_values_side_.JxW(k);
152  k++;
153  }
154  }
155  else if (bc_type == EqFields::bc_type_displacement_normal)
156  {
157  unsigned int k=0;
158  for (auto p : this->boundary_points(cell_side) ) {
159  for (unsigned int i=0; i<n_dofs_; i++)
160  for (unsigned int j=0; j<n_dofs_; j++)
161  local_matrix_[i*n_dofs_+j] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
162  arma::dot(vec_view_side_->value(i,k), fe_values_side_.normal_vector(k)) *
164  k++;
165  }
166  }
167 
169  }
170 
171 
172  /// Assembles between elements of different dimensions.
173  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
174  if (dim == 1) return;
175  ASSERT_EQ_DBG(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
176 
177  cell_lower_dim.get_dof_indices(side_dof_indices_[0]);
178  ElementAccessor<3> cell_sub = cell_lower_dim.elm();
179  fe_values_sub_.reinit(cell_sub);
180 
181  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
182  cell_higher_dim.get_dof_indices(side_dof_indices_[1]);
183  fe_values_side_.reinit(neighb_side.side());
184 
185  // Element id's for testing if they belong to local partition.
186  bool own_element_id[2];
187  own_element_id[0] = cell_lower_dim.is_own();
188  own_element_id[1] = cell_higher_dim.is_own();
189 
190  for (unsigned int n=0; n<2; ++n)
191  for (unsigned int i=0; i<n_dofs_; i++)
192  for (unsigned int m=0; m<2; ++m)
193  for (unsigned int j=0; j<n_dofs_; j++)
194  local_matrix_ngh_[n][m][i*(n_dofs_)+j] = 0;
195 
196  // set transmission conditions
197  unsigned int k=0;
198  for (auto p_high : this->coupling_points(neighb_side) )
199  {
200  auto p_low = p_high.lower_dim(cell_lower_dim);
202 
203  for (int n=0; n<2; n++)
204  {
205  if (!own_element_id[n]) continue;
206 
207  for (unsigned int i=0; i<n_dofs_ngh_[n]; i++)
208  {
209  arma::vec3 vi = (n==0) ? arma::zeros(3) : vec_view_side_->value(i,k);
210  arma::vec3 vf = (n==1) ? arma::zeros(3) : vec_view_sub_->value(i,k);
211  arma::mat33 gvft = (n==0) ? vec_view_sub_->grad(i,k) : arma::zeros(3,3);
212 
213  for (int m=0; m<2; m++)
214  for (unsigned int j=0; j<n_dofs_ngh_[m]; j++) {
215  arma::vec3 ui = (m==0) ? arma::zeros(3) : vec_view_side_->value(j,k);
216  arma::vec3 uf = (m==1) ? arma::zeros(3) : vec_view_sub_->value(j,k);
217  arma::mat33 guit = (m==1) ? mat_t(vec_view_side_->grad(j,k),nv) : arma::zeros(3,3);
218  double divuit = (m==1) ? arma::trace(guit) : 0;
219 
220  local_matrix_ngh_[n][m][i*n_dofs_ngh_[m] + j] +=
221  eq_fields_->fracture_sigma(p_low)*(
222  arma::dot(vf-vi,
223  2/eq_fields_->cross_section(p_low)*(eq_fields_->lame_mu(p_low)*(uf-ui)+(eq_fields_->lame_mu(p_low)+eq_fields_->lame_lambda(p_low))*(arma::dot(uf-ui,nv)*nv))
224  + eq_fields_->lame_mu(p_low)*arma::trans(guit)*nv
225  + eq_fields_->lame_lambda(p_low)*divuit*nv
226  )
227  - arma::dot(gvft, eq_fields_->lame_mu(p_low)*arma::kron(nv,ui.t()) + eq_fields_->lame_lambda(p_low)*arma::dot(ui,nv)*arma::eye(3,3))
228  )*fe_values_sub_.JxW(k);
229  }
230 
231  }
232  }
233  k++;
234  }
235 
236  for (unsigned int n=0; n<2; ++n)
237  for (unsigned int m=0; m<2; ++m)
239  }
240 
241 
242 
243 private:
244  inline arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
245  {
246  arma::mat33 mt = m - m*arma::kron(n,n.t());
247  return mt;
248  }
249 
250 
251 
252  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
253  shared_ptr<FiniteElement<dim-1>> fe_low_; ///< Finite element for the solution of the advection-diffusion equation (dim-1).
254 
255  /// Data objects shared with Elasticity
256  EqFields *eq_fields_;
257  EqData *eq_data_;
258 
259  /// Sub field set contains fields used in calculation.
261 
262  unsigned int n_dofs_; ///< Number of dofs
263  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
264  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
265  FEValues<3> fe_values_; ///< FEValues of cell object (FESystem of P disc finite element type)
266  FEValues<3> fe_values_side_; ///< FEValues of side object
267  FEValues<3> fe_values_sub_; ///< FEValues of lower dimension cell object
268 
269  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
270  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
271  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
272  vector<vector<vector<PetscScalar>>> local_matrix_ngh_; ///< Auxiliary vectors for assemble ngh integral
273  const FEValuesViews::Vector<3> * vec_view_; ///< Vector view in cell integral calculation.
274  const FEValuesViews::Vector<3> * vec_view_side_; ///< Vector view in boundary / neighbour calculation.
275  const FEValuesViews::Vector<3> * vec_view_sub_; ///< Vector view of low dim element in neighbour calculation.
276 
277  template < template<IntDim...> class DimAssembly>
278  friend class GenericAssembly;
279 
280 };
281 
282 
283 template <unsigned int dim>
285 {
286 public:
287  typedef typename Elasticity::EqFields EqFields;
288  typedef typename Elasticity::EqData EqData;
289 
290  static constexpr const char * name() { return "RhsAssemblyElasticity"; }
291 
292  /// Constructor.
293  RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
294  : AssemblyBase<dim>(1), eq_fields_(eq_fields), eq_data_(eq_data) {
297  this->used_fields_ += eq_fields_->load;
301  this->used_fields_ += eq_fields_->bc_type;
304  }
305 
306  /// Destructor.
308 
309  /// Initialize auxiliary vectors and other data members
310  void initialize(ElementCacheMap *element_cache_map) {
311  //this->balance_ = eq_data_->balance_;
312  this->element_cache_map_ = element_cache_map;
313 
314  shared_ptr<FE_P<dim>> fe_p = std::make_shared< FE_P<dim> >(1);
315  shared_ptr<FE_P<dim-1>> fe_p_low = std::make_shared< FE_P<dim-1> >(1);
316  fe_ = std::make_shared<FESystem<dim>>(fe_p, FEVector, 3);
317  fe_low_ = std::make_shared<FESystem<dim-1>>(fe_p_low, FEVector, 3);
318  fe_values_.initialize(*this->quad_, *fe_,
320  fe_values_bdr_side_.initialize(*this->quad_low_, *fe_,
322  fe_values_side_.initialize(*this->quad_low_, *fe_,
324  fe_values_sub_.initialize(*this->quad_low_, *fe_low_,
326  n_dofs_ = fe_->n_dofs();
327  n_dofs_sub_ = fe_low_->n_dofs();
329  dof_indices_.resize(n_dofs_);
330  side_dof_indices_.resize(2);
331  side_dof_indices_[0].resize(n_dofs_sub_); // index 0 = element with lower dimension,
332  side_dof_indices_[1].resize(n_dofs_); // index 1 = side of element with higher dimension
333  local_rhs_.resize(n_dofs_);
334  local_rhs_ngh_.resize(2);
335  for (uint n=0; n<2; ++n) local_rhs_ngh_[n].resize(n_dofs_);
337  vec_view_bdr_ = &fe_values_bdr_side_.vector_view(0);
339  if (dim>1) vec_view_sub_ = &fe_values_sub_.vector_view(0);
340  }
341 
342 
343  /// Assemble integral over element
344  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
345  {
346  if (cell.dim() != dim) return;
347  if (!cell.is_own()) return;
348 
349  ElementAccessor<3> elm_acc = cell.elm();
350 
351  fe_values_.reinit(elm_acc);
353 
354  // assemble the local stiffness matrix
355  fill_n(&(local_rhs_[0]), n_dofs_, 0);
356  //local_source_balance_vector.assign(n_dofs_, 0);
357  //local_source_balance_rhs.assign(n_dofs_, 0);
358 
359  // compute sources
360  unsigned int k=0;
361  for (auto p : this->bulk_points(element_patch_idx) )
362  {
363  for (unsigned int i=0; i<n_dofs_; i++)
364  local_rhs_[i] += (
365  arma::dot(eq_fields_->load(p), vec_view_->value(i,k))
368  ++k;
369  }
370  eq_data_->ls->rhs_set_values(n_dofs_, dof_indices_.data(), &(local_rhs_[0]));
371 
372 // for (unsigned int i=0; i<n_dofs_; i++)
373 // {
374 // for (unsigned int k=0; k<qsize_; k++) // point range
375 // local_source_balance_vector[i] -= 0;//sources_sigma[k]*fe_values_[vec_view_].value(i,k)*fe_values_.JxW(k);
376 //
377 // local_source_balance_rhs[i] += local_rhs_[i];
378 // }
379 // balance_->add_source_matrix_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_vector);
380 // balance_->add_source_vec_values(subst_idx, elm_acc.region().bulk_idx(), dof_indices_, local_source_balance_rhs);
381  }
382 
383  /// Assembles boundary integral.
384  inline void boundary_side_integral(DHCellSide cell_side)
385  {
386  ASSERT_EQ_DBG(cell_side.dim(), dim).error("Dimension of element mismatch!");
387  if (!cell_side.cell().is_own()) return;
388 
389  Side side = cell_side.side();
390  const DHCellAccessor &dh_cell = cell_side.cell();
391  dh_cell.get_dof_indices(dof_indices_);
392  fe_values_bdr_side_.reinit(side);
393 
394  auto p_side = *( this->boundary_points(cell_side).begin() );
395  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
396  unsigned int bc_type = eq_fields_->bc_type(p_bdr);
397 
398  fill_n(&(local_rhs_[0]), n_dofs_, 0);
399  // local_flux_balance_vector.assign(n_dofs_, 0);
400  // local_flux_balance_rhs = 0;
401 
402  unsigned int k=0;
403  if (bc_type == EqFields::bc_type_displacement)
404  {
405  double side_measure = cell_side.measure();
406  for (auto p : this->boundary_points(cell_side) )
407  {
408  for (unsigned int i=0; i<n_dofs_; i++)
409  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
410  arma::dot(eq_fields_->bc_displacement(p), vec_view_bdr_->value(i,k)) *
411  fe_values_bdr_side_.JxW(k);
412  ++k;
413  }
414  }
415  else if (bc_type == EqFields::bc_type_displacement_normal)
416  {
417  double side_measure = cell_side.measure();
418  for (auto p : this->boundary_points(cell_side) )
419  {
420  for (unsigned int i=0; i<n_dofs_; i++)
421  local_rhs_[i] += (eq_fields_->dirichlet_penalty(p) / side_measure) *
422  arma::dot(eq_fields_->bc_displacement(p), fe_values_bdr_side_.normal_vector(k)) *
423  arma::dot(vec_view_bdr_->value(i,k), fe_values_bdr_side_.normal_vector(k)) *
424  fe_values_bdr_side_.JxW(k);
425  ++k;
426  }
427  }
428  else if (bc_type == EqFields::bc_type_traction)
429  {
430  for (auto p : this->boundary_points(cell_side) )
431  {
432  for (unsigned int i=0; i<n_dofs_; i++)
433  local_rhs_[i] += eq_fields_->cross_section(p) *
434  arma::dot(vec_view_bdr_->value(i,k), eq_fields_->bc_traction(p) + eq_fields_->potential_load(p) * fe_values_bdr_side_.normal_vector(k)) *
435  fe_values_bdr_side_.JxW(k);
436  ++k;
437  }
438  }
439  eq_data_->ls->rhs_set_values(n_dofs_, dof_indices_.data(), &(local_rhs_[0]));
440 
441 
442 // balance_->add_flux_matrix_values(subst_idx, loc_b, side_dof_indices, local_flux_balance_vector);
443 // balance_->add_flux_vec_value(subst_idx, loc_b, local_flux_balance_rhs);
444  // ++loc_b;
445  }
446 
447 
448  /// Assembles between elements of different dimensions.
449  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
450  if (dim == 1) return;
451  ASSERT_EQ_DBG(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
452 
453  cell_lower_dim.get_dof_indices(side_dof_indices_[0]);
454  ElementAccessor<3> cell_sub = cell_lower_dim.elm();
455  fe_values_sub_.reinit(cell_sub);
456 
457  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
458  cell_higher_dim.get_dof_indices(side_dof_indices_[1]);
459  fe_values_side_.reinit(neighb_side.side());
460 
461  // Element id's for testing if they belong to local partition.
462  bool own_element_id[2];
463  own_element_id[0] = cell_lower_dim.is_own();
464  own_element_id[1] = cell_higher_dim.is_own();
465 
466  for (unsigned int n=0; n<2; ++n)
467  for (unsigned int i=0; i<n_dofs_; i++)
468  local_rhs_ngh_[n][i] = 0;
469 
470  // set transmission conditions
471  unsigned int k=0;
472  for (auto p_high : this->coupling_points(neighb_side) )
473  {
474  auto p_low = p_high.lower_dim(cell_lower_dim);
476 
477  for (int n=0; n<2; n++)
478  {
479  if (!own_element_id[n]) continue;
480 
481  for (unsigned int i=0; i<n_dofs_ngh_[n]; i++)
482  {
483  arma::vec3 vi = (n==0) ? arma::zeros(3) : vec_view_side_->value(i,k);
484  arma::vec3 vf = (n==1) ? arma::zeros(3) : vec_view_sub_->value(i,k);
485 
486  local_rhs_ngh_[n][i] -= eq_fields_->fracture_sigma(p_low) * eq_fields_->cross_section(p_high) *
487  arma::dot(vf-vi, eq_fields_->potential_load(p_high) * nv) * fe_values_sub_.JxW(k);
488  }
489  }
490  ++k;
491  }
492 
493  for (unsigned int n=0; n<2; ++n)
494  eq_data_->ls->rhs_set_values(n_dofs_ngh_[n], side_dof_indices_[n].data(), &(local_rhs_ngh_[n][0]));
495  }
496 
497 
498 
499 private:
500  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
501  shared_ptr<FiniteElement<dim-1>> fe_low_; ///< Finite element for the solution of the advection-diffusion equation (dim-1).
502 
503  /// Data objects shared with Elasticity
504  EqFields *eq_fields_;
505  EqData *eq_data_;
506 
507  /// Sub field set contains fields used in calculation.
509 
510  unsigned int n_dofs_; ///< Number of dofs
511  unsigned int n_dofs_sub_; ///< Number of dofs (on lower dim element)
512  std::vector<unsigned int> n_dofs_ngh_; ///< Number of dofs on lower and higher dimension element (vector of 2 items)
513  FEValues<3> fe_values_; ///< FEValues of cell object (FESystem of P disc finite element type)
514  FEValues<3> fe_values_bdr_side_; ///< FEValues of side (boundary integral) object
515  FEValues<3> fe_values_side_; ///< FEValues of side (neighbour integral) object
516  FEValues<3> fe_values_sub_; ///< FEValues of lower dimension cell object
517 
518  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
519  vector<vector<LongIdx> > side_dof_indices_; ///< 2 items vector of DOF indices in neighbour calculation.
520  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for assemble methods
521  vector<vector<PetscScalar>> local_rhs_ngh_; ///< Auxiliary vectors for assemble ngh integral
522  const FEValuesViews::Vector<3> * vec_view_; ///< Vector view in cell integral calculation.
523  const FEValuesViews::Vector<3> * vec_view_bdr_; ///< Vector view in boundary calculation.
524  const FEValuesViews::Vector<3> * vec_view_side_; ///< Vector view in neighbour calculation.
525  const FEValuesViews::Vector<3> * vec_view_sub_; ///< Vector view of low dim element in neighbour calculation.
526 
527 
528  template < template<IntDim...> class DimAssembly>
529  friend class GenericAssembly;
530 
531 };
532 
533 template <unsigned int dim>
535 {
536 public:
537  typedef typename Elasticity::EqFields EqFields;
538  typedef typename Elasticity::EqData EqData;
539 
540  static constexpr const char * name() { return "OutpuFieldsAssemblyElasticity"; }
541 
542  /// Constructor.
543  OutpuFieldsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
544  : AssemblyBase<dim>(0), eq_fields_(eq_fields), eq_data_(eq_data) {
547  this->used_fields_ += eq_fields_->lame_mu;
549  }
550 
551  /// Destructor.
553 
554  /// Initialize auxiliary vectors and other data members
555  void initialize(ElementCacheMap *element_cache_map) {
556  //this->balance_ = eq_data_->balance_;
557  this->element_cache_map_ = element_cache_map;
558 
559  shared_ptr<FE_P<dim>> fe_p = std::make_shared< FE_P<dim> >(1);
560  fe_ = std::make_shared<FESystem<dim>>(fe_p, FEVector, 3);
561  fv_.initialize(*this->quad_, *fe_,
563  fsv_.initialize(*this->quad_low_, *fe_,
565  n_dofs_ = fe_->n_dofs();
566  vec_view_ = &fv_.vector_view(0);
567  // if (dim>1) ??
568  vec_view_side_ = &fsv_.vector_view(0);
569 
570  output_vec_ = eq_fields_->output_field_ptr->vec();
571  output_stress_vec_ = eq_fields_->output_stress_ptr->vec();
572  output_von_mises_stress_vec_ = eq_fields_->output_von_mises_stress_ptr->vec();
573  output_cross_sec_vec_ = eq_fields_->output_cross_section_ptr->vec();
574  output_div_vec_ = eq_fields_->output_div_ptr->vec();
575  }
576 
577 
578  /// Assemble integral over element
579  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
580  {
581  if (cell.dim() != dim) return;
582  if (!cell.is_own()) return;
583  DHCellAccessor cell_tensor = cell.cell_with_other_dh(eq_data_->dh_tensor_.get());
584  DHCellAccessor cell_scalar = cell.cell_with_other_dh(eq_data_->dh_scalar_.get());
585 
586  auto elm = cell.elm();
587 
588  fv_.reinit(elm);
590  dof_indices_scalar_ = cell_scalar.get_loc_dof_indices();
591  dof_indices_tensor_ = cell_tensor.get_loc_dof_indices();
592 
593  auto p = *( this->bulk_points(element_patch_idx).begin() );
594 
595  arma::mat33 stress = arma::zeros(3,3);
596  double div = 0;
597  for (unsigned int i=0; i<n_dofs_; i++)
598  {
599  stress += (2*eq_fields_->lame_mu(p)*vec_view_->sym_grad(i,0) + eq_fields_->lame_lambda(p)*vec_view_->divergence(i,0)*arma::eye(3,3))
600  * output_vec_.get(dof_indices_[i]);
601  div += vec_view_->divergence(i,0)*output_vec_.get(dof_indices_[i]);
602  }
603 
604  arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
605  double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
606  output_div_vec_.add(dof_indices_scalar_[0], div);
607 
608  for (unsigned int i=0; i<3; i++)
609  for (unsigned int j=0; j<3; j++)
610  output_stress_vec_.add( dof_indices_tensor_[i*3+j], stress(i,j) );
611  output_von_mises_stress_vec_.set( dof_indices_scalar_[0], von_mises_stress );
612 
613  output_cross_sec_vec_.add( dof_indices_scalar_[0], eq_fields_->cross_section(p) );
614  }
615 
616 
617  /// Assembles between elements of different dimensions.
618  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
619  if (dim == 1) return;
620  ASSERT_EQ_DBG(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
621 
622  normal_displacement_ = 0;
623  normal_stress_.zeros();
624 
625  DHCellAccessor cell_higher_dim = neighb_side.cell();
626  DHCellAccessor cell_tensor = cell_lower_dim.cell_with_other_dh(eq_data_->dh_tensor_.get());
627  DHCellAccessor cell_scalar = cell_lower_dim.cell_with_other_dh(eq_data_->dh_scalar_.get());
628  fsv_.reinit(neighb_side.side());
629 
630  dof_indices_ = cell_higher_dim.get_loc_dof_indices();
631  auto p_high = *( this->coupling_points(neighb_side).begin() );
632  auto p_low = p_high.lower_dim(cell_lower_dim);
633 
634  for (unsigned int i=0; i<n_dofs_; i++)
635  {
636  normal_displacement_ -= arma::dot(vec_view_side_->value(i,0)*output_vec_.get(dof_indices_[i]), fsv_.normal_vector(0));
637  arma::mat33 grad = -arma::kron(vec_view_side_->value(i,0)*output_vec_.get(dof_indices_[i]), fsv_.normal_vector(0).t()) / eq_fields_->cross_section(p_low);
638  normal_stress_ += eq_fields_->lame_mu(p_low)*(grad+grad.t()) + eq_fields_->lame_lambda(p_low)*arma::trace(grad)*arma::eye(3,3);
639  }
640 
641  LocDofVec dof_indices_scalar_ = cell_scalar.get_loc_dof_indices();
642  LocDofVec dof_indices_tensor_ = cell_tensor.get_loc_dof_indices();
643  for (unsigned int i=0; i<3; i++)
644  for (unsigned int j=0; j<3; j++)
645  output_stress_vec_.add( dof_indices_tensor_[i*3+j], normal_stress_(i,j) );
646  output_cross_sec_vec_.add( dof_indices_scalar_[0], normal_displacement_ );
647  output_div_vec_.add( dof_indices_scalar_[0], normal_displacement_ / eq_fields_->cross_section(p_low) );
648  }
649 
650 
651 
652 private:
653  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
654 
655  /// Data objects shared with Elasticity
656  EqFields *eq_fields_;
657  EqData *eq_data_;
658 
659  /// Sub field set contains fields used in calculation.
661 
662  unsigned int n_dofs_; ///< Number of dofs
663  FEValues<3> fv_; ///< FEValues of cell object (FESystem of P disc finite element type)
664  FEValues<3> fsv_; ///< FEValues of side (neighbour integral) object
665 
666  LocDofVec dof_indices_; ///< Vector of local DOF indices of vector fields
667  LocDofVec dof_indices_scalar_; ///< Vector of local DOF indices of scalar fields
668  LocDofVec dof_indices_tensor_; ///< Vector of local DOF indices of tensor fields
669  const FEValuesViews::Vector<3> * vec_view_; ///< Vector view in cell integral calculation.
670  const FEValuesViews::Vector<3> * vec_view_side_; ///< Vector view in neighbour calculation.
671 
672  double normal_displacement_; ///< Holds constributions of normal displacement.
673  arma::mat33 normal_stress_; ///< Holds constributions of normal stress.
674 
675  /// Data vectors of output fields (FieldFE).
681 
682  template < template<IntDim...> class DimAssembly>
683  friend class GenericAssembly;
684 
685 };
686 
687 
688 #endif /* ASSEMBLY_ELASTICITY_HH_ */
689 
Shape function values.
Definition: update_flags.hh:87
FEValues< 3 > fe_values_side_
FEValues of side (neighbour integral) object.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
LocDofVec dof_indices_
Vector of local DOF indices of vector fields.
EqFields * eq_fields_
Data objects shared with Elasticity.
const FEValuesViews::Vector< 3 > * vec_view_sub_
Vector view of low dim element in neighbour calculation.
Transformed quadrature weight for cell sides.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:277
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
Field< 3, FieldValue< 3 >::Scalar > lame_mu
Definition: elasticity.hh:87
vector< vector< LongIdx > > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
unsigned int uint
int active_integrals_
Holds mask of active integrals.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
unsigned int dim() const
Return dimension of element appropriate to the side.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
Definition: elasticity.hh:93
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:73
LinSys * ls
Linear algebra system for the transport equation.
Definition: elasticity.hh:131
unsigned int n_dofs_
Number of dofs.
~RhsAssemblyElasticity()
Destructor.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
unsigned int n_dofs_
Number of dofs.
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
static constexpr const char * name()
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
Field< 3, FieldValue< 3 >::Scalar > dirichlet_penalty
Definition: elasticity.hh:89
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
unsigned int dim() const
Return dimension of element appropriate to cell.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
Directing class of FieldValueCache.
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:66
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in neighbour calculation.
FEValues< 3 > fe_values_sub_
FEValues of lower dimension cell object.
Field< 3, FieldValue< 3 >::Scalar > lame_lambda
Definition: elasticity.hh:88
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
double measure() const
vector< vector< LongIdx > > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:67
bool is_own() const
Return true if accessor represents own element (false for ghost element)
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
Definition: elasticity.hh:124
Side side() const
Return Side of given cell and side_idx.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:541
LocDofVec dof_indices_tensor_
Vector of local DOF indices of tensor fields.
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
EqFields * eq_fields_
Data objects shared with Elasticity.
FEM for linear elasticity.
double normal_displacement_
Holds constributions of normal displacement.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
OutpuFieldsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
Definition: elasticity.hh:96
Transformed quadrature points.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
EqFields * eq_fields_
Data objects shared with Elasticity.
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
Definition: elasticity.hh:125
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
FieldSet used_fields_
Sub field set contains fields used in calculation.
static constexpr const char * name()
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
Definition: elasticity.hh:74
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
ElementAccessor< 3 > element() const
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:137
FEValues< 3 > fe_values_bdr_side_
FEValues of side (boundary integral) object.
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:451
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
Definition: elasticity.hh:94
vector< vector< vector< PetscScalar > > > local_matrix_ngh_
Auxiliary vectors for assemble ngh integral.
Elasticity::EqData EqData
FEValues< 3 > fv_
FEValues of cell object (FESystem of P disc finite element type)
arma::mat33 normal_stress_
Holds constributions of normal stress.
Shape function gradients.
Definition: update_flags.hh:95
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
const FEValuesViews::Vector< 3 > * vec_view_bdr_
Vector view in boundary calculation.
Elasticity::EqFields EqFields
Normal vectors.
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in neighbour calculation.
FEValues< 3 > fe_values_
FEValues of cell object (FESystem of P disc finite element type)
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:64
const FEValuesViews::Vector< 3 > * vec_view_sub_
Vector view of low dim element in neighbour calculation.
static constexpr const char * name()
unsigned int n_dofs_
Number of dofs.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
VectorMPI output_vec_
Data vectors of output fields (FieldFE).
FieldSet used_fields_
Sub field set contains fields used in calculation.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
FieldSet used_fields_
Sub field set contains fields used in calculation.
std::shared_ptr< DOFHandlerMultiDim > dh_
Objects for distribution of dofs.
Definition: elasticity.hh:123
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
vector< vector< PetscScalar > > local_rhs_ngh_
Auxiliary vectors for assemble ngh integral.
FEValues< 3 > fe_values_side_
FEValues of side object.
Conforming Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:86
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Generic class of assemblation.
Definitions of particular quadrature rules on simplices.
LocDofVec dof_indices_scalar_
Vector of local DOF indices of scalar fields.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
Definition: elasticity.hh:97
vector< PetscScalar > local_rhs_
Auxiliary vector for assemble methods.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
FEValues< 3 > fsv_
FEValues of side (neighbour integral) object.
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...
FEValues< 3 > fe_values_sub_
FEValues of lower dimension cell object.
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:223
Boundary cond() const
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:65
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
arma::mat::fixed< spacedim, spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of vector-valued shape function.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
Quadrature * quad_
Quadrature used in assembling methods.
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in boundary / neighbour calculation.
Side accessor allows to iterate over sides of DOF handler cell.
ElementAccessor< 3 > element_accessor()
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Definition: elasticity.hh:70
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
Definition: elasticity.hh:95
FEValues< 3 > fe_values_
FEValues of cell object (FESystem of P disc finite element type)
Transformed quadrature weights.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:257