Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
arruda_boyce.f90
Go to the documentation of this file.
1 !*********************************************************************
2 !* Illinois Open Source License *
3 !* *
4 !* University of Illinois/NCSA *
5 !* Open Source License *
6 !* *
7 !* Copyright@2008, University of Illinois. All rights reserved. *
8 !* *
9 !* Developed by: *
10 !* *
11 !* Center for Simulation of Advanced Rockets *
12 !* *
13 !* University of Illinois *
14 !* *
15 !* www.csar.uiuc.edu *
16 !* *
17 !* Permission is hereby granted, free of charge, to any person *
18 !* obtaining a copy of this software and associated documentation *
19 !* files (the "Software"), to deal with the Software without *
20 !* restriction, including without limitation the rights to use, *
21 !* copy, modify, merge, publish, distribute, sublicense, and/or *
22 !* sell copies of the Software, and to permit persons to whom the *
23 !* Software is furnished to do so, subject to the following *
24 !* conditions: *
25 !* *
26 !* *
27 !* @ Redistributions of source code must retain the above copyright *
28 !* notice, this list of conditions and the following disclaimers. *
29 !* *
30 !* @ Redistributions in binary form must reproduce the above *
31 !* copyright notice, this list of conditions and the following *
32 !* disclaimers in the documentation and/or other materials *
33 !* provided with the distribution. *
34 !* *
35 !* @ Neither the names of the Center for Simulation of Advanced *
36 !* Rockets, the University of Illinois, nor the names of its *
37 !* contributors may be used to endorse or promote products derived *
38 !* from this Software without specific prior written permission. *
39 !* *
40 !* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
41 !* EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
42 !* OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
43 !* NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
44 !* COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
45 !* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
46 !* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
47 !* USE OR OTHER DEALINGS WITH THE SOFTWARE. *
48 !*********************************************************************
49 !* Please acknowledge The University of Illinois Center for *
50 !* Simulation of Advanced Rockets in works and publications *
51 !* resulting from this software or its derivatives. *
52 !*********************************************************************
53 SUBROUTINE arruda_boyce(Cij, &
54  s11,s22,s33,s12,s23,s13,ielem,mu,kappa)
55 
56 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90
57 !!
58 !! NAME
59 !! ARRUDA_BOYCE
60 !!
61 !! FUNCTION
62 !! Arruda-Boyce constitutive model, returns
63 !! 2nd Piola-Kircheoff Stresses
64 !!
65 !! INPUTS
66 !! Cij = right Cauchy-Green deformation tensor
67 !! mu, kappa -- material parameters
68 !! ielem -- element id number
69 !!
70 !! OUTPUT
71 !! S11,S22,S33,S12,S23,S13 -- componets of
72 !! the 2nd Piola-Kircheoff Stresses
73 !!
74 !! USES
75 !! rs, Solve_x
76 !!
77 !!
78 !!***
79 
80  IMPLICIT NONE
81 
82 !-- Variables
83  INTEGER :: j,k,l,m,ierr,ielem
84  INTEGER :: istep
85  REAL*8 :: s11,s22,s33,s12,s13,s23
86  REAL*8 :: shear_modulus
87  REAL*8 :: delta
88  REAL*8 :: cij(3,3), bulk, sqrt_n, cr, n
89  REAL*8 :: e_vec(3,3), sigma_a(3)
90  REAL*8 :: e_val(3), e_chain, xi3, xmu, xxx, xmu_max, xjay
91  REAL*8 :: fv1(3), fv2(3), stretch(3)
92  REAL*8 :: fact
93  REAL*8 :: mu, kappa
94 
95 !-- (0) initial parameters
96 
97 ! CR : (initial shear modulus)/3
98 ! Bulk : bulk modulus, kappa
99 ! N : chain locking stretch ,for example, N=8
100 
101  cr = mu/3.d0
102 
103  bulk = 1.d0*kappa
104 
105  n = 4.d0
106 
107  sqrt_n = sqrt(n)
108 
109 !
110 ! obtains the tensor b, left cauchy-green tensor using equation ???
111 !
112 !$$$ btens(1,1) = F11**2+F12**2+F13**2
113 !$$$ btens(1,2) = F21*F11+F12*F22+F13*F23
114 !$$$ btens(1,3) = F31*F11+F12*F32+F13*F33
115 !$$$ btens(2,1) = btens(1,2)
116 !$$$ btens(2,2) = F21**2+F22**2+F23**2
117 !$$$ btens(2,3) = F21*F31+F32*F22+F23*F33
118 !$$$ btens(3,1) = btens(1,3)
119 !$$$ btens(3,2) = btens(2,3)
120 !$$$ btens(3,3) = F31**2+F32**2+F33**2
121 !$$$
122 !$$$ CALL jacobi(btens,stret,princ)
123 
124 
125 
126 !-- (3) compute eigen values (e_val) and eigen vectors (e_vec) of Cij
127 
128  DO j=1,3
129  stretch(j)=0.0d0
130  e_val(j)=0.0d0
131  fv1(j)=0.0d0
132  fv2(j)=0.0d0
133  DO k=1,3
134  e_vec(j,k)=0.0d0
135  END DO
136  END DO
137 !
138  CALL rs(3,3,cij,e_val,1,e_vec,fv1,fv2,ierr)
139  IF(ierr .NE. 0) THEN
140  print *,' error occurs at element ', ielem
141  END IF
142 
143 !-- (4) calculate stretch
144 
145  DO j=1,3
146  stretch(j)=sqrt(e_val(j))
147  END DO
148 
149 !-- (5) compute ramda chain
150 
151  e_chain=sqrt(stretch(1)**2+stretch(2)**2+stretch(3)**2)
152  e_chain=1.0d0/sqrt(3.0d0)*(e_chain)
153 
154 !$$$ e_vec(:,:) = princ(:,:)
155 !$$$
156 !$$$ stretch(1:3) = stret(1:3)
157 
158 
159 !-- (6) compute I3 and Jay
160 
161  xjay= stretch(1)*stretch(2)*stretch(3)
162  xi3 = (xjay)**2
163 
164 !-- (7) compute xmu
165 
166  xmu = e_chain/sqrt_n
167 
168 !-- (8) Solve for x
169 
170  CALL solve_x(xxx,xmu,ielem,stretch)
171 
172 !-- (9) Compute sigma_A
173 
174  DO k=1,3
175  sigma_a(k)=cr*sqrt_n*(stretch(k)**2-e_chain**2)/ &
176  e_chain*xxx + bulk*log(sqrt(xi3))
177  END DO
178 
179 !-- (10) Compute 2nd Piola-Kircheoff Stresses
180 
181  s11 = 0.0d0
182  s22 = 0.0d0
183  s33 = 0.0d0
184  s12 = 0.0d0
185  s23 = 0.0d0
186  s13 = 0.0d0
187 
188  DO k = 1, 3
189  fact = xjay * sigma_a(k) / (stretch(k)**2)
190 
191  s11 = s11 + e_vec(1,k) * e_vec(1,k) * fact
192  s22 = s22 + e_vec(2,k) * e_vec(2,k) * fact
193  s33 = s33 + e_vec(3,k) * e_vec(3,k) * fact
194  s12 = s12 + e_vec(1,k) * e_vec(2,k) * fact
195  s23 = s23 + e_vec(2,k) * e_vec(3,k) * fact
196  s13 = s13 + e_vec(1,k) * e_vec(3,k) * fact
197  END DO
198 
199  RETURN
200 END SUBROUTINE arruda_boyce
201 
202 
203 SUBROUTINE solve_x(x,xmu,i,stretch)
204 
205 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/Solve_x
206 !!
207 !! NAME
208 !! Solve_x
209 !!
210 !! FUNCTION
211 !! Solve for the x
212 !!
213 !! INPUTS
214 !! xmu -- mu
215 !! i -- element id
216 !! stretch -- principle stretch
217 !!
218 !! OUTPUT
219 !! x -- x
220 !!
221 !!***
222 
223  IMPLICIT NONE
224 
225 
226  REAL*8 x, xmu, f, fp, dx, x0
227  REAL*8 sinh, cosh, coth
228  REAL*8 stretch(3)
229  INTEGER iter, i, k
230 
231  x0 = 0.1d0
232  x = x0
233  iter = 1
234 
235 10 sinh = (exp(x) - exp(-x))*0.5d0
236  cosh = (exp(x) + exp(-x))*0.5d0
237  coth = cosh/sinh
238  f = coth -1.d0/x - xmu
239  IF(abs(f) .LE. 1.0d-10) go to 20
240  fp = -1.d0/(sinh**2) + 1.d0/(x**2)
241  dx = -f/fp
242  x = x + dx
243 
244  IF(iter .GT. 10000) THEN
245  print *,' iteration exceeds 10000'
246  print *,' program stop ! '
247  print *,' at element ',i
248  WRITE(*,*) (stretch(k),k=1,3)
249  stop
250  END IF
251 
252  iter=iter+1
253  go to 10
254 
255 20 CONTINUE
256  RETURN
257 END SUBROUTINE solve_x
258 
259 SUBROUTINE rs(nm,n,a,w,matz,z,fv1,fv2,ierr)
260 
261  IMPLICIT NONE
262 
263  INTEGER n,nm,ierr,matz
264  REAL*8 a(nm,n),w(n),z(nm,n),fv1(n),fv2(n)
265 
266 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/rs
267 !!
268 !! NAME
269 !! rs
270 !!
271 !! FUNCTION
272 !! calls the recommended SEQUENCE of
273 !! subroutines from the eigensystem SUBROUTINE package (eispack)
274 !! to find the eigenvalues and eigenvectors (IF desired)
275 !! of a REAL symmetric matrix.
276 !!
277 !!
278 !! INPUTS
279 !!
280 !! nm must be set to the row DIMENSION of the two-dimensional
281 !! array parameters as declared in the calling PROGRAM
282 !! DIMENSION statement.
283 !!
284 !! n is the order of the matrix a.
285 !!
286 !! a CONTAINS the REAL symmetric matrix.
287 !!
288 !! matz is an INTEGER variable set equal to zero IF
289 !! ONLY eigenvalues are desired. otherwise it is set to
290 !! any non-zero INTEGER for both eigenvalues and eigenvectors.
291 !!
292 !! OUTPUT
293 !!
294 !! w CONTAINS the eigenvalues in ascending order.
295 !!
296 !! z CONTAINS the eigenvectors IF matz is not zero.
297 !!
298 !! ierr is an INTEGER output variable set equal to an error
299 !! completion code described in the documentation for tqlrat
300 !! and tql2. the normal completion code is zero.
301 !!
302 !! fv1 and fv2 are temporary storage arrays.
303 !!
304 !! NOTES
305 !!
306 !! questions and comments should be directed to burton s. garbow,
307 !! mathematics and computer science div, argonne national laboratory
308 !!
309 !! this version dated august 1983.
310 !!
311 !!***
312 
313  IF (n .LE. nm) go to 10
314  ierr = 10 * n
315  go to 50
316 !
317 10 IF (matz .NE. 0) go to 20
318 ! .......... find eigenvalues ONLY ..........
319  CALL tred1(nm,n,a,w,fv1,fv2)
320 ! tqlrat encounters catastrophic underflow on the Vax
321 ! CALL tqlrat(n,w,fv2,ierr)
322 CALL tql1(n,w,fv1,ierr)
323 go to 50
324 ! .......... find both eigenvalues and eigenvectors ..........
325 20 CALL tred2(nm,n,a,w,fv1,z)
326 CALL tql2(nm,n,w,fv1,z,ierr)
327 50 RETURN
328 END SUBROUTINE rs
329 
330 REAL*8 FUNCTION pythag(a,b)
331 
332 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/pythag
333 !!
334 !! NAME
335 !! pythag
336 !!
337 !! FUNCTION
338 !! finds SQRT(a**2+b**2) without overflow or destructive underflow
339 !!*****
340 
341  IMPLICIT NONE
342  REAL*8 a,b
343 
344  REAL*8 p,r,s,t,u
345 ! DOUBLE PRECISION p,r,s,t,u
346  p = dmax1(abs(a),abs(b))
347  IF (p .EQ. 0.0d0) go to 20
348  r = (dmin1(abs(a),abs(b))/p)**2
349 10 CONTINUE
350  t = 4.0d0 + r
351  IF (t .EQ. 4.0d0) go to 20
352  s = r/t
353  u = 1.0d0 + 2.0d0*s
354  p = u*p
355  r = (s/u)**2 * r
356  go to 10
357 20 pythag = p
358  RETURN
359 END FUNCTION pythag
360 
361 SUBROUTINE tql1(n,d,e,ierr)
362 
363 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/tql1
364 !!
365 !! NAME
366 !! tql1
367 !!
368 !! FUNCTION
369 !! finds the eigenvalues of a symmetric
370 !! tridiagonal matrix by the ql method.
371 !!
372 !!
373 !! INPUTS
374 !!
375 !! n is the order of the matrix.
376 !!
377 !! d CONTAINS the diagonal elements of the input matrix.
378 !!
379 !! e CONTAINS the subdiagonal elements of the input matrix
380 !! in its last n-1 positions. e(1) is arbitrary.
381 !!
382 !! OUTPUT
383 !!
384 !! d CONTAINS the eigenvalues in ascending order. IF an
385 !! error EXIT is made, the eigenvalues are correct and
386 !! ordered for indices 1,2,...ierr-1, but may not be
387 !! the smallest eigenvalues.
388 !!
389 !! e has been destroyed.
390 !!
391 !! ierr is set to
392 !! zero for normal RETURN,
393 !! j IF the j-th eigenvalue has not been
394 !! determined after 30 iterations.
395 !! USES
396 !!
397 !! pythag for SQRT(a*a + b*b) .
398 !!
399 !! NOTES
400 !!
401 !! this SUBROUTINE is a translation of the algol PROCEDURE tql1,
402 !! num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
403 !! wilkinson.
404 !! handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
405 !!
406 !!
407 !! questions and comments should be directed to burton s. garbow,
408 !! mathematics and computer science div, argonne national laboratory
409 !!
410 !! this version dated august 1983.
411 !!
412 !!*****
413  IMPLICIT NONE
414 
415  INTEGER i,j,l,m,n,ii,l1,l2,mml,ierr
416  REAL*8 d(n),e(n)
417  REAL*8 c,c2,c3,dl1,el1,f,g,h,p,r,s,s2,tst1,tst2,pythag
418 
419 
420 
421 ! ------------------------------------------------------------------
422 !
423  ierr = 0
424  IF (n .EQ. 1) go to 1001
425 
426  DO 100 i = 2, n
427 100 e(i-1) = e(i)
428 
429  f = 0.0d0
430  tst1 = 0.0d0
431  e(n) = 0.0d0
432 
433  DO 290 l = 1, n
434  j = 0
435  h = abs(d(l)) + abs(e(l))
436  IF (tst1 .LT. h) tst1 = h
437 ! .......... look for small sub-diagonal element ..........
438  DO 110 m = l, n
439  tst2 = tst1 + abs(e(m))
440  IF (tst2 .EQ. tst1) go to 120
441 ! .......... e(n) is always zero, so there is no EXIT
442 ! through the bottom of the loop ..........
443 110 CONTINUE
444 !
445 120 IF (m .EQ. l) go to 210
446 130 IF (j .EQ. 30) go to 1000
447  j = j + 1
448 ! .......... form shift ..........
449  l1 = l + 1
450  l2 = l1 + 1
451  g = d(l)
452  p = (d(l1) - g) / (2.0d0 * e(l))
453  r = pythag(p,1.0d0)
454  d(l) = e(l) / (p + sign(r,p))
455  d(l1) = e(l) * (p + sign(r,p))
456  dl1 = d(l1)
457  h = g - d(l)
458  IF (l2 .GT. n) go to 145
459 !
460  DO 140 i = l2, n
461  140 d(i) = d(i) - h
462 !
463  145 f = f + h
464 ! .......... ql transformation ..........
465  p = d(m)
466  c = 1.0d0
467  c2 = c
468  el1 = e(l1)
469  s = 0.0d0
470  mml = m - l
471 ! .......... for i=m-1 step -1 until l DO -- ..........
472  DO 200 ii = 1, mml
473  c3 = c2
474  c2 = c
475  s2 = s
476  i = m - ii
477  g = c * e(i)
478  h = c * p
479  r = pythag(p,e(i))
480  e(i+1) = s * r
481  s = e(i) / r
482  c = p / r
483  p = c * d(i) - s * g
484  d(i+1) = h + s * (c * g + s * d(i))
485  200 CONTINUE
486 !
487  p = -s * s2 * c3 * el1 * e(l) / dl1
488  e(l) = s * p
489  d(l) = c * p
490  tst2 = tst1 + abs(e(l))
491  IF (tst2 .GT. tst1) go to 130
492  210 p = d(l) + f
493 ! .......... order eigenvalues ..........
494  IF (l .EQ. 1) go to 250
495 ! .......... for i=l step -1 until 2 DO -- ..........
496  DO 230 ii = 2, l
497  i = l + 2 - ii
498  IF (p .GE. d(i-1)) go to 270
499  d(i) = d(i-1)
500  230 CONTINUE
501 !
502  250 i = 1
503  270 d(i) = p
504  290 CONTINUE
505 !
506  go to 1001
507 ! .......... set error -- no convergence to an
508 ! eigenvalue after 30 iterations ..........
509  1000 ierr = l
510  1001 RETURN
511  END
512  SUBROUTINE tql2(nm,n,d,e,z,ierr)
513 
514 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/tql2
515 !!
516 !! NAME
517 !! tql2
518 !!
519 !! FUNCTION
520 !!
521 !! Finds the eigenvalues and eigenvectors
522 !! of a symmetric tridiagonal matrix by the ql method.
523 !! the eigenvectors of a full symmetric matrix can also
524 !! be found IF tred2 has been used to reduce this
525 !! full matrix to tridiagonal form.
526 
527 !! NOTES
528 !! this SUBROUTINE is a translation of the algol PROCEDURE tql2,
529 !! num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
530 !! wilkinson.
531 !! handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
532 !!
533 !!
534 !! questions and comments should be directed to burton s. garbow,
535 !! mathematics and computer science div, argonne national laboratory
536 !!
537 !! this version dated august 1983.
538 !!
539 !!
540 !! INPUTS
541 !!
542 !! nm must be set to the row DIMENSION of two-dimensional
543 !! array parameters as declared in the calling PROGRAM
544 !! DIMENSION statement.
545 !!
546 !! n is the order of the matrix.
547 !!
548 !! d CONTAINS the diagonal elements of the input matrix.
549 !!
550 !! e CONTAINS the subdiagonal elements of the input matrix
551 !! in its last n-1 positions. e(1) is arbitrary.
552 !!
553 !! z CONTAINS the transformation matrix produced in the
554 !! reduction by tred2, IF performed. IF the eigenvectors
555 !! of the tridiagonal matrix are desired, z must contain
556 !! the identity matrix.
557 !!
558 !! OUTPUT
559 !!
560 !! d CONTAINS the eigenvalues in ascending order. IF an
561 !! error EXIT is made, the eigenvalues are correct but
562 !! unordered for indices 1,2,...,ierr-1.
563 !!
564 !! e has been destroyed.
565 !!
566 !! z CONTAINS orthonormal eigenvectors of the symmetric
567 !! tridiagonal (or full) matrix. IF an error EXIT is made,
568 !! z CONTAINS the eigenvectors associated WITH the stored
569 !! eigenvalues.
570 !!
571 !! ierr is set to
572 !! zero for normal RETURN,
573 !! j IF the j-th eigenvalue has not been
574 !! determined after 30 iterations.
575 !!
576 !! USES
577 !! pythag for SQRT(a*a + b*b) .
578 !!
579 !!***
580 
581  IMPLICIT NONE
582 !
583  INTEGER i,j,k,l,m,n,ii,l1,l2,nm,mml,ierr
584  REAL*8 d(n),e(n),z(nm,n)
585  REAL*8 c,c2,c3,dl1,el1,f,g,h,p,r,s,s2,tst1,tst2,pythag
586 ! DOUBLE PRECISION d(n),e(n),z(nm,n)
587 ! DOUBLE PRECISION c,c2,c3,dl1,el1,f,g,h,p,r,s,s2,tst1,tst2,pythag
588 !
589  ierr = 0
590  IF (n .EQ. 1) go to 1001
591 
592  DO 100 i = 2, n
593  100 e(i-1) = e(i)
594 
595  f = 0.0d0
596  tst1 = 0.0d0
597  e(n) = 0.0d0
598 
599  DO 240 l = 1, n
600  j = 0
601  h = abs(d(l)) + abs(e(l))
602  IF (tst1 .LT. h) tst1 = h
603 ! .......... look for small sub-diagonal element ..........
604  DO 110 m = l, n
605  tst2 = tst1 + abs(e(m))
606  IF (tst2 .EQ. tst1) go to 120
607 ! .......... e(n) is always zero, so there is no EXIT
608 ! through the bottom of the loop ..........
609  110 CONTINUE
610 !
611  120 IF (m .EQ. l) go to 220
612  130 IF (j .EQ. 30) go to 1000
613  j = j + 1
614 ! .......... form shift ..........
615  l1 = l + 1
616  l2 = l1 + 1
617  g = d(l)
618  p = (d(l1) - g) / (2.0d0 * e(l))
619  r = pythag(p,1.0d0)
620  d(l) = e(l) / (p + sign(r,p))
621  d(l1) = e(l) * (p + sign(r,p))
622  dl1 = d(l1)
623  h = g - d(l)
624  IF (l2 .GT. n) go to 145
625 !
626  DO 140 i = l2, n
627  140 d(i) = d(i) - h
628 !
629  145 f = f + h
630 ! .......... ql transformation ..........
631  p = d(m)
632  c = 1.0d0
633  c2 = c
634  el1 = e(l1)
635  s = 0.0d0
636  mml = m - l
637 ! .......... for i=m-1 step -1 until l DO -- ..........
638  DO 200 ii = 1, mml
639  c3 = c2
640  c2 = c
641  s2 = s
642  i = m - ii
643  g = c * e(i)
644  h = c * p
645  r = pythag(p,e(i))
646  e(i+1) = s * r
647  s = e(i) / r
648  c = p / r
649  p = c * d(i) - s * g
650  d(i+1) = h + s * (c * g + s * d(i))
651 ! .......... form vector ..........
652  DO 180 k = 1, n
653  h = z(k,i+1)
654  z(k,i+1) = s * z(k,i) + c * h
655  z(k,i) = c * z(k,i) - s * h
656  180 CONTINUE
657 !
658  200 CONTINUE
659 !
660  p = -s * s2 * c3 * el1 * e(l) / dl1
661  e(l) = s * p
662  d(l) = c * p
663  tst2 = tst1 + abs(e(l))
664  IF (tst2 .GT. tst1) go to 130
665  220 d(l) = d(l) + f
666  240 CONTINUE
667 ! .......... order eigenvalues and eigenvectors ..........
668  DO 300 ii = 2, n
669  i = ii - 1
670  k = i
671  p = d(i)
672 !
673  DO 260 j = ii, n
674  IF (d(j) .GE. p) go to 260
675  k = j
676  p = d(j)
677  260 CONTINUE
678 !
679  IF (k .EQ. i) go to 300
680  d(k) = d(i)
681  d(i) = p
682 !
683  DO 280 j = 1, n
684  p = z(j,i)
685  z(j,i) = z(j,k)
686  z(j,k) = p
687  280 CONTINUE
688 !
689  300 CONTINUE
690 !
691  go to 1001
692 ! .......... set error -- no convergence to an
693 ! eigenvalue after 30 iterations ..........
694  1000 ierr = l
695  1001 RETURN
696  END
697  SUBROUTINE tred1(nm,n,a,d,e,e2)
698 
699  IMPLICIT NONE
700 !
701  INTEGER i,j,k,l,n,ii,nm,jp1
702  REAL*8 a(nm,n),d(n),e(n),e2(n)
703  REAL*8 f,g,h,scale
704 
705 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/tred1
706 !!
707 !! NAME
708 !! tred1
709 !!
710 !! FUNCTION
711 !! Reduces a REAL symmetric matrix
712 !! to a symmetric tridiagonal matrix using
713 !! orthogonal similarity transformations.
714 !!
715 !! INPUTS
716 !!
717 !! nm must be set to the row DIMENSION of two-dimensional
718 !! array parameters as declared in the calling PROGRAM
719 !! DIMENSION statement.
720 !!
721 !! n is the order of the matrix.
722 !!
723 !! a CONTAINS the REAL symmetric input matrix. ONLY the
724 !! lower triangle of the matrix need be supplied.
725 !!
726 !! OUTPUT
727 !!
728 !! a CONTAINS information about the orthogonal trans-
729 !! formations used in the reduction in its strict lower
730 !! triangle. the full upper triangle of a is unaltered.
731 !!
732 !! d CONTAINS the diagonal elements of the tridiagonal matrix.
733 !!
734 !! e CONTAINS the subdiagonal elements of the tridiagonal
735 !! matrix in its last n-1 positions. e(1) is set to zero.
736 !!
737 !! e2 CONTAINS the squares of the corresponding elements of e.
738 !! e2 may coincide WITH e IF the squares are not needed.
739 !!
740 !! NOTES
741 !!
742 !! this SUBROUTINE is a translation of the algol PROCEDURE tred1,
743 !! num. math. 11, 181-195(1968) by martin, reinsch, and wilkinson.
744 !! handbook for auto. comp., vol.ii-linear algebra, 212-226(1971).
745 !!
746 !! questions and comments should be directed to burton s. garbow,
747 !! mathematics and computer science div, argonne national laboratory
748 !!
749 !! this version dated august 1983.
750 !!
751 !!****
752 
753  DO 100 i = 1, n
754  d(i) = a(n,i)
755  a(n,i) = a(i,i)
756  100 CONTINUE
757 ! .......... for i=n step -1 until 1 DO -- ..........
758  DO 300 ii = 1, n
759  i = n + 1 - ii
760  l = i - 1
761  h = 0.0d0
762  scale = 0.0d0
763  IF (l .LT. 1) go to 130
764 ! .......... scale row (algol tol THEN not needed) ..........
765  DO 120 k = 1, l
766  120 scale = scale + abs(d(k))
767 
768  IF (scale .NE. 0.0d0) go to 140
769 
770  DO 125 j = 1, l
771  d(j) = a(l,j)
772  a(l,j) = a(i,j)
773  a(i,j) = 0.0d0
774  125 CONTINUE
775 
776  130 e(i) = 0.0d0
777  e2(i) = 0.0d0
778  go to 300
779 
780  140 DO 150 k = 1, l
781  d(k) = d(k) / scale
782  h = h + d(k) * d(k)
783  150 CONTINUE
784 
785  e2(i) = scale * scale * h
786  f = d(l)
787  g = -sign(sqrt(h),f)
788  e(i) = scale * g
789  h = h - f * g
790  d(l) = f - g
791  IF (l .EQ. 1) go to 285
792 ! .......... form a*u ..........
793  DO 170 j = 1, l
794  170 e(j) = 0.0d0
795 
796  DO 240 j = 1, l
797  f = d(j)
798  g = e(j) + a(j,j) * f
799  jp1 = j + 1
800  IF (l .LT. jp1) go to 220
801 
802  DO 200 k = jp1, l
803  g = g + a(k,j) * d(k)
804  e(k) = e(k) + a(k,j) * f
805  200 CONTINUE
806 
807  220 e(j) = g
808  240 CONTINUE
809 ! .......... form p ..........
810  f = 0.0d0
811 
812  DO 245 j = 1, l
813  e(j) = e(j) / h
814  f = f + e(j) * d(j)
815  245 CONTINUE
816 
817  h = f / (h + h)
818 ! .......... form q ..........
819  DO 250 j = 1, l
820  250 e(j) = e(j) - h * d(j)
821 ! .......... form reduced a ..........
822  DO 280 j = 1, l
823  f = d(j)
824  g = e(j)
825 
826  DO 260 k = j, l
827  260 a(k,j) = a(k,j) - f * e(k) - g * d(k)
828 
829  280 CONTINUE
830 
831  285 DO 290 j = 1, l
832  f = d(j)
833  d(j) = a(l,j)
834  a(l,j) = a(i,j)
835  a(i,j) = f * scale
836  290 CONTINUE
837 
838  300 CONTINUE
839 
840  RETURN
841  END
842 
843  SUBROUTINE tred2(nm,n,a,d,e,z)
844 
845 !!****f* Rocfrac/Rocfrac/Source/arruda_boyce.f90/tred2
846 !!
847 !! NAME
848 !! tred2
849 !!
850 !! FUNCTION
851 !! reduces a REAL symmetric matrix to a
852 !! symmetric tridiagonal matrix using and accumulating
853 !! orthogonal similarity transformations.
854 !!
855 !! INPUTS
856 !!
857 !! nm must be set to the row DIMENSION of two-dimensional
858 !! array parameters as declared in the calling PROGRAM
859 !! DIMENSION statement.
860 !!
861 !! n is the order of the matrix.
862 !!
863 !! a CONTAINS the REAL symmetric input matrix. ONLY the
864 !! lower triangle of the matrix need be supplied.
865 !!
866 !! OUTPUT
867 !!
868 !! d CONTAINS the diagonal elements of the tridiagonal matrix.
869 !!
870 !! e CONTAINS the subdiagonal elements of the tridiagonal
871 !! matrix in its last n-1 positions. e(1) is set to zero.
872 !!
873 !! z CONTAINS the orthogonal transformation matrix
874 !! produced in the reduction.
875 !!
876 !! a and z may coincide. IF distinct, a is unaltered.
877 !!
878 !! NOTES
879 !! this SUBROUTINE is a translation of the algol PROCEDURE tred2,
880 !! num. math. 11, 181-195(1968) by martin, reinsch, and wilkinson.
881 !! handbook for auto. comp., vol.ii-linear algebra, 212-226(1971).
882 !! questions and comments should be directed to burton s. garbow,
883 !! mathematics and computer science div, argonne national laboratory
884 !!
885 !! this version dated august 1983.
886 !!****
887 
888 
889  IMPLICIT NONE
890 
891  INTEGER i,j,k,l,n,ii,nm,jp1
892  REAL*8 a(nm,n),d(n),e(n),z(nm,n)
893  REAL*8 f,g,h,hh,scale
894 
895  DO 100 i = 1, n
896 
897  DO 80 j = i, n
898  80 z(j,i) = a(j,i)
899 
900  d(i) = a(n,i)
901  100 CONTINUE
902 
903  IF (n .EQ. 1) go to 510
904 ! .......... for i=n step -1 until 2 DO -- ..........
905  DO 300 ii = 2, n
906  i = n + 2 - ii
907  l = i - 1
908  h = 0.0d0
909  scale = 0.0d0
910  IF (l .LT. 2) go to 130
911 ! .......... scale row (algol tol THEN not needed) ..........
912  DO 120 k = 1, l
913  120 scale = scale + abs(d(k))
914 
915  IF (scale .NE. 0.0d0) go to 140
916  130 e(i) = d(l)
917 
918  DO 135 j = 1, l
919  d(j) = z(l,j)
920  z(i,j) = 0.0d0
921  z(j,i) = 0.0d0
922  135 CONTINUE
923 
924  go to 290
925 
926  140 DO 150 k = 1, l
927  d(k) = d(k) / scale
928  h = h + d(k) * d(k)
929  150 CONTINUE
930 
931  f = d(l)
932  g = -sign(sqrt(h),f)
933  e(i) = scale * g
934  h = h - f * g
935  d(l) = f - g
936 ! .......... form a*u ..........
937  DO 170 j = 1, l
938  170 e(j) = 0.0d0
939 
940  DO 240 j = 1, l
941  f = d(j)
942  z(j,i) = f
943  g = e(j) + z(j,j) * f
944  jp1 = j + 1
945  IF (l .LT. jp1) go to 220
946 
947  DO 200 k = jp1, l
948  g = g + z(k,j) * d(k)
949  e(k) = e(k) + z(k,j) * f
950  200 CONTINUE
951 
952  220 e(j) = g
953  240 CONTINUE
954 ! .......... form p ..........
955  f = 0.0d0
956 
957  DO 245 j = 1, l
958  e(j) = e(j) / h
959  f = f + e(j) * d(j)
960  245 CONTINUE
961 
962  hh = f / (h + h)
963 ! .......... form q ..........
964  DO 250 j = 1, l
965  250 e(j) = e(j) - hh * d(j)
966 ! .......... form reduced a ..........
967  DO 280 j = 1, l
968  f = d(j)
969  g = e(j)
970 
971  DO 260 k = j, l
972  260 z(k,j) = z(k,j) - f * e(k) - g * d(k)
973 
974  d(j) = z(l,j)
975  z(i,j) = 0.0d0
976  280 CONTINUE
977 
978  290 d(i) = h
979  300 CONTINUE
980 ! .......... accumulation of transformation matrices ..........
981  DO 500 i = 2, n
982  l = i - 1
983  z(n,l) = z(l,l)
984  z(l,l) = 1.0d0
985  h = d(i)
986  IF (h .EQ. 0.0d0) go to 380
987 
988  DO 330 k = 1, l
989  330 d(k) = z(k,i) / h
990 
991  DO 360 j = 1, l
992  g = 0.0d0
993 
994  DO 340 k = 1, l
995  340 g = g + z(k,i) * z(k,j)
996 
997  DO 360 k = 1, l
998  z(k,j) = z(k,j) - g * d(k)
999  360 CONTINUE
1000 
1001  380 DO 400 k = 1, l
1002  400 z(k,i) = 0.0d0
1003 
1004  500 CONTINUE
1005 
1006  510 DO 520 i = 1, n
1007  d(i) = z(n,i)
1008  z(n,i) = 0.0d0
1009  520 CONTINUE
1010 
1011  z(n,n) = 1.0d0
1012  e(1) = 0.0d0
1013  RETURN
1014  END
1015 
subroutine tql1(n, d, e, ierr)
unsigned char r() const
Definition: Color.h:68
FT m(int i, int j) const
subroutine rs(nm, n, a, w, matz, z, fv1, fv2, ierr)
static SURF_BEGIN_NAMESPACE double sign(double x)
const NT & d
j indices k indices k
Definition: Indexing.h:6
NT dx
double s
Definition: blastest.C:80
subroutine tred2(nm, n, a, d, e, z)
unsigned char b() const
Definition: Color.h:70
subroutine tql2(nm, n, d, e, z, ierr)
double sqrt(double d)
Definition: double.h:73
RT c() const
Definition: Line_2.h:150
subroutine tred1(nm, n, a, d, e, e2)
void scale(const Real &a, Nodal_data &x)
void int int int REAL REAL REAL * z
Definition: write.cpp:76
blockLoc i
Definition: read.cpp:79
void int int REAL * x
Definition: read.cpp:74
const NT & n
RT delta(int i) const
Definition: Direction_2.h:159
virtual std::ostream & print(std::ostream &os) const
j indices j
Definition: Indexing.h:6
REAL *8 function pythag(a, b)
subroutine arruda_boyce(Cij, S11, S22, S33, S12, S23, S13, ielem, mu, kappa)
subroutine solve_x(x, xmu, i, stretch)
RT a() const
Definition: Line_2.h:140
unsigned char g() const
Definition: Color.h:69