NLSL
sscale.f90
Go to the documentation of this file.
1 c NLSL Version 1.3 7/17/93
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine SSCALE
5 c =========================
6 c
7 c Given the array <data> with <npt> data points and a matrix of <nsite>
8 c calculated spectra, each with <npt> points, calculate a least-squares
9 c set of <nsite> scaling factors. This routine uses the QR factorization
10 c routine that belongs to MINPACK together with the associated routines
11 c for solving the least-squares problem using QR factoriaation.
12 c
13 c Inputs:
14 c data(npt) Data values to be fit by scaled spectra
15 c spectr(ldspct,nsi) Spectrum/spectra to be scaled
16 c ldspct Leading dimension of spectr array
17 c nsite # spectra contained in spectr array
18 c npt # points in data/spectra arrays
19 c ctol Tolerance for linearly independent components
20 c iscal Flags =1 for automatic scaling of each site
21 c (otherwise used fixed value passed in sfac)
22 c noneg Flag =1: zero any negative coeffients
23 c work(npt) Temporary work array used in QR solution
24 c
25 c Output:
26 c sfac Array of scaling factors
27 c Scale factors corresponding to spectra that
28 c were eliminated on the basis of linear dependence
29 c or were negative (if noneg=1)
30 c
31 c resid Residuals (data minus scaled calculated spectra)
32 c
33 c Includes:
34 c nlsdim.inc
35 c
36 c Uses:
37 c qrfac Calculate QR factorization of design matrix (MINPACK)
38 c qtbvec Calculate Q(transpose)*b from Householder transformations
39 c rsolve Solve R*x = Q(transpose)*b
40 c dchex Update QR factors for column permutation of design matrix
41 c
42 c----------------------------------------------------------------------
43  subroutine sscale( data,spct,wspct,ldspct,work,nsite,npt,ctol,
44  # noneg,iscal,sfac,resid )
45 c
46  use nlsdim
47 c
48  implicit none
49 c
50  integer ldspct,noneg,npt,nsite
51  integer iscal(nsite)
52  double precision data(npt),work(npt),spct(ldspct,nsite),
53  # wspct(ldspct,nsite),resid(npt),sfac(nsite),ctol
54 c
55  integer ixspc(mxsite),jpvt(mxsite)
56  double precision qtd(mxsite),qraux(mxsite),rdiag(mxsite),
57  # tmpfac(mxsite),wa1(mxsite),wa2(mxsite)
58 c
59  integer i,info,j,jtmp,k,m,mneg,nscl
60  double precision smin,sumc2,sumdc,tmp
61 c
62  double precision ZERO
63  parameter(zero=0.0d0)
64 c
65 c
66 c######################################################################
67 c
68  do j=1,npt
69  resid(j)=data(j)
70  end do
71 c
72 c --------------------------------------------------------------
73 c Copy any spectra with the autoscale flag set to the work array
74 c Subtract spectra with fixed scale factors from the residuals
75 c --------------------------------------------------------------
76  nscl=0
77  do i=1,nsite
78  if (iscal(i).ne.0) then
79  nscl=nscl+1
80  ixspc(nscl)=i
81  do j=1,npt
82  wspct(j,nscl)=spct(j,i)
83  end do
84  else
85  do j=1,npt
86  resid(j)=resid(j)-sfac(i)*spct(j,i)
87  end do
88  end if
89  end do
90 c
91  if (nscl.eq.0) return
92 c
93 c----------------------------------------------------------------------
94 c If there is only one site, calculate the least-squares
95 c scaling factor directly
96 c----------------------------------------------------------------------
97  if (nscl.eq.1) then
98 c
99  sumdc=zero
100  sumc2=zero
101  do i=1,npt
102  sumdc=sumdc+resid(i)*wspct(i,1)
103  sumc2=sumc2+wspct(i,1)*wspct(i,1)
104  end do
105 c
106  if (sumc2.ne.zero) then
107  tmpfac(1)=sumdc/sumc2
108  else
109  tmpfac(1)=zero
110  end if
111 c
112 c---------------------------------------------------------------
113 c For multiple sites, use QR decomposition to find linear
114 c least-squares scaling coefficients and rms deviation
115 c---------------------------------------------------------------
116  else
117 c
118 c --------------------------------------------
119 c Compute the QR factorization of the spectra
120 c --------------------------------------------
121  call qrfac(npt,nscl,wspct,ldspct,.true.,jpvt,nscl,rdiag,
122  # wa1,wa2)
123 c
124  do i=1,nscl
125  qraux(i)=wspct(i,i)
126  wspct(i,i)=rdiag(i)
127  end do
128 c
129 c --------------------------------------------------------------
130 c Determine which spectra are linearly dependent and remove them
131 c from the problem
132 c --------------------------------------------------------------
133  k=0
134  do i=1,nscl
135  if (dabs(rdiag(i)) .le. ctol*dabs(rdiag(1)) ) goto 14
136  k=i
137  end do
138 c
139  14 call qtbvec(npt,k,wspct,ldspct,qraux,data,work)
140 c
141 c ---------------------------------------------------------
142 c Loop to here if QR factorization has been updated
143 c to eliminatespectra entering with negative coefficients
144 c ---------------------------------------------------------
145 c
146  15 call rsolve( npt,k,wspct,ldspct,qraux,work,tmpfac,
147  # .false.,tmp )
148 c
149 c ------------------------------------------------------------
150 c Optional check for negative coefficients: permute spectra
151 c with negative coefficients into the truncated part of the QR
152 c factorization and re-solve the problem
153 c (But don't truncate the last spectrum!)
154 c ------------------------------------------------------------
155  if (noneg.ne.0.and.k.gt.1) then
156  smin=zero
157  mneg=0
158  do i=1,k
159  if (tmpfac(i).lt.smin) then
160  mneg=i
161  smin=tmpfac(i)
162  end if
163  end do
164 c
165  if (mneg.ne.0) then
166  if (mneg.lt.k) then
167  call dchex( wspct,ldspct,nscl,mneg,k,
168  # work,ldspct,1,wa1,wa2,2 )
169 c
170  jtmp=jpvt(mneg)
171  do j=mneg,k-1
172  jpvt(j)=jpvt(j+1)
173  end do
174  jpvt(k)=jtmp
175  end if
176  k=k-1
177  go to 15
178  end if
179  end if
180 c
181 c ---------------------------------------------------------------
182 c Set the unused components of tmpfac to zero and initialize jpvt
183 c for unscrambling scaling coefficients from pivoted order
184 c ---------------------------------------------------------------
185  do i=1,nscl
186  jpvt(i)=-jpvt(i)
187  if (i.gt.k) tmpfac(i)=zero
188  end do
189 c
190 c -----------------------------------------------------
191 c Unscramble the solution from pivoted order using jpvt
192 c -----------------------------------------------------
193  do 70 i=1,nscl
194  if (jpvt(i).le.0) then
195  k=-jpvt(i)
196  jpvt(i)=k
197 c
198  50 if (k.eq.i) goto 70
199  tmp=tmpfac(i)
200  tmpfac(i)=tmpfac(k)
201  tmpfac(k)=tmp
202  tmp=rdiag(i)
203  rdiag(i)=rdiag(k)
204  rdiag(k)=tmp
205  jpvt(k)=-jpvt(k)
206  k=jpvt(k)
207  go to 50
208  end if
209  70 continue
210 c
211 c *** End if (nscl.eq.1) ... else ...
212  end if
213 c
214 c ------------------------------------------------------
215 c Calculate residuals using least-squares scale factors
216 c ------------------------------------------------------
217  do i=1,nscl
218  k=ixspc(i)
219  sfac(k)=tmpfac(i)
220  do j=1,npt
221  resid(j)=resid(j)-sfac(k)*spct(j,k)
222  end do
223  end do
224 c
225  return
226  end
subroutine qrfac(m, n, a, lda, pivot, ipvt, lipvt, rdiag, acnorm, wa)
Definition: qrfac.f90:2
subroutine dchex(r, ldr, p, k, l, z, ldz, nz, c, s, job)
Definition: dchex.f90:2
integer, parameter mxsite
Definition: nlsdim.f90:39
subroutine qtbvec(m, n, q, ldq, qraux, b, qtb)
Definition: qrutil.f90:39
subroutine rsolve(m, n, q, ldq, qraux, qtb, x, rcalc, rsd)
Definition: qrutil.f90:117
subroutine sscale(data, spct, wspct, ldspct, work, nsite, npt, ctol, noneg, iscal, sfac, resid)
Definition: sscale.f90:45