dune-composites (unstable)

getStress.hh
1// -*- tab-width: 2; indent-tabs-mode: nil -*-
2// vi: set et ts=2 sw=2 sts=2:
3#ifndef GET_STRESS_HH
4#define GET_STRESS_HH
5
6#include <dune/geometry/type.hh>
7#include <dune/geometry/referenceelements.hh>
8#include <dune/geometry/quadraturerules.hh>
9
10namespace Dune {
11 namespace PDELab {
12
13 template<typename GV, typename PROBLEM, typename DGF, int dofel>
14 class getStressLinear :
23 public NumericalJacobianApplyVolume<getStressLinear<GV,PROBLEM,DGF,dofel> >,
24 public FullVolumePattern,
25 public LocalOperatorDefaultFlags,
26 public InstationaryLocalOperatorDefaultMethods<double>,
27 public NumericalJacobianVolume<getStressLinear<GV,PROBLEM,DGF, dofel> >
28 {
29 public:
30 // pattern assembly flags
31 enum { doPatternVolume = false };
32
33 // residual assembly flags
34 enum { doAlphaVolume = true };
35
36 //Constructor
37 getStressLinear (const GV& gv_, const PROBLEM& prob_, const DGF& dgf_, double& maxDisp_)
38 : gv(gv_)
39 , problem(prob_)
40 , dgf(dgf_)
41 , maxDisp(maxDisp_)
42 {}
43 //
44 // alpha_volume for getStress
45 //
46 // Volume Integral Depending on Test and Ansatz Functions
47 template<typename EG, typename LFSU, typename X, typename LFSV, typename R>
48 void alpha_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv, R& r) const
49 {
50 const int dim = 3;
51
52 //Unwrap function spaces
53 typedef typename LFSU::template Child<0>::Type LFS1;
54 typedef typename LFS1::Traits::FiniteElementType::Traits::LocalBasisType::Traits::JacobianType JacobianType;
55 //Unwrap displacements
56 const auto& lfs1 = lfsu.template child<0>();
57 const auto& lfs2 = lfsu.template child<1>();
58 const auto& lfs3 = lfsu.template child<2>();
59 //Unwrap stress
60 const auto& lfsv1 = lfsv.template child<0>();
61 const auto& lfsv2 = lfsv.template child<1>();
62 const auto& lfsv3 = lfsv.template child<2>();
63 const auto& lfsv4 = lfsv.template child<3>();
64 const auto& lfsv5 = lfsv.template child<4>();
65 const auto& lfsv6 = lfsv.template child<5>();
66
67 //Evaluate elasticity tensor
68 Dune::FieldMatrix<double,6,6> C(0.0),Cthermal(0.0);
69 problem.evaluate(eg.entity(),C);
70 problem.evaluate_strain(eg.entity(),Cthermal);
71
72 //Residual heat stresses
73 Dune::FieldVector<double,6> delta(0.0);
74 if(problem.getMaterialType(eg.entity()) == 0)//resin
75 problem.evaluateHeat(delta,true);
76 else
77 problem.evaluateHeat(delta,false); //ply
78
79 //Calculate input heat stress
80 Dune::FieldVector<double,6> appliedStress(0.0);
81 Cthermal.mv(delta, appliedStress);
82
83 int nodes_per_element = lfsv1.size();
84 int dof = lfs1.size();
85
86 // select quadrature rule
87 for (int iter = 0; iter < nodes_per_element; iter++)
88 {
89 Dune::FieldVector<double,1> num_elem(1.0);
90 Dune::FieldVector<double,dim> pt = eg.geometry().local(eg.geometry().corner(iter));
91 Dune::FieldVector<double,dim> pt_cen = eg.geometry().local(eg.geometry().center());
92
93 auto theta = problem.getThetaFct(eg.entity(),pt);
94 auto Cijkl = getMatrix(C,theta);
95
96 // Evaluate Jacobian
97 std::vector<JacobianType> js(dof);
98 lfs1.finiteElement().localBasis().evaluateJacobian(pt,js);
99
100 // Transform gradient to real element and evaluate at integration points
101 const auto jac = eg.geometry().jacobianInverseTransposed(pt);
102 std::vector<Dune::FieldVector<double,dim> > gradphi(dof);
103 for (int i=0; i < dof; i++){
104 gradphi[i] = 0.0;
105 jac.umv(js[i][0],gradphi[i]);
106 }
107 Dune::FieldMatrix<double,6,dofel> B(0.0);
108 for (int i = 0; i < dof; i++)
109 {
110 B[0][i] = gradphi[i][0]; // E11
111 B[1][i + dof] = gradphi[i][1]; // E22
112 B[2][i + 2 * dof] = gradphi[i][2]; // E33
113 B[3][i + dof] = gradphi[i][2]; B[3][i + 2 * dof] = gradphi[i][1]; // E23
114 B[4][i] = gradphi[i][2]; B[4][i + 2 * dof] = gradphi[i][0]; // E13
115 B[5][i] = gradphi[i][1]; B[5][i + dof] = gradphi[i][0]; // E12
116 }
117
118 Dune::FieldVector<double,6> e(0.0),s(0.0);
119 Dune::FieldVector<double,dofel> xlc(0.0);
120 for(int i = 0; i < dof; i++)
121 {
122 xlc[i] = x(lfs1,i);
123 xlc[i+dof] = x(lfs2,i);
124 xlc[i+2*dof] = x(lfs3,i);
125 }
126
127 B.mv(xlc,e); // e = B * x
128 Cijkl.mv(e,s); // s = Cijkl * e
129
130 s -= appliedStress;
131
132 // Evaluate Jacobian
133 std::vector<JacobianType> js_cen(dof);
134 lfs1.finiteElement().localBasis().evaluateJacobian(pt_cen,js_cen);
135
136 // Transform gradient to real element and evaluate at integration points
137 const auto jac_cen = eg.geometry().jacobianInverseTransposed(pt_cen);
138 std::vector<Dune::FieldVector<double,dim> > gradphi_cen(dof);
139 for (int i=0; i < dof; i++)
140 {
141 gradphi_cen[i] = 0.0;
142 jac_cen.umv(js_cen[i][0],gradphi_cen[i]);
143 }
144 Dune::FieldMatrix<double,6,dofel> B_cen(0.0);
145 for (int i = 0; i < dof; i++)
146 {
147 B_cen[0][i] = gradphi_cen[i][0]; // E11
148 B_cen[1][i + dof] = gradphi_cen[i][1]; // E22
149 B_cen[2][i + 2 * dof] = gradphi_cen[i][2]; // E33
150 B_cen[3][i + dof] = gradphi_cen[i][2]; B_cen[3][i + 2 * dof] = gradphi_cen[i][1]; // E23
151 B_cen[4][i] = gradphi_cen[i][2]; B_cen[4][i + 2 * dof] = gradphi_cen[i][0]; // E13
152 B_cen[5][i] = gradphi_cen[i][1]; B_cen[5][i + dof] = gradphi_cen[i][0]; // E12
153 }
154
155 Dune::FieldVector<double,6> e_cen(0.0),s_cen(0.0);
156 B_cen.mv(xlc,e_cen); // e = B * x
157 Cijkl.mv(e_cen,s_cen); // s = Cijkl * e
158
159 s_cen -= appliedStress;
160
161 for(int i = 0; i < dof; i++)
162 {
163 double maxUpdate = sqrt(xlc[i]*xlc[i]+xlc[i+dof]*xlc[i+dof]+xlc[i+2*dof]*xlc[i+2*dof]);
164 if(maxUpdate > maxDisp){ maxDisp = maxUpdate; }
165 }
166
167 //Stresses
168 dgf.evaluate(eg.entity(),pt,num_elem);
169 r.accumulate(lfsv1, iter, s[0]);
170 r.accumulate(lfsv2, iter, s[1] /num_elem[0]);
171 r.accumulate(lfsv3, iter, s_cen[2]/num_elem[0]);
172 r.accumulate(lfsv4, iter, s[3]);
173 r.accumulate(lfsv5, iter, s[4] /num_elem[0]);
174 r.accumulate(lfsv6, iter, s[5]);
175 }
176
177
178 }
179 private:
180
181 Dune::FieldMatrix<double,6,6> getR(Dune::FieldMatrix<double,3,3> A) const
182 {
183 Dune::FieldMatrix<double,6,6> R;
184 R[0][0] = A[0][0]*A[0][0]; R[0][1] = A[0][1]*A[0][1]; R[0][2] = A[0][2]*A[0][2];
185 R[1][0] = A[1][0]*A[1][0]; R[1][1] = A[1][1]*A[1][1]; R[1][2] = A[1][2]*A[1][2];
186 R[2][0] = A[2][0]*A[2][0]; R[2][1] = A[2][1]*A[2][1]; R[2][2] = A[2][2]*A[2][2];
187
188 R[0][3] = 2.0*A[0][1]*A[0][2]; R[0][4] = 2.0*A[0][0]*A[0][2]; R[0][5] = 2.0*A[0][0]*A[0][1];
189 R[1][3] = 2.0*A[1][1]*A[1][2]; R[1][4] = 2.0*A[1][0]*A[1][2]; R[1][5] = 2.0*A[1][0]*A[1][1];
190 R[2][3] = 2.0*A[2][1]*A[2][2]; R[2][4] = 2.0*A[2][0]*A[2][2]; R[2][5] = 2.0*A[2][0]*A[2][1];
191
192 R[3][0] = A[1][0]*A[2][0]; R[3][1] = A[1][1]*A[2][1]; R[3][2] = A[1][2]*A[2][2];
193 R[4][0] = A[0][0]*A[2][0]; R[4][1] = A[0][1]*A[2][1]; R[4][2] = A[0][2]*A[2][2];
194 R[5][0] = A[0][0]*A[1][0]; R[5][1] = A[0][1]*A[1][1]; R[5][2] = A[0][2]*A[1][2];
195
196 R[3][3] = A[1][1]*A[2][2]+A[1][2]*A[2][1]; R[3][4] = A[1][0]*A[2][2]+A[1][2]*A[2][0]; R[3][5] = A[1][0]*A[2][1]+A[1][1]*A[2][0];
197 R[4][3] = A[0][1]*A[2][2]+A[2][1]*A[0][2]; R[4][4] = A[0][0]*A[2][2]+A[2][0]*A[0][2]; R[4][5] = A[2][0]*A[0][1]+A[2][1]*A[0][0];
198 R[5][3] = A[0][1]*A[1][2]+A[0][2]*A[1][1]; R[5][4] = A[0][0]*A[1][2]+A[0][2]*A[1][0]; R[5][5] = A[0][0]*A[1][1]+A[0][1]*A[1][0];
199 return R;
200 }
201
202 Dune::FieldMatrix<double,6,6> getMatrix(Dune::FieldMatrix<double,6,6> C, std::vector<double> theta) const
203 {
204 double c = cos(theta[1]*M_PI/180.0);
205 double s = sin(theta[1]*M_PI/180.0);
206 Dune::FieldMatrix<double,3,3> A = {{c,0,-s},{0,1,0},{s,0,c}};
207
208 Dune::FieldMatrix<double,6,6> R = getR(A);
209
210 auto Ch = C;
211 Dune::FieldMatrix<double,6,6> RT(0.0);
212 Composites::transpose<Dune::FieldMatrix<double,6,6>,6>(R,RT);
213 Ch.rightmultiply(RT);
214 return Ch;
215 }
216
217 const GV& gv;
218 const PROBLEM& problem;
219 const DGF& dgf;
220 double& maxDisp;
221 };
222
223
224 template<typename GV, typename MODEL, int dofel>
225 class getStress :
226 public NumericalJacobianApplyVolume<getStress<GV,MODEL,dofel> >,
227 public FullVolumePattern,
228 public LocalOperatorDefaultFlags,
229 public InstationaryLocalOperatorDefaultMethods<double>,
230 public NumericalJacobianVolume<getStress<GV,MODEL, dofel> >
231 {
232 public:
233 // pattern assembly flags
234 enum { doPatternVolume = false };
235
236 // residual assembly flags
237 enum { doAlphaVolume = true };
238
239 //Constructor
240 getStress (const GV& gv_, const MODEL& model_) : gv(gv_), model(model_) {}
241 //
242 // alpha_volume for getStress
243 //
244 // Volume Integral Depending on Test and Ansatz Functions
245 template<typename EG, typename LFSU, typename X, typename LFSV, typename R>
246 void alpha_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv, R& r) const
247 {
248 const int dim = 3;
249
250 //Unwrap function spaces
251 typedef typename LFSU::template Child<0>::Type LFS1;
252 typedef typename LFS1::Traits::FiniteElementType::Traits::LocalBasisType::Traits::DomainFieldType DF;
253 typedef typename LFS1::Traits::FiniteElementType::Traits::LocalBasisType::Traits::JacobianType JacobianType;
254 //Unwrap displacements
255 const auto& lfs1 = lfsu.template child<0>();
256 const auto& lfs2 = lfsu.template child<1>();
257 const auto& lfs3 = lfsu.template child<2>();
258 //Unwrap stress
259 const auto& lfsv1 = lfsv.template child<0>();
260 const auto& lfsv2 = lfsv.template child<1>();
261 const auto& lfsv3 = lfsv.template child<2>();
262 const auto& lfsv4 = lfsv.template child<3>();
263 const auto& lfsv5 = lfsv.template child<4>();
264 const auto& lfsv6 = lfsv.template child<5>();
265
266 // Get element id
267 int id = gv.indexSet().index(eg.entity());
268
269 //Evaluate elasticity tensor
270 Dune::FieldMatrix<double,6,6> C = model.getElasticTensor(id);
271 auto theta = model.returnTheta(eg.geometry().center());
272 auto Cthermal = model.TensorRotationRight(C,-theta[2],2);
273
274 int nodes_per_element = lfsv1.size();
275 int dof = lfs1.size();
276
277 //Residual heat stresses
278 Dune::FieldVector<double,6> delta(0.0);
279 model.evaluateHeat(delta, id);
280
281 //Calculate input heat stress
282 Dune::FieldVector<double,6> appliedStress(0.0);
283 Cthermal.mv(delta, appliedStress);
284
285 // select quadrature rule
286 GeometryType gt = eg.geometry().type();
287 const QuadratureRule<DF,dim>& rule = QuadratureRules<DF,dim>::rule(gt,1);
288 for (typename Dune::QuadratureRule<DF,dim>::const_iterator it = rule.begin(),endit = rule.end(); it != endit; ++it)
289 {
290 Dune::FieldVector<double,dim> pt = it->position();
291 Dune::FieldVector<DF,dim> local = eg.geometry().global(pt);
292
293 auto theta = model.returnTheta(local);
294 auto Cijkl = model.TensorRotationRight(C,theta[1],1);
295
296 // Evaluate Jacobian
297 std::vector<JacobianType> js(dof);
298 lfs1.finiteElement().localBasis().evaluateJacobian(pt,js);
299
300 // Transform gradient to real element and evaluate at element center
301 const auto jac = eg.geometry().jacobianInverseTransposed(pt);
302 std::vector<Dune::FieldVector<double,dim> > gradphi(dof);
303 for (int i=0; i < dof; i++)
304 {
305 gradphi[i] = 0.0;
306 jac.umv(js[i][0],gradphi[i]);
307 }
308 Dune::FieldMatrix<double,6,dofel> B(0.0);
309 for (int i = 0; i < dof; i++)
310 {
311 B[0][i] = gradphi[i][0]; // E11
312 B[1][i + dof] = gradphi[i][1]; // E22
313 B[2][i + 2 * dof] = gradphi[i][2]; // E33
314 B[3][i + dof] = gradphi[i][2]; B[3][i + 2 * dof] = gradphi[i][1]; // E23
315 B[4][i] = gradphi[i][2]; B[4][i + 2 * dof] = gradphi[i][0]; // E13
316 B[5][i] = gradphi[i][1]; B[5][i + dof] = gradphi[i][0]; // E12
317 }
318
319 Dune::FieldVector<double,6> e(0.0),s(0.0);
320 Dune::FieldVector<double,dofel> xlc(0.0);
321 for(int i = 0; i < dof; i++)
322 {
323 xlc[i] = x(lfs1,i);
324 xlc[i+dof] = x(lfs2,i);
325 xlc[i+2*dof] = x(lfs3,i);
326 }
327
328 B.mv(xlc,e); // e = B * x
329 Cijkl.mv(e,s); // s = Cijkl * e
330
331 s-=appliedStress;
332
333 for (int i = 0; i < nodes_per_element; i++)
334 {
335 r.accumulate(lfsv1, i, s[0] * it->weight());
336 r.accumulate(lfsv2, i, s[1] * it->weight());
337 r.accumulate(lfsv3, i, s[2] * it->weight());
338 r.accumulate(lfsv4, i, s[3] * it->weight());
339 r.accumulate(lfsv5, i, s[4] * it->weight());
340 r.accumulate(lfsv6, i, s[5] * it->weight());
341 }
342 }
343 }
344 const GV& gv;
345 const MODEL& model;
346
347 };
348
349 class countElem :
350 public NumericalJacobianApplyVolume<countElem >,
351 public FullVolumePattern,
352 public LocalOperatorDefaultFlags,
353 public InstationaryLocalOperatorDefaultMethods<double>,
354 public NumericalJacobianVolume<countElem >
355 {
356 public:
357 // pattern assembly flags
358 enum { doPatternVolume = false };
359
360 // residual assembly flags
361 enum { doLambdaVolume = true };
362
363 //Constructor
364 countElem() = default;
365
366 // alpha_volume for countElem
367 template<typename EG, typename LFSV, typename R>
368 void lambda_volume (const EG& eg, const LFSV& lfsv, R& r) const
369 {
370 for (typename LFSV::Traits::SizeType i = 0; i < lfsv.size(); i++)
371 {
372 r.accumulate(lfsv, i, 1.);
373 }
374 }
375 };
376
377
378 template<typename GV, typename PROBLEM, typename DGF, int dofel>
379 class getStressUG :
380 public NumericalJacobianApplyVolume<getStressUG<GV,PROBLEM,DGF,dofel> >,
381 public FullVolumePattern,
382 public LocalOperatorDefaultFlags,
383 public InstationaryLocalOperatorDefaultMethods<double>,
384 public NumericalJacobianVolume<getStressUG<GV,PROBLEM,DGF, dofel> >
385 {
386 public:
387 // pattern assembly flags
388 enum { doPatternVolume = false };
389
390 // residual assembly flags
391 enum { doAlphaVolume = true };
392
393 //Constructor
394 getStressUG (const GV& gv_, const PROBLEM& prob_, const DGF& dgf_, double& maxDisp_)
395 : gv(gv_)
396 , model(prob_)
397 , dgf(dgf_)
398 , maxDisp(maxDisp_)
399 {}
400 //
401 // alpha_volume for getStress
402 //
403 // Volume Integral Depending on Test and Ansatz Functions
404 template<typename EG, typename LFSU, typename X, typename LFSV, typename R>
405 void alpha_volume (const EG& eg, const LFSU& lfsu, const X& x, const LFSV& lfsv, R& r) const
406 {
407 const int dim = 3;
408
409 //Unwrap function spaces
410 typedef typename LFSU::template Child<0>::Type LFS1;
411 typedef typename LFS1::Traits::FiniteElementType::Traits::LocalBasisType::Traits::DomainFieldType DF;
412 typedef typename LFS1::Traits::FiniteElementType::Traits::LocalBasisType::Traits::JacobianType JacobianType;
413 //Unwrap displacements
414 const auto& lfs1 = lfsu.template child<0>();
415 const auto& lfs2 = lfsu.template child<1>();
416 const auto& lfs3 = lfsu.template child<2>();
417 //Unwrap stress
418 const auto& lfsv1 = lfsv.template child<0>();
419 const auto& lfsv2 = lfsv.template child<1>();
420 const auto& lfsv3 = lfsv.template child<2>();
421 const auto& lfsv4 = lfsv.template child<3>();
422 const auto& lfsv5 = lfsv.template child<4>();
423 const auto& lfsv6 = lfsv.template child<5>();
424
425 // Get element id
426 int id = gv.indexSet().index(eg.entity());
427
428 //Evaluate elasticity tensor
429 Dune::FieldMatrix<double,6,6> C = model.getElasticTensor(id);
430 auto theta = model.returnTheta(eg.geometry().center());
431 auto Cthermal = model.TensorRotationRight(C,-theta[2],2);
432
433 int nodes_per_element = lfsv1.size();
434 int dof = lfs1.size();
435
436 //Residual heat stresses
437 Dune::FieldVector<double,6> delta(0.0);
438 model.evaluateHeat(delta, id);
439
440 //Calculate input heat stress
441 Dune::FieldVector<double,6> appliedStress(0.0);
442 Cthermal.mv(delta, appliedStress);
443
444 // select quadrature rule
445 GeometryType gt = eg.geometry().type();
446 const QuadratureRule<DF,dim>& rule = QuadratureRules<DF,dim>::rule(gt,1);
447 for (typename Dune::QuadratureRule<DF,dim>::const_iterator it = rule.begin(),endit = rule.end(); it != endit; ++it)
448 {
449 Dune::FieldVector<double,dim> pt = it->position();
450 Dune::FieldVector<DF,dim> local = eg.geometry().global(pt);
451
452 auto theta = model.returnTheta(local);
453 auto Cijkl = model.TensorRotationRight(C,theta[1],1);
454
455 // Evaluate Jacobian
456 std::vector<JacobianType> js(dof);
457 lfs1.finiteElement().localBasis().evaluateJacobian(pt,js);
458
459 // Transform gradient to real element and evaluate at element center
460 const auto jac = eg.geometry().jacobianInverseTransposed(pt);
461 std::vector<Dune::FieldVector<double,dim> > gradphi(dof);
462 for (int i=0; i < dof; i++)
463 {
464 gradphi[i] = 0.0;
465 jac.umv(js[i][0],gradphi[i]);
466 }
467 Dune::FieldMatrix<double,6,dofel> B(0.0);
468 for (int i = 0; i < dof; i++)
469 {
470 B[0][i] = gradphi[i][0]; // E11
471 B[1][i + dof] = gradphi[i][1]; // E22
472 B[2][i + 2 * dof] = gradphi[i][2]; // E33
473 B[3][i + dof] = gradphi[i][2]; B[3][i + 2 * dof] = gradphi[i][1]; // E23
474 B[4][i] = gradphi[i][2]; B[4][i + 2 * dof] = gradphi[i][0]; // E13
475 B[5][i] = gradphi[i][1]; B[5][i + dof] = gradphi[i][0]; // E12
476 }
477
478 Dune::FieldVector<double,6> e(0.0),s(0.0);
479 Dune::FieldVector<double,dofel> xlc(0.0);
480 for(int i = 0; i < dof; i++)
481 {
482 xlc[i] = x(lfs1,i);
483 xlc[i+dof] = x(lfs2,i);
484 xlc[i+2*dof] = x(lfs3,i);
485 }
486
487 B.mv(xlc,e); // e = B * x
488 Cijkl.mv(e,s); // s = Cijkl * e
489
490 s-=appliedStress;
491
492 //Stresses
493 for(int i = 0; i < nodes_per_element; i++)
494 {
495 Dune::FieldVector<double,3> pt = eg.geometry().local(eg.geometry().corner(i));
496 Dune::FieldVector<double,1> num_elem(1.0);
497 dgf.evaluate(eg.entity(),pt,num_elem);
498 r.accumulate(lfsv1, i, s[0]/num_elem[0]*it->weight());
499 r.accumulate(lfsv2, i, s[1]/num_elem[0]*it->weight());
500 r.accumulate(lfsv3, i, s[2]/num_elem[0]*it->weight());
501 r.accumulate(lfsv4, i, s[3]/num_elem[0]*it->weight());
502 r.accumulate(lfsv5, i, s[4]/num_elem[0]*it->weight());
503 r.accumulate(lfsv6, i, s[5]/num_elem[0]*it->weight());
504 }
505
506 }
507
508 }
509 private:
510 const GV& gv;
511 const PROBLEM& model;
512 const DGF& dgf;
513 double& maxDisp;
514 };
515
516 }
517}
518
519#endif
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 5, 23:02, 2025)