API Reference
Power spectral densities
Models and functions
Pioran.approx
— Functionapprox(psd_model, f_min, f_max, n_components=20, norm=1.0,S_low::Real=20., S_high::Real=20. ; is_integrated_power::Bool = true, basis_function="SHO")
Approximate the PSD with a sum of basis functions to form a covariance function. The PSD model is approximated between f0=f_min/S_low
and fM=f_min*S_high
. By default it is normalised by its integral from f_min
to f_max
but it can also be normalised by its integral from 0 to infinity using the variance
argument.
Arguments
psd_model::PowerSpectralDensity
: model of the PSDf_min::Real
: the minimum frequency in the time seriesf_max::Real
: the maximum frequency in the time seriesn_components::Integer=20
: the number of basis functions to usenorm::Real=1.0
: normalisation of the PSD.S_low::Real=20.0
: scaling factor for the lowest frequency in the approximation.S_high::Real=20.0
: scaling factor for the highest frequency in the approximation.is_integrated_power::Bool = true
: if the norm corresponds to integral of the PSD betweenf_min
andf_max
, if not it is the variance of the process, integral of the PSD from 0 to +inf.basis_function::String="SHO"
: the basis function to use, either "SHO" or "DRWCelerite"
Return
covariance::SumOfSemiSeparable
: the covariance function
Example
using Pioran
𝓟 = SingleBendingPowerLaw(1.0, 1.0, 2.0)
𝓡 = approx(𝓟, 1e-4, 1e-1, 30, 2.31,basis_function="SHO")
Pioran.get_covariance_from_psd
— Method get_covariance_from_psd(psd_features)
Get the covariance function from the PSD features
Tonari.DoubleBendingPowerLaw
— Type DoubleBendingPowerLaw(A, α₁, f₁, α₂, f₂, α₃)
Double bending power law model for the power spectral density
A
: the amplitudeα₁
: the first power law indexf₁
: the first bend frequencyα₂
: the second power law indexf₂
: the second bend frequencyα₃
: the third power law index
\[\mathcal{P}(f) = A\frac{(f/f₁)^{-α₁}}{1 + (f / f₁)^{α₂ - α₁}}\frac{1}{1 + (f / f₂)^{α₃ - α₂}}\]
Tonari.Lorentzian
— TypeLorentzian(A, γ, f₀)
Lorentzian model for the power spectral density
A
: the amplitudeγ
: the width of the peakf₀
: the central frequency
\[\mathcal{P}(f) = \frac{A}{4\pi^2 (f - f₀)^2 + γ^2}\]
Tonari.PowerLaw
— Type PowerLaw(α)
Power law model for the power spectral density
α
: the power law index
\[\mathcal{P}(f) = A f^{-α}\]
Tonari.QPO
— TypeQPO(S₀, f₀,A Q)
QPO model
S₀
: the amplitude at the peakf₀
: the central frequencyQ
: quality factor
\[\mathcal{P}(f) = \frac{S_0 {f_0}^4 } {\left(f^ 2 -{f_0}^2\right)^ 2 + f^2 {f_0}^2 / Q^2 }\]
Tonari.SingleBendingPowerLaw
— Type SingleBendingPowerLaw(A, α₁, f₁, α₂)
Single bending power law model for the power spectral density
A
: the amplitudeα₁
: the first power law indexf₁
: the first bend frequencyα₂
: the second power law index
\[\mathcal{P}(f) = A \frac{(f/f₁)^{-α₁}}{1 + (f / f₁)^{α₂ - α₁}}\]
Helper functions
Pioran.approximated_psd
— Method approximated_psd(f, psd_model, f0, fM; n_components=20, norm=1.0, basis_function="SHO")
Return the approximated PSD. This is essentially to check that the model and the approximation are consistent.
Arguments
f::AbstractVector{<:Real}
: the frequencies at which to evaluate the PSDpsd_model::PowerSpectralDensity
: model of the PSDf0::Real
: the lowest frequencyfM::Real
: the highest frequencyn_components::Integer=20
: the number of basis functions to usenorm::Real=1.0
: normalisation of the PSDbasis_function::String="SHO"
: the basis function to use, either "SHO" or "DRWCelerite"individual::Bool=false
: return the individual components
Pioran.build_approx
— Method build_approx(J, f0, fM; basis_function="SHO")
Prepare the approximation of a PSD
Arguments
J::Integer
: the number of basis functionsf0::Real
: the lowest frequencyfM::Real
: the highest frequencybasis_function::String="SHO"
: the basis function to use, either "SHO" or "DRWCelerite"
Return
spectral_points::Vector{Real}
: the spectral pointsspectral_matrix::Matrix{Real}
: the spectral matrix
Pioran.convert_feature
— Method convert_feature(psd_feature)
Convert a PSD feature to a Celerite covariance function Only QPO is implemented
Arguments
psd_feature::PowerSpectralDensity
: the PSD feature
Return
covariance::SemiSeparable
: the covariance function
Pioran.get_approx_coefficients
— Method get_approx_coefficients(psd_model, f0, fM; n_components=20, basis_function="SHO")
Get the coefficients of the approximated PSD
Arguments
psd_model::PowerSpectralDensity
: model of the PSDf0::Real
: the lowest frequencyfM::Real
: the highest frequencyn_components::Integer=20
: the number of basis functions to usebasis_function::String="SHO"
: the basis function to use, either "SHO" or "DRWCelerite"
Return
amplitudes::Vector{Real}
: the amplitudes of the basis functions
Pioran.get_norm_psd
— Functionget_norm_psd(amplitudes,spectral_points,f_min,f_max,basis_function,is_integrated_power, cov_features=nothing)
Get the normalisation of the sum of basis functions.
Arguments
amplitudes
: amplitude of the basis functionspectral_points
: spectral points of the basis functionf_min::Real
: the minimum frequency in the time seriesf_max::Real
: the maximum frequency in the time seriesbasis_function::String="SHO"
: the basis function to use, either "SHO" or "DRWCelerite"is_integrated_power::Bool=true
: if the norm corresponds to integral of the PSD betweenf_min
andf_max
or if it is the integral from 0 to infinity.cov_features
: PSD features of the PSD, is nothing if there are no features.
Pioran.get_normalised_psd
— Method get_normalised_psd(psd_model::PowerSpectralDensity, spectral_points::AbstractVector{<:Real})
Get the PSD normalised at the lowest frequency
Pioran.integral_celerite
— Methodintegral_celerite(a, b, c, d, x)
Computes the integral of the Celerite power spectrum:
Pioran.integral_drwcelerite
— Methodintegral_drwcelerite(a, c, x)
Computes the integral of the DRWCelerite basis function of amplitude a
and width c
for a given x
.
The DRWCelerite basis function is defined as:
\[ \int \dfrac{a\, {d}x}{(x/c)^6+1} =\dfrac{ac}{3} \left[ \arctan{(x/c)} +\dfrac{\sqrt3}{4}\ln{\left(\dfrac{x^2+xc\sqrt3+c^2}{x^2-xc\sqrt3+c^2}\right)}+\dfrac{1}{2}\arctan{\left(\dfrac{x^2-c^2}{xc}\right)}\right]\]
Pioran.integral_sho
— Methodintegral_sho(a, c, x)
Computes the integral of the SHO basis function of amplitude a and width c for a given x.
This integral is obtained using Equation: 4.2.7.1.3 from the "Handbook of Mathematical Formulas and Integrals" 2009
\[ \int \dfrac{a\, {d}x}{(x/c)^4+1} =\dfrac{ac}{4\sqrt2} \left[\ln{\left(\dfrac{x^2+cx\sqrt2+c^2}{x^2-cx\sqrt2+c^2}\right)}+2\arctan{\left(\dfrac{cx\sqrt2}{c^2-x^2}\right)}\right]\]
Pioran.integrate_basis_function
— Methodintegrate_basis_function(a,c,x₁,x₂,basis_function)
Computes the integral of the basis function between x₁ and x₂ for a given amplitude a and width c.
Pioran.integrate_psd_feature
— Methodintegrate_psd_feature(a, b, c, d, x₁, x₂)
Computes the integral of a celerite power spectral density with coefficients (a,b,c,d) between x₁ and x₂.
Pioran.psd_decomp
— Method psd_decomp(psd_normalised, spectral_matrix)
Get amplitudes of the basis functions by solving the linear system of the approximation
Tonari.separate_psd
— Method separate_psd(psd::PowerSpectralDensity)
Separate the PSD into its BendingPowerLaw components and other components if it is a sum of PSDs
Arguments
psd::PowerSpectralDensity
: power spectral density or sum of PowerSpectralDensity objects
Return
psd_continuum::Union{SumOfPowerSpectralDensity,PowerSpectralDensity,nothing}
: continuum part of the psdpsd_line::Union{PowerSpectralDensity,nothing,Vector{PowerSpectralDensity}}
: non-continuum part of the psd
Covariance functions
Models and functions
Pioran.Exp
— TypeExp(A, α)
Exponential covariance Function
A
: the amplitude of the covariance functionα
: the decay rate of the covariance function
\[k(τ) = A \exp(-α τ)\]
Example
Exp(1.0, 0.25)
Pioran.evaluate
— Methodevaluate(R::Exp, f)
This is the right formula but it disagrees with the Celerite implementation...
Evaluate the power spectral density at frequency f
Pioran.SHO
— Type SHO(A, ω₀, Q)
Simple Harmonic Oscillator covariance Function
A
: the amplitude of the covariance functionω₀
: the angular frequency of the simple harmonic oscillatorQ
: the quality factor of the simple harmonic oscillator
\[k(τ) = A \exp(-ω₀ τ / Q / 2) \left\{\begin{matrix} 2(1 + ω₀ τ) & Q = 1/2 \\ \cos(η ω₀ τ) + \frac{\sin(η ω₀ τ)}{2η Q} & Q < 1/2 \\ \cosh(η ω₀ τ) + \frac{\sinh(η ω₀ τ)}{2η Q} & Q \geq 1/2 \end{matrix}\right.\\ η = \sqrt{\left|1 - \frac{1}{4 Q^2}\right|}\]
See Foreman-Mackey et al. (2017) for more details.
Pioran.Celerite
— Type Celerite(a, b, c, d)
Celerite covariance Function
a
: the amplitude of the first termb
: the amplitude of the second termc
: the decay rate of the covariance functiond
: theperiod
of the covariance function
\[k(τ) = \exp(-c τ) (a \cos(d τ) + b \sin(d τ))\]
See Foreman-Mackey et al. (2017) for more details.
Pioran.evaluate
— Methodevaluate(f, C::Celerite)
evaluate the power spectral density at frequency f
Helper functions
Base.:+
— Method +(::SemiSeparable, ::SemiSeparable)
Sum of semi-separable covariance functions
Pioran.celerite_coefs
— Method celerite_coefs(covariance)
Get the celerite coefficients
Gaussian processes
Pioran.ScalableGP
— TypeA scalable Gaussian process has a covariance function formed of semi-separable kernels
Example
using Pioran
𝓟 = SingleBendingPowerLaw(1.0, 1.0, 2.0)
𝓡 = approx(𝓟, 1e-4, 1e-1, 30, 2.31,basis_function="SHO")
μ = 1.2
f = ScalableGP(𝓡) # zero-mean GP
f = ScalableGP(μ, 𝓡) # with mean μ
See Foreman-Mackey et al. (2017) for more details.
Base.rand
— Functionrand(rng::AbstractRNG, fp::PosteriorGP, N::Int=1)
rand(rng::AbstractRNG, fp::PosteriorGP, τ::AbstractVecOrMat{<:Real}, N::Int=1)
rand(fp::PosteriorGP, N::Int=1)
Sample N
realisations from the posterior GP fp
at the points τ
.
Base.rand
— Methodrand(rng::AbstractRNG, f::ScalableGP)
rand(rng::AbstractRNG, f::ScalableGP, t::AbstractVecOrMat{<:Real})
rand(f::ScalableGP)
rand(f::ScalableGP, t::AbstractVecOrMat{<:Real})
Draw a realisation from the GP f
at the points t
.
Distributions.logpdf
— Methodlogpdf(f::ScalableGP, Y::AbstractVecOrMat{<:Real})
Compute the log-likelihood of the data Y given the GP f.
Pioran._predict_cov
— Method_predict_cov(fp::PosteriorGP, τ::AbstractVecOrMat{<:Real})
Compute the posterior covariance of the GP at the points τ.
Pioran._predict_mean
— Method_predict_mean(fp::PosteriorGP, τ::AbstractVecOrMat{<:Real})
Compute the Posterior mean of the GP at the points τ.
Pioran.posterior
— Methodposterior(f::ScalableGP, y::AbstractVecOrMat{<:Real})
Compute the posterior Gaussian process fp
given the GP f
and the data y
.
Statistics.cov
— Methodcov(fp::PosteriorGP, τ::AbstractVecOrMat{<:Real})
cov(fp::PosteriorGP)
Compute the covariance of the posterior GP at the points τ.
Statistics.mean
— Methodmean(fp::PosteriorGP, τ::AbstractVecOrMat{<:Real})
mean(fp::PosteriorGP)
Compute the mean of the posterior GP at the points τ.
Statistics.std
— Methodstd(fp::PosteriorGP, τ::AbstractVecOrMat{<:Real})
std(fp::PosteriorGP)
Compute the standard deviation of the posterior GP at the points τ.
Solvers
Celerite solver
Pioran.compute_nll
— Methodcompute_nll(t, y, σ², a, b, c, d)
Compute the likelihood using the vectorised implementation in `get_values!`.
Still experimental.
Pioran.get_values!
— Methodget_values!(a, b, c, d, zp, U, V, P, D, t)
Compute the values of the matrices and vectors needed for the celerite algorithm. This is a vectorised version of the init_semi_separable!
and solve_prec!
functions. This function appears to be faster than the two previous functions when J > 16 but it also uses more memory.
More study of this implementation is needed.
Pioran.log_likelihood
— Methodlog_likelihood(cov, τ, y, σ2)
Compute the log-likelihood of a semi-separable covariance function using the celerite algorithm.
Arguments
cov::SumOfSemiSeparable
orcov::CARMA
orcov::SemiSeparable
: the covariance functionτ::Vector
: the time pointsy::Vector
: the dataσ2::Vector
: the measurement variances
Pioran.logl
— Methodlogl(a, b, c, d, τ, y, σ2)
Compute the log-likelihood of a GP with a semi-separable covariance function using the celerite algorithm.
Arguments
a::Vector
b::Vector
c::Vector
d::Vector
τ::Vector
: the time pointsy::Vector
: the dataσ2::Vector
: the measurement variances
See Foreman-Mackey et al. (2017) for more details.
Pioran.predict
— Methodpredict(cov, τ, t, y, σ2)
Compute the posterior mean of the GP at the points τ given the data y and time t.
Arguments
cov::SumOfSemiSeparable
orcov::CARMA
orcov::SemiSeparable
: the covariance functionτ::Vector
: the time pointst::Vector
: the data time pointsy::Vector
: the dataσ2::Vector
: the measurement variances
Pioran.simulate
— Methodsimulate(rng, cov, τ, σ2)
simulate(cov, τ, σ2)
Draw a realisation from the GP with the covariance function cov at the points τ with the variances σ2.
Arguments
rng::AbstractRNG
: the random number generatorcov::SumOfSemiSeparable
orcov::CARMA
orcov::SemiSeparable
: the covariance functionτ::Vector
: the time pointsσ2::Vector
: the measurement variances
Pioran.solve_prec!
— Methodsolve_prec!(z, y, U, W, D, ϕ)
Forward and backward substitution of the celerite algorithm.
See Foreman-Mackey et al. (2017) for more details.
Direct solver
Pioran.log_likelihood_direct
— Methodlog_likelihood_direct(cov::KernelFunctions.SimpleKernel, t::Vector, y::Vector, σ²::Vector)
Compute the log-likelihood of the data Y given the GP f with the direct solver.
Pioran.predict_cov
— Methodpredict_direct(cov::KernelFunctions.SimpleKernel, τ::AbstractVector, t::AbstractVector, σ²::AbstractVector)
Compute the posterior covariance of the GP at the points τ given the times t and the noise variance σ².
Pioran.predict_direct
— Functionpredict_direct(cov::KernelFunctions.SimpleKernel, τ::AbstractVector, t::AbstractVector, y::AbstractVector, σ²::AbstractVector, with_covariance::Bool=false)
Compute the posterior mean of the GP at the points τ given the data (t, y) and the noise variance σ².
Plotting
Diagnostics
Pioran.run_diagnostics
— Functionrun_diagnostics(prior_samples, norm_samples, f_min, f_max, model, S_low=20.0, S_high=20.0; path="", basis_function="SHO", n_components=20)
Run the prior predictive checks for the model and the approximation of the PSD
Arguments
prior_samples::Array{Float64, 2}
: Model samples from the prior distributionnorm_samples::Array{Float64, 1}
: The samples of the normalisation of the PSDf_min::Float64
: The minimum frequency of the observed dataf_max::Float64
: The maximum frequency of the observed datamodel::Function
: The modelS_low::Float64=20.0
: the scaling factor for the appoximation at low frequenciesS_high::Float64=20.0
: the scaling factor for the appoximation at high frequenciespath::String=""
: The path to save the plotsbasis_function::String="SHO"
: The basis function for the approximation of the PSDn_components::Int=20
: The number of components to use for the approximation of the PSD
Pioran.run_posterior_predict_checks
— Methodrun_posterior_predict_checks(samples, paramnames, t, y, yerr, model, with_log_transform; S_low = 20, S_high = 20, is_integrated_power = true, plots = "all", n_samples = 100, path = "", basis_function = "SHO", n_frequencies = 1000, plot_f_P = false, n_components = 20)
Run the posterior predictive checks for the model and the approximation of the PSD
Arguments
samples::Array{Float64, 2}
: The samples from the posterior distributionparamnames::Array{String, 1}
: The names of the parameterst::Array{Float64, 1}
: The time seriesy::Array{Float64, 1}
: The values of the time seriesyerr::Array{Float64, 1}
: The errors of the time seriesmodel::Function
: The model or a string representing the modelwith_log_transform::Bool
: If true, the flux is log-transformedS_low::Float64=20.0
: the scaling factor for the appoximation at low frequenciesS_high::Float64=20.0
: the scaling factor for the appoximation at high frequenciesplots::String or Array{String, 1}
: The type of plots to make. It can be "all", "psd", "lsp", "timeseries" or a combination of them in an arrayn_samples::Int=200
: The number of samples to draw from the posterior predictive distributionpath::String="all"
: The path to save the plotsbasis_function::String="SHO"
: The basis function for the approximation of the PSDis_integrated_power = true
: if the norm corresponds to integral of the PSD betweenf_min
andf_max
or if it is the integral from 0 to infinity.n_frequencies::Int=1000
: The number of frequencies to use for the approximation of the PSDplot_f_P::Bool=false
: If true, the plots are made in terms of f * PSDn_components::Int=20
: The number of components to use for the approximation of the PSD
Pioran.sample_approx_model
— Methodsample_approx_model(samples, norm_samples, f0, fM, model; n_frequencies=1_000, basis_function="SHO", n_components=20)
Check the approximation of the PSD by computing the residuals and the ratios of the PSD and the approximated PSD
Arguments
samples::Array{Float64, n}
: The model samplesnorm_samples::Array{Float64, 1}
: The normalisation samplesf0::Float64
: The minimum frequency for the approximation of the PSDfM::Float64
: The maximum frequency for the approximation of the PSDmodel::Function
: The modeln_frequencies::Int=1_000
: The number of frequencies to use for the approximation of the PSDbasis_function::String="SHO"
: The basis function for the approximation of the PSDn_components::Int=20
: The number of components to use for the approximation of the PSDis_integrated_power=true
: Does the normalisation correspond to the integral of the PSD between two frequencies? If false, it corresponds to the true variance of the process.
Return
psd::Array{Float64, 2}
: The PSDpsd_approx::Array{Float64, 2}
: The approximated PSDresiduals::Array{Float64, 2}
: The residuals (psd-approx_psd)ratios::Array{Float64, 2}
: The ratios (approx_psd/psd)
Individual plotting functions
Pioran.get_ppc_timeseries
— Methodget_ppc_timeseries(samples, t, y, yerr, GP_model, with_log_transform; t_pred = nothing, n_samples = 1000)
Get the random posterior predictive time series from the model and samples
Pioran.plot_boxplot_psd_approx
— Methodplot_boxplot_psd_approx(residuals, ratios; path="")
Plot the boxplot of the residuals and ratios for the PSD approximation
Pioran.plot_lsp_ppc
— Methodplot_lsp_ppc(samples, t, y, yerr, GP_model; plot_f_P=false, n_frequencies=1000, n_samples=1000, is_integrated_power = true, n_components=20, bin_fact=10, path="", basis_function="SHO",with_log_transform = true)
Plot the posterior predictive Lomb-Scargle periodogram of random time series from the model to compare with the one of the data.
Arguments
samples::Array{Float64, n}
: Posterior samplest::Array{Float64, 1}
: Time seriesy::Array{Float64, 1}
: The values of the time seriesyerr::Array{Float64, 1}
: The errors of the time seriesGP_model::Function
: The GP model is a function that returns a ScalableGP conditionned on t and yerr, its arguments are t,y,yerr,paramsS_low = 20.0
: scaling factor for the lowest frequency in the approximation.S_high = 20.0
: scaling factor for the highest frequency in the approximation.plot_f_P::Bool=false
: If true, the plot is made in terms of f * PSDn_frequencies::Int=1000
: The number of frequencies to use for the approximation of the PSDn_samples::Int=1000
: The number of samples to draw from the posterior predictive distributionis_integrated_power::Bool = true
: if the norm corresponds to integral of the PSD betweenf_min
andf_max
or if it is the integral from 0 to infinity.n_components::Int=20
: The number of components to use for the approximation of the PSDbin_fact::Int=10
: The binning factor for the LSPpath::String=""
: The path to save the plotsbasis_function::String="SHO"
: The basis function for the approximation of the PSDwith_log_transform = true
: Apply a log transform to the data.
Pioran.plot_mean_approx
— Methodplot_mean_approx(f, residuals, ratios; path="")
Plot the frequency-averaged residuals and ratios
Pioran.plot_ppc_timeseries
— Methodplot_ppc_timeseries(samples, t, y, yerr, GP_model, with_log_transform; t_pred = nothing, n_samples = 100, path = "")
Plot the posterior predictive time series and the residuals
Arguments
samples::Matrix{Float64}
: The samples of the model parameterst::Vector{Float64}
: The time seriesy::Vector{Float64}
: The values of the time seriesyerr::Vector{Float64}
: The errors of the time seriesGP_model::PowerSpectralDensity
: The modelwith_log_transform::Bool
: If true, the flux was log-transformed for the inferencet_pred::Vector{Float64}=nothing
: The prediction timesn_samples::Int=100
: The number of samples to draw from the posterior predictive distribution
Pioran.plot_psd_ppc
— Methodplot_psd_ppc(samples_𝓟, samples_norm, samples_ν, t, y, yerr, model; S_low = 20.0, S_high = 20.0, plot_f_P = false, n_frequencies = 1000, path = "", n_components = 20, basis_function = "SHO", is_integrated_power = true, with_log_transform = false, save_samples = false)
Plot the posterior predictive power spectral density and the noise levels
Arguments
samples_𝓟::Array{Float64, n}
: The samples of the model parameterssamples_norm::Array{Float64, 1}
: The normalisation samples, either variance or integrated power betweenf_min
orf_max
.samples_ν::Array{Float64, 1}
: The ν samplest::Array{Float64, 1}
: The time seriesy::Array{Float64, 1}
: The values of the time seriesyerr::Array{Float64, 1}
: The errors of the time seriesmodel::Function
: The modelS_low = 20.0
: scaling factor for the lowest frequency in the approximation.S_high = 20.0
: scaling factor for the highest frequency in the approximation.plot_f_P::Bool = false
: If true, the plot is made in terms of f * PSDn_frequencies::Int = 1000
: The number of frequencies to use for the approximation of the PSDpath::String = ""
: The path to save the plotsn_components::Int = 20
: The number of components to use for the approximation of the PSDbasis_function::String = "SHO"
: The basis function for the approximation of the PSDis_integrated_power::Bool = true
: if the norm corresponds to integral of the PSD betweenf_min
andf_max
or if it is the integral from 0 to infinity.with_log_transform::Bool = false
: If true, the flux is log-transformed.save_samples::Bool = false
: Save samples of the psd and approx
Pioran.plot_quantiles_approx
— Methodplot_quantiles_approx(f, f_min, f_max, residuals, ratios; path="")
Plot the quantiles of the residuals and ratios (with respect to the approximated PSD) of the PSD
Pioran.plot_residuals_diagnostics
— Methodplot_residuals_diagnostics(t, mean_res, res_quantiles; confidence_intervals=[95, 99], path="")
Plot the residuals and the autocorrelation function of the residuals
Arguments
t::Vector{Float64}
: The time seriesmean_res::Vector{Float64}
: The mean of the residualsres_quantiles::Array{Float64, 2}
: The quantiles of the residualsconfidence_intervals::Array{Int, 1}
: The confidence intervalspath::String=""
: The path to save the plot
Pioran.plot_simu_ppc_timeseries
— Methodplot_simu_ppc_timeseries(t_pred, ts_quantiles, t, y, yerr; path="")
Plot the posterior predictive simulated time series
Arguments
t_pred::Vector{Float64}
: The prediction timests_quantiles::Array{Float64, 2}
: The quantiles of the posterior predictive time seriest::Vector{Float64}
: The time seriesy::Vector{Float64}
: The values of the time seriesyerr::Vector{Float64}
: The errors of the time seriespath::String=""
: The path to save the plot
Utilities
Pioran.check_conjugate_pair
— Methodcheck_conjugate_pair(r::Vector{Complex})
Check if the roots are complex conjugate pairs and negative real parts Returns true if the roots are complex conjugate pairs and false otherwise
Pioran.check_order_imag_roots
— Methodcheck_order_imag_roots(r::Vector{Complex})
Check if the imaginary parts of the roots are in ascending order
Pioran.check_roots_bounds
— Methodcheck_roots_bounds(r::Vector{Complex},f_min::Float64,f_max::Float64)
Check if the roots are within the bounds of the frequency range
Pioran.extract_subset
— Methodextract_subset(rng, prefix, t, y, yerr; n_perc=0.03, take_log=true, suffix="")
extract_subset(seed, prefix, t, y, yerr; n_perc=0.03, take_log=true)
Extract a subset of the data for the analysis and return initial guesses for the mean and variance. Either a random number generator or a seed can be provided.
Arguments
seed::Int64
: Seed for the random number generator.rng::AbstractRNG
: Random number generator.prefix::String
: Prefix for the output files.t::Array{Float64,1}
: Time array.y::Array{Float64,1}
: Time series array.yerr::Array{Float64,1}
: Time series error array.n_perc::Float64
: Percentage of the time series to extract.take_log::Bool
: If true, log transform the time series for the estimation of the mean and variance.suffix::String
: Suffix for the output files.
Return
t_subset::Array{Float64,1}
: Time array of the subset.y_subset::Array{Float64,1}
: Time series array of the subset.yerr_subset::Array{Float64,1}
: Time series error array of the subset.x̄::Float64
: Mean of the normal distribution for μ.va::Float64
: Variance of the normal distribution for μ.
Pioran.makelist_namessplit
— Methodsimply put entries of paramnames_split into lists even when it's a string
Pioran.separate_samples
— Methodseparate_samples(samples, paramnames, with_log_transform::Bool)
Separate the samples into the parameters of the model and the parameters of the power spectral density.
Prior distributions
Pioran.ThreeUniformDependent
— TypeThreeUniformDependent(a, b, c, ϵ)
ThreeUniformDependent(a, b, c) (constructor with default ϵ = 1e-10)
Multivariate distribution to model three random variables where the first one x1 is given by U[a,b] and the second one x2 is given by U[x1,c] and the third one x3 is given by U[x2,c]. where a<b<c.
a
: lower bound of the first distributionb
: upper bound of the first distributionc
: upper bound of the second and third distributionϵ
: small value to make sure that the lower and upper bounds of each distribution are different
This means that the lower bound of the second distribution is dependent on the value of the first distribution and so on... This is implemented to overcome the limitations of the current Turing's implementation for dependent priors with dynamic support. See the following issues for more details: [1],[2],[3].
Pioran.TwoLogUniformDependent
— TypeTwoLogUniformDependent(a, b, ϵ)
TwoLogUniformDependent(a, b) (constructor with default ϵ = 1e-10
Multivariate distribution to model three random variables where the first one x1 is given by log-U[a,b] and the second one x2 is given by log-U[x1,b].
a
: lower bound of the first distributionb
: upper bound of the first distributionϵ
: small value to make sure that the lower and upper bounds of each distribution are different
This means that the lower bound of the second distribution is dependent on the value of the first distribution. This is implemented to overcome the limitations of the current Turing's implementation for dependent priors with dynamic support. See the following issues for more details: [1],[2],[3].
Pioran.TwoUniformDependent
— TypeTwoUniformDependent(a, b, c, ϵ)
TwoUniformDependent(a, b, c) (constructor with default ϵ = 1e-10)
Multivariate distribution to model two random variables where the first one is given by U[a,b] and the second one is given by U[x,c], where x is a random variable sampled from the first distribution.
a
: lower bound of the first distributionb
: upper bound of the first distributionc
: upper bound of the second distributionϵ
: small value to make sure that the lower and upper bounds of each distribution are different
This means that the lower bound of the second distribution is dependent on the value of the first distribution.This is implemented to overcome the limitations of the current Turing's implementation for dependent priors with dynamic support. See the following issues for more details: [1],[2],[3].
Example
#```jldoctest julia> using Pioran, Distributions julia> d = TwoUniformDependent(0, 1, 2) TwoUniformDependent(0.0, 1.0, 2.0)
julia> rand(d) 2-element Array{Float64,1}: 0.123 1.234 #```
Bijectors.bijector
— MethodBijectors.bijector(d::TwoUniformDependent)
Create a bijector for the TwoUniformDependent distribution. This is used to sample from the distribution using the Bijectors package. Adapted from the following issues [1],[2],[3].
CARMA
Pioran.CARMA
— TypeCARMA(p, q, rα, β, σ²)
Continuous-time AutoRegressive Moving Average (CARMA) model for the power spectral density
p
: the order of the autoregressive polynomialq
: the order of the moving average polynomialrα
: roots of the autoregressive polynomial length p+1β
: the moving average coefficients length q+1σ²
: the variance of the process
The power spectral density of the CARMA model is given by:
\[\mathcal{P}(f) = \sigma^2 \left|\dfrac{\sum\limits_{k=0}^q \beta_k \left(2\pi\mathrm{i}f\right)^k }{\sum\limits_{l=0}^p \alpha_l \left(2\pi\mathrm{i}f\right)^l}\right|^2\]
Pioran.evaluate
— Methodevaluate(model::CARMA, f)
evaluate the power spectral density of the CARMA model at frequency f.
Pioran.quad2roots
— Methodquad2roots(quad)
Convert the coefficients of a quadratic polynomial to its roots.
Arguments
quad::Vector{Real}
: Coefficients of the quadratic polynomial.
Returns
r::Vector{Complex}
: Roots of the polynomial.
Pioran.roots2coeffs
— Methodroots2coeffs(r)
Convert the roots of a polynomial to its coefficients.
Arguments
r::Vector{Complex}
: Roots of the polynomial.
Returns
c::Vector{Complex}
: Coefficients of the polynomial.