DiskIntegratedPhaseFunc¶
- class sbpy.photometry.DiskIntegratedPhaseFunc(*args, radius=None, wfb=None, **kwargs)[source]¶
Bases:
Fittable1DModel
Base class for disk-integrated phase function model
Examples
Define a linear phase function with phase slope 0.04 mag/deg, and study its properties:
>>> # Define a disk-integrated phase function model >>> import numpy as np >>> import astropy.units as u >>> from astropy.modeling import Parameter >>> from sbpy.calib import solar_fluxd >>> from sbpy.photometry import DiskIntegratedPhaseFunc >>> >>> class LinearPhaseFunc(DiskIntegratedPhaseFunc): ... ... _unit = 'mag' ... H = Parameter() ... S = Parameter() ... ... @staticmethod ... def evaluate(a, H, S): ... return H + S * a ... >>> linear_phasefunc = LinearPhaseFunc(5 * u.mag, 0.04 * u.mag/u.deg, ... radius = 300 * u.km, wfb = 'V') >>> pha = np.linspace(0, 180, 200) * u.deg >>> with solar_fluxd.set({'V': -26.77 * u.mag}): ... mag = linear_phasefunc.to_mag(pha) ... ref = linear_phasefunc.to_ref(pha) ... geomalb = linear_phasefunc.geomalb ... phaseint = linear_phasefunc.phaseint ... bondalb = linear_phasefunc.bondalb >>> print('Geometric albedo is {0:.3}'.format(geomalb)) Geometric albedo is 0.0487 >>> print('Bond albedo is {0:.3}'.format(bondalb)) Bond albedo is 0.0179 >>> print('Phase integral is {0:.3}'.format(phaseint)) Phase integral is 0.367
Initialization with subclass of
DataClass
:The subclassed models can either be initialized by model parameters, or by subclass of
DataClass
. Below example uses theHG
model class.>>> from sbpy.photometry import HG >>> from sbpy.data import Phys, Orbit, Ephem >>> >>> # Initialize from physical parameters pulled from JPL SBDB >>> phys = Phys.from_sbdb('Ceres') >>> print(phys['targetname','H','G']) <QTable length=1> targetname H G mag str17 float64 float64 ----------------- ------- ------- 1 Ceres (A801 AA) 3.31 0.12 >>> m = HG.from_phys(phys) INFO: Model initialized for 1 Ceres (A801 AA). [sbpy.photometry.core] >>> print(m) Model: HG Inputs: ('x',) Outputs: ('y',) Model set size: 1 Parameters: H G mag ---- ---- 3.31 0.12 >>> print(m.meta['targetname']) 1 Ceres (A801 AA) >>> print(m.radius) 469.7 km >>> >>> # Initialize from orbital elements pulled from JPL Horizons that also >>> # contain the H and G parameters >>> elem = Orbit.from_horizons('Ceres') >>> print(elem['targetname','H','G']) <QTable length=1> targetname H G mag str17 float64 float64 ----------------- ------- ------- 1 Ceres (A801 AA) 3.33 0.12 >>> m = HG.from_phys(elem) INFO: Model initialized for 1 Ceres (A801 AA). [sbpy.photometry.core] >>> >>> # Failed initialization due to the lack of field 'G' >>> phys = Phys.from_sbdb('12893') >>> print('G' in phys.field_names) False >>> m = HG(data=phys) Traceback (most recent call last): File "<stdin>", line 1, in <module> KeyError: 'field G not available.'
Initialize DiskIntegratedPhaseFunc
- Parameters:
- radiusastropy.units.Quantity, optional
Radius of object. Required if conversion between magnitude and reflectance is involved.
- wfb
Quantity
,SpectralElement
, string Wavelengths, frequencies, or bandpasses. Bandpasses may be a filter name (string). Required if conversion between magnitude and reflectance is involved.
- **kwargsoptional parameters accepted by
astropy.modeling.Model.__init__()
Attributes Summary
Bond albedo
Geometric albedo
Phase integral
Methods Summary
from_obs
(obs, fitter[, fields, init])Instantiate a photometric model class object from data
from_phys
(phys, **kwargs)Initialize an object from
Phys
objectto_mag
(eph[, unit, append_results])Calculate phase function in magnitude
to_phys
()Wrap the model into a
sbpy.data.Phys
objectto_ref
(eph[, normalized, append_results])Calculate phase function in average bidirectional reflectance
Attributes Documentation
- bondalb¶
Bond albedo
- geomalb¶
Geometric albedo
- input_units = {'x': Unit("rad")}¶
- input_units_allow_dimensionless = {'x': True}¶
- phaseint¶
Phase integral
Methods Documentation
- classmethod from_obs(obs, fitter, fields='mag', init=None, **kwargs)[source]¶
Instantiate a photometric model class object from data
- Parameters:
- obs
DataClass
, dict_like If
DataClass
or dict_like, must contain'phaseangle'
or the equivalent names (seeDataClass
). If any distance (heliocentric and geocentric) is provided, then they will be used to correct magnitude to 1 au before fitting.- fitter
Fitter
The fitter to be used for fitting.
- fieldsstr or array_like of str
The field name or names in
obs
to be fitted. If an array_like str, then multiple fields will be fitted one by one and a model set will be returned. In this case,.meta['fields']
of the returned object contains the names of fields fitted.- initnumpy array,
Quantity
, optional The initial parameters for model fitting. Its first dimension has the length of the model parameters, and its second dimension has the length of
n_model
if multiple models are fitted.- **kwargsoptional parameters accepted by
fitter()
. Note that the magnitude uncertainty can also be supplied to the fit via
weights
keyword for all fitters provided byfitting
.
- obs
- Returns:
- Object of
DiskIntegratedPhaseFunc
subclass The best-fit model class object.
- Object of
Examples
>>> from sbpy.photometry import HG >>> from sbpy.data import Misc >>> from astropy.modeling.fitting import SLSQPLSQFitter >>> fitter = SLSQPLSQFitter() >>> obs = Misc.mpc_observations('Bennu') >>> hg = HG() >>> best_hg = hg.from_obs(obs, eph['mag'], fitter)
- classmethod from_phys(phys, **kwargs)[source]¶
Initialize an object from
Phys
object- Parameters:
- Returns:
- Object of
DiskIntegratedPhaseFunc
subclass The phase function model object
- Object of
Examples
Initialization with
Phys
. This example uses theHG
model class.>>> from sbpy.photometry import HG >>> from sbpy.data import Phys >>> >>> # Initialize from physical parameters pulled from JPL SBDB >>> phys = Phys.from_sbdb('Ceres') >>> print(phys['targetname','H','G']) <QTable length=1> targetname H G mag str17 float64 float64 ----------------- ------- ------- 1 Ceres (A801 AA) 3.31 0.12 >>> m = HG.from_phys(phys) INFO: Model initialized for 1 Ceres (A801 AA). [sbpy.photometry.core] >>> print(m) Model: HG Inputs: ('x',) Outputs: ('y',) Model set size: 1 Parameters: H G mag ---- ---- 3.31 0.12 >>> print(m.meta['targetname']) 1 Ceres (A801 AA) >>> print(m.radius) 469.7 km >>> >>> # Failed initialization due to the lack of field 'G' >>> phys = Phys.from_sbdb('12893') >>> print('G' in phys.field_names) False >>> m = HG.from_phys(phys) Traceback (most recent call last): File "<stdin>", line 1, in <module> KeyError: 'field G not available.'
- to_mag(eph, unit=None, append_results=False, **kwargs)[source]¶
Calculate phase function in magnitude
- Parameters:
- eph
Ephem
, numbers, iterables of numbers, or Quantity
IfEphem
or dict_like, ephemerides of the object that can include phase angle, heliocentric and geocentric distances via keywordsphase
,r
anddelta
. If float or array_like, then the phase angle of object. If any distance (heliocentric and geocentric) is not provided, then it will be assumed to be 1 au. If no unit is provided via typeQuantity
, then radians is assumed for phase angle, and au is assumed for distances.- unit
astropy.units.mag
,astropy.units.MagUnit
, optional The unit of output magnitude. The corresponding solar magnitude must be available either through
sun
module or set byset
.- append_resultsbool, optional
Controls the return of this method.
- **kwargsoptional parameters accepted by
- eph
- Returns:
Quantity
, array ifappend_results == False
Ephem
ifappend_results == True
- When
append_results == False
: The calculated magnitude will be - returned.
- When
append_results == True
: Ifeph
is aEphem
- object, then the calculated magnitude will be appended to
eph
as - a new column. Otherwise a new
Ephem
object is created to - contain the input
eph
and the calculated magnitude in two columns.
Examples
>>> import numpy as np >>> from astropy import units as u >>> from sbpy.photometry import HG >>> from sbpy.data import Ephem >>> ceres_hg = HG(3.34 * u.mag, 0.12) >>> # parameter `eph` as `~sbpy.data.Ephem` type >>> eph = Ephem.from_dict({'alpha': np.linspace( ... 0, np.pi * 0.9, 200) * u.rad, ... 'r': np.repeat(2.7 * u.au, 200), ... 'delta': np.repeat(1.8 * u.au, 200)}) >>> mag1 = ceres_hg.to_mag(eph) >>> # parameter `eph` as numpy array >>> pha = np.linspace(0, 170, 200) * u.deg >>> mag2 = ceres_hg.to_mag(pha)
- to_phys()[source]¶
Wrap the model into a
sbpy.data.Phys
object- Returns:
Phys
object
Examples
>>> import astropy.units as u >>> from sbpy.calib import solar_fluxd >>> from sbpy.photometry import HG >>> from sbpy.data import Phys >>> >>> # Initialize from physical parameters pulled from JPL SBDB >>> phys = Phys.from_sbdb('Ceres') >>> print(phys['targetname','radius','H','G']) <QTable length=1> targetname radius H G km mag str17 float64 float64 float64 ----------------- ------- ------- ------- 1 Ceres (A801 AA) 469.7 3.31 0.12 >>> m = HG.from_phys(phys) INFO: Model initialized for 1 Ceres (A801 AA). [sbpy.photometry.core] >>> m.wfb = 'V' >>> with solar_fluxd.set({'V': -26.77 * u.mag}): ... p = m.to_phys() >>> print(type(p)) <class 'sbpy.data.phys.Phys'> >>> p.table.pprint(max_width=-1) targetname diameter H G pv A km mag ----------------- -------- ---- ---- ------------------- -------------------- 1 Ceres (A801 AA) 939.4 3.31 0.12 0.07624470768627523 0.027779803126557152
- to_ref(eph, normalized=None, append_results=False, **kwargs)[source]¶
Calculate phase function in average bidirectional reflectance
- Parameters:
- eph
Ephem
, numbers, iterables of numbers, or Quantity
IfEphem
or dict_like, ephemerides of the object that can include phase angle, heliocentric and geocentric distances via keywordsphase
,r
anddelta
. If float or array_like, then the phase angle of object. If any distance (heliocentric and geocentric) is not provided, then it will be assumed to be 1 au. If no unit is provided via typeQuantity
, then radians is assumed for phase angle, and au is assumed for distances.- normalizednumber,
Quantity
The angle to which the reflectance is normalized.
- append_resultsbool
Controls the return of this method.
- **kwargsoptional parameters accepted by
- eph
- Returns:
Quantity
, array ifappend_results == False
Ephem
ifappend_results == True
- When
append_results == False
: The calculated reflectance will be - returned.
- When
append_results == True
: Ifeph
is aEphem
- object, then the calculated reflectance will be appended to
eph
as - a new column. Otherwise a new
Ephem
object is created to - contain the input
eph
and the calculated reflectance in two - columns.
Examples
>>> import numpy as np >>> from astropy import units as u >>> from sbpy.calib import solar_fluxd >>> from sbpy.photometry import HG >>> from sbpy.data import Ephem >>> ceres_hg = HG(3.34 * u.mag, 0.12, radius = 480 * u.km, wfb= 'V') >>> # parameter `eph` as `~sbpy.data.Ephem` type >>> eph = Ephem.from_dict({'alpha': np.linspace( ... 0, np.pi * 0.9, 200) * u.rad, ... 'r': np.repeat(2.7 * u.au, 200), ... 'delta': np.repeat(1.8 * u.au, 200)}) >>> with solar_fluxd.set({'V': -26.77 * u.mag}): ... ref1 = ceres_hg.to_ref(eph) ... # parameter `eph` as numpy array ... pha = np.linspace(0, 170, 200) * u.deg ... ref2 = ceres_hg.to_ref(pha)