1 #ifndef LINEARELASTICITY_H 2 #define LINEARELASTICITY_H 1
49 #include <lifev/core/LifeV.hpp> 51 #include <lifev/core/array/VectorEpetra.hpp> 53 #include <lifev/core/array/MatrixEpetra.hpp> 55 #include <lifev/core/fem/BCManage.hpp> 57 #include <lifev/eta/fem/ETFESpace.hpp> 59 #include <lifev/core/fem/FESpace.hpp> 61 #include <lifev/eta/expression/Integrate.hpp> 66 class LinearElasticity
69 typedef Epetra_Comm comm_Type;
71 typedef std::shared_ptr< comm_Type > commPtr_Type;
73 typedef VectorEpetra vector_Type;
75 typedef std::shared_ptr<vector_Type> vectorPtr_Type;
79 typedef std::shared_ptr<matrix_Type> matrixPtr_Type;
81 typedef std::shared_ptr<BCHandler> bcPtr_Type;
83 typedef RegionMesh<LinearTetra> mesh_Type;
85 typedef std::shared_ptr<mesh_Type> meshPtr_Type;
89 typedef ETFESpace< mesh_Type, map_Type, 3, 3 > ETFESpace_displacement;
97 LinearElasticity(
const commPtr_Type& communicator );
108 void setCoefficients (
const Real density,
const Real young,
const Real poisson);
118 void setCoefficientsThinLayer (
const Real density,
const Real young,
const Real poisson,
const Real thickness,
const UInt interface );
125 void setup(
const meshPtr_Type& mesh,
const std::string dOrder );
133 void assemble_matrices (
const Real timestep,
const Real beta, bcPtr_Type & bc,
bool useBDF =
false );
139 matrixPtr_Type
const& mass_matrix_no_bc ( )
const {
return M_mass_no_bc; };
145 matrixPtr_Type
const& stiffness_matrix_no_bc ( )
const {
return M_stiffness_no_bc; };
151 matrixPtr_Type
const& jacobian ( )
const {
return M_jacobian; };
158 const std::shared_ptr<FESpace<mesh_Type, map_Type> >& fespace ( )
const {
return M_displacementFESpace; };
165 const std::shared_ptr<ETFESpace_displacement >& et_fespace ( )
const {
return M_displacementFESpace_ETA; };
182 std::string M_dOrder;
185 std::shared_ptr<FESpace<mesh_Type, map_Type> > M_displacementFESpace;
188 std::shared_ptr<ETFESpace_displacement > M_displacementFESpace_ETA;
191 matrixPtr_Type M_mass_no_bc;
192 matrixPtr_Type M_stiffness_no_bc;
195 matrixPtr_Type M_mass_no_bc_thin;
196 matrixPtr_Type M_stiffness_no_bc_thin;
199 matrixPtr_Type M_jacobian;
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
Epetra_Import const & importer()
Getter for the Epetra_Import.
double Real
Generic real data.
class ETFESpace A light, templated version of the FESpace
uint32_type UInt
generic unsigned integer (used mainly for addressing)