The peakbagging module#
Fit functions#
- apollinaire.peakbagging.stellar_framework(freq, psd, r=1, m=1, teff=5770, nwalkers=64, dnu=None, back=None, wdw=None, numax=None, Hmax=None, Wenv=None, epsilon=None, width=None, alpha=None, d02=None, d01=None, d13=None, n_harvey=2, low_cut=10.0, guess_back=None, low_bounds_back=None, up_bounds_back=None, frozen_param_back=None, power_law=False, high_cut_plaw=20.0, spectro=False, filename_back='background.png', back_mcmc=True, back_mle=False, filemcmc_back='mcmc_sampler_background.h5', nsteps_mcmc_back=1000, guess_pattern=None, low_bounds_pattern=None, up_bounds_pattern=None, pattern_mcmc=True, pattern_mle=False, n_order=3, n_order_peakbagging=3, h0_screening=False, format_cornerplot='pdf', filename_pattern='pattern.png', filemcmc_pattern='mcmc_sampler_pattern.h5', nsteps_mcmc_pattern=1000, parallelise=False, fit_l1=True, fit_l3=False, use_numax_back=False, progress=False, a2z_file='modes_param.a2z', a2z_in=None, nsteps_mcmc_peakbagging=1000, mcmcDir='.', instr='geometric', amp_l=[1.0, 1.5, 0.7, 0.2], filename_peakbagging='summary_peakbagging.png', nopeakbagging=False, discard_back=200, discard_pattern=200, discard_pkb=200, thin_back=1, thin_pattern=1, thin_pkb=1, quickfit=False, num=500, num_numax=500, reboxing_behaviour='advanced_reboxing', reboxing_strategy='median', save_only_after_sampling=False, apodisation=False, fit_angle=False, guess_angle=90, fit_splittings=False, fit_splittings_per_order_in_global=False, guess_split=0, extended_splittings=False, fit_amplitude_ratio=False, frozen_harvey_exponent=False, fit_amp=False, extended=False, author=None, pkb_fmt=None, projected_splittings=False, bins=100, existing_chains_back='read', existing_chains_pattern='read', existing_chains_pkb='read', strategy='order', common_width=False, save_background_vector=True, dpi=300, show_corner=True, estimate_autocorrelation=False, plot_datapoints=True, show_gaussian=True)#
Framework for stellar peakbagging considering only a few input parameters. Background, global pattern and individual mode parameters are successively fitted.
- Parameters:
freq (ndarray) – frequency vector
psd (ndarray) – real power or power spectral density vector
r (float) – stellar radius. Given in solar radius. Optional, default 1.
m (float) – stellar mass. Given in solar masses. Optional, default 1.
teff (float) – stellar effective temperature. Given in Kelvins. Optional, default 5770.
nwalkers (int) – number of walkers for the MCMC process. Optional, default 64.
dnu (float) – large frequency separation. Must be given in µHz. If
None, will be computed through scaling laws. Optional, defaultNone.back (ndarray) – array of fitted background, of same length as
freqandpsd. If given, no fit will be performed on background and the pattern step will directly take place.Hmax,Wenvandnumaxmust in this case also be given. Optional, defaultNone.wdw (ndarray) – observational window of the light curve to take into account the side lobes pattern in the modes. Must be 1 when an actuak observation is acquired and zero otherwise. Optional, default
None.numax (float) – maximum of power in p-mode envelope. Must be given in µHz. Optional, default
None.Hmax (float) – amplitude of p-mode envelope. Optional, default
None.Wenv (float) – width of p-mode envelope. Optional, default
None.epsilon (float) – epsilon guess value of mode patterns. If given, will override the value computed by the automatic guess. Optional, default
None.width (float) – width guess value of mode patterns. If given, will override the value computed by the automatic guess. Optional, default
None.alpha (float) – alpha value for global pattern guess. If specified, this value will override guess value computed by the automatic guess. Optional, default
None.d02 (float) – d02 value for global pattern guess. If specified, this value will override guess value computed by the automatic guess. Optional, default
None.d01 (float) – d01 value for global pattern guess. If specified, this value will override guess value computed by the automatic guess. Optional, default
None.d13 (float) – d13 value for global pattern guess. If specified, this value will override guess value computed by the automatic guess. Optional, default
None.n_harvey (int) – number of Harvey laws to use to build the background model. Optional, default 2. With more than two Harvey laws, it is strongly recommended to manually feed the ‘guess’ parameter.
guess_back (array-like.) – first guess directly passed by the user. If guess is
None, the function will automatically infer a first guess. Optional, defaultNone. Backgrounds parameters given in the following order: param_harvey (3n_harvey), param_plaw (2), param_gaussian (3), white noise (1). You can usecreate_background_guess_arraysto create guess and bounds arrays with proper structure.low_bounds_back (ndarray) – lower bounds to consider in the background parameters space exploration. Must have the same dimension than
guess_back.up_bounds_back (ndarray) – upper bounds to consider in the background parameters space exploration. Must have the same dimension than
guess_back.frozen_param_back (boolean array) – boolean array of the same size as
guess_back. Components set toTruewill not be fitted.power_law (bool) – set to
Trueto fit a power law on the background. Optional, defaultFalse.high_cut_plaw (float) – if
power_law=True, the function will consider the PSD segment betweenlow_cutandhigh_cut_plawto compute the power law parameters initial guess.spectro (bool) – if set to
True, make the plot with unit consistent with radial velocity, else with photometry. Optional, defaultFalse.filename_back (str) – name of the of the dat file of the fitted background, the background parameters pdf file where the plot will be stored. Optional, default
background.back_mcmc (bool) – If set to
False, no MCMC sampling will be performed for the background. Optional, defaultTrue.back_mle (bool) – If set to
False, no MLE optimisation will be performed for the background. Optional, defaultTrue.filemcmc_back (str) – Name of the hdf5 file where the MCMC background chain will be stored. Optional, default
mcmc_sampler_background.h5.nsteps_mcmc_back (int) – number of MCMC steps for the background parameters exploration. Optional, default 1000.
guess_pattern (array-like.) – first guess directly passed by the user. If guess is
None, the function will automatically infer a first guess. Optional, defaultNone. Pattern parameters have to be given in the following order:eps,alpha,Dnu,numax,Hmax,Wenv,w,d02,b02,d01,b01,d13,b03. A first version of the guess and bounds array can be created through thecreate_pattern_guess_arrays.low_bounds_pattern (ndarray) – lower bounds to consider in the background parameters space exploration. Must have the same dimension than
guess_pattern.up_bounds_pattern (ndarray) – upper bounds to consider in the background parameters space exploration. Must have the same dimension than
guess_pattern.pattern_mle (bool) – If set to
False, no MLE optimisation will be performed for the pattern. Optional, defaultTrue.pattern_mcmc (bool) – If set to
False, no MCMC sampling will be performed for the pattern. Optional, defaultTrue.n_order (int) – number of orders to consider for the global pattern fit. Optional, default 3.
n_order_peakbagging (int) – number of orders to fit at the individual mode parameters step. Optional, default 3.
h0_screening (bool) – if set to
True, a fast frequency screening of the p-mode region will be performed in order to determine the mode that it possible to fit.n_order_peakbaggingvalue will be ignored in this case. Optional, defaultFalse.format_cornerplot (str) – set cornerplot file format. Optional, default
pdf.filename_pattern (str) – name of the pdf file where the fitted global pattern parameters and the plot of the fitted global pattern will be stored. Optional, default
pattern.filemcmc_pattern (str) – Name of the hdf5 file where the MCMC pattern chain will be stored. Optional, default
mcmc_sampler_pattern.h5.nsteps_mcmc_pattern (int) – number of MCMC steps for the pattern parameters exploration. Optional, default 1000.
parallelise (bool) – if set to
True, multiprocessing will be used to sample the MCMC. Optional, defaultFalse.fit_l1 (bool) – set to
Trueto fit the d01 and b03 param and create guess for l=1 modes. Optional, defaultTrue.fit_l3 (bool) – set to
Trueto fit the d03 and b03 param and create guess for l=3 modes. Optional, defaultFalse.fit_l1must be set toTrue.use_numax_back (bool) – if set to
True, thenumaxvalues computed with the background fit will be used as input guess for the pattern fit. Otherwise, the initialnumaxfor the pattern fit will be taken from thenumaxargument or computed with the scaling laws.progress (bool) – show progress bar in terminal. Optional, default
False.a2z_file (str) – name of the a2z file where the individual parameters will be stored. This will also be used to name the output pkb file. Optional, default
modes_param.a2za2z_in (str) – a2z file with guess input for individual parameters. If provided, the global pattern step will be ignored and the function will directly sample individual modes parameters (after sampling the background if it has not been provided too). Optional, default
None.nsteps_mcmc_peakbagging (int) – number of MCMC steps for the peakbagging parameters exploration. Optional, default 1000.
mcmcDir (str) – Name of the directory where the MCMC of the individual parameters should be stored. The current directory will be used if no path is provided. Optional, default
None.instr (str) – name of the instrument for which m-amplitude ratios will be computed.
geometric,golfandvirgoare available.amp_l (array-like) – amplitude ratio between harmonic degrees (relative to l=0 mode). Optional, default [1., 1.5, 0.7, 0.2].
filename_peakbagging (str) – name of the file where the plot summary of the individual mode parameters peakbagging will be stored. Optional, default
summary_peakbagging.pdf.nopeakbagging (bool) – if set to
True, individual modes parameters will not be fitted. Optional, defaultFalse.discard_back (int) – number of step to discard for the background sampling. Optional, default 200.
discard_pattern (int) – number of step to discard for the pattern sampling. Optional, default 200.
discard_pkb (int) – number of step to discard for the peakbagging sampling. Optional, default 200.
thin_back (int) – take only every
thinsteps from the chain sampled for backgrounds parameters. Optional, default 10.thin_pattern (int) – take only every
thinsteps from the chain sampled for pattern parameters. Optional, default 10.quickfit (bool) – if set to
True, the fit will be performed over a smoothed and logarithmically resampled background. Optional, defaultFalse.num (int) – Only considered if
quickfit=True. Number of point in the resampled PSD, or target number of point for region away fromnumaxifadvanced_reboxingis selected. Optional, default5000.num_numax (int) – Only considered if
quickfit=Trueandreboxing_behaviour=advanced_reboxing. Number of resampling points in thenumaxregion. Optional, default1000.reboxing_behaviour (str) – behaviour for
quickfitnew power vector computation. It can besmoothing,reboxingoradvanced_reboxing. Optional, defaultadvanced_reboxing.reboxing_strategy (str) – reboxing strategy to take box values when using
quickfitandreboxing_behaviour=reboxing. Can bemedianormean. Optional, defaultmedian.save_only_after_sampling – if set to True, hdf5 file with chains information will only be saved at the end of the sampling process. If set to False, the file will be saved step by step (see
emceedocumentation).apodisation (bool) – if set to
True, distort the model to take the apodisation into account. Optional, defaultFalse.fit_angle (bool) – if set to
True, the peakbagging process will include inclination angle of the star in the parameters to fit. Optional, defaultFalse.guess_angle (float) – initial guess value for angle. Will be taken into account only if
fit_angleisTrue. Optional, default 90.fit_splittings (bool) – if set to
True, the peakbagging process will include mode splittings in the parameters to fit. Optional, defaultFalse.fit_splittings_per_order_in_global (bool) – if set to
Trueandfit_splittings=Truewithstrategy='global', the peakbagging process will include one mode splittings per order in the parameters to fit. Optional, defaultFalse.guess_split (float) – initial guess value for splittings. Will be taken into account only if
fit_splittingsisTrue. Optional, default 0.extended_splittings (bool) – if set to
True, the splittings will be fitted according to the formula provided by Eq.(S9) and (S15) of Benomar et al. (2018).fit_amplitude_ratio (bool) – if the selected strategy is
global, the globalamp_lamplitude ratio between degrees parameters will be fitted only if this parameter is set toTrue. Optional, defaultFalse.frozen_harvey_exponent (bool) – set to True to have the Harvey models exponent fixed to 4. Optional, default
False.fit_amp (bool) – if set to
True, amplitude of the modes will be fitted instead of heights (input heights will be transformed inside the code and retransformed as outputs).extended (bool) – if set to
True, the returned pkb array will be extended to contain asymmetric uncertainties over each parameter.author (str) – identifier of the peakbagger, if given it will be specified in the header of the output pkb file. Optional, default
None.pkb_fmt – float formatting of the pkb file that will be created by the function. Optional, default
None.projected_splittings (bool) – if set to
True, the function will consider that thesplitparameters of the input are projected splittings and will build the model consequently. Optional, defaultFalse.bins (int) – number of bins for each cornerplot panel. Optional, default 100.
existing_chains_back – controls the behaviour of the function concerning a background parameters existing hdf5 file. If,
read, existing files will be read without sampling and function output will be updated consequently, ifreset, the backend will be cleared and the chain will be sampled from scratch, ifsamplethe function will sample the chain from where the previous exploration was stopped, ifignore, the chain is not read. Optional, defaultread. The optionignorewill apply only to individual mode chains, it will be changed toreadfor background and pattern steps.existing_chains_pattern (str) – Same as
existing_chains_backbut for a global pattern parameters existing hdf5 file.existing_chains_pkb (str) – Same as
existing_chains_backbut for individual mode parameters existing hdf5 files.strategy (str) – fitting strategy for individual modes. Can be
orderorglobal. Optional, defaultorder. Note that theglobalbehaviour is experimental.common_width (bool) – if set to
True, only a global width parameter will be set in the a2z DataFrame. Optional, defaultFalse.save_background_vector (bool) – the full background vector will be saved only if this option is set to
True. Optional, defaultTrue.dpi (int) – dot per inch value for figures.
show_corner (bool) – if set to
True, corner plots will be shown and saved. Optional, defaultTrue.estimate_autocorrelation (bool) – if set to
True, the autocorrelation time of the sampled chains will be estimated byemcee. Optional, defaultFalse.plot_datapoints (bool) – data points outside contours will be drawn if set to
True. Optional, defaultTrue.show_gaussian (bool) – If set to
True, show p-mode acoustic hump modelled by a Gaussian function in the background summary plot. Optional, defaultTrue.
- Returns:
None
- apollinaire.peakbagging.peakbagging(a2z_file, freq, psd, back=None, wdw=None, dnu=None, spectro=True, nsteps_mcmc=1000, discard=200, show_corner=False, store_chains=False, mcmcDir=None, order_to_fit=None, parallelise=False, progress=False, strategy='order', fit_02=True, fit_13=True, nwalkers=64, normalise=False, instr='geometric', show_summary=False, filename_summary=None, show_prior=False, use_sinc=False, asym_profile='nigam-kosovichev', thin=1, restrict=False, remove_asymmetry_outside_window=True, do_not_use_dnu=False, save_only_after_sampling=False, size_window=None, fit_amp=False, extended=False, projected_splittings=False, bins=100, existing_chains='read', fit_splittings=True, extended_splittings=False, fit_angle=True, fit_amplitude_ratio=False, dpi=300, format_cornerplot='pdf', estimate_autocorrelation=False, plot_datapoints=True)#
Perform the peakbagging process for a given set of modes.
- Parameters:
a2z_file (string) – name of the file to read the parameters.
freq (ndarray) – frequency vector
psd (ndarray) – real power vector
back (ndarray.) – activity background vector that will be used to complete the model to fit. Optional default
None. Must have the same length than freq and psd.wdw (ndarray.) – observation window (0 and 1 array of the same lenght as the original timeseries) to analyse in order to predict the sidelob pattern. Optional, default
None.dnu (float.) – large separation of the modes. In this function, it will only be used to adapt the fitting window for fitted modes above 1200 muHz. If not given and if at least two l=0 modes are specified in the a2z file, dnu will be automatically computed.
nsteps_mcmc (int) – number of steps to process into each MCMC exploration.
discard (int) – number of steps to discard into each MCMC exploration.
show_corner (bool) – if set to
True, show and save the corner plot for the posterior distribution sampling.store_chains (bool) – if set to
True, each MCMC sampler will be stored as an hdf5 files. Filename will be autogenerated with modes parameters. Optional, defaultFalsemcmcDir (string) – directory where to save the MCMC sampler. If not path is provided, the current directory will be used. Optional, default
None.order_to_fit (array-like) – list of order to fit if the input a2z file contains order that are supposed not to be fitted. Optional, default
None.parallelise (bool) – If set to
True, use Python multiprocessing tool to parallelise process. Optional, defaultFalse.strategy (str) – strategy to use for the fit,
pair,order,global. Optional, defaultorder. If ‘pair’ is used, a2z input must contain individual heights, widths and splittings for each degree. Theglobalbehaviour is experimental.fit_02 (bool) – if strategy is set to
pair, l=0 and 2 modes will only be fitted if this parameter is set toTrue. Optional, defaultTrue.fit_13 (bool) – if strategy is set to
pair, l=1 and 3 modes will only be fitted if this parameter is set toTrue. Optional, defaultTrue.nwalkers (int) – number of walkers to use for the sampling.
instr (str) – instrument to consider (amplitude ratio inside degrees depend on geometry AND instrument and should be adaptated). Possible argument :
geometric,kepler,golf,virgo. Optional, defaultgeometric.show_summary (bool) – show summary plot at the end of the fit. Optional, default
False.use_sinc (bool) – if set to
True, mode profiles will be computed using cardinal sinus and not Lorentzians. No asymmetry term will be used if it is the case. Optional, defaultFalse.asym_profile (str) – depending on the chosen argument, asymmetric profiles will be computed following Korzennik 2005 (
korzennik) or Nigam & Kosovichev 1998 (nigam-kosovichev). Defaultnigam-kosovichev.restrict (bool) – if set to True, will only use the modes with unfrozen parameters in the model.
remove_asymmetry_outside_window (bool) – if set to True, asymmetry of modes outside of the fitting window will be set to 0. Optional, default
True.do_not_use_dnu (bool) – if set to
True, fitting window will be computed without using dnu value.save_only_after_sampling – if set to True, hdf5 file with chains information will only be saved at the end of the sampling process. If set to False, the file will be saved step by step (see
emceedocumentation).size_window (float) – size of the fitting window, in muHz. If not given, the size of the fitting window will be automatically set, using dnu or the input frequencies of the parameter to fit. Optional, default
None.fit_amp (bool) – if set to
True, amplitude of the modes will be fitted instead of heights (input heights will be transformed inside the code and retransformed as outputs).extended (bool) – if set to
True, the returned pkb array will be extended to contain asymmetric uncertainties over each parameter.projected_splittings (bool) – if set to
True, the function will consider that thesplitparameters of the input are projected splittings and will build the model consequently. Optional, defaultFalse.bins (int) – number of bins for each cornerplot panel. Optional, default 100.
existing_chains (str) – controls the behaviour of the function concerning existing hdf5 files. If,
read, existing files will be read without sampling and function output will be updated consequently, ifreset, the backend will be cleared and the chain will be sampled from scratch, ifsamplethe function will sample the chain from where the previous exploration was stopped, ifignore, the chain is not read. Optional, defaultread.fit_splittings (bool) – if the selected strategy is
global, the globalsplitparameter will be fitted only if this parameter is set toTrue. Optional, defaultTrue.extended_splittings (bool) – if set to
True, the splittings will be fitted according to the formula provided by Eq.(S9) and (S15) of Benomar et al. (2018).fit_angle (bool) – if the selected strategy is
global, the globalangleparameter will be fitted only if this parameter is set toTrue. Optional, defaultTrue.fit_amplitude_ratio (bool) – if the selected strategy is
global, the globalamp_lamplitude ratio between degrees parameters will be fitted only if this parameter is set toTrue. Optional, defaultTrue.dpi (int) – dot per inch value for figures.
estimate_autocorrelation (bool) – if set to
True, the autocorrelation time of the sampled chains will be estimated byemcee. Optional, defaultFalse.plot_datapoints (bool) – data points outside contours will be drawn if set to
True. Optional, defaultTrue.
- Returns:
a2z fitted modes parameters as a DataFrame, and corresponding pkb array (extended pkb array if extended is set to
True)- Return type:
tuple of objects
- apollinaire.peakbagging.perform_mle_background(freq, psd, n_harvey=2, r=1.0, m=1.0, teff=5770.0, dnu=None, numax=None, guess=None, frozen_param=None, power_law=False, gaussian=True, frozen_harvey_exponent=False, low_cut=100.0, method='Powell', fit_log=False, low_bounds=None, up_bounds=None, no_bounds=False, show=True, show_guess=False, filename=None, spectro=True, quickfit=False, num=5000, num_numax=500, two_gaussian=False, asy_gaussian=False, reboxing_behaviour='advanced_reboxing', reboxing_strategy='median', apodisation=False, high_cut_plaw=20.0, show_gaussian=True, **kwargs)#
Perform MLE over background model.
- Parameters:
freq (ndarray) – frequency vector in µHz.
psd (ndarray) – power density vector in ppm^2/µHz or (m/s)^2/muHz.
n_harvey (int) – number of Harvey laws to use to build the background model. Optional, default 2. With more than two Harvey laws, it is strongly recommended to manually feed the ‘guest’ parameter.
r (float) – stellar radius. Given in solar radius. Optional, default 1.
m (float) – stellar mass. Given in solar mass. Optional, default 1.
teff (float) – stellar effective temperature. Given in K. Optional, default 5770.
dnu (float) – large separation of the p-mode. If given, it will superseed the scaling law guess using
r,mandteff.numax (float) – maximum p-mode bump power frequency. If given, it will superseed the scaling law guess using
r,mandteff.guess (array-like.) – first guess directly passed by the users. If guess is
None, the function will automatically infer a first guess. Optional, defaultNone. Backgrounds parameters given in the following order: param_harvey (3n_harvey), param_plaw (2), param_gaussian (3), white noise (1)frozen_param (boolean array) – boolean array of the same size as guess. Components set to
Truewill not be fitted.power_law (bool) – set to
Trueto fit a power law on the background. Optional, defaultFalse.gaussian (bool) – set to
Trueto fit the p-mode gaussian on the background. Optional, defaultTrue.frozen_harvey_exponent (bool) – set to
Trueto freeze and not fit Harvey laws exponent to 4. Optional, defaultFalse.low_cut (float) – Spectrum below this frequency will be ignored for the fit. The frequency value must be given in µHz. Optional, default 100.
method – minimization method used by the scipy minimize function. Optional, default
"Powell"fit_log (bool) – if set to
True, fit natural logarithm of the parameters. Optional, defaultFalse.low_bounds (ndarray) – lower bounds to consider in the parameter space exploration. Must have the same structure than guess.
up_bounds (ndarray) – upper bounds to consider in the parameter space exploration. Must have the same structure than guess.
no_bounds (bool) – If set to
True, no bounds will be considered for the parameter space exploration. Optional, defaultFalse.show (bool) – if set to
True, will show at the end a plot summarising the fit. Optional, defaultTrue.show_guess (bool) – if set to
True, will show at the beginning a plot summarising the guess. Optional, defaultFalse.filename (str) – if given, the summary plot will be saved under this filename. Optional, default
None.showargument must have been set toTrue.spectro (bool) – if set to
True, make the plot with unit consistent with radial velocity, else with photometry. Automated guess will also be computed consistently with spectroscopic measurements. Optional, defaultTrue.quickfit (bool) – if set to
True, the fit will be performed over a smoothed and logarithmically resampled background. Optional, defaultFalse.num (int) – Only considered if
quickfit=True. Number of point in the resampled PSD, or target number of point for region away fromnumaxifadvanced_reboxingis selected. Optional, default5000.num_numax (int) – Only considered if
quickfit=Trueandreboxing_behaviour=advanced_reboxing. Number of resampling points in thenumaxregion. Optional, default500.reboxing_behaviour (str) – behaviour for
quickfitnew power vector computation. It can besmoothing,reboxingoradvanced_reboxing. Optional, defaultadvanced_reboxing.reboxing_strategy (str) – reboxing strategy to take box values when using
quickfitandreboxing_behaviour=reboxing. Can bemedianormean. Optional, defaultmedian.apodisation (bool) – if set to
True, distort the model to take the apodisation into account. Optional, defaultFalse.high_cut_plaw (float) – high cut in frequency to compute the power law guess. Optional, default 20.
show_gaussian (bool) – If set to
True, show p-mode acoustic hump modelled by a Gaussian function in the summary plot. Optional, defaultTrue.
- Returns:
fitted background and fitted background parameters.
- Return type:
tuple of array
- apollinaire.peakbagging.explore_distribution_background(freq, psd, n_harvey=2, r=1.0, m=1.0, teff=5770.0, dnu=None, numax=None, guess=None, frozen_param=None, power_law=False, gaussian=True, frozen_harvey_exponent=False, low_cut=100.0, fit_log=False, low_bounds=None, up_bounds=None, spectro=True, show=True, show_guess=False, show_corner=True, nsteps=1000, discard=200, filename=None, parallelise=False, progress=False, nwalkers=64, filemcmc=None, thin=1, quickfit=False, num=5000, num_numax=500, estimate_autocorrelation=False, reboxing_behaviour='advanced_reboxing', reboxing_strategy='median', format_cornerplot='pdf', save_only_after_sampling=False, apodisation=False, high_cut_plaw=20.0, bins=100, two_gaussian=False, asy_gaussian=False, existing_chains='read', norm=None, plot_datapoints=True, xlim=None, ylim=None, show_gaussian=True, **kwargs)#
Use a MCMC framework to fit the background.
- Parameters:
freq (ndarray) – frequency vector
psd – real power vector
n_harvey (int) – number of Harvey laws to use to build the background model. Optional, default 2. With more than two Harvey laws, it is strongly recommended to manually feed the
guessparameter.r (float) – stellar radius. Given in solar radius. Optional, default 1.
m (float) – stellar mass. Given in solar mass. Optional, default 1.
teff (float) – stellar effective temperature. Given in K. Optional, default 5770.
dnu (float) – large separation of the p-mode. If given, it will superseed the scaling law guess using
r,mandteff.numax (float) – maximum p-mode bump power frequency. If given, it will superseed the scaling law guess using
r,mandteff.guess (array-like.) – first guess directly passed by the user. If guess is
None, the function will automatically infer a first guess. Optional, defaultNone. Backgrounds parameters given in the following order: param_harvey (3n_harvey), param_plaw (2), param_gaussian (3), white noise (1)frozen_param (boolean array) – boolean array of the same size as guess. Components set to
Truewill not be fitted.power_law (bool) – set to
Trueto fit a power law on the background. Optional, defaultFalse.gaussian (bool) – set to
Trueto fit the p-mode gaussian on the background. Optional, defaultTrue.frozen_harvey_exponent (bool) – set to
Trueto freeze and not fit Harvey laws exponent to 4. Optional, defaultFalse.low_cut (float) – Spectrum below this frequency will be ignored for the fit. The frequency value must be given in µHz. Optional, default 100.
fit_log (bool) – if set to
True, fit natural logarithm of the parameters. Optional, defaultFalse.low_bounds (ndarray) – lower bounds to consider in the parameter space exploration. Must have the same structure than
guess.up_bounds (ndarray) – upper bounds to consider in the parameter space exploration. Must have the same structure than
guess.spectro (bool) – if set to
True, make the plot with unit consistent with radial velocity, else with photometry. Automated guess will also be computed consistently with spectroscopic measurements. Optional, defaultTrue.filename (string) – Name of the file to save the summary plot. If filename is
None, the name will not be stored. Optional, defaultNone.filemcmc – name of the hdf5 where to store the chain. If filename is
None, the name will not be stored. Optional, defaultNone.parallelise (bool) – If set to
True, use Python multiprocessing tool to parallelise process. Optional, defaultFalse.show (bool) – if set to
True, will show at the end a plot summarising the fit. Optional, defaultTrue.show_corner – if set to
True, will show the corner plot summarising the MCMC process. Plot will be saved as a pdf isfilemcmcis also specified. Optional, defaultTrue.nsteps (int) – number of MCMC steps to perform.
discard (int) – number of step to discard in the sampling. Optional, default 200.
show_guess (bool) – if set to
True, will show at the beginning a plot summarising the guess. Optional, defaultFalse.thin (int) – take only every
thinsteps from the chain. Optional, default 1.quickfit (bool) – if set to
True, the fit will be performed over a smoothed and logarithmically resampled background. Optional, defaultFalse.num (int) – Only considered if
quickfit=True. Number of point in the resampled array, or target number of point for region away fromnumaxifadvanced_reboxingis selected. Optional, default5000.num_numax (int) – Only considered if
quickfit=Trueandreboxing_behaviour=advanced_reboxing. Number of resampling points in thenumaxregion. Optional, default1000.reboxing_behaviour (str) – behaviour for
quickfitnew power vector computation. It can besmoothing,reboxingoradvanced_reboxing. Optional, defaultadvanced_reboxing.reboxing_strategy (str) – reboxing strategy to take box values when using
quickfitandreboxing_behaviour=reboxing. Can bemedianormean. Optional, defaultmedian.save_only_after_sampling – if set to True, hdf5 file with chains information will only be saved at the end of the sampling process. If set to False, the file will be saved step by step (see
emceedocumentation).apodisation (bool) – if set to
True, distort the model to take the apodisation into account. Optional, defaultFalse.high_cut_plaw (float) – high cut in frequency to compute the power law guess. Optional, default 20.
bins (int) – number of bins for each cornerplot panel. Optional, default 100.
existing_chains (str) – controls the behaviour of the function concerning existing hdf5 files. If,
read, existing files will be read without sampling and function output will be updated consequently, ifreset, the backend will be cleared and the chain will be sampled from scratch, ifsamplethe function will sample the chain from where the previous exploration was stopped. Optional, defaultread.estimate_autocorrelation (bool) – if set to
True, the autocorrelation time of the sampled chains will be estimated byemcee. Optional, defaultFalse.plot_datapoints (bool) – data points outside contours will be drawn if set to
True. Optional, defaultTrue.show_gaussian (bool) – If set to
True, show p-mode acoustic hump modelled by a Gaussian function in the summary plot. Optional, defaultTrue.
- Returns:
the fitted background, its param and their sigma values as determined by the MCMC walk.
- Return type:
tuple of array
Chain management functions#
- apollinaire.peakbagging.read_chain(filename, thin=1, discard=0, chain_type='peakbagging', fit_amp=False, projected_splittings=False)#
Read a chain created by
apollinaireand return flattened sampled chain with auxiliary informations.- Parameters:
filename (str) – name of the hdf5 file with the stored chain.
thin (int) – the returned chain will be thinned according to this factor. Optional, default 1.
discard (int) – number of iteration step to discard. Optional, default 0.
chain_type (str) – can be
peakbagging,patternorbackground. Optional, defaultpeakbagging.fit_amp (bool) – set to
Trueif amplitudes were fitted instead of heights.projected_splittings (bool) – set to
Trueif projected splittings were fitted instead of splittings.
- Returns:
tuple with flatchain, param labels, degrees of param and order of modes if
chain_type=peakbagging, flatchain and labels otherwise.- Return type:
tuple of ndarrays
- apollinaire.peakbagging.make_cornerplot(flatchain, labels=None, filename=None, bins=10, figsize=(24, 24), reverse=False, quantiles=None, title_fmt='.4f', **kwargs)#
Make a cornerplot from a flattened chain. This is simply a wrapper for
corner.corner.- Parameters:
flatchain (ndarray) – the flattened chain.
labels (array-like) – parameters labels
filename (str) – name of the file where the cornerplot will be saved.
bins (int) – number of bins in the histograms. Optional, default
10.quantiles (list or tuple of three elements.) – quantiles to show on the 1d distributions. Optional, default
None(in this case[0.16, 05. 0.84]will be used).
- apollinaire.peakbagging.hdf5_to_a2z(workDir='.', a2zname='summary_fit.a2z', discard=0, thin=10, instr='geometric', centile=False, verbose=False)#
Create a a2z file with the hdf5 files stored in a given folder.
- Parameters:
workDir (str) – the directory where the hdf5 file will be read.
a2zname – the name of the output file. Set to
Noneif you do not want to save any file.discard (int) – the number of elements to ignore at the beginning of the chain.
thin (int) – one element of the chain every
thinelements will be considered.centile (bool) – if set to True, change the structure of the output to write a centile frame. Optional, default
False.instr (str) – instrument to consider. Possible argument :
geometric,kepler,golf,virgo. Optional, defaultgeometric.verbose (bool) – if set to
True, print the name of the considered files. Optional, defaultFalse.
- Returns:
a2z DataFrame
- Return type:
pandas DataFrame
- apollinaire.peakbagging.hdf5_to_pkb(workDir='.', pkbname='summary_fit.pkb', discard=0, thin=10, instr='geometric', extended=False, verbose=False)#
Create a pkb file with the hdf5 files stored in a given folder.
- Parameters:
workDir (str) – the directory where the hdf5 file will be read.
pkbname (str) – the name of the output pkb. Set to
Noneif you do not want to save a file.discard (int) – the number of elements to ignore at the beginning of the chain.
thin (int) – one element of the chain every
thinelements will be considered.extended (bool) – if set to True, change the structure of the output to write a pkb extended. Optional, default
False.instr (str) – instrument to consider. Possible argument :
geometric,kepler,golf,virgo. Optional, defaultgeometric.verbose (bool) – if set to
True, print the name of the considered files. Optional, defaultFalse.
- Returns:
pkb array
- Return type:
ndarray
- apollinaire.peakbagging.chain_to_a2z(filename, thin=1, discard=0, add_ampl=False, instr='geometric')#
Function to create an a2z dataframe from a sampled chain.
- Parameters:
filename (str) – name of the chain. The name format is the following
mcmc_sampler_order_[n_order]_degrees_[n_degrees].h5withn_degreesbeing13or02.instr (str) – instrument to consider. Possible argument :
geometric,kepler,golf,virgo. Optional, defaultgeometric.
- Returns:
a2z style dataframe.
- Return type:
pandas DataFrame
Input/output files management functions#
- apollinaire.peakbagging.read_a2z(a2z_file)#
Read a file with a a2z standard syntax (doc to be written) and return a2z-style parameters (sorted by orders and degrees).
- Parameters:
a2z_file (string) – name of the file to read the parameters.
- Returns:
input parameters as a pandas DataFrame with the a2z syntax.
- Return type:
pandas DataFrame
- apollinaire.peakbagging.save_a2z(filename, df)#
Write with a a2z standard syntax (doc to be written).
- Parameters:
filename (string) – name of the file where to write the parameters.
- apollinaire.peakbagging.save_pkb(filename, pkb, author=None, spectro=False, extended=False, fmt=None, projected_splittings=False, nwalkers=None, nsteps=None, discard=None, fit_amp=False)#
Save pkb file with dedicated header.
- apollinaire.peakbagging.a2z_to_pkb(df_a2z, nopandas=True)#
Take a2z dataframe and return pkb_style parameters. Frequency are given in muHz.
- Parameters:
df_a2z (pandas DataFrame) – input parameters as a pandas DataFrame with the a2z syntax.
- Returns:
array under pkb format.
- Return type:
ndarray
- apollinaire.peakbagging.test_a2z(a2z_file, verbose=True)#
Test a a2z file to check that it is valid. The function checks the bounds set for the parameters, convert the a2z DataFrame to pkb
- Parameters:
a2z_file (str) – path of the a2z file to test.
- Returns:
a2z DataFrame and pkb array.
- Return type:
tuple
Quality assurance functions#
- apollinaire.peakbagging.bayes_factor(freq, psd, back, df_a2z, n, strategy='pair', l02=True, dnu=None, do_not_use_dnu=True, size_window=None, thin=10, discard=0, instr='geometric', use_sinc=False, asym_profile='nigam-kosovichev', hdf5Dir='./', wdw=None, add_ampl=False, parallelise=False, fit_amp=False, projected_splittings=False, nodes=4)#
- Parameters:
freq (ndarray) – frequency vector
psd (ndarray) – real power vector
back (ndarray) – model background power vector.
df_a2z (pandas DataFrame) – a2z input DataFrame that was used for the fit. The DataFrame is only used to compute the window inside which the fit was performed.
n (int) – mode order for which the statistical test will evaluate the model.
strategy (str) – strategy that was used for the fit,
pairororder. Optional, defaultpair.l02 (bool) – set to True if strategy is
pairand if the pair to consider is even. Set to False if the fitted mode was odd.dnu (float.) – large separation of the modes. In this function, it will only be used to adapt the fitting window for fitted modes above 1200 muHz. If not given, it will be computed (if possible) with the input df_a2z. Optional, default
None.do_not_use_dnu (bool) – if set to
True, fitting window will be computed without using dnu value. Optional, defaultTrue.size_window (float) – size of the window that was used for the fit, it it was given manually to the
peakbaggingfunction. Optional, defaultNone.thin (int) – one element of the chain every
thinelements will be considered.discard (int) – the number of elements to ignore at the beginning of the chain.
use_sinc (bool) – if set to
True, mode profiles will be computed using cardinal sinus and not Lorentzians. No asymmetry term will be used if it is the case. Optional, defaultFalse.asym_profile (str) – depending on the chosen argument, asymmetric profiles will be computed following Korzennik 2005 (
korzennik) or Nigam & Kosovichev 1998 (nigam-kosovichev). Defaultnigam-kosovichev.hdf5Dir (str) – directory where the hdf5 file storing the MCMC should be read. Optional, default ‘./’.
wdw (ndarray) – Observation window of the considered timeseries. Optional, default
None.add_ampl (bool) – if set to True, standard values will be used to compute the mode height ratios (if heights have not been specified individually for each degree). Optional, default False.
parallelise (bool) – set to True to parallelise process. Optional, default False.
- Returns:
tuple of float, in the following order: probability of model with two modes, probability of model with only the strong mode, probaility of the H0 hypothesis, number of models on which the probability has been estimated.
- Return type:
ndarray
Plot functions#
- apollinaire.peakbagging.plot_from_param(param_pkb, freq, psd, back=None, wdw=None, smoothing=50, spectro=True, correct_width=1.0, show=False, filename=None, instr='geometric', use_sinc=False, asym_profile='korzennik', projected_splittings=False, extended_splittings=False, **kwargs)#
Plot the results of a fit according to an input given with a pkb format.
- Parameters:
param_pkb (ndarray) – parameters contained in the pkb files.
freq (ndarray) – frequency vector, must be given in muHz.
psd (ndarray) – real power vector of the observed data.
back (ndarray) – real power vector of the fitted background.
wdw (bool) – set to
Trueif the mode have been fitted using the sidelobes fitting method, defaultFalse.smoothing (int) – size of the rolling window used to smooth the psd in the plot.
spectro (bool) – set to
Trueif the instruments uses spectroscopy, set the units in m/s instead of ppm, defaultTrue.correct_width (float) – param to adjust the width of the Lorentzian if it has been manually modified during the fitting
show (bool) – automatically show the plot, default
False.instr (str) – instrument to consider (amplitude ratio inside degrees depend on geometry AND instrument and should be adaptated). Possible argument :
geometric,kepler,golf,virgo. Optional, defaultgeometric.use_sinc (bool) – if set to
True, mode profiles will be computed using cardinal sinus and not Lorentzians. No asymmetry term will be used if it is the case. Optional, defaultFalse.asym_profile (str) – depending on the chosen argument, asymmetric profiles will be computed following Korzennik 2005 (
korzennik) or Nigam & Kosovichev 1998 (nigam-kosovichev).
- Returns:
None
- apollinaire.peakbagging.banana_diagram(freq, psd, back, pkb, grid=None, n=30, k=30, angle_min=0, angle_max=90, split_min=0, split_max=2, figsize=(9, 6), shading='auto', marker_color='black', cmap='plasma', contour_color='black', marker='o', annotate_max=True, show_contour=True, param_wdw=None, instr='geometric', normalise_likelihood=True, use_sinc=False, asym_profile='nigam-kosovichev', fit_amp=False, projected_splittings=False, add_colorbar=False, levels=None, centile_max=95, centile_min=75, vmin=None, vmax=None, pcolormesh_options={}, contour_options={})#
Compute banana diagram (see Ballot et al 2006) for a given set of fitted parameters. Parameters are fixed, except for angles and splittings which are set to vary between 0 and 90 degrees and 0 and 2 muHz, respectively.
- Parameters:
freq (ndarray) – frequency vector.
psd (ndarray) – psd vector.
back (ndarray) – background vector.
pkb (ndarray) – pkb array that will be used to compute the model and the likelihood.
grid (ndarray) – if
gridis not None, no computation will be performed and the plot will directly be made with the array passed with this argument.n (int) – number of elements along angle axis. Optional, default 30.
k (int) – number of elements along splittings axis. Optional, default 30.
angle_min – minimal value for angle on the diagram. Optional, default 0.
angle_max – maximal value for angle on the diagram. Optional, default 90.
split_min – minimal value for splittings on the diagram. Optional, default 0.
figsize (tuple) – size of the matplotlib Figure. Optional, default (9,6)
split_max – maximal value for splittings on the diagram. Optional, default 2.
shading (str) – shading parameter. Optional, default
auto.marker_color (str) – color of markers. Optional, default
black.cmap (str or colormap instance.) – colormap. Optional, default
plasma.contour_color (str) – color of contours. Optional, default
black.marker (str) – marker type. Optional, default
o.annotate_max (bool) – if
True, coordinates of the maximum in the computed grid will be given as a legend. Optional, defaultTrue.show_contour (bool) – whether to show contour or not. Optional, default
True.param_wdw (ndarray.) – window parameters to apply on the mode model.
instr (str) – instrument used to collect the time series. Optional, default
geometric.normalise_likelihood (bool) – if set to
True, log likelihood will be divided by the maximum of its absolute value. Optional, defaultTrue.pcolormesh_options (dict) – ditionary that will be passed as the
kwargsfor matplotlib.axes.Axes.pcolormesh.contour_options (dict) – ditionary that will be passed as the
kwargsfor matplotlib.axes.Axes.contour.
- Returns:
tuple with likelihood-grid and figure
- Return type:
tuple
- apollinaire.peakbagging.echelle_diagram(freq, PSD, dnu, twice=False, fig=None, index=111, figsize=(16, 16), title=None, smooth=10, cmap='cividis', cmap_scale='linear', mode_freq=None, mode_freq_err=None, vmin=None, vmax=None, scatter_color='white', fmt='+', ylim=None, shading='gouraud', mfc='none', ms=20, index_offset=None, mec=None, xlabel=None, ylabel=None, **kwargs)#
Build the echelle diagram of a given PSD.
- Parameters:
freq (ndarray) – input vector of frequencies.
PSD (ndarray) – input vector of power. Must be of same size than freq.
dnu (float) – the large frequency separation use to cut slices into the diagram.
twice (bool) – slice using 2 x dnu instead of dnu, default False.
fig (matplotlib Figure) – figure on which the echelle diagram will be plotted. If
None, a new figure instance will be created. Optional, defaultNone.index (int) – position of the echelle diagram Axe in the figure. Optional, default
111.figsize (tuple) – size of the echelle diagram to plot.
title (str) – title of the figure. Optional, default
(16, 16)smooth (int) – size of the rolling window used to smooth the PSD. Default 10.
cmap (str) – select one available color map provided by matplotlib, default
cividiscmap_scale (str) – scale use for the colormap. Can be ‘linear’ or ‘logarithmic’. Optional, default ‘linear’.
mode_freq (ndarray or tuple of array) – frequency array of the modes to represent on the diagram. It can be single array or a tuple of array.
mode_freq_err (ndarray or tuple of array) – frequency uncertainty of the modes to represent on the diagram. It can be a single array or a tuple of array.
vmin (float) – minimum value for the colormap.
vmax (float) – maximum value for the colormap.
scatter_color (str) – color of the scatter point of the mode frequencies. Optional, default
white.fmt (str or tuple) – the format of the errorbar to plot. Can be a single string or a tuple of string with the same dimension that
mode_freq.ylim (tuple) – the y-bounds of the echelle diagram.
mew (float) – marker edge width. Optional, default 1.
markersize (float) – size of the markers used for the errorbar plot. Optional, default 10.
capsize (float) – length of the error bar caps. Optional, default 2.
- Returns:
the matplotlib Figure with the echelle diagram.
- Return type:
matplotlib Figure
Model building functions#
- apollinaire.peakbagging.compute_model(freq, param_pkb, param_wdw=None, correct_width=1.0, instr='kepler', use_sinc=False, asym_profile='nigam-kosovichev', fit_amp=False, projected_splittings=False, splitting_array=None)#
Compute a p-mode model from a given set of parameters.
- Parameters:
freq (ndarray) – frequency vector.
param_pkb (ndarray) – parameters contained in the pkb files.
param_wdw (ndarray) – parameters given by the analysis of the window.
correct_width (float) – param to adjust the width of the Lorentzian if it has been manually modified during the fitting
instr (str) – instrument to consider (amplitude ratio inside degrees depend on geometry AND instrument and should be adaptated). Possible argument :
geometric,kepler,golf,virgo.use_sinc (bool) – if set to
True, mode profiles will be computed using cardinal sinus and not Lorentzians. No asymmetry term will be used if it is the case. Optional, defaultFalse.asym_profile (str) – depending on the chosen argument, asymmetric profiles will be computed following Korzennik 2005 (
korzennik) or Nigam & Kosovichev 1998 (nigam-kosovichev).fit_amp (bool) – if set to
True, the function consider that it got amplitudes and not heights as input parameters. Optional, defaultFalse.projected_splittings (bool) – if set to
True, the function will consider that thesplitparameters of the input are projected splittings and will build the model consequently. Optional, defaultFalse.splitting_array (ndarray) – A 3d array storing the splitting information at the m-component level. It is generated from a dictionary by
construct_splitting_array. To be used for generation of synthetic spectra.
- Returns:
computed model
- Return type:
ndarray
- apollinaire.peakbagging.compute_model_from_pattern(freq, pattern, back, n_order=3, amp_l=[1.0, 1.5, 0.7, 2])#
Compute a p-mode model from global pattern parameters.
- Parameters:
freq (ndarray) – frequency vector.
pattern (ndarray) – pattern parameters. Pattern parameters have to be given in the following order:
eps,alpha,Dnu,numax,Hmax,Wenv,w,d02,b02,d01,b01,d13,b03. A first version of the guess and bounds array can be created through thecreate_pattern_guess_arrays.back (ndarray) – background vector to use for the model.
n_order (int) – number of orders to consider around
numaxto build the model. Optional, default3.
- Returns:
computed model
- Return type:
ndarray
- apollinaire.peakbagging.build_background(freq, param, n_harvey=2, apodisation=False, remove_gaussian=True, two_gaussian=False, asy_gaussian=False)#
Rebuild background vector from complete set of parameters.
- Parameters:
freq (ndarray) – frequency vector (in muHz).
param (array like) – 1D vector with backgrounds parameters given in the following order: param_harvey (3x``n_harvey``), param_plaw (2), param_gaussian (3), white noise (1)
n_harvey (int) – number of Harvey laws to use to build the background model. Optional, default 2.
apodisation (bool) – if set to
True, distort the model to take the apodisation into account. Optional, defaultFalse.remove_gaussian – if set to
True, the p-mode gaussian bump will not be taken into account when building the vector. Optional, defaultTrue.remove_gaussian – bool
- Returns:
background array
- Return type:
ndarray
Convenience functions#
- apollinaire.peakbagging.create_background_guess_arrays(freq, psd, r=1.0, m=1.0, teff=5770.0, dnu=None, numax=None, n_harvey=2, spectro=False, power_law=False, high_cut_plaw=20.0, return_labels=False, two_gaussian=False, asy_gaussian=False)#
Create guess and bound arrays for background fit.
- Parameters:
freq (ndarray) – frequency vector
psd – real power vector
n_harvey (int) – number of Harvey laws to use to build the background model. Optional, default 2. With more than two Harvey laws, it is strongly recommended to manually feed the
guessparameter.r (float) – stellar radius. Given in solar radius. Optional, default 1.
m (float) – stellar mass. Given in solar mass. Optional, default 1.
teff (float) – stellar effective temperature. Given in K. Optional, default 5770.
dnu (float) – large separation of the p-mode. If given, it will superseed the scaling law guess using
r,mandteff.numax (float) – maximum p-mode bump power frequency. If given, it will superseed the scaling law guess using
r,mandteff.power_law (bool) – set to
Trueto fit a power law on the background. Optional, defaultFalse.spectro (bool) – if set to
True, make the plot with unit consistent with radial velocity, else with photometry. Automated guess will also be computed consistently with spectroscopic measurements. Optional, defaultTrue.high_cut_plaw (float) – high cut in frequency to compute the power law guess. Optional, default 20.
return_labels (bool) – if set to
True, the label array will be returned among the returned tuple of arrays.
- Returns:
guess, low bounds, up bounds (and labels if
return_labels=True) arrays.- Return type:
tuple of arrays
- apollinaire.peakbagging.create_pattern_guess_arrays(Hmax, Wenv, r=1.0, m=1.0, teff=5770.0, dnu=None, numax=None, mcmc=True)#
Create guess and bound arrays for pattern fit.
- Parameters:
Hmax (float) – height at numax.
Wenv (float) – width of the p-mode gaussian.
r (float) – stellar radius. Given in solar radius. Optional, default 1.
m (float) – stellar mass. Given in solar mass. Optional, default 1.
teff (float) – stellar effective temperature. Given in K. Optional, default 5770.
dnu (float) – large separation of the p-mode. If given, it will superseed the scaling law guess using
r,mandteff.numax (float) – maximum p-mode bump power frequency. If given, it will superseed the scaling law guess using
r,mandteff.mcmc (bool) – if set to
True, will compute the guess values used for the MCMC exploration. Otherwise, will compute the guess values for the MLE optimisation. Optional, defaultTrue.
- Returns:
guess, low bounds and up bounds array
- Return type:
tuple of arrays
- apollinaire.peakbagging.hide_modes(freq, psd, pkb, back, l=1, n_width=4, remove_left=True, remove_zero_bin=True)#
Hide PSD region around modes of given l.
- Parameters:
freq (ndarray) – frequency vector
psd (ndarray) – power vector
back (ndarray) – background power vector
pkb (ndarray) – pkb array with modes parameters.
l (int) – mode degrees to remove. Optional, default
1.n_width (float) – frequency bins will be removed
n_widthfrom the central frequency of the remaining modes. Optional, default4.remove_left (bool) – if set to
True, remove a region Dnu/2 wide before the first mode given in the pkb. Optional, defaultTrue.remove_zero_bin (bool) – if set to
True, bins with zero values will be removed from returned vector. Optional, defaultTrue.
- apollinaire.peakbagging.remove_pattern(freq, psd, back, pkb, l=[0, 2])#
Remove a given pattern from the PSD by dividing the PSD by the computed model.
- Parameters:
freq (ndarray) – frequency vector
psd (ndarray) – power vector
back (ndarray) – background power vector
pkb (ndarray) – pkb array with modes parameters.
l (int or array-like) – mode degrees to remove. Optional, default
[0,2].