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