share/tensor/itensor.lisp: make X and D shared lexical variables for the functions...
[maxima.git] / share / minpack / fortran / qform.f
blob087b2478b9c5c673a13e39eb0c3d366bfca426c9
1 subroutine qform(m,n,q,ldq,wa)
2 integer m,n,ldq
3 double precision q(ldq,m),wa(m)
4 c **********
6 c subroutine qform
8 c this subroutine proceeds from the computed qr factorization of
9 c an m by n matrix a to accumulate the m by m orthogonal matrix
10 c q from its factored form.
12 c the subroutine statement is
14 c subroutine qform(m,n,q,ldq,wa)
16 c where
18 c m is a positive integer input variable set to the number
19 c of rows of a and the order of q.
21 c n is a positive integer input variable set to the number
22 c of columns of a.
24 c q is an m by m array. on input the full lower trapezoid in
25 c the first min(m,n) columns of q contains the factored form.
26 c on output q has been accumulated into a square matrix.
28 c ldq is a positive integer input variable not less than m
29 c which specifies the leading dimension of the array q.
31 c wa is a work array of length m.
33 c subprograms called
35 c fortran-supplied ... min0
37 c argonne national laboratory. minpack project. march 1980.
38 c burton s. garbow, kenneth e. hillstrom, jorge j. more
40 c **********
41 integer i,j,jm1,k,l,minmn,np1
42 double precision one,sum,temp,zero
43 data one,zero /1.0d0,0.0d0/
45 c zero out upper triangle of q in the first min(m,n) columns.
47 minmn = min0(m,n)
48 if (minmn .lt. 2) go to 30
49 do 20 j = 2, minmn
50 jm1 = j - 1
51 do 10 i = 1, jm1
52 q(i,j) = zero
53 10 continue
54 20 continue
55 30 continue
57 c initialize remaining columns to those of the identity matrix.
59 np1 = n + 1
60 if (m .lt. np1) go to 60
61 do 50 j = np1, m
62 do 40 i = 1, m
63 q(i,j) = zero
64 40 continue
65 q(j,j) = one
66 50 continue
67 60 continue
69 c accumulate q from its factored form.
71 do 120 l = 1, minmn
72 k = minmn - l + 1
73 do 70 i = k, m
74 wa(i) = q(i,k)
75 q(i,k) = zero
76 70 continue
77 q(k,k) = one
78 if (wa(k) .eq. zero) go to 110
79 do 100 j = k, m
80 sum = zero
81 do 80 i = k, m
82 sum = sum + q(i,j)*wa(i)
83 80 continue
84 temp = sum/wa(k)
85 do 90 i = k, m
86 q(i,j) = q(i,j) - temp*wa(i)
87 90 continue
88 100 continue
89 110 continue
90 120 continue
91 return
93 c last card of subroutine qform.
95 end