pyqz package

Submodules

pyqz.pyqz module

pyqz: a Python module to derive the ionization parameter and oxygen abundance of HII regions, by comparing their strong emission line ratios with MAPPINGS photoionization models.

Copyright (C) 2014-2016, F.P.A. Vogt


This file contains the master pyqz routines. See the dedicated website for more info:

http://fpavogt.github.io/pyqz

Updated April 2016, F.P.A. Vogt - frederic.vogt@alumni.anu.edu.au


This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

pyqz.pyqz.check_grid(ratios, coeffs=[[1, 0], [0, 1]], Pk=5.0, kappa=inf, struct=’pp’, sampling=1)

Creates the diagram of a given line ratio grid for rapid inspection, and identify the location of “wraps” (= fold-over ambiguous regions).

Args:
ratios: string

The line ratios defining the grid, e.g. ‘[NII]/[SII]+;[OIII]/Hb’

coeffs: list of list [default: [[1,0],[0,1]] ]

The different coefficients with which to mix the line ratios. The size of each sub-list must be equal to the number of line ratios involved. Used for projected 3D diagnostic grids.

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled grid ?

Returns:
out: list

A list of bad segments with the associated node coordinates.

pyqz.pyqz.get_global_qz(data, data_cols, which_grids, ids=None, qzs=[‘LogQ’, ‘Tot[O]+12’], Pk=5.0, kappa=inf, struct=’pp’, sampling=1, error_pdf=’normal’, srs=400, flag_level=2.0, KDE_method=’gauss’, KDE_qz_sampling=101j, KDE_do_singles=True, KDE_pickle_loc=None, verbose=True, nproc=1)

Derives the LogQ and [O]+12 values for a given set of line fluxes based on a given set of diagnostic grids.

This function only uses “valid” grids defined in pyqz_metadata.py.

For each set of line ratios, this function combines the estimates from invdividual diagnostics into a combined best estimate of LogQ and [O]+12. Observational errors are propagated via the joint probability density function, reconstructed via a Kernel Density Estimation from “srs” realizations of the measured line fluxes, constructed randomly according to the observational errors and probability distribution function.

Args:
data: numpy array of size Nx2*M

N sets of M line fluxes and errors. An error =-1 implies that the associated line flux is an upper limit.

data_cols: list of Mx2 strings

Description of the coumns in “data”, e.g. [‘[NII]+’,’stdHa’,…] Mind the MAPPINGS + pyqz conventions.

which_grids: list of strings

The list of the diagnostic grids to use for the estimations, e.g. [‘[NII]+/Ha;[OII]/Hb’, …]

ids: list [default: None]

An optional length-N list of numbers or strings to identify each data point (e.g. spaxel number, source ID, etc …)

qzs: list of strings [default: [‘LogQ’,’Tot[O]+12’]]

list of Q/Z values to compute: [‘LogQ’,’Tot[O]+12’] or [‘LogQ’,’gas[O]+12’]

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf]

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled MAPPINGS grid ?

error_pdf: string [default = ‘normal’]

The shape of the error function for the line fluxes. Currently, only ‘normal’ is supported.

srs: int [default: 400]

The number of random line fluxes generated to discretize (and reconstruct) the joint probability function. If srs = 0, no KDE is computed.

flag_level: float [default: 2]

A ‘flag’ is raised (in the output array) when the direct q/z estimates and the KDE q/z estimates (q_rs and z_rs) are offset by more than flag_level * standard_deviation. Might indicate trouble. A flag is always raised when a point is outside all of the grids (8), when the KDE PDF is multipeaked (9) or when srs = 0 (-1, no KDE computed).

KDE_method: string [default: ‘gauss’]

Whether to use scipy.stats.gaussian_kde (‘gauss’) or sm.nonparametric.KDEMultivariate (‘multiv’) to perform the 2-D Kernel Density Estimation.

KDE_qz_sampling: complex [default: 101j]

The sampling of the QZ plane for the KDE reconstruction.

KDE_do_singles: bool [default: True]

Weather to compute KDE for single diagnostics (in addition to the KDE for all the diagnostics together).

KDE_pickle_loc: string [default: None]

If specified, the location to save the pickle files generated for each input points. Useful for subsequent plotting of the output via pyqz_plots.plot_global_qz()

verbose: bool [default: True]

Do you want to read anything ?

nproc: int [default = 1]

Defines how many process to use. Set to -1 for using as many as possible.

Returns:
out: [P,r]

A list where P contains all the estimated Log(Q) and [O]+12 values, and r contains the columns name

pyqz.pyqz.get_global_qz_ff(fn, which_grids, qzs=[‘LogQ’, ‘Tot[O]+12’], Pk=5.0, kappa=inf, struct=’pp’, sampling=1, error_pdf=’normal’, srs=400, flag_level=2.0, KDE_method=’gauss’, KDE_qz_sampling=101j, KDE_do_singles=True, KDE_pickle_loc=False, decimals=5, missing_values=’$$$’, suffix_out=’_out’, verbose=True, nproc=1)

The get_global_qz-‘from a file’ function. Gets the line ratios from a csv file, and not directly as Python arrays.

Save to file the log(q) and 12+log(O/H) values for a given set of line fluxes based on a given set of diagnostic grids. Requires a specific input file, and can batch process several measurements automatically.

Args:
fn: string

The path+filename of the input file, in CSV format.

which_grids: list of strings

The list of the diagnostic grids to use for the estimations, e.g. [‘[NII]+/Ha;[OII]/Hb’, …]

qzs: list of strings [default: [‘LogQ’,’Tot[O]+12’]]

list of Q/Z values to compute: [‘LogQ’,’Tot[O]+12’] or [‘LogQ’,’gas[O]+12’]

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf]

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled MAPPINGS grid ?

error_pdf: string [default = ‘normal’]

The shape of the error function for the line fluxes. Currently, only ‘normal’ is supported.

srs: int [default: 400]

The number of random line fluxes generated to discretize (and reconstruct) the joint probability function. If srs = 0, no KDE is computed.

flag_level: float [default: 2]

A ‘flag’ is raised (in the output array) when the direct q/z estimates and the KDE q/z estimates (q_rs and z_rs) are offset by more than flag_level * standard_deviation. Might indicate trouble. A flag is always raised when a point is outside all of the grids (8), when the KDE PDF is multipeaked (9) or when srs = 0 (-1, no KDE computed).

KDE_method: string [default: ‘gauss’]

Whether to use scipy.stats.gaussian_kde (‘gauss’) or sm.nonparametric.KDEMultivariate (‘multiv’) to perform the 2-D Kernel Density Estimation.

KDE_qz_sampling: complex [default: 101j]

The sampling of the QZ plane for the KDE reconstruction.

KDE_do_singles: bool [default: True]

Weather to compute KDE for single diagnostics (in addition to the KDE for all the diagnostics together).

KDE_pickle_loc: string [default: None]

If specified, the location to save the pickle files generated for each input points. Useful for subsequent plotting of the output via pyqz_plots.plot_global_qz()

decimals: int [default: 5]

The number of decimals to print in the final CSV file.

missing_values: string [default: ‘$$$’]

The symbol used to mark missing values in the CSV file.

suffix_out: string [default = ‘_ggQZff-out’}

The string to add to the input filename to create the output filename.

verbose: bool [default: True]

Do you want to read anything ?

nproc: int [default = 1]

Defines how many process to use. Set to -1 for using as many as possible.

pyqz.pyqz.get_global_qz_singlespec((j, final_cols, data, data_cols, which_grids, ids, qzs, Pk, kappa, struct, sampling, error_pdf, srs, flag_level, KDE_method, KDE_qz_sampling, KDE_do_singles, KDE_pickle_loc, verbose))

The single-process function associated with get_global_qz().

This function is used to spawn as many process as input spectra sent by the user to get_global_qz(). See http://stackoverflow.com/a/3843313 for more details.

Notes:The data is fed as one big list. See get_global_qz() for details on each item.
pyqz.pyqz.get_global_qz_singlespec_init(q)

Some function required to get the multiprocessing working fine. See http://stackoverflow.com/a/3843313 for details.

pyqz.pyqz.get_grid(ratios, coeffs=[[1, 0], [0, 1]], Pk=5.0, kappa=inf, struct=’pp’, sampling=1)

Returns a given line ratio diagnostic grid generated using MAPPINGS for a given set of parameters.

Args:
ratios: string

The line ratios defining the grid, e.g. ‘[NII]/[SII]+;[OIII]/Hb’

coeffs: list of list [default: [[1,0],[0,1]] ]

The different coefficients with which to mix the line ratios. The size of each sub-list must be equal to the number of line ratios involved. Used for projected 3D diagnostic grids.

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled grid ?

Returns:
out: list

Returns a list [grid,grid_cols,metadata], where: grid [numpy array]: the diagnostic grid as a numpy array. grid_cols [list]: labels of each column inside grid. metadata [list]: basic info from the associated MAPPINGS simulations.

pyqz.pyqz.interp_qz(qz, data, ratios, coeffs=[[1, 0], [0, 1]], Pk=5.0, kappa=inf, struct=’pp’, sampling=1, method=’linear’)

The core function of pyqz.

Returns the ‘q’ or ‘z’ value for a given set of line ratios based on a given 2-D diagnostic grid.

Args:
qz: string

Which estimate to return; ‘LogQ’, ‘Tot[O]+12’ or ‘gas[O]+12’

data: list of numpy array

List of Arrays of the line ratio values. One array per line ratio.

ratios: string

The line ratios defining the grid, e.g. ‘[NII]/[SII]+;[OIII]/Hb’

coeffs: list of list [default: [[1,0],[0,1]] ]

The different coefficients with which to mix the line ratios. The size of each sub-list must be equal to the number of line ratios involved. Used for projected 3D diagnostic grids.

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled grid ?

method: string [default: ‘linear’]

‘linear’ or ‘cubic’ interpolation method of scipy.interpolate.griddata.

Returns:
out: numpy array

The computed estimates, with nans when the points land outside the grid.

Notes:

WARNING: using method = ‘cubic’ can lead to estimates slightly “outside” the grid area. Proceed with caution. Usually, it is better to use a resampled grid with method = ‘linear’.

pyqz.pyqz_metadata module

pyqz: a Python module to derive the ionization parameter and oxygen abundance of HII regions, by comparing their strong emission line ratios with MAPPINGS photoionization models.

Copyright (C) 2014-2017, F.P.A. Vogt


This file contains some generic metadata used throughout the pyqz module, including the version number, etc …

Updated May 2017, F.P.A. Vogt - frederic.vogt@alumni.anu.edu.au

pyqz.pyqz_plots module

pyqz: a Python module to derive the ionization parameter and oxygen abundance of HII regions, by comparing their strong emission line ratios with MAPPINGS photoionization models.

Copyright (C) 2014-2016, F.P.A. Vogt


This file contains tools for the pyqz routines to create pretty plots seamlessly.

Updated April 2016, F.P.A. Vogt - frederic.vogt@alumni.anu.edu.au

pyqz.pyqz_plots.get_plot_labels(rats, coeffs)

Create the x and y label strings for the line ratio diagrams. Make sure they look as pretty as can be !

Args:
rats: list of string

The list of line ratios involved.

coeffs: list of list

The mixing coeffificents of the line ratios.

Returns:
labelx,labely: string, string

The pretty labels

pyqz.pyqz_plots.plot_global_qz(pickle_fn, do_all_diags=True, show_plots=True, save_loc=None)

A plotting function designed to plot the pickle output from pyqz.get_global_qz().

Args:
pickel_fn: string

the path+filename of the pickle filed to process.

do_all_diags: bool [default: True]

Whether to construct the diagnostic grid plots for all individual diagnostics or not.

show_plots: bool [default: True]

Whether to display the plots or not.

save_loc: string [default: None][

If set, the location where the plots will be saved. If not set, no plots will be saved.

Notes:

You must set KDE_pickle_loc to the location of your choice when running pyqz.get_global_qz() and pyqz.get_global_qz_ff(), if you want the associated pickle file to be generated.

pyqz.pyqz_plots.plot_grid(ratios, coeffs=[[1, 0], [0, 1]], Pk=5.0, kappa=inf, struct=’pp’, sampling=1, color_mode=’Tot[O]+12’, figname=None, show_plot=True, data=None, interp_data=None)

Creates the diagram of a given line ratio grid for rapid inspection, including wraps (fold-over regions), and of certain line ratios, if “data” is specified.

Args:
ratios: string

The line ratios defining the grid, e.g. ‘[NII]/[SII]+;[OIII]/Hb’

coeffs: list of list [default: [[1,0],[0,1]] ]

The different coefficients with which to mix the line ratios. The size of each sub-list must be equal to the number of line ratios involved. Used for projected 3D diagnostic grids.

Pk: float [default: 5.0]

MAPPINGS model pressure. This value must match an existing grid file.

kappa: float [default: np.inf

The kappa value. This value must match an existing grid file.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions. This value must match an existing reference grid file.

sampling: int [default: 1]

Use a resampled grid ?

color_mode: string [default: ‘Tot[O]+12’]

Color the grid according to ‘Tot[O]+12’, ‘gas[O]+12’ or ‘LogQ’.

figname: string [default: None]

‘path+name+format’ to save the Figure to.

show_plot: bool [default: True]

Do you want to display the Figure ?

data: list of numpy array [default: None]

List of Arrays of the line ratio values. One array per line ratio.

interp_data: numpy array [default: ‘linear’]

interpolated line ratios (via pyqz.interp_qz)

pyqz.pyqz_tools module

pyqz: a Python module to derive the ionization parameter and oxygen abundance of HII regions, by comparing their strong emission line ratios with MAPPINGS photoionization models.

Copyright (C) 2014-2017, F.P.A. Vogt


This file contains general tools for the pyqz routines.

Updated April 2016, F.P.A. Vogt - frederic.vogt@alumni.anu.edu.au Updated June 2017, F.P.A. Vogt - frederic.vogt@alumni.anu.edu.au

pyqz.pyqz_tools.ccw(A, B, C)

Checks whether 3 points are listed in a counter-clockwise series or not.

Args:
A: 2-D numpy array

Defining point 1

B: 2-D numpy array

Defining point 2

C: 2-D numpy array

Defining point 3

Returns:
out: bool

True or False ?

Notes:

Adapted from Bryce Boe: http://bryceboe.com/2006/10/23/line-segment-intersection-algorithm/

pyqz.pyqz_tools.get_KDE_picklefn(ids, Pk, kappa, struct, sampling, KDE_method)

Construct the filename of the pickle file associated with a given run of pyqz.get_global_qz.

Args:
ids: string

a string to identify the data point (e.g. spaxel number, source ID, …)

Pk: float

MAPPINGS model pressure.

kappa: float

The kappa value.

struct: string

spherical (‘sph’) or plane-parallel (‘pp’) HII regions.

sampling: int

Use a resampled grid ?

KDE_method: string

Whether to use scipy.stats.gaussian_kde (‘gauss’) or sm.nonparametric.KDEMultivariate (‘multiv’) to perform the 2-D Kernel Density Estimation.

Returns:
out: string

The pickle filename.

pyqz.pyqz_tools.get_MVphotogrid_fn(Pk=5.0, calibs=’GCZO’, kappa=inf, struct=’pp’, sampling=1)

Returns the filename of a MAPPINGS V grid, given a set of parameters.

Args:
Pk: float [default: 5.0]

MAPPINGS model pressure.

calibs: string [default: ‘GCZO’]

Input models to the MAPPINGS simulations.

kappa: float [default: np.inf]

The kappa value.

struct: string [default: ‘pp’]

spherical (‘sph’) or plane-parallel (‘pp’) HII regions.

sampling: int [default: 1]

Grid resampling factor.

Returns:
out: string

The MV filename, incl. the full path.

pyqz.pyqz_tools.get_MVphotogrid_metadata(fn)

Returns the first 4 lines of a MAPPINGS V file (compiled via the Hawk script).

Args:
fn: string

The filename (inculding path if required) as string.

Returns:
out: dictionnary of strings

The various parameters of the MAPPINGS simulation, including date & id.

pyqz.pyqz_tools.resample_MVphotogrid(fn, sampling=2)

A function to interpolate the MAPPINGS grid into a finer grid. Uses Akyma splines (1-D) in the Q/Z planes (one axis at a time).

Args:
fn: string

filename (inculding path if required) as string of the MAPPINGS models to resample.

sampling: int [default: 2]

The refining factor, i.e. number of interpolated points between each native grid node.

Notes:

The resampled grid is directly saved to a dedicated file.

pyqz.pyqz_tools.run_awk_loop(Pks=[5.0], kappas=[‘inf’], structs=[‘pp’], ncpu=1, awk_loc=’.’)

Runs the awk script provided with MAPPINGS in a loop, to easily generate all the default MAPPINGS grid required by pyqz.

Args:
Pks: list [deafult: [‘5.0’]]

The values of Pk to compute the grids at.

kappas: list [default: [‘inf’]]

The values of kappa to compute the grids at.

struct: list [default: [‘pp’]]

The values of struct to compute the grids at.

ncpu: int [default: 1]

The number of cpus to use.

awk_loc: string [default: ‘.’]

The path to the awk script provided with MAPPINGS.

Notes:

This function is not intended for general use - proceed with caution !

pyqz.pyqz_tools.seg_intersect(A, B, C, D)

Checks whether 2 segments defined by 4 points intersect.

Args:
A: 2-D numpy array

Defining point 1

B: 2-D numpy array

Defining point 2

C: 2-D numpy array

Defining point 3

D: 2-D numpy array

Defining point 4

Returns:
out: bool

True or False ?

Notes:

Adapted from Bryce Boe: http://bryceboe.com/2006/10/23/line-segment-intersection-algorithm/

Module contents