Module aovper !! Purpose: <br> !! State-of-art routines for period search in unevenly sampled time series. !! Package contains routines for phase folding and binning, Fourier series fit !! power spectrum and search of planetary transits, for weighted and unweighted !! observations. !! !! Calling:<br> !! Demo calling example is contained in `AovDrv` routine in this file. !! Because of numerical considerations times of observations need to be shifted !! by subtraction of `EPOCH` time so that they span 0. This is not done within !! ` aovper`routines to enable the user use of one value of EPOCH for !! independent runs.<br> !! All routines come in two flavours: for !! weighted and unweighted observations. Their calling sequences are uniform, !! though depending on routine some parameters are ignored.<br> !! CHAR sptype - determines type of frequency spectrum returned from routine; !! for `RAW`, `AOV`, `MODEL`, `RESID` the spectrum obeys: !! none, F, beta and another beta probability distributions<br> !! INT no - number of observations<br> !! TIMES t(No) - times of observations<br> !! FLOAT v(no) - values of observations<br> !! FLOAT w(no) - weights of observations <br> !! INT npar - input number of parameters of mathematical model in use;<br> !! INT ncov - input number of different phase origins employed to smooth !! effects of time discretization<br> !! SPEC_T spec - on input contains spectrum array `th(no)` and parameters of its !! frequency grid `nfr`, `fr0`, `frs`, where `freq=ifr*frs+fr0, ifr=0,...,nfr-1`; !! on output `th` contains resultant frequency spectrum.<br> !! !! Methods:<br> !! The Aov package contains routines for period search and !! planetary transit search by phase folding and binning of light curves !! and/or by multi-harmonic Fourier series fit. For a particular set up !! these routines reproduce popular box-least squares and generalised Lomb !! Scargle algorithm, ableit in a more efficient implementations involving sub bins !! or orthogonal projections. <br> !! !! The period routines split into two groups: !! (i) `AOV/TR & AOVMH` calculate quadratic norm of a fitted model. !! As such their results obey Parceval theorem. Depending on `SPTYPE` parameter !! the result spectrum may be converted to classical statistics: either aov !! (Fisher F) spectrum or the corresponding beta statistics spectra of normalized !! 'model power' or 'residual chi2' flavour.<br> !! (i.A) `AOV/TR`: Depending on sign of npar: for npar>0 `aov/tr` routine fits !! folded light curve with a histogram of npar bins, for `npar<0` `aov/tr` routine !! fits light curve with just two phase bins: in- and out-of suspect transit. !! The transit width is fixed by npar while its position is fine-tuned !! depending on `ncov`. (Reference: A.Schwarzenberg-Czerny, 1989, MNRAS 241, 153; !! A.Schwarzenberg-Czerny, & J.-Ph. Beaulieu, 2006, MNRAS 365, 165.<br> !! (i.B) `AOVMH` fits Fourier series of npar terms by orthogonal projection on !! Szego polynomials. Reference: A.Schwarzenberg-Czerny, 1996, ApJ 460L, 107. !! In particular for `npar=3` `AOVMH` returns Generalised !! Lomb-Scargle spectrum of Ferraz-Mello, 1981 and Zechmeister & Kurster, 2009.<br> !! (ii) The other two routines, `PWSP` & `PWSPW`, return squared amplitude of a fitted !! sinusoid. The underlying algorithms differ from each other even for unary !! weights. For uneven sampling the result power spectra obey no Parceval theorem !! hence power spectra obey no classical probability distribution. Provided all !! input values are 1 `PWSP` returns window function (`PWSPW` does not).<br> !! Several general utility routines are provided:<br> !! `GRID`, `PEAK` and `PEAKRM` to set up frequency spectrum and analyse its peaks. !! `NORM_TIMES` from module `AovObs` may be used to subtract `EPOCH` from times. <br> !! Copyrights: <br> !! This package is subject to copyrights by its author, !! (C) A. Schwarzenberg-Czerny 1998-2018, alex@camk.edu.pl !! Its distribution is free, except that distribution of modifications is !! prohibited. Please acknowlege use by citation of papers listed in the !! corresponding routines. Use aovconst Use mprms, Only: EPOCH, ERRIN, VERBOSE, CTMIN Use aovsub Use aovobs ! Implicit None !subroutine spectrum(sptype, no, npar, vf, spec) !function aovmh(sptype,t,vin,w,npar,ncov,spec) !function aovmhw(sptype,t,vin,win,npar,ncov,spec) !function aov(sptype,t,vin,w,np,ncov,spec) !function aovw(sptype,t,vin,win,np,ncov,spec) !function powsp(sptype, t,vin,w, npar, ncov, spec) !function powspw(sptype,t,vin,win,npar,ncov,spec) !function aovdrv(method,sptype,t,v,w,npar,ncov,spec) ! Private !public aovdrv,spectrum,,aovmh,aovmhw,aov,aovw,powsp,powspw Public aovdrv ! Contains ! Subroutine spectrum (sptype, no, npar, vf, spec) !! Purpose: <br> !! Converts quadratic norm of fitted model into !! one of frequency spectrum classic Fisher statistics Integer, Intent (In) :: no !! number of observations Integer, Intent (In) :: npar !! number of model parameters Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Real (SP), Intent (In) :: vf !! sum of squares Type (SPEC_T) :: spec !! frequency and spectrum grid ! Real (SP) :: fac Select Case (trim(sptype))! SP_RAW,SP_AOV,SP_MODEL,SP_RESID Case ('RAW')! do nothing Continue Case ('AOV')! convert to F(npar,no-npar;th) fac = (no-npar) / (npar-1._SP) spec%th = spec%th / Max (vf-spec%th, TINY(vf)) * fac Case ('MODEL')! convert to ibeta(npar/2,(no-npar)/2;th) spec%th = spec%th / vf Case ('RESID')! convert to ibeta((no-npar)/2,npar/2;th) spec%th = 1._SP - spec%th / vf Case Default If (VERBOSE) Print *, "Spectrum: wrong SPTYPE" End Select ! End Subroutine spectrum ! ! ========================================================================== ! AOVMH & AOVMHW routines for multiharmonic Fourier spectrum ! Function aovmh (sptype, t, vin, w, npar, ncov, spec) !! Routine for LSQ fit of Fourier harmonics via projection on orthogonal !! polynomialsfor. No weighting of observations ! !!Purpose: !!Respectively for no weights and weights present AOVMH & AOVMHW fit !!unevenlysampled time series with a multiharmonic Fourier series and !!return spectrum for a range of Fourier frequencies. On input NPAR !!sets number of Fourier terrms (NPAR/2 harmonics). For NPAR=3 and !!SPTYPE=SP_MODEL, AOVMH & AOVMHW results are identical to !!Ferraz-Mello(1981) and Zechmeister & Kurster (2009) generalized !!Lomb-Scargle (GLS) spectrum. NCOV is ignored. !! !!Method: !!By Szego recurrence AOVMH & AOVMHW generate orthogonal trigonometric !!polynomials and next by orthogonal projection expand uneven sampled time !!series in terms of these polynomials. Such an algorithm constitutes the !!fastest implementation of Fourier series least squares fit. For NPAR>>3 !!these algorithms are particularly suitable to detect smooth non-sinusoidal !!oscillations involving sharp features (e.g. eclipses and narrow pulses). !!For detection of sine signals in statistical sense (sensitivity) AOVMH & !!AOVMHW never perform worse than either Lomb-Scargle or Power Spectrum. !!For non-sinusoidal signals AOVMH & AOVMHW advantage becomes significant. !!For identical NPAR, AOVMH & AOVMHW are about 10 times slower than more crude !!AOV/AOVW routines. ! !!Copyrights and Distribution: !!This routines are subject to copyrights by its author, !!(C) A. Schwarzenberg-Czerny 1998-2005, alex@camk.edu.pl !!Its distribution is free, except that distribution of modifications !!is prohibited. !!Please quote A.Schwarzenberg-Czerny, 1996, Astrophys. J.,460, L107 !!while using AOVMH & AOVMHW. Integer :: aovmh !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: w (:)!! dummy Integer, Intent (In) :: npar !! number of parameters (2*nh+1) Integer, Intent (In) :: ncov !! dummy Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: l, n, nn2, nh, no Real (TIME) :: ph (size(t)) Real (SP) :: sn, sav, v (size(t)) Complex (CP), Dimension (size(t)) :: cf, p, z, zn Complex (CP) :: al, sc ! aovmh = ncov + size (w)! dummy use aovmh = 1 no = size (t) nh = Max (1, npar/2) nn2 = nh + nh ! ! normalize observations sav = sum (vin) / no v = vin - sav sav = sum (v**2) If (no .Lt. nn2+2) Then If (VERBOSE) Print *, 'AOVMH:error: too few observations' Return End If ! spec%th = 0. Do l = 1, spec%nfr ph = (spec%frs*(l-1)+spec%fr0) * t ! ph = PI2 * (ph-floor(ph)) z = cmplx (Cos(ph), Sin(ph), kind=CP) ph = ph * nh cf = v * cmplx (Cos(ph), Sin(ph), kind=CP) zn = 1._SP p = 1. ! Do n = 0, nn2 sn = sum (Abs(p)**2) al = sum (z*p) ! NOTE: dot_product(a,b)=sum(conjg(a)*b) sc = dot_product (p, cf) ! ! sn = Max (sn, 10._SP*tiny(sn)) sn = Max (sn, epsilon(sn)) al = al / sn spec%th (l) = spec%th(l) + Abs (sc) ** 2 / sn p = p * z - al * zn * conjg (p) zn = zn * z End Do End Do Call spectrum (sptype, no, nn2+1, sav, spec) aovmh = 0 End Function aovmh ! Function aovmhw (sptype, t, vin, win, npar, ncov, spec) !! Routine for LSQ fit of Fourier harmonics via projection on orthogonal !! polynomialsfor. Weights are applied to observations. ! !!Purpose:<br> !!Respectively for no weights and weights present AOVMH & AOVMHW fit !!unevenlysampled time series with a multiharmonic Fourier series and !!return spectrum for a range of Fourier frequencies. On input NPAR !!sets number of Fourier terrms (NPAR/2 harmonics). For NPAR=3 and !!SPTYPE=SP_MODEL, AOVMH & AOVMHW results are identical to !!Ferraz-Mello(1981) and Zechmeister & Kurster (2009) generalized !!Lomb-Scargle (GLS) spectrum. NCOV is ignored. !! !!Method: <br> !!By Szego recurrence AOVMH & AOVMHW generate orthogonal trigonometric !!polynomials and next by orthogonal projection expand uneven sampled time !!series in terms of these polynomials. Such an algorithm constitutes the !!fastest implementation of Fourier series least squares fit. For NPAR>>3 !!these algorithms are particularly suitable to detect smooth non-sinusoidal !!oscillations involving sharp features (e.g. eclipses and narrow pulses). !!For detection of sine signals in statistical sense (sensitivity) AOVMH & !!AOVMHW never perform worse than either Lomb-Scargle or Power Spectrum. !!For non-sinusoidal signals AOVMH & AOVMHW advantage becomes significant. !!For identical NPAR, AOVMH & AOVMHW are about 10 times slower than more crude !!AOV/AOVW routines. ! !!Copyrights and Distribution:<br> !!This routines are subject to copyrights by its author, !!(C) A. Schwarzenberg-Czerny 1998-2005, alex@camk.edu.pl !!Its distribution is free, except that distribution of modifications !!is prohibited. !!Please quote A.Schwarzenberg-Czerny, 1996, Astrophys. J.,460, L107 !!while using AOVMH & AOVMHW. Integer :: aovmhw !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: win (:)!! weights of observations Integer, Intent (In) :: npar !! number of parameters (2*nh+1) Integer, Intent (In) :: ncov !! dummy Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: l, n, nn2, nh, no ! Real (SP) :: sn, sav, rw (size(t)), vr (size(t)) Real (TIME) :: ph (size(t)) Complex (CP), Dimension (size(t)) :: cf, p, z, zn Complex (CP) :: al, sc ! aovmhw = ncov ! dummy use aovmhw = 1 nh = Max (1, npar/2) nn2 = nh + nh ! ! normalize observations rw = 0._SP Where (win > 0._SP) rw = Sqrt (win) no = count (rw > 0._SP) vr = vin * rw sav = sum (vr*rw) / sum (rw**2) vr = vr - sav * rw sav = sum (vr**2) If (no .Lt. nn2+2) Then If (VERBOSE) Print *, 'AOVMHW:error: too few observations' Return End If spec%th = 0. Do l = 1, spec%nfr ph = (spec%frs*(l-1)+spec%fr0) * t ! ! (f,g)=sum(w*conjg(f)*g) -definition ph = PI2 * (ph-floor(ph)) z = cmplx (Cos(ph), Sin(ph), kind=CP) ph = ph * nh cf = vr * cmplx (Cos(ph), Sin(ph), kind=CP) zn = 1._SP p = rw ! Do n = 0, nn2 sn = sum (Abs(p)**2) al = sum (rw*z*p) ! NOTE: dot_product(a,b)=sum(conjg(a)*b) sc = dot_product (p, cf) ! sn = Max (sn, epsilon(sn)) al = al / sn spec%th (l) = spec%th(l) + Abs (sc) ** 2 / sn p = p * z - al * zn * conjg (p) zn = zn * z End Do End Do Call spectrum (sptype, no, nn2+1, sav, spec) aovmhw = 0 End Function aovmhw ! End of AOVMH & AOVMHW routines for multiharmonic Fourier spectrum ! ========================================================================== ! Function aov (sptype, t, vin, w, np, ncov, spec) !! Routine for phase folding and binning of observations using no weights<br> ! !! Purpose:<br> !! The state-of-art routines for either period search by phase folding and !! binning of observations (for NP>0) or for search for periodic transits/pulses !! of short duty cycle, with constat signal elsewhere (for NP<0). In the latter !! case data are fitted with just two unequal phase bins: in- and out-of transit. ! !! Method:<br> !! In any case observations are phase folded and binned into NP equal phase !! bins starting at NCOV sub-bin phases. Subsequently the algorithm !! either sums squares of bin averages (for NP>0). For transit search !! (NP<0) algorithm searches for the maximum value phase bin (considered !! in-transit for magnitude units). Its squared value and identities stemming !! from Parceval theorem are employed to evaluate sum of squares for two !! phase-bins in- and out-of transit. Note: here transits are defined as MAXIMA !! (i.e. in magnitude units). For flux units (minima) reverse sign of observations. !! If algorithm detects poor phase coverage (bins with <CTMIN observations) !! it rebins observations into bins of near equal occupancy, by sorting !! observations in phase. This is a rare situation hence amounts to no !! significant slowing down. Thanks to this feature, for null frequency !! these routines return useful variability indicator as long as the sort !! routine is stable (e.g. merge sort). !! !! Advantage:<br> !! By checking and repairing phase coverage these routines perform robust and !! near-optimal period search, returnind spectrum of known classical !! probability distribution, equally well for small and large data. Slight !! statistical inefficiency stems from step-like function yielding slightly !! worse fit that equal size Fourier fit. This may be remedied by calling !! AOVMH/AOVMHV instead, at the speed cost. In terms of speed the present !! AOV/AOVW implementation is state-of-art as it employs pre-binning into !! NP*NCOV sub-bins and only afterwords combines subbins into suitably !! phased bins. See B. Enoch et al., 2012, A&A 548, A48 !! for independent evaluation of the present sub-bin algorithm. ! !! Copyrights and Distribution:<br> !! This package is subject to copyrights by its author, !! (C) A. Schwarzenberg-Czerny 2003-2018, alex@camk.edu.pl !! Its distribution is free, except that distribution of modifications !! is prohibited. !! Please quote A. Schwarzenberg-Czerny & J.-Ph. Beaulieu, 2005, MNRAS 365, 165, !! while using AOV or AOVW. Integer :: aov !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: w (:)!! dummy Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Integer, Intent (In) :: np !! number of phase bins Integer, Intent (In) :: ncov !! number of phase coverages Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: nct ((Abs(np)+1)*ncov), ind (size(t)), i, ibin, ip, & nbc, ifr, iflex, npar, no Real (SP) :: ph (size(t)), v (size(t)), ave((Abs(np)+1)*ncov), sav, v2 Real (TIME) :: fr, dph Logical :: transit ! aov = size (w)!! dummy use aov = 1 no = size (t) npar = Abs (np) transit = np < 0 ! ! normalize observations sav = sum (vin) / no v = vin - sav v2 = sum (v**2) If (no <= npar+npar) Then If (VERBOSE) Print *, 'AOV: error: too few observations' Return End If ! nbc = npar * ncov ! Loop over all frequencies iflex = 0 Do ifr = 1, spec%nfr ! Loop over frequencies fr = real (ifr-1, kind=TIME) * spec%frs + spec%fr0 Do ip = 0, 1 ave = 0._SP nct = 0 If (ip == 0) Then ! Try default fixed bins ... Do i = 1, no ! MOST LABOR HERE dph = t (i) * fr ! real(TIME) :: dph, t, fr dph = dph - floor (dph) ph (i) = real (dph, SP) ibin = Int (dph*nbc) + 1 ave (ibin) = ave (ibin) + v (i) nct (ibin) = nct (ibin) + 1 End Do Else !... elastic bins, if necesseary iflex = iflex + 1 ! sort index ind using key ph Call sortx (ph, ind) ! NRf77 indexx would do Do i = 1, no ibin = i * nbc / no + 1 ave (ibin) = ave (ibin) + v (ind(i)) nct (ibin) = nct (ibin) + 1 End Do End If ! counts: sub-bins=>bins nct (1+nbc:ncov+nbc) = nct (1:ncov) i = 0 Do ibin = ncov + nbc, 1, - 1 i = i + nct (ibin) nct (ibin) = i End Do Do ibin = 1, nbc nct (ibin) = nct (ibin) - nct (ibin+ncov) End Do If (minval(nct(1:nbc)) >= CTMIN) Exit End Do ! data: sub-bins=>bins */ ave (1+nbc:ncov+nbc) = ave (1:ncov) sav = 0._SP Do ibin = ncov + nbc, 1, - 1 sav = sav + ave (ibin) ave (ibin) = sav End Do Do ibin = 1, nbc ave (ibin) = ave (ibin) - ave (ibin+ncov) End Do ! If (transit) Then ! statistics for transits ... ave (1:nbc) = ave (1:nbc) / nct (1:nbc) ibin = maxloc (ave(1:nbc), dim=1) sav = ave (ibin) sav = sav * sav * nct (ibin) * no / real (no-nct(ibin), kind=SP) Else ! the same for the classical AoV statistics: sav = sum (ave(1:nbc)**2/nct(1:nbc)) / ncov End If spec%th (ifr) = sav End Do i = npar If (transit) i = 2 ! If (VERBOSE .And. iflex > 0) Print *, & 'AOV:warning: poor phase coverage at ', iflex, ' frequencies' Call spectrum (sptype, no, i, v2, spec) aov = 0 End Function aov ! Function aovw (sptype, t, vin, win, np, ncov, spec) !! Routine for phase folding and binning of observations using weights ! !! Purpose:<br> !! The state-of-art routines for either period search by phase folding and !! binning of observations (for NP>0) or for search for periodic transits/pulses !! of short duty cycle, with constat signal elsewhere (for NP<0). In the latter !! case data are fitted with just two unequal phase bins: in- and out-of transit. ! !! Method:<br> !! In any case observations are phase folded and binned into NP equal phase !! bins starting at NCOV sub-bin phases. Subsequently the algorithm !! either sums squares of bin averages (for NP>0). For transit search !! (NP<0) algorithm searches for the maximum value phase bin (considered !! in-transit for magnitude units). Its squared value and identities stemming !! from Parceval theorem are employed to evaluate sum of squares for two !! phase-bins in- and out-of transit. Note: here transits are defined as MAXIMA !! (i.e. in magnitude units). For flux units (minima) reverse sign of observations. !! If algorithm detects poor phase coverage (bins with <CTMIN observations) !! it rebins observations into bins of near equal occupancy, by sorting !! observations in phase. This is a rare situation hence amounts to no !! significant slowing down. Thanks to this feature, for null frequency !! these routines return useful variability indicator as long as the sort !! routine is stable (e.g. merge sort). !! !! Advantage:<br> !! By checking and repairing phase coverage these routines perform robust and !! near-optimal period search, returnind spectrum of known classical !! probability distribution, equally well for small and large data. Slight !! statistical inefficiency stems from step-like function yielding slightly !! worse fit that equal size Fourier fit. This may be remedied by calling !! AOVMH/AOVMHV instead, at the speed cost. In terms of speed the present !! AOV/AOVW implementation is state-of-art as it employs pre-binning into !! NP*NCOV sub-bins and only afterwords combines subbins into suitably !! phased bins. See B. Enoch et al., 2012, A&A 548, A48 !! for independent evaluation of the present sub-bin algorithm. ! !! Copyrights and Distribution:<br> !! This package is subject to copyrights by its author, !! (C) A. Schwarzenberg-Czerny 2003-2018, alex@camk.edu.pl !! Its distribution is free, except that distribution of modifications !! is prohibited. !! Please quote A. Schwarzenberg-Czerny & J.-Ph. Beaulieu, 2005, MNRAS 365, 165, !! while using AOV or AOVW. Integer :: aovw !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: win (:)!! weights of observations Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Integer, Intent (In) :: np !! number of phase bins Integer, Intent (In) :: ncov !! number of phase coverages Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: ind (size(t)), i, ibin, ip, nbc, ifr, iflex, npar, & no, noeff Real (SP) :: ph (size(t)), v (size(t)), w (size(t)), sav, ctw, sw, & nct((Abs(np)+1)*ncov), ave ((Abs(np)+1)*ncov), v2 Real (TIME) :: fr, dph Logical :: transit ! aovw = 1 no = size (t) npar = Abs (np) transit = np < 0 ! ! normalize observations w = 0._SP Where (win > 0._SP) w = win noeff = count (w > 0.) sw = sum (w) v = vin - sum (vin*w) / sw v2 = sum (w*v**2) v = v * w If (noeff <= npar+npar .Or. noeff < CTMIN*npar) Then If (VERBOSE) Print *, 'AOVw: error: too few observations' Return End If ! nbc = npar * ncov iflex = 0 ctw = CTMIN * sw / noeff Do ifr = 1, spec%nfr ! Loop over frequencies fr = real (ifr-1, kind=TIME) * spec%frs + spec%fr0 Do ip = 0, 1 ave = 0._SP nct = 0._SP If (ip == 0) Then ! Try default fixed bins ... Do i = 1, no ! MOST LABOR HERE dph = t (i) * fr ! real(TIME) :: dph, t, fr dph = dph - floor (dph) ph (i) = real (dph, SP) ! note: mod accounts for rounding up due to double-->real conversion ibin = Int (dph*nbc) + 1 ave (ibin) = ave (ibin) + v (i) nct (ibin) = nct (ibin) + w (i) End Do Else !... elastic bins, if necesseary ! iflex = iflex + 1 ! sort index ind using key ph Call sortx (ph, ind)! NRf77 indexx would do Do i = 1, no ibin = i * nbc / no + 1 ave (ibin) = ave (ibin) + v (ind(i)) nct (ibin) = nct (ibin) + w (ind(i)) End Do End If ! counts: sub-bins=>bins nct (1+nbc:ncov+nbc) = nct (1:ncov) sav = 0._SP Do ibin = ncov + nbc, 1, - 1 sav = sav + nct (ibin) nct (ibin) = sav End Do Do ibin = 1, nbc nct (ibin) = nct (ibin) - nct (ibin+ncov) End Do If (minval(nct(1:nbc)) >= ctw) Exit End Do ! data: sub-bins=>bins */ ave (1+nbc:ncov+nbc) = ave (1:ncov) sav = 0._SP Do ibin = ncov + nbc, 1, - 1 sav = sav + ave (ibin) ave (ibin) = sav End Do Do ibin = 1, nbc ave (ibin) = ave (ibin) - ave (ibin+ncov) End Do ! If (transit) Then ! AoV statistics for transits ... ave (1:nbc) = ave (1:nbc) / nct (1:nbc) ibin = maxloc (ave(1:nbc), dim=1) sav = ave (ibin) sav = sav * sav * nct (ibin) * sw / (sw-nct(ibin)) Else ! the same for the classical AoV statistics: sav = sum (ave(1:nbc)**2/nct(1:nbc)) / ncov End If spec%th (ifr) = sav End Do i = npar If (transit) i = 2 ! Call spectrum (sptype, noeff, i, v2, spec) ! If (VERBOSE .And. iflex > 0) Print *, & 'AOVw:warning: poor phase coverage at ', iflex, ' frequencies' aovw = 0 End Function aovw ! End of AOV & AOVW routines !========================================================================== ! POWSP & POWSPW routines returning two different versions of power spectrum ! Function powsp (sptype, t, vin, w, npar, ncov, spec) ! !!Purpose:<br> !!POWSP - returns discrete power spectrum (DPS) for unevenly sampled !!observations. For obs.v==1 calculates window function (WF). !!SPTYPE, NPAR and NCOV are ignored. The normalization factor is 1/NOBS**2, !!suitable for WF. Power is defined in a usual way as sum of squares of !!projections onto sine and cosine functions. This concept does not extend !!on weighted observations, hence POWSPW employs different algorithm.<br> !!Reference: T.Deeming, 1975, Ap&SS 36, 137. !! !! Copyrights and Distribution:<br> !! This routine is subject to copyrights by its author, !! (C) A. Schwarzenberg-Czerny 1998-2018, alex@camk.edu.pl !! Its distribution is free, except that distribution of modifications is !! prohibited. Integer :: powsp !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: w (:)!! dummy Integer, Intent (In) :: npar !! dummy Integer, Intent (In) :: ncov !! dummy Character (Len=*), Intent (In) :: sptype !! dummy Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: l, no Real (SP) :: fac, v (size(t)) Real (TIME) :: ph (size(t)) ! powsp = ncov + npar + len (sptype) + size (w)! dummy use powsp = 1 no = size (t) ! normalize observations v = vin - sum (vin) / no If (no .Lt. 4) Then If (VERBOSE) Print *, 'POWSP:error: wrong size of arrays' Return End If ! fac = 1._SP / no ** 2 Do l = 1, spec%nfr ph = (spec%frs*(l-1)+spec%fr0) * t ! (f,g)=sum(w*conjg(f)*g) -definition ph = PI2 * (ph-floor(ph)) spec%th (l) = Abs(sum(v*cmplx(Cos(ph), Sin(ph), kind=CP)))**2 * fac End Do powsp = 0 End Function powsp ! Function powspw (sptype, t, vin, win, npar, ncov, spec) ! !!Purpose: !!POWSPW - returns discrete power spectrum (DPS) for weighted unevenly sampled !!observations defined differently than in POWSP. No window function (WF) is !!available and NPAR and NCOV are ignored. The algorithm by weighted least !!squares fits data with Asin(om*t)+Bcos(om*t) and for each frequency returns !!A**2+B**2. The normalization is peculiar to this implementation and depends !!on time sampling. !!Reference: S. T. Fletcher et al., 2011, MNRAS 415, 1310 (see Appendix). !! !! Copyrights and Distribution: !! This routine is subject to copyrights by its author, !! (C) A. Schwarzenberg-Czerny 2018, alex@camk.edu.pl !! Its distribution is free, except that distribution of modifications is !! prohibited. Integer :: powspw !! status Real (TIME), Intent (In) :: t (:)!! times of observations Real (SP), Intent (In) :: vin (:)!! values of observations Real (SP), Intent (In) :: win (:)!! weights of observations Integer, Intent (In) :: npar !! dummy Integer, Intent (In) :: ncov !! dummy Character (Len=*), Intent (In) :: sptype !! dummy Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: l, no Real (SP), Dimension (size(t)) :: c, s, v, w Real (SP) :: sx, cx, ss, cc, sc Real (TIME) :: dph (size(t)) ! powspw = Int (ncov+npar+len(sptype))! dummy use powspw = 1 ! normalize observations w = 0._SP Where (win > 0._SP) w = win no = count (w > 0._SP) v = (vin-sum(vin*w)/sum(w)) * w If (no .Lt. 4) Then If (VERBOSE) Print *, 'POWSPW:error: wrong size of arrays' Return End If ! Do l = 1, spec%nfr dph = (spec%frs*(l-1)+spec%fr0) * t s = real (PI2*(dph-floor(dph)), SP) c = Cos (s) s = Sin (s) sx = sum (s*v) cx = sum (c*v) ss = sum (w*s**2) cc = sum (w*c**2) sc = sum (w*c*s) spec%th (l) = 0.25_SP * ((sx*cc-cx*sc)**2+(cx*ss-sx*sc)**2) / & Max ((ss*cc-sc*sc)**2, TINY(ss)) End Do powspw = 0 End Function powspw ! End of POWSP & POWSPW routines !======================================================== Function aovdrv (method, sptype, t, v, w, npar, ncov, spec) !! AOVDRV driver routine for `AovPer` frequency spectrum routines !! Purpose: example calling of spectrum routines. Integer :: aovdrv !! status Character (Len=*), Intent (In) :: method !! periodogram method Character (Len=*), Intent (In) :: sptype !! type of returned spectrum Real (TIME), Intent (In) :: t (:)!! times Real (SP), Intent (In) :: v (:)!! values Real (SP), Intent (In) :: w (:)!! weights or errors, depending on ERRIN Integer, Intent (In) :: npar !! number of parameters Integer, Intent (In) :: ncov !! number of phase coverages Type (SPEC_T) :: spec !! frequency grid and result spectrum ! Integer :: i ! get spectrum/periodogram Select Case (trim(method)) Case ('AOV') i = aov (sptype, t, v, w, Max(2, npar), ncov, spec) Case ('ATR') i = aov (sptype, t, v, w,-Max(2, npar), ncov, spec) Case ('AMH') i = aovmh (sptype, t, v, w, Max(2, npar), ncov, spec) Case ('PWS') i = powsp (sptype, t, v, w, Max(2, npar), ncov, spec) Case ('AOVW') i = aovw (sptype, t, v, w, Max(2, npar), ncov, spec) Case ('ATRW') i = aovw (sptype, t, v, w,-Max(2, npar), ncov, spec) Case ('AMHW') i = aovmhw (sptype, t, v, w, Max(2, npar), ncov, spec) Case ('PWSW') i = powspw (sptype, t, v, w, Max(2, npar), ncov, spec) Case Default If (VERBOSE) Print *, 'aovdrv: error: Wrong mode' Stop End Select aovdrv = i End Function aovdrv ! End Module aovper