FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
Hyperelastic.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7  use hecmw_util
8  use mmaterial
9  implicit none
10 
11 contains
12 
14  subroutine cderiv( matl, sectType, ctn, itn, &
15  inv1b, inv2b, inv3b, dibdc, d2ibdc2, strain, direction, inv4b, dibdc_ani, d2ibdc2_ani )
16  type( tmaterial ), intent(in) :: matl
17  integer, intent(in) :: sectType
18  real(kind=kreal), intent(out) :: inv1b
19  real(kind=kreal), intent(out) :: inv2b
20  real(kind=kreal), intent(out) :: inv3b
21  real(kind=kreal), intent(out) :: dibdc(3,3,3)
22  real(kind=kreal), intent(out) :: d2ibdc2(3,3,3,3,3)
23  real(kind=kreal), intent(in) :: strain(6)
24  real(kind=kreal), intent(out) :: ctn(3,3)
25  real(kind=kreal), intent(out) :: itn(3,3)
26  real(kind=kreal), intent(in), optional :: direction(3)
27  real(kind=kreal), intent(out), optional :: inv4b
28  real(kind=kreal), intent(out), optional :: dibdc_ani(3,3)
29  real(kind=kreal), intent(out), optional :: d2ibdc2_ani(3,3,3,3)
30 
31  integer :: i, j, k, l, m, n
32 
33  real(kind=kreal) :: inv1, inv2, inv3, inv33
34  real(kind=kreal) :: delta(3,3)
35  real(kind=kreal) :: didc(3,3,3), ctninv(3,3)
36  real(kind=kreal) :: d2idc2(3,3,3,3,3)
37 
38  ! Kronecker's delta
39  delta(:,:) = 0.d0
40  delta(1,1) = 1.d0
41  delta(2,2) = 1.d0
42  delta(3,3) = 1.d0
43 
44  ! Green-Lagrange strain to right Cauchy-Green deformation tensor
45  ctn(1,1)=strain(1)*2.d0+1.d0
46  ctn(2,2)=strain(2)*2.d0+1.d0
47  ctn(3,3)=strain(3)*2.d0+1.d0
48  ctn(1,2)=strain(4); ctn(2,1)=ctn(1,2)
49  ctn(2,3)=strain(5); ctn(3,2)=ctn(2,3)
50  ctn(3,1)=strain(6); ctn(1,3)=ctn(3,1)
51 
52  itn(:,:) = delta(:,:)
53 
54  ! ----- calculate the invariant of C
55  inv1 = ctn(1,1)+ctn(2,2)+ctn(3,3)
56  inv2 = ctn(2,2)*ctn(3,3)+ctn(1,1)*ctn(3,3)+ctn(1,1)*ctn(2,2) &
57  -ctn(2,3)*ctn(2,3)-ctn(1,3)*ctn(1,3)-ctn(1,2)*ctn(1,2)
58  inv3 = ctn(1,1)*ctn(2,2)*ctn(3,3) &
59  +ctn(2,1)*ctn(3,2)*ctn(1,3) &
60  +ctn(3,1)*ctn(1,2)*ctn(2,3) &
61  -ctn(3,1)*ctn(2,2)*ctn(1,3) &
62  -ctn(2,1)*ctn(1,2)*ctn(3,3) &
63  -ctn(1,1)*ctn(3,2)*ctn(2,3)
64  inv33 = inv3**(-1.d0/3.d0)
65 
66  ! ----- calculate the inverse of C
67  ctninv(1,1)=(ctn(2,2)*ctn(3,3)-ctn(2,3)*ctn(2,3))/inv3
68  ctninv(2,2)=(ctn(1,1)*ctn(3,3)-ctn(1,3)*ctn(1,3))/inv3
69  ctninv(3,3)=(ctn(1,1)*ctn(2,2)-ctn(1,2)*ctn(1,2))/inv3
70  ctninv(1,2)=(ctn(1,3)*ctn(2,3)-ctn(1,2)*ctn(3,3))/inv3
71  ctninv(1,3)=(ctn(1,2)*ctn(2,3)-ctn(2,2)*ctn(1,3))/inv3
72  ctninv(2,3)=(ctn(1,2)*ctn(1,3)-ctn(1,1)*ctn(2,3))/inv3
73  ctninv(2,1)=ctninv(1,2)
74  ctninv(3,1)=ctninv(1,3)
75  ctninv(3,2)=ctninv(2,3)
76 
77  ! ----- calculate the derivative of the C-invariant with respect to c(i,j)
78  do i=1,3
79  do j=1,3
80  didc(i,j,1) = delta(i,j)
81  didc(i,j,2) = inv1*delta(i,j)-ctn(i,j)
82  didc(i,j,3) = inv3*ctninv(i,j)
83  enddo
84  enddo
85 
86  ! ----- calculate the second derivative of the C-invariant
87  do k=1,3
88  do l=1,3
89  do m=1,3
90  do n=1,3
91  d2idc2(k,l,m,n,1)=0.d0
92  d2idc2(k,l,m,n,2)=delta(k,l)*delta(m,n)- &
93  (delta(k,m)*delta(l,n)+delta(k,n)*delta(l,m))/2.d0
94  d2idc2(k,l,m,n,3)=inv3*(ctninv(m,n)*ctninv(k,l)- &
95  (ctninv(k,m)*ctninv(n,l)+ctninv(k,n)*ctninv(m,l))/2.d0)
96  enddo
97  enddo
98  enddo
99  enddo
100 
101  ! ----- derivatives for the reduced invariants
102  inv1b = inv1*inv33
103  inv2b = inv2*inv33*inv33
104  inv3b = dsqrt(inv3)
105 
106  ! --- first derivative the reduced c-invarians w.r.t. c(i,j)
107  do i=1,3
108  do j=1,3
109  dibdc(i,j,1) = -inv33**4*inv1*didc(i,j,3)/3.d0 + inv33*didc(i,j,1)
110  dibdc(i,j,2) = -2.d0*inv33**5*inv2*didc(i,j,3)/3.d0 + inv33**2*didc(i,j,2)
111  dibdc(i,j,3) = didc(i,j,3)/(2.d0*dsqrt(inv3))
112  enddo
113  enddo
114 
115  ! --- second derivative of the reduced c-invariants w.r.t. c(i,j)
116  do i=1,3
117  do j=1,3
118  do k=1,3
119  do l=1,3
120  d2ibdc2(i,j,k,l,1) = 4.d0/9.d0*inv33**7*inv1*didc(i,j,3)*didc(k,l,3) &
121  -inv33**4/3.d0*(didc(k,l,1)*didc(i,j,3)+didc(i,j,1)*didc(k,l,3)) &
122  -inv33**4/3.d0*inv1*d2idc2(i,j,k,l,3) &
123  +inv33*d2idc2(i,j,k,l,1)
124  d2ibdc2(i,j,k,l,2) = 10.d0/9.d0*inv33**8*inv2*didc(i,j,3)*didc(k,l,3) &
125  -2.d0/3.d0*inv33**5*(didc(k,l,2)*didc(i,j,3)+didc(i,j,2)*didc(k,l,3)) &
126  -2.d0/3.d0*inv33**5*inv2*d2idc2(i,j,k,l,3) &
127  +inv33**2*d2idc2(i,j,k,l,2)
128  d2ibdc2(i,j,k,l,3) = -didc(i,j,3)*didc(k,l,3)/(4.d0*inv3**1.5d0) &
129  +d2idc2(i,j,k,l,3)/(2.d0*dsqrt(inv3))
130  enddo
131  enddo
132  enddo
133  enddo
134 
135  if( present(direction) ) call cderiv_aniso( ctn, inv3, didc(:,:,3), d2idc2(:,:,:,:,3), &
136  direction, inv4b, dibdc_ani, d2ibdc2_ani )
137 
138  end subroutine cderiv
139 
142  subroutine cderiv_aniso( ctn, inv3, didc_3, d2idc2_3, direction, inv4b, dibdc_ani, d2ibdc2_ani )
143  real(kind=kreal), intent(in) :: ctn(3,3)
144  real(kind=kreal), intent(in) :: inv3
145  real(kind=kreal), intent(in) :: didc_3(3,3)
146  real(kind=kreal), intent(in) :: d2idc2_3(3,3,3,3)
147  real(kind=kreal), intent(in) :: direction(3)
148  real(kind=kreal), intent(out) :: inv4b
149  real(kind=kreal), intent(out) :: dibdc_ani(3,3)
150  real(kind=kreal), intent(out) :: d2ibdc2_ani(3,3,3,3)
151 
152  integer :: i, j, k, l, m, n
153  real(kind=kreal) :: inv4, inv33, inv3m43, inv4d3
154  real(kind=kreal) :: d2ibdc24
155 
156  inv33 = inv3**(-1.d0/3.d0) ! = I_3^{-1/3}
157  inv3m43 = inv33/inv3 ! = I_3^{-4/3}
158 
159  inv4 = 0.d0
160  do i=1,3
161  do j=1,3
162  inv4 = inv4 + direction(j)*ctn(j,i)*direction(i)
163  enddo
164  enddo
165  inv4b = inv33*inv4 ! I4
166  inv4d3 = inv4/inv3 ! I4*I_3^{-1}
167 
168  ! --- first derivative the reduced 4th c-invarians w.r.t. c(i,j)
169  do i=1,3
170  do j=1,3
171  dibdc_ani(i,j) = inv33*( (-1.d0/3.d0)*didc_3(i,j)*inv4d3+direction(i)*direction(j) )
172  enddo
173  enddo
174 
175  ! --- second derivative of the reduced c-invariants w.r.t. c(i,j)
176  do k=1,3
177  do l=1,3
178  do m=1,3
179  do n=1,3
180  d2ibdc24 = (2.d0/3.d0)*didc_3(k,l)*didc_3(m,n)*inv4d3
181  d2ibdc24 = d2ibdc24 - d2idc2_3(k,l,m,n)*inv4
182  d2ibdc24 = d2ibdc24 - direction(k)*direction(l)*didc_3(m,n)
183  d2ibdc24 = d2ibdc24 - didc_3(k,l)*direction(m)*direction(n)
184  d2ibdc2_ani(k,l,m,n) = (1.d0/3.d0)*inv3m43*d2ibdc24
185  enddo
186  enddo
187  enddo
188  enddo
189 
190  end subroutine
191 
192 
193  !-------------------------------------------------------------------------------
195  !
196  !-------------------------------------------------------------------------------
197  subroutine calelasticarrudaboyce( matl, sectType, cijkl, strain )
198  type( tmaterial ), intent(in) :: matl
199  integer, intent(in) :: secttype
200  real(kind=kreal), intent(out) :: cijkl(3,3,3,3)
201  real(kind=kreal), intent(in) :: strain(6)
202 
203  integer :: i, j, k, l
204  real(kind=kreal) :: ctn(3,3), itn(3,3)
205  real(kind=kreal) :: inv1b, inv2b, inv3b
206  real(kind=kreal) :: dibdc(3,3,3)
207  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
208  real(kind=kreal) :: constant(3), coef
209 
210  call cderiv( matl, secttype, ctn, itn, inv1b, inv2b, inv3b, &
211  dibdc, d2ibdc2, strain )
212  constant(1:3)=matl%variables(m_plconst1:m_plconst3)
213  coef = constant(2)
214 
215  do l=1,3
216  do k=1,3
217  do j=1,3
218  do i=1,3
219  cijkl(i,j,k,l) = constant(1)*(1.d0/(10.d0*coef**2) &
220  +66.d0*inv1b/(1050.d0*coef**4) &
221  +228.d0*inv1b**2/(7000.d0*coef**6) &
222  +10380.d0*inv1b**3/(673750.d0*coef**8)) &
223  *dibdc(i,j,1)*dibdc(k,l,1) &
224  +constant(1)*(0.5d0+inv1b/(10.d0*coef**2) &
225  +33.d0*inv1b**2/(1050.d0*coef**4) &
226  +76.d0*inv1b**3/(7000.d0*coef**6) &
227  +2595.d0*inv1b**4/(673750.d0*coef**8)) &
228  *d2ibdc2(i,j,k,l,1) &
229  +(1.d0+1.d0/inv3b**2)*dibdc(i,j,3)*dibdc(k,l,3)/constant(3) &
230  +(inv3b-1.d0/inv3b)*d2ibdc2(i,j,k,l,3)/constant(3)
231  enddo
232  enddo
233  enddo
234  enddo
235  cijkl(:,:,:,:) = 4.d0*cijkl(:,:,:,:)
236 
237  end subroutine calelasticarrudaboyce
238 
239  !-------------------------------------------------------------------------------
241  subroutine calupdateelasticarrudaboyce( matl, sectType, dstrain, dstress )
242  type( tmaterial ), intent(in) :: matl
243  integer, intent(in) :: secttype
244  real(kind=kreal), intent(out) :: dstress(6)
245  real(kind=kreal), intent(in) :: dstrain(6)
246 
247  integer :: i, j
248  real(kind=kreal) :: ctn(3,3), itn(3,3)
249  real(kind=kreal) :: inv1b, inv2b, inv3b
250  real(kind=kreal) :: dibdc(3,3,3)
251  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
252  real(kind=kreal) :: constant(3), coef
253  real(kind=kreal) :: pkstress(3,3)
254 
255 
256  constant(1:3)=matl%variables(m_plconst1:m_plconst3)
257  coef = constant(2)
258  call cderiv( matl, secttype, ctn, itn, &
259  inv1b, inv2b, inv3b, &
260  dibdc, d2ibdc2, dstrain )
261 
262 
263  ! ----- calculate the stress
264  do i=1,3
265  do j=1,3
266  pkstress(i,j) = constant(1)*( 0.5d0+inv1b/(10.d0*coef**2) &
267  +33.d0*inv1b*inv1b/(1050.d0*coef**4) &
268  +76.d0*inv1b**3/(7000.d0*coef**6) &
269  +2595.d0*inv1b**4/(673750.d0*coef**8))*dibdc(i,j,1) &
270  +(inv3b-1.d0/inv3b)*dibdc(i,j,3)/constant(3)
271  enddo
272  enddo
273 
274  dstress(1) = 2.d0*pkstress(1,1)
275  dstress(2) = 2.d0*pkstress(2,2)
276  dstress(3) = 2.d0*pkstress(3,3)
277  dstress(4) = pkstress(1,2) + pkstress(2,1)
278  dstress(5) = pkstress(2,3) + pkstress(3,2)
279  dstress(6) = pkstress(1,3) + pkstress(3,1)
280 
281  end subroutine calupdateelasticarrudaboyce
282 
283  !-------------------------------------------------------------------------------
285  !-------------------------------------------------------------------------------
286  subroutine calelasticmooneyrivlin( matl, sectType, cijkl, strain, hdflag )
287  type( tmaterial ), intent(in) :: matl
288  integer, intent(in) :: secttype
289  real(kind=kreal), intent(out) :: cijkl(3,3,3,3)
290  real(kind=kreal), intent(in) :: strain(6)
291  integer(kind=kint), intent(in), optional :: hdflag
292 
293  integer :: k, l, m, n
294  real(kind=kreal) :: ctn(3,3), itn(3,3)
295  real(kind=kreal) :: inv1b, inv2b, inv3b
296  real(kind=kreal) :: dibdc(3,3,3)
297  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
298  real(kind=kreal) :: constant(3)
299 
300 
301  constant(1:3)=matl%variables(m_plconst1:m_plconst3)
302  if( present(hdflag) ) then
303  if( hdflag == 1 ) then
304  constant(3)=0.d0
305  else if( hdflag == 2 ) then
306  constant(1:2)=0.d0
307  end if
308  end if
309  if( dabs(constant(3))>1.d-12 ) constant(3) = 1.d0/constant(3)
310 
311  call cderiv( matl, secttype, ctn, itn, inv1b, inv2b, inv3b, &
312  dibdc, d2ibdc2, strain )
313 
314  do n=1,3
315  do m=1,3
316  do l=1,3
317  do k=1,3
318  cijkl(k,l,m,n) = d2ibdc2(k,l,m,n,1)*constant(1) + &
319  d2ibdc2(k,l,m,n,2)*constant(2) + &
320  2.d0*(dibdc(k,l,3)*dibdc(m,n,3)+ &
321  (inv3b-1.d0)*d2ibdc2(k,l,m,n,3))*constant(3)
322  enddo
323  enddo
324  enddo
325  enddo
326  cijkl(:,:,:,:)=4.d0*cijkl(:,:,:,:)
327 
328  end subroutine calelasticmooneyrivlin
329 
330  !-------------------------------------------------------------------------------
332  !-------------------------------------------------------------------------------
333  subroutine calupdateelasticmooneyrivlin( matl, sectType, strain, stress, strain_energy, hdflag )
334  type( tmaterial ), intent(in) :: matl
335  integer, intent(in) :: secttype
336  real(kind=kreal), intent(out) :: stress(6)
337  real(kind=kreal), intent(in) :: strain(6)
338  real(kind=kreal), intent(out) :: strain_energy
339  integer(kind=kint), intent(in), optional :: hdflag
340 
341  integer :: k, l
342  real(kind=kreal) :: ctn(3,3), itn(3,3)
343  real(kind=kreal) :: inv1b, inv2b, inv3b
344  real(kind=kreal) :: dibdc(3,3,3)
345  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
346  real(kind=kreal) :: constant(3)
347  real(kind=kreal) :: dudc(3,3)
348 
349  constant(1:3)=matl%variables(m_plconst1:m_plconst3)
350  if( present(hdflag) ) then
351  if( hdflag == 1 ) then
352  constant(3)=0.d0
353  else if( hdflag == 2 ) then
354  constant(1:2)=0.d0
355  end if
356  end if
357  if( dabs(constant(3))>1.d-12 ) constant(3) = 1.d0/constant(3)
358 
359  call cderiv( matl, secttype, ctn, itn, inv1b, inv2b, inv3b, &
360  dibdc, d2ibdc2, strain )
361 
362 
363  ! ----- stress
364  do l=1,3
365  do k=1,3
366  dudc(k,l) = dibdc(k,l,1)*constant(1)+dibdc(k,l,2)*constant(2) &
367  +2.d0*(inv3b-1.d0)*dibdc(k,l,3)*constant(3)
368  enddo
369  enddo
370 
371  stress(1)=2.d0*dudc(1,1)
372  stress(2)=2.d0*dudc(2,2)
373  stress(3)=2.d0*dudc(3,3)
374  stress(4)=2.d0*dudc(1,2)
375  stress(5)=2.d0*dudc(2,3)
376  stress(6)=2.d0*dudc(1,3)
377 
378  strain_energy = (inv1b-3.d0)*constant(1)+(inv2b-3.d0)*constant(2) &
379  +(inv3b-1.d0)*(inv3b-1.d0)*constant(3)
380 
381  end subroutine calupdateelasticmooneyrivlin
382 
383  !-------------------------------------------------------------------------------
385  !-------------------------------------------------------------------------------
386  subroutine calelasticmooneyrivlinaniso( matl, sectType, cijkl, strain, cdsys, hdflag )
387  type( tmaterial ), intent(in) :: matl
388  integer, intent(in) :: secttype
389  real(kind=kreal), intent(out) :: cijkl(3,3,3,3)
390  real(kind=kreal), intent(in) :: strain(6)
391  real(kind=kreal), intent(in) :: cdsys(3,3)
392  integer(kind=kint), intent(in), optional :: hdflag
393 
394  integer :: i, j, k, l, m, n, jj
395  real(kind=kreal) :: ctn(3,3), itn(3,3)
396  real(kind=kreal) :: inv1b, inv2b, inv3b, inv4b
397  real(kind=kreal) :: dibdc(3,3,3)
398  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
399  real(kind=kreal) :: constant(5)
400  real(kind=kreal) :: dibdc_ani(3,3)
401  real(kind=kreal) :: d2ibdc2_ani(3,3,3,3)
402 
403 
404  constant(1:5)=matl%variables(m_plconst1:m_plconst5)
405  if( present(hdflag) ) then
406  if( hdflag == 1 ) then
407  constant(3)=0.d0
408  else if( hdflag == 2 ) then
409  constant(1:2)=0.d0
410  constant(4:5)=0.d0
411  end if
412  end if
413  if( dabs(constant(3))>1.d-12 ) constant(3) = 1.d0/constant(3)
414 
415  call cderiv( matl, secttype, ctn, itn, inv1b, inv2b, inv3b, &
416  dibdc, d2ibdc2, strain, cdsys(1,1:3), inv4b, dibdc_ani, d2ibdc2_ani )
417 
418  do n=1,3
419  do m=1,3
420  do l=1,3
421  do k=1,3
422  cijkl(k,l,m,n) = d2ibdc2(k,l,m,n,1)*constant(1) + &
423  d2ibdc2(k,l,m,n,2)*constant(2) + &
424  2.d0*(dibdc(k,l,3)*dibdc(m,n,3)+ &
425  (inv3b-1.d0)*d2ibdc2(k,l,m,n,3))*constant(3)+ &
426  (2.d0*constant(4)+6.d0*(inv4b-1.d0)*constant(5))*dibdc_ani(k,l)*dibdc_ani(m,n)+ &
427  (inv4b-1.d0)*(2.d0*constant(4)+3.d0*(inv4b-1.d0)*constant(5))*d2ibdc2_ani(k,l,m,n)
428  enddo
429  enddo
430  enddo
431  enddo
432  cijkl(:,:,:,:)=4.d0*cijkl(:,:,:,:)
433 
434  end subroutine calelasticmooneyrivlinaniso
435 
436  !-------------------------------------------------------------------------------
438  !-------------------------------------------------------------------------------
439  subroutine calupdateelasticmooneyrivlinaniso( matl, sectType, strain, stress, cdsys, hdflag )
440  type( tmaterial ), intent(in) :: matl
441  integer, intent(in) :: secttype
442  real(kind=kreal), intent(out) :: stress(6)
443  real(kind=kreal), intent(in) :: strain(6)
444  real(kind=kreal), intent(in) :: cdsys(3,3)
445  integer(kind=kint), intent(in), optional :: hdflag
446 
447  integer :: k, l
448  real(kind=kreal) :: ctn(3,3), itn(3,3)
449  real(kind=kreal) :: inv1b, inv2b, inv3b, inv4b
450  real(kind=kreal) :: dibdc(3,3,3)
451  real(kind=kreal) :: d2ibdc2(3,3,3,3,3)
452  real(kind=kreal) :: constant(5)
453  real(kind=kreal) :: dudc(3,3)
454  real(kind=kreal) :: dibdc_ani(3,3)
455  real(kind=kreal) :: d2ibdc2_ani(3,3,3,3)
456 
457  constant(1:5)=matl%variables(m_plconst1:m_plconst5)
458  if( present(hdflag) ) then
459  if( hdflag == 1 ) then
460  constant(3)=0.d0
461  else if( hdflag == 2 ) then
462  constant(1:2)=0.d0
463  constant(4:5)=0.d0
464  end if
465  end if
466  if( dabs(constant(3))>1.d-12 ) constant(3) = 1.d0/constant(3)
467 
468  call cderiv( matl, secttype, ctn, itn, inv1b, inv2b, inv3b, &
469  dibdc, d2ibdc2, strain, cdsys(1,1:3), inv4b, dibdc_ani, d2ibdc2_ani )
470 
471  ! ----- stress
472  do l=1,3
473  do k=1,3
474  dudc(k,l) = dibdc(k,l,1)*constant(1)+dibdc(k,l,2)*constant(2) &
475  +2.d0*(inv3b-1.d0)*dibdc(k,l,3)*constant(3) &
476  +(inv4b-1.d0)*(2.d0*constant(4)+3.d0*(inv4b-1.d0)*constant(5))*dibdc_ani(k,l)
477  enddo
478  enddo
479 
480  stress(1)=2.d0*dudc(1,1)
481  stress(2)=2.d0*dudc(2,2)
482  stress(3)=2.d0*dudc(3,3)
483  stress(4)=2.d0*dudc(1,2)
484  stress(5)=2.d0*dudc(2,3)
485  stress(6)=2.d0*dudc(1,3)
486 
487  end subroutine calupdateelasticmooneyrivlinaniso
488 
489 end module
mhyperelastic::calelasticmooneyrivlinaniso
subroutine calelasticmooneyrivlinaniso(matl, sectType, cijkl, strain, cdsys, hdflag)
This subroutine provides elastic tangent coefficient for anisotropic Mooney-Rivlin hyperelastic mater...
Definition: Hyperelastic.f90:387
mhyperelastic
This module provides functions for hyperelastic calculation.
Definition: Hyperelastic.f90:6
mmaterial::m_plconst5
integer(kind=kint), parameter m_plconst5
Definition: material.f90:102
mmaterial::m_plconst3
integer(kind=kint), parameter m_plconst3
Definition: material.f90:100
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
mhyperelastic::calupdateelasticmooneyrivlinaniso
subroutine calupdateelasticmooneyrivlinaniso(matl, sectType, strain, stress, cdsys, hdflag)
This subroutine provides to update stress and strain for anisotropic Mooney-Rivlin material.
Definition: Hyperelastic.f90:440
mhyperelastic::calupdateelasticarrudaboyce
subroutine calupdateelasticarrudaboyce(matl, sectType, dstrain, dstress)
This subroutine provides to update stress and strain for Arrude-Royce material.
Definition: Hyperelastic.f90:242
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
mhyperelastic::cderiv
subroutine cderiv(matl, sectType, ctn, itn, inv1b, inv2b, inv3b, dibdc, d2ibdc2, strain, direction, inv4b, dibdc_ani, d2ibdc2_ani)
This subroutine calculates derivative of the invariant with respect to Cauchy-Green tensor.
Definition: Hyperelastic.f90:16
mmaterial::m_plconst1
integer(kind=kint), parameter m_plconst1
Definition: material.f90:98
mhyperelastic::calupdateelasticmooneyrivlin
subroutine calupdateelasticmooneyrivlin(matl, sectType, strain, stress, strain_energy, hdflag)
This subroutine provides to update stress and strain for Mooney-Rivlin material.
Definition: Hyperelastic.f90:334
mmaterial::tmaterial
Structure to manage all material related data.
Definition: material.f90:167
mhyperelastic::cderiv_aniso
subroutine cderiv_aniso(ctn, inv3, didc_3, d2idc2_3, direction, inv4b, dibdc_ani, d2ibdc2_ani)
This subroutine calculates derivative of the 4th reduced invariant I4 with respect to Cauchy-Green te...
Definition: Hyperelastic.f90:143
mhyperelastic::calelasticmooneyrivlin
subroutine calelasticmooneyrivlin(matl, sectType, cijkl, strain, hdflag)
This subroutine provides elastic tangent coefficient for Mooney-Rivlin hyperelastic material.
Definition: Hyperelastic.f90:287
mmaterial
This module summarizes all information of material properties.
Definition: material.f90:6
mhyperelastic::calelasticarrudaboyce
subroutine calelasticarrudaboyce(matl, sectType, cijkl, strain)
This subroutine provides elastic tangent coefficient for Arruda-Boyce hyperelastic material.
Definition: Hyperelastic.f90:198