fieldstars: module of field stars tools

Tools for field stars

author:
  1. Bednarski

license:

GNU GPL v3.0 (https://github.com/danmoser/pyhdust/blob/master/LICENSE)

pyhdust.fieldstars.binData(objarr, xarr, yarr, zarr=None, ignx=False, igny=False, prevent=False)[source]

Bin data

objarr is a list [[objects], [filters], [flags]]

xarr is a list [[x values], [sx values]] (idem for yarr and zarr)

If zarr is void, bin just the y values and calculate error by propagation only (no stddev). Otherwise, bin only the z values.

  • prevent: prevent to bin when flags contain ‘no-std’? If True,

    case there are more than one point able to be binned, bin only the data which don’t have ‘no-std’ flag. The output lists will have a copy of every line from those in input lists with ‘no-std’ flag, as like as the binned lines among those that just don’t have ‘no-std’ flag.

CAUTION: only bins when the contents of objarr AND xarr are the same (AND yarr also, case zarr != None)! Only if ignx or/and igny parameter have value ‘True’ the routine ignores when xarr or/and yarr have not the same values. In this case, the x or/and y value to be returned are taken from the first occurrence.

CONSIDERATIONS:

  • The error of binned data is just the propagated error, sqrt(sum(sigma_i^2))/n, and doesn’t consider the stddev

pyhdust.fieldstars.fitSerk(larr, parr, sarr, star='', law='w82', n_burnin=300, n_mcmc=800, n_walkers=120, extens='pdf')[source]

Fit Serkowski law using Markov Chain Monte Carlo from emcee code. The likelihood function (L) supposes Gaussian errors around the Pmax values:

log(L) = -0.5*chi2 -0.5*sum(ln(2*pi*sy^2))

INPUT:

larr: array/list with lambda values parr: array/list with P values sarr: array/list with the sigma_P values star: star name to be printed in the graph and

its filename. If it’s a void str ‘’, this routine give a random number to prevent overwriting data.

law: what K value in Serkowski’s law use?

(see polt.serkowski).

n_burnin: number of iterations for burning-in

n_mcmc: number of iterations to run emcee

n_walkers: number of walkers to map the posterior

probabilities.

extens: extension for the graph file

OUTPUT: sorted like “pmax_fit, lmax_fit, chi2”

pmax_fit: [Pmax, sPmax_+, sPmax_-], the Pmax value and

its errors (at right and left from it). Pmax is the median of distribution probability and sPmax_+, sPmax_- are the range within which there are 68.3% of the points in such distribution.

lmax_fit: Idem, for lmax.

chi2: Reduced chi2.

pyhdust.fieldstars.fixName(star)[source]

Fix the name of star ‘star’, returning the printable name.

‘hd’ -> ‘HD ‘ ‘2MASS-’ -> ‘’ ‘hip’ -> ‘HIP ‘ ‘TYC-’ -> ‘TYC ‘

Return identical ‘star’ name if the name was not identified

pyhdust.fieldstars.genInt(be, path=None, vfilter=['no-std'], extens='pdf')[source]

Call polt.graf_qu() for Be star ‘be’ and save the intrinsic angles inside thet_int.csv.

path: path inside which the logfile ‘be’.log (out from

polt.genTarget) is located. The code will try to open thet_int.csv inside the current directory (.). If it was not found, it will be created. Otherwise, the routine will append a new line inside it, asking about to overwrite an eventual existing line for star ‘be’.

CONSIDERATIONS:

  • Propagates errors from standard star.

pyhdust.fieldstars.gen_color(csvfile, be, obj, onlyY=False)[source]

Receives a string ‘csvfile’ pointing to the csvfile, a Be star ‘be’ and the field star named ‘obj’. Returns a smart computed np array for the color of such field star.

Note that if you run this task with distinct values for ‘onlyY’, the color returned may be distinct also.

‘be’ parameter is necessary because some field stars are used for more than one Be star.

pyhdust.fieldstars.gencsv(csvin, path=None, skipdth=False, delta=3.5, epssig=2.0)[source]

Generate a csvfile with every observations for the standard stars listed in pyhdust/refs/pol_hip.txt.

Compute lambda_eff using color index and an mean airmass X=1.35. The error is the (lambda_eff(X=1.7)-lambda_eff(X=1.))/2

INPUT
csvin The csv file with the informations

about the standard stars. The columns are indentified through idx1[] dictionary.

path The path to the red directory, inside

which are the nights with the observations

skipdth print all target values anyway, even when

there are no standard star?

epssig sigP/P max for unpolarized target (sigP/P

up to epssig doesn’t need standard star case skipdth==False)

delta tolerance for the angle between the two beams

of calcite.

pyhdust.fieldstars.getTable(data, x, y, z=None, sx=None, sy=None, sz=None, vfilter=['no-std'], bin_data=True, onlyY=False, unbias='wk')[source]

Receive the list data with many collumns and return the columns concerning to the x, y, z quantities (and their errors sx, sy, sz). Returns also ‘objarr’, a list with object names, filters and validation status.

IMPORTANT: 1) Propagates error of delta theta factor for

‘q’, ‘u’, ‘p’ and ‘thet’ labels

  1. Apply unbias correction for P and theta values when specified through ‘unbias’ variable.

  2. Case some label is ‘q’, ‘u’, ‘p’ or ‘thet’, ‘no-std’ is not in vfilter list and bin_data==True, don’t bin the data of those which have ‘no-std’ flag to prevent to compute values WITH NO MEANING. Return various lines for every unbinnable data and a line for the others else binned data.

INPUT

x,y,z,sx,sy,sz Label (key) of idx2 dictionary. bin_data Bin the data in the observations for

the same object+filter+x variable? Note: Case y and/or z are ‘p’ or ‘thet’, compute the bins over the Stokes parameters and then return the values to ‘p’ or ‘thet’!

onlyY True = only selects the list with

a “Y” marked

unbias Estimator to unbias the data when ‘y’ or ‘z’
is equals to ‘p’ or ‘thet’. Three options:
  1. ‘ml’: Maximum Likelihood (K=1.41)

  2. ‘wk’: Wardle & Kronberg (K=1.0)

  3. ‘’ : None (K=0.)

where K is the factor to unbias (read the description of routines ‘unbiasData’ and ‘meanAngle’.

vfilter List whose elements are the labels (tags)

to be filtered from the output. It can be ‘full’, ‘prob’ or ‘comp’ for these pre-defined lists in polt.vfil dictionary.

To mount your own list vfilter, select the current tags:

# General tags for target/standard observation - ‘bad-mod’ bad modulation - ‘very-bad-mod’ very bad modulation - ‘incomp-mods’ some modulations are incompatible - ‘obs-prob’ some observational problem/error - ‘iagpol-prob’ polarimeter problem suspected - ‘other-prob’ another relevant problem - ‘obs!=pub’ std incompatible with the published

# Tags for standard status - ‘no-std’ no standard in the night - ‘oth-day-std’ standard from another day - ‘oth-dth’ delta theta estimated from another filter

pyhdust.fieldstars.graf_field(csvfile, be, vfilter=['no-std'], squared=True, save=False, bin_data=True, onlyY=False, extens='pdf', unbias='wk')[source]

Plot a field graph with polarization directions for Be star be through csvfile data table.

squared: use the same scale in both x and y axes?

pyhdust.fieldstars.graf_inst(logfile, mode=1, vfilter=['no-std'], save=False, extens='pdf')[source]

Plot a QU diagram for unpolarized standard stars in the logfile .log file (the outfile from polt.genTarget, appended by the name of the object at each line). Propagates error of standard star.

mode=1 plot BVRI graphs in the same figure; mode=2 plot UBVRI graphs in separeted figures.

Return 3 lists:

  1. qarr = [[mean Q], [propagated Q error], [Q stddev]]

  2. uarr = [[mean U], [propagated U error], [U stddev]]

  3. narr = [n] (the number of data used to compute the averages)

Where [mean Q], [propagated Q error], …, [n] are lists with the values for each filter.

pyhdust.fieldstars.graf_p(csvfile, be, thetfile=None, path=None, vfilter=[], save=False, bin_data=True, onlyY=False, every=False, rotate=False, fit=True, unbias='wk', law='w82', extens='pdf')[source]

Plot P x wavelength for star ‘be’ and operate over a /’be’_is.csv file. The field stars are read from ‘csvfile’.

Fix the polarization bias using Wardle & Kronberg formula.

‘csvfile’ : location of dados.csv (field stars data).

‘path’the path where is located the log file for star ‘be’

(out from polt.genTarget). If None, it is supposed inside the current directory (.).

‘bin_data’: bin data in graphs?

‘onlyY’use only the field stars with ‘Y’ marks (field

stars originally selected for Be star ‘be’)?

‘thetfile’file with the intrinsic angles (oufile from

fs.genInt). If ‘thetfile’ != None, plot the ISP component parallel to the disk of the Be star. It is done by rotating the data in QU diagram at the angle corresponding to the disk inclination.

‘unbias’Estimator to unbias the data when ‘y’ or ‘z’
is equals to ‘p’ or ‘thet’. Three options:
  1. ‘ml’: Maximum Likelihood (K=1.41)

  2. ‘wk’: Wardle & Kronberg (K=1.0)

  3. ‘’ : None (K=0.)

where K is the factor to unbias (read the description of routines ‘unbiasData’ and ‘meanAngle’.

‘every’use one intrinsic angle for each one filter to

obtain the // component? If every=False makes all data to use a mean value at the -4:-2 collums (22th to 24th) from ‘thetfile’

‘fit’: fit the ISP using MCMC? Case True, generate the

graphs and a file ./’be’_is.csv with the best values. Write the mean polarization angles inside this file also. Case fit==True and there exists this file, this routine will ask to the user if he wants to read the previous fitted parameters and plot the Serkowski curves using them, or if he wants to run MCMC again, overwriting this file. If there exists a ./’be’_is.csv file and some standard star was not fitted yet, this routine will do that and append a line to the csv file.

‘rotate’DON’T WORKING PROPERLY. Rotate the polarization of

field stars in QU diagram by the mean angle? A option to be explored to replace the unbias procedure.

CONSIDERATIONS:

  • Error for parallel component is the combination between

the propagated error and stddev. - Error of standard stars is propagated to the field stars - The polarization angles (PA) inside ‘be’_is.csv are binned

over Stokes parameters ALWAYS. The errors are computed by the simple propagation.

  • The mean PA for each star is computed using QU parameters in all filters, even when p/sigma_p<1 (when the individual PA error is equals to 51.96), because we are actually operating over QU space!

  • p/sigma_p inside ‘be’_is.csv is computed by the ratio of 28.65 and the PA error.

pyhdust.fieldstars.graf_pradial(csvfile, be, filt='pmax', vfilter=[], isfile=None, fit=False, bin_data=True, onlyY=True, save=False, extens='pdf', unbias='wk')[source]

Plot a field graph with polarization versus distance.

filtis the filter to plot in y axes - ‘pmax’,’u’,’b’,

‘v’,’r’,’i’. If ‘pmax’ use the Pmax values from ./’be’_is.csv file to plot.

isfile‘be’_is.csv file location (out from fs.graf_p).

If None, it is suposed ./’be’_is.csv

fit(only for filt==’pmax’) fit a line in graph? This

routine will not use in fitting the points whose values in 22th column inside isfile have a ‘0’ character (this column defines the points to be used in the adjust). The points not used will be marked with a ‘x’ in the graph.

csvfile : location of dados.csv.

Considerations:

  • If filt==’pmax’, don’t plot the data whose values in 21th

column inside isfile have a ‘0’ character.

  • If filt==’pmax’, don’t use in fitting the points whose

values in 22th column inside isfile have a ‘0’ character when fit=True. The points not used will be marked with a ‘x’ in the graph.

  • Skip the lines inside isfile commented with a ‘#’, or with

a void first element

pyhdust.fieldstars.meanAngle(q, u, sq, su, estim='wk')[source]

Return the mean angle and the error propagated

The computation is over QU parameters.

Unbias theta error using ‘estim’ estimator:

if p/sp <= K, s_theta = psi otherwise, s_theta = propagated error

where K is given by the estimator related to the ‘estim’ variable:

  1. ‘ml’ : Maximum Likelihood (K=1.41, psi=51.96)

  2. ‘wk’ : Wardle & Kronberg (K=1.0, psi=51.96)

  3. ‘’ : None (K=0, psi=51.96)

  4. ‘mts’: Maier, Tenzer & Santangelo (estimates

    from Bayesian analysis, psi=61.14)

pyhdust.fieldstars.meanAngle_star(csvfile, be, obj, filts='ubvri', vfilter=['no-std'], onlyY=False, estim='wk')[source]

Return the meanAngle for star ‘obj’ in field of Be ‘be’ computed in all filters specified in string format in the variable ‘filts’.

For example, if filts=’ubv’, this routine will compute the mean angle among UBV filters.

Unbias theta error using ‘estim’ estimator:

if p/sp <= K, s_theta = psi otherwise, s_theta = propagated error

where K is given by the estimator related to the ‘estim’ variable:

  1. ‘ml’ : Maximum Likelihood (K=1.41, psi=51.96)

  2. ‘wk’ : Wardle & Kronberg (K=1.0, psi=51.96)

  3. ‘’ : None (K=0, psi=51.96)

  4. ‘mts’: Maier, Tenzer & Santangelo (estimates

    from Bayesian analysis, psi=61.14)

CONSIDERATIONS:

  • Operations over QU parameters

  • Error from standard star is propagated

  • If estim != ‘’ and the final polarization p is smaller than its uncertainty multiplied by the estimator factor K, the error of theta has the value 51.96 or 61.14, depending of the model

pyhdust.fieldstars.readcsv(csvfile, be)[source]

Read lines of Be star be from cvsfile and return a list with all components:

[
[[star 1, filter 1],
[star 1, filter 2]],
[[star 2, filter 1],
[star 2, filter 2]],

]

Ignores the observations with ‘E’ flag.

pyhdust.fieldstars.rotQU(q, u, sq, su, ang, sang)[source]

Rotates lists/arrays q and u in QU diagram at an angle 2*(ang +- sang) (clockwise).

Look that if ‘ang’ is the mean polarization angle, this rotation will transpose all polarization to Q parameter. U parameter will have residual variations with respect to the 0.

Return q_rot, u_rot, sq_rot, su_rot

Todo: sometimes it’s better don’t use any error value for ang?

pyhdust.fieldstars.rotQUBe(be, thetfile, path=None, every=False, vfilter=['no-std'])[source]

Rotates the QU values for Be star ‘be’ in the intrinsic angle specified inside thetfile and computes the <Q’> and <U’> of parameters rotated, returning four lists: lbd (the wavelength values), <U’>, sigma U’ and stddev U’, with the values in each one of the filters UBVRI.

‘thetfile’the location (path+filename) of thet_int.csv

file with the intrinsic angles (out from fs.gen).

‘path’the path where is located the log file for star

‘be’. If None, is supposed inside ‘.’.

‘every’use one intrinsic angle for each one filter to

rotate them? If every=False makes all data to use a mean value at the -4:-2 collums (22th to 24th) from ‘thetfile’.

CONSIDERATIONS:

  • Propagates errors from standard star.

The rotated parameters are just the polarization components perpendicular and parallel to the orientation of the disk.

If every==True, uses the intrinsic angle of each filter to rotate its QU values (if some angle==0, skip this filter!); otherwise, use the same value specified in last 4 columns for every one of the 5 filters.

pyhdust.fieldstars.unbiasData(p, s, estim='wk')[source]

Unbias P data, appling the operation over input lists/arrays:

p^2 -> p^2 - K^2 s^2

INPUT

p: % polarization s: pol error

estim: estimative model:
  1. ‘ml’: Maximum Likelihood (K=1.41)

  2. ‘wk’: Wardle & Kronberg (K=1.0)