NLSL
writr.f90
Go to the documentation of this file.
1 c Version 1.5 beta 11/23/95
2 c----------------------------------------------------------------------
3 c =========================
4 c WRITR subroutines
5 c =========================
6 c
7 c Outputs various information regarding intermediate spectral
8 c calculations for a specific iteration during a nonlinear least-squares
9 c fit.
10 c
11 c Contains the routines:
12 c wrfun Write function (calculated spectra) into <nlsfun.nnn>
13 c wrspc Write spectra (experimental) and current fit into <dataid>.spc
14 c wrjac Write Jacobian
15 c wrtrdg Write tridiagonal matrices
16 c
17 c----------------------------------------------------------------------
18 c =========================
19 c subroutine WRFUN
20 c =========================
21 c
22 c Outputs the current calculated spectrum matrix (spectr in commmon
23 c /mspctr/ into an ASCII file named <nlsfun.nnn> where nnn is a three-
24 c digit iteration number (passed as an argument)
25 c----------------------------------------------------------------------
26 c
27  subroutine wrfun(iter)
28 c
29  use nlsdim
30  use expdat
31  use parcom
32  use lmcom
33  use mspctr
34  use stdio
35 c
36  integer iter
37  integer i,isp,ixs,j,k
38  character tmpnam*40
39 c
40  integer itrim
41  external itrim
42 c
43 c --------------------------------------------------------
44 c Open file for printing iterates (function evaluations)
45 c --------------------------------------------------------
46  write(tmpnam,2001) iter
47  open(unit=ludisk,file=tmpnam,status='unknown',
48  # access='sequential',form='formatted')
49 c
50 c --------------------
51 c Loop over spectra
52 c --------------------
53  do isp=1,nser
54  ixs=ixsp(isp)
55 c
56 c ----------------------------------------------------------
57 c Output calculation information and spectra for all sites
58 c ----------------------------------------------------------
59 c
60  write (ludisk,3005) (tag(i)(:itrim(tag(i))),x(i),i=1,nprm)
61  write (ludisk,3006) isp,sbi(isp)-shft(isp)-tmpshft(isp),
62  # sdb(isp),(sfac(j,isp),j=1,nsite)
63  do i=1,npts(isp)
64  k=ixs+i-1
65  write (ludisk,2003) (spectr(k,j),j=1,nsite)
66  end do
67 c
68  end do
69 c
70  close(ludisk)
71  return
72 c
73 c######################################################################
74 c
75  2001 format('nlsfun.',i3.3)
76  2003 format(10(g15.8,1x))
77  3005 format('c ',a9,' = ',f10.5)
78  3006 format('c Spectrum ',i2,', fieldi=',f9.3,' step=',f7.4/
79  # 'c Scales:'/
80  # 'c',10(g15.8,1x))
81  end
82 
83 c----------------------------------------------------------------------
84 c =========================
85 c subroutine WRSPC
86 c =========================
87 c
88 c Outputs the current experimental spectra and calculated function
89 c into the files <dataID1>.spc <dataID2>.spc, ... for each member
90 c of the series being fit.
91 c----------------------------------------------------------------------
92  subroutine wrspc()
93 c
94  use nlsdim
95  use expdat
96  use parcom
97  use lmcom
98  use iterat
99  use nlsnam
100  use mspctr
101  use stdio
102 c
103  implicit none
104  integer i,j,k,l
105  double precision field
106 c
107  do i=1,nspc
108  call setdat( dataid(i) )
109 c
110  open(ludisk,file=spname(:lthdnm),status='unknown',
111  # access='sequential',form='formatted')
112  field=sbi(i)
113  do j=1,npts(i)
114  k=ixsp(i)+j-1
115 c *** Multiple sites
116  if (nsite.gt.1) then
117 c *** weighted residuals
118  if (iwflag.ne.0) then
119  write(ludisk,1045) field,data(k),
120  # data(k)-rmsn(i)*fvec(k),
121  # (sfac(l,i)*spectr(k,l),l=1,nsite)
122 
123 c *** unweighted residuals
124  else
125  write(ludisk,1045) field,data(k),
126  # data(k)-fvec(k),
127  # (sfac(l,i)*spectr(k,l),l=1,nsite)
128  end if
129 
130 c *** Single site
131  else
132 c *** weighted residuals
133  if (iwflag.ne.0) then
134  write(ludisk,1045) field,data(k),
135  # data(k)-rmsn(i)*fvec(k)
136 
137 c *** unweighted residuals
138  else
139 
140  write(ludisk,1045) field,data(k),
141  # data(k)-fvec(k)
142  end if
143 
144  end if
145  field=field+sdb(i)
146  end do
147  close (ludisk)
148 c
149  end do
150  written=1
151  return
152 c
153  1045 format(f10.3,6(' ',g14.7))
154 c
155  end
156 
157 c----------------------------------------------------------------------
158 c =========================
159 c subroutine WRJAC
160 c =========================
161 c
162 c Outputs the current Jacobian matrix in array fjac (commmon /lmcom/
163 c into an ASCII file named <nlsjac.nnn> where nnn is a three-digit
164 c iteration number (passed as an argument)
165 c----------------------------------------------------------------------
166  subroutine wrjac(iter)
167 c
168  use nlsdim
169  use expdat
170  use lmcom
171  use parcom
172  use stdio
173 c
174  implicit none
175  integer iter,i,j
176  character tmpnam*40
177 c
178  write(tmpnam,1000) iter
179  open(unit=ludisk,file=tmpnam,status='unknown',
180  # access='sequential',form='formatted')
181 c
182  write (ludisk,1002) (tag(j),j=1,njcol)
183  write (ludisk,1003) (xfdstp(j),j=1,njcol)
184 c
185  do i=1,ndatot
186  write (ludisk,1001) (fjac(i,j),j=1,njcol)
187  end do
188 c
189  close(ludisk)
190  return
191 c
192 c######################################################################
193 c
194  1000 format('nlsjac.',i3.3)
195  1001 format(10(g15.8,1x))
196  1002 format('c ',10(a9,6x))
197  1003 format('c ',10(g15.7))
198  end
199 
200 
201 
202 c----------------------------------------------------------------------
203 c =========================
204 c subroutine WRTRDG
205 c =========================
206 c
207 c Outputs the current set of tridiagonal matrices into an ASCII
208 c file named <nlstri.nnn> where nnn is a three-digit iteration
209 c number (passed as an argument)
210 c----------------------------------------------------------------------
211  subroutine wrtrdg(iter)
212 c
213  use nlsdim
214  use expdat
215  use tridag
216  use parcom
217  use stdio
218 c
219  implicit none
220  integer iter,isp,isi,ixt,j
221  character tmpnam*40
222 c
223  write(tmpnam,1000) iter
224  open(unit=ludisk,file=tmpnam,status='unknown',
225  # access='sequential',form='formatted')
226 c
227  do isp=1,nser
228  do isi=1,nsite
229  ixt=ixtd(isi,isp)
230  write (ludisk,1001) isp,isi
231  write (ludisk,1002) (j,alpha(ixt+j-1),
232  # beta(ixt+j-1),j=1,ltd(isi,isp))
233  end do
234  end do
235 c
236  close(ludisk)
237  return
238 c
239 c######################################################################
240 c
241  1000 format('nlstri.',i3.3)
242  1001 format('c Spectrum ',i2,', site ',i2)
243  1002 format(i4,' (',g12.5,',',g12.5,') (',g12.5,',',g12.5,')')
244  end
245 
double precision, dimension(mxspc), save sbi
Definition: expdat.f90:40
double complex, dimension(mxtdg), save alpha
Definition: tridag.f90:32
integer, save ndatot
Definition: expdat.f90:45
double precision, dimension(mxspc), save sdb
Definition: expdat.f90:40
subroutine wrtrdg(iter)
Definition: writr.f90:212
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
double precision, dimension(mxspc), save tmpshft
Definition: expdat.f90:40
double precision, dimension(mxsite, mxspc), save sfac
Definition: mspctr.f90:33
double precision, dimension(mxspc), save shft
Definition: expdat.f90:40
subroutine wrjac(iter)
Definition: writr.f90:167
Definition: stdio.f90:26
character *30, dimension(mxspc), save dataid
Definition: expdat.f90:51
subroutine wrspc()
Definition: writr.f90:93
integer, save lthdnm
Definition: nlsnam.f90:27
integer, save iwflag
Definition: iterat.f90:14
double precision, dimension(mxvar), save xfdstp
Definition: parcom.f90:56
integer, save written
Definition: expdat.f90:45
integer, save nprm
Definition: parcom.f90:62
integer, save nspc
Definition: expdat.f90:45
integer, save nsite
Definition: parcom.f90:62
integer, dimension(mxsite, mxspc), save ltd
Definition: tridag.f90:28
subroutine setdat(dataid)
Definition: setnm.f90:67
double precision, dimension(mxpt, mxjcol), save fjac
Definition: lmcom.f90:17
double precision, dimension(mxpt), save fvec
Definition: lmcom.f90:17
integer, save nser
Definition: parcom.f90:62
subroutine wrfun(iter)
Definition: writr.f90:28
integer, dimension(mxsite, mxspc), save ixtd
Definition: tridag.f90:28
Definition: lmcom.f90:13
integer, parameter ludisk
Definition: stdio.f90:29
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
double complex, dimension(mxtdg), save beta
Definition: tridag.f90:32
character *30, save spname
Definition: nlsnam.f90:28
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
double precision, dimension(mxpt, mxsite), save spectr
Definition: mspctr.f90:33
double precision, dimension(mxjcol), save x
Definition: lmcom.f90:17
integer, save njcol
Definition: parcom.f90:62