LSST Applications g063fba187b+cac8b7c890,g0f08755f38+6aee506743,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g1a2382251a+b4475c5878,g1dcb35cd9c+8f9bc1652e,g20f6ffc8e0+6aee506743,g217e2c1bcf+73dee94bd0,g28da252d5a+1f19c529b9,g2bbee38e9b+3f2625acfc,g2bc492864f+3f2625acfc,g3156d2b45e+6e55a43351,g32e5bea42b+1bb94961c2,g347aa1857d+3f2625acfc,g35bb328faa+a8ce1bb630,g3a166c0a6a+3f2625acfc,g3e281a1b8c+c5dd892a6c,g3e8969e208+a8ce1bb630,g414038480c+5927e1bc1e,g41af890bb2+8a9e676b2a,g7af13505b9+809c143d88,g80478fca09+6ef8b1810f,g82479be7b0+f568feb641,g858d7b2824+6aee506743,g89c8672015+f4add4ffd5,g9125e01d80+a8ce1bb630,ga5288a1d22+2903d499ea,gb58c049af0+d64f4d3760,gc28159a63d+3f2625acfc,gcab2d0539d+b12535109e,gcf0d15dbbd+46a3f46ba9,gda6a2b7d83+46a3f46ba9,gdaeeff99f8+1711a396fd,ge79ae78c31+3f2625acfc,gef2f8181fd+0a71e47438,gf0baf85859+c1f95f4921,gfa517265be+6aee506743,gfa999e8aa5+17cd334064,w.2024.51
LSST Data Management Base Package
Loading...
Searching...
No Matches
makeTransmissionCurves.py
Go to the documentation of this file.
1# This file is part of obs_decam.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (http://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <http://www.gnu.org/licenses/>.
21
22import os
23import astropy.io.fits
24import numpy as np
25
26from lsst.afw.image import TransmissionCurve
27from lsst.utils import getPackageDir
28from .decamFilters import DECAM_FILTER_DEFINITIONS
29
30__all__ = ("getDESSystemTransmission", "getDESAtmosphereTransmission")
31
32DATA_DIR = os.path.join(getPackageDir("obs_decam_data"), "decam", "transmission_curves")
33
34DECAM_BEGIN = "2012-09-12" # Initial date for DECam first light.
35
36
38 """Return a dictionary of TransmissionCurves describing the atmospheric
39 throughput at CTIO.
40
41 Dictionary keys are string dates (YYYY-MM-DD) indicating the beginning of
42 the validity period for the curve stored as the associated dictionary
43 value. The curve is guaranteed not to be spatially-varying.
44 """
45 table = astropy.io.fits.getdata(os.path.join(DATA_DIR, "des", "des_atm_std.fits"))
46
47 atm = TransmissionCurve.makeSpatiallyConstant(
48 throughput=table["throughput_atm"].astype(np.float64),
49 wavelengths=table["lambda"].astype(np.float64),
50 throughputAtMin=table["throughput_atm"][0],
51 throughputAtMax=table["throughput_atm"][-1],
52 )
53
54 return {DECAM_BEGIN: atm}
55
56
58 """Return a nested dictionary of TransmissionCurves describing the
59 system throughput (optics + filter + detector) at the location of
60 each detector.
61
62 Outer dictionary keys are string dates (YYYY-MM-DD). The next level
63 dictionary maps the physical filter name to another dict. The inner
64 dict is keyed by detector number.
65 """
66 # We have DES detector throughputs for the grizy bands.
67
68 bands = ["g", "r", "i", "z", "Y"]
69
70 filter_dict = {}
71 for band in bands:
72 for filter_def in DECAM_FILTER_DEFINITIONS:
73 if band == filter_def.band:
74 physical_filter = filter_def.physical_filter
75 break
76
77 table = astropy.io.fits.getdata(
78 os.path.join(DATA_DIR, "des", f"{band}_band_per_detector_throughput.fits"),
79 )
80
81 detector_dict = {}
82 for index in range(table['throughput_ccd'].shape[1]):
83 # The DECam detector starts at 1.
84 detector = index + 1
85
86 tput = TransmissionCurve.makeSpatiallyConstant(
87 throughput=table["throughput_ccd"][:, index].astype(np.float64),
88 wavelengths=table["lambda"].astype(np.float64),
89 throughputAtMin=0.0,
90 throughputAtMax=0.0,
91 )
92
93 detector_dict[detector] = tput
94
95 filter_dict[physical_filter] = detector_dict
96
97 return {DECAM_BEGIN: filter_dict}