1 /** @file gsMaterialMatrixLinear.h
2 
3     @brief Provides linear material matrices
4 
5     This file is part of the G+Smo library.
6 
7     This Source Code Form is subject to the terms of the Mozilla Public
8     License, v. 2.0. If a copy of the MPL was not distributed with this
9     file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11     Author(s):
12         H.M. Verhelst   (2019-..., TU Delft)
13         A. Mantzaflaris (2019-..., Inria)
14 */
15 
16 #pragma once
17 
18 #include <gsKLShell/gsMaterialMatrixBase.h>
19 #include <gsKLShell/gsMaterialMatrixBaseDim.h>
20 #include <gsKLShell/gsMaterialMatrixUtils.h>
21 #include <gsIO/gsOptionList.h>
22 #include <gsCore/gsFuncData.h>
23 
24 namespace gismo
25 {
26 
27 
28 /**
29  * @brief      This class defines a linear material
30  *
31  * @tparam     dim   The dimension of the problem (2 = planar, 3 = surface)
32  * @tparam     T     Real type
33  *
34  * @ingroup    KLShell
35  *
36  */
37 template <  short_t dim,
38             class T
39          >
40 class gsMaterialMatrixLinear : public gsMaterialMatrixBaseDim<dim,T>
41 {
42 public:
43     using Base = gsMaterialMatrixBaseDim<dim,T>;
44 
45     /**
46      * @brief      Constructor without deformed multipatch and density
47      *
48      * @param[in]  mp         Original geometry
49      * @param[in]  thickness  Thickness function
50      * @param[in]  pars       Vector with parameters (E, nu)
51      */
52     gsMaterialMatrixLinear(   const gsFunctionSet<T> & mp,
53                         const gsFunction<T> & thickness,
54                         const std::vector<gsFunction<T> *> &pars);
55 
56     /**
57      * @brief      Full constructor
58      *
59      * @param[in]  mp         Original geometry
60      * @param[in]  mp_def     Deformed geometry
61      * @param[in]  thickness  Thickness function
62      * @param[in]  pars       Vector with parameters (E, nu)
63      * @param[in]  Density    Density function
64      */
65     gsMaterialMatrixLinear(   const gsFunctionSet<T> & mp,
66                         const gsFunction<T> & thickness,
67                         const std::vector<gsFunction<T> *> &pars,
68                         const gsFunction<T> & Density);
69 
70     /// Destructor
gsMaterialMatrixLinear()71     gsMaterialMatrixLinear() { }
72 
73     /// See \ref gsMaterialMatrixBase for details
isMatIntegrated()74     inline enum MatIntegration isMatIntegrated() const {return MatIntegration::Constant; }
75 
76     /// See \ref gsMaterialMatrixBase for details
isVecIntegrated()77     inline enum MatIntegration isVecIntegrated() const {return MatIntegration::Constant; }
78 
79     /// See \ref gsMaterialMatrixBase for details
options()80     gsOptionList & options() {return m_options;}
81 
82     /// See \ref gsMaterialMatrixBase for details
setOptions(gsOptionList opt)83     void setOptions(gsOptionList opt) {m_options.update(opt,gsOptionList::addIfUnknown); }
84 
85     /// See \ref gsMaterialMatrixBase for details
86     void density_into(const index_t patch, const gsMatrix<T>& u, gsMatrix<T>& result) const;
87 
88     /// See \ref gsMaterialMatrixBase for details
89     void stretch_into(const index_t patch, const gsMatrix<T>& u, gsMatrix<T>& result) const;
90 
91     /// See \ref gsMaterialMatrixBase for details
92     void stretchDir_into(const index_t patch, const gsMatrix<T>& u, gsMatrix<T>& result) const;
93 
94     /// See \ref gsMaterialMatrixBase for details
95     void thickness_into(const index_t patch, const gsMatrix<T> & u, gsMatrix<T>& result) const;
96 
97     /// See \ref gsMaterialMatrixBase for details
98     void transform_into(const index_t patch, const gsMatrix<T> & u, gsMatrix<T>& result) const;
99 
100     /// See \ref gsMaterialMatrixBase for details
101     gsMatrix<T> eval3D_matrix(const index_t patch, const gsMatrix<T> & u, const gsMatrix<T>& z, enum MaterialOutput out = MaterialOutput::Generic) const;
102 
103     /// See \ref gsMaterialMatrixBase for details
104     gsMatrix<T> eval3D_vector(const index_t patch, const gsMatrix<T> & u, const gsMatrix<T>& z, enum MaterialOutput out = MaterialOutput::Generic) const;
105 
106     /// See \ref gsMaterialMatrixBase for details
107     gsMatrix<T> eval3D_pstress(const index_t patch, const gsMatrix<T> & u, const gsMatrix<T>& z, enum MaterialOutput out = MaterialOutput::Generic) const;
108 
109     /// See \ref gsMaterialMatrixBase for details
setParameters(const std::vector<gsFunction<T> * > & pars)110     void setParameters(const std::vector<gsFunction<T>*> &pars)
111     {
112         m_pars = pars;
113         m_numPars = m_pars.size();
114     }
115 
116     /// See \ref gsMaterialMatrixBase for details
117     void info() const;
118 
119 public:
120     /// Shared pointer for gsMaterialMatrixLinear
121     typedef memory::shared_ptr< gsMaterialMatrixLinear > Ptr;
122 
123     /// Unique pointer for gsMaterialMatrixLinear
124     typedef memory::unique_ptr< gsMaterialMatrixLinear > uPtr;
125 
126 protected:
127     /**
128      * @brief      Initializes the object.
129      *
130      * Initializes options, flags and defines the number of parameters
131      *
132      */
133     void _initialize();
134 
135     /**
136      * @brief      Sets default options
137      */
138     void _defaultOptions();
139 
140 protected:
141     /**
142      * @brief      Computes the linear material matrix entry with indices \a i \a j \a k \a l
143      *
144      * The entry is computed by \f$ \mathcal{C}^{ijkl} = \frac{2\lambda\mu}{\lambda+2\mu}a^{ij}a^{kl} + \mu (a^{ik}a^{jl} + a^{il}a^{jk})\f$ where \f$\lambda\f$ and \f$\mu\f$
145      * are the Lamé parameters and \f$a^{ij}=\mathbf{a}^i\cdot \mathbf{a}^j$ with \f$\mathbf{a}^i\f$ is the contravariant vector \f$ i \f$
146      *
147      * @param[in]  i,j,k,l  Indices
148      *
149      * @return     Cijkl
150      */
151     T _Cijkl  (const index_t i, const index_t j, const index_t k, const index_t l) const;
152 
153     /**
154      * @brief      Computes the linear force/moment entry with indices \a i \a j \a k \a l at height z
155      *
156      * Computes the thickness-integrated stress tensor, i.e. the normal force (0th thickness-moment) or the bending moment (1st thickness-moment).
157      * Sij is computed as \f$ \mathcal{C}^{ijkl} : \mathbf{E}_{ij} \f$ where \f$\mathbf{E}_{ij} = a_{ij} z b_{ij}\f$ with \f$ a_{ij}\f$ the in-plane metric and \f$b_{ij}\f$ the curvature.
158      * According to this definition, MaterialOutput out==VectorN would return the 0th moment, hence the \f$a_{ij}\f$-part would be relevant and for VectorM the \f$b_{ij}$ part is relevant.
159      * The term could also be integrated (out==Generalized).
160      *
161      * @param[in]  i,j   Indices
162      * @param[in]  z     Through-thickness coordinate
163      * @param[in]  out   Output specification
164      *
165      * @return     Sij
166      */
167     T _Sij    (const index_t i, const index_t j, const T z, enum MaterialOutput out) const;
168 
169     /**
170      * @brief      Computes the map, the metric quantities and the parameters on
171      *             specified points.
172      *
173      * @param[in]  patch  The patch index
174      * @param[in]  u      The in-plane point coordinates
175      */
176     void _computePoints(const index_t patch, const gsMatrix<T> & u, bool basis = true) const;
177 
178     /// Computes the stretch given deformation tensor C, into class members m_stretches and m_stretchDirs
179     void _computePStress(const gsMatrix<T> & C ) const;
180 
181     /// Computes the stretch given deformation tensor C, into a pair
182     std::pair<gsVector<T>,gsMatrix<T>> _evalPStress(const gsMatrix<T> & C ) const;
183 
184 protected:
185     // general
186     index_t m_numPars; // how many parameters for the material model?
187 
188     // constructor
189     using Base::m_patches;
190     using Base::m_defpatches;
191     const gsFunction<T> * m_thickness;
192     std::vector<gsFunction<T>* > m_pars;
193     const gsFunction<T> * m_density;
194 
195     mutable gsMatrix<T> m_Emat,m_Nmat,m_Tmat,m_rhomat;
196     mutable real_t m_lambda, m_mu, m_Cconstant;
197 
198     mutable gsMatrix<T>                 m_parmat;
199     mutable gsVector<T>                 m_parvals;
200 
201     mutable gsMatrix<T> m_pstress, m_pstressvec;
202 
203     // Geometric data point
204     using Base::m_map;
205     using Base::m_map_def;
206 
207     using Base::m_Acov_ori;
208     using Base::m_Acon_ori;
209     using Base::m_Acov_def;
210     using Base::m_Acon_def;
211     using Base::m_Bcov_ori;
212     using Base::m_Bcon_ori;
213     using Base::m_Bcov_def;
214     using Base::m_Bcon_def;
215     using Base::m_acov_ori;
216     using Base::m_acon_ori;
217     using Base::m_acov_def;
218     using Base::m_acon_def;
219     using Base::m_ncov_ori;
220     using Base::m_ncov_def;
221     using Base::m_Gcov_ori;
222     using Base::m_Gcon_ori;
223     using Base::m_Gcov_def;
224     using Base::m_Gcon_def;
225     using Base::m_Gcov_ori_L;
226     using Base::m_Gcov_def_L;
227     using Base::m_gcov_ori;
228     using Base::m_gcov_def;
229     using Base::m_gcon_ori;
230     using Base::m_gcon_def;
231     using Base::m_Acov_ori_mat;
232     using Base::m_Acon_ori_mat;
233     using Base::m_Acov_def_mat;
234     using Base::m_Acon_def_mat;
235     using Base::m_Bcov_ori_mat;
236     using Base::m_Bcov_def_mat;
237     using Base::m_acov_ori_mat;
238     using Base::m_acon_ori_mat;
239     using Base::m_acov_def_mat;
240     using Base::m_acon_def_mat;
241     using Base::m_ncov_ori_mat;
242     using Base::m_ncov_def_mat;
243 
244     using Base::m_stretches;
245     using Base::m_stretchvec;
246 
247     using Base::m_J0_sq;
248     using Base::m_J_sq;
249 
250 
251 
252     gsOptionList m_options;
253 
254 };
255 
256 } // namespace
257 
258 
259 #ifndef GISMO_BUILD_LIB
260 #include GISMO_HPP_HEADER(gsMaterialMatrixLinear.hpp)
261 #endif
262