Parameter law added to use Rayleigh-Jeans approximation

This commit is contained in:
Lukas Klass 2020-10-06 11:37:05 +02:00
parent e210cf0fa2
commit f8f582f254
2 changed files with 38 additions and 5 deletions

View File

@ -21,7 +21,7 @@ Model a target as a black body of a given temperature and apparent magnitude.
.. code-block:: xml .. code-block:: xml
<target type="BlackBodyTarget" temp="5778" temp_unit="K" mag="10" mag_unit="mag" band="B" size="point"/> <target type="BlackBodyTarget" temp="5778" temp_unit="K" mag="10" mag_unit="mag" band="B" size="point" law="Planck"/>
Attributes: Attributes:
* | **temp:** float * | **temp:** float
@ -32,8 +32,10 @@ Attributes:
| The apparent magnitude of the black body in magnitudes. In case of None or magnitude per solid angle, an extended target is assumed. | The apparent magnitude of the black body in magnitudes. In case of None or magnitude per solid angle, an extended target is assumed.
* | **mag_unit:** str, *optional* = "mag" * | **mag_unit:** str, *optional* = "mag"
| The unit of the black body's magnitude. This has to be [``mag``, ``mag / arcsec**2``, ``mag / sr``]. The default is ``mag``. | The unit of the black body's magnitude. This has to be [``mag``, ``mag / arcsec**2``, ``mag / sr``]. The default is ``mag``.
* | **band:** str * | **band:** str, *optional*
| The band used for fitting the black body's flux density to Vega's flux density. This has to be one of [``U``, ``B``, ``V``, ``R``, ``I``, ``J``, ``H``, ``K``, ``L``, ``M``, ``N``]. | The band used for fitting the black body's flux density to Vega's flux density. This has to be one of [``U``, ``B``, ``V``, ``R``, ``I``, ``J``, ``H``, ``K``, ``L``, ``M``, ``N``].
* | **law:** str, *optional*
| The law used for the black body emission. This can be either ``Planck`` for using Planck's law or ``RJ`` for using the Rayleigh-Jeans approximation of Planck's law.
.. _filetarget: .. _filetarget:

View File

@ -2,6 +2,7 @@ from ..target.ATarget import ATarget
from ..SpectralQty import SpectralQty from ..SpectralQty import SpectralQty
import astropy.units as u import astropy.units as u
from astropy.modeling.models import BlackBody from astropy.modeling.models import BlackBody
from astropy.constants import c, k_B
from ...lib.logger import logger from ...lib.logger import logger
from ..Entry import Entry from ..Entry import Entry
from typing import Union from typing import Union
@ -26,7 +27,7 @@ class BlackBodyTarget(ATarget):
@u.quantity_input(wl_bins='length', temp=[u.Kelvin, u.Celsius], mag=[u.mag, u.mag / u.sr]) @u.quantity_input(wl_bins='length', temp=[u.Kelvin, u.Celsius], mag=[u.mag, u.mag / u.sr])
def __init__(self, wl_bins: u.Quantity, temp: u.Quantity = 5778 * u.K, mag: u.Quantity = None, def __init__(self, wl_bins: u.Quantity, temp: u.Quantity = 5778 * u.K, mag: u.Quantity = None,
band: str = "V"): band: str = "V", law: str = "Planck"):
""" """
Initialize a new black body point source Initialize a new black body point source
@ -41,6 +42,9 @@ class BlackBodyTarget(ATarget):
unit, an extended source will be assumed. unit, an extended source will be assumed.
band : str band : str
Band used for fitting the planck curve to a star of 0th magnitude. Can be one of [U, B, V, R, I, J, H, K]. Band used for fitting the planck curve to a star of 0th magnitude. Can be one of [U, B, V, R, I, J, H, K].
law : str
Which law to use for the calculation of the flux values. Can be either 'Planck' for using Planck's law or
'RJ' to use the Rayleigh-Jeans approximation.
Returns Returns
------- -------
@ -48,8 +52,13 @@ class BlackBodyTarget(ATarget):
if band.upper() not in self._band.keys(): if band.upper() not in self._band.keys():
logger.error("Band has to be one of '[" + ", ".join(list(self._band.keys())) + "]'") logger.error("Band has to be one of '[" + ", ".join(list(self._band.keys())) + "]'")
# Create blackbody model with given temperature # Create blackbody model with given temperature
bb = BlackBody(temperature=temp, scale=1 * u.W / (u.m ** 2 * u.nm * u.sr)) bb = None
if law.lower() == "planck":
bb = BlackBody(temperature=temp, scale=1 * u.W / (u.m ** 2 * u.nm * u.sr))
elif law.upper() == "RJ":
bb = self.__rayleigh_jeans_factory(temp)
else:
logger.error("Unknown law '" + law + "' for target type BlackBody.")
if mag is not None: if mag is not None:
# Calculate the correction factor for a star of 0th magnitude using the spectral flux density # Calculate the correction factor for a star of 0th magnitude using the spectral flux density
# for the central wavelength of the given band # for the central wavelength of the given band
@ -67,6 +76,24 @@ class BlackBodyTarget(ATarget):
# Initialize super class # Initialize super class
super().__init__(SpectralQty(wl_bins, sfd), wl_bins) super().__init__(SpectralQty(wl_bins, sfd), wl_bins)
@staticmethod
@u.quantity_input(temp=[u.Kelvin, u.Celsius])
def __rayleigh_jeans_factory(temp: u.Quantity):
"""
Create a lambda function for the Rayleigh-Jeans law
Parameters
----------
temp : u.Quantity
The temperature in Kelvins
Returns
-------
res : lambda
A lambda function for the Rayleigh-Jeans law with the variable lambda wavelength
"""
return lambda wl: (2 * c * k_B * temp / wl ** 4 / u.sr).to(u.W / (u.m ** 2 * u.nm * u.sr))
@staticmethod @staticmethod
def check_config(conf: Entry) -> Union[None, str]: def check_config(conf: Entry) -> Union[None, str]:
""" """
@ -94,3 +121,7 @@ class BlackBodyTarget(ATarget):
mes = conf.check_selection("band", ["U", "B", "V", "R", "I", "J", "H", "K", "L", "M", "N"]) mes = conf.check_selection("band", ["U", "B", "V", "R", "I", "J", "H", "K", "L", "M", "N"])
if mes is not None: if mes is not None:
return mes return mes
if hasattr(conf, "law"):
mes = conf.check_selection("law", ["Planck", "RJ"])
if mes is not None:
return mes