1 /** @file gsElasticityFunctions.hpp
2
3 @brief Provides useful classes derived from gsFunction which can be used
4 for visualization or coupling.
5
6 This file is part of the G+Smo library.
7
8 This Source Code Form is subject to the terms of the Mozilla Public
9 License, v. 2.0. If a copy of the MPL was not distributed with this
10 file, You can obtain one at http://mozilla.org/MPL/2.0/.
11
12 Author(s):
13 A.Shamanskiy (2016 - ...., TU Kaiserslautern)
14 */
15
16 #pragma once
17
18 #include <gsElasticity/gsElasticityFunctions.h>
19 #include <gsCore/gsFuncData.h>
20 #include <gsAssembler/gsAssembler.h>
21
22 namespace gismo
23 {
24
25 template <class T>
linearElastic(const gsMatrix<T> & u,gsMatrix<T> & result) const26 void gsCauchyStressFunction<T>::linearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const
27 {
28 result.setZero(targetDim(),outputCols(u.cols()));
29 // evaluating the fields
30 gsMapData<T> mdGeo(NEED_GRAD_TRANSFORM);
31 mdGeo.points = u;
32 m_geometry->patch(m_patch).computeMap(mdGeo);
33 gsMapData<T> mdDisp(NEED_DERIV);
34 mdDisp.points = u;
35 m_displacement->patch(m_patch).computeMap(mdDisp);
36 // define temporary matrices here for efficieny
37 gsMatrix<T> I = gsMatrix<T>::Identity(m_dim,m_dim);
38 gsMatrix<T> sigma,eps,dispGrad;
39 // material parameters
40 T YM = m_options.getReal("YoungsModulus");
41 T PR = m_options.getReal("PoissonsRatio");
42 T lambda = YM * PR / ( ( 1. + PR ) * ( 1. - 2. * PR ) );
43 T mu = YM / ( 2. * ( 1. + PR ) );
44
45 for (index_t q = 0; q < u.cols(); ++q)
46 {
47 // linear strain tensor eps = (gradU+gradU^T)/2
48 if (mdGeo.jacobian(q).determinant() <= 0)
49 gsInfo << "Invalid domain parametrization: J = " << mdGeo.jacobian(q).determinant() <<
50 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
51 if (abs(mdGeo.jacobian(q).determinant()) > 1e-20)
52 dispGrad = mdDisp.jacobian(q)*(mdGeo.jacobian(q).cramerInverse());
53 else
54 dispGrad = gsMatrix<T>::Zero(m_dim,m_dim);
55 eps = (dispGrad + dispGrad.transpose())/2;
56 // Cauchy stress tensor
57 sigma = lambda*eps.trace()*I + 2*mu*eps;
58 saveStress(sigma,result,q);
59 }
60 }
61
62 template <class T>
nonLinearElastic(const gsMatrix<T> & u,gsMatrix<T> & result) const63 void gsCauchyStressFunction<T>::nonLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const
64 {
65 result.setZero(targetDim(),outputCols(u.cols()));
66 // evaluating the fields
67 gsMapData<T> mdGeo(NEED_GRAD_TRANSFORM);
68 mdGeo.points = u;
69 m_geometry->patch(m_patch).computeMap(mdGeo);
70 gsMapData<T> mdDisp(NEED_DERIV);
71 mdDisp.points = u;
72 m_displacement->patch(m_patch).computeMap(mdDisp);
73 // define temporary matrices here for efficieny
74 gsMatrix<T> I = gsMatrix<T>::Identity(m_dim,m_dim);
75 gsMatrix<T> S,sigma,F,C,E;
76 // material parameters
77 T YM = m_options.getReal("YoungsModulus");
78 T PR = m_options.getReal("PoissonsRatio");
79 T lambda = YM * PR / ( ( 1. + PR ) * ( 1. - 2. * PR ) );
80 T mu = YM / ( 2. * ( 1. + PR ) );
81
82 for (index_t q = 0; q < u.cols(); ++q)
83 {
84 // deformation gradient F = I + gradU*gradGeo^-1
85 if (mdGeo.jacobian(q).determinant() <= 0)
86 gsInfo << "Invalid domain parametrization: J = " << mdGeo.jacobian(q).determinant() <<
87 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
88 if (abs(mdGeo.jacobian(q).determinant()) > 1e-20)
89 F = I + mdDisp.jacobian(q)*(mdGeo.jacobian(q).cramerInverse());
90 else
91 F = I;
92 T J = F.determinant();
93 if (J <= 0)
94 gsInfo << "Invalid displacement field: J = " << J <<
95 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
96 // Second Piola-Kirchhoff stress tensor
97 if (material_law::law(m_options.getInt("MaterialLaw")) == material_law::saint_venant_kirchhoff)
98 {
99 // Green-Lagrange strain tensor E = 0.5(F^T*F-I)
100 E = (F.transpose() * F - I)/2;
101 S = lambda*E.trace()*I + 2*mu*E;
102 }
103 if (material_law::law(m_options.getInt("MaterialLaw")) == material_law::neo_hooke_ln)
104 {
105 // Right Cauchy Green strain, C = F'*F
106 C = F.transpose() * F;
107 S = (lambda*log(J)-mu)*(C.cramerInverse()) + mu*I;
108 }
109 if (material_law::law(m_options.getInt("MaterialLaw")) == material_law::neo_hooke_quad)
110 {
111 // Right Cauchy Green strain, C = F'*F
112 C = F.transpose() * F;
113 S = (lambda*(J*J-1)/2-mu)*(C.cramerInverse()) + mu*I;
114 }
115 // transformation to Cauchy stress
116 sigma = F*S*F.transpose()/J;
117 saveStress(sigma,result,q);
118 }
119 }
120
121 template <class T>
mixedLinearElastic(const gsMatrix<T> & u,gsMatrix<T> & result) const122 void gsCauchyStressFunction<T>::mixedLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const
123 {
124 result.setZero(targetDim(),outputCols(u.cols()));
125 // evaluating the fields
126 gsMapData<T> mdGeo(NEED_GRAD_TRANSFORM);
127 mdGeo.points = u;
128 m_geometry->patch(m_patch).computeMap(mdGeo);
129 gsMapData<T> mdDisp(NEED_DERIV);
130 mdDisp.points = u;
131 m_displacement->patch(m_patch).computeMap(mdDisp);
132 gsMatrix<T> presVals;
133 m_pressure->patch(m_patch).eval_into(u,presVals);
134 // define temporary matrices here for efficieny
135 gsMatrix<T> I = gsMatrix<T>::Identity(m_dim,m_dim);
136 gsMatrix<T> sigma,eps,dispGrad;
137 // material parameters
138 T YM = m_options.getReal("YoungsModulus");
139 T PR = m_options.getReal("PoissonsRatio");
140 T mu = YM / ( 2. * ( 1. + PR ) );
141
142 for (index_t q = 0; q < u.cols(); ++q)
143 {
144 // linear strain tensor eps = (gradU+gradU^T)/2
145 if (mdGeo.jacobian(q).determinant() <= 0)
146 gsInfo << "Invalid domain parametrization: J = " << mdGeo.jacobian(q).determinant() <<
147 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
148 if (abs(mdGeo.jacobian(q).determinant()) > 1e-20)
149 dispGrad = mdDisp.jacobian(q)*(mdGeo.jacobian(q).cramerInverse());
150 else
151 dispGrad = gsMatrix<T>::Zero(m_dim,m_dim);
152 eps = (dispGrad + dispGrad.transpose())/2;
153 // Cauchy stress tensor
154 sigma = presVals.at(q)*I + 2*mu*eps;
155 saveStress(sigma,result,q);
156 }
157 }
158
159 template <class T>
mixedNonLinearElastic(const gsMatrix<T> & u,gsMatrix<T> & result) const160 void gsCauchyStressFunction<T>::mixedNonLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const
161 {
162 result.setZero(targetDim(),outputCols(u.cols()));
163 // evaluating the fields
164 gsMapData<T> mdGeo(NEED_GRAD_TRANSFORM);
165 mdGeo.points = u;
166 m_geometry->patch(m_patch).computeMap(mdGeo);
167 gsMapData<T> mdDisp(NEED_DERIV);
168 mdDisp.points = u;
169 m_displacement->patch(m_patch).computeMap(mdDisp);
170 gsMatrix<T> presVals;
171 m_pressure->patch(m_patch).eval_into(u,presVals);
172 // define temporary matrices here for efficieny
173 gsMatrix<T> I = gsMatrix<T>::Identity(m_dim,m_dim);
174 gsMatrix<T> S,sigma,F,C,E;
175 // material parameters
176 T YM = m_options.getReal("YoungsModulus");
177 T PR = m_options.getReal("PoissonsRatio");
178 T mu = YM / ( 2. * ( 1. + PR ) );
179
180 for (index_t q = 0; q < u.cols(); ++q)
181 {
182 if (mdGeo.jacobian(q).determinant() <= 0)
183 gsInfo << "Invalid domain parametrization: J = " << mdGeo.jacobian(q).determinant() <<
184 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
185 // deformation gradient F = I + gradU*gradGeo^-1
186 if (abs(mdGeo.jacobian(q).determinant()) > 1e-20)
187 F = I + mdDisp.jacobian(q)*(mdGeo.jacobian(q).cramerInverse());
188 else
189 F = I;
190 T J = F.determinant();
191 if (J <= 0)
192 gsInfo << "Invalid displacement field: J = " << J <<
193 " at point (" << u.col(q).transpose() << ") of patch " << m_patch << std::endl;
194 // Second Piola-Kirchhoff stress tensor
195 C = F.transpose() * F;
196 S = (presVals.at(q)-mu)*(C.cramerInverse()) + mu*I;
197 // transformation to Cauchy stress
198 sigma = F*S*F.transpose()/J;
199 saveStress(sigma,result,q);
200 }
201 }
202
203 template <class T>
saveStress(const gsMatrix<T> & S,gsMatrix<T> & result,index_t q) const204 void gsCauchyStressFunction<T>::saveStress(const gsMatrix<T> & S, gsMatrix<T> & result, index_t q) const
205 {
206 switch (m_type)
207 {
208 case stress_components::von_mises :
209 {
210 if (m_geometry->parDim() == 2)
211 result(0,q) = sqrt( S(0,0)*S(0,0) + S(1,1)*S(1,1) - S(0,0)*S(1,1) + 3*S(0,1)*S(0,1) );
212 if (m_geometry->parDim() == 3)
213 result(0,q) = sqrt(0.5*( pow(S(0,0)-S(1,1),2) + pow(S(0,0)-S(2,2),2) + pow(S(1,1)-S(2,2),2) +
214 6 * (pow(S(0,1),2) + pow(S(0,2),2) + pow(S(1,2),2) ) ) );
215 return;
216 }
217 case stress_components::all_2D_vector : result(0,q) = S(0,0); result(1,q) = S(1,1); result(2,q) = S(0,1); return;
218 case stress_components::all_2D_matrix : result.middleCols(q*2,2) = S; return;
219 case stress_components::normal_3D_vector : result(0,q) = S(0,0); result(1,q) = S(1,1); result(2,q) = S(2,2); return;
220 case stress_components::shear_3D_vector : result(0,q) = S(0,1); result(1,q) = S(0,2); result(2,q) = S(1,2); return;
221 case stress_components::all_3D_matrix : result.middleCols(q*3,3) = S; return;
222 }
223 }
224
225
226 template <class T>
eval_into(const gsMatrix<T> & u,gsMatrix<T> & result) const227 void gsDetFunction<T>::eval_into(const gsMatrix<T> & u, gsMatrix<T> & result) const
228 {
229 result.resize(1,u.cols());
230
231 gsMapData<T> mappingData;
232 mappingData.points = u;
233 mappingData.flags = NEED_DERIV;
234 m_geo.patch(m_patch).computeMap(mappingData);
235
236 for (index_t i = 0; i < u.cols(); ++i)
237 result(0,i) = mappingData.jacobian(i).determinant();
238 }
239
240 template <class T>
eval_into(const gsMatrix<T> & u,gsMatrix<T> & result) const241 void gsFsiLoad<T>::eval_into(const gsMatrix<T> & u, gsMatrix<T> & result) const
242 {
243 result.setZero(targetDim(),u.cols());
244 // mapping points back to the parameter space via the reference configuration
245 gsMatrix<T> paramPoints;
246 m_geo.patch(m_patchGeo).invertPoints(u,paramPoints);
247 // evaluate reference geometry mapping at the param points
248 // NEED_GRAD_TRANSFORM for velocity gradients transformation from parametric to reference domain
249 gsMapData<T> mdGeo(NEED_GRAD_TRANSFORM);
250 mdGeo.points = paramPoints;
251 m_geo.patch(m_patchGeo).computeMap(mdGeo);
252 // evaluate velocity at the param points
253 // NEED_DERIV for velocity gradients
254 gsMapData<T> mdVel(NEED_DERIV);
255 mdVel.points = paramPoints;
256 m_vel.patch(m_patchVP).computeMap(mdVel);
257 // evaluate pressure at the quad points
258 gsMatrix<T> pressureValues;
259 m_pres.patch(m_patchVP).eval_into(paramPoints,pressureValues);
260 // evaluate ALE dispacement at the param points
261 // NEED_DERIV for gradients
262 gsMapData<T> mdALE(NEED_DERIV);
263 mdALE.points = paramPoints;
264 m_ale.patch(m_patchGeo).computeMap(mdALE);
265
266 gsMatrix<T> I = gsMatrix<T>::Identity(targetDim(),targetDim());
267 for (index_t p = 0; p < paramPoints.cols(); ++p)
268 {
269 // transform velocity gradients from parametric to reference
270 gsMatrix<T> physGradVel = mdVel.jacobian(p)*(mdGeo.jacobian(p).cramerInverse());
271 // ALE jacobian (identity + physical displacement gradient)
272 gsMatrix<T> physJacALE = I + mdALE.jacobian(p)*(mdGeo.jacobian(p).cramerInverse());
273 // inverse ALE jacobian
274 gsMatrix<T> invJacALE = physJacALE.cramerInverse();
275 // ALE stress tensor
276 gsMatrix<T> sigma = pressureValues.at(p)*I
277 - m_density*m_viscosity*(physGradVel*invJacALE +
278 invJacALE.transpose()*physGradVel.transpose());
279 // stress tensor pull back
280 gsMatrix<T> sigmaALE = physJacALE.determinant()*sigma*(invJacALE.transpose());
281
282 // normal length is the local measure
283 gsVector<T> normal;
284 outerNormal(mdGeo,p,m_sideGeo,normal);
285 result.col(p) = sigmaALE * normal / normal.norm();
286 }
287 }
288
289 } // namespace gismo ends
290