GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
dorth.f
Go to the documentation of this file.
1C Work performed under the auspices of the U.S. Department of Energy
2C by Lawrence Livermore National Laboratory under contract number
3C W-7405-Eng-48.
4C
5 SUBROUTINE dorth (VNEW, V, HES, N, LL, LDHES, KMP, SNORMW)
6C
7C***BEGIN PROLOGUE DORTH
8C***DATE WRITTEN 890101 (YYMMDD)
9C***REVISION DATE 900926 (YYMMDD)
10C
11C
12C-----------------------------------------------------------------------
13C***DESCRIPTION
14C
15C This routine orthogonalizes the vector VNEW against the previous
16C KMP vectors in the V array. It uses a modified Gram-Schmidt
17C orthogonalization procedure with conditional reorthogonalization.
18C
19C On entry
20C
21C VNEW = The vector of length N containing a scaled product
22C OF The Jacobian and the vector V(*,LL).
23C
24C V = The N x LL array containing the previous LL
25C orthogonal vectors V(*,1) to V(*,LL).
26C
27C HES = An LL x LL upper Hessenberg matrix containing,
28C in HES(I,K), K.LT.LL, scaled inner products of
29C A*V(*,K) and V(*,I).
30C
31C LDHES = The leading dimension of the HES array.
32C
33C N = The order of the matrix A, and the length of VNEW.
34C
35C LL = The current order of the matrix HES.
36C
37C KMP = The number of previous vectors the new vector VNEW
38C must be made orthogonal to (KMP .LE. MAXL).
39C
40C
41C On return
42C
43C VNEW = The new vector orthogonal to V(*,I0),
44C where I0 = MAX(1, LL-KMP+1).
45C
46C HES = Upper Hessenberg matrix with column LL filled in with
47C scaled inner products of A*V(*,LL) and V(*,I).
48C
49C SNORMW = L-2 norm of VNEW.
50C
51C-----------------------------------------------------------------------
52C***ROUTINES CALLED
53C DDOT, DNRM2, DAXPY
54C
55C***END PROLOGUE DORTH
56C
57 INTEGER N, LL, LDHES, KMP
58 DOUBLE PRECISION VNEW, V, HES, SNORMW
59 dimension vnew(*), v(n,*), hes(ldhes,*)
60 INTEGER I, I0
61 DOUBLE PRECISION ARG, DDOT, DNRM2, SUMDSQ, TEM, VNRM
62C
63C-----------------------------------------------------------------------
64C Get norm of unaltered VNEW for later use.
65C-----------------------------------------------------------------------
66 vnrm = dnrm2(n, vnew, 1)
67C-----------------------------------------------------------------------
68C Do Modified Gram-Schmidt on VNEW = A*V(LL).
69C Scaled inner products give new column of HES.
70C Projections of earlier vectors are subtracted from VNEW.
71C-----------------------------------------------------------------------
72 i0 = max0(1,ll-kmp+1)
73 DO 10 i = i0,ll
74 hes(i,ll) = ddot(n, v(1,i), 1, vnew, 1)
75 tem = -hes(i,ll)
76 CALL daxpy (n, tem, v(1,i), 1, vnew, 1)
77 10 CONTINUE
78C-----------------------------------------------------------------------
79C Compute SNORMW = norm of VNEW.
80C If VNEW is small compared to its input value (in norm), then
81C Reorthogonalize VNEW to V(*,1) through V(*,LL).
82C Correct if relative correction exceeds 1000*(unit roundoff).
83C Finally, correct SNORMW using the dot products involved.
84C-----------------------------------------------------------------------
85 snormw = dnrm2(n, vnew, 1)
86 IF (vnrm + 0.001d0*snormw .NE. vnrm) RETURN
87 sumdsq = 0.0d0
88 DO 30 i = i0,ll
89 tem = -ddot(n, v(1,i), 1, vnew, 1)
90 IF (hes(i,ll) + 0.001d0*tem .EQ. hes(i,ll)) GO TO 30
91 hes(i,ll) = hes(i,ll) - tem
92 CALL daxpy (n, tem, v(1,i), 1, vnew, 1)
93 sumdsq = sumdsq + tem**2
94 30 CONTINUE
95 IF (sumdsq .EQ. 0.0d0) RETURN
96 arg = max(0.0d0,snormw**2 - sumdsq)
97 snormw = sqrt(arg)
98 RETURN
99C
100C------END OF SUBROUTINE DORTH------------------------------------------
101 END
charNDArray max(char d, const charNDArray &m)
Definition chNDArray.cc:230
subroutine dorth(vnew, v, hes, n, ll, ldhes, kmp, snormw)
Definition dorth.f:6