Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use nlamba=1 for PSF generation in test to be a bit faster. #105

Merged
merged 1 commit into from
Mar 14, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
14 changes: 10 additions & 4 deletions romanisim/image.py
Original file line number Diff line number Diff line change
Expand Up @@ -475,7 +475,8 @@ def simulate_counts(metadata, objlist,
rng=None, seed=None,
ignore_distant_sources=10, usecrds=True,
webbpsf=True,
darkrate=None, flat=None):
darkrate=None, flat=None,
psf_keywords=dict()):
"""Simulate total counts in a single SCA.

This gives the total counts in an idealized instrument with no systematics;
Expand All @@ -499,6 +500,8 @@ def simulate_counts(metadata, objlist,
dark rate image to use (electrons / s)
flat : float or np.ndarray[float]
flat field to use
psf_keywords : dict
keywords passed to PSF generation routine

Returns
-------
Expand Down Expand Up @@ -538,7 +541,7 @@ def simulate_counts(metadata, objlist,
chromatic = True
psf = romanisim.psf.make_psf(sca, filter_name, wcs=imwcs,
chromatic=chromatic, webbpsf=webbpsf,
variable=True)
variable=True, **psf_keywords)
image = galsim.ImageF(roman.n_pix, roman.n_pix, wcs=imwcs, xmin=0, ymin=0)
SCA_cent_pos = imwcs.toWorld(image.true_center)
sky_level = roman.getSkyLevel(bandpass, world_pos=SCA_cent_pos,
Expand Down Expand Up @@ -667,7 +670,8 @@ def gather_reference_data(image_mod, usecrds=False):

def simulate(metadata, objlist,
usecrds=True, webbpsf=True, level=2, crparam=dict(),
persistence=None, seed=None, rng=None, **kwargs
persistence=None, seed=None, rng=None,
psf_keywords=dict(), **kwargs
):
"""Simulate a sequence of observations on a field in different bandpasses.

Expand Down Expand Up @@ -701,6 +705,8 @@ def simulate(metadata, objlist,
Random number generator to use
seed : int
Seed for populating RNG. Only used if rng is None.
psf_keywords : dict
Keywords passed to the PSF generation routine

Returns
-------
Expand Down Expand Up @@ -762,7 +768,7 @@ def simulate(metadata, objlist,
log.info('Simulating filter {0}...'.format(filter_name))
counts, simcatobj = simulate_counts(
image_mod.meta, objlist, rng=rng, usecrds=usecrds, darkrate=darkrate,
webbpsf=webbpsf, flat=flat)
webbpsf=webbpsf, flat=flat, psf_keywords=psf_keywords)
if level == 0:
im = dict(data=counts.array, meta=dict(image_mod.meta.items()))
else:
Expand Down
12 changes: 7 additions & 5 deletions romanisim/psf.py
Original file line number Diff line number Diff line change
Expand Up @@ -24,7 +24,7 @@


def make_one_psf(sca, filter_name, wcs=None, webbpsf=True, pix=None,
chromatic=False, **kw):
chromatic=False, oversample=4, **kw):
"""Make a PSF profile for Roman at a specific detector location.

Can construct both PSFs using galsim's built-in galsim.roman.roman_psfs
Expand All @@ -41,8 +41,11 @@ def make_one_psf(sca, filter_name, wcs=None, webbpsf=True, pix=None,
scale of image for webbpsf PSFs
pix : tuple (float, float)
pixel location of PSF on focal plane
oversample : int
oversampling with which to sample WebbPSF PSF
**kw : dict
Additional keywords passed to galsim.roman.getPSF
Additional keywords passed to galsim.roman.getPSF or webbpsf.calc_psf,
depending on whether webbpsf is set.

Returns
-------
Expand Down Expand Up @@ -75,8 +78,7 @@ def make_one_psf(sca, filter_name, wcs=None, webbpsf=True, pix=None,
wfi.detector = f'SCA{sca:02d}'
wfi.filter = filter_name
wfi.detector_position = pix
oversample = kw.get('oversample', 4)
psf = wfi.calc_psf(oversample=oversample)
psf = wfi.calc_psf(oversample=oversample, **kw)
# webbpsf doesn't do distortion
# calc_psf gives something aligned with the pixels, but with
# a constant pixel scale equal to wfi.pixelscale / oversample.
Expand Down Expand Up @@ -137,7 +139,7 @@ def make_psf(sca, filter_name, wcs=None, webbpsf=True, pix=None,
variable : bool
True if a variable PSF object is desired
**kw : dict
Additional keywords passed to galsim.roman.getPSF
Additional keywords passed to make_one_psf

Returns
-------
Expand Down
34 changes: 21 additions & 13 deletions romanisim/tests/test_image.py
Original file line number Diff line number Diff line change
Expand Up @@ -111,8 +111,10 @@ def test_make_l2():
def set_up_image_rendering_things():
im = galsim.Image(100, 100, scale=0.1, xmin=0, ymin=0)
filter_name = 'F158'
impsfgray = psf.make_psf(1, filter_name, webbpsf=True, chromatic=False)
impsfchromatic = psf.make_psf(1, filter_name, webbpsf=False, chromatic=True)
impsfgray = psf.make_psf(1, filter_name, webbpsf=True, chromatic=False,
nlambda=1) # nlambda = 1 speeds tests
impsfchromatic = psf.make_psf(1, filter_name, webbpsf=False,
chromatic=True)
bandpass = roman.getBandpasses(AB_zeropoint=True)['H158']
counts = 1000
fluxdict = {filter_name: counts}
Expand Down Expand Up @@ -170,7 +172,7 @@ def test_image_rendering():
sca = 1
pos = [50, 50]
impsfgray = psf.make_psf(sca, filter_name, webbpsf=True, chromatic=False,
pix=pos, oversample=oversample)
pix=pos, oversample=oversample, nlambda=1)
counts = 100000
fluxdict = {filter_name: counts}
psfcatalog = [catalog.CatalogObject(None, galsim.DeltaFunction(), fluxdict)]
Expand All @@ -180,7 +182,7 @@ def test_image_rendering():
wfi.detector_position = pos
# oversample = kw.get('oversample', 4)
# webbpsf doesn't do distortion
psfob = wfi.calc_psf(oversample=oversample)
psfob = wfi.calc_psf(oversample=oversample, nlambda=1)
psfim = psfob[1].data * counts
# PSF from WebbPSF
im = galsim.Image(101, 101, scale=wfi.pixelscale, xmin=0, ymin=0)
Expand Down Expand Up @@ -242,7 +244,7 @@ def test_image_rendering():
modim = mod(xx, yy)
from scipy.signal import fftconvolve
# convolve with the PSF
psfob = wfi.calc_psf(oversample=oversample, fov_pixels=45)
psfob = wfi.calc_psf(oversample=oversample, fov_pixels=45, nlambda=1)
psfim = psfob[0].data
modim = fftconvolve(modim, psfim, mode='same')
modim = np.sum(modim.reshape(-1, imsz, oversample), axis=-1)
Expand Down Expand Up @@ -439,7 +441,8 @@ def test_simulate_counts():
webbpsf=False, ignore_distant_sources=100)
im2 = image.simulate_counts(meta, graycat,
usecrds=False, webbpsf=True,
ignore_distant_sources=100)
ignore_distant_sources=100,
psf_keywords=dict(nlambda=1))
im1 = im1[0].array
im2 = im2[0].array
maxim = np.where(im1 > im2, im1, im2)
Expand Down Expand Up @@ -497,9 +500,10 @@ def test_simulate():
imdict['tabcatalog'][filter_name] = (
imdict['tabcatalog'][filter_name] / abfluxdict[filter_name])
l0 = image.simulate(meta, graycat, webbpsf=True, level=0,
usecrds=False)
usecrds=False, psf_keywords=dict(nlambda=1))
l0tab = image.simulate(
meta, imdict['tabcatalog'], webbpsf=True, level=0, usecrds=False)
meta, imdict['tabcatalog'], webbpsf=True, level=0, usecrds=False,
psf_keywords=dict(nlambda=1))
# seed = 0 is special and means "don't actually use a seed." Any other
# choice of seed gives deterministic behavior
# note that we have scaled down the size of the image to 100x100 pix
Expand All @@ -509,7 +513,8 @@ def test_simulate():
# i.e., there are 1600x too many CRs. Fine for unit tests?
rng = galsim.BaseDeviate(1)
l1 = image.simulate(meta, graycat, webbpsf=True, level=1,
crparam=dict(), usecrds=False, rng=rng)
crparam=dict(), usecrds=False, rng=rng,
psf_keywords=dict(nlambda=1))
peakloc = np.nonzero(l0[0]['data'] == np.max(l0[0]['data']))

# check that the location with the most flux is the location where the
Expand All @@ -528,11 +533,13 @@ def test_simulate():

rng = galsim.BaseDeviate(1)
l1_nocr = image.simulate(meta, graycat, webbpsf=True, level=1,
usecrds=False, crparam=None, rng=rng)
usecrds=False, crparam=None, rng=rng,
psf_keywords=dict(nlambda=1))
assert np.all(l1[0].data >= l1_nocr[0].data)
log.info('DMS221: Successfully added cosmic rays to an L1 image.')
l2 = image.simulate(meta, graycat, webbpsf=True, level=2,
usecrds=False, crparam=dict())
usecrds=False, crparam=dict(),
psf_keywords=dict(nlambda=1))
# throw in some CRs for fun
l2c = image.simulate(meta, chromcat, webbpsf=False, level=2,
usecrds=False)
Expand All @@ -542,7 +549,8 @@ def test_simulate():
# zap the whole frame, 100 seconds ago.
rng = galsim.BaseDeviate(1)
l1p = image.simulate(meta, graycat, webbpsf=True, level=1, usecrds=False,
persistence=persist, crparam=None, rng=rng)
persistence=persist, crparam=None, rng=rng,
psf_keywords=dict(nlambda=1))
# the random number gets instatiated from the same seed, but the order in
# which the numbers are generated is different so we can't guarantee, e.g.,
# that all of the new values are strictly greater than the old ones.
Expand Down Expand Up @@ -629,7 +637,7 @@ def test_reference_file_crds_match(level):
im, simcatobj = image.simulate(
metadata, cat, usecrds=True,
webbpsf=True, level=level,
rng=rng)
rng=rng, psf_keywords=dict(nlambda=1))

# Confirm that CRDS keyword was updated
assert im.meta.ref_file.crds.sw_version != '12.3.1'
Expand Down
12 changes: 7 additions & 5 deletions romanisim/tests/test_psf.py
Original file line number Diff line number Diff line change
Expand Up @@ -22,18 +22,20 @@ def local(self, *args, **kwargs):

def test_make_psf():
psfs = []
psfs.append(psf.make_psf(1, 'F087'))
psfs.append(psf.make_psf(2, 'F184', chromatic=False))
pkw = dict(nlambda=1)
psfs.append(psf.make_psf(1, 'F087', **pkw))
psfs.append(psf.make_psf(2, 'F184', chromatic=False, **pkw))
psfs.append(psf.make_psf(3, 'F184', webbpsf=False))
psfs.append(psf.make_psf(4, 'H158', webbpsf=False, chromatic=True))
psfs.append(psf.make_psf(5, 'F184', pix=(1000, 1000), webbpsf=False))
psfs.append(psf.make_psf(6, 'F184', pix=(1000, 1000), webbpsf=True))
psfs.append(psf.make_psf(7, 'F129', wcs=FakeWCS()))
psfs.append(psf.make_psf(6, 'F184', pix=(1000, 1000), webbpsf=True,
**pkw))
psfs.append(psf.make_psf(7, 'F129', wcs=FakeWCS(), **pkw))
chromatic = [False] * 7
chromatic[3] = True
bandpass = galsim.roman.getBandpasses(AB_zeropoint=True)['H158']
vega_sed = galsim.SED('vega.txt', 'nm', 'flambda')
varpsf = psf.make_psf(8, 'F087', webbpsf=True, variable=True)
varpsf = psf.make_psf(8, 'F087', webbpsf=True, variable=True, **pkw)
psfs.append(varpsf.at_position(100, 100))
for p, c in zip(psfs, chromatic):
if not c:
Expand Down