G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsMaterialMatrixLinear.hpp
Go to the documentation of this file.
1
16/*
17 To Do [updated 16-06-2020]:
18 - Make beta (compressible materials) and material parameters universal for all integration points over the thickness. So get them out of the dPsi functions etc and move them into the integration loops as global variables.
19
20*/
21
22#pragma once
23
25
26namespace gismo
27{
28
29template <short_t dim, class T>
36
37template <short_t dim, class T>
39 const gsFunctionSet<T> & mp,
40 const gsFunctionSet<T> & thickness
41 )
42 :
43 Base(&mp,&thickness,nullptr)
44{
46}
47
48template <short_t dim, class T >
50 const gsFunctionSet<T> & mp,
51 const gsFunctionSet<T> & thickness,
52 const gsFunctionSet<T> & YoungsModulus,
53 const gsFunctionSet<T> & PoissonsRatio
54 )
55 :
56 gsMaterialMatrixLinear(&mp,&thickness,YoungsModulus,PoissonsRatio,nullptr)
57{
58}
59
60template <short_t dim, class T >
62 const gsFunctionSet<T> & mp,
63 const gsFunctionSet<T> & thickness,
64 const gsFunctionSet<T> & YoungsModulus,
65 const gsFunctionSet<T> & PoissonsRatio,
66 const gsFunctionSet<T> & Density
67 )
68 :
69 gsMaterialMatrixLinear(&mp,&thickness,YoungsModulus,PoissonsRatio,&Density)
70{
71}
72
73template <short_t dim, class T >
75 const gsFunctionSet<T> * mp,
76 const gsFunctionSet<T> * thickness,
77 const gsFunctionSet<T> & YoungsModulus,
78 const gsFunctionSet<T> & PoissonsRatio,
79 const gsFunctionSet<T> * Density
80 )
81 :
82 Base(mp,thickness,Density)
83{
84 m_pars.resize(2);
85 this->setYoungsModulus(YoungsModulus);
86 this->setPoissonsRatio(PoissonsRatio);
88}
89
90template <short_t dim, class T >
92 const gsFunctionSet<T> & mp,
93 const gsFunctionSet<T> & thickness,
94 const std::vector<gsFunctionSet<T>*> &pars
95 )
96 :
97 gsMaterialMatrixLinear(&mp,&thickness,pars,nullptr)
98{
99}
100
101template <short_t dim, class T>
103 const gsFunctionSet<T> & mp,
104 const gsFunctionSet<T> & thickness,
105 const std::vector<gsFunctionSet<T>*> &pars,
106 const gsFunctionSet<T> & density
107 )
108 :
109 gsMaterialMatrixLinear(&mp,&thickness,pars,&density)
110{
111}
112
113template <short_t dim, class T >
115 const gsFunctionSet<T> & thickness,
116 const std::vector<gsFunctionSet<T>*> &pars
117 )
118 :
119 gsMaterialMatrixLinear(nullptr,&thickness,pars,nullptr)
120{
121}
122
123template <short_t dim, class T>
125 const gsFunctionSet<T> & thickness,
126 const std::vector<gsFunctionSet<T>*> &pars,
127 const gsFunctionSet<T> & density
128 )
129 :
130 gsMaterialMatrixLinear(nullptr,&thickness,pars,&density)
131{
132}
133
134template <short_t dim, class T>
136 const gsFunctionSet<T> * mp,
137 const gsFunctionSet<T> * thickness,
138 const std::vector<gsFunctionSet<T>*> &pars,
139 const gsFunctionSet<T> * density
140 )
141 :
142 Base(mp,thickness,density)
143{
144 GISMO_ASSERT(pars.size()==2,"Two material parameters should be assigned!");
145 this->setParameters(pars);
146 _initialize();
147}
148
149template <short_t dim, class T>
150std::ostream & gsMaterialMatrixLinear<dim,T>::print(std::ostream &os) const
151{
152 os <<"---------------------------------------------------------------------\n"
153 <<"---------------------Elastic Material Info---------------------------\n"
154 <<"---------------------------------------------------------------------\n\n";
155
156 os <<"Material model: \t";
157 os <<"Saint-Venant Kirchhoff";
158 os <<"\n";
159 os <<"---------------------------------------------------------------------\n\n";
160 return os;
161}
162
163template <short_t dim, class T>
165{
166 Base::defaultOptions();
167 m_options.addInt("NumGauss","Number of Gaussian points through thickness",4);
168}
169
170template <short_t dim, class T>
172{
173 // Set default options
174 this->defaultOptions();
175}
176
177template <short_t dim, class T>
179{
180 // gsInfo<<"TO DO: evaluate moments using thickness";
181 // Input: u in-plane points
182 // z matrix with, per point, a column with z integration points
183 // Output: (n=u.cols(), m=z.rows())
184 // [(u1,z1) (u2,z1) .. (un,z1), (u1,z2) .. (un,z2), .., (u1,zm) .. (un,zm)]
185
186 this->_computePoints(patch,u);
187 gsMatrix<T> result(9, u.cols() * z.rows());
188 result.setZero();
189
190 for (index_t k=0; k!=u.cols(); k++)
191 {
192 // Evaluate material properties on the quadrature point
193 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
194 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,k);
195
196 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
197 {
198 this->_getMetric(k, z(j, k) * m_data.mine().m_Tmat(0, k)); // on point i, on height z(0,j)
199
200 gsAsMatrix<T, Dynamic, Dynamic> C = result.reshapeCol(j*u.cols()+k,3,3);
201 /*
202 C = C1111, C1122, C1112
203 symm, C2222, C2212
204 symm, symm, C1212
205 */
206 C(0,0) = _Cijkl(0,0,0,0); // C1111
207 C(1,1) = _Cijkl(1,1,1,1); // C2222
208 C(2,2) = _Cijkl(0,1,0,1); // C1212
209 C(1,0) = C(0,1) = _Cijkl(0,0,1,1); // C1122
210 C(2,0) = C(0,2) = _Cijkl(0,0,0,1); // C1112
211 C(2,1) = C(1,2) = _Cijkl(1,1,0,1); // C2212
212 }
213 }
214 return result;
215}
216
217template <short_t dim, class T>
218gsMatrix<T> gsMaterialMatrixLinear<dim,T>::eval3D_matrix_C(const gsMatrix<T> & Cmat, const index_t patch, const gsVector<T> & u, const T z, enum MaterialOutput /*out*/) const
219{
220 // gsInfo<<"TO DO: evaluate moments using thickness";
221 // Input: u in-plane points
222 // z matrix with, per point, a column with z integration points
223 // Output: (n=u.cols(), m=z.rows())
224 // [(u1,z1) (u2,z1) .. (un,z1), (u1,z2) .. (un,z2), .., (u1,zm) .. (un,zm)]
225
226 this->_computePoints(patch,u);
227
228 gsMatrix<T> result(9, 1);
229 result.setZero();
230
231 // Evaluate material properties on the quadrature point
232 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
233 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,0);
234
235 // Get metric
236 this->_getMetricUndeformed(0, z * m_data.mine().m_Tmat(0, 0)); // on point i, on height z(0,j)
237 this->_getMetricDeformed(Cmat);
238
240 /*
241 C = C1111, C1122, C1112
242 symm, C2222, C2212
243 symm, symm, C1212
244 */
245 C(0,0) = _Cijkl(0,0,0,0); // C1111
246 C(1,1) = _Cijkl(1,1,1,1); // C2222
247 C(2,2) = _Cijkl(0,1,0,1); // C1212
248 C(1,0) = C(0,1) = _Cijkl(0,0,1,1); // C1122
249 C(2,0) = C(0,2) = _Cijkl(0,0,0,1); // C1112
250 C(2,1) = C(1,2) = _Cijkl(1,1,0,1); // C2212
251 return result;
252}
253
254template <short_t dim, class T>
256{
257 return gsMatrix<T>::Zero(27,u.cols()*z.rows());
258}
259
260// template <short_t dim, class T>
261// gsMatrix<T> gsMaterialMatrixLinear<dim,T>::eval3D_matrix(const index_t patch, const gsMatrix<T> & u, const gsMatrix<T> & z, enum MaterialOutput out) const
262// {
263// // gsInfo<<"TO DO: evaluate moments using thickness";
264// // Input: u in-plane points
265// // z matrix with, per point, a column with z integration points
266// // Output: (n=u.cols(), m=z.cols())
267// // [(u1,z1) (u2,z1) .. (un,z1), (u1,z2) .. (un,z2), .., (u1,zm) .. (un,zm)]
268// typedef std::function<gsMatrix<T> (gsMatrix<T> const &)> Stress_t;
269// Stress_t Sfun = [this]( const gsMatrix<T> & E )
270// {
271// gsMatrix<T> S(3,1);
272// S(0,0) = _Sij(0, 0, E);
273// S(1,0) = _Sij(1, 1, E);
274// S(2,0) = _Sij(0, 1, E);
275// return S;
276// };
277
278// this->_computePoints(patch,u);
279
280// gsMatrix<T> result(3, u.cols() * z.rows());
281// result.setZero();
282
283// gsMatrix<T> E;
284// for (index_t k=0; k!=u.cols(); k++)
285// {
286// // Evaluate material properties on the quadrature point
287// for (index_t v=0; v!=m_parmat.rows(); v++)
288// m_parvals.at(v) = m_parmat(v,k);
289
290// for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
291// {
292// this->_getMetric(k, z(j, k) * m_Tmat(0, k)); // on point i, on height z(0,j)
293
294// E = _E(z(j, k) * m_Tmat(0, k),out);
295
296// result(0, j * u.cols() + k) = ;
297// result(1, j * u.cols() + k) = _Sij(1, 1, E);
298// result(2, j * u.cols() + k) = _Sij(0, 1, E);
299// }
300// }
301
302// return result;
303// }
304
305template <short_t dim, class T>
307{
308 return eval3D_stress(patch,u,z,out);
309}
310
311template <short_t dim, class T>
313{
314 return eval3D_CauchyStress(patch,u,z,out);
315}
316
317template <short_t dim, class T>
318gsMatrix<T> gsMaterialMatrixLinear<dim,T>::eval3D_vector_C(const gsMatrix<T> & Cmat, const index_t patch, const gsVector<T> & u, const T z, enum MaterialOutput out) const
319{
320 // gsInfo<<"TO DO: evaluate moments using thickness";
321 // Input: u in-plane points
322 // z matrix with, per point, a column with z integration points
323 // Output: (n=u.cols(), m=z.cols())
324 // [(u1,z1) (u2,z1) .. (un,z1), (u1,z2) .. (un,z2), .., (u1,zm) .. (un,zm)]
325
326 this->_computePoints(patch,u);
327
328 gsMatrix<T> result(3, 1);
329 result.setZero();
330
331 gsMatrix<T> E;
332 // Evaluate material properties on the quadrature point
333 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
334 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,0);
335
336 // Get metric
337 this->_getMetricUndeformed(0, z * m_data.mine().m_Tmat(0, 0)); // on point i, on height z(0,j)
338 this->_getMetricDeformed(Cmat);
339
340 // GISMO_ASSERT( out == MaterialOutput::VectorN ||
341 // out == MaterialOutput::PStrainN ||
342 // out == MaterialOutput::PStressN ||
343 // out == MaterialOutput::TensionField ||
344 // out == MaterialOutput::StrainN,
345 // "MaterialOutput is wrong!");
346 GISMO_ASSERT(z==0.0,"This only works for z=0");
347
348 m_data.mine().m_Acov_def = m_data.mine().m_Gcov_def.block(0,0,2,2);
349 m_data.mine().m_Acon_def = m_data.mine().m_Gcon_def.block(0,0,2,2);
350 m_data.mine().m_Bcov_def.setZero();
351 m_data.mine().m_Bcon_def.setZero();
352
353 E = _E(z * m_data.mine().m_Tmat(0, 0),out);
354
355 result(0, 0) = _Sij(0, 0, E);
356 result(1, 0) = _Sij(1, 1, E);
357 result(2, 0) = _Sij(0, 1, E);
358 return result;
359}
360
361template <short_t dim, class T>
363{
364 gsMatrix<T> Smat = eval3D_stress(patch,u,z,out);
365
366 gsMatrix<T> result(2, u.cols() * z.rows());
367 result.setZero();
369 index_t colIdx;
370 std::pair<gsVector<T>,gsMatrix<T>> res;
371 for (index_t k=0; k!=u.cols(); k++)
372 {
373 // Evaluate material properties on the quadrature point
374 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
375 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,k);
376
377 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
378 {
379 colIdx = j * u.cols() + k;
380 S.setZero();
381 S(0,0) = Smat(0,colIdx);
382 S(1,1) = Smat(1,colIdx);
383 S(0,1) = S(1,0) = Smat(2,colIdx);
384 S(2,2) = 0;
385 res = this->_evalPStress(S);
386 result.col(j * u.cols() + k) = res.first;
387 }
388 }
389
390 return result;
391
392}
393
394template <short_t dim, class T>
396{
397 gsMatrix<T> Smat = eval3D_stress(patch,u,z,out);
398
399 gsMatrix<T> result(9, u.cols() * z.rows());
400 result.setZero();
402 index_t colIdx;
403 std::pair<gsVector<T>,gsMatrix<T>> res;
404 for (index_t k=0; k!=u.cols(); k++)
405 {
406 // Evaluate material properties on the quadrature point
407 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
408 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,k);
409
410 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
411 {
412 colIdx = j * u.cols() + k;
413 S.setZero();
414 S(0,0) = Smat(0,colIdx);
415 S(1,1) = Smat(1,colIdx);
416 S(0,1) = S(1,0) = Smat(2,colIdx);
417 S(2,2) = 0;
418 res = this->_evalPStress(S);
419 result.col(j * u.cols() + k) = res.second.reshape(9,1);
420 }
421 }
422
423 return result;
424
425}
426
427template <short_t dim, class T>
429{
430 gsMatrix<T> Smat = eval3D_CauchyStress(patch,u,z,out);
431
432 gsMatrix<T> result(2, u.cols() * z.rows());
433 result.setZero();
435 index_t colIdx;
436 std::pair<gsVector<T>,gsMatrix<T>> res;
437 for (index_t k=0; k!=u.cols(); k++)
438 {
439 // Evaluate material properties on the quadrature point
440 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
441 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,k);
442
443 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
444 {
445 colIdx = j * u.cols() + k;
446 S.setZero();
447 S(0,0) = Smat(0,colIdx);
448 S(1,1) = Smat(1,colIdx);
449 S(0,1) = S(1,0) = Smat(2,colIdx);
450 S(2,2) = 0;
451 res = this->_evalPStress(S);
452 result.col(j * u.cols() + k) = res.first;
453 }
454 }
455
456 return result;
457
458}
459
460template <short_t dim, class T>
462{
463 // Input: u in-plane points
464 // z matrix with, per point, a column with z integration points
465 // Output: (n=u.cols(), m=z.cols())
466 // [(u1,z1) (u2,z1) .. (un,z1), (u1,z2) .. (un,z2), .., (u1,zm) .. (un,zm)]
467
468 this->_computePoints(patch,u);
469
470 gsMatrix<T> result(3, u.cols() * z.rows());
471 result.setZero();
473 for (index_t k=0; k!=u.cols(); k++)
474 {
475 // Evaluate material properties on the quadrature point
476 for (index_t v=0; v!=m_data.mine().m_parmat.rows(); v++)
477 m_data.mine().m_parvals.at(v) = m_data.mine().m_parmat(v,k);
478
479 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
480 {
481 this->_getMetric(k, z(j, k) * m_data.mine().m_Tmat(0, k)); // on point i, on height z(0,j)
482
483 // E.setZero();
484 E.block(0,0,2,2) = _E(0,out);
485 result(0, j * u.cols() + k) = _Sij(0, 0, E);
486 result(1, j * u.cols() + k) = _Sij(1, 1, E);
487 result(2, j * u.cols() + k) = _Sij(0, 1, E);
488 }
489 }
490
491 return result;
492
493}
494
495template <short_t dim, class T>
497{
498 this->_computePoints(patch,u);
499 gsMatrix<T> result(1, u.cols() * z.rows());
500 result.setZero();
501 for (index_t k=0; k!=u.cols(); k++)
502 {
503 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
504 {
505 this->_getMetric(k, z(j, k) * m_data.mine().m_Tmat(0, k)); // on point i, on height z(0,j)
506 result(0,j * u.cols() + k) = math::sqrt(m_data.mine().m_J0_sq*1.0);
507 }
508 }
509 return result;
510}
511
512template <short_t dim, class T>
514{
515 this->_computePoints(patch,u);
516 gsMatrix<T> result = eval3D_stress(patch,u,z,out);
517 index_t colIdx;
518 T detF;
519 for (index_t k=0; k!=u.cols(); k++)
520 {
521 for( index_t j=0; j < z.rows(); ++j ) // through-thickness points
522 {
523 colIdx = j * u.cols() + k;
524 this->_getMetric(k, z(j, k) * m_data.mine().m_Tmat(0, k)); // on point i, on height z(0,j)
525 detF = math::sqrt(m_data.mine().m_J0_sq*1.0);
526 result.col(colIdx) /= detF;
527 }
528 }
529 return result;
530}
531
532/*
533 Available class members:
534 - m_data.mine().m_parvals
535 - m_metric
536 - m_metric_def
537*/
538template <short_t dim, class T>
539T gsMaterialMatrixLinear<dim,T>::_Cijkl(const index_t i, const index_t j, const index_t k, const index_t l) const
540{
541 GISMO_ENSURE( ( (i < 2) && (j < 2) && (k < 2) && (l < 2) ) , "Index out of range. i="<<i<<", j="<<j<<", k="<<k<<", l="<<l);
542
543
544 // --------------------------
545 // Saint Venant Kirchhoff
546 // --------------------------
547 T lambda, mu, Cconstant;
548
549 mu = m_data.mine().m_parvals.at(0) / (2.*(1. + m_data.mine().m_parvals.at(1)));
550 GISMO_ENSURE((1.-2.*m_data.mine().m_parvals.at(1)) != 0, "Division by zero in construction of SvK material parameters! (1.-2.*nu) = "<<(1.-2.*m_data.mine().m_parvals.at(1))<<"; m_data.mine().m_parvals.at(1) = "<<m_data.mine().m_parvals.at(1));
551 lambda = m_data.mine().m_parvals.at(0) * m_data.mine().m_parvals.at(1) / ( (1. + m_data.mine().m_parvals.at(1))*(1.-2.*m_data.mine().m_parvals.at(1))) ;
552 Cconstant = 2*lambda*mu/(lambda+2*mu);
553
554 return Cconstant*m_data.mine().m_Acon_ori(i,j)*m_data.mine().m_Acon_ori(k,l) + mu*(m_data.mine().m_Acon_ori(i,k)*m_data.mine().m_Acon_ori(j,l) + m_data.mine().m_Acon_ori(i,l)*m_data.mine().m_Acon_ori(j,k));
555}
556
557template <short_t dim, class T>
558T gsMaterialMatrixLinear<dim,T>::_Sij(const index_t i, const index_t j, const gsMatrix<T> & strain) const
559{
560 T result = _Cijkl(i,j,0,0) * strain(0,0) + _Cijkl(i,j,0,1) * strain(0,1)
561 + _Cijkl(i,j,1,0) * strain(1,0) + _Cijkl(i,j,1,1) * strain(1,1);
562 return result;
563}
564
565template <short_t dim, class T>
567{
568 gsWarn<<"This is dangerous, since it does not compute material parameters in advance!\n";
569 GISMO_ASSERT(strain.rows()==strain.cols() && strain.rows()==2, "Strain tensor must be 2x2!");
570 gsMatrix<T> result(2,2);
571 for (index_t i = 0; i!=2; i++)
572 for (index_t j = 0; j!=2; j++)
573 result(i,j) = _Cijkl(i,j,0,0) * strain(0,0) + _Cijkl(i,j,0,1) * strain(0,1)
574 + _Cijkl(i,j,1,0) * strain(1,0) + _Cijkl(i,j,1,1) * strain(1,1);
575
576 return result;
577}
578
579template <short_t dim, class T>
580gsMatrix<T> gsMaterialMatrixLinear<dim,T>::C(const gsMatrix<T> &) const
581{
582 gsMatrix<T> result(2,2);
583 result(0,0) = _Cijkl(0,0,0,0); // C1111
584 result(1,1) = _Cijkl(1,1,1,1); // C2222
585 result(2,2) = _Cijkl(0,1,0,1); // C1212
586 result(1,0) = result(0,1) = _Cijkl(0,0,1,1); // C1122
587 result(2,0) = result(0,2) = _Cijkl(0,0,0,1); // C1112
588 result(2,1) = result(1,2) = _Cijkl(1,1,0,1); // C2212
589 return result;
590}
591
592template <short_t dim, class T>
594{
595 gsMatrix<T> strain;
596 if ( out == MaterialOutput::VectorN ||
597 out == MaterialOutput::CauchyVectorN ||
598 out == MaterialOutput::StrainN ||
599 out == MaterialOutput::StressN ||
600 out == MaterialOutput::CauchyStressN ||
601 out == MaterialOutput::PStrainN ||
602 out == MaterialOutput::PStressN ||
603 out == MaterialOutput::PCauchyStressN ||
604 out == MaterialOutput::TensionField ||
605 out == MaterialOutput::StrainN ) // To be used with multiplyZ_into
606 strain = 0.5*(m_data.mine().m_Acov_def - m_data.mine().m_Acov_ori);
607 else if ( out == MaterialOutput::VectorM ||
608 out == MaterialOutput::CauchyVectorM ||
609 out == MaterialOutput::StressM ||
610 out == MaterialOutput::CauchyStressM ||
611 out == MaterialOutput::StrainM ||
612 out == MaterialOutput::PStrainM ||
613 out == MaterialOutput::PStressM ||
614 out == MaterialOutput::PCauchyStressM ||
615 out == MaterialOutput::StrainM )
616 strain = (m_data.mine().m_Bcov_ori - m_data.mine().m_Bcov_def);
617 else if ( out == MaterialOutput::Generic ||
618 out == MaterialOutput::Strain ||
619 out == MaterialOutput::Stress ||
620 out == MaterialOutput::PStress ) // To be used with multiplyLinZ_into or integrateZ_into
621 strain = 0.5*(m_data.mine().m_Gcov_def.block(0,0,2,2) - m_data.mine().m_Gcov_ori.block(0,0,2,2));
622 else
623 GISMO_ERROR("Output type MaterialOutput::" + std::to_string((short_t)(out)) + " not understood. See gsMaterialMatrixUtils.h");
624
625 return strain;
626}
627
628// template <short_t dim, class T>
629// gsMatrix<T> gsMaterialMatrixLinear<dim,T>::_E(const gsMatrix<T> & C, const T z, enum MaterialOutput out) const
630// {
631// gsMatrix<T> strain;
632// if (out == MaterialOutput::VectorN || out == MaterialOutput::PStrainN || out == MaterialOutput::PStressN || out == MaterialOutput::TensionField || out == MaterialOutput::StrainN) // To be used with multiplyZ_into
633// strain = 0.5*(C - m_data.mine().m_Acov_ori);
634// else if (out == MaterialOutput::VectorM || out == MaterialOutput::PStrainM || out == MaterialOutput::PStressM || out == MaterialOutput::TensionField || out == MaterialOutput::StrainM) // To be used with multiplyZ_into
635// strain = (C - m_data.mine().m_Bcov_def);
636// else if (out == MaterialOutput::Generic || out == MaterialOutput::Strain) // To be used with multiplyLinZ_into or integrateZ_into
637// strain = 0.5*(C - m_data.mine().m_Gcov_ori);
638// else
639// GISMO_ERROR("Output type is not VectorN, PStrainN, PStressN, VectorM, PStrainM. PStressM, TensionField or Generic!");
640
641// return strain;
642// }
643
644//--------------------------------------------------------------------------------------------------------------------------------------
645
646namespace internal
647{
648
652template<short_t d, class T>
653class gsXml< gsMaterialMatrixLinear<d,T> >
654{
655private:
656 gsXml() { }
657 typedef gsMaterialMatrixLinear<d,T> Object;
658
659public:
660 GSXML_COMMON_FUNCTIONS(gsMaterialMatrixLinear<TMPLA2(d,T)>);
661 static std::string tag () { return "MaterialMatrix"; }
662 static std::string type () { return "Linear" + to_string(d); }
663
664 GSXML_GET_POINTER(Object);
665
666 // static Object * get(gsXmlNode * node)
667 // {
668 // Object result;
669 // get_into(node, result);
670 // return result.clone().release();
671 // }
672
673 static void get_into(gsXmlNode * node,Object & obj)
674 {
675 obj = getMaterialMatrixFromXml< Object >( node );
676 }
677
678 static gsXmlNode * put (const Object & obj,
679 gsXmlTree & data)
680 {
681 return putMaterialMatrixToXml< Object >( obj,data );
682 // GISMO_NO_IMPLEMENTATION;
683 // return putGeometryToXml(obj,data);
684 }
685};
686
687}// namespace internal
688
689} // end namespace
Creates a mapped object or data pointer to a matrix without copying data.
Definition gsAsMatrix.h:32
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition gsFunctionSet.h:219
This class defines the base class for material matrices.
Definition gsMaterialMatrixBaseDim.h:37
virtual void setParameters(const std::vector< function_ptr > &pars)
Sets the material parameters.
Definition gsMaterialMatrixBase.h:684
This class defines a linear material.
Definition gsMaterialMatrixLinear.h:38
gsMatrix< T > eval3D_vector(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:306
gsMatrix< T > eval3D_matrix_C(const gsMatrix< T > &Cmat, const index_t patch, const gsVector< T > &u, const T z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:218
virtual gsMatrix< T > eval3D_CauchyStress(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:513
void _initialize()
Initializes the object.
Definition gsMaterialMatrixLinear.hpp:171
void setYoungsModulus(const gsFunctionSet< T > &YoungsModulus) override
Sets the YoungsModulus.
Definition gsMaterialMatrixLinear.h:218
T _Sij(const index_t i, const index_t j, const gsMatrix< T > &z) const
Computes the linear material matrix entry with indices i j k l.
Definition gsMaterialMatrixLinear.hpp:558
gsMatrix< T > eval3D_pstressDir(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:395
void defaultOptions() override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:164
gsMatrix< T > eval3D_detF(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:496
std::ostream & print(std::ostream &os) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:150
gsMatrix< T > eval3D_CauchyVector(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:312
gsMatrix< T > _E(const T z, enum MaterialOutput out) const
Computes the strain tensor.
Definition gsMaterialMatrixLinear.hpp:593
gsMaterialMatrixLinear()
Empty constructor.
Definition gsMaterialMatrixLinear.hpp:30
void setPoissonsRatio(const gsFunctionSet< T > &PoissonsRatio) override
Sets the Poisson's Ratio.
Definition gsMaterialMatrixLinear.h:224
gsMatrix< T > eval3D_stress(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:461
gsMatrix< T > eval3D_vector_C(const gsMatrix< T > &Cmat, const index_t patch, const gsVector< T > &u, const T z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:318
gsMatrix< T > eval3D_pstress(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:362
gsMatrix< T > eval3D_dmatrix(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:255
gsMatrix< T > eval3D_matrix(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:178
gsMatrix< T > eval3D_CauchyPStress(const index_t patch, const gsMatrix< T > &u, const gsMatrix< T > &z, enum MaterialOutput out=MaterialOutput::Generic) const override
See gsMaterialMatrixBase for details.
Definition gsMaterialMatrixLinear.hpp:428
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
gsAsMatrix< T, Dynamic, Dynamic > reshapeCol(index_t c, index_t n, index_t m)
Returns column c of the matrix resized to n x m matrix This function assumes that the matrix is size ...
Definition gsMatrix.h:231
gsAsMatrix< T, Dynamic, Dynamic > reshape(index_t n, index_t m)
Returns the matrix resized to n x m matrix (data is not copied) This function assumes that the matrix...
Definition gsMatrix.h:221
T at(index_t i) const
Returns the i-th element of the vectorization of the matrix.
Definition gsMatrix.h:211
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
MaterialOutput
This class describes the output type.
Definition gsMaterialMatrixUtils.h:99
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define gsWarn
Definition gsDebug.h:50
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
std::string to_string(const unsigned &i)
Helper to convert small unsigned to string.
Definition gsXml.cpp:74
The G+Smo namespace, containing all definitions for the library.