├── README.md └── umat_foam.for /README.md: -------------------------------------------------------------------------------- 1 | # ElastomericFoam 2 | Abaqus UMAT and sample input files for modeling of open-cell elastomeric foams 3 | 4 | # Cite 5 | If used, please cite Alexander K. Landauer, Xiuqi Li, Christian Franck, and David L. Henann. Experimental characterization and hyperelastic constitutive modeling of open-cell elastomeric foams. Journal of the Mechanics and Physics of Solids 133 (2019) 103701. 6 | -------------------------------------------------------------------------------- /umat_foam.for: -------------------------------------------------------------------------------- 1 | !************************************************************************ 2 | ! 3 | ! User material subroutine (UMAT) for the large-deformation behavior of 4 | ! compressible, elastomeric foams, using a logarithmic-strain-based 5 | ! hyperelastic model. This UMAT is not for use in plane stress or in any 6 | ! other situation in which there are more strain terms than stress terms. 7 | ! 8 | ! Xiuqi Li, March 2019 9 | ! 10 | !************************************************************************ 11 | ! Usage: 12 | !************************************************************************ 13 | ! 14 | ! Material Properties Vector (*user material, constants = 14) 15 | ! -------------------------------------------------------------- 16 | ! G0 = props(1) ! Shear modulus 17 | ! B = props(2) ! Bulk modulus 18 | ! Jmin = props(3) ! Locking value of J in f-fuction 19 | ! C1 = props(4) ! C_1 in L-function 20 | ! K10 = props(5) ! K1_0 in X-function 21 | ! delta_K = props(6) ! delta_K in X-function 22 | ! X1 = props(7) ! X1' in X-function 23 | ! X2 = props(8) ! X2' in X-function 24 | ! C0 = props(9) ! C_0 in L-function 25 | ! p = props(10) ! power p in L-function 26 | ! q = props(11) ! power q in L-function 27 | ! C2 = props(12) ! C_2 in f-function 28 | ! C3 = props(13) ! C_3 in f-function 29 | ! r = props(14) ! power r in f-function 30 | ! 31 | ! State Variables (*depvar 3) (for visualization purposes only) 32 | ! -------------------------------------------------------------- 33 | ! statev(1) = K1 ---- ln(J), range(-inf, inf) 34 | ! statev(2) = K2 ---- amount of distortional deformation, range[0, inf) 35 | ! statev(3) = K3 ---- deformation mode, range[-1, 1], 36 | ! -1 is compression, 0 is shear, 1 is tension 37 | ! 38 | !************************************************************************ 39 | subroutine umat(stress,statev,ddsdde,sse,spd,scd, 40 | + rpl,ddsddt,drplde,drpldt, 41 | + stran,dstran,time,dtime,temp,dtemp,predef,dpred,cmname, 42 | + ndi,nshr,ntens,nstatv,props,nprops,coords,drot,pnewdt, 43 | + celent,dfgrd0,dfgrd1,noel,npt,layer,kspt,kstep,kinc) 44 | ! 45 | include 'aba_param.inc' 46 | ! 47 | dimension stress(ntens),statev(nstatv), 48 | + ddsdde(ntens,ntens),ddsddt(ntens),drplde(ntens), 49 | + stran(ntens),dstran(ntens),time(2),predef(1),dpred(1), 50 | + props(nprops),coords(3),drot(3,3),dfgrd0(3,3),dfgrd1(3,3) 51 | ! 52 | character*8 cmname 53 | ! 54 | ! Variables defined and used in the UMAT 55 | ! 56 | integer i,j,stat 57 | ! 58 | real*8 F_tau(3,3),T_tau(3,3),K1,K2,K3,psi,SpTanMod(3,3,3,3) 59 | ! 60 | real*8 zero,one,two,three,fourth,third,half,six,four 61 | parameter(zero=0.d0,one=1.d0,two=2.d0,three=3.d0,four=4.d0, 62 | + six=6.d0,third=1.d0/3.d0,half=1.d0/2.d0,fourth=1.d0/4.d0) 63 | 64 | 65 | ! Do nothing if a "dummy" step 66 | ! 67 | if(dtime.eq.zero) return 68 | 69 | 70 | ! Get the deformation gradient at the end of the increment. The 71 | ! subscript tau denotes the time at the end of the increment. 72 | ! 73 | do i=1,3 74 | do j=1,3 75 | F_tau(i,j) = dfgrd1(i,j) 76 | end do 77 | end do 78 | 79 | 80 | ! At the start of an Abaqus calculation, the state 81 | ! variables are passed into UMAT with zero values. 82 | ! Initialize the state variables. At this point, 83 | ! the time total_time and step_time both have a value 84 | ! equal to zero and the step counter, kstep, is 85 | ! equal to 1. 86 | ! 87 | if ((time(1).eq.zero).and.(kstep.eq.1)) then 88 | statev(1) = zero 89 | statev(2) = zero 90 | statev(3) = zero 91 | end if 92 | 93 | 94 | !@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ 95 | ! 96 | ! Perform the constitutive update 97 | ! 98 | call foam(props,nprops,dtime,F_tau, 99 | + T_tau,K1,K2,K3,SpTanMod,psi,stat) 100 | ! 101 | if (stat.eq.0) then 102 | pnewdt = 0.5 103 | return 104 | endif 105 | ! 106 | !@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ 107 | 108 | 109 | ! Update the state variables 110 | ! 111 | statev(1) = K1 112 | statev(2) = K2 113 | statev(3) = K3 114 | 115 | 116 | ! Update the stress 117 | ! 118 | stress = zero 119 | ! 120 | do i=1,ndi 121 | stress(i) = T_tau(i,i) 122 | end do 123 | ! 124 | if (nshr.ne.0) then 125 | stress(ndi+1) = T_tau(1,2) 126 | if (nshr.ne.1) then 127 | stress(ndi+2) = T_tau(1,3) 128 | if (nshr.ne.2) then 129 | stress(ndi+3) = T_tau(2,3) 130 | endif 131 | endif 132 | endif 133 | 134 | 135 | ! Update the tangents 136 | ! 137 | ddsdde = zero 138 | ! 139 | do i=1,ndi 140 | do j=1,ndi 141 | ddsdde(i,j) = SpTanMod(i,i,j,j) 142 | end do 143 | end do 144 | ! 145 | if (nshr.ne.0) then 146 | do i=1,ndi 147 | ddsdde(i,ndi+1) = SpTanMod(i,i,1,2) 148 | ddsdde(ndi+1,i) = SpTanMod(1,2,i,i) 149 | end do 150 | ddsdde(ndi+1,ndi+1) = SpTanMod(1,2,1,2) 151 | if (nshr.ne.1) then 152 | do i=1,ndi 153 | ddsdde(i,ndi+2) = SpTanMod(i,i,1,3) 154 | ddsdde(ndi+2,i) = SpTanMod(1,3,i,i) 155 | end do 156 | ddsdde(ndi+2,ndi+2) = SpTanMod(1,3,1,3) 157 | ddsdde(ndi+1,ndi+2) = SpTanMod(1,2,1,3) 158 | ddsdde(ndi+2,ndi+1) = SpTanMod(1,3,1,2) 159 | if (nshr.ne.2) then 160 | do i=1,ndi 161 | ddsdde(i,ndi+3) = SpTanMod(i,i,2,3) 162 | ddsdde(ndi+3,i) = SpTanMod(2,3,i,i) 163 | end do 164 | ddsdde(ndi+3,ndi+3) = SpTanMod(2,3,2,3) 165 | ddsdde(ndi+1,ndi+3) = SpTanMod(1,2,2,3) 166 | ddsdde(ndi+3,ndi+1) = SpTanMod(2,3,1,2) 167 | ddsdde(ndi+2,ndi+3) = SpTanMod(1,3,2,3) 168 | ddsdde(ndi+3,ndi+2) = SpTanMod(2,3,1,3) 169 | endif 170 | endif 171 | endif 172 | 173 | 174 | ! Update the free energy 175 | ! 176 | sse = psi 177 | 178 | 179 | return 180 | end subroutine umat 181 | 182 | !************************************************************************ 183 | ! Material subroutines 184 | !************************************************************************ 185 | 186 | subroutine foam( 187 | + ! Inputs 188 | + props,nprops,dtime,F_tau, 189 | + ! Outputs 190 | + T_tau,K1,K2,K3,SpTanMod,psi,stat) 191 | 192 | implicit none 193 | ! 194 | integer i,j,k,l,m,w,s,t,nprops,stat 195 | ! 196 | real*8 props(nprops),dtime,F_tau(3,3),T_tau(3,3),K1,K2,K3, 197 | + SpTanMod(3,3,3,3),psi,Iden(3,3),G0,B,Jmin,C1,K10,delta_K,X1,X2, 198 | + C0,p,q,C2,C3,r,detF,Rot(3,3),U_tau(3,3),E_tau(3,3),B_tau(3,3), 199 | + dev_E(3,3),N(3,3),Y(3,3),f,dfdK1,dfdK1_2,X,dXdK1,dXdK1_2,Lf, 200 | + dLdK2,dLdK3,dLdK2dK2,dLdK2dK3,dpsidK1,dpsidK2,dpsidK3,dpsidK1_2, 201 | + dpsidK2_2,dpsidK2dK1,dpsidK2dK3,DTkDE(3,3,3,3),DlnBDB(3,3,3,3) 202 | ! 203 | real*8 zero,one,two,three,fourth,third,half,six,four 204 | parameter(zero=0.d0,one=1.d0,two=2.d0,three=3.d0,four=4.d0, 205 | + six=6.d0,third=1.d0/3.d0,half=1.d0/2.d0,fourth=1.d0/4.d0) 206 | 207 | 208 | ! Identity matrix 209 | ! 210 | call onem(Iden) 211 | 212 | 213 | ! Obtain material properties 214 | ! 215 | G0 = props(1) ! Shear modulus 216 | B = props(2) ! Bulk modulus 217 | Jmin = props(3) ! Locking value of J in f-fuction 218 | C1 = props(4) ! C_1 in L-function 219 | K10 = props(5) ! K1_0 in X-function 220 | delta_K = props(6) ! delta_K in X-function 221 | X1 = props(7) ! X1' in X-function 222 | X2 = props(8) ! X2' in X-function 223 | C0 = props(9) ! C_0 in L-function 224 | p = props(10) ! power p in L-function 225 | q = props(11) ! power q in L-function 226 | C2 = props(12) ! C_2 in f-function 227 | C3 = props(13)! C_3 in f-function 228 | r = props(14) ! power r in f-function 229 | 230 | 231 | ! Calulate kinematic quantities 232 | ! 233 | call mdet(F_tau,detF) 234 | call skinem(F_tau,Rot,U_tau,E_tau,stat) 235 | E_tau = matmul(matmul(Rot,E_tau),transpose(Rot)) 236 | B_tau = matmul(F_tau,transpose(F_tau)) 237 | 238 | 239 | ! Calulate the invariants K1, K2, K3 240 | ! 241 | K1 = E_tau(1,1) + E_tau(2,2) + E_tau(3,3) 242 | dev_E = E_tau - third*K1*Iden 243 | K2 = dsqrt(sum(dev_E*dev_E)) 244 | if (K2.eq.zero) then 245 | ! if K2 == 0, assign an arbitrary value to N, K3, and Y 246 | N = zero 247 | K3 = zero 248 | Y = zero 249 | else 250 | ! else compute N, K3, and Y 251 | N = dev_E/K2 252 | call mdet(N,K3) 253 | K3 = K3*three*dsqrt(six) 254 | Y = three*dsqrt(six)*matmul(N,N) 255 | + - dsqrt(six)*Iden - three*K3*N 256 | end if 257 | 258 | 259 | ! Compute the free energy function and its derivatives 260 | ! 261 | call f_fun(f,dfdK1,dfdK1_2,K1,Jmin,C2,C3,r) 262 | call X_fun(X,dXdK1,dXdK1_2,K1,X1,X2,K10,delta_K) 263 | call L_fun(Lf,dLdK2,dLdK3,dLdK2dK2,dLdK2dK3,K2,K3,C0,C1,p,q) 264 | dpsidK1 = G0*dXdK1*K2**two + B*dfdK1 265 | dpsidK2 = G0*(two*X*K2 + dLdK2) 266 | dpsidK3 = G0*dLdK3 267 | dpsidK1_2 = G0*dXdK1_2*K2**two + B*dfdK1_2 268 | dpsidK2_2 = G0*(two*X + dLdK2dK2) 269 | dpsidK2dK1 = G0*two*dXdK1*K2 270 | dpsidK2dK3 = G0*dLdK2dK3 271 | psi = G0*(X*(K2**two) + Lf) + B*f 272 | 273 | 274 | ! Compute the Cauchy stress and the material contribution 275 | ! to the tangent (derivative of Kirchhoff stress wrt log strain) 276 | ! 277 | if (K2.eq.zero) then 278 | ! 279 | T_tau = (one/detF)*(dPSIdK1*Iden) 280 | DTkDE = zero 281 | do i=1,3 282 | do j=1,3 283 | do k=1,3 284 | do l=1,3 285 | DTkDE(i,j,k,l) = DTkDE(i,j,k,l) 286 | + + two*B*dfdK1_2*Iden(i,j)*Iden(k,l) 287 | + + two*G0*(half*Iden(i,k)*Iden(j,l) 288 | + + half*Iden(i,l)*Iden(j,k) 289 | + - third*Iden(i,j)*Iden(k,l)) 290 | end do 291 | end do 292 | end do 293 | end do 294 | ! 295 | else 296 | ! 297 | T_tau = (one/detF)*(dpsidK1*Iden + dpsidK2*N + dpsidK3/K2*Y) 298 | DTkDE = zero 299 | do i=1,3 300 | do j=1,3 301 | do k=1,3 302 | do l=1,3 303 | DTkDE(i,j,k,l) = DTkDE(i,j,k,l) + 304 | + (dev_E(k,l)*Iden(i,j) + Iden(k,l)*dev_E(i,j))* 305 | + (-two*dsqrt(six)/K2**three*dpsidK3 + one/K2*dpsidK2dK1)+ 306 | + (dev_E(k,l)*Y(i,j) + Y(k,l)*dev_E(i,j))* 307 | + (one/K2**two*dpsidK2dK3 - three/K2**three*dpsidK3)+ 308 | + dev_E(k,l)*dev_E(i,j)* 309 | + (-one/K2**three*dpsidK2 + one/K2**two*dpsidK2_2 - 310 | + three*K3/K2**four*dpsidK3)+ 311 | + (Iden(i,k)*dev_E(l,j) + Iden(j,l)*dev_E(i,k) + 312 | + Iden(i,l)*dev_E(k,j)+Iden(k,j)*dev_E(i,l)) 313 | + *three/two*dsqrt(six)/K2**three*dpsidK3 314 | + +(half*Iden(i,k)*Iden(j,l) + half*Iden(i,l)*Iden(j,k) - 315 | + third*Iden(i,j)*Iden(k,l))* 316 | + (one/K2*dpsidK2 - three*K3/K2**two*dpsidK3)+ 317 | + Iden(i,j)*Iden(k,l)*dpsidK1_2 318 | end do 319 | end do 320 | end do 321 | end do 322 | ! 323 | end if 324 | 325 | 326 | ! Compute the fourth order tensor dlnBDB 327 | ! 328 | call dlnxdx(B_tau,DlnBDB) 329 | ! 330 | ! Compute the tangent required by Abaqus 331 | ! 332 | SpTanMod = zero 333 | do i=1,3 334 | do j=1,3 335 | do k=1,3 336 | do l=1,3 337 | do m=1,3 338 | do w=1,3 339 | do t=1,3 340 | do s=1,3 341 | SpTanMod(i,j,k,l) = SpTanMod(i,j,k,l) + 342 | + fourth/detF*(B_tau(t,l)*Iden(k,s)+B_tau(s,l)*Iden(t,k)+ 343 | + B_tau(t,k)*Iden(l,s)+B_tau(s,k)*Iden(t,l)) 344 | + *DTkDE(i,j,m,w)*DlnBDB(m,w,t,s) 345 | end do 346 | end do 347 | end do 348 | end do 349 | end do 350 | end do 351 | end do 352 | end do 353 | 354 | 355 | return 356 | end subroutine foam 357 | 358 | !**************************************************************************** 359 | 360 | subroutine f_fun(f,dfdK1,dfdK1_2,K1,Jmin,C2,C3,r) 361 | ! 362 | ! This subroutine calculates the f-function and its derivatives 363 | ! 364 | implicit none 365 | ! 366 | real*8 f,dfdK1,dfdK1_2,K1,Jmin,B,C2,C3,r,J 367 | ! 368 | real*8 zero,one,two,three,fourth,third,half,six,four 369 | parameter(zero=0.d0,one=1.d0,two=2.d0,three=3.d0,four=4.d0, 370 | + six=6.d0,third=1.d0/3.d0,half=1.d0/2.d0,fourth=1.d0/4.d0) 371 | 372 | 373 | J = dexp(K1) 374 | 375 | 376 | ! Compute f 377 | ! 378 | f = (dexp(C2*(K1))-C2*K1 - one)/(C2**two) + (C3/(r-one))* 379 | + (-J**(one-r) + (J-Jmin)**(one-r)*(one-Jmin)**r + Jmin) 380 | 381 | 382 | ! Compute first derivative of f wrt K1 383 | ! 384 | dfdK1 = (dexp(C2*(K1))-one)/C2 + 385 | + C3*J*(J**(-r)-(one-Jmin)**r/(J-Jmin)**r) 386 | 387 | 388 | ! Compute second derivative of f wrt K1 389 | ! 390 | dfdK1_2 = dexp(C2*K1) + C3*(J**(-r)-(one-Jmin)**(r) 391 | + /(J-Jmin)**r)*dexp(K1) 392 | + + C3*(-r*J**(-r-one) + r*(one-Jmin)**(r)/(J-Jmin)**(r+one)) 393 | + *dexp(K1)*dexp(K1) 394 | 395 | 396 | return 397 | end subroutine f_fun 398 | 399 | !**************************************************************************** 400 | 401 | subroutine X_fun(X,dXdK1,dXdK1_2,K1,X1,X2,K10,delta_K) 402 | ! 403 | ! This subroutine calculates the X-function and its derivatives 404 | ! 405 | implicit none 406 | ! 407 | real*8 K1,X,dXdK1,X1,X2,delta_K,K10,dXdK1_2 408 | ! 409 | real*8 zero,one,two,three,fourth,third,half,pi 410 | parameter(zero=0.d0,one=1.d0,two=2.d0,three=3.d0,fourth=1.d0/4.d0, 411 | + third=1.d0/3.d0,half=1.d0/2.d0, pi=3.141592653d0) 412 | 413 | 414 | ! Compute X 415 | X = (X1+X2)/two*K1 + (X1-X2)/two*delta_K* 416 | + dlog(dcosh((K1-K10)/delta_K)/dcosh(K10/delta_K)) + one 417 | 418 | 419 | ! Compute first derivative of X wrt K1 420 | ! 421 | dXdK1 = (X1+X2)/two + (X1-X2)/two*dtanh((K1 - K10)/delta_K) 422 | 423 | 424 | ! Compute second derivative of X wrt K1 425 | ! 426 | dXdK1_2 = (X1-X2)/two/dcosh((K10-K1)/delta_K)**two/delta_K 427 | 428 | 429 | return 430 | end subroutine X_fun 431 | 432 | !**************************************************************************** 433 | 434 | subroutine L_fun(Lf,dLdK2,dLdK3,dLdK2dK2,dLdK2dK3,K2,K3,C0,C1,p,q) 435 | ! 436 | ! This subroutine calculates the L-function and its derivatives 437 | ! 438 | implicit none 439 | ! 440 | real*8 Lf,dLdK2,dLdK3,dLdK2dK2,dLdK2dK3,K2,K3,C1,C0,p,q 441 | ! 442 | real*8 zero,one,two,three,fourth,third,half,pi 443 | parameter(zero=0.d0,one=1.d0,two=2.d0,three=3.d0,fourth=1.d0/4.d0, 444 | + third=1.d0/3.d0,half=1.d0/2.d0, pi=3.141592653d0) 445 | 446 | 447 | 448 | ! Compute L 449 | ! 450 | Lf = C0*K2**p + C1*(one+K3)*K2**q 451 | 452 | 453 | ! Compute first derivative of L wrt K2 454 | ! 455 | dLdK2 = C0*p*K2**(p-one) + C1*(one+K3)*q*K2**(q-one) 456 | 457 | 458 | ! Compute first derivative of L wrt K3 459 | ! 460 | dLdK3 = C1*K2**q 461 | 462 | 463 | ! Compute second derivative of L wrt K2 464 | ! 465 | dLdK2dK2 = C0*p*(p-one)*K2**(p-two) + 466 | + C1*(one+K3)*q*(q-one)*K2**(q-two) 467 | 468 | 469 | ! Compute mixed derivative of L wrt K2 and K3 470 | ! 471 | dLdK2dK3 = C1*q*K2**(q-one) 472 | 473 | 474 | return 475 | end subroutine L_fun 476 | 477 | !**************************************************************************** 478 | ! The next subroutine calculates various kinematical quantities 479 | ! associated with the deformation gradient 480 | !**************************************************************************** 481 | 482 | subroutine skinem(F,R,U,E,istat) 483 | ! 484 | ! This subroutine performs the right polar decomposition 485 | ! F = RU of the deformation gradient F into a rotation 486 | ! R and the right stretch tensor U. The logarithmic 487 | ! strain E = ln(U) is also returned. 488 | ! 489 | ! F(3,3): the deformation gradient; input 490 | ! detF: the determinant of F; detF > 0 491 | ! R(3,3): the rotation matrix; output 492 | ! U(3,3): the right stretch tensor; output 493 | ! Uinv(3,3): the inverse of U 494 | ! C(3,3): the right Cauchy-Green tensor 495 | ! omega(3): the squares of the principal stretches 496 | ! Ueigval(3): the principal stretches 497 | ! eigvec(3,3): matrix of eigenvectors of U 498 | ! E(3,3): the logarithmic strain tensor; output 499 | ! istat: success flag, istat=0 for a failed attempt; output 500 | ! 501 | implicit none 502 | ! 503 | integer istat 504 | ! 505 | real*8 F(3,3),C(3,3),omega(3),Ueigval(3),eigvec(3,3), 506 | + U(3,3),E(3,3),Uinv(3,3),R(3,3),detF 507 | 508 | 509 | ! Store the identity matrix in R, U, and Uinv 510 | ! 511 | call onem(R) 512 | call onem(U) 513 | call onem(Uinv) 514 | 515 | 516 | ! Store the zero matrix in E 517 | ! 518 | E = 0.d0 519 | 520 | 521 | ! Check if the determinant of F is greater than zero. 522 | ! If not, then print a diagnostic and cut back the 523 | ! time increment. 524 | ! 525 | call mdet(F,detF) 526 | if (detF.le.0.d0) then 527 | write(*,'(/5X,A/)') '--problem in kinematics-- the', 528 | + ' determinant of F is not greater than 0' 529 | istat = 0 530 | return 531 | end if 532 | 533 | 534 | ! Calculate the right Cauchy-Green tensor C 535 | ! 536 | C = matmul(transpose(F),F) 537 | 538 | 539 | ! Calculate the eigenvalues and eigenvectors of C 540 | ! 541 | call spectral(C,omega,eigvec,istat) 542 | 543 | 544 | ! Calculate the principal values of U and E 545 | ! 546 | Ueigval(1) = dsqrt(omega(1)) 547 | Ueigval(2) = dsqrt(omega(2)) 548 | Ueigval(3) = dsqrt(omega(3)) 549 | ! 550 | U(1,1) = Ueigval(1) 551 | U(2,2) = Ueigval(2) 552 | U(3,3) = Ueigval(3) 553 | ! 554 | E(1,1) = dlog(Ueigval(1)) 555 | E(2,2) = dlog(Ueigval(2)) 556 | E(3,3) = dlog(Ueigval(3)) 557 | 558 | 559 | ! Calculate the complete tensors U and E 560 | ! 561 | U = matmul(matmul(eigvec,U),transpose(eigvec)) 562 | E = matmul(matmul(eigvec,E),transpose(eigvec)) 563 | 564 | 565 | ! Calculate Uinv 566 | ! 567 | call matInv3D(U,Uinv,detF,istat) 568 | 569 | 570 | ! calculate R 571 | ! 572 | R = matmul(F,Uinv) 573 | 574 | 575 | return 576 | end subroutine skinem 577 | 578 | !**************************************************************************** 579 | ! The following subroutines calculate the spectral 580 | ! decomposition of a symmetric 3 by 3 matrix 581 | !**************************************************************************** 582 | 583 | subroutine spectral(A,D,V,istat) 584 | ! 585 | ! This subroutine calculates the eigenvalues and eigenvectors of 586 | ! a symmetric 3 by 3 matrix A. 587 | ! 588 | ! The output consists of a vector D containing the three 589 | ! eigenvalues in ascending order, and a matrix V whose 590 | ! columns contain the corresponding eigenvectors. 591 | ! 592 | implicit none 593 | ! 594 | integer np,nrot,i,j,istat 595 | parameter(np=3) 596 | ! 597 | real*8 D(3),V(3,3),A(3,3),E(3,3) 598 | 599 | 600 | E = A 601 | ! 602 | call jacobi(E,3,np,D,V,nrot,istat) 603 | call eigsrt(D,V,3,np) 604 | 605 | 606 | return 607 | end subroutine spectral 608 | 609 | !**************************************************************************** 610 | 611 | subroutine jacobi(A,n,np,D,V,nrot,istat) 612 | ! 613 | ! Computes all eigenvalues and eigenvectors of a real symmetric 614 | ! matrix A, which is of size n by n, stored in a physical 615 | ! np by np array. On output, elements of A above the diagonal 616 | ! are destroyed, but the diagonal and sub-diagonal are unchanged 617 | ! and give full information about the original symmetric matrix. 618 | ! Vector D returns the eigenvalues of A in its first n elements. 619 | ! V is a matrix with the same logical and physical dimensions as 620 | ! A whose columns contain, upon output, the normalized 621 | ! eigenvectors of A. nrot returns the number of Jacobi rotation 622 | ! which were required. 623 | ! 624 | ! This subroutine is taken from 'Numerical Recipes.' 625 | ! 626 | implicit none 627 | ! 628 | integer ip,iq,n,nmax,np,nrot,i,j,istat 629 | parameter (nmax=100) 630 | ! 631 | real*8 A(np,np),D(np),V(np,np),B(nmax),Z(nmax), 632 | + sm,tresh,G,T,H,theta,S,C,tau 633 | 634 | 635 | ! Initialize V to the identity matrix 636 | ! 637 | call onem(V) 638 | 639 | 640 | ! Initialize B and D to the diagonal of A, and Z to zero. 641 | ! The vector Z will accumulate terms of the form T*A_PQ as 642 | ! in equation (11.1.14) 643 | ! 644 | do ip = 1,n 645 | B(ip) = A(ip,ip) 646 | D(ip) = B(ip) 647 | Z(ip) = 0.d0 648 | end do 649 | 650 | 651 | ! Begin iteration 652 | ! 653 | nrot = 0 654 | do i=1,50 655 | ! 656 | ! Sum off-diagonal elements 657 | ! 658 | sm = 0.d0 659 | do ip=1,n-1 660 | do iq=ip+1,n 661 | sm = sm + dabs(A(ip,iq)) 662 | end do 663 | end do 664 | ! 665 | ! If sm = 0., then return. This is the normal return, 666 | ! which relies on quadratic convergence to machine 667 | ! underflow. 668 | ! 669 | if (sm.eq.0.d0) return 670 | ! 671 | ! In the first three sweeps carry out the PQ rotation only if 672 | ! |A_PQ| > tresh, where tresh is some threshold value, 673 | ! see equation (11.1.25). Thereafter tresh = 0. 674 | ! 675 | if (i.lt.4) then 676 | tresh = 0.2d0*sm/n**2 677 | else 678 | tresh = 0.d0 679 | end if 680 | ! 681 | do ip=1,n-1 682 | do iq=ip+1,n 683 | G = 100.d0*dabs(A(ip,iq)) 684 | ! 685 | ! After four sweeps, skip the rotation if the 686 | ! off-diagonal element is small. 687 | ! 688 | if ((i.gt.4).and.(dabs(D(ip))+G.eq.dabs(D(ip))) 689 | + .and.(dabs(D(iq))+G.eq.dabs(D(iq)))) then 690 | A(ip,iq) = 0.d0 691 | else if (dabs(A(ip,iq)).gt.tresh) then 692 | H = D(iq) - D(ip) 693 | if (dabs(H)+G.eq.dabs(H)) then 694 | ! 695 | ! T = 1./(2.*theta), equation (11.1.10) 696 | ! 697 | T =A(ip,iq)/H 698 | else 699 | theta = 0.5d0*H/A(ip,iq) 700 | T =1.d0/(dabs(theta)+dsqrt(1.d0+theta**2.d0)) 701 | if (theta.lt.0.d0) T = -T 702 | end if 703 | C = 1.d0/dsqrt(1.d0 + T**2.d0) 704 | S = T*C 705 | tau = S/(1.d0 + C) 706 | H = T*A(ip,iq) 707 | Z(ip) = Z(ip) - H 708 | Z(iq) = Z(iq) + H 709 | D(ip) = D(ip) - H 710 | D(iq) = D(iq) + H 711 | A(ip,iq) = 0.d0 712 | ! 713 | ! Case of rotations 1 <= J < P 714 | ! 715 | do j=1,ip-1 716 | G = A(j,ip) 717 | H = A(j,iq) 718 | A(j,ip) = G - S*(H + G*tau) 719 | A(j,iq) = H + S*(G - H*tau) 720 | end do 721 | ! 722 | ! Case of rotations P < J < Q 723 | ! 724 | do j=ip+1,iq-1 725 | G = A(ip,j) 726 | H = A(j,iq) 727 | A(ip,j) = G - S*(H + G*tau) 728 | A(j,iq) = H + S*(G - H*tau) 729 | end do 730 | ! 731 | ! Case of rotations Q < J <= N 732 | ! 733 | do j=iq+1,n 734 | G = A(ip,j) 735 | H = A(iq,j) 736 | A(ip,j) = G - S*(H + G*tau) 737 | A(iq,j) = H + S*(G - H*tau) 738 | end do 739 | do j = 1,n 740 | G = V(j,ip) 741 | H = V(j,iq) 742 | V(j,ip) = G - S*(H + G*tau) 743 | V(j,iq) = H + S*(G - H*tau) 744 | end do 745 | nrot = nrot + 1 746 | end if 747 | end do 748 | end do 749 | ! 750 | ! Update D with the sum of T*A_PQ, and reinitialize Z 751 | ! 752 | do ip=1,n 753 | B(ip) = B(ip) + Z(ip) 754 | D(ip) = B(ip) 755 | Z(ip) = 0.d0 756 | end do 757 | end do 758 | 759 | 760 | ! If the algorithm has reached this stage, then there 761 | ! are too many sweeps. Print a diagnostic and cut the 762 | ! time increment. 763 | ! 764 | write (*,'(/1X,A/)') '50 iterations in jacobi should never happen' 765 | istat = 0 766 | 767 | 768 | return 769 | end subroutine jacobi 770 | 771 | !**************************************************************************** 772 | 773 | subroutine eigsrt(D,V,n,np) 774 | ! 775 | ! Given the eigenvalues D and eigenvectors V as output from 776 | ! jacobi, this subroutine sorts the eigenvales into ascending 777 | ! order and rearranges the colmns of V accordingly. 778 | ! 779 | ! The subroutine is taken from 'Numerical Recipes.' 780 | ! 781 | implicit none 782 | ! 783 | integer n,np,i,j,k 784 | ! 785 | real*8 D(np),V(np,np),P 786 | 787 | 788 | do i=1,n-1 789 | k = i 790 | P = D(i) 791 | do j=i+1,n 792 | if (D(j).ge.P) then 793 | k = j 794 | P = D(j) 795 | end if 796 | end do 797 | if (k.ne.i) then 798 | D(k) = D(i) 799 | D(i) = P 800 | do j=1,n 801 | P = V(j,i) 802 | V(j,i) = V(j,k) 803 | V(j,k) = P 804 | end do 805 | end if 806 | end do 807 | 808 | 809 | return 810 | end subroutine eigsrt 811 | 812 | !**************************************************************************** 813 | ! Utility subroutines 814 | !**************************************************************************** 815 | 816 | subroutine matInv3D(A,A_inv,det_A,istat) 817 | ! 818 | ! Returns A_inv, the inverse and det_A, the determinant 819 | ! Note that the det is of the original matrix, not the 820 | ! inverse 821 | ! 822 | implicit none 823 | ! 824 | integer istat 825 | ! 826 | real*8 A(3,3),A_inv(3,3),det_A,det_A_inv 827 | 828 | 829 | istat = 1 830 | 831 | det_A = A(1,1)*(A(2,2)*A(3,3) - A(3,2)*A(2,3)) - 832 | + A(2,1)*(A(1,2)*A(3,3) - A(3,2)*A(1,3)) + 833 | + A(3,1)*(A(1,2)*A(2,3) - A(2,2)*A(1,3)) 834 | 835 | if (det_A .le. 0.d0) then 836 | write(*,*) 'WARNING: SUBROUTINE matInv3:' 837 | write(*,*) 'WARNING: DET of MAT=',DET_A 838 | istat = 0 839 | return 840 | end if 841 | 842 | det_A_inv = 1.d0/det_A 843 | 844 | A_inv(1,1) = det_A_inv*(A(2,2)*A(3,3)-A(3,2)*A(2,3)) 845 | A_inv(1,2) = det_A_inv*(A(3,2)*A(1,3)-A(1,2)*A(3,3)) 846 | A_inv(1,3) = det_A_inv*(A(1,2)*A(2,3)-A(2,2)*A(1,3)) 847 | A_inv(2,1) = det_A_inv*(A(3,1)*A(2,3)-A(2,1)*A(3,3)) 848 | A_inv(2,2) = det_A_inv*(A(1,1)*A(3,3)-A(3,1)*A(1,3)) 849 | A_inv(2,3) = det_A_inv*(A(2,1)*A(1,3)-A(1,1)*A(2,3)) 850 | A_inv(3,1) = det_A_inv*(A(2,1)*A(3,2)-A(3,1)*A(2,2)) 851 | A_inv(3,2) = det_A_inv*(A(3,1)*A(1,2)-A(1,1)*A(3,2)) 852 | A_inv(3,3) = det_A_inv*(A(1,1)*A(2,2)-A(2,1)*A(1,2)) 853 | 854 | 855 | return 856 | end subroutine matInv3D 857 | 858 | !**************************************************************************** 859 | 860 | subroutine mdet(A,det) 861 | ! 862 | ! This subroutine calculates the determinant 863 | ! of a 3 by 3 matrix [A] 864 | ! 865 | implicit none 866 | ! 867 | real*8 A(3,3),det 868 | 869 | 870 | det = A(1,1)*A(2,2)*A(3,3) 871 | + + A(1,2)*A(2,3)*A(3,1) 872 | + + A(1,3)*A(2,1)*A(3,2) 873 | + - A(3,1)*A(2,2)*A(1,3) 874 | + - A(3,2)*A(2,3)*A(1,1) 875 | + - A(3,3)*A(2,1)*A(1,2) 876 | 877 | 878 | return 879 | end subroutine mdet 880 | 881 | !**************************************************************************** 882 | 883 | subroutine onem(A) 884 | ! 885 | ! This subroutine stores the identity matrix in the 886 | ! 3 by 3 matrix [A] 887 | ! 888 | implicit none 889 | ! 890 | integer i,j 891 | ! 892 | real*8 A(3,3) 893 | 894 | 895 | do i=1,3 896 | do J=1,3 897 | if (i .eq. j) then 898 | A(i,j) = 1.0 899 | else 900 | A(i,j) = 0.0 901 | end if 902 | end do 903 | end do 904 | 905 | 906 | return 907 | end subroutine onem 908 | 909 | !**************************************************************************** 910 | 911 | subroutine dlnxdx(X,DYDX) 912 | ! 913 | ! This subroutine calculates the derivative of the logarithm 914 | ! of a symmetric tensor with respect to that tensor 915 | ! 916 | implicit none 917 | ! 918 | integer i,j,k,l 919 | ! 920 | real*8 X(3,3),DYDX(3,3,3,3),Iden(3,3),Iden4(3,3,3,3),eigval(3), 921 | + eigvec(3,3),ehat1(3),ehat2(3),ehat3(3),E1(3,3),E2(3,3),E3(3,3), 922 | + y(3),DX2DX(3,3,3,3),s1,s2,s3,s4,s5,s6 923 | ! 924 | real*8 zero,one,two,half,three,third,small 925 | parameter(zero=0.d0,one=1.d0,two=2.d0,half=0.5d0, 926 | + three=3.d0,third=1.d0/3.d0,small=1.d-12) 927 | 928 | 929 | ! Initialize 930 | ! 931 | ! Second order identity tensor 932 | ! 933 | call onem(Iden) 934 | ! 935 | ! Fourth order symmetric identity tensor 936 | ! 937 | Iden4 = zero 938 | do i=1,3 939 | do j=1,3 940 | do k=1,3 941 | do l=1,3 942 | Iden4(i,j,k,l) = half*(Iden(i,k)*Iden(j,l) + 943 | + Iden(i,l)*Iden(j,k)) 944 | end do 945 | end do 946 | end do 947 | end do 948 | 949 | 950 | ! Calculate the eigenvalues and eigenvectors of X 951 | ! 952 | call spectral(X,eigval,eigvec) 953 | ! 954 | ! Extract the eigenvectors 955 | ! 956 | do i=1,3 957 | ehat1(i) = eigvec(i,1) 958 | ehat2(i) = eigvec(i,2) 959 | ehat3(i) = eigvec(i,3) 960 | end do 961 | ! 962 | ! Assemble the eigenprojections 963 | ! 964 | do i=1,3 965 | do j=1,3 966 | E1(i,j) = ehat1(i)*ehat1(j) 967 | E2(i,j) = ehat2(i)*ehat2(j) 968 | E3(i,j) = ehat3(i)*ehat3(j) 969 | end do 970 | end do 971 | 972 | 973 | ! Calculate the eigenvalues of Y = ln(X) 974 | ! 975 | y = dlog(eigval) 976 | 977 | 978 | ! Calculate the derivative of X^2 with respect to X 979 | ! 980 | DX2DX = zero 981 | do i=1,3 982 | do j=1,3 983 | do k=1,3 984 | do l=1,3 985 | DX2DX(i,j,k,l) = half*(Iden(i,k)*X(j,l) + 986 | + Iden(i,l)*X(j,k) + 987 | + X(i,k)*Iden(j,l) + 988 | + X(i,l)*Iden(j,k)) 989 | end do 990 | end do 991 | end do 992 | end do 993 | 994 | 995 | ! Calculate DYDX 996 | ! 997 | DYDX = zero 998 | if (dabs(eigval(1)-eigval(3)).le.small) then 999 | ! 1000 | ! Three repeated eigenvalues 1001 | ! 1002 | do i=1,3 1003 | do j=1,3 1004 | do k=1,3 1005 | do l=1,3 1006 | DYDX(i,j,k,l) = (one/eigval(1))*Iden4(i,j,k,l) 1007 | end do 1008 | end do 1009 | end do 1010 | end do 1011 | ! 1012 | elseif (dabs(eigval(2)-eigval(3)).le.small) then 1013 | ! 1014 | ! The eigenvalues 2 and 3 are repeated. Eigenvalue 1 is distinct. 1015 | ! 1016 | s1 = (y(1) - y(2))/((eigval(1)-eigval(2))**two) - 1017 | + (one/eigval(2))/(eigval(1)-eigval(2)) 1018 | s2 = two*eigval(2)*(y(1)-y(2))/((eigval(1)-eigval(2))**two) - 1019 | + (one/eigval(2))*(eigval(1)+eigval(2))/(eigval(1)-eigval(2)) 1020 | s3 = two*(y(1)-y(2))/((eigval(1)-eigval(2))**three) - 1021 | + ((one/eigval(1))+(one/eigval(2)))/((eigval(1)-eigval(2))**two) 1022 | s4 = eigval(2)*s3 1023 | s5 = s4 1024 | s6 = (eigval(2)**two)*s3 1025 | ! 1026 | do i=1,3 1027 | do j=1,3 1028 | do k=1,3 1029 | do l=1,3 1030 | DYDX(i,j,k,l) = s1*DX2DX(i,j,k,l) - s2*Iden4(i,j,k,l) - 1031 | + s3*X(i,j)*X(k,l) + s4*X(i,j)*Iden(k,l) + 1032 | + s5*Iden(i,j)*X(k,l) - s6*Iden(i,j)*Iden(k,l) 1033 | end do 1034 | end do 1035 | end do 1036 | end do 1037 | ! 1038 | elseif (dabs(eigval(1)-eigval(2)).le.small) then 1039 | ! 1040 | ! The eigenvalues 1 and 2 are repeated. Eigenvalue 3 is distinct. 1041 | ! 1042 | s1 = (y(3) - y(2))/((eigval(3)-eigval(2))**two) - 1043 | + (one/eigval(2))/(eigval(3)-eigval(2)) 1044 | s2 = two*eigval(2)*(y(3)-y(2))/((eigval(3)-eigval(2))**two) - 1045 | + (one/eigval(2))*(eigval(3)+eigval(2))/(eigval(3)-eigval(2)) 1046 | s3 = two*(y(3)-y(2))/((eigval(3)-eigval(2))**three) - 1047 | + ((one/eigval(3))+(one/eigval(2)))/((eigval(3)-eigval(2))**two) 1048 | s4 = eigval(2)*s3 1049 | s5 = s4 1050 | s6 = (eigval(2)**two)*s3 1051 | ! 1052 | do i=1,3 1053 | do j=1,3 1054 | do k=1,3 1055 | do l=1,3 1056 | DYDX(i,j,k,l) = s1*DX2DX(i,j,k,l) - s2*Iden4(i,j,k,l) - 1057 | + s3*X(i,j)*X(k,l) + s4*X(i,j)*Iden(k,l) + 1058 | + s5*Iden(i,j)*X(k,l) - s6*Iden(i,j)*Iden(k,l) 1059 | end do 1060 | end do 1061 | end do 1062 | end do 1063 | ! 1064 | else 1065 | ! 1066 | ! Eigenvalues are distinct. 1067 | ! 1068 | do i=1,3 1069 | do j=1,3 1070 | do k=1,3 1071 | do l=1,3 1072 | DYDX(i,j,k,l) = (y(1)/((eigval(1)-eigval(2))* 1073 | + (eigval(1)-eigval(3))))* 1074 | + (DX2DX(i,j,k,l) - (eigval(2)+eigval(3))*Iden4(i,j,k,l) - 1075 | + ((eigval(1)-eigval(2))+(eigval(1)-eigval(3)))*E1(i,j)*E1(k,l) - 1076 | + (eigval(2)-eigval(3))*(E2(i,j)*E2(k,l)-E3(i,j)*E3(k,l))) + 1077 | + (one/eigval(1))*E1(i,j)*E1(k,l) + 1078 | + (y(2)/((eigval(2)-eigval(1))* 1079 | + (eigval(2)-eigval(3))))* 1080 | + (DX2DX(i,j,k,l) - (eigval(1)+eigval(3))*Iden4(i,j,k,l) - 1081 | + ((eigval(2)-eigval(1))+(eigval(2)-eigval(3)))*E2(i,j)*E2(k,l) - 1082 | + (eigval(1)-eigval(3))*(E1(i,j)*E1(k,l)-E3(i,j)*E3(k,l))) + 1083 | + (one/eigval(2))*E2(i,j)*E2(k,l) + 1084 | + (y(3)/((eigval(3)-eigval(1))* 1085 | + (eigval(3)-eigval(2))))* 1086 | + (DX2DX(i,j,k,l) - (eigval(1)+eigval(2))*Iden4(i,j,k,l) - 1087 | + ((eigval(3)-eigval(1))+(eigval(3)-eigval(2)))*E3(i,j)*E3(k,l) - 1088 | + (eigval(1)-eigval(2))*(E1(i,j)*E1(k,l)-E2(i,j)*E2(k,l))) + 1089 | + (one/eigval(3))*E3(i,j)*E3(k,l) 1090 | end do 1091 | end do 1092 | end do 1093 | end do 1094 | ! 1095 | end if 1096 | 1097 | return 1098 | end subroutine dlnxdx 1099 | 1100 | !**************************************************************************** 1101 | 1102 | subroutine mprod4(A,B,C) 1103 | ! 1104 | ! This subroutine calculates the product of two fourth order tensors, 1105 | ! A and B, and places the product in C: 1106 | ! C_{ijkl} = A_{ijmn}B_{mnkl}. 1107 | ! 1108 | implicit none 1109 | ! 1110 | integer i,j,k,l,m,n 1111 | ! 1112 | real*8 A(3,3,3,3),B(3,3,3,3),C(3,3,3,3) 1113 | 1114 | 1115 | C = 0.d0 1116 | do i=1,3 1117 | do j=1,3 1118 | do k=1,3 1119 | do l=1,3 1120 | do m=1,3 1121 | do n=1,3 1122 | C(i,j,k,l) = C(i,j,k,l) + A(i,j,m,n)*B(m,n,k,l) 1123 | end do 1124 | end do 1125 | end do 1126 | end do 1127 | end do 1128 | end do 1129 | 1130 | 1131 | return 1132 | end subroutine mprod4 1133 | 1134 | !**************************************************************************** --------------------------------------------------------------------------------