csst_cpic_sim package
Submodules
csst_cpic_sim.main module
- csst_cpic_sim.main.deduplicate_names_add_count(names: list)
remove duplicate names and add count
- csst_cpic_sim.main.main(argv=None)
Command line interface of csst_cpic_sim
- Parameters
argv (list) – input arguments. Default is None. If None, sys.argv is used.
- csst_cpic_sim.main.observation_simulation_from_config(obs_file, config_file)
Run observation simulation from config file
- Parameters
obs_file (str) – observation file or folder.
config_file (str) – config file.
Examples – see examples in example folder.
- csst_cpic_sim.main.quick_run_v2(target_str: str, band: str, expt: float, emgain: float, nframe: int, skybg: Optional[float] = None, rotation: float = 0, shift: list = [0, 0], emset_input: bool = False, cr_frame: bool = True, camera_effect: bool = True, prograss_bar=False, output='./') ndarray
A quick run function for CPIC.
- Parameters
target_str (str) – target string. See the input of target.target_file_load for more details.
band (str) – band name
expt (float) – exposure time in unit of second
emgain (float) – EM gain value. Note that emgain is not the same as emset, controled by emset_input keyword.
nframe (int) – number of frames
skybg (float or None) – sky background in unit of magnitude. If None, sky background will not be set.
rotation (float) – rotation of the telescope. in unit of degree. 0 means North is up.
shift (list) – target shift in unit of arcsec
emset_input (bool) – if True, emgain paramter is emset. Else, emgain is set using emgain parameter.
cr_frame (bool) – if True, cosmic ray frame will be added to the image.
camera_effect (bool) – if True, camera effect will be added to the image.
prograss_bar (bool) – if True, a prograss_bar will be shown.
output (str) – the output directory.
- Returns
image_cube
- Return type
np.ndarray
- csst_cpic_sim.main.vis_observation(target: dict, skybg: float, expt: float, nframe: int, band: str, emset: int, obsid: int = 41000000000, rotation: float = 0, shift: list = [0, 0], gnc_info: dict = {}, csst_format: bool = True, camera=<csst_cpic_sim.camera.CpicVisEmccd object>, crmaker=<csst_cpic_sim.camera.CosmicRayFrameMaker object>, nsample: int = 1, emgain=None, prograss_bar=None, output='./') ndarray
Observation simulation main process on visable band using EMCCD.
- Parameters
target (dict) – target dictionary. See the input of target.spectrum_generator for more details.
skybg (float) – sky background in unit of magnitude/arcsec^2
expt (float) – exposure time in unit of second
nframe (int) – number of frames
band (str) – band name
emset (int) – EM gain setting value. 1023(0x3FF) for ~1.0× EM gain.
obsid (int) – observation ID. Start from 4 for CPIC, 01 for science observation. See the input of io.obsid_parser for more details.
rotation (float) – rotation of the telescope. in unit of degree. 0 means North is up.
shift (list) – target shift in unit of arcsec
gnc_info (dict) – gnc_info dictionary. See the input of io.primary_hdu for more details.
csst_format (bool) – if True, the output fits file will be in the csst format.
crmaker (CosmicRayFrameMaker) – CosmicRayFrameMaker object. See the input of camera.CosmicRayFrameMaker for more details.
nsample (int) – number of samples for wide bandpass.
emgain (float or None) – if None, emgain are set using emset parameter. Else, emgain are set using emgain parameter.
prograss_bar (bool) – if True, a prograss_bar will be shown.
output (str) – the output directory.
- Returns
image_cube
- Return type
np.ndarray
csst_cpic_sim.target module
- class csst_cpic_sim.target.AlbedoCat(phase: float, metallicity: float, f_sed: float)
Bases:
TabularSpectralElement
,Icat
Generate albedo spectrum from the planet reflection model in Batalha et al. 2018
- Parameters
phase (float) – [degree], The phase angle of the planet, range from 0 to 180.
metallicity (float) – The metallicity of the planet. log(Z/Zsun), range from 0 to 2.
f_sed (float) – The sedimentation efficiency of the cloud. log(f_sed), range from -2 to 2. 2 is for cloud free case.
Notes
f_sed is only reliable between -2 and log10(6) and the cloud free case (2). Values between log10(6) and 2 are interpolated from the cloud free case (2) and log10(6).
Color Classification of Extrasolar Giant Planets: Prospects and Cautions Natasha E. Batalha et al 2018 AJ 156 158
- class csst_cpic_sim.target.TargetOjbect(info, cstar=None)
Bases:
object
A helper class to generate target spectrum and albedo spectrum
- x
x of target in arcsec
- Type
float
- y
y of target in arcsec
- Type
fload
- ra
ra string for center star, such as ‘15d23m05s’
- Type
str
- dec
dec string for center star
- Type
str
- distance
distance of center star in pc
- Type
float
- image
image of the target
- Type
2d np.array
- spectrum
spectrum of the target
- Type
pysynphot.spectrum.Spectrum
- csst_cpic_sim.target.bcc_spectrum(coe_cloud: float, coe_metalicity: float) ArraySpectralElement
Albedo spectrum of planet using BCC model (Batalha et al. 2018),
- Parameters
coe_cloud (float) – The sedimentation efficiency of the cloud. log(f_sed), range from -2 to 2. 2 is for cloud free case.
coe_metalicity (float) – The metallicity of the planet. log(Z/Zsun), range from 0 to 2.
- Returns
albedo spectrum of the planet
- Return type
pysynphot.spectrum.ArrayBandpass
- csst_cpic_sim.target.detect_template_path(template: str) str
Find template file in catalog folder or current folder.
- Parameters
template (str) – template file name
- Returns
absolute path of template file
- Return type
str
- csst_cpic_sim.target.extract_target_x_y(target: dict, ra0: Optional[str] = None, dec0: Optional[str] = None) tuple
extract x, y of target from target dict
- Parameters
target (dict) – target dict. must contain either (ra, dec) or (pangle, spearation)
ra0 (str) – ra of center star. must be provided if (ra, dec) of target is used
dec0 (str) – dec of center star. must be provided if (ra, dec) of target is used
- Returns
x, y – x, y of target in arcsec
- Return type
float
- Raises
ValueError – if (ra, dec) of target is used but (ra, dec) of center star is not provided.
ValueError – one of (ra, dec) or (pangle, spearation) is not provided.
- csst_cpic_sim.target.hybrid_albedo_spectrum(coe_b: float, coe_r: float) ArraySpectralElement
Albedo spectrum of planet using hybrid-jupiter-neptune model (Lacy et al. 2018) jupiter and neptune spectrum is from Karkoschka’s 1994
- Parameters
coe_b (float) – coefficient of blue spectrum, 1 for jupiter, 0 for neptune
coe_r (float) – coefficient of red spectrum, 1 for jupiter, 0 for neptune
- Returns
albedo spectrum of the planet
- Return type
pysynphot.spectrum.ArrayBandpass
- csst_cpic_sim.target.planet_contrast(planet_x_au: float, planet_y_au: float, phase_angle: float, radius: float) float
calculate the contrast of a planet
- Parameters
planet_x_au (float) – x position of the planet in au
planet_y_au (float) – y position of the planet in au
phase_angle (float) – phase angle of the planet in degree
radius (float) – radius of the planet in jupiter radius
- Returns
contrast – contrast of the planet
- Return type
float
- csst_cpic_sim.target.spectrum_generator(targets: dict) list
generate the spectrum due to the input target list
- Parameters
targets (dict) –
target dictionary which contains keyword ‘cstar’ (necessary), ‘stars’(optional), ‘planets’ (optional). The values are: - cstar: dict
center star information. must contain keywords ra, dec, distance, magnitude, sptype
- objects: list of dict, optional, recommended! V2.0 new!
list of targets. each dict must contain ra, dec, magnitude, sptype
- Returns
obj_sp_list – list of [x, y, spectrum, image] of each target
- Return type
list
- csst_cpic_sim.target.standard_spectrum(magnitude: float) ArraySourceSpectrum
Standard spectrum of magnitude system. For example, the standard_spectrum of vega megnitude is vega spectrum. The standard spectrum of ab magnitude is a flat spectrum.
- Parameters
magnitude (float) – magnitude of the standard spectrum
- Returns
star_sp
- Return type
S.spectrum.Spectrum
- csst_cpic_sim.target.star_photlam(magnitude: float, sptype: str, is_blackbody: bool = False, mag_input_band: str = 'f661') ArraySourceSpectrum
genrate flux spectrum of a star by its observal magnitude and spectral type
- Parameters
magnitude (float) – magnitude of the star
sptype (str) – spectral type of the star
is_blackbody (bool) – if True, use blackbody spectrum instead of ck04model
mag_input_band (str) – bandpass of the input magnitude, default is f661
- Returns
spectrum of the star in photlam unit
- Return type
pysynphot.spectrum.ArraySpectrum
- csst_cpic_sim.target.target_file_load(target: Union[dict, str]) dict
Generate target dict from file, string or dict.
- Parameters
target (dict) – target file name, formated string or target dict.
Outputs –
-------- –
target – dictionary of target. Format of input
Note
If target is a string start with *, it will be treated as a formatted string. e.g. “*5.1/25.3(1.3,1.5)/22.1(2.3,-4.5)” which means a central object with magnitude 5.1, and two substellar with magnitude 25.3 and 22.1, respectively. The spectrum of each object is standard refernece spectrum of ABmag. The first number in the parenthesis is the x position in arcsec, and the second is the y position.
If target is a string without *, it will be treated as file name. And find the file in catalog folder. The file need to be in yaml format. And end with .yaml (note .yml not work). If not .yaml will be added.
If target is a dict, it will be returned directly.
If all the above conditions are not met, an empty dict will be returned.
csst_cpic_sim.optics module
- csst_cpic_sim.optics.filter_throughput(filter_name)
Totally throughput of the each CPIC band. Including the throughput of the filter, telescope, cpic, and camera QE. If the filter_name is not supported, return the throughput of the default filter(f661).
- Parameters
filter_name (str) – The name of the filter. One of [‘f565’, ‘f661’(default), ‘f743’, ‘f883’, ‘f940’, ‘f1265’, ‘f1425’, ‘f1542’]
- Returns
The throughput of the filter.
- Return type
synphot.Bandpass
- csst_cpic_sim.optics.focal_convolve(band: str, targets: list, init_shifts: list = [0, 0], rotation: float = 0, nsample: int = 5, error: float = 0, platesize: list = [1024, 1024]) ndarray
PSF convolution of the ideal focus image.
- Parameters
band (str) – The name of the band.
target (list) – The list of thetargets. See the output of spectrum_generator for details.
init_shifts (list) – The shifts of the targets to simulate the miss alignment. Unit: arcsec
rotation (float) – The rotation of the image. Unit: degree
error (float) – The error of the DM acceleration. Unit: nm
platesize (list) – The size of the image. Unit: pixel
- Return type
np.ndarray
- csst_cpic_sim.optics.focal_mask(image, iwa, platescale, throughtput=1e-06)
Mask the image outside the inner working angle.
- Parameters
image (np.ndarray) – The image to be masked.
iwa (float) – The inner working angle. Unit: arcsec.
platescale (float) – The platescale of the image. Unit: arcsec/pixel.
throughtput (float) – The throughtput of the mask. The default is 1e-6.
- Returns
The masked image.
- Return type
np.ndarray
- csst_cpic_sim.optics.ideal_focus_image(bandpass: SpectralElement, targets: list, platescale, platesize: list = [1024, 1024], init_shifts: list = [0, 0], rotation: float = 0) ndarray
Ideal focus image of the targets. Each star is a little point of 1pixel.
- Parameters
bandpass (synphot.SpectralElement) – The bandpass of the filter.
targets (list) – The list of the targets. See the output of spectrum_generator for details.
platescale (float) – The platescale of the camera. Unit: arcsec/pixel
platesize (list) – The size of the image. Unit: pixel
init_shifts (list) – The shifts of the targets to simulate the miss alignment. Unit: arcsec
rotation (float) – The rotation of the image. Unit: degree
- Returns
The ideal focus image.
- Return type
np.ndarray
csst_cpic_sim.camera module
- class csst_cpic_sim.camera.CRobj(flux, angle, sigma, length)
Bases:
object
Cosmic ray object.
- flux
The flux of the cosmic ray.
- Type
float
- angle
The angle of the cosmic ray.
- Type
float
- sigma
The width of the cosmic ray.
- Type
float
- length
The length of the cosmic ray.
- Type
int
- class csst_cpic_sim.camera.CosmicRayFrameMaker
Bases:
object
Cosmic ray frame maker.
- Parameters
depth (float) – The depth of the camera pixel in um.
pitch (float) – The pitch of the camera pixel in um.
cr_rate (float) – The cosmic ray rate per second per cm2.
- make_CR(length, sigma, seed=- 1)
make a image of cosmic ray with given length and sigma.
- Parameters
length (int) – The length of the cosmic ray in pixel.
sigma (float) – The width of the cosmic ray in pixel.
- Returns
output – The image of cosmic ray.
- Return type
numpy.ndarray
- make_cr_frame(shape, expt, seed=- 1)
make a cosmic ray frame.
- Parameters
shape (list) – The size of the image in pixel.
expt (float) – The exposure time in second.
seed (int, optional) – The random seed. The default is -1. If seed is -1, the seed will be randomly selected.
- Returns
image – The cosmic ray frame.
- Return type
numpy.ndarray
- random_CR_parameter(expt, pixsize)
randomly generate cosmic ray parameters, including number, length, flux, sigma and angle.
- Parameters
expt (float) – The exposure time in second.
pixsize (list) – The size of the image in pixel.
- Returns
CRs – A list of cosmic ray objects.
- Return type
list
- class csst_cpic_sim.camera.CpicVisEmccd(config_dict=None)
Bases:
object
The class for Cpic EMCCD camera.
- switch
A dictionary to switch on/off the camera effects.
- Type
dict
- bias_frame()
Generate bias frame The bias frame contains vertical, horizontal, peper-salt noise, bias drift effect. Can be configurable using self.switch.
- Returns
bias frame
- Return type
np.ndarray
- config_init()
initialize the camera. If the config is set, call this function to update the config.
- em_fix_fuc_fit(emgain)
Calculate the emgain fix coeficient to fix the gamma distribution. The coeficient is from fixing of ideal emgain distribution.
- Parameters
emgain (float) – The emgain.
- Returns
The coeficient.
- Return type
float
- emgain_set(em_set, ccd_temp=None, self_update=True)
Set emgain from em set value.
- Parameters
em_set (int) – em set value. 3FF is about 1.17×, 200 is about 1000×.
ccd_temp (float, optional) – CCD temperature. If not given, use the current ccd temperature.
self_update (bool, optional) – if True, update the emgain and emset. Default is True. if False, only return the emgain.
- emregester_blooming(image, max_iteration=5)
emregester blooming effect
- nonlinear_effect(image)
nonlinear effect
- readout(image_focal, em_set, expt_set, image_cosmic_ray=False, emgain=None)
From focal planet image to ccd output. Interface function for emccd. Simulate the readout process.
- Parameters
image_focal (np.ndarray) – focal planet image. Unit: photon-electron/s/pixel
em_set (int) – emgain set value.
expt_set (float) – exposure time set value. Unit: s (0 mains 1ms)
image_cosmic_ray (np.ndarray, optional) – cosmic ray image. Unit: photon-electron/pixel
emgain (float, optional) – if not None, use the given emgain. Else, calculate the emgain using em_set
- Returns
ccd output image. Unit: ADU
- Return type
np.ndarray
- time_syn(t, readout=True, initial=False)
time synchronization and update the system time and ccd temperature
- Parameters
t (float) – relative time
readout (bool, optional) – If the camera is readout before time synchronization, set readout to True if True, the ccd temperature will increase, otherwise, it will decrease
initial (bool, optional) – If inital is True, the ccd will be intialized to the cooler temperature
- vertical_blooming(image)
vertical blooming effect
- csst_cpic_sim.camera.sky_frame_maker(band, skybg, platescale, shape)
generate a sky background frame.
- Parameters
band (str) – The band of the sky background.
skybg (str) – The sky background file name.
platescale (float) – The platescale of the camera in arcsec/pixel.
shape (tuple) – The shape of the output frame. (y, x)
- Returns
sky_bkg_frame – The sky background frame.
- Return type
numpy.ndarray
csst_cpic_sim.io module
- csst_cpic_sim.io.check_and_update_fits_header(header)
Check the header keywords and update the description according to the data model.
- Parameters
header (astropy.io.fits.header.Header) – The header to be checked.
- Return type
None
- csst_cpic_sim.io.datetime_obj_to_mjd(time_obj)
transfer datetime object to mean julian date (MJD).
- Parameters
time_obj (datetime.datetime) – The datetime object.
- Returns
The mean julian date (MJD).
- Return type
float
- csst_cpic_sim.io.frame_header(obs_info, index, primary_header, camera_dict)
Generate the header for a single frame.
- Parameters
obs_info (dict) – Dictionary of parameters. See save_fits function.
index (int) – Frame index.
primary_header (dict) – Primary header
- Return type
astropy.io.fits.Header
- csst_cpic_sim.io.obsid_parser(obsid: int)
Parse the obsid to get the obstype.
- Parameters
obsid (str) – The obsid of the observation. Obsid must be 11 digits and start with 5 for CPIC.
- Returns
The obstype of the observation.
- Return type
str
- Raises
ValueError – If the obsid is not 11 digits or does not start with 5.
- csst_cpic_sim.io.primary_hdu(obs_info: dict, gnc_info: dict, filename_output=False)
Generate the primary hdu of the fits file.
- Parameters
obs_info (dict) – The parameters of the observation. See save_fits function.
gnc_info (dict) – The gnc information of the observation.
filename_output (bool (optional)) – If True, the folder and the filename will be returned.
- Returns
fits.PrimaryHDU – The primary hdu of the fits file.
str, str – The folder and filename of the fits file. Only returned if filename_output is True.
Notes
The gnc_info dict should contain the information of orbit and observation. these informations are used to genrated the hdu header. Refer to the data model for more information.
- csst_cpic_sim.io.save_fits(images, obs_info, gnc_info, camera_dict={}, csst_format=True, output_folder='./')
Save the image to a fits file.
- Parameters
images (numpy.ndarray) – Image array to be saved.
obs_info (dict) –
Dictionary of observation informations. Must contain the following keys
- band: str
Band of the image.
- expt: float
Exposure time of the each image.
- nframe: int
Number of frames in the image.
- emgain: int
EM gain of the camera.
- obsid: str
Observation ID. Obsid must be 11 digits and start with 5 for CPIC. See pharse_obsid() for details.
- rotation: float
Rotation angle of the image.
- shift: list
Shift of the image.
gnc_info (dict) – Dictionary of GNC information. Contains the keywords in the primary header. See primary_hdu() for details.
csst_format (bool, optional) – Whether to save the fits file in CSST format, by default True.
- csst_cpic_sim.io.save_fits_simple(images, obs_info, output_folder='./')
Save the image to a fits file with a simple header to TMP directory.
- Parameters
images (numpy.ndarray) – Image array to be written.
obs_info (dict) – Dictionary of observation informations. See save_fits function.
- Return type
Filename of the saved fits file.
- csst_cpic_sim.io.set_up_logger(log_dir)
csst_cpic_sim.config module
- csst_cpic_sim.config.iso_time(time)
Transfer relative time to iso time format
- Parameters
time (str or float) – The relative time in seconds.
- Returns
The iso time format.
- Return type
str
- csst_cpic_sim.config.load_refdata_path(config_aim)
Load refdata path. The refdata path is stored in config_aim file. If not found, set it.
- Parameters
config_aim (str) – config_aim file path
- csst_cpic_sim.config.relative_time(time)
Transfer iso time format to relative time in seconds
- Parameters
time (str or float) – The iso time format.
- Returns
The relative time in seconds.
- Return type
float
- csst_cpic_sim.config.replace_cpism_refdata(config: dict, output: str = '$') None
Replace the cpism_refdata in the config. In the config file, we use ${cpism_refdata} to indicate the cpism_refdata. This function is used to replace the cpism_refdata in the config, or replace back.
- Parameters
config (dict) – config dict.
output (str) – ‘$’ or ‘other’. If output is ‘$’, then replace the cpism_refdata in the config with ${cpism_refdata}. If output is ‘other’, then replace the ${cpism_refdata} in the config file with the real path. ‘$’ is used meanly to generate a demo config file.
- csst_cpic_sim.config.setup_config(new_config)
Set up config from a dict. Some of the configs need to calcuate.
- Parameters
new_config (dict) – new config dict.
- Return type
None
- csst_cpic_sim.config.which_focalplane(band)
Return the name of the focalplane which the band belongs to. right now only support vis
- Parameters
band (str) – The name of the band.
- Returns
The name of the focalplane. ‘vis’ or ‘nir’ or ‘wfs’
- Return type
str
- Raises
ValueError – If the band is not in [‘f565’, ‘f661’, ‘f743’, ‘f883’, ‘f940’, ‘f1265’, ‘f1425’, ‘f1542’, ‘wfs’]
csst_cpic_sim.utils module
- class csst_cpic_sim.utils.Logger(filename, level='INFO')
Bases:
object
- csst_cpic_sim.utils.psf_imshow(psf, vmin=1e-08, vmax=0.1, log=True, region=1)
- csst_cpic_sim.utils.random_seed_select(seed=- 1)
Select a random seed for numpy.random and return it.
- csst_cpic_sim.utils.region_replace(background: ndarray, front: ndarray, shift: list, bmask: float = 1.0, fmask: float = 1.0, padded_in: bool = False, padded_out: bool = False, subpix: bool = False)
replace a region of the background with the front image.
- Parameters
background (np.ndarray) – The background image.
front (np.ndarray) – The front image.
shift (list) – The [x, y] shift of the front image. Unit: pixel. Relative to the lower-left corner of the background image. [0, 0] means the lower-left corner of the front image is at the lower-left corner of the background image.
bmask (float) – The mask of the background image. Default: 1.0 0.0 means the background image is masked. 1.0 means the background image is fully added.
fmask (float) – The mask of the front image. Default: 1.0 0.0 means the front image is masked (not added). 1.0 means the front image is fully added.
padded_in (bool) – Whether the input background image is padded. Default: False In the function, the background image is padded by the size of the front image. If True, means the background image is padded.
padded_out (bool) – Whether the output image is padded. Default: False In the function, the background image is padded by the size of the front image. If True, means the output image is padded. padded_in and padded_out are designed for the case that replace_region fuction is called multiple times.
subpix (bool) – Whether the shift is subpixel. Default: False If True, the shift is subpixel, using scipy.ndimage.shift to shift the front image. If False, the shift is integer, using numpy slicing to shift the front image.
- Returns
The output image. shape = background.shape if padded_out = False shape = background.shape + 2 * front.shape if padded_out = True
- Return type
np.ndarray
Module contents
- class csst_cpic_sim.CosmicRayFrameMaker
Bases:
object
Cosmic ray frame maker.
- Parameters
depth (float) – The depth of the camera pixel in um.
pitch (float) – The pitch of the camera pixel in um.
cr_rate (float) – The cosmic ray rate per second per cm2.
- make_CR(length, sigma, seed=- 1)
make a image of cosmic ray with given length and sigma.
- Parameters
length (int) – The length of the cosmic ray in pixel.
sigma (float) – The width of the cosmic ray in pixel.
- Returns
output – The image of cosmic ray.
- Return type
numpy.ndarray
- make_cr_frame(shape, expt, seed=- 1)
make a cosmic ray frame.
- Parameters
shape (list) – The size of the image in pixel.
expt (float) – The exposure time in second.
seed (int, optional) – The random seed. The default is -1. If seed is -1, the seed will be randomly selected.
- Returns
image – The cosmic ray frame.
- Return type
numpy.ndarray
- random_CR_parameter(expt, pixsize)
randomly generate cosmic ray parameters, including number, length, flux, sigma and angle.
- Parameters
expt (float) – The exposure time in second.
pixsize (list) – The size of the image in pixel.
- Returns
CRs – A list of cosmic ray objects.
- Return type
list
- csst_cpic_sim.extract_target_x_y(target: dict, ra0: Optional[str] = None, dec0: Optional[str] = None) tuple
extract x, y of target from target dict
- Parameters
target (dict) – target dict. must contain either (ra, dec) or (pangle, spearation)
ra0 (str) – ra of center star. must be provided if (ra, dec) of target is used
dec0 (str) – dec of center star. must be provided if (ra, dec) of target is used
- Returns
x, y – x, y of target in arcsec
- Return type
float
- Raises
ValueError – if (ra, dec) of target is used but (ra, dec) of center star is not provided.
ValueError – one of (ra, dec) or (pangle, spearation) is not provided.
- csst_cpic_sim.focal_mask(image, iwa, platescale, throughtput=1e-06)
Mask the image outside the inner working angle.
- Parameters
image (np.ndarray) – The image to be masked.
iwa (float) – The inner working angle. Unit: arcsec.
platescale (float) – The platescale of the image. Unit: arcsec/pixel.
throughtput (float) – The throughtput of the mask. The default is 1e-6.
- Returns
The masked image.
- Return type
np.ndarray
- csst_cpic_sim.planet_contrast(planet_x_au: float, planet_y_au: float, phase_angle: float, radius: float) float
calculate the contrast of a planet
- Parameters
planet_x_au (float) – x position of the planet in au
planet_y_au (float) – y position of the planet in au
phase_angle (float) – phase angle of the planet in degree
radius (float) – radius of the planet in jupiter radius
- Returns
contrast – contrast of the planet
- Return type
float
- csst_cpic_sim.quick_run_v2(target_str: str, band: str, expt: float, emgain: float, nframe: int, skybg: Optional[float] = None, rotation: float = 0, shift: list = [0, 0], emset_input: bool = False, cr_frame: bool = True, camera_effect: bool = True, prograss_bar=False, output='./') ndarray
A quick run function for CPIC.
- Parameters
target_str (str) – target string. See the input of target.target_file_load for more details.
band (str) – band name
expt (float) – exposure time in unit of second
emgain (float) – EM gain value. Note that emgain is not the same as emset, controled by emset_input keyword.
nframe (int) – number of frames
skybg (float or None) – sky background in unit of magnitude. If None, sky background will not be set.
rotation (float) – rotation of the telescope. in unit of degree. 0 means North is up.
shift (list) – target shift in unit of arcsec
emset_input (bool) – if True, emgain paramter is emset. Else, emgain is set using emgain parameter.
cr_frame (bool) – if True, cosmic ray frame will be added to the image.
camera_effect (bool) – if True, camera effect will be added to the image.
prograss_bar (bool) – if True, a prograss_bar will be shown.
output (str) – the output directory.
- Returns
image_cube
- Return type
np.ndarray
- csst_cpic_sim.sky_frame_maker(band, skybg, platescale, shape)
generate a sky background frame.
- Parameters
band (str) – The band of the sky background.
skybg (str) – The sky background file name.
platescale (float) – The platescale of the camera in arcsec/pixel.
shape (tuple) – The shape of the output frame. (y, x)
- Returns
sky_bkg_frame – The sky background frame.
- Return type
numpy.ndarray
- csst_cpic_sim.spectrum_generator(targets: dict) list
generate the spectrum due to the input target list
- Parameters
targets (dict) –
target dictionary which contains keyword ‘cstar’ (necessary), ‘stars’(optional), ‘planets’ (optional). The values are: - cstar: dict
center star information. must contain keywords ra, dec, distance, magnitude, sptype
- objects: list of dict, optional, recommended! V2.0 new!
list of targets. each dict must contain ra, dec, magnitude, sptype
- Returns
obj_sp_list – list of [x, y, spectrum, image] of each target
- Return type
list
- csst_cpic_sim.star_photlam(magnitude: float, sptype: str, is_blackbody: bool = False, mag_input_band: str = 'f661') ArraySourceSpectrum
genrate flux spectrum of a star by its observal magnitude and spectral type
- Parameters
magnitude (float) – magnitude of the star
sptype (str) – spectral type of the star
is_blackbody (bool) – if True, use blackbody spectrum instead of ck04model
mag_input_band (str) – bandpass of the input magnitude, default is f661
- Returns
spectrum of the star in photlam unit
- Return type
pysynphot.spectrum.ArraySpectrum
- csst_cpic_sim.vis_observation(target: dict, skybg: float, expt: float, nframe: int, band: str, emset: int, obsid: int = 41000000000, rotation: float = 0, shift: list = [0, 0], gnc_info: dict = {}, csst_format: bool = True, camera=<csst_cpic_sim.camera.CpicVisEmccd object>, crmaker=<csst_cpic_sim.camera.CosmicRayFrameMaker object>, nsample: int = 1, emgain=None, prograss_bar=None, output='./') ndarray
Observation simulation main process on visable band using EMCCD.
- Parameters
target (dict) – target dictionary. See the input of target.spectrum_generator for more details.
skybg (float) – sky background in unit of magnitude/arcsec^2
expt (float) – exposure time in unit of second
nframe (int) – number of frames
band (str) – band name
emset (int) – EM gain setting value. 1023(0x3FF) for ~1.0× EM gain.
obsid (int) – observation ID. Start from 4 for CPIC, 01 for science observation. See the input of io.obsid_parser for more details.
rotation (float) – rotation of the telescope. in unit of degree. 0 means North is up.
shift (list) – target shift in unit of arcsec
gnc_info (dict) – gnc_info dictionary. See the input of io.primary_hdu for more details.
csst_format (bool) – if True, the output fits file will be in the csst format.
crmaker (CosmicRayFrameMaker) – CosmicRayFrameMaker object. See the input of camera.CosmicRayFrameMaker for more details.
nsample (int) – number of samples for wide bandpass.
emgain (float or None) – if None, emgain are set using emset parameter. Else, emgain are set using emgain parameter.
prograss_bar (bool) – if True, a prograss_bar will be shown.
output (str) – the output directory.
- Returns
image_cube
- Return type
np.ndarray