aovsub.f90 Source File


Contents

Source Code


Source Code

Module aovsub
!! Utility routines for Aov.<br>
!! Note that only `public` routines are reliable and tested.
!! Other ones are kept here for past/future development.
!
      Use aovconst
!
      Implicit None
!
      Private
      Public sortx, givensacc, givens, tolower, toupper, dict_val
      Interface sortx
         Module Procedure sortx_m_r4, sortx_m_r8
      End Interface sortx
!
      Character (Len=*), Parameter :: Lowcase = 'abcdefghijklmnopqrstuvwxyz'
      Character (Len=*), Parameter :: Upcase = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
!
Contains
!
      Function dict_val (dict, isize, key)
!! Find key position in a keyword dictionary.
!! Restriction: dictionary entries are of equal size,
!! all kept in one string.
         Integer :: dict_val !! Key position
         Character (Len=*), Intent (In) :: dict !! Dictionary string
         Character (Len=*), Intent (In) :: key !! key string
         Integer, Intent (In) :: isize !! size of keys
         dict_val = index (dict, key)
         If (dict_val > 0) dict_val = dict_val / isize + 1
      End Function dict_val
!
!
      Function tolower (str)!! Change string case to lower
         Character (Len=*), Intent (In) :: str !! input string
         Character (Len=Len(str)) :: tolower !! result output string
!
         Integer l, i, j
         l = len (str)
         tolower = str
         Do i = 1, l
            j = index (Upcase, tolower(i:i))
            If (j > 0) tolower (i:i) = Lowcase (j:j)
         End Do
      End Function tolower
!
      Function toupper (str)!! Change string case to upper
         Character (Len=*), Intent (In) :: str !! input string
         Character (Len=Len(str)) :: toupper !! result output string
!
         Integer l, i, j
         l = len (str)
         toupper = str
         Do i = 1, l
            j = index (Lowcase, toupper(i:i))
            If (j > 0) toupper (i:i) = Upcase (j:j)
         End Do
      End Function toupper
!
      Subroutine givensacc (a, r)
!! Purpose: GIVENSACC prepares the least squares solution (LSQ) of an
!! overdetermined linear system by Givens rotations, a state-of-art
!! algorithm. Routine takes current observation equation from a and
!! accumulates it into r:
!! (C) Alex Schwarzenberg-Czerny,2000      alex@camk.edu.pl
         Implicit None
         Real (SP), Intent (In) :: a (:)
!! a(n1)    - an equation with its RHS as a(n1) element,
!!            multiplied by sqrt of weight.
!!            Note: sum of all weights must equall no.observations,
!!            lest GIVENS returns wrong covariance estimates.
         Real (SP), Intent (Inout) :: r (:, :)
!! Input: r(n1,n1) - triangle system R and its R.H.S
!!            set r = 0 before te first call of GIVENSACC.
!! Output: r(n1,n1) - updated system accounting of a
         Integer :: n1, n, i, j
         Real (SP) :: e (size(a)), p, q, s, rii, ei
         n1 = size (a)
         n = n1 - 1
         e = a
         r (n1, n1) = r (n1, n1) + e (n1) * e (n1)
         Do i = 1, n
            rii = r (i, i)
            ei = e (i)
            s = rii * rii + ei * ei
            If (s <= 0._SP) Then
               p = 1._SP
               q = 0._SP
            Else
               s = sign (Sqrt(s), rii)
               p = rii / s
               q = ei / s
            End If
            r (i, i) = s
            e (i) = 0._SP
            Do j = i + 1, n1
               s = q * e (j) + p * r (i, j)
               e (j) = p * e (j) - q * r (i, j)
               r (i, j) = s
            End Do
         End Do
!
      End Subroutine givensacc
!
      Function givens (r, idf, detnth)
!! Purpose: GIVENS calculates the least squares solution of an
!! overdetermined linear system by accurate Givens rotations.
!! Prior to GIVENS call m-times GIVENSA in order to build the
!! triangle system R
!! (C) Alex Schwarzenberg-Czerny,2000      alex@camk.edu.pl
         Implicit None
         Real (SP) :: givens !! result condition ratio (min eigenv/max eigenv),
!! 0 signalls singular system
         Integer, Intent (In) :: idf
!! idf - >0 solve normal equations for idf degrees of freedom
!!       <0 and r(n1,n1)=1: inverse triangle matrix R
         Real (SP), Intent (Inout) :: r (:, :)
!! Input: r(n1,n1) - triangle system R and its R.H.S
!! Output: For idf>0:
!! r(n1,:)   - unknowns
!! (r(k,k),k=1,n) variances
!! (r(:k,k),k=1,n) half of covariance matrix
!! (r(k+1:,k),k=2,n) half of inverse triangle root Q
!! r(n1,n1) - chi^2/idf (=RHS variance=std.dev^2)
!!         For idf<0:
!! ((r(j,k),j=1,k),k=1,n) inverse triangle root  Q
         Real (SP), Intent (Out) :: detnth !! =Det(R)^(1/n) - n-th root of determinant
!
         Integer :: n1, n, i, j
         Real (SP), Parameter :: tol = 30._SP * epsilon (1._SP)
         Real (SP) :: rii, riimx, riimn
!
         n1 = size (r, dim=1)
         n = n1 - 1
         r (n1, n1) = (r(n1, n1)-dot_product(r(1:n, n1), r(1:n, n1))) / &
            Max (1, idf)
         Do i = 1, n
            rii = 0._SP
            If (Abs(r(i, i)) > tol) rii = 1._SP / r (i, i)
            r (i, i) = - 1._SP
            Do j = 1, i
               r (i, j) = r (i, j) * rii
               r (i+1, j) = - dot_product (r(j:i, j), r(j:i, i+1))
!en(j,i+1)=-dot_product(en(j,k:i),en(i+1,k:i))
            End Do
         End Do
!Already solved, square Q to get  inverse/covariance
         detnth = 0.
         riimx = Abs (r(1, 1))
         riimn = riimx
         Do i = 1, n
            rii = Abs (r(i, i))
            If (idf .Gt. 0) Then
               Do j = 1, i
                  r (j, i) = dot_product (r(i:n, i), r(i:n, j)) * r (n1, n1)
               End Do
            End If
            If (rii > 0._SP) Then
               detnth = detnth - Log (rii)
            End If
            riimx = Max (rii, riimx)
            riimn = Min (rii, riimn)
         End Do
         givens = riimn / riimx
         detnth = Exp (2._SP*detnth/n)
!
      End Function givens
!
      Subroutine sortx_h (dat, ind)
!! (C) Alex Schwarzenberg-Czerny, 2001,2005 alex@camk.edu.pl
!! Unstable index sorting. Method based on heap sort routine in
!! Numerical Recipes, 1992, Press et al., Cambridge U.P.
         Implicit None
         Real (SP), Intent (In) :: dat (:)!! unsorted data
         Integer, Intent (Out) :: ind (:)!! sort index
!
         Integer :: i, j, n
         n = size (dat)
         If (size(ind) < n) Then
            Write (*,*) 'sortx_h: wrong dimension'
            Stop
         End If
!
         Do i = 1, n
            ind (i) = i
         End Do
!
         Do i = n / 2, 1, - 1
            Call put (i, n)
         End Do
         Do i = n, 2, - 1
            j = ind (1)
            ind (1) = ind (i)
            ind (i) = j
            Call put (1, i-1)
         End Do
!
      Contains
!
         Subroutine put (low, up)
! Swap index values when needed
! Auxillary routine for hep sort
! ind,dat are global
            Integer, Intent (In) :: low ! lower index
            Integer, Intent (In) :: up ! upper index
!
            Integer :: i, j, jsav
            Real (SP) :: sav
            jsav = ind (low)
            sav = dat (jsav)
            j = low
            i = low + low
            Do
               If (i > up) Exit
               If (i < up) Then
                  If (dat(ind(i)) < dat(ind(i+1))) i = i + 1
               End If
               If (sav >= dat(ind(i))) Exit
               ind (j) = ind (i)
               j = i
               i = i + i
            End Do
            ind (j) = jsav
!
         End Subroutine put
!
      End Subroutine sortx_h
!
      Subroutine sortx_m_r4 (dat, ina)
!! (C) Alex Schwarzenberg-Czerny, 2005      alex@camk.edu.pl
!! Stable index sorting by merge sort. Requires extra memory (inb)
         Implicit None
         Real (4), Intent (In) :: dat (:)!! unsorted data
         Integer, Intent (Out) :: ina (:)!! sort index
!
         Integer :: inb ((size(dat)+1)/2), i, n
!
         n = size (dat)
         If (size(ina) < n) Then
            Write (*,*) 'sortx_m: wrong dimension'
            Stop
         End If
!
         Forall (i=1:n) ina (i) = i
!
         Call sort_part (1, n)
!
      Contains
!
         Recursive Subroutine sort_part (low, up)
! Sort part of data between incdices
! Auxillary recursive subroutine for merge sort
            Integer, Intent (In) :: low ! lower sort index
            Integer, Intent (In) :: up ! upper sort index
!
            Integer :: med
!
            If (low < up) Then
               med = (low+up) / 2
               Call sort_part (low, med)
               Call sort_part (med+1, up)
               Call merge_part (low, med, up)
            End If
         End Subroutine sort_part
!
         Subroutine merge_part (low, med, up)
!! merge 2 sorted parts, divided by med
            Integer, Intent (In) :: low !! lower index
            Integer, Intent (In) :: med !! middle index
            Integer, Intent (In) :: up !! upper index
!
            Integer :: i, j, r
!
            inb (1:med-low+1) = ina (low:med)
!
            i = 1
            r = low
            j = med + 1
            Do
               If (r >= j .Or. j > up) Exit
               If (dat(inb(i)) <= dat(ina(j))) Then
                  ina (r) = inb (i)
                  i = i + 1
               Else
                  ina (r) = ina (j)
                  j = j + 1
               End If
               r = r + 1
            End Do
!
            ina (r:j-1) = inb (i:i+j-1-r)
         End Subroutine merge_part
!
      End Subroutine sortx_m_r4
!
      Subroutine sortx_m_r8 (dat, ina)
!! (C) Alex Schwarzenberg-Czerny, 2005      alex@camk.edu.pl
!! Stable index sorting by merge sort. Requires extra memory (inb)
         Implicit None
         Real (8), Intent (In) :: dat (:)!! unsorted data
         Integer, Intent (Out) :: ina (:)!! sort index
!
         Integer :: inb ((size(dat)+1)/2), i, n
!
         n = size (dat)
         If (size(ina) < n) Then
            Write (*,*) 'sortx_m: wrong dimension'
            Stop
         End If
!
         Forall (i=1:n) ina (i) = i
!
         Call sort_part (1, n)
!
      Contains
!
         Recursive Subroutine sort_part (low, up)
! Sort part of data between incdices
! Auxillary recursive subroutine for merge sort
            Integer, Intent (In) :: low ! lower sort index
            Integer, Intent (In) :: up ! upper sort index
!
            Integer :: med
!
            If (low < up) Then
               med = (low+up) / 2
               Call sort_part (low, med)
               Call sort_part (med+1, up)
               Call merge_part (low, med, up)
            End If
         End Subroutine sort_part
!
         Subroutine merge_part (low, med, up)
! merge 2 sorted parts, divided by med
! Auxillary subroutine for merge sort
            Integer, Intent (In) :: low ! lower index
            Integer, Intent (In) :: med ! middle index
            Integer, Intent (In) :: up ! upper index
!
            Integer :: i, j, r
!
            inb (1:med-low+1) = ina (low:med)
!
            i = 1
            r = low
            j = med + 1
            Do
               If (r >= j .Or. j > up) Exit
               If (dat(inb(i)) <= dat(ina(j))) Then
                  ina (r) = inb (i)
                  i = i + 1
               Else
                  ina (r) = ina (j)
                  j = j + 1
               End If
               r = r + 1
            End Do
!
            ina (r:j-1) = inb (i:i+j-1-r)
         End Subroutine merge_part
!
      End Subroutine sortx_m_r8
!
      Subroutine covar (t1, v1, e1, t2, v2, e2, nlag, eps, iscale, &
     & ifunct, n1, n2, lav, lmi, lmx, cc, cmi, cmx)
!! TO BE DONE:
!! -more flexibility in choosing the lag grid
!! Purpose:
!! Calculates cross-correlation function(CCF) of two time series sampled unevenly
!! in time. Note that for good result each hump in the time series should be covered
!! by several observations (over-sampling).
!!
!! Method/Reference:
!! \bibitem[Alexander(1997)]{1997ASSL..218..163A}
!! Alexander, T.\ 1997, Astronomical Time Series,
!! Eds. D. Maoz, A. Sternberg, and E.M. Leibowitz, 1997 (Dordrecht: Kluwer),218, 163
!! Title: Is AGN Variability Correlated with Other AGN Properties? ZDCF Analysis of
!! Small Samples of Sparse Light Curves
!
!!  Copyrights and Distribution:
!!  This package is subject to copyrights by its author,
!!  (C) A. Schwarzenberg-Czerny 2011,  alex@camk.edu.pl
!
         Integer, Intent (In) :: iscale !!   (ignored)output scale: /=1 -linear; ==1 -logarythmic
         Integer, Intent (In) :: ifunct !!   (ignored)output function: /=1 -correlation; ==1 -structure
         Real (SP), Intent (In) :: eps !!   max separation within lag
         Integer, Intent (In) :: n1 !!    number of observations in 1st set
         Real (TIME), Intent (In) :: t1 (n1)!!   1st set times
         Real (SP), Intent (In) :: v1 (n1)!!   1st set values
         Real (SP), Intent (In) :: e1 (n1)!!   1st set errors
         Integer, Intent (In) :: n2 !!    number of observations in 2nd set
         Real (TIME), Intent (In) :: t2 (n2)!!   2nd set times
         Real (SP), Intent (In) :: v2 (n2)!!   2nd set values
         Real (SP), Intent (In) :: e2 (n2)!!   2nd set errors
         Integer, Intent (In) :: nlag !!   number of lags
         Real (SP), Intent (Out) :: cc (nlag)!! average of correlations
         Real (SP), Intent (Out) :: cmi (nlag)!! min of correlations
         Real (SP), Intent (Out) :: cmx (nlag)!! max of correlations
         Real (SP), Intent (Out) :: lav (nlag)!! average of lags
         Real (SP), Intent (Out) :: lmi (nlag)!! min of lags
         Real (SP), Intent (Out) :: lmx (nlag)!! max of lags
! f2py Integer,Intent(IN),optional   :: iscale=0,ifunct=0
! f2py Real(SP),Intent(IN),optional  :: eps=0.
!
         Integer :: ind (n1*n2), mdt, mbin, mct, idt, ibin, ict, i, i1, &
        & i2, nct, ib1 (2*n1*n2/nlag), ib2 (2*n1*n2/nlag)
         Real (SP) :: dt (n1*n2), z, s, r, r2, ad1, ad2, sd1, sd2, sw ! n1*n2 could be large
         Real (SP), Allocatable :: db1 (:), db2 (:), wb1 (:), wb2 (:), &
        & wb (:), cb (:), lags (:)
!
         mdt = iscale + ifunct ! use dummy variables to prevent warnings
         mdt = n1 * n2
         nct = mdt / nlag
         mbin = mdt / nct
         mct = 2 * nct
         lav = 0._SP
         lmi = 0._SP
         lmx = 0._SP
         cc = 0._SP
         cmi = 0._SP
         cmx = 0._SP
!
         Forall (i1=1:n1, i2=1:n2) dt ((i1-1)*n2+i2) = real &
        & (t2(i2)-t1(i1), SP)
         Call sortx (dt, ind)
!
         idt = 0
         ibin = 0
         Do ! increment bins
            If (ibin >= mbin) Exit
            ibin = ibin + 1
            ict = 0
!
            Do ! increment pairs
               If (idt >= mdt) Exit
               idt = idt + 1
! check If lag bin finished
               If (ict > nct .And. dt(ind(idt))-dt(ind(Max(1, idt-1))) &
              & > eps .And. mdt-idt > nct/2) Exit
               i1 = ind (idt) / n2 + 1
               i2 = ind (idt) - (i1-1) * n2
!
               Do i = 1, ict ! check dependence
                  If (ib1(i) == i1 .Or. ib2(i) == i2) Exit
               End Do
               If (i <= ict) Then ! dependent
                  If (e1(i1)*e2(i2) < e1(ib1(i))*e2(ib2(i))) Then
                     ib1 (i) = i1
                     ib2 (i) = i2
                  End If
               Else ! independent
                  If (ict >= mct) Exit! this should never happen
                  ict = ict + 1
                  ib1 (ict) = i1
                  ib2 (ict) = i2
               End If
            End Do
!
! evaluate correlation within a bin
            Allocate (db1(ict), db2(ict), wb1(ict), wb2(ict), wb(ict), &
               cb(ict), lags(ict))
            db1 = v1 (ib1(:ict))
            wb1 = 0._SP
            Where (e1(ib1(:ict)) > 0._SP) wb1 = 1._SP / e1 (ib1(:ict))**2
            db2 = v2 (ib2(:ict))
            wb2 = 0._SP
            Where (e2(ib2(:ict)) > 0._SP) wb2 = 1._SP / e2 (ib2(:ict))**2
            wb = Sqrt (wb1*wb2)
            sw = sum (wb)
!
            ad1 = sum (db1*wb1) / sum (wb1)
            db1 = db1 - ad1 ! data variance and averages
            ad2 = sum (db2*wb2) / sum (wb2)
            db2 = db2 - ad2
            sd1 = sum (db1**2*wb1)
            sd2 = sum (db2**2*wb2)
!
 ! correlation and its variance in z-variable
            r = sum (db1*db2*wb) / Sqrt (sd1*sd2)
            cc (ibin) = r
            r2 = r * r
            z = (((((r2*3._SP+2._SP)*r2+11._SP)*0.5_SP/(ict-1)+r2+5._SP)* &
               0.25/(ict-1)+1._SP)*r/(ict-1)+Log(Max(1._SP-r, tiny(z))/ &
               (1._SP+r))) * 0.5_SP
            s = Sqrt (((((-r2*3._SP-6._SP)*r2+22._SP)/3._SP/(ict-1)+4._SP-r2) &
               *0.5/(ict-1)+1._SP)/(ict-1))
            cmi (ibin) = r - Abs (Tanh(z-s)-r)
            cmx (ibin) = r + Abs (Tanh(z+s)+r)
!
            lags = real (t2(ib2(:ict))-t1(ib1(:ict)), SP)! average & extreme intervals
            lav (ibin) = sum (lags*wb) / sw
            lmx (ibin) = maxval (lags)
            lmi (ibin) = minval (lags)
            Deallocate (db1, db2, wb1, wb2, wb, cb, lags)
!
            If (idt >= mdt) Exit
         End Do
      End Subroutine covar
!
      Subroutine covar1 (t1, d1, v1, t2, d2, v2, nlag, lag0, lstep, &
     & iscale, ifunct, n1, n2, lgs, cc, cv, av1, var1, av2, var2)
!!   Compute discrete covariance function for unequaly sampled signals
!!   Reference: Edelson R.A. and Krolik, J.H., 1988, ApJ 333, 646
!! this routine is incorrect (weighting not implementrd fully/correctly)
!!  Copyrights and Distribution:
!!  This package is subject to copyrights by its author,
!!  (C) A. Schwarzenberg-Czerny 1998-2011,  alex@camk.edu.pl
!
         Integer, Intent (In) :: iscale !!   (ignored)output scale: /=1 -linear; ==1 -logarythmic
         Integer, Intent (In) :: ifunct !!   (ignored)output function: /=1 -correlation; ==1 -structure
         Integer, Intent (In) :: n1 !!   number of observations in 1st set
         Real (TIME), Intent (In) :: t1 (n1)!!   1st set times
         Real (SP), Intent (In) :: d1 (n1)!!   1st set values
         Real (SP), Intent (In) :: v1 (n1)!!   1st set variances
         Integer, Intent (In) :: n2 !!   number of observations in 2nd set
         Real (TIME), Intent (In) :: t2 (n2)!!   2nd set times
         Real (SP), Intent (In) :: d2 (n2)!!   2nd set values
         Real (SP), Intent (In) :: v2 (n2)!!   2nd set variances
         Integer, Intent (In) :: nlag !!   number of lags
         Real (TIME), Intent (In) :: lag0 !!   first time lag
         Real (TIME), Intent (In) :: lstep !!   lag increment
         Real (SP), Intent (Out) :: lgs (nlag)!! time lag
         Real (SP), Intent (Out) :: cc (nlag)!! covariance
         Real (SP), Intent (Out) :: cv (nlag)!! its variance
         Real (SP), Intent (Out) :: av1 !! average 1-st set
         Real (SP), Intent (Out) :: var1 !! variance 1-st set
         Real (SP), Intent (Out) :: av2 !! average 2-nd set
         Real (SP), Intent (Out) :: var2 !! variance 2-nd set
!
         Integer :: iobs1, iobs2, ilag
!  Real(SP)   :: s1(n1),s2(n2)
         Real (SP) :: ct1, ct2, cnt, oldaver, aver, v, vi, ei, var12, t, &
            tlow, tup, tlag0
         Real (TIME) :: ti
!
         If (lstep <= 0.) Then
            Print *, 'COVAR: lstep must be positive'
            Return
         End If
         tlow = lag0
         tup = lag0 + (nlag-1) * lstep ! lag scale range
         If (iscale == 1) Then
            If (lag0 <= 0.) Then
               Print *, 'COVAR: lag0 must be positive for log scale'
               Return
            End If
            tlow = 10._SP ** tlow
            tup = 10._SP ** tup ! de-log lag scale range
         End If
!   Compute weighted averages and variances
!
         ct1 = sum (1._SP/v1)
         ct2 = sum (1._SP/v2)
         av1 = sum (d1/v1) / ct1
         av2 = sum (d2/v2) / ct2
         var1 = sum ((d1-av1)**2/v1) / (n1-1) * n1 / ct1
         var2 = sum ((d2-av2)**2/v2) / (n2-1) * n2 / ct2
         var12 = Sqrt (var1*var2)
!
!   Compute binned discrete covariance function
         cc = 0._SP
         cv = 0._SP
         lgs = 0._SP
         tlag0 = lag0
         Do iobs1 = 1, n1 ! loop over data 1
            ti = t1 (iobs1)
            vi = d1 (iobs1) - av1
            ei = v1 (iobs1)
            Do iobs2 = 1, n2 ! loop over all data 1 & data 2 pairs
               t = real (ti-t2(iobs2), Kind=SP)
               If ((t-tlow)*(tup-t) >= 0._SP) Then
                  If (iscale == 1) t = Log10 (t)! log scale
                  ilag = Nint ((t-tlag0)/lstep) + 1
                  If (ilag > 0 .And. ilag <= nlag) Then
                     v = av2 - d2 (iobs2) + vi
                     v = v * v - ei - v2 (iobs2)
!   Use rounding proof algorithm for accumulation of means and variances
                     cnt = lgs (ilag) + 1._SP
                     oldaver = cc (ilag)
                     aver = oldaver + (v-oldaver) / cnt
                     cc (ilag) = aver
                     cv (ilag) = cv (ilag) + (v-oldaver) * (v-aver)
                     lgs (ilag) = cnt
                  End If
               End If
            End Do
         End Do
!
!   Get errors and time lag grid
!
         Do ilag = 1, nlag
            If (lgs(ilag) > 1.5_SP) Then
               cv (ilag) = Sqrt (cv(ilag)/(lgs(ilag)-1._SP))
               If (ifunct /= 1) cv (ilag) = var12 - cv (ilag)! calculate CCF variance
            End If
            If (ifunct /= 1) cc (ilag) = var12 - cc (ilag)! calculate CCF
            lgs (ilag) = lag0 + lstep * (ilag-1)
         End Do
      End Subroutine covar1
!
      Subroutine fouw (t, valin, er, frin, nh2, frout, valout, cof, &
     & dcof, nobs)
!
!! FOUW - Fit Fourier series (cof) and Return vector of residuals
!!          (fout), adjusting frequency (fr) by NLSQ
!!
!!  Purpose:
!!   Twofold: (i) Find Fourier coefficients and improve frequency
!!   (ii) Pre-whiten data with a given frequency
!
!!  Method:
!!   Fits Fourier series of nh2/2 harmonics with frequency adjustment
!!   NLSQ by Newton-Raphson iterations:
!!   fin(ph)=cof(1)+sum_{n=1}^{nh2/2}(cof(2n-1)Cos(ph*n)+cof(2n)*Sin(ph*n))
!!   where ph=2pi*frout*(t-t0) and t0=(max(t)+min(t))/2
!!  Input:
!!   t(:),fin(:),er(:) -times, values and weights of observations;
!!   nh2- number of model parameters: nh2/2-no. of harmonics
!!   frin- initial value of frequency is abs(fr)
!!       for frin<0 on Return frout=frin
!!       for frin>0 value of frout is adjusted by NLSQ iterations
!!  Output:
!!   frout-final value of frequency (see frin)
!!   cof(:),dcof(:)-Fourier coefficients and errors
!!      cof(nh2+2)-approximate t0 (for full accuracy use t0=(max(t)+min(t))/2)
!!      dcof(nh2+2)-frequency error (0 for fr<0)
!!
!!  Copyrights and Distribution:
!!  This package is subject to copyrights by its author,
!!  (C) A. Schwarzenberg-Czerny 1998-2011,  alex@camk.edu.pl
!!  Please quote A.Schwarzenberg-Czerny, 1995, Astr. & Astroph. Suppl, 110,405 ;
!
         Integer, Parameter :: ITMAX = 100 !! max number of Newton-Raphson iterations
         Integer, Intent (In) :: nobs !! number of observations
         Real (TIME), Intent (In) :: t (nobs)!! observation times
         Real (SP), Intent (In) :: valin (nobs)!! observation values
         Real (SP), Intent (In) :: er (nobs)!! observation weights
         Integer, Intent (In) :: nh2 !! number of Fourier coef's (nh2=2*nh+1)
         Real (TIME), Intent (In) :: frin !! initial frequency approximmation
         Real (TIME), Intent (Out) :: frout !! improved frequency
         Real (SP), Intent (Out) :: valout (nobs)!! residuals from fit/prewhitened data
         Real (SP), Intent (Out) :: cof (nh2/2*2+2)!!   Fourier coefficients
         Real (SP), Intent (Out) :: dcof (nh2/2*2+2)!!   coefficient errors
!
         Logical :: finish
         Integer :: l, n, nn2, nh, nx, nx1, it
         Real (SP), Allocatable :: e (:), a (:, :)
         Real (SP) :: cond, detnth, rhs, sig, avstep
         Real (TIME) :: ph, t0, dfr
!
         nh = Max (1, nh2/2)
         nn2 = nh + nh
!
!
         If (nobs < nh+nh+2 .Or. size(valin) /= nobs .Or. size(er) /= nobs &
               .Or. size(valout) /= nobs) Then
            Write (*,*) 'FOUW:error: wrong Size of arrays'
            Return
         End If
!
         t0 = (maxval(t)+minval(t)) * 0.5_TIME
         dfr = 1. / (maxval(t)-minval(t))
!
         cof = 0._SP
         dcof = 0._SP
         frout = Abs (frin)
         Print *, '  it    sigma   step    cond'
         Do it = 1, ITMAX
            nx = nn2 + 2
            nx1 = nx ! adjust only coefficients
            If (it > 1) nx1 = nx + 1 ! adjust also frequency
            Allocate (e(nx1), a(nx1, nx1))
            a = 0._SP
!
            Do l = 1, nobs
               ph = frout * (t(l)-t0)
               ph = PI2 * (ph-floor(ph)) * 0.5_SP
               rhs = valin (l) - cof (1)
               e = 0._SP
               e (1) = 1._SP
!
               Do n = 2, nn2, 2
                  e (n) = real (Cos(ph*n), SP)
                  e (n+1) = real (Sin(ph*n), SP)
                  e (nx) = e (nx) + (-e(n+1)*cof(n)+e(n)*cof(n+1)) * n
                  rhs = rhs - e (n) * cof (n) - e (n+1) * cof (n+1)
               End Do
!
               e (nx) = e (nx) * real (PI2*0.5_SP*(t(l)-t0), SP)
               e (nx1) = rhs
               valout (l) = rhs
               If (er(l) > 0._SP) Then
                  e = e / er (l)
                  Call givensacc (e, a)
               End If
            End Do
!
            finish = (frin < 0. .And. it >= 2)
            If ( .Not. finish) Then
               cond = givens (a, nobs-nx1+1, detnth)
               sig = Sqrt (a(nx1, nx1))
               avstep = Sqrt (sum(a(nx1, :nx1-1)**2)/(nx1-1))
               cof (nx) = 0._SP
               cof (:nx1-1) = cof (:nx1-1) + a (nx1, :nx1-1)
               frout = frout + cof (nx)
               Do l = 1, nx1 - 1
                  dcof (l) = Sqrt (Max(a(l, l), tiny(a)))
               End Do
               finish = Abs (cof(nx)) > dfr.Or.Abs (cof(nx)) < dcof(nx) * 0.01
            End If
            If (Nint(Exp(Nint(2.5*Log(real(it)))/2.5)) == it .Or. finish) &
               Print '(i5,f10.5,2(1pe8.1))', it, sig, avstep, cond
            If (finish) Exit
            Deallocate (e, a)
         End Do
         cof (nx) = real (t0, SP)
!
         Print *, 'frout=', frout
         Do l = 1, nx1 - 1
            Print '(2(f10.5,a))', cof (l), ' +/-', dcof (l)
         End Do
      End Subroutine fouw
!
End Module aovsub