Flow123d
JS_before_hm-1711-g43b724b37
|
Go to the documentation of this file.
18 #ifndef ASSEMBLY_ELASTICITY_HH_
19 #define ASSEMBLY_ELASTICITY_HH_
34 template <
unsigned int dim>
41 static constexpr
const char *
name() {
return "StiffnessAssemblyElasticity"; }
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);
83 for (
uint m=0; m<2; ++m) {
85 for (
uint n=0; n<2; ++n)
97 if (cell.
dim() != dim)
return;
105 for (
unsigned int i=0; i<
n_dofs_; i++)
106 for (
unsigned int j=0; j<
n_dofs_; j++)
110 for (
auto p : this->
bulk_points(element_patch_idx) )
112 for (
unsigned int i=0; i<
n_dofs_; i++)
114 for (
unsigned int j=0; j<
n_dofs_; j++)
128 ASSERT_EQ_DBG(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
136 for (
unsigned int i=0; i<
n_dofs_; i++)
137 for (
unsigned int j=0; j<
n_dofs_; j++)
143 double side_measure = cell_side.
measure();
148 for (
unsigned int i=0; i<
n_dofs_; i++)
149 for (
unsigned int j=0; j<
n_dofs_; j++)
159 for (
unsigned int i=0; i<
n_dofs_; i++)
160 for (
unsigned int j=0; j<
n_dofs_; j++)
174 if (dim == 1)
return;
175 ASSERT_EQ_DBG(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
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();
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++)
200 auto p_low = p_high.lower_dim(cell_lower_dim);
203 for (
int n=0; n<2; n++)
205 if (!own_element_id[n])
continue;
213 for (
int m=0; m<2; m++)
218 double divuit = (m==1) ? arma::trace(guit) : 0;
236 for (
unsigned int n=0; n<2; ++n)
237 for (
unsigned int m=0; m<2; ++m)
252 shared_ptr<FiniteElement<dim>>
fe_;
277 template <
template<
IntDim...>
class DimAssembly>
283 template <
unsigned int dim>
290 static constexpr
const char *
name() {
return "RhsAssemblyElasticity"; }
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);
346 if (cell.
dim() != dim)
return;
347 if (!cell.
is_own())
return;
361 for (
auto p : this->
bulk_points(element_patch_idx) )
363 for (
unsigned int i=0; i<
n_dofs_; i++)
386 ASSERT_EQ_DBG(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
405 double side_measure = cell_side.
measure();
409 for (
unsigned int i=0; i<
n_dofs_; i++)
418 double side_measure = cell_side.
measure();
422 for (
unsigned int i=0; i<
n_dofs_; i++)
435 for (
unsigned int i=0; i<
n_dofs_; i++)
453 if (dim == 1)
return;
454 ASSERT_EQ_DBG(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
465 bool own_element_id[2];
466 own_element_id[0] = cell_lower_dim.
is_own();
467 own_element_id[1] = cell_higher_dim.
is_own();
469 for (
unsigned int n=0; n<2; ++n)
470 for (
unsigned int i=0; i<
n_dofs_; i++)
477 auto p_low = p_high.lower_dim(cell_lower_dim);
480 for (
int n=0; n<2; n++)
482 if (!own_element_id[n])
continue;
496 for (
unsigned int n=0; n<2; ++n)
503 shared_ptr<FiniteElement<dim>>
fe_;
531 template <
template<
IntDim...>
class DimAssembly>
536 template <
unsigned int dim>
543 static constexpr
const char *
name() {
return "OutpuFieldsAssemblyElasticity"; }
562 shared_ptr<FE_P<dim>> fe_p = std::make_shared< FE_P<dim> >(1);
563 fe_ = std::make_shared<FESystem<dim>>(fe_p,
FEVector, 3);
584 if (cell.
dim() != dim)
return;
585 if (!cell.
is_own())
return;
589 auto elm = cell.
elm();
596 auto p = *( this->
bulk_points(element_patch_idx).begin() );
600 for (
unsigned int i=0; i<
n_dofs_; i++)
607 arma::mat33 stress_dev = stress - arma::trace(stress)/3*arma::eye(3,3);
608 double von_mises_stress = sqrt(1.5*arma::dot(stress_dev, stress_dev));
611 for (
unsigned int i=0; i<3; i++)
612 for (
unsigned int j=0; j<3; j++)
622 if (dim == 1)
return;
623 ASSERT_EQ_DBG(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
635 auto p_low = p_high.lower_dim(cell_lower_dim);
637 for (
unsigned int i=0; i<
n_dofs_; i++)
646 for (
unsigned int i=0; i<3; i++)
647 for (
unsigned int j=0; j<3; j++)
656 shared_ptr<FiniteElement<dim>>
fe_;
685 template <
template<
IntDim...>
class DimAssembly>
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Elasticity::EqData EqData
arma::Col< IntIdx > LocDofVec
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::TensorFixed > > output_stress_ptr
const FEValuesViews::Vector< 3 > * vec_view_sub_
Vector view of low dim element in neighbour calculation.
@ update_gradients
Shape function gradients.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Side side() const
Return Side of given cell and side_idx.
Elasticity::EqFields EqFields
double divergence(unsigned int function_no, unsigned int point_no) const
Return divergence of vector-valued shape function.
static constexpr const char * name()
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
FEValues< 3 > fe_values_bdr_side_
FEValues of side (boundary integral) object.
FieldSet used_fields_
Sub field set contains fields used in calculation.
FEValues< 3 > fe_values_sub_
FEValues of lower dimension cell object.
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in boundary / neighbour calculation.
vector< vector< PetscScalar > > local_rhs_ngh_
Auxiliary vectors for assemble ngh integral.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
LinSys * ls
Linear algebraic system.
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Directing class of FieldValueCache.
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
FEValues< 3 > fe_values_sub_
FEValues of lower dimension cell object.
BCField< 3, FieldValue< 3 >::Enum > bc_type
EqFields * eq_fields_
Data objects shared with Elasticity.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
bool is_own() const
Return true if accessor represents own element (false for ghost element)
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
unsigned int n_dofs_
Number of dofs.
~StiffnessAssemblyElasticity()
Destructor.
@ update_values
Shape function values.
~RhsAssemblyElasticity()
Destructor.
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
double normal_displacement_
Holds constributions of normal displacement.
vector< vector< LongIdx > > side_dof_indices_
2 items vector of DOF indices in neighbour calculation.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
@ update_quadrature_points
Transformed quadrature points.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
unsigned int dim() const
Return dimension of element appropriate to cell.
@ update_normal_vectors
Normal vectors.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_div_ptr
unsigned int dim() const
Return dimension of element appropriate to the side.
FEValues< 3 > fe_values_
FEValues of cell object (FESystem of P disc finite element type)
LocDofVec dof_indices_
Vector of local DOF indices of vector fields.
std::shared_ptr< DOFHandlerMultiDim > dh_
Objects for distribution of dofs.
Elasticity::EqFields EqFields
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
void set(unsigned int pos, double val)
Set value on given position.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Side accessor allows to iterate over sides of DOF handler cell.
arma::mat::fixed< spacedim, spacedim > sym_grad(unsigned int function_no, unsigned int point_no) const
Return symmetric gradient of vector-valued shape function.
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Elasticity::EqData EqData
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_von_mises_stress_ptr
FEValues< 3 > fe_values_
FEValues of cell object (FESystem of P disc finite element type)
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
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...
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in neighbour calculation.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
VectorMPI output_vec_
Data vectors of output fields (FieldFE).
Conforming Lagrangean finite element on dim dimensional simplex.
arma::mat33 normal_stress_
Holds constributions of normal stress.
FEM for linear elasticity.
Definitions of particular quadrature rules on simplices.
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
EqFields * eq_fields_
Data objects shared with Elasticity.
const FEValuesViews::Vector< 3 > * vec_view_side_
Vector view in neighbour calculation.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
StiffnessAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
Field< 3, FieldValue< 3 >::Scalar > lame_mu
VectorMPI output_von_mises_stress_vec_
arma::mat::fixed< spacedim, spacedim > grad(unsigned int function_no, unsigned int point_no) const
Return gradient of vector-valued shape function.
arma::mat33 mat_t(const arma::mat33 &m, const arma::vec3 &n)
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > output_field_ptr
RhsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
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.
Container for various descendants of FieldCommonBase.
unsigned int n_dofs_sub_
Number of dofs (on lower dim element)
const FEValuesViews::Vector< 3 > * vec_view_sub_
Vector view of low dim element in neighbour calculation.
ElementAccessor< 3 > element_accessor()
LocDofVec dof_indices_tensor_
Vector of local DOF indices of tensor fields.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
VectorMPI output_div_vec_
void add(unsigned int pos, double val)
Add value to item on given position.
ElementAccessor< 3 > element() const
VectorMPI output_cross_sec_vec_
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
OutpuFieldsAssemblyElasticity(EqFields *eq_fields, EqData *eq_data)
Constructor.
FEValues< 3 > fsv_
FEValues of side (neighbour integral) object.
std::shared_ptr< DOFHandlerMultiDim > dh_tensor_
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Elasticity::EqFields EqFields
Cell accessor allow iterate over DOF handler cells.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
FEValues< 3 > fe_values_side_
FEValues of side object.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Field< 3, FieldValue< 3 >::Scalar > dirichlet_penalty
vector< vector< vector< PetscScalar > > > local_matrix_ngh_
Auxiliary vectors for assemble ngh integral.
Field< 3, FieldValue< 3 >::Scalar > lame_lambda
Field< 3, FieldValue< 3 >::VectorFixed > load
static constexpr const char * name()
int active_integrals_
Holds mask of active integrals.
@ update_JxW_values
Transformed quadrature weights.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
arma::vec::fixed< spacedim > value(unsigned int function_no, unsigned int point_no) const
Return value of vector-valued shape function.
EqFields * eq_fields_
Data objects shared with Elasticity.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
double get(unsigned int pos) const
Return value on given position.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
VectorMPI output_stress_vec_
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
FEValues< 3 > fv_
FEValues of cell object (FESystem of P disc finite element type)
const FEValuesViews::Vector< 3 > * vec_view_
Vector view in cell integral calculation.
FieldSet used_fields_
Sub field set contains fields used in calculation.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles between elements of different dimensions.
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > output_cross_section_ptr
unsigned int n_dofs_
Number of dofs.
FieldSet used_fields_
Sub field set contains fields used in calculation.
std::shared_ptr< DOFHandlerMultiDim > dh_scalar_
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
void boundary_side_integral(DHCellSide cell_side)
Assembles boundary integral.
Generic class of assemblation.
Compound finite element on dim dimensional simplex.
static constexpr const char * name()
FEValues< 3 > fe_values_side_
FEValues of side (neighbour integral) object.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
Elasticity::EqData EqData
unsigned int n_dofs_
Number of dofs.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
std::vector< unsigned int > n_dofs_ngh_
Number of dofs on lower and higher dimension element (vector of 2 items)
LocDofVec dof_indices_scalar_
Vector of local DOF indices of scalar fields.
Field< 3, FieldValue< 3 >::Scalar > potential_load
Potential of an additional (external) load.
const FEValuesViews::Vector< 3 > * vec_view_bdr_
Vector view in boundary calculation.
Quadrature * quad_
Quadrature used in assembling methods.
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
vector< PetscScalar > local_rhs_
Auxiliary vector for assemble methods.
@ bc_type_displacement_normal
~OutpuFieldsAssemblyElasticity()
Destructor.
unsigned int IntDim
A dimension index type.