NLSL
cgltri.f90
Go to the documentation of this file.
1 c NLS VERSION 20 May 1992
2 c**********************************************************************
3 c
4 c SUBROUTINE : CGLTRI
5 c -------------------
6 c
7 c This subroutine will generate the "Lanczos" tridiagonal
8 c from some of the quantities which arise naturally in the
9 c course of the conjugate gradients algorithm (in the absence
10 c of preconditioning). This is accomplished by multiplying
11 c out the factors of the Cholesky decomposition of the Lanczos
12 c tridiagonal matrix which are implicitly constructed during
13 c the conjugate gradients procedure. See pp. 370-1 in "Matrix
14 c Computations", G. Golub and C. Van Loan, Johns Hopkins Univ.
15 c Press, 1983 and pp. in " ", J. Cullum and R. Willoughby,
16 c Springer-Verlag, 1985 for more information.
17 c
18 c written by DJS 20-SEP-87
19 c
20 c Includes:
21 c nlsdim.inc
22 c rndoff.inc
23 c
24 c Uses:
25 c
26 c**********************************************************************
27 c
28  subroutine cgltri(ndone,a,b,shiftr,shifti)
29 c
30  use nlsdim
31  use rnddbl
32 c
33  integer ndone
34  double precision shiftr,shifti
35  double precision a,b
36  dimension a(2,mxstep),b(2,mxstep)
37 c
38  integer i,j
39  double precision amp,phase,tr,ti,trhalf,tihalf
40  double precision c
41  dimension c(2,mxstep)
42 c
43 c######################################################################
44 c
45  do 10 i=1,ndone
46  c(1,i)=a(1,i)
47  c(2,i)=a(2,i)
48  10 continue
49 c
50  do 20 i=2,ndone
51  j=i-1
52  tr=b(1,i)
53  ti=b(2,i)
54 c
55  amp=sqrt(sqrt(tr*tr+ti*ti))
56  if (amp.gt.rndoff) then
57  phase=0.5d0*datan2(ti,tr)
58  trhalf=amp*cos(phase)
59  tihalf=amp*sin(phase)
60  if (abs(trhalf).lt.rndoff) trhalf=0.0d0
61  if (abs(tihalf).lt.rndoff) tihalf=0.0d0
62  else
63  trhalf=0.0d0
64  tihalf=0.0d0
65  end if
66 c
67  b(1,j)=c(1,j)*trhalf-c(2,j)*tihalf
68  b(2,j)=c(1,j)*tihalf+c(2,j)*trhalf
69  if (b(1,j).lt.0.0d0) then
70  b(1,j)=-b(1,j)
71  b(2,j)=-b(2,j)
72  end if
73  a(1,i)=c(1,i)+c(1,j)*tr-c(2,j)*ti
74  a(2,i)=c(2,i)+c(1,j)*ti+c(2,j)*tr
75  20 continue
76 c
77  do 30 i=1,ndone
78  a(1,i)=a(1,i)-shiftr
79  a(2,i)=a(2,i)-shifti
80  30 continue
81 c
82  b(1,ndone)=0.0d0
83  b(2,ndone)=0.0d0
84 c
85 c----------------------------------------------------------------------
86 c return to caller
87 c----------------------------------------------------------------------
88 c
89  return
90  end
subroutine cgltri(ndone, a, b, shiftr, shifti)
Definition: cgltri.f90:29
integer, parameter mxstep
Definition: nlsdim.f90:39
double precision, parameter rndoff
Definition: rnddbl.f90:86