mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
speckleAmpPSD.hpp File Reference

Calculates the PSD of speckle intensity given a modified Fourier mode amplitude PSD. More...

Calculates the PSD of speckle intensity given a modified Fourier mode amplitude PSD.

Author
Jared R. Males (jared.nosp@m.male.nosp@m.s@gma.nosp@m.il.c.nosp@m.om)

Definition in file speckleAmpPSD.hpp.

Go to the source code of this file.

Namespaces

 mx
 The mxlib c++ namespace.
 

Functions

template<typename realT >
int mx::AO::analysis::speckleAmpPSD (std::vector< realT > &spFreq, std::vector< realT > &spPSD, const std::vector< realT > &freq, const std::vector< realT > &fmPSD, const std::vector< std::complex< realT >> &fmXferFxn, const std::vector< realT > &nPSD, const std::vector< std::complex< realT >> &nXferFxn, int N, std::vector< realT > *vars=nullptr, std::vector< realT > *bins=nullptr, bool noPSD=false)
 Calculate the PSD of the speckle intensity given the PSD of Fourier mode amplitude. More...
 
template<typename realT >
int mx::AO::analysis::speckleAmpVarMean (std::vector< realT > &vars, std::vector< realT > &bins, std::vector< realT > &freq, std::vector< realT > &fmPSD, std::vector< std::complex< realT >> &fmXferFxn, std::vector< realT > &nPSD, std::vector< std::complex< realT >> &nXferFxn, int N)
 Calculate the variance of the mean vs bin size in a speckle time-series. More...
 

Function Documentation

◆ speckleAmpPSD()

template<typename realT >
int mx::AO::analysis::speckleAmpPSD ( std::vector< realT > &  spFreq,
std::vector< realT > &  spPSD,
const std::vector< realT > &  freq,
const std::vector< realT > &  fmPSD,
const std::vector< std::complex< realT >> &  fmXferFxn,
const std::vector< realT > &  nPSD,
const std::vector< std::complex< realT >> &  nXferFxn,
int  N,
std::vector< realT > *  vars = nullptr,
std::vector< realT > *  bins = nullptr,
bool  noPSD = false 
)

Calculate the PSD of the speckle intensity given the PSD of Fourier mode amplitude.

Does so by generating a time-series from the PSD using Fourier-domain convolution, and calculating the average periodogram.

A Hann window is used.

Statistics are not currently normalized. You need to normalize to match expected contrast if desired.

Will also calculate the binned-variances in the generated time-series, if the arguments vars and bins are not null pointers.

Todo:

Figure out what's going on with psdFilter normalization!

probably need to not do overlapped averaging in periodogram. Could drop mean-sub then.

Parameters
[out]spFreqThe frequency grid of the output PSD
[out]spPSDThe speckle amplitude PSD corresponding to the freq coordinates. Will be resized.
[in]freqThe Frequency grid of the input PSD
[in]fmPSDThe Fourier mode PSD.
[in]fmXferFxnThe complex error transfer function, as a function of freq
[in]nPSDThe open-loop noise PSD
[in]nXferFxnThe noise transfer function, as a function of freq
[in]NThe number of trials to use in calculating the amplitude PSD.
[out]vars[optional] The binned variances of the time series generated.
[in]bins[optional] The bin sizes to use in calculating vars.
[in]noPSD[optional] if true then the PSD is not actually calculated, only the binned variances

Definition at line 54 of file speckleAmpPSD.hpp.

References mx::sigproc::augment1SidedPSD(), mx::sigproc::psdVar(), mx::sigproc::averagePeriodogram< realT >::size(), and mx::AO::analysis::speckleAmpPSD().

Referenced by mx::AO::analysis::speckleAmpPSD(), and mx::AO::analysis::speckleAmpVarMean().

◆ speckleAmpVarMean()

template<typename realT >
int mx::AO::analysis::speckleAmpVarMean ( std::vector< realT > &  vars,
std::vector< realT > &  bins,
std::vector< realT > &  freq,
std::vector< realT > &  fmPSD,
std::vector< std::complex< realT >> &  fmXferFxn,
std::vector< realT > &  nPSD,
std::vector< std::complex< realT >> &  nXferFxn,
int  N 
)

Calculate the variance of the mean vs bin size in a speckle time-series.

Does so by generating a time-series from the PSD using Fourier-domain convolution, then calculates the binned-variances in the generated time-series.

This is just a wrapper for speckleAmpPSD, with no actual PSD calculation.

Parameters
[out]varsThe binned variances of the time series generated.
[in]binsThe bin sizes to use to calculate the variances.
[in]freqThe Frequency grid of the input PSD
[in]fmPSDThe open-loop Fourier mode PSD.
[in]fmXferFxnThe complex error transfer function, as a function of freq
[in]nPSDThe open-loop noise PSD
[in]nXferFxnThe noise transfer function, as a function of freq
[in]NThe number of trials to use in calculating the amplitude time-series.

Definition at line 273 of file speckleAmpPSD.hpp.

References mx::AO::analysis::speckleAmpPSD(), and mx::AO::analysis::speckleAmpVarMean().

Referenced by mx::AO::analysis::speckleAmpVarMean().