FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
utilities.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 !-------------------------------------------------------------------------------
6 module m_utilities
7  use hecmw
8  implicit none
9 
10  real(kind=kreal), parameter, private :: pi=3.14159265358979d0
11 
12 contains
13 
15  subroutine memget(var,dimn,syze)
16  integer :: var,dimn,syze,bite
17  parameter(bite=1)
18  var = var + dimn*syze*bite
19  end subroutine memget
20 
21 
23  subroutine append_int2name( n, fname, n1 )
24  integer, intent(in) :: n
25  integer, intent(in), optional :: n1
26  character(len=*), intent(inout) :: fname
27  integer :: npos, nlen
28  character(len=128) :: tmpname, tmp
29 
30  npos = scan( fname, '.')
31  nlen = len_trim( fname )
32  if( nlen>128 ) stop "String too long(>128) in append_int2name"
33  if( n>100000 ) stop "Integer too big>100000 in append_int2name"
34  tmpname = fname
35  if( npos==0 ) then
36  write( fname, '(a,i6)') fname(1:nlen),n
37  else
38  write( tmp, '(i6,a)') n,tmpname(npos:nlen)
39  fname = tmpname(1:npos-1) // adjustl(tmp)
40  endif
41  if(present(n1).and.n1/=0)then
42  write(tmp,'(i8)')n1
43  fname = fname(1:len_trim(fname))//'.'//adjustl(tmp)
44  endif
45  end subroutine
46 
48  subroutine insert_int2array( iin, carray )
49  integer, intent(in) :: iin
50  integer, pointer :: carray(:)
51 
52  integer :: i, oldsize
53  integer, pointer :: dumarray(:) => null()
54  if( .not. associated(carray) ) then
55  allocate( carray(1) )
56  carray(1) = iin
57  else
58  oldsize = size( carray )
59  allocate( dumarray(oldsize) )
60  do i=1,oldsize
61  dumarray(i) = carray(i)
62  enddo
63  deallocate( carray )
64  allocate( carray(oldsize+1) )
65  do i=1,oldsize
66  carray(i) = dumarray(i)
67  enddo
68  carray(oldsize+1) = iin
69  endif
70  if( associated(dumarray) ) deallocate( dumarray )
71  end subroutine
72 
74  subroutine tensor_eigen3( tensor, eigval )
75  real(kind=kreal), intent(in) :: tensor(6)
76  real(kind=kreal), intent(out) :: eigval(3)
77 
78  real(kind=kreal) :: i1,i2,i3,r,theta,q, x(3,3), xx(3,3)
79  real(kind=kreal), parameter :: small = 1.0d-10
80 
81  x(1,1)=tensor(1); x(2,2)=tensor(2); x(3,3)=tensor(3)
82  x(1,2)=tensor(4); x(2,1)=x(1,2)
83  x(2,3)=tensor(5); x(3,2)=x(2,3)
84  x(3,1)=tensor(6); x(1,3)=x(3,1)
85 
86  xx= matmul( x,x )
87  i1= x(1,1)+x(2,2)+x(3,3)
88  i2= 0.5d0*( i1*i1 - (xx(1,1)+xx(2,2)+xx(3,3)) )
89  i3= x(1,1)*x(2,2)*x(3,3)+x(2,1)*x(3,2)*x(1,3)+x(3,1)*x(1,2)*x(2,3) &
90  -x(3,1)*x(2,2)*x(1,3)-x(2,1)*x(1,2)*x(3,3)-x(1,1)*x(3,2)*x(2,3)
91 
92  r=(-2.d0*i1*i1*i1+9.d0*i1*i2-27.d0*i3)/54.d0
93  q=(i1*i1-3.d0*i2)/9.d0
94 
95  ! Check discriminant to ensure all eigenvalues are real
96  ! For a 3x3 symmetric matrix, Delta = R^2 - Q^3 must be negative or zero
97  if(r*r - q*q*q > small) then
98  write(*,*) 'ERROR in tensor_eigen3: Discriminant R^2 - Q^3 > 0'
99  write(*,*) 'This indicates the algorithm has numerical issues.'
100  write(*,*) 'Discriminant =', r*r - q*q*q
101  stop 'Eigenvalue calculation failed - possible numerical error'
102  endif
103 
104  ! Add handling for special cases (when Q is close to 0)
105  if(abs(q) < small) then
106  ! Matrix has nearly equal diagonal elements (e.g. close to identity matrix)
107  eigval(1) = i1/3.d0
108  eigval(2) = i1/3.d0
109  eigval(3) = i1/3.d0
110  else
111  ! For Q^3 sufficiently small, address numerical stability
112  if(abs(q*q*q) < small) then
113  ! Q^3 is very small, but not zero - special handling
114  if(abs(r) < small) then
115  ! Both R and Q^3 small means nearly triple root
116  eigval(1) = i1/3.d0
117  eigval(2) = i1/3.d0
118  eigval(3) = i1/3.d0
119  else
120  ! When Q^3 is very small, determine theta based on sign of R only
121  ! This avoids division by a very small number which could cause instability
122  if(r >= 0.0d0) then
123  theta = 0.0d0 ! R > 0 case
124  else
125  theta = pi ! R < 0 case
126  endif
127  endif
128  else
129  ! Ensure R/sqrt(Q^3) is within [-1,1] for numerical stability
130  theta = acos(max(-1.0d0, min(1.0d0, r/dsqrt(q*q*q))))
131  endif
132 
133  eigval(1) = -2.d0*dsqrt(q)*cos(theta/3.d0)+i1/3.d0
134  eigval(2) = -2.d0*dsqrt(q)*cos((theta+2.d0*pi)/3.d0)+i1/3.d0
135  eigval(3) = -2.d0*dsqrt(q)*cos((theta-2.d0*pi)/3.d0)+i1/3.d0
136  endif
137 
138  end subroutine
139 
142  subroutine eigen3 (tensor, eigval, princ)
143  real(kind=kreal), intent(in) :: tensor(6)
144  real(kind=kreal), intent(out) :: eigval(3)
145  real(kind=kreal), intent(out) :: princ(3, 3)
146 
147  integer, parameter :: msweep = 50
148  integer :: i,j, is, ip, iq, ir
149  real(kind=kreal) :: fsum, od, theta, t, c, s, tau, g, h, hd, btens(3,3), factor
150 
151  btens(1,1)=tensor(1); btens(2,2)=tensor(2); btens(3,3)=tensor(3)
152  btens(1,2)=tensor(4); btens(2,1)=btens(1,2)
153  btens(2,3)=tensor(5); btens(3,2)=btens(2,3)
154  btens(3,1)=tensor(6); btens(1,3)=btens(3,1)
155  !
156  ! Initialise princ to the identity
157  !
158  factor = 0.d0
159  do i = 1, 3
160  do j = 1, 3
161  princ(i, j) = 0.d0
162  end do
163  princ(i, i) = 1.d0
164  eigval(i) = btens(i, i)
165  factor = factor + dabs(btens(i,i))
166  end do
167  ! Scaling and iszero/isnan exception
168  if( factor == 0.d0 .or. factor /= factor ) then
169  return
170  else
171  eigval(1:3) = eigval(1:3)/factor
172  btens(1:3,1:3) = btens(1:3,1:3)/factor
173  end if
174 
175  !
176  ! Starts sweeping.
177  !
178  do is = 1, msweep
179  fsum = 0.d0
180  do ip = 1, 2
181  do iq = ip + 1, 3
182  fsum = fsum + abs( btens(ip, iq) )
183  end do
184  end do
185  !
186  ! If the fsum of off-diagonal terms is zero returns
187  !
188  if ( fsum < 1.d-10 ) then
189  eigval(1:3) = eigval(1:3)*factor
190  return
191  endif
192  !
193  ! Performs the sweep in three rotations. One per off diagonal term
194  !
195  do ip = 1, 2
196  do iq = ip + 1, 3
197  od = 100.d0 * abs(btens(ip, iq) )
198  if ( (od+abs(eigval(ip) ) /= abs(eigval(ip) )) &
199  .and. (od+abs(eigval(iq) ) /= abs(eigval(iq) ))) then
200  hd = eigval(iq) - eigval(ip)
201  !
202  ! Evaluates the rotation angle
203  !
204  if ( abs(hd) + od == abs(hd) ) then
205  t = btens(ip, iq) / hd
206  else
207  theta = 0.5d0 * hd / btens(ip, iq)
208  t = 1.d0 / (abs(theta) + sqrt(1.d0 + theta**2) )
209  if ( theta < 0.d0 ) t = - t
210  end if
211  !
212  ! Re-evaluates the diagonal terms
213  !
214  c = 1.d0 / sqrt(1.d0 + t**2)
215  s = t * c
216  tau = s / (1.d0 + c)
217  h = t * btens(ip, iq)
218  eigval(ip) = eigval(ip) - h
219  eigval(iq) = eigval(iq) + h
220  !
221  ! Re-evaluates the remaining off-diagonal terms
222  !
223  ir = 6 - ip - iq
224  g = btens(min(ir, ip), max(ir, ip) )
225  h = btens(min(ir, iq), max(ir, iq) )
226  btens(min(ir, ip), max(ir, ip) ) = g &
227  - s * (h + g * tau)
228  btens(min(ir, iq), max(ir, iq) ) = h &
229  + s * (g - h * tau)
230  !
231  ! Rotates the eigenvectors
232  !
233  do ir = 1, 3
234  g = princ(ir, ip)
235  h = princ(ir, iq)
236  princ(ir, ip) = g - s * (h + g * tau)
237  princ(ir, iq) = h + s * (g - h * tau)
238  end do
239  end if
240  btens(ip, iq) = 0.d0
241  end do
242  end do
243  end do ! over sweeps
244  !
245  ! If convergence is not achieved stops
246  !
247  stop ' Jacobi iteration unable to converge'
248  end subroutine eigen3
249 
251  real(kind=kreal) function determinant( mat )
252  real(kind=kreal) :: mat(6)
253  real(kind=kreal) :: xj(3,3)
254 
255  xj(1,1)=mat(1); xj(2,2)=mat(2); xj(3,3)=mat(3)
256  xj(1,2)=mat(4); xj(2,1)=xj(1,2)
257  xj(2,3)=mat(5); xj(3,2)=xj(2,3)
258  xj(3,1)=mat(6); xj(1,3)=xj(3,1)
259 
260  determinant=xj(1,1)*xj(2,2)*xj(3,3) &
261  +xj(2,1)*xj(3,2)*xj(1,3) &
262  +xj(3,1)*xj(1,2)*xj(2,3) &
263  -xj(3,1)*xj(2,2)*xj(1,3) &
264  -xj(2,1)*xj(1,2)*xj(3,3) &
265  -xj(1,1)*xj(3,2)*xj(2,3)
266  end function determinant
267 
269  real(kind=kreal) function determinant33( XJ )
270  real(kind=kreal) :: xj(3,3)
271 
272  determinant33=xj(1,1)*xj(2,2)*xj(3,3) &
273  +xj(2,1)*xj(3,2)*xj(1,3) &
274  +xj(3,1)*xj(1,2)*xj(2,3) &
275  -xj(3,1)*xj(2,2)*xj(1,3) &
276  -xj(2,1)*xj(1,2)*xj(3,3) &
277  -xj(1,1)*xj(3,2)*xj(2,3)
278  end function determinant33
279 
280  subroutine fstr_chk_alloc( imsg, sub_name, ierr )
281  use hecmw
282  character(*) :: sub_name
283  integer(kind=kint) :: imsg
284  integer(kind=kint) :: ierr
285 
286  if( ierr /= 0 ) then
287  write(imsg,*) 'Memory overflow at ', sub_name
288  write(*,*) 'Memory overflow at ', sub_name
289  call hecmw_abort( hecmw_comm_get_comm( ) )
290  endif
291  end subroutine fstr_chk_alloc
292 
294  subroutine calinverse(NN, A)
295  integer, intent(in) :: NN
296  real(kind=kreal), intent(inout) :: a(nn,nn)
297 
298  integer :: I, J,K,IW,LR,IP(NN)
299  real(kind=kreal) :: w,wmax,pivot,api,eps,det
300  data eps/1.0e-35/
301  det=1.d0
302  lr = 0.0d0
303  do i=1,nn
304  ip(i)=i
305  enddo
306  do k=1,nn
307  wmax=0.d0
308  do i=k,nn
309  w=dabs(a(i,k))
310  if (w.GT.wmax) then
311  wmax=w
312  lr=i
313  endif
314  enddo
315  pivot=a(lr,k)
316  api=abs(pivot)
317  if(api.LE.eps) then
318  write(*,'(''PIVOT ERROR AT'',I5)') k
319  stop
320  end if
321  det=det*pivot
322  if (lr.NE.k) then
323  det=-det
324  iw=ip(k)
325  ip(k)=ip(lr)
326  ip(lr)=iw
327  do j=1,nn
328  w=a(k,j)
329  a(k,j)=a(lr,j)
330  a(lr,j)=w
331  enddo
332  endif
333  do i=1,nn
334  a(k,i)=a(k,i)/pivot
335  enddo
336  do i=1,nn
337  if (i.NE.k) then
338  w=a(i,k)
339  if (w.NE.0.) then
340  do j=1,nn
341  if (j.NE.k) a(i,j)=a(i,j)-w*a(k,j)
342  enddo
343  a(i,k)=-w/pivot
344  endif
345  endif
346  enddo
347  a(k,k)=1.d0/pivot
348  enddo
349 
350  do i=1,nn
351  k=ip(i)
352  if (k.NE.i) then
353  iw=ip(k)
354  ip(k)=ip(i)
355  ip(i)=iw
356  do j=1,nn
357  w=a(j,i)
358  a(j,i)=a(j,k)
359  a(j,k)=w
360  enddo
361  endif
362  enddo
363 
364  end subroutine calinverse
365 
366  subroutine cross_product(v1,v2,vn)
367  real(kind=kreal),intent(in) :: v1(3),v2(3)
368  real(kind=kreal),intent(out) :: vn(3)
369 
370  vn(1) = v1(2)*v2(3) - v1(3)*v2(2)
371  vn(2) = v1(3)*v2(1) - v1(1)*v2(3)
372  vn(3) = v1(1)*v2(2) - v1(2)*v2(1)
373  end subroutine cross_product
374 
375  subroutine transformation(jacob, tm)
376  real(kind=kreal),intent(in) :: jacob(3,3)
377  real(kind=kreal),intent(out) :: tm(6,6)
378 
379  integer :: i,j
380 
381  do i=1,3
382  do j=1,3
383  tm(i,j)= jacob(i,j)*jacob(i,j)
384  enddo
385  tm(i,4) = jacob(i,1)*jacob(i,2)
386  tm(i,5) = jacob(i,2)*jacob(i,3)
387  tm(i,6) = jacob(i,3)*jacob(i,1)
388  enddo
389  tm(4,1) = 2.d0*jacob(1,1)*jacob(2,1)
390  tm(5,1) = 2.d0*jacob(2,1)*jacob(3,1)
391  tm(6,1) = 2.d0*jacob(3,1)*jacob(1,1)
392  tm(4,2) = 2.d0*jacob(1,2)*jacob(2,2)
393  tm(5,2) = 2.d0*jacob(2,2)*jacob(3,2)
394  tm(6,2) = 2.d0*jacob(3,2)*jacob(1,2)
395  tm(4,3) = 2.d0*jacob(1,3)*jacob(2,3)
396  tm(5,3) = 2.d0*jacob(2,3)*jacob(3,3)
397  tm(6,3) = 2.d0*jacob(3,3)*jacob(1,3)
398  tm(4,4) = jacob(1,1)*jacob(2,2) + jacob(1,2)*jacob(2,1)
399  tm(5,4) = jacob(2,1)*jacob(3,2) + jacob(2,2)*jacob(3,1)
400  tm(6,4) = jacob(3,1)*jacob(1,2) + jacob(3,2)*jacob(1,1)
401  tm(4,5) = jacob(1,2)*jacob(2,3) + jacob(1,3)*jacob(2,2)
402  tm(5,5) = jacob(2,2)*jacob(3,3) + jacob(2,3)*jacob(3,2)
403  tm(6,5) = jacob(3,2)*jacob(1,3) + jacob(3,3)*jacob(1,2)
404  tm(4,6) = jacob(1,3)*jacob(2,1) + jacob(1,1)*jacob(2,3)
405  tm(5,6) = jacob(2,3)*jacob(3,1) + jacob(2,1)*jacob(3,3)
406  tm(6,6) = jacob(3,3)*jacob(1,1) + jacob(3,1)*jacob(1,3)
407 
408  end subroutine transformation
409 
410  subroutine get_principal (tensor, eigval, princmatrix)
411 
412  implicit none
413  integer i,j
414  real(kind=kreal) :: tensor(1:6)
415  real(kind=kreal) :: eigval(3)
416  real(kind=kreal) :: princmatrix(3,3)
417  real(kind=kreal) :: princnormal(3,3)
418  real(kind=kreal) :: tempv(3)
419  real(kind=kreal) :: temps
420 
421  call eigen3(tensor,eigval,princnormal)
422 
423  if (eigval(1)<eigval(2)) then
424  temps=eigval(1)
425  eigval(1)=eigval(2)
426  eigval(2)=temps
427  tempv(:)=princnormal(:,1)
428  princnormal(:,1)=princnormal(:,2)
429  princnormal(:,2)=tempv(:)
430  end if
431  if (eigval(1)<eigval(3)) then
432  temps=eigval(1)
433  eigval(1)=eigval(3)
434  eigval(3)=temps
435  tempv(:)=princnormal(:,1)
436  princnormal(:,1)=princnormal(:,3)
437  princnormal(:,3)=tempv(:)
438  end if
439  if (eigval(2)<eigval(3)) then
440  temps=eigval(2)
441  eigval(2)=eigval(3)
442  eigval(3)=temps
443  tempv(:)=princnormal(:,2)
444  princnormal(:,2)=princnormal(:,3)
445  princnormal(:,3)=tempv(:)
446  end if
447 
448  do j=1,3
449  do i=1,3
450  princmatrix(i,j) = princnormal(i,j) * eigval(j)
451  end do
452  end do
453 
454  end subroutine get_principal
455 
456  subroutine eigen3d (tensor, eigval, princ)
457  implicit none
458 
459  real(kind=kreal) :: tensor(6)
460  real(kind=kreal) :: eigval(3)
461  real(kind=kreal) :: princ(3,3)
462 
463  real(kind=kreal) :: s11, s22, s33, s12, s23, s13, j1, j2, j3
464  real(kind=kreal) :: ml,nl
465  complex(kind=kreal):: x1,x2,x3
466  real(kind=kreal):: rtemp
467  integer :: i
468  s11 = tensor(1)
469  s22 = tensor(2)
470  s33 = tensor(3)
471  s12 = tensor(4)
472  s23 = tensor(5)
473  s13 = tensor(6)
474 
475  ! invariants of stress tensor
476  j1 = s11 + s22 + s33
477  j2 = -s11*s22 - s22*s33 - s33*s11 + s12**2 + s23**2 + s13**2
478  j3 = s11*s22*s33 + 2*s12*s23*s13 - s11*s23**2 - s22*s13**2 - s33*s12**2
479  ! Cardano's method
480  ! x^3+ ax^2 + bx +c =0
481  ! s^3 - J1*s^2 -J2s -J3 =0
482  call cardano(-j1, -j2, -j3, x1, x2, x3)
483  eigval(1)= real(x1)
484  eigval(2)= real(x2)
485  eigval(3)= real(x3)
486  if (eigval(1)<eigval(2)) then
487  rtemp=eigval(1)
488  eigval(1)=eigval(2)
489  eigval(2)=rtemp
490  end if
491  if (eigval(1)<eigval(3)) then
492  rtemp=eigval(1)
493  eigval(1)=eigval(3)
494  eigval(3)=rtemp
495  end if
496  if (eigval(2)<eigval(3)) then
497  rtemp=eigval(2)
498  eigval(2)=eigval(3)
499  eigval(3)=rtemp
500  end if
501 
502  do i=1,3
503  if (eigval(i)/(eigval(1)+eigval(2)+eigval(3)) < 1.0d-10 )then
504  eigval(i) = 0.0d0
505  princ(i,:) = 0.0d0
506  exit
507  end if
508  ml = ( s23*s13 - s12*(s33-eigval(i)) ) / ( -s23**2 + (s22-eigval(i))*(s33-eigval(i)) )
509  nl = ( s12**2 - (s22-eigval(i))*(s11-eigval(i)) ) / ( s12*s23 - s13*(s22-eigval(i)) )
510  if (abs(ml) >= huge(ml)) then
511  ml=0.0d0
512  end if
513  if (abs(nl) >= huge(nl)) then
514  nl=0.0d0
515  end if
516  princ(i,1) = eigval(i)/sqrt( 1 + ml**2 + nl**2)
517  princ(i,2) = ml * princ(i,1)
518  princ(i,3) = nl * princ(i,1)
519  end do
520 
521  write(*,*)
522  end subroutine eigen3d
523 
524  subroutine cardano(a,b,c,x1,x2,x3)
525  real(kind=kreal):: a,b,c
526  real(kind=kreal):: p,q,d
527  complex(kind=kreal):: w
528  complex(kind=kreal):: u,v,y
529  complex(kind=kreal):: x1,x2,x3
530  w = (-1.0d0 + sqrt(dcmplx(-3.0d0)))/2.0d0
531  p = -a**2/9.0d0 + b/3.0d0
532  q = 2.0d0/2.7d1*a**3 - a*b/3.0d0 + c
533  d = q**2 + 4.0d0*p**3
534 
535  u = ((-dcmplx(q) + sqrt(dcmplx(d)))/2.0d0)**(1.0d0/3.0d0)
536 
537  if(u.ne.0.0d0) then
538  v = -dcmplx(p)/u
539  x1 = u + v -dcmplx(a)/3.0d0
540  x2 = u*w + v*w**2 -dcmplx(a)/3.0d0
541  x3 = u*w**2 + v*w -dcmplx(a)/3.0d0
542  else
543  y = (-dcmplx(q))**(1.0d0/3.0d0)
544  x1 = y -dcmplx(a)/3.0d0
545  x2 = y*w -dcmplx(a)/3.0d0
546  x3 = y*w**2 -dcmplx(a)/3.0d0
547  end if
548 
549  end subroutine cardano
550 
552  real(kind=kreal),intent(in) :: r(3)
553  real(kind=kreal),intent(inout) :: v(3)
554 
555  real(kind=kreal) :: rotv(3), rv
556  real(kind=kreal) :: cosx, sinc(2)
557  real(kind=kreal) :: x, x2, x4, x6
558  real(kind=kreal), parameter :: c0 = 0.5d0
559  real(kind=kreal), parameter :: c2 = -4.166666666666666d-002
560  real(kind=kreal), parameter :: c4 = 1.388888888888889d-003
561  real(kind=kreal), parameter :: c6 = -2.480158730158730d-005
562 
563  x2 = dot_product(r,r)
564  x = dsqrt(x2)
565  cosx = dcos(x)
566  if( x < 1.d-4 ) then
567  x4 = x2*x2
568  x6 = x4*x2
569  sinc(1) = 1.d0-x2/6.d0+x4/120.d0
570  sinc(2) = c0+c2*x2+c4*x4+c6*x6
571  else
572  sinc(1) = dsin(x)/x
573  sinc(2) = (1.d0-cosx)/x2
574  endif
575 
576  ! calc Rot*v
577  rv = dot_product(r,v)
578  rotv(1:3) = cosx*v(1:3)
579  rotv(1:3) = rotv(1:3)+rv*sinc(2)*r(1:3)
580  rotv(1) = rotv(1) + (-v(2)*r(3)+v(3)*r(2))*sinc(1)
581  rotv(2) = rotv(2) + (-v(3)*r(1)+v(1)*r(3))*sinc(1)
582  rotv(3) = rotv(3) + (-v(1)*r(2)+v(2)*r(1))*sinc(1)
583  v = rotv
584 
585  end subroutine
586 
588  subroutine deriv_general_iso_tensor_func_3d(dpydpx, dydx, eigv, px, py)
589  real(kind=kreal), intent(in) :: dpydpx(3,3)
590  real(kind=kreal), intent(out) :: dydx(6,6)
591  real(kind=kreal), intent(in) :: eigv(3,3)
592  real(kind=kreal), intent(in) :: px(3), py(3)
593 
594  real(kind=kreal) :: x(6)
595  real(kind=kreal) :: ep(6,3), dx2dx(6,6)
596  integer(kind=kint) :: i, j, m1, m2, m3, ia, ib, ic, ip, jp, k
597  real(kind=kreal) :: xmax, dif12, dif23, c1, c2, c3, c4, d1, d2, d3
598  real(kind=kreal) :: xa, xc, ya, yc, daa, dac, dca, dcb, dcc
599  real(kind=kreal) :: xa_xc, xa_xc2, xa_xc3, ya_yc, xc2
600  real(kind=kreal) :: s1, s2, s3, s4, s5, s6
601  real(kind=kreal), parameter :: i2(6) = [1.d0, 1.d0, 1.d0, 0.d0, 0.d0, 0.d0]
602  real(kind=kreal), parameter :: is(6,6) = reshape([&
603  1.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0,&
604  0.d0, 1.d0, 0.d0, 0.d0, 0.d0, 0.d0,&
605  0.d0, 0.d0, 1.d0, 0.d0, 0.d0, 0.d0,&
606  0.d0, 0.d0, 0.d0, 0.5d0, 0.d0, 0.d0,&
607  0.d0, 0.d0, 0.d0, 0.d0, 0.5d0, 0.d0,&
608  0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.5d0],&
609  [6, 6])
610  real(kind=kreal), parameter :: eps=1.d-6
611 
612  do i=1,3
613  ep(1,i) = eigv(1,i)*eigv(1,i)
614  ep(2,i) = eigv(2,i)*eigv(2,i)
615  ep(3,i) = eigv(3,i)*eigv(3,i)
616  ep(4,i) = eigv(1,i)*eigv(2,i)
617  ep(5,i) = eigv(2,i)*eigv(3,i)
618  ep(6,i) = eigv(3,i)*eigv(1,i)
619  enddo
620  x(:) = 0.d0
621  do i=1,3
622  do j=1,6
623  x(j) = x(j)+px(i)*ep(j,i)
624  enddo
625  enddo
626  dx2dx = reshape([2.d0*x(1), 0.d0, 0.d0, x(4), 0.d0, x(6),&
627  & 0.d0, 2.d0*x(2), 0.d0, x(4), x(5), 0.d0,&
628  & 0.d0, 0.d0, 2.d0*x(3), 0.d0, x(5), x(6),&
629  & x(4), x(4), 0.d0, (x(1)+x(2))/2.d0, x(6)/2.d0, x(5)/2.d0,&
630  & 0.d0, x(5), x(6), x(6)/2.d0, (x(2)+x(3))/2.d0, x(4)/2.d0,&
631  & x(6), 0.d0, x(6), x(5)/2.d0, x(4)/2.d0, (x(1)+x(3))/2.d0],&
632  [6, 6])
633  m1 = maxloc(px, 1)
634  m3 = minloc(px, 1)
635  if( m1 == m3 ) then
636  m1 = 1; m2 = 2; m3 = 3
637  else
638  m2 = 6 - (m1 + m3)
639  endif
640  xmax = maxval(abs(px), 1)
641  if( xmax < eps ) xmax = 1.d0
642  dif12 = abs(px(m1) - px(m2)) / xmax
643  dif23 = abs(px(m2) - px(m3)) / xmax
644  if( dif12 < eps .and. dif23 < eps ) then
645  c1 = dpydpx(1,1)-dpydpx(1,2)
646  c2 = dpydpx(1,2)
647  do j=1,6
648  do i=1,6
649  dydx(i,j) = c1*is(i,j)+c2*i2(i)*i2(j)
650  enddo
651  enddo
652  else if( dif12 < eps .or. dif23 < eps ) then
653  if( dif12 < eps ) then
654  ia = m3; ib = m1; ic = m2
655  else
656  ia = m1; ib = m2; ic = m3
657  endif
658  ya = py(ia); yc = py(ic)
659  xa = px(ia); xc = px(ic)
660  daa = dpydpx(ia,ia)
661  dac = dpydpx(ia,ic)
662  dca = dpydpx(ic,ia)
663  dcb = dpydpx(ic,ib)
664  dcc = dpydpx(ic,ic)
665  xa_xc = xa-xc
666  xa_xc2 = xa_xc*xa_xc
667  xa_xc3 = xa_xc2*xa_xc
668  ya_yc = ya-yc
669  xc2 = xc*xc
670  c1 = ya_yc/xa_xc2
671  c2 = ya_yc/xa_xc3
672  c3 = xc/xa_xc2
673  c4 = (xa+xc)/xa_xc
674  d1 = dac+dca
675  d2 = daa+dcc
676  d3 = d1-d2
677  s1 = c1+(dcb-dcc)/xa_xc
678  s2 = 2.d0*xc*c1+(dcb-dcc)*c4
679  s3 = 2.d0*c2+d3/xa_xc2
680  s4 = 2.d0*xc*c2+(dac-dcb)/xa_xc+d3*c3
681  s5 = 2.d0*xc*c2+(dca-dcb)/xa_xc+d3*c3
682  s6 = 2.d0*xc2*c2+(d1*xa-d2*xc)*c3-dcb*c4
683  do j=1,6
684  do i=1,6
685  dydx(i,j) = s1*dx2dx(i,j)-s2*is(i,j)-s3*x(i)*x(j)+s4*x(i)*i2(j)+s5*i2(i)*x(j)-s6*i2(i)*i2(j)
686  enddo
687  enddo
688  else
689  do k=1,3
690  select case(k)
691  case(1)
692  ia = 1; ib = 2; ic = 3
693  case(2)
694  ia = 2; ib = 3; ic = 1
695  case(3)
696  ia = 3; ib = 1; ic = 2
697  end select
698  c1 = py(ia)/((px(ia)-px(ib))*(px(ia)-px(ic)))
699  c2 = px(ib)+px(ic)
700  c3 = (px(ia)-px(ib))+(px(ia)-px(ic))
701  c4 = px(ib)-px(ic)
702  do j=1,6
703  do i=1,6
704  dydx(i,j) = dydx(i,j)+c1*(dx2dx(i,j)-c2*is(i,j)-c3*ep(i,ia)*ep(j,ia)-c4*(ep(i,ib)*ep(j,ib)-ep(i,ic)*ep(j,ic)))
705  enddo
706  enddo
707  enddo
708  do jp=1,3
709  do ip=1,3
710  c1 = dpydpx(ip,jp)
711  do j=1,6
712  do i=1,6
713  dydx(i,j) = dydx(i,j)+c1*ep(i,ip)*ep(j,jp)
714  enddo
715  enddo
716  enddo
717  enddo
718  endif
719  end subroutine deriv_general_iso_tensor_func_3d
720 
721 end module
Definition: hecmw.f90:6
This module provides aux functions.
Definition: utilities.f90:6
subroutine eigen3(tensor, eigval, princ)
Compute eigenvalue and eigenvetor for symmetric 3*3 tensor using Jacobi iteration adapted from numeri...
Definition: utilities.f90:143
subroutine cross_product(v1, v2, vn)
Definition: utilities.f90:367
subroutine eigen3d(tensor, eigval, princ)
Definition: utilities.f90:457
real(kind=kreal) function determinant33(XJ)
Compute determinant for 3*3 matrix.
Definition: utilities.f90:270
subroutine insert_int2array(iin, carray)
Insert an integer into a integer array.
Definition: utilities.f90:49
real(kind=kreal) function determinant(mat)
Compute determinant for symmetric 3*3 matrix.
Definition: utilities.f90:252
subroutine append_int2name(n, fname, n1)
Insert an integer at end of a file name.
Definition: utilities.f90:24
subroutine tensor_eigen3(tensor, eigval)
Given symmetric 3x3 matrix M, compute the eigenvalues.
Definition: utilities.f90:75
subroutine get_principal(tensor, eigval, princmatrix)
Definition: utilities.f90:411
subroutine transformation(jacob, tm)
Definition: utilities.f90:376
subroutine fstr_chk_alloc(imsg, sub_name, ierr)
Definition: utilities.f90:281
subroutine deriv_general_iso_tensor_func_3d(dpydpx, dydx, eigv, px, py)
Compute derivative of a general isotropic tensor function of one tensor.
Definition: utilities.f90:589
subroutine cardano(a, b, c, x1, x2, x3)
Definition: utilities.f90:525
subroutine memget(var, dimn, syze)
Record used memory.
Definition: utilities.f90:16
subroutine rotate_3dvector_by_rodrigues_formula(r, v)
Definition: utilities.f90:552
subroutine calinverse(NN, A)
calculate inverse of matrix a
Definition: utilities.f90:295