NLSL
tdsqz.f90
Go to the documentation of this file.
1 c NLSL Version 1.3.2 2/22/94
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine TDSQZ
5 c =========================
6 c
7 c Loops through blocks stored in common /tridag/ and moves tridiagonal
8 c matrices that may be re-used to the bottom of the storage space
9 c
10 c----------------------------------------------------------------------
11 c
12  subroutine tdsqz()
13 c
14  use nlsdim
15  use expdat
16  use tridag
17  use stdio
18 c
19  implicit none
20  integer i,j,next,isp,isi,ixt,k,n
21 c
22  next=1
23  n=0
24  do i=1,ntd
25  isi=tdsite(i)
26  isp=tdspec(i)
27 c
28 c Check whether next tridiagonal matrix may be kept
29 c
30  if (modtd(isi,isp).eq.0) then
31 c
32 c Move the block to lowest available storage location
33 c
34  ixt=ixtd(isi,isp)
35  if (ixt.gt.next) then
36  ixtd(isi,isp)=next
37  do j=1,ltd(isi,isp)
38  k=ixt+j-1
39  alpha(next)=alpha(k)
40  beta(next)=beta(k)
41  next=next+1
42  end do
43  else if (ixt.eq.next) then
44  next=next+ltd(isi,isp)
45  else
46  write (luout,1005)
47  end if
48  n=n+1
49  tdsite(n)=isi
50  tdspec(n)=isp
51  end if
52 c
53  end do
54  nexttd=next
55  ntd=n
56 c
57  1005 format('*** Error in tridiagonal matrix storage ***')
58 c
59  return
60  end
61 
62 
63 
64 
65 
66 
67 
68 
69 
70 
integer, save luout
Definition: stdio.f90:32
double complex, dimension(mxtdg), save alpha
Definition: tridag.f90:32
integer, dimension(mxtdm), save tdspec
Definition: tridag.f90:28
integer, save nexttd
Definition: tridag.f90:28
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
Definition: stdio.f90:26
subroutine tdsqz()
Definition: tdsqz.f90:13
integer, dimension(mxtdm), save tdsite
Definition: tridag.f90:28
integer, dimension(mxsite, mxspc), save ltd
Definition: tridag.f90:28
integer, dimension(mxsite, mxspc), save ixtd
Definition: tridag.f90:28
double complex, dimension(mxtdg), save beta
Definition: tridag.f90:32
integer, save ntd
Definition: tridag.f90:28