NLSL
sshift.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 7/1/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine SSHIFT
5 c =========================
6 c
7 c Routine to calculate the shift needed to obtain optimal overlap
8 c between a spectrum (or set of spectra) and experimental data
9 c
10 c Inputs:
11 c data(npt) Data values to be fit by shifting/scaling
12 c spectr(ldspct,nsi) Spectrum/spectra to be shifted/scaled
13 c ldspct Leading dimension of spectr array
14 c nsi # spectra contained in spectr array
15 c npt # points in data/spectra arrays
16 c nft Array length for FFT (next power of 2 > npt)
17 c ideriv Flag =1 for 1st derivative data
18 c srange Allowed range for shifting (fraction of npt)
19 c ctol Tolerance for linearly independent components
20 c noneg Flag =1: zero any negative coeffients
21 c tmpdat(2*nft) Temporary storage for 0-padded data (for FFT)
22 c tmpclc(2*nft) Temporary storage for 0-padded spectrum (for FFT)
23 c wspec(ldspct,nsi) Work array used to store correlation functions
24 c work(2*nft) Temporary work array used in FFT and QR solution
25 c
26 c Output:
27 c sfac Scale factors for each site
28 c
29 c Includes:
30 c nlsdim.inc
31 c
32 c Uses:
33 c correl Calculate correlation function (from Numerical Recipes)
34 c qrfac QR factorization of normal equations (part of MINPACK)
35 c rsolve Solves the equation R*x = Q(transpose)*b
36 c qtbvec Calculates the vector Q(transpose)*b from the Householder
37 c factors of Q
38 c dchex Updates QR factorization (from LINPACK)
39 c fmomnt Returns first moment of an array (in index units)
40 c
41 c----------------------------------------------------------------------
42  function sshift( data,spct,ldspct,nsi,npt,nft,ideriv,srange,
43  # ctol,noneg,tmpdat,tmpclc,wspec,work,sfac )
44 c
45  use nlsdim
46 c
47  implicit none
48  double precision sshift
49 c
50  integer iarr,ideriv,ldspct,nsi,npt,nft,noneg
51  double precision spct(ldspct,nsi),data(npt),tmpdat(2*nft),
52  # tmpclc(2*nft),wspec(ldspct,nsi),work(2*nft),
53  # sfac(nsi),srange,ctol,err
54 c
55  integer i,info,irng,isi,ixmx,ixw1,ixw2,ixw3,j,jtmp,k,mneg,
56  # mnrng,mxrng
57  integer jpvt(mxsite)
58  double precision a,approx,asum,b,c,c1m,d1m,dummy,ovmax,
59  # scl,shift,smax,smin,temp,xfrac
60  double precision amat(mxsite,mxsite),qraux(mxsite),
61  # rdiag(mxsite),tmpscl(mxsite)
62 c
63  double precision ZERO
64  parameter(zero=0.0d0)
65 c
66 c double precision fmomnt,lgrint
67  double precision fmomnt
68  double complex lgrint
69  external fmomnt,lgrint
70 c
71 c######################################################################
72 c
73 c----------------------------------------------------------------------
74 c Determine the allowed range for shifting. First, find an
75 c approximate shift using the first moments of the experimental
76 c and calculated spectra
77 c----------------------------------------------------------------------
78  smax = dabs( srange*npt )
79  d1m=dabs( fmomnt( data,npt,ideriv ) )
80  c1m=zero
81  do i=1,nsi
82  c1m=c1m+dabs( fmomnt( spct(1,i),npt,ideriv ) )
83  end do
84 c
85  approx = d1m-c1m/dble(nsi)
86  if (approx.gt.smax .or. approx.lt.-smax) approx=zero
87  mnrng = max0( int(approx-smax)+nft/2, 1 )
88  mxrng = min0( int(approx+smax)+nft/2, npt )
89 c
90 c----------------------------------------------------------------------
91 c *** Make a copy of the data and zero-pad it for FFT
92 c (A temporary array is needed to accommodate zero-padding)
93 c----------------------------------------------------------------------
94  do j=1,npt
95  tmpdat(j)=data(j)
96  end do
97  do j=npt+1,nft
98  tmpdat(j)=zero
99  end do
100 c
101 c----------------------------------------------------------------------
102 c *** Loop over each site for this spectrum
103 c----------------------------------------------------------------------
104 c
105 c *** Copy each calculated spectrum to temporary array and zero-pad it
106 c (A temporary array is needed to accommodate zero-padding)
107 c
108  do isi=1,nsi
109  do j=1,npt
110  tmpclc(j)=spct(j,isi)
111  end do
112  do j=npt+1,nft
113  tmpclc(j)=zero
114  end do
115 c
116  call correl( tmpdat,tmpclc,nft,wspec(1,isi),work )
117 c
118  end do
119 c
120 c----------------------------------------------------------------------
121 c *** Calculate the normal equations matrix for linear least-squares
122 c----------------------------------------------------------------------
123  do j=1,nsi
124  do k=j,nsi
125  asum=zero
126  do i=1,npt
127  asum=asum+spct(i,j)*spct(i,k)
128  end do
129  amat(j,k)=asum
130  if (k.ne.j) amat(k,j)=asum
131  end do
132  end do
133 c
134 c----------------------------------------------------------------------
135 c Calculate the QR decomposition of the normal equation matrix
136 c----------------------------------------------------------------------
137  ixw1=1+nsi
138  ixw2=ixw1+nsi
139  ixw3=ixw2+nsi
140  call qrfac(nsi,nsi,amat,mxsite,.true.,jpvt,nsi,rdiag,
141  # work(ixw1),work(ixw2) )
142 c
143 c ---------------------------------------------------------------
144 c Store diagonal of Q (returned in diagonal of A) in work vector
145 c Replace diagonal of A with diagonal of R
146 c (this is the arrangement expected by qtbvec and rsolve)
147 c ---------------------------------------------------------------
148  do i=1,nsi
149  work(i)=amat(i,i)
150  amat(i,i)=rdiag(i)
151  end do
152 c
153 c -----------------------------------------------------------------
154 c Determine which spectra are linearly dependent and truncate them
155 c -----------------------------------------------------------------
156  k=0
157  do i=1,nsi
158  if (dabs(amat(i,i)) .le. ctol*dabs(amat(1,1)) ) go to 21
159  k=i
160  end do
161 c
162 c----------------------------------------------------------------------
163 c Calculate the correlation function of the data with the
164 c least-squares sum of spectra by solving the (possibly truncated)
165 c linear least-squares problem for each shift value in the
166 c allowed shifting range. These will form the vector part of
167 c the normal equations at the given shift value.
168 c Store the correlation functions in the tmpclc array.
169 c----------------------------------------------------------------------
170 c
171 c ----------------------------------------------------------------
172 c Loop through all possible shift values in discrete correlation
173 c function(s)
174 c ----------------------------------------------------------------
175  21 do irng=mnrng,mxrng
176 c
177 c ------------------------------------------------------------
178 c Store in tmpdat the values of the correlation functions of
179 c each individual site spectrum with the data for the current
180 c shift value
181 c ------------------------------------------------------------
182  do j=1,nsi
183  tmpdat(j)=wspec(irng,j)
184  end do
185 c
186 c -----------------------------------------------------------
187 c Calculate Q(trans)*vector (returned in work(1) array) and
188 c solve for scaling coefficients
189 c -----------------------------------------------------------
190  call qtbvec( nsi,k,amat,mxsite,work,tmpdat,work(ixw1) )
191  call rsolve( nsi,k,amat,mxsite,work,work(ixw1),tmpscl,
192  # .false.,work(ixw1) )
193 c
194 c ------------------------------------------------------------
195 c Calculate correlation function using optimal scale factors
196 c for the current shift value
197 c ------------------------------------------------------------
198  tmpclc(irng)=zero
199  do j=1,nsi
200  scl=tmpscl( jpvt(j) )
201  if (j.gt.k) scl=zero
202  tmpclc(irng)=tmpclc(irng)+scl*wspec(irng,j)
203  end do
204 c
205  end do
206 c
207 c--------------------------------------------------------------------------
208 c Search correlation function stored in tmpclc for the maximum value
209 c within the allowed shifting range.
210 c--------------------------------------------------------------------------
211  ixmx=mnrng
212  ovmax=zero
213  do irng=mnrng,mxrng
214  if (tmpclc(irng).gt.ovmax) then
215  ixmx=irng
216  ovmax=tmpclc(irng)
217  end if
218  end do
219 c
220 c ------------------------------------------------------------------
221 c Find the ordinate for the extremum (in this case a known maximum)
222 c ------------------------------------------------------------------
223 c if (ixmx.gt.mnrng .and. ixmx.lt.mxrng) then
224 c xfrac = -(tmpclc(ixmx+1)-tmpclc(ixmx-1))/
225 c # (tmpclc(ixmx+1)+tmpclc(ixmx-1)-2.0d0*tmpclc(ixmx))
226 c else
227 c xfrac=ZERO
228 c end if
229 c
230  iarr=max(ixmx-2,mnrng)
231  iarr=min(iarr,mxrng-4)
232  xfrac=lgrint(tmpclc(iarr),err)
233 c
234  shift=dble(ixmx-(nft/2)-1)+xfrac
235 c
236 c ----------------------------------------------------------------------
237 c Find the values of the scaling coefficients for the calculated
238 c fractional shift: place interpolated values of correlation functions
239 c in tmpdat
240 c ----------------------------------------------------------------------
241  if (ixmx.gt.mnrng.and. ixmx.lt.mxrng) then
242  do j=1,nsi
243  a = wspec(ixmx,j)
244  b = 0.5d0*(wspec(ixmx+1,j)-wspec(ixmx-1,j))
245  c = 0.5d0*(wspec(ixmx+1,j)+wspec(ixmx-1,j))-a
246  tmpdat(j)= a + b*xfrac + c*xfrac*xfrac
247  end do
248 c
249  else
250  do j=1,nsi
251  tmpdat(j)=wspec(ixmx,j)
252  end do
253  end if
254 c
255 c --------------------------------------------------
256 c Now solve the linear least-squares one more time
257 c to find the spectral scaling coefficients at the
258 c correlation function maximum
259 c --------------------------------------------------
260  call qtbvec( nsi,k,amat,mxsite,work,tmpdat,work(ixw1) )
261  call rsolve( nsi,k,amat,mxsite,work,work(ixw1),tmpscl,.false.,
262  # dummy )
263 c
264 c ------------------------------------------------------------
265 c Optional check for negative coefficients: permute spectra
266 c with negative coefficients into the truncated part of the QR
267 c factorization and re-solve the problem
268 c
269 c (But don't truncate the last spectrum!)
270 c ------------------------------------------------------------
271 c
272  if (noneg.ne.0.and.k.gt.1) then
273  smin=zero
274  mneg=0
275  do i=1,k
276  if (tmpscl(i).lt.smin) then
277  mneg=i
278  smin=tmpscl(i)
279  end if
280  end do
281 c
282  if (mneg.ne.0) then
283  if (mneg.lt.k) then
284  call dchex(amat,mxsite,nsi,mneg,k,work(ixw1),
285  # mxsite,1,work(ixw2),work(ixw3),2)
286 c
287  jtmp=jpvt(mneg)
288  do j=mneg,k-1
289  jpvt(j)=jpvt(j+1)
290  end do
291  jpvt(k)=jtmp
292  end if
293  k=k-1
294  go to 21
295  end if
296  end if
297 c
298 c ------------------------------------------------------------
299 c Zero the unused components of sfac and initialize the
300 c permutation array jpvt for unscrambling from pivoted order
301 c ------------------------------------------------------------
302  do i=1,nsi
303  jpvt(i) = -jpvt(i)
304  sfac(i) = tmpscl(i)
305  if (i.gt.k) sfac(i)=zero
306  end do
307 c
308 c ------------------------------------------------------------
309 c Unscramble the solution from its pivoted order using jpvt
310 c ------------------------------------------------------------
311  do 70 i=1,nsi
312  if (jpvt(i).le.0) then
313  k=-jpvt(i)
314  jpvt(i)=k
315 c
316  50 if (k.eq.i) goto 70
317  temp=sfac(i)
318  sfac(i)=sfac(k)
319  sfac(k)=temp
320  temp=rdiag(i)
321  rdiag(i)=rdiag(k)
322  rdiag(k)=temp
323  jpvt(k)=-jpvt(k)
324  k=jpvt(k)
325  go to 50
326  end if
327  70 continue
328 c
329  sshift=shift
330  return
331  end
332 
333 
334 
335 
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
double precision function sshift(data, spct, ldspct, nsi, npt, nft, ideriv, srange, ctol, noneg, tmpdat, tmpclc, wspec, work, sfac)
Definition: sshift.f90:44
subroutine correl(data, spctr, npt, ans, fft)
Definition: correl.f90:28
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