NLSL
tdchek.f90
Go to the documentation of this file.
1 c Version 1.4 10/10/94
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine TDCHEK
5 c =========================
6 c
7 c This routine checks whether a tridiagonal matrix is available for
8 c a given site and spectrum and allocates one if necessary.
9 c
10 c----------------------------------------------------------------------
11  subroutine tdchek(isite,ispc,ierr)
12 c
13  use nlsdim
14  use eprprm
15  use parcom
16  use tridag
17  use basis
18  use mtsdef
19  use errmsg
20 c
21  implicit none
22  integer isite,ispc,ierr
23 c
24  integer idummy,lthtd,mtstmp(mxmts)
25  character*30 ndummy
26 c
27  integer setmts
28  external setmts
29 c
30 c------------------------------------------------------------------------
31 c
32  if (ispc.le.0 .or. ispc.gt.mxspc
33  # .or. isite.le.0 .or. isite.gt.mxsite) return
34 c
35 c -------------------------------------------------------------------
36 c Check whether new tridiagonal matrix is needed for the calculation
37 c -------------------------------------------------------------------
38 c
39 c ------------------------------------------------------------
40 c SPECIAL CHECK: Avoid repeating a MOMD calculation.
41 c If a MOMD calculation is being performed for later spectra
42 c in a PSI series, check whether it is possible to use the
43 c tridiagonal matrix calculated from the first spectrum in the
44 c series (i.e. make sure the first matrix is unmodified).
45 c ------------------------------------------------------------
46  if (iparm(inort,isite).gt.1 .and. iser.eq.ipsi
47  # .and. ispc.gt.1 .and. modtd(isite,1).eq.0)
48  # then
49  ixtd(isite,ispc)=ixtd(isite,1)
50  ltd(isite,ispc)=ltd(isite,1)
51  modtd(isite,ispc)=0
52  end if
53 c
54  if (modtd(isite,ispc).ne.0) then
55 c
56 c ------------------------------------------------
57 c if nstep has been set, lthtd <- nstep
58 c else if pruned basis set in use lthtd <- ltbas
59 c else if ndim has been set lthtd <- ndim
60 c else determine ndim and set lthtd, nstep <- ndim
61 c ------------------------------------------------
62 c
63  lthtd=iparm(instep,isite)
64  if (lthtd.le.0) then
65  if (basno(isite,ispc).eq.0) then
66  lthtd=iparm(indim,isite)
67 c
68 c set NDIM if necessary
69 c
70  if (lthtd.le.0) then
71  ierr=setmts(fparm(1,isite),iparm(1,isite),mtstmp)
72  ndummy=' '
73  idummy=0
74  call lbasix(ndummy,idummy,mtstmp,lthtd,0,1,ierr)
75  iparm(indim,isite)=lthtd
76  end if
77  else
78  lthtd=ltbas(basno(isite,ispc))
79  end if
80  end if
81 c
82  if (iparm(inort,isite).gt.1) lthtd=lthtd*iparm(inort,isite)
83 c
84 c -------------------------------------------------------
85 c Check that there is room for the new tridiagonal matrix
86 c -------------------------------------------------------
87 c
88  if (nexttd+lthtd.gt.mxtdg) then
89  ierr=tdgbig
90  return
91  end if
92 c
93 c ----------------------------------
94 c Allocate a new tridiagonal matrix
95 c ----------------------------------
96 c
97  modtd(isite,ispc)=1
98  ltd(isite,ispc)=lthtd
99  ixtd(isite,ispc)=nexttd
100  nexttd=nexttd+lthtd
101  ntd=ntd+1
102  tdsite(ntd)=isite
103  tdspec(ntd)=ispc
104  end if
105  return
106  end
subroutine tdchek(isite, ispc, ierr)
Definition: tdchek.f90:12
integer, dimension(mxtdm), save tdspec
Definition: tridag.f90:28
integer, save nexttd
Definition: tridag.f90:28
integer, dimension(mxsite, mxspc), save basno
Definition: basis.f90:23
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
integer, parameter mxmts
Definition: nlsdim.f90:39
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
Definition: basis.f90:19
integer, parameter ipsi
Definition: eprprm.f90:92
integer, dimension(mxtdm), save tdsite
Definition: tridag.f90:28
integer, parameter mxspc
Definition: nlsdim.f90:39
integer, parameter indim
Definition: eprprm.f90:101
integer, dimension(mxsite, mxspc), save ltd
Definition: tridag.f90:28
integer, parameter instep
Definition: eprprm.f90:101
integer, dimension(niprm, mxsite), target, save iparm
Definition: parcom.f90:60
integer, parameter mxsite
Definition: nlsdim.f90:39
integer, parameter tdgbig
Definition: errmsg.f90:51
integer, dimension(mxsite, mxspc), save ixtd
Definition: tridag.f90:28
integer, save iser
Definition: parcom.f90:62
integer, parameter inort
Definition: eprprm.f90:101
subroutine lbasix(ixname, bss, mts, lthb, maxb, new, ierr)
Definition: lbasix.f90:39
integer, parameter mxtdg
Definition: nlsdim.f90:39
integer, save ntd
Definition: tridag.f90:28
integer, dimension(mxtdm), save ltbas
Definition: basis.f90:23