| |
(14) |
Quadratic convergence means that the square of the error
at one iteration is proportional to the error at the next iteration
| |
(15) |
Another interesting feature of the Newton iteration is that all iterations (except possibly the initial guess) are above the ultimate square root. This is obvious from equation (15).
We can insert spectral functions in the Newton square-root iteration,
for example
and
.Where the first guess a0 happens to match
,it will match
at all iterations.
Something inspires Wilson to write the Newton iteration as
| (16) |
| |
(17) |
Now we are ready for the algorithm: Compute the right side of (17) by polynomial division forwards and backwards and then add 1. Then abandon negative lags and take half of the zero lag. Now you have At+1(Z) / At(Z). Multiply out (convolve) the denominator At(Z), and you have the desired result At+1(Z). Iterate as long as you wish.
(Parenthetically, for those people familiar with the idea of minimum phase (if not, see FGDP or PVI), we show that At+1(Z) is minimum phase: Both sides of (17) are positive, as noted earlier. Both terms on the right are positive. Since the Newton iteration always overestimates, the 1 dominates the rightmost term. After masking off the negative powers of Z (and half the zero power), the right side of (17) adds two wavelets. The 1/2 is wholly real, and hence its real part always dominates the real part of the rightmost term. Thus (after masking negative powers) the wavelet on the right side of (17) has a positive real part, so the phase cannot loop about the origin. This wavelet multiplies At(Z) to give the final wavelet At+1(Z) and the product of two minimum-phase wavelets is minimum phase.)
The input of the program is the spectrum S(Z)
and the output is the factor A(Z),
a function with the spectrum S(Z).
I mention here that in later chapters of this book,
the factor A(Z) is known as the inverse Prediction-Error Filter (PEF).
In the Wilson-Burg code below,
S(Z) and A(Z) are Z-transform polynomials
but their lead coefficients are extracted off,
so for example,
module wilson { # Wilson's factorization
use helicon
use polydiv
integer, private :: n
real, dimension( :), allocatable, private :: auto, bb, cc, b, c
contains
subroutine wilson_init( nmax) {
integer, intent (in) :: nmax; n = nmax
allocate ( auto( 2*n-1), bb( 2*n-1), cc( 2*n-1), b(n), c(n))
}
subroutine wilson_factor( niter, s0, ss, a0, aa, verb) {
integer, intent( in) :: niter # Newton iterations
real, intent( in) :: s0 # autocorrelation zero lag
type( filter), intent( in) :: ss # autocorrelation, other lags
real, intent( out) :: a0 # factor, zero lag
type( filter) :: aa # factor, other lags
logical, intent( in) :: verb
optional :: verb
real :: eps
integer :: i, stat
auto = 0.; auto( n) = s0; b( 1) =1. # initialize
auto( n+ss%lag) = ss%flt # autocorrelation
auto( n-ss%lag) = ss%flt # symmetrize input auto.
call helicon_init( aa) # multiply polynoms
call polydiv_init( 2*n-1, aa) # divide polynoms
do i = 1, niter {
stat= polydiv_lop(.false.,.false., auto, bb) # bb = S/A
stat= polydiv_lop( .true.,.false., cc, bb) # cc = S/(AA')
b( 2:n) = 0.5*( cc( n+1:2*n-1 ) +
cc( n-1:1 :-1)) / cc( n) # b = plusside(1+cc)
eps = maxval( abs( b(2:n))); # "L1 norm"
if (present (verb)) { if (verb) write (0,*) i, eps}
stat= helicon_lop( .false., .false., b, c) # c = A b
aa%flt = c( 1+aa%lag) # put on helix
if( eps < epsilon( a0)) break # convergence
}
a0 = sqrt( cc( n))
}
subroutine wilson_close () {
deallocate( auto, bb, cc, b, c)
call polydiv_close()
}
}
is broken into the two parts a0 and aa.
#$=head1 NAME
#$
#$wilson - wilson's factorization
#$
#$=head1 SYNOPSIS
#$
#$C<call wilson_init(nmax)>
#$
#$C<call wilson_factor(niter,s0,ss,a0,aa)>
#$
#$C<call wilson_close()>
#$
#$=head1 PARAMETERS
#$
#$=over 4
#$
#$=item nmax - integer
#$
#$ maximum number of space needed (n1*10 good number)
#$
#$=item s0 - real
#$
#$ zero lag value of input
#$
#$=item ss - type(filter)
#$
#$ auto correlation
#$
#$=item a0 - real
#$
#$ Output zero lag value
#$
#$=item aa - type(filter)
#$
#$ Minimum phase filter
#$
#$=back
#$
#$=head1 DESCRIPTION
#$
#$Find a minumum phase filter with the given auto-correaltion function
#$
#$=head1 SEE ALSO
#$
#$L<cross_wilson>,L<lapfac>
#$
#$=head1 LIBRARY
#$
#$B<geef90>
#$
#$=cut
