LifeV
LinearElasticity.hpp
Go to the documentation of this file.
1 #ifndef LINEARELASTICITY_H
2 #define LINEARELASTICITY_H 1
3 
4 //@HEADER
5 /*
6  *******************************************************************************
7 
8  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
9  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
10 
11  This file is part of LifeV.
12 
13  LifeV is free software; you can redistribute it and/or modify
14  it under the terms of the GNU Lesser General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  LifeV is distributed in the hope that it will be useful,
19  but WITHOUT ANY WARRANTY; without even the implied warranty of
20  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
21  Lesser General Public License for more details.
22 
23  You should have received a copy of the GNU Lesser General Public License
24  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
25 
26  *******************************************************************************
27  */
28 //@HEADER
29 /*!
30  @file
31  @brief Implementation of Linear elastic structural model
32  @author Davide Forti <davide.forti@epfl.ch>
33  @maintainer Antonello Gerbi <antonello.gerbi@epfl.ch>
34  @date 25-05-2016
35 
36  This file contains an ETA implementation of a linear elastic structural model. This class
37  temporary as it will disappear when the new structure model will be released.
38 
39  This class contains also an implementation of a thin layer structural model which can be
40  used together with the 3D structure to form a multilayered structure. As a reference for applications
41  of such multilayered computational model in the context of fluid-structure interaction, please refer to:
42 
43  * D. Forti, M. Bukac, A. Quaini, S. Čanić, and S. Deparis, "A Monolithic Approach to Fluid-Composite Structure
44  Interaction" available as NA & SC Preprint Series No. 47 - February, 2016. Accepted in JSC.
45  (http://www.uh.edu/nsm/_docs/math/NASC-preprint-series/2015_2016/Preprint_No16-47.pdf)
46 
47 */
48 
49 #include <lifev/core/LifeV.hpp>
50 
51 #include <lifev/core/array/VectorEpetra.hpp>
52 
53 #include <lifev/core/array/MatrixEpetra.hpp>
54 
55 #include <lifev/core/fem/BCManage.hpp>
56 
57 #include <lifev/eta/fem/ETFESpace.hpp>
58 
59 #include <lifev/core/fem/FESpace.hpp>
60 
61 #include <lifev/eta/expression/Integrate.hpp>
62 
63 namespace LifeV
64 {
65 
66 class LinearElasticity
67 {
68 
69  typedef Epetra_Comm comm_Type;
70 
71  typedef std::shared_ptr< comm_Type > commPtr_Type;
72 
73  typedef VectorEpetra vector_Type;
74 
75  typedef std::shared_ptr<vector_Type> vectorPtr_Type;
76 
77  typedef MatrixEpetra<Real> matrix_Type;
78 
79  typedef std::shared_ptr<matrix_Type> matrixPtr_Type;
80 
81  typedef std::shared_ptr<BCHandler> bcPtr_Type;
82 
83  typedef RegionMesh<LinearTetra> mesh_Type;
84 
85  typedef std::shared_ptr<mesh_Type> meshPtr_Type;
86 
87  typedef MapEpetra map_Type;
88 
89  typedef ETFESpace< mesh_Type, map_Type, 3, 3 > ETFESpace_displacement;
90 
91 public:
92 
93  // Constructor
94  /*!
95  * @param communicator Epetra communicator
96  */
97  LinearElasticity( const commPtr_Type& communicator );
98 
99  // empty destructor
100  ~LinearElasticity();
101 
102  //! Set the physical parameters of the structure
103  /*!
104  * @param density Density of the structure
105  * @param young Young modulus of the structure
106  * @param poisson Poisson ratio of the structure
107  */
108  void setCoefficients ( const Real density, const Real young, const Real poisson);
109 
110  //! Set the physical parameters of the thin layered structure
111  /*!
112  * @param density Density of the thin layered structure
113  * @param young Young modulus of the thin layered structure
114  * @param poisson Poisson ratio of the thin layered structure
115  * @param thickness Thickness of the thin layered structure
116  * @param interface Flag of the thin layered structure
117  */
118  void setCoefficientsThinLayer ( const Real density, const Real young, const Real poisson, const Real thickness, const UInt interface );
119 
120  //! Setup: build the FE space and the ETA fespace
121  /*!
122  * @param mesh Computational mesh
123  * @param dOrder Degree of the finite element used
124  */
125  void setup( const meshPtr_Type& mesh, const std::string dOrder );
126 
127  //! Assembling matrices (mass and stiffness) which for this model are constant in a time dependent simulation
128  /*!
129  * @param timestep Value of the time step used
130  * @param bc Boundary conditions
131  * @param useBDF Boolean variable which specifies if BDF is used for the structure. By default it is false (therefore, Newmark is used).
132  */
133  void assemble_matrices ( const Real timestep, const Real beta, bcPtr_Type & bc, bool useBDF = false );
134 
135  //! Getter of the mass matrix whithout boundary conditions applied
136  /*!
137  * @return M_mass_no_bc mass matrix whithout boundary conditions applied
138  */
139  matrixPtr_Type const& mass_matrix_no_bc ( ) const { return M_mass_no_bc; };
140 
141  //! Getter of the stiffness matrix whithout boundary conditions applied
142  /*!
143  * @return M_stiffness_no_bc stiffness matrix whithout boundary conditions applied
144  */
145  matrixPtr_Type const& stiffness_matrix_no_bc ( ) const { return M_stiffness_no_bc; };
146 
147  //! Getter of the Jacobian matrix
148  /*!
149  * @return M_jacobian stiffness matrix
150  */
151  matrixPtr_Type const& jacobian ( ) const { return M_jacobian; };
152 
153  //! Getter of the standard FE space
154  /*!
155  * @return M_displacementFESpace FE space used for the displacement
156  *
157  */
158  const std::shared_ptr<FESpace<mesh_Type, map_Type> >& fespace ( ) const { return M_displacementFESpace; };
159 
160  //! Getter of the ETA FE space
161  /*!
162  * @return M_displacementFESpace_ETA ETA FE space used for the displacement
163  *
164  */
165  const std::shared_ptr<ETFESpace_displacement >& et_fespace ( ) const { return M_displacementFESpace_ETA; };
166 
167 private:
168 
169  // communicator
170  commPtr_Type M_comm;
171 
172  // Coefficients for the structure
173  Real M_density;
174  Real M_young;
175  Real M_poisson;
176 
177  // Lame coefficients computed on M_young and M_poisson
178  Real M_lambda;
179  Real M_mu;
180 
181  // order FE space displacement
182  std::string M_dOrder;
183 
184  // FE space
185  std::shared_ptr<FESpace<mesh_Type, map_Type> > M_displacementFESpace;
186 
187  // ET FE Space
188  std::shared_ptr<ETFESpace_displacement > M_displacementFESpace_ETA;
189 
190  // Matrices without bc applied
191  matrixPtr_Type M_mass_no_bc;
192  matrixPtr_Type M_stiffness_no_bc;
193 
194  // Matrices without bc applied thin layer
195  matrixPtr_Type M_mass_no_bc_thin;
196  matrixPtr_Type M_stiffness_no_bc_thin;
197 
198  // jacobian matrix
199  matrixPtr_Type M_jacobian;
200 
201  // Parameters used to deal with the thin layer
202  bool M_thinLayer;
203  Real M_thinLayerThickness;
204  Real M_thinLayerDensity;
205  Real M_thinLayerLameI;
206  Real M_thinLayerLameII;
207  UInt M_interfaceFlag;
208 
209 };
210 
211 } // end namespace LifeV
212 
213 #endif
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.
Definition: MapEpetra.cpp:394
double Real
Generic real data.
Definition: LifeV.hpp:175
class ETFESpace A light, templated version of the FESpace
Definition: ETFESpace.hpp:93
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191