G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsMaterialMatrixEval.hpp
Go to the documentation of this file.
1
16#pragma once
17
20
21namespace gismo
22{
23
24template <class T, enum MaterialOutput out>
25gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( const gsMaterialMatrixContainer<T> & materialMatrices,
26 const gsFunctionSet<T> * deformed,
27 const gsMatrix<T> z)
28:
29m_materialMatrices(materialMatrices),
30m_deformed(deformed),
31m_z(z)
32{
33 for (index_t p = 0; p!=deformed->nPieces(); ++p)
34 GISMO_ASSERT(materialMatrices.piece(p)->initialized(),"Material matrix "<<p<<" is incomplete!");
35
36 this->_makePieces();
37}
38
39template <class T, enum MaterialOutput out>
40gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( gsMaterialMatrixBase<T> * materialMatrix,
41 const gsFunctionSet<T> * deformed,
42 const gsMatrix<T> z)
43:
44m_materialMatrices(deformed->nPieces()),
45m_deformed(deformed),
46m_z(z)
47{
48 GISMO_ASSERT(materialMatrix->initialized(),"Material matrix is incomplete!");
49 for (index_t p = 0; p!=deformed->nPieces(); ++p)
50 m_materialMatrices.set(p,materialMatrix);
51 this->_makePieces();
52}
53
54template <class T, enum MaterialOutput out>
55gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( typename gsMaterialMatrixBase<T>::uPtr & materialMatrix,
56 const gsFunctionSet<T> * deformed,
57 const gsMatrix<T> z)
58
59:
60gsMaterialMatrixEval(materialMatrix.get(),deformed,z)
61{}
62
63template <class T, enum MaterialOutput out>
64gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( const gsMaterialMatrixContainer<T> & materialMatrices,
65 const gsFunctionSet<T> * undeformed,
66 const gsFunctionSet<T> * deformed,
67 const gsMatrix<T> z)
68:
69m_materialMatrices(materialMatrices),
70m_deformed(deformed),
71m_z(z)
72{
73 for (index_t p = 0; p!=deformed->nPieces(); ++p)
74 GISMO_ASSERT(materialMatrices.piece(p)->initialized(),"Material matrix "<<p<<" is incomplete!");
75
76 this->_makePieces(undeformed);
77}
78
79template <class T, enum MaterialOutput out>
80gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( gsMaterialMatrixBase<T> * materialMatrix,
81 const gsFunctionSet<T> * undeformed,
82 const gsFunctionSet<T> * deformed,
83 const gsMatrix<T> z)
84:
85m_materialMatrices(deformed->nPieces()),
86m_deformed(deformed),
87m_z(z)
88{
89 GISMO_ASSERT(materialMatrix->initialized(),"Material matrix is incomplete!");
90 for (index_t p = 0; p!=deformed->nPieces(); ++p)
91 m_materialMatrices.set(p,materialMatrix);
92 this->_makePieces(undeformed);
93}
94
95template <class T, enum MaterialOutput out>
96gsMaterialMatrixEval<T,out>::gsMaterialMatrixEval( typename gsMaterialMatrixBase<T>::uPtr & materialMatrix,
97 const gsFunctionSet<T> * undeformed,
98 const gsFunctionSet<T> * deformed,
99 const gsMatrix<T> z)
100:
101gsMaterialMatrixEval(materialMatrix.get(),undeformed,deformed,z)
102{}
103
104template <class T, enum MaterialOutput out>
105gsMaterialMatrixEval<T,out>::~gsMaterialMatrixEval()
106{
107 freeAll(m_pieces);
108}
109
110template <class T, enum MaterialOutput out>
111void gsMaterialMatrixEval<T, out>::_makePieces()
112{
113 m_pieces.resize(m_deformed->nPieces());
114 for (size_t p = 0; p != m_pieces.size(); ++p)
115 m_pieces.at(p) = new gsMaterialMatrixEvalSingle<T, out>(p, m_materialMatrices.piece(p), m_deformed, m_z);
116}
117
118template <class T, enum MaterialOutput out>
119void gsMaterialMatrixEval<T, out>::_makePieces(const gsFunctionSet<T> * undeformed)
120{
121 m_pieces.resize(m_deformed->nPieces());
122 for (size_t p = 0; p != m_pieces.size(); ++p)
123 m_pieces.at(p) = new gsMaterialMatrixEvalSingle<T, out>(p, m_materialMatrices.piece(p), undeformed, m_deformed, m_z);
124}
125
126template <class T, enum MaterialOutput out>
128 gsMaterialMatrixBase<T> * materialMatrix,
129 const gsFunctionSet<T> * deformed,
130 const gsMatrix<T> z
131 )
132:
133m_pIndex(patch),
134m_materialMat(materialMatrix),
135m_z(z)
136{
137 m_materialMat->setDeformed(deformed);
138
139 if (m_z.cols()==0 || m_z.rows()==0)
140 {
141 m_z.resize(1,1);
142 m_z.setZero();
143 }
144 GISMO_ASSERT(z.cols()==1,"Z coordinates should be provided row-wise in one column");
145
146 // m_materialMat = new gsMaterialMatrix(materialMatrix);
147}
148
149template <class T, enum MaterialOutput out>
151 gsMaterialMatrixBase<T> * materialMatrix,
152 const gsFunctionSet<T> * /*undeformed*/,
153 const gsFunctionSet<T> * deformed,
154 const gsMatrix<T> z
155 )
156:
157m_pIndex(patch),
158m_materialMat(materialMatrix),
159m_z(z)
160{
161 m_materialMat->setUndeformed(deformed);
162 m_materialMat->setDeformed(deformed);
163
164 if (m_z.cols()==0 || m_z.rows()==0)
165 {
166 m_z.resize(1,1);
167 m_z.setZero();
168 }
169 GISMO_ASSERT(z.cols()==1,"Z coordinates should be provided row-wise in one column");
170
171 // m_materialMat = new gsMaterialMatrix(materialMatrix);
172}
173
174template <class T, enum MaterialOutput out>
176{
178// #pragma omp critical (gsMaterialMatrixEvalSingle_eval_into)
179 this->eval_into_impl<out>(u,result);
180}
181
182template <class T, enum MaterialOutput out>
183template <enum MaterialOutput _out>
184typename std::enable_if<_out==MaterialOutput::Density, void>::type
186{
187 m_materialMat->density_into(m_pIndex,u,result);
188}
189
190template <class T, enum MaterialOutput out>
191template <enum MaterialOutput _out>
192typename std::enable_if<_out==MaterialOutput::VectorN ||
193 _out==MaterialOutput::VectorM ||
194 _out==MaterialOutput::Generic, void>::type
196{
197 result = m_materialMat->eval3D_vector(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
198}
199
200template <class T, enum MaterialOutput out>
201template <enum MaterialOutput _out>
202typename std::enable_if<_out==MaterialOutput::CauchyVectorN ||
203 _out==MaterialOutput::CauchyVectorM, void>::type
204gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
205{
206 result = m_materialMat->eval3D_CauchyVector(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
207}
208
209
210template <class T, enum MaterialOutput out>
211template <enum MaterialOutput _out>
212typename std::enable_if< _out==MaterialOutput::MatrixA || _out==MaterialOutput::MatrixB
213 || _out==MaterialOutput::MatrixC || _out==MaterialOutput::MatrixD, void>::type
214gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
215{
216 result = m_materialMat->eval3D_matrix(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
217}
218
219template <class T, enum MaterialOutput out>
220template <enum MaterialOutput _out>
221typename std::enable_if<_out==MaterialOutput::Stretch, void>::type
222gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
223{
224 result = m_materialMat->eval3D_pstretch(m_pIndex,u,m_z.replicate(1,u.cols()));
225}
226
227template <class T, enum MaterialOutput out>
228template <enum MaterialOutput _out>
229typename std::enable_if<_out==MaterialOutput::PStress ||
230 _out==MaterialOutput::PStressN ||
231 _out==MaterialOutput::PStressM , void>::type
232gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
233{
234 result = m_materialMat->eval3D_pstress(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
235}
236
237template <class T, enum MaterialOutput out>
238template <enum MaterialOutput _out>
239typename std::enable_if<_out==MaterialOutput::PCauchyStressN ||
240 _out==MaterialOutput::PCauchyStressM, void>::type
241gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
242{
243 result = m_materialMat->eval3D_CauchyPStress(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
244}
245
246template <class T, enum MaterialOutput out>
247template <enum MaterialOutput _out>
248typename std::enable_if<_out==MaterialOutput::PStrainN || _out==MaterialOutput::PStrainM, void>::type
249gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
250{
251 result = m_materialMat->eval3D_pstrain(m_pIndex,u,m_z.replicate(1,u.cols()));
252}
253
254template <class T, enum MaterialOutput out>
255template <enum MaterialOutput _out>
256typename std::enable_if<_out==MaterialOutput::StretchDir, void>::type
257gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
258{
259 result = m_materialMat->eval3D_pstretchDir(m_pIndex,u,m_z.replicate(1,u.cols()));
260}
261
262template <class T, enum MaterialOutput out>
263template <enum MaterialOutput _out>
264typename std::enable_if<_out==MaterialOutput::PStressDir, void>::type
265gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
266{
267 result = m_materialMat->eval3D_pstressDir(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
268}
269
270template <class T, enum MaterialOutput out>
271template <enum MaterialOutput _out>
272typename std::enable_if<_out==MaterialOutput::StretchTransform, void>::type
273gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
274{
275 result = m_materialMat->eval3D_pstretchTransform(m_pIndex,u,m_z.replicate(1,u.cols()));
276}
277
278template <class T, enum MaterialOutput out>
279template <enum MaterialOutput _out>
280typename std::enable_if<_out==MaterialOutput::PStressTransform, void>::type
281gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
282{
283 result = m_materialMat->eval3D_pstressTransform(m_pIndex,u,m_z.replicate(1,u.cols()));
284}
285
286template <class T, enum MaterialOutput out>
287template <enum MaterialOutput _out>
288typename std::enable_if<_out==MaterialOutput::Spec2CovTransform, void>::type
289gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
290{
291 result = m_materialMat->eval3D_spec2cov(m_pIndex,u,m_z.replicate(1,u.cols()));
292}
293
294template <class T, enum MaterialOutput out>
295template <enum MaterialOutput _out>
296typename std::enable_if<_out==MaterialOutput::Spec2ConTransform, void>::type
297gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
298{
299 result = m_materialMat->eval3D_spec2con(m_pIndex,u,m_z.replicate(1,u.cols()));
300}
301
302template <class T, enum MaterialOutput out>
303template <enum MaterialOutput _out>
304typename std::enable_if<_out==MaterialOutput::Cov2CartTransform, void>::type
305gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
306{
307 result = m_materialMat->eval3D_cov2cart(m_pIndex,u,m_z.replicate(1,u.cols()));
308}
309
310template <class T, enum MaterialOutput out>
311template <enum MaterialOutput _out>
312typename std::enable_if<_out==MaterialOutput::Con2CartTransform, void>::type
313gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
314{
315 result = m_materialMat->eval3D_con2cart(m_pIndex,u,m_z.replicate(1,u.cols()));
316}
317
318template <class T, enum MaterialOutput out>
319template <enum MaterialOutput _out>
320typename std::enable_if<_out==MaterialOutput::TensionField, void>::type
321gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
322{
323 result = m_materialMat->eval3D_tensionfield(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
324}
325
326template <class T, enum MaterialOutput out>
327template <enum MaterialOutput _out>
328typename std::enable_if<_out==MaterialOutput::Theta, void>::type
329gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
330{
331 result = m_materialMat->eval3D_theta(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
332}
333
334template <class T, enum MaterialOutput out>
335template <enum MaterialOutput _out>
336typename std::enable_if<_out==MaterialOutput::Gamma, void>::type
337gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
338{
339 result = m_materialMat->eval3D_gamma(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
340}
341
342template <class T, enum MaterialOutput out>
343template<enum MaterialOutput _out>
344typename std::enable_if<_out==MaterialOutput::Strain ||
345 _out==MaterialOutput::StrainN ||
346 _out==MaterialOutput::StrainM , void>::type
347gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
348{
349 result = m_materialMat->eval3D_strain(m_pIndex,u,m_z.replicate(1,u.cols()));
350}
351
352template <class T, enum MaterialOutput out>
353template<enum MaterialOutput _out>
354typename std::enable_if<_out==MaterialOutput::Stress ||
355 _out==MaterialOutput::StressN ||
356 _out==MaterialOutput::StressM , void>::type
357gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
358{
359 result = m_materialMat->eval3D_stress(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
360}
361
362template <class T, enum MaterialOutput out>
363template<enum MaterialOutput _out>
364typename std::enable_if<_out==MaterialOutput::CauchyStress ||
365 _out==MaterialOutput::CauchyStressN ||
366 _out==MaterialOutput::CauchyStressM , void>::type
367gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
368{
369 result = m_materialMat->eval3D_CauchyStress(m_pIndex,u,m_z.replicate(1,u.cols()),_out);
370}
371
372
373template <class T, enum MaterialOutput out>
374template <enum MaterialOutput _out>
375typename std::enable_if<_out==MaterialOutput::Thickness, void>::type
376gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
377{
378 m_materialMat->thickness_into(m_pIndex,u,result);
379}
380
381template <class T, enum MaterialOutput out>
382template <enum MaterialOutput _out>
383typename std::enable_if<_out==MaterialOutput::Parameters, void>::type
384gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
385{
386 m_materialMat->parameters_into(m_pIndex,u,result);
387}
388
389template <class T, enum MaterialOutput out>
390template <enum MaterialOutput _out>
391typename std::enable_if<_out==MaterialOutput::Deformation, void>::type
392gsMaterialMatrixEvalSingle<T,out>::eval_into_impl(const gsMatrix<T>& u, gsMatrix<T>& result) const
393{
394 result = m_materialMat->eval3D_deformation(m_pIndex,u,m_z.replicate(1,u.cols()));
395}
396
397
398} // end namespace
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 gsMaterialMatrixBase.h:33
This class serves as the evaluator of material matrices, based on gsMaterialMatrixBase.
Definition gsMaterialMatrixEval.h:107
void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Implementation of eval_into, see gsFunction.
Definition gsMaterialMatrixEval.hpp:175
std::enable_if< _out==MaterialOutput::Density, void >::type eval_into_impl(const gsMatrix< T > &u, gsMatrix< T > &result) const
Specialisation of eval_into for densities.
Definition gsMaterialMatrixEval.hpp:185
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
#define index_t
Definition gsConfig.h:32
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides the Gauss-Legendre quadrature rule.
Provides an evaluator for material matrices for thin shells.
The G+Smo namespace, containing all definitions for the library.
void freeAll(It begin, It end)
Frees all pointers in the range [begin end)
Definition gsMemory.h:312