18 #ifndef ASSEMBLY_ELASTICITY_HH_
19 #define ASSEMBLY_ELASTICITY_HH_
34 template <
unsigned int dim>
41 static constexpr
const char *
name() {
return "StiffnessAssemblyElasticity"; }
72 shared_ptr<
FE_P<dim-1>> fe_p_low = std::make_shared<
FE_P<dim-1> >(1);
89 if (cell.
dim() != dim)
return;
94 for (
unsigned int i=0; i<
n_dofs_; i++)
95 for (
unsigned int j=0; j<
n_dofs_; j++)
98 for (
auto p : this->
bulk_points(element_patch_idx) )
100 for (
unsigned int i=0; i<
n_dofs_; i++)
102 for (
unsigned int j=0; j<
n_dofs_; j++)
115 ASSERT_EQ(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
122 for (
unsigned int i=0; i<
n_dofs_; i++)
123 for (
unsigned int j=0; j<
n_dofs_; j++)
129 double side_measure = cell_side.
measure();
133 for (
unsigned int i=0; i<
n_dofs_; i++)
134 for (
unsigned int j=0; j<
n_dofs_; j++)
142 for (
unsigned int i=0; i<
n_dofs_; i++)
143 for (
unsigned int j=0; j<
n_dofs_; j++)
156 if (dim == 1)
return;
157 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
160 for(
unsigned int i=0; i<n_indices; ++i) {
166 for(
unsigned int i=0; i<n_indices; ++i) {
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();
182 auto p_low = p_high.lower_dim(cell_lower_dim);
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);
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);
201 double div_deform_j = arma::trace(grad_deform_j);
205 arma::dot(diff_deform_i,
255 template <
template<
IntDim...>
class DimAssembly>
261 template <
unsigned int dim>
268 static constexpr
const char *
name() {
return "RhsAssemblyElasticity"; }
303 shared_ptr<
FE_P<dim-1>> fe_p_low = std::make_shared<
FE_P<dim-1> >(1);
320 if (cell.
dim() != dim)
return;
321 if (!cell.
is_own())
return;
331 for (
auto p : this->
bulk_points(element_patch_idx) )
333 for (
unsigned int i=0; i<
n_dofs_; i++)
356 ASSERT_EQ(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
373 for (
unsigned int i=0; i<
n_dofs_; i++)
382 double side_measure = cell_side.
measure();
386 for (
unsigned int i=0; i<
n_dofs_; i++)
394 double side_measure = cell_side.
measure();
398 for (
unsigned int i=0; i<
n_dofs_; i++)
410 for (
unsigned int i=0; i<
n_dofs_; i++)
421 for (
unsigned int i=0; i<
n_dofs_; i++)
440 if (dim == 1)
return;
441 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
444 for(
unsigned int i=0; i<n_indices; ++i) {
450 for(
unsigned int i=0; i<n_indices; ++i) {
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();
459 for (
unsigned int i=0; i<2*
n_dofs_; i++)
465 auto p_low = p_high.lower_dim(cell_lower_dim);
469 uint is_high_i = join_shape_i.is_high_dim();
470 if (!own_element_id[is_high_i])
continue;
512 template <
template<
IntDim...>
class DimAssembly>
517 template <
unsigned int dim>
524 static constexpr
const char *
name() {
return "OutpuFieldsAssemblyElasticity"; }
566 if (cell.
dim() != dim)
return;
567 if (!cell.
is_own())
return;
575 auto p = *( this->
bulk_points(element_patch_idx).begin() );
579 for (
unsigned int i=0; i<
n_dofs_; i++)
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;
591 for (
unsigned int i=0; i<3; i++)
592 for (
unsigned int j=0; j<3; j++)
603 if (dim == 1)
return;
604 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
615 auto p_low = p_high.lower_dim(cell_lower_dim);
617 for (
unsigned int i=0; i<
n_dofs_; i++)
626 for (
unsigned int i=0; i<3; i++)
627 for (
unsigned int j=0; j<3; j++)
666 template <
template<
IntDim...>
class DimAssembly>
676 template <
unsigned int dim>
683 static constexpr
const char *
name() {
return "ConstraintAssemblyElasticity"; }
714 if (dim == 1)
return;
715 if (!cell_lower_dim.
is_own())
return;
717 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
722 for (
unsigned int i=0; i<
n_dofs_; i++)
728 double local_vector = 0;
731 auto p_low = p_high.lower_dim(cell_lower_dim);
736 for (
unsigned int i=0; i<
n_dofs_; i++)
770 template <
template<
IntDim...>
class DimAssembly>
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
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()
~ConstraintAssemblyElasticity()
Destructor.
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.
FeQ< Vector > deform_side_
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.
Elasticity::EqFields EqFields
Elasticity::EqData EqData
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.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
ElementAccessor< 3 > element() const
Data of equation parameters.
LinSys * ls
Linear algebraic system.
std::shared_ptr< DOFHandlerMultiDim > dh_
Objects for distribution of dofs.
std::map< LongIdx, LongIdx > constraint_idx
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
@ bc_type_displacement_normal
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
Field< 3, FieldValue< 3 >::Scalar > cross_section_min
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
BCField< 3, FieldValue< 3 >::Enum > bc_type
Field< 3, FieldValue< 3 >::Scalar > lame_lambda
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Field< 3, FieldValue< 3 >::TensorFixed > initial_stress
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Field< 3, FieldValue< 3 >::Scalar > ref_potential_load
Potential of reference external load on boundary. TODO: Switch to BCField when possible.
Field< 3, FieldValue< 3 >::VectorFixed > load
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_mean_stress_ptr
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Field< 3, FieldValue< 3 >::Scalar > dirichlet_penalty
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
BCField< 3, FieldValue< 3 >::TensorFixed > bc_stress
Field< 3, FieldValue< 3 >::Scalar > lame_mu
Data of output parameters.
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
Objects for distribution of dofs.
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
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...
Directing class of FieldValueCache.
unsigned int n_neighs_vb() const
Return number of neighbours.
Compound finite element on dim dimensional simplex.
Conforming Lagrangean finite element on dim dimensional simplex.
Container for various descendants of FieldCommonBase.
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
FeQ< Vector > deform_side_
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.
VectorMPI output_div_vec_
VectorMPI output_stress_vec_
FeQ< Tensor > gras_deform_
VectorMPI output_von_mises_stress_vec_
FeQ< Tensor > sym_grad_deform_
LocDofVec dof_indices_tensor_
Vector of local DOF indices of tensor fields.
EqFields * eq_fields_
Data objects shared with Elasticity.
FeQ< Scalar > div_deform_
arma::mat33 normal_stress_
Holds constributions of normal stress.
Elasticity::EqFields EqFields
VectorMPI output_cross_sec_vec_
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.
VectorMPI output_mean_stress_vec_
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.
~OutpuFieldsAssemblyElasticity()
Destructor.
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.
ElQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
~RhsAssemblyElasticity()
Destructor.
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
FeQ< Scalar > div_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()
FeQ< Vector > deform_side_
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
FieldSet used_fields_
Sub field set contains fields used in calculation.
FeQ< Tensor > gras_deform_
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.
Elasticity::EqFields EqFields
Elasticity::EqData EqData
vector< LongIdx > side_dof_indices_
vector of DOF indices in neighbour calculation.
static constexpr const char * name()
FeQ< Vector > deform_side_
ElQ< Scalar > JxW_
Following data members represent Element quantities and FE quantities.
FeQ< Tensor > gras_deform_
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
~StiffnessAssemblyElasticity()
Destructor.
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.
FeQ< Tensor > sym_grad_deform_
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)
FeQ< Scalar > div_deform_
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.
void set(unsigned int pos, double val)
Set value on given position.
void add(unsigned int pos, double val)
Add value to item on given position.
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,...
arma::Col< IntIdx > LocDofVec
unsigned int IntDim
A dimension index type.
Definitions of particular quadrature rules on simplices.