Coverage for python/lsst/sims/GalSimInterface/galSimInterpreter.py : 9%

Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1"""
2This file defines the following classes:
4GalSimInterpreter -- a class which takes objects passed by a GalSim Instance Catalog
5(see galSimCatalogs.py) and uses GalSim to write them to FITS images.
7GalSimDetector -- a class which stored information about a detector in a way that
8GalSimInterpreter expects.
9"""
10from __future__ import print_function
12import math
13from builtins import object
14import os
15import pickle
16import tempfile
17import gzip
18import numpy as np
19import astropy
20import galsim
21from lsst.sims.utils import radiansFromArcsec, observedFromPupilCoords
22from lsst.sims.GalSimInterface import make_galsim_detector, SNRdocumentPSF, \
23 Kolmogorov_and_Gaussian_PSF, LsstObservatory
25__all__ = ["make_gs_interpreter", "GalSimInterpreter", "GalSimSiliconInterpreter",
26 "ObjectFlags"]
29def make_gs_interpreter(obs_md, detectors, bandpassDict, noiseWrapper,
30 epoch=None, seed=None, apply_sensor_model=False,
31 bf_strength=1):
32 if apply_sensor_model:
33 return GalSimSiliconInterpreter(obs_metadata=obs_md, detectors=detectors,
34 bandpassDict=bandpassDict, noiseWrapper=noiseWrapper,
35 epoch=epoch, seed=seed, bf_strength=bf_strength)
37 return GalSimInterpreter(obs_metadata=obs_md, detectors=detectors,
38 bandpassDict=bandpassDict, noiseWrapper=noiseWrapper,
39 epoch=epoch, seed=seed)
41class GalSimInterpreter(object):
42 """
43 This is the class which actually takes the objects contained in the GalSim
44 InstanceCatalog and converts them into FITS images.
45 """
46 _observatory = LsstObservatory()
48 def __init__(self, obs_metadata=None, detectors=None,
49 bandpassDict=None, noiseWrapper=None,
50 epoch=None, seed=None):
52 """
53 @param [in] obs_metadata is an instantiation of the ObservationMetaData class which
54 carries data about this particular observation (telescope site and pointing information)
56 @param [in] detectors is a list of GalSimDetectors for which we are drawing FITS images
58 @param [in] bandpassDict is a BandpassDict containing all of the bandpasses for which we are
59 generating images
61 @param [in] noiseWrapper is an instantiation of a NoiseAndBackgroundBase
62 class which tells the interpreter how to add sky noise to its images.
64 @param [in] seed is an integer that will use to seed the random number generator
65 used when drawing images (if None, GalSim will automatically create a random number
66 generator seeded with the system clock)
67 """
69 self.obs_metadata = obs_metadata
70 self.epoch = epoch
71 self.PSF = None
72 self.noiseWrapper = noiseWrapper
74 if seed is not None:
75 self._rng = galsim.UniformDeviate(seed)
76 else:
77 self._rng = None
79 if detectors is None:
80 raise RuntimeError("Will not create images; you passed no detectors to the GalSimInterpreter")
82 self.detectors = detectors
84 self.detectorImages = {} # this dict will contain the FITS images (as GalSim images)
85 self.bandpassDict = bandpassDict
86 self.blankImageCache = {} # this dict will cache blank images associated with specific detectors.
87 # It turns out that calling the image's constructor is more
88 # time-consuming than returning a deep copy
89 self.checkpoint_file = None
90 self.drawn_objects = set()
91 self.nobj_checkpoint = 1000
93 self.centroid_base_name = None
94 self.centroid_handles = {} # This dict will contain the file handles for each
95 # centroid file where sources are found.
96 self.centroid_list = [] # This is a list of the centroid objects which
97 # will be written to the file.
99 def setPSF(self, PSF=None):
100 """
101 Set the PSF wrapper for this GalSimInterpreter
103 @param [in] PSF is an instantiation of a class which inherits from PSFbase and defines _getPSF()
104 """
105 self.PSF = PSF
107 def _getFileName(self, detector=None, bandpassName=None):
108 """
109 Given a detector and a bandpass name, return the name of the FITS file to be written
111 @param [in] detector is an instantiation of GalSimDetector
113 @param [in] bandpassName is a string i.e. 'u' denoting the filter being drawn
115 The resulting filename will be detectorName_bandpassName.fits
116 """
117 return detector.fileName+'_'+bandpassName+'.fits'
119 def _doesObjectImpingeOnDetector(self, xPupil=None, yPupil=None, detector=None,
120 imgScale=None, nonZeroPixels=None):
121 """
122 Compare an astronomical object to a detector and determine whether or not that object will cast any
123 light on that detector (in case the object is near the edge of a detector and will cast some
124 incidental light onto it).
126 This method is called by the method findAllDetectors. findAllDetectors will generate a test image
127 of an astronomical object. It will find all of the pixels in that test image with flux above
128 a certain threshold and pass that list of pixels into this method along with data characterizing
129 the detector in question. This method compares the pupil coordinates of those pixels with the pupil
130 coordinate domain of the detector. If some of those pixels fall inside the detector, then this method
131 returns True (signifying that the astronomical object does cast light on the detector). If not, this
132 method returns False.
134 @param [in] xPupil the x pupil coordinate of the image's origin in arcseconds
136 @param [in] yPupil the y pupil coordinate of the image's origin in arcseconds
138 @param [in] detector an instantiation of GalSimDetector. This is the detector against
139 which we will compare the object.
141 @param [in] nonZeroPixels is a numpy array of non-zero pixels from the test image referenced
142 above. nonZeroPixels[0] is their x coordinate (in pixel value). nonZeroPixels[1] is
143 ther y coordinate.
145 @param [in] imgScale is the platescale of the test image in arcseconds per pixel
146 """
148 if detector is None:
149 return False
151 xPupilList = radiansFromArcsec(np.array([xPupil + ix*imgScale for ix in nonZeroPixels[0]]))
152 yPupilList = radiansFromArcsec(np.array([yPupil + iy*imgScale for iy in nonZeroPixels[1]]))
154 answer = detector.containsPupilCoordinates(xPupilList, yPupilList)
156 if True in answer:
157 return True
158 else:
159 return False
161 def findAllDetectors(self, gsObject, conservative_factor=10.):
162 """
163 Find all of the detectors on which a given astronomical object might cast light.
165 Note: This is a bit conservative. Later, once we actually have the real flux, we
166 can figure out a better estimate for the stamp size to use, at which point some objects
167 might not get drawn. For now, we just use the nominal stamp size from GalSim, scaled
168 up by the factor `conservative_factor` (default=10).
170 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
171 carrying information about the object whose image is to be drawn
173 @param [in] conservative_factor is a factor (should be > 1) by which we scale up the
174 nominal stamp size. Brighter objects should use larger factors, but the default value
175 of 10 should be fairly conservative and not waste too many cycles on things off the
176 edges of detectors.
178 @param [out] outputString is a string indicating which chips the object illumines
179 (suitable for the GalSim InstanceCatalog classes)
181 @param [out] outputList is a list of detector instantiations indicating which
182 detectors the object illumines
184 @param [out] centeredObj is a GalSim GSObject centered on the chip
186 Note: parameters that only apply to Sersic profiles will be ignored in the case of
187 pointSources, etc.
188 """
190 # create a GalSim Object centered on the chip.
191 centeredObj = self.createCenteredObject(gsObject)
193 sizeArcsec = centeredObj.getGoodImageSize(1.0) # pixel_scale = 1.0 means size is in arcsec.
194 sizeArcsec *= conservative_factor
195 xmax = gsObject.xPupilArcsec + sizeArcsec/2.
196 xmin = gsObject.xPupilArcsec - sizeArcsec/2.
197 ymax = gsObject.yPupilArcsec + sizeArcsec/2.
198 ymin = gsObject.yPupilArcsec - sizeArcsec/2.
200 outputString = ''
201 outputList = []
203 # first assemble a list of detectors which have any hope
204 # of overlapping the test image
205 viableDetectors = []
206 for dd in self.detectors:
207 xOverLaps = min(xmax, dd.xMaxArcsec) > max(xmin, dd.xMinArcsec)
208 yOverLaps = min(ymax, dd.yMaxArcsec) > max(ymin, dd.yMinArcsec)
210 if xOverLaps and yOverLaps:
211 if outputString != '':
212 outputString += '//'
213 outputString += dd.name
214 outputList.append(dd)
216 if outputString == '':
217 outputString = None
219 return outputString, outputList, centeredObj
221 def blankImage(self, detector=None):
222 """
223 Draw a blank image associated with a specific detector. The image will have the correct size
224 for the given detector.
226 param [in] detector is an instantiation of GalSimDetector
227 """
229 # in order to speed up the code (by a factor of ~2), this method
230 # only draws a new blank image the first time it is called on a
231 # given detector. It then caches the blank images it has drawn and
232 # uses GalSim's copy() method to return copies of cached blank images
233 # whenever they are called for again.
235 if detector.name in self.blankImageCache:
236 return self.blankImageCache[detector.name].copy()
237 else:
238 image = galsim.Image(detector.xMaxPix-detector.xMinPix+1, detector.yMaxPix-detector.yMinPix+1,
239 wcs=detector.wcs)
241 self.blankImageCache[detector.name] = image
242 return image.copy()
244 def drawObject(self, gsObject, max_flux_simple=0, sensor_limit=0, fft_sb_thresh=None):
245 """
246 Draw an astronomical object on all of the relevant FITS files.
248 @param [in] gsObject is an instantiation of the GalSimCelestialObject
249 class carrying all of the information for the object whose image
250 is to be drawn
252 @param [in] max_flux_simple is ignored here. (Used by GalSimSiliconInterpreter)
254 @param [in] sensor_limit is ignored here. (Used by GalSimSiliconInterpreter)
256 @param [in] fft_sb_thresh is ignored here. (Used by GalSimSiliconInterpreter)
258 @param [out] outputString is a string denoting which detectors the astronomical
259 object illumines, suitable for output in the GalSim InstanceCatalog
260 """
261 object_flags = ObjectFlags()
262 object_flags.set_flag('no_silicon')
264 # find the detectors which the astronomical object illumines
265 outputString, \
266 detectorList, \
267 centeredObj = self.findAllDetectors(gsObject)
269 # Make sure this object is marked as "drawn" since we only
270 # care that this method has been called for this object.
271 self.drawn_objects.add(gsObject.uniqueId)
273 # Compute the realized object fluxes for each band and return
274 # if all values are zero in order to save compute.
275 fluxes = [gsObject.flux(bandpassName) for bandpassName in self.bandpassDict]
276 realized_fluxes = [galsim.PoissonDeviate(self._rng, mean=f)() for f in fluxes]
277 if all([f == 0 for f in realized_fluxes]):
278 object_flags.set_flag('skipped')
279 self._store_zero_flux_centroid_info(detectorList, fluxes, gsObject,
280 object_flags.value)
281 return outputString
283 if len(detectorList) == 0:
284 # there is nothing to draw
285 return outputString
287 self._addNoiseAndBackground(detectorList)
289 for bandpassName, realized_flux, flux in zip(self.bandpassDict, realized_fluxes, fluxes):
290 for detector in detectorList:
292 name = self._getFileName(detector=detector, bandpassName=bandpassName)
294 xPix, yPix = detector.camera_wrapper.pixelCoordsFromPupilCoords(gsObject.xPupilRadians,
295 gsObject.yPupilRadians,
296 detector.name,
297 self.obs_metadata)
299 # Set the object flux to the value realized from the
300 # Poisson distribution.
301 obj = centeredObj.withFlux(realized_flux)
303 obj.drawImage(method='phot',
304 gain=detector.photParams.gain,
305 offset=galsim.PositionD(xPix-detector.xCenterPix,
306 yPix-detector.yCenterPix),
307 rng=self._rng,
308 maxN=int(1e6),
309 image=self.detectorImages[name],
310 poisson_flux=False,
311 add_to_image=True)
313 # If we are writing centroid files, store the entry.
314 if self.centroid_base_name is not None:
315 centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId,
316 flux, realized_flux, xPix, yPix, object_flags.value,
317 gsObject.galSimType)
318 self.centroid_list.append(centroid_tuple)
320 # Because rendering FitsImage object types can take a long
321 # time for bright objects (>1e4 photons takes longer than ~30s
322 # on cori-haswell), force a checkpoint after each object is
323 # drawn.
324 force_checkpoint = ((gsObject.galSimType == 'FitsImage') and
325 realized_flux > 1e4)
326 self.write_checkpoint(force=force_checkpoint)
327 return outputString
329 def _store_zero_flux_centroid_info(self, detectorList, fluxes, gsObject, obj_flags_value):
330 if self.centroid_base_name is None:
331 return
332 realized_flux = 0
333 for bandpassName, flux in zip(self.bandpassDict, fluxes):
334 for detector in detectorList:
335 xPix, yPix = detector.camera_wrapper.pixelCoordsFromPupilCoords(gsObject.xPupilRadians,
336 gsObject.yPupilRadians,
337 detector.name,
338 self.obs_metadata)
339 centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId,
340 flux, realized_flux, xPix, yPix, obj_flags_value,
341 gsObject.galSimType)
342 self.centroid_list.append(centroid_tuple)
344 def _addNoiseAndBackground(self, detectorList):
345 """
346 Go through the list of detector/bandpass combinations and
347 initialize all of the FITS files we will need (if they have
348 not already been initialized)
349 """
350 for detector in detectorList:
351 for bandpassName in self.bandpassDict:
352 name = self._getFileName(detector=detector, bandpassName=bandpassName)
353 if name not in self.detectorImages:
354 self.detectorImages[name] = self.blankImage(detector=detector)
355 if self.noiseWrapper is not None:
356 # Add sky background and noise to the image
357 self.detectorImages[name] = \
358 self.noiseWrapper.addNoiseAndBackground(self.detectorImages[name],
359 bandpass=self.bandpassDict[bandpassName],
360 m5=self.obs_metadata.m5[bandpassName],
361 FWHMeff=self.
362 obs_metadata.seeing[bandpassName],
363 photParams=detector.photParams,
364 detector=detector)
366 self.write_checkpoint(force=True, object_list=set())
368 def drawPointSource(self, gsObject, psf=None):
369 """
370 Draw an image of a point source.
372 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
373 carrying information about the object whose image is to be drawn
375 @param [in] psf PSF to use for the convolution. If None, then use self.PSF.
376 """
377 if psf is None:
378 psf = self.PSF
379 return self._drawPointSource(gsObject, psf=psf)
381 def _drawPointSource(self, gsObject, psf=None):
382 if psf is None:
383 raise RuntimeError("Cannot draw a point source in GalSim "
384 "without a PSF")
385 return psf.applyPSF(xPupil=gsObject.xPupilArcsec, yPupil=gsObject.yPupilArcsec)
387 def drawSersic(self, gsObject, psf=None):
388 """
389 Draw the image of a Sersic profile.
391 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
392 carrying information about the object whose image is to be drawn
394 @param [in] psf PSF to use for the convolution. If None, then use self.PSF.
395 """
397 if psf is None:
398 psf = self.PSF
399 return self._drawSersic(gsObject, psf=psf)
401 def _drawSersic(self, gsObject, psf=None):
402 # create a Sersic profile
403 centeredObj = galsim.Sersic(n=float(gsObject.sindex),
404 half_light_radius=float(gsObject.halfLightRadiusArcsec))
406 # Turn the Sersic profile into an ellipse
407 centeredObj = centeredObj.shear(q=gsObject.minorAxisRadians/gsObject.majorAxisRadians,
408 beta=(0.5*np.pi+gsObject.positionAngleRadians)*galsim.radians)
410 # Apply weak lensing distortion.
411 centeredObj = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu)
413 # Apply the PSF.
414 if psf is not None:
415 centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec,
416 yPupil=gsObject.yPupilArcsec,
417 obj=centeredObj)
419 return centeredObj
421 def drawRandomWalk(self, gsObject, psf=None):
422 """
423 Draw the image of a RandomWalk light profile. In orider to allow for
424 reproducibility, the specific realisation of the random walk is seeded
425 by the object unique identifier, if provided.
427 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
428 carrying information about the object whose image is to be drawn
430 @param [in] psf PSF to use for the convolution. If None, then use self.PSF.
431 """
432 if psf is None:
433 psf = self.PSF
434 return self._drawRandomWalk(gsObject, psf=psf)
436 def _drawRandomWalk(self, gsObject, psf=None):
437 # Seeds the random walk with the object id if available
438 if gsObject.uniqueId is None:
439 rng = None
440 else:
441 rng = galsim.BaseDeviate(int(gsObject.uniqueId))
443 # Create the RandomWalk profile
444 centeredObj = galsim.RandomKnots(npoints=int(gsObject.npoints),
445 half_light_radius=float(gsObject.halfLightRadiusArcsec),
446 rng=rng)
448 # Apply intrinsic ellipticity to the profile
449 centeredObj = centeredObj.shear(q=gsObject.minorAxisRadians/gsObject.majorAxisRadians,
450 beta=(0.5*np.pi+gsObject.positionAngleRadians)*galsim.radians)
452 # Apply weak lensing distortion.
453 centeredObj = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu)
455 # Apply the PSF.
456 if psf is not None:
457 centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec,
458 yPupil=gsObject.yPupilArcsec,
459 obj=centeredObj)
461 return centeredObj
463 def drawFitsImage(self, gsObject, psf=None):
464 """
465 Draw the image of a FitsImage light profile.
467 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
468 carrying information about the object whose image is to be drawn
470 @param [in] psf PSF to use for the convolution. If None, then use self.PSF.
471 """
472 if psf is None:
473 psf = self.PSF
474 return self._drawFitsImage(gsObject, psf=psf)
476 def _drawFitsImage(self, gsObject, psf=None):
477 # Create the galsim.InterpolatedImage profile from the FITS image.
478 centeredObj = galsim.InterpolatedImage(gsObject.fits_image_file,
479 scale=gsObject.pixel_scale)
480 if gsObject.rotation_angle != 0:
481 centeredObj = centeredObj.rotate(gsObject.rotation_angle*galsim.degrees)
483 # Apply weak lensing distortion.
484 centeredObj = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu)
486 # Apply the PSF
487 if psf is not None:
488 centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec,
489 yPupil=gsObject.yPupilArcsec,
490 obj=centeredObj)
492 return centeredObj
494 def createCenteredObject(self, gsObject, psf=None):
495 """
496 Create a centered GalSim Object (i.e. if we were just to draw this object as an image,
497 the object would be centered on the frame)
499 @param [in] gsObject is an instantiation of the GalSimCelestialObject class
500 carrying information about the object whose image is to be drawn
502 Note: parameters that obviously only apply to Sersic profiles will be ignored in the case
503 of point sources
504 """
505 if psf is None:
506 psf = self.PSF
507 return self._createCenteredObject(gsObject, psf=psf)
509 def _createCenteredObject(self, gsObject, psf=None):
510 if gsObject.galSimType == 'sersic':
511 centeredObj = self._drawSersic(gsObject, psf=psf)
513 elif gsObject.galSimType == 'pointSource':
514 centeredObj = self._drawPointSource(gsObject, psf=psf)
516 elif gsObject.galSimType == 'RandomWalk':
517 centeredObj = self._drawRandomWalk(gsObject, psf=psf)
519 elif gsObject.galSimType == 'FitsImage':
520 centeredObj = self._drawFitsImage(gsObject, psf=psf)
522 else:
523 raise RuntimeError("Apologies: the GalSimInterpreter does not yet have a method to draw " +
524 gsObject.galSimType + " objects")
526 return centeredObj
528 def writeImages(self, nameRoot=None):
529 """
530 Write the FITS files to disk.
532 @param [in] nameRoot is a string that will be prepended to the names of the output
533 FITS files. The files will be named like
535 @param [out] namesWritten is a list of the names of the FITS files written
537 nameRoot_detectorName_bandpassName.fits
539 myImages_R_0_0_S_1_1_y.fits is an example of an image for an LSST-like camera with
540 nameRoot = 'myImages'
541 """
542 namesWritten = []
543 for name in self.detectorImages:
544 if nameRoot is not None:
545 fileName = nameRoot+'_'+name
546 else:
547 fileName = name
548 self.detectorImages[name].write(file_name=fileName)
549 namesWritten.append(fileName)
551 return namesWritten
553 def open_centroid_file(self, centroid_name):
554 """
555 Open a centroid file. This file will have one line per-object and the
556 it will be labeled with the objectID and then followed by the average X
557 Y position of the photons from the object. Either the true photon
558 position or the average of the pixelated electrons collected on a finite
559 sensor can be chosen.
560 """
562 visitID = self.obs_metadata.OpsimMetaData['obshistID']
563 file_name = self.centroid_base_name + str(visitID) + '_' + centroid_name + '.txt.gz'
565 # Open the centroid file for this sensor with the gzip module to write
566 # the centroid files in gzipped format. Note the 'wt' which writes in
567 # text mode which you must explicitly specify with gzip.
568 self.centroid_handles[centroid_name] = gzip.open(file_name, 'wt')
569 self.centroid_handles[centroid_name].write('{:15} {:>15} {:>15} {:>10} {:>10} {:>11} {:>15}\n'.
570 format('SourceID', 'Flux', 'Realized flux',
571 'xPix', 'yPix', 'flags', 'GalSimType'))
573 def _writeObjectToCentroidFile(self, detector_name, bandpass_name, uniqueId,
574 flux, realized_flux, xPix, yPix,
575 obj_flags_value, object_type):
576 """
577 Write the flux and the the object position on the sensor for this object
578 into a centroid file. First check if a centroid file exists for this
579 detector and, if it doesn't create it.
581 @param [in] detector_name is the name of the sensor the gsObject falls on.
583 @param [in] bandpass_name is the name of the filter used in this exposure.
585 @param [in] uniqueId is the unique ID of the gsObject.
587 @param [in] flux is the calculated flux for the gsObject in the given bandpass.
589 @param [in] realized_flux is the Poisson realization of the object flux.
591 @param [in] xPix x-pixel coordinate of object.
593 @param [in] yPix y-pixel coordinate of object.
595 @param [in] obj_flags_value is the bit flags for the object handling composed
596 as an integer.
598 @param [in] object_type is the gsObject.galSimType
599 """
601 centroid_name = detector_name + '_' + bandpass_name
603 # If we haven't seen this sensor before open a centroid file for it.
604 if centroid_name not in self.centroid_handles:
605 self.open_centroid_file(centroid_name)
607 # Write the object to the file
608 self.centroid_handles[centroid_name].write('{:<15} {:15.5f} {:15.5f} {:10.2f} {:10.2f} {:11d} {:>15}\n'.
609 format(uniqueId, flux, realized_flux, xPix, yPix,
610 obj_flags_value, object_type))
612 def write_centroid_files(self):
613 """
614 Write the centroid data structure out to the files.
616 This function loops over the entries in the centroid list and
617 then sends them each to be writen to a file. The
618 _writeObjectToCentroidFile will decide how to put them in files.
620 After writing the files are closed.
621 """
622 # Loop over entries
623 for centroid_tuple in self.centroid_list:
624 self._writeObjectToCentroidFile(*centroid_tuple)
626 # Now close the centroid files.
627 for name in self.centroid_handles:
628 self.centroid_handles[name].close()
630 def write_checkpoint(self, force=False, object_list=None):
631 """
632 Write a pickle file of detector images packaged with the
633 objects that have been drawn. By default, write the checkpoint
634 every self.nobj_checkpoint objects.
635 """
636 if self.checkpoint_file is None:
637 return
638 if force or len(self.drawn_objects) % self.nobj_checkpoint == 0:
639 # The galsim.Images in self.detectorImages cannot be
640 # pickled because they contain references to unpickleable
641 # afw objects, so just save the array data and rebuild
642 # the galsim.Images from scratch, given the detector name.
643 images = {key: value.array for key, value
644 in self.detectorImages.items()}
645 drawn_objects = self.drawn_objects if object_list is None \
646 else object_list
647 image_state = dict(images=images,
648 rng=self._rng,
649 drawn_objects=drawn_objects,
650 centroid_objects=self.centroid_list)
651 with tempfile.NamedTemporaryFile(mode='wb', delete=False,
652 dir='.') as tmp:
653 pickle.dump(image_state, tmp)
654 tmp.flush()
655 os.fsync(tmp.fileno())
656 os.chmod(tmp.name, 0o660)
657 os.rename(tmp.name, self.checkpoint_file)
659 def restore_checkpoint(self, camera_wrapper, phot_params, obs_metadata,
660 epoch=2000.0):
661 """
662 Restore self.detectorImages, self._rng, and self.drawn_objects states
663 from the checkpoint file.
665 Parameters
666 ----------
667 camera_wrapper: lsst.sims.GalSimInterface.GalSimCameraWrapper
668 An object representing the camera being simulated
670 phot_params: lsst.sims.photUtils.PhotometricParameters
671 An object containing the physical parameters representing
672 the photometric properties of the system
674 obs_metadata: lsst.sims.utils.ObservationMetaData
675 Characterizing the pointing of the telescope
677 epoch: float
678 Representing the Julian epoch against which RA, Dec are
679 reckoned (default = 2000)
680 """
681 if (self.checkpoint_file is None
682 or not os.path.isfile(self.checkpoint_file)):
683 return
684 with open(self.checkpoint_file, 'rb') as input_:
685 image_state = pickle.load(input_)
686 images = image_state['images']
687 for key in images:
688 # Unmangle the detector name.
689 detname = "R:{},{} S:{},{}".format(*tuple(key[1:3] + key[5:7]))
690 # Create the galsim.Image from scratch as a blank image and
691 # set the pixel data from the persisted image data array.
692 detector = make_galsim_detector(camera_wrapper, detname,
693 phot_params, obs_metadata,
694 epoch=epoch)
695 self.detectorImages[key] = self.blankImage(detector=detector)
696 self.detectorImages[key] += image_state['images'][key]
697 self._rng = image_state['rng']
698 self.drawn_objects = image_state['drawn_objects']
699 self.centroid_list = image_state['centroid_objects']
701 def getHourAngle(self, mjd, ra):
702 """
703 Compute the local hour angle of an object for the specified
704 MJD and RA.
706 Parameters
707 ----------
708 mjd: float
709 Modified Julian Date of the observation.
710 ra: float
711 Right Ascension (in degrees) of the object.
713 Returns
714 -------
715 float: hour angle in degrees
716 """
717 time = astropy.time.Time(mjd, format='mjd',
718 location=self.observatory.getLocation())
719 # Get the local apparent sidereal time.
720 last = time.sidereal_time('apparent').degree
721 ha = last - ra
722 return ha
724 @property
725 def observatory(self):
726 return self._observatory
729class GalSimSiliconInterpreter(GalSimInterpreter):
730 """
731 This subclass of GalSimInterpreter applies the Silicon sensor
732 model to the drawn objects.
733 """
734 def __init__(self, obs_metadata=None, detectors=None, bandpassDict=None,
735 noiseWrapper=None, epoch=None, seed=None, bf_strength=1):
736 super(GalSimSiliconInterpreter, self)\
737 .__init__(obs_metadata=obs_metadata, detectors=detectors,
738 bandpassDict=bandpassDict, noiseWrapper=noiseWrapper,
739 epoch=epoch, seed=seed)
741 self.gs_bandpass_dict = {}
742 for bandpassName in bandpassDict:
743 bandpass = bandpassDict[bandpassName]
744 index = np.where(bandpass.sb != 0)
745 bp_lut = galsim.LookupTable(x=bandpass.wavelen[index],
746 f=bandpass.sb[index])
747 self.gs_bandpass_dict[bandpassName] \
748 = galsim.Bandpass(bp_lut, wave_type='nm')
750 self.sky_bg_per_pixel = None
752 # Create a PSF that's fast to evaluate for the postage stamp
753 # size calculation for extended objects in .getStampBounds.
754 FWHMgeom = obs_metadata.OpsimMetaData['FWHMgeom']
755 self._double_gaussian_psf = SNRdocumentPSF(FWHMgeom)
757 # Save the parameters needed to create a Kolmogorov PSF for a
758 # custom value of gsparams.folding_threshold. That PSF will
759 # to be used in the .getStampBounds function for bright stars.
760 altRad = np.radians(obs_metadata.OpsimMetaData['altitude'])
761 self._airmass = 1.0/np.sqrt(1.0-0.96*(np.sin(0.5*np.pi-altRad))**2)
762 self._rawSeeing = obs_metadata.OpsimMetaData['rawSeeing']
763 self._band = obs_metadata.bandpass
765 # Save the default folding threshold for determining when to recompute
766 # the PSF for bright point sources.
767 self._ft_default = galsim.GSParams().folding_threshold
769 # Save these, which are needed for DCR
770 self.local_hour_angle \
771 = self.getHourAngle(self.obs_metadata.mjd.TAI,
772 self.obs_metadata.pointingRA)*galsim.degrees
773 self.obs_latitude = self.observatory.getLatitude().asDegrees()*galsim.degrees
775 # Make a trivial SED to use for faint things.
776 blue_limit = np.min([bp.blue_limit for bp in self.gs_bandpass_dict.values()])
777 red_limit = np.max([bp.red_limit for bp in self.gs_bandpass_dict.values()])
778 constant_func = galsim.LookupTable([blue_limit, red_limit], [1,1], interpolant='linear')
779 self.trivial_sed = galsim.SED(constant_func, wave_type='nm', flux_type='fphotons')
781 # Create SiliconSensor objects for each detector.
782 self.sensor = dict()
783 for det in detectors:
784 self.sensor[det.name] \
785 = galsim.SiliconSensor(strength=bf_strength,
786 treering_center=det.tree_rings.center,
787 treering_func=det.tree_rings.func,
788 transpose=True)
790 def drawObject(self, gsObject, max_flux_simple=0, sensor_limit=0, fft_sb_thresh=None):
791 """
792 Draw an astronomical object on all of the relevant FITS files.
794 @param [in] gsObject is an instantiation of the GalSimCelestialObject
795 class carrying all of the information for the object whose image
796 is to be drawn
798 @param [in] max_flux_simple is the maximum flux at which various simplifying
799 approximations are used. These include using a flat SED and possibly omitting
800 the realistic sensor effects. (default = 0, which means always use the full SED)
802 @param [in] sensor_limit is the limiting value of the existing flux in the
803 postage stamp image, above which the use of a SiliconSensor model is forced.
804 For faint things, if there is not already flux at this level, then a simple
805 sensor model will be used instead. (default = 0, which means the SiliconSensor
806 is always used, even for the faint things)
808 @param [in] fft_sb_thresh is a surface brightness (photons/pixel) where we will
809 switch from photon shooting to drawing with fft if any pixel is above this.
810 Should be at least the saturation level, if not higher. (default = None, which means
811 never switch to fft.)
813 @param [out] outputString is a string denoting which detectors the astronomical
814 object illumines, suitable for output in the GalSim InstanceCatalog
815 """
816 object_flags = ObjectFlags()
818 # find the detectors which the astronomical object illumines
819 outputString, \
820 detectorList, \
821 centeredObj = self.findAllDetectors(gsObject)
823 # Make sure this object is marked as "drawn" since we only
824 # care that this method has been called for this object.
825 self.drawn_objects.add(gsObject.uniqueId)
827 # Compute the realized object fluxes (as drawn from the
828 # corresponding Poisson distribution) for each band and return
829 # right away if all values are zero in order to save compute.
830 fluxes = [gsObject.flux(bandpassName) for bandpassName in self.bandpassDict]
831 realized_fluxes = [galsim.PoissonDeviate(self._rng, mean=f)() for f in fluxes]
832 if all([f == 0 for f in realized_fluxes]):
833 # All fluxes are 0, so no photons will be shot.
834 object_flags.set_flag('skipped')
835 self._store_zero_flux_centroid_info(detectorList, fluxes, gsObject,
836 object_flags.value)
837 return outputString
839 if len(detectorList) == 0:
840 # There is nothing to draw
841 return outputString
843 self._addNoiseAndBackground(detectorList)
845 # Create a surface operation to sample incident angles and a
846 # galsim.SED object for sampling the wavelengths of the
847 # incident photons.
848 fratio = 1.234 # From https://www.lsst.org/scientists/keynumbers
849 obscuration = 0.606 # (8.4**2 - 6.68**2)**0.5 / 8.4
850 angles = galsim.FRatioAngles(fratio, obscuration, self._rng)
852 faint = all([f < max_flux_simple for f in realized_fluxes])
854 if faint:
855 # For faint things, use a very simple SED, since we don't really care about getting
856 # the exact right distribution of wavelengths here. (Impacts DCR and electron
857 # conversion depth in silicon)
858 gs_sed = self.trivial_sed
859 object_flags.set_flag('simple_sed')
860 else:
861 sed_lut = galsim.LookupTable(x=gsObject.sed.wavelen,
862 f=gsObject.sed.flambda)
863 gs_sed = galsim.SED(sed_lut, wave_type='nm', flux_type='flambda',
864 redshift=0.)
866 ra_obs, dec_obs = observedFromPupilCoords(gsObject.xPupilRadians,
867 gsObject.yPupilRadians,
868 obs_metadata=self.obs_metadata)
869 obj_coord = galsim.CelestialCoord(ra_obs*galsim.degrees,
870 dec_obs*galsim.degrees)
872 for bandpassName, realized_flux, flux in zip(self.bandpassDict, realized_fluxes, fluxes):
873 gs_bandpass = self.gs_bandpass_dict[bandpassName]
874 waves = galsim.WavelengthSampler(sed=gs_sed, bandpass=gs_bandpass,
875 rng=self._rng)
876 dcr = galsim.PhotonDCR(base_wavelength=gs_bandpass.effective_wavelength,
877 HA=self.local_hour_angle,
878 latitude=self.obs_latitude,
879 obj_coord=obj_coord)
881 # Set the object flux to the value realized from the
882 # Poisson distribution.
883 obj = centeredObj.withFlux(realized_flux)
885 use_fft = False
886 if realized_flux > 1.e6 and fft_sb_thresh is not None and realized_flux > fft_sb_thresh:
887 # Note: Don't bother with this check unless the total flux is > thresh.
888 # Otherwise, there is no chance that the flux in 1 pixel is > thresh.
889 # Also, the cross-over point for time to where the fft becomes faster is
890 # emprically around 1.e6 photons, so also don't bother unless the flux
891 # is more than this.
892 obj, use_fft = self.maybeSwitchPSF(gsObject, obj, fft_sb_thresh)
894 if use_fft:
895 object_flags.set_flag('fft_rendered')
896 object_flags.set_flag('no_silicon')
898 for detector in detectorList:
900 name = self._getFileName(detector=detector,
901 bandpassName=bandpassName)
903 xPix, yPix = detector.camera_wrapper\
904 .pixelCoordsFromPupilCoords(gsObject.xPupilRadians,
905 gsObject.yPupilRadians,
906 chipName=detector.name,
907 obs_metadata=self.obs_metadata)
909 # Desired position to draw the object.
910 image_pos = galsim.PositionD(xPix, yPix)
912 # Find a postage stamp region to draw onto. Use (sky
913 # noise)/3. as the nominal minimum surface brightness
914 # for rendering an extended object.
915 keep_sb_level = np.sqrt(self.sky_bg_per_pixel)/3.
916 full_bounds = self.getStampBounds(gsObject, realized_flux, image_pos,
917 keep_sb_level, 3*keep_sb_level)
919 # Ensure the bounds of the postage stamp lie within the image.
920 bounds = full_bounds & self.detectorImages[name].bounds
922 if not bounds.isDefined():
923 continue
925 # Offset is relative to the "true" center of the postage stamp.
926 offset = image_pos - bounds.true_center
928 image = self.detectorImages[name][bounds]
930 if faint:
931 # For faint things, only use the silicon sensor if there is already
932 # some significant flux on the image near the object.
933 # Brighter-fatter doesn't start having any measurable effect until at least
934 # around 1000 e-/pixel. So a limit of 200 is conservative by a factor of 5.
935 # Do the calculation relative to the median, since a perfectly flat sky level
936 # will not have any B/F effect. (But noise fluctuations due to the sky will
937 # be properly included here if the sky is drawn first.)
938 if np.max(image.array) > np.median(image.array) + sensor_limit:
939 sensor = self.sensor[detector.name]
940 else:
941 sensor = None
942 object_flags.set_flag('no_silicon')
943 else:
944 sensor = self.sensor[detector.name]
946 if sensor:
947 # Ensure the rng used by the sensor object is set to the desired state.
948 self.sensor[detector.name].rng.reset(self._rng)
949 surface_ops = [waves, dcr, angles]
950 else:
951 # Don't need angles if not doing silicon sensor.
952 surface_ops = [waves, dcr]
954 if use_fft:
955 # When drawing with FFTs, large offsets can be a problem, since they
956 # can blow up the required FFT size. We'll guard for that below with
957 # a try block, but we can minimize how often this happens by making sure
958 # the offset is close to 0,0.
959 if abs(offset.x) > 2 or abs(offset.y) > 2:
960 # Make a larger image that has the object near the center.
961 fft_image = galsim.Image(full_bounds, dtype=image.dtype, wcs=image.wcs)
962 fft_image[bounds] = image
963 fft_offset = image_pos - full_bounds.true_center
964 else:
965 fft_image = image.copy()
966 fft_offset = offset
968 try:
969 obj.drawImage(method='fft',
970 offset=fft_offset,
971 image=fft_image,
972 gain=detector.photParams.gain)
973 except galsim.errors.GalSimFFTSizeError:
974 use_fft = False
975 object_flags.unset_flag('fft_rendered')
976 if sensor is not None:
977 object_flags.unset_flag('no_silicon')
978 else:
979 # Some pixels can end up negative from FFT numerics. Just set them to 0.
980 fft_image.array[fft_image.array < 0] = 0.
981 fft_image.addNoise(galsim.PoissonNoise(rng=self._rng))
982 # In case we had to make a bigger image, just copy the part we need.
983 image += fft_image[bounds]
984 if not use_fft:
985 obj.drawImage(method='phot',
986 offset=offset,
987 rng=self._rng,
988 maxN=int(1e6),
989 image=image,
990 sensor=sensor,
991 surface_ops=surface_ops,
992 add_to_image=True,
993 poisson_flux=False,
994 gain=detector.photParams.gain)
996 # If we are writing centroid files,store the entry.
997 if self.centroid_base_name is not None:
998 centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId,
999 flux, realized_flux, xPix, yPix, object_flags.value,
1000 gsObject.galSimType)
1001 self.centroid_list.append(centroid_tuple)
1003 # Because rendering FitsImage object types can take a long
1004 # time for bright objects (>1e4 photons takes longer than ~30s
1005 # on cori-haswell), force a checkpoint after each object is
1006 # drawn.
1007 force_checkpoint = ((gsObject.galSimType == 'FitsImage') and
1008 realized_flux > 1e4)
1009 self.write_checkpoint(force=force_checkpoint)
1010 return outputString
1012 @staticmethod
1013 def maybeSwitchPSF(gsObject, obj, fft_sb_thresh, pixel_scale=0.2):
1014 """
1015 Check if the maximum surface brightness of the object is high enough that we should
1016 switch to using an fft method rather than photon shooting.
1018 When we do this, we also switch the PSF model to something slightly simpler with
1019 roughly the same wings, but not as complicated in the central core. Thus, this
1020 should only be done when the core is going to be saturated anyway, so we only really
1021 care about the wings of the PSF.
1023 Note: This function assumes that obj at this point is a convolution with the PSF at the
1024 end, and that it has had its flux set to a new value with `withFlux()`.
1025 If this is not the case, an AttributeError will be raised.
1027 Parameters
1028 ----------
1029 gsObject: GalSimCelestialObject
1030 This contains the information needed to construct a
1031 galsim.GSObject convolved with the desired PSF.
1032 obj: galsim.GSObject
1033 The current GSObject to draw, which might need to be modified
1034 if we decide to switch to fft drawing.
1035 fft_sb_thresh: float
1036 The surface brightness (photons/pixel) where we will switch from
1037 photon shooting to drawing with fft if any pixel is above this.
1038 Should be at least the saturation level, if not higher.
1039 pixel_scale: float [0.2]
1040 The CCD pixel scale in arcsec.
1042 Returns
1043 -------
1044 galsim.GSObj, bool: obj = the object to actually use
1045 use_fft = whether to use fft drawing
1046 """
1047 if not fft_sb_thresh:
1048 return obj, False
1050 # obj.original should be a Convolution with the PSF at the end. Extract it.
1051 geom_psf = obj.original.obj_list[-1]
1052 all_but_psf = obj.original.obj_list[:-1]
1053 try:
1054 screen_list = geom_psf.screen_list
1055 except AttributeError:
1056 # If it's not a galsim.PhaseScreenPSF, just use whatever it is.
1057 fft_psf = [geom_psf]
1058 else:
1059 # If geom_psf is a PhaseScreenPSF, then make a simpler one the just convolves
1060 # a Kolmogorov profile with an OpticalPSF.
1061 opt_screens = [s for s in geom_psf.screen_list if isinstance(s, galsim.OpticalScreen)]
1062 if len(opt_screens) >= 1:
1063 # Should never be more than 1, but it there weirdly is, just use the first.
1064 opt_screen = opt_screens[0]
1065 optical_psf = galsim.OpticalPSF(
1066 lam=geom_psf.lam,
1067 diam=opt_screen.diam,
1068 aberrations=opt_screen.aberrations,
1069 annular_zernike=opt_screen.annular_zernike,
1070 obscuration=opt_screen.obscuration,
1071 gsparams=geom_psf.gsparams)
1072 fft_psf = [optical_psf]
1073 else:
1074 fft_psf = []
1075 r0_500 = screen_list.r0_500_effective
1076 atm_psf = galsim.Kolmogorov(lam=geom_psf.lam, r0_500=r0_500,
1077 gsparams=geom_psf.gsparams)
1078 fft_psf.append(atm_psf)
1080 fft_obj = galsim.Convolve(all_but_psf + fft_psf).withFlux(obj.flux)
1082 # Now this object should have a much better estimate of the real maximum surface brightness
1083 # than the original geom_psf did.
1084 # However, the max_sb feature gives an over-estimate, whereas to be conservative, we would
1085 # rather an under-estimate. For this kind of profile, dividing by 2 does a good job
1086 # of giving us an underestimate of the max surface brightness.
1087 # Also note that `max_sb` is in photons/arcsec^2, so multiply by pixel_scale**2
1088 # to get photons/pixel, which we compare to fft_sb_thresh.
1089 if fft_obj.max_sb/2. * pixel_scale**2 > fft_sb_thresh:
1090 return fft_obj, True
1091 else:
1092 return obj, False
1094 def getStampBounds(self, gsObject, flux, image_pos, keep_sb_level,
1095 large_object_sb_level, Nmax=1400, pixel_scale=0.2):
1096 """
1097 Get the postage stamp bounds for drawing an object within the stamp
1098 to include the specified minimum surface brightness. Use the
1099 folding_threshold criterion for point source objects. For
1100 extended objects, use the getGoodPhotImageSize function, where
1101 if the initial stamp is too large (> Nmax**2 ~ 1GB of RSS
1102 memory for a 72 vertex/pixel sensor model), use the relaxed
1103 surface brightness level for large objects.
1105 Parameters
1106 ----------
1107 gsObject: GalSimCelestialObject
1108 This contains the information needed to construct a
1109 galsim.GSObject convolved with the desired PSF.
1110 flux: float
1111 The flux of the object in e-.
1112 keep_sb_level: float
1113 The minimum surface brightness (photons/pixel) out to which to
1114 extend the postage stamp, e.g., a value of
1115 sqrt(sky_bg_per_pixel)/3 would be 1/3 the Poisson noise
1116 per pixel from the sky background.
1117 large_object_sb_level: float
1118 Surface brightness level to use for large/bright objects that
1119 would otherwise yield stamps with more than Nmax**2 pixels.
1120 Nmax: int [1400]
1121 The largest stamp size to consider at the nominal keep_sb_level.
1122 1400**2*72*8/1024**3 = 1GB.
1123 pixel_scale: float [0.2]
1124 The CCD pixel scale in arcsec.
1126 Returns
1127 -------
1128 galsim.BoundsI: The postage stamp bounds.
1130 """
1131 if flux < 10:
1132 # For really faint things, don't try too hard. Just use 32x32.
1133 image_size = 32
1134 elif gsObject.galSimType.lower() == "pointsource":
1135 # For bright stars, set the folding threshold for the
1136 # stamp size calculation. Use a
1137 # Kolmogorov_and_Gaussian_PSF since it is faster to
1138 # evaluate than an AtmosphericPSF.
1139 folding_threshold = self.sky_bg_per_pixel/flux
1140 if folding_threshold >= self._ft_default:
1141 gsparams = None
1142 else:
1143 gsparams = galsim.GSParams(folding_threshold=folding_threshold)
1144 psf = Kolmogorov_and_Gaussian_PSF(airmass=self._airmass,
1145 rawSeeing=self._rawSeeing,
1146 band=self._band,
1147 gsparams=gsparams)
1148 obj = self.drawPointSource(gsObject, psf=psf)
1149 image_size = obj.getGoodImageSize(pixel_scale)
1150 else:
1151 # For extended objects, recreate the object to draw, but
1152 # convolved with the faster DoubleGaussian PSF.
1153 obj = self.createCenteredObject(gsObject,
1154 psf=self._double_gaussian_psf)
1155 obj = obj.withFlux(flux)
1157 # Start with GalSim's estimate of a good box size.
1158 image_size = obj.getGoodImageSize(pixel_scale)
1160 # For bright things, defined as having an average of at least 10 photons per
1161 # pixel on average, try to be careful about not truncating the surface brightness
1162 # at the edge of the box.
1163 if flux > 10 * image_size**2:
1164 image_size = self._getGoodPhotImageSize(gsObject, flux, keep_sb_level,
1165 pixel_scale=pixel_scale)
1167 # If the above size comes out really huge, scale back to what you get for
1168 # a somewhat brighter surface brightness limit.
1169 if image_size > Nmax:
1170 image_size = self._getGoodPhotImageSize(gsObject, flux, large_object_sb_level,
1171 pixel_scale=pixel_scale)
1172 image_size = max(image_size, Nmax)
1174 # Create the bounds object centered on the desired location.
1175 xmin = int(math.floor(image_pos.x) - image_size/2)
1176 xmax = int(math.ceil(image_pos.x) + image_size/2)
1177 ymin = int(math.floor(image_pos.y) - image_size/2)
1178 ymax = int(math.ceil(image_pos.y) + image_size/2)
1180 return galsim.BoundsI(xmin, xmax, ymin, ymax)
1182 def _getGoodPhotImageSize(self, gsObject, flux, keep_sb_level,
1183 pixel_scale=0.2):
1184 point_source = self.drawPointSource(gsObject, self._double_gaussian_psf)
1185 point_source = point_source.withFlux(flux)
1186 ps_size = getGoodPhotImageSize(point_source, keep_sb_level,
1187 pixel_scale=pixel_scale)
1188 unconvolved_obj = self._createCenteredObject(gsObject, psf=None)
1189 unconvolved_obj = unconvolved_obj.withFlux(flux)
1190 obj_size = getGoodPhotImageSize(unconvolved_obj, keep_sb_level,
1191 pixel_scale=pixel_scale)
1192 return int(np.sqrt(ps_size**2 + obj_size**2))
1195def getGoodPhotImageSize(obj, keep_sb_level, pixel_scale=0.2):
1196 """
1197 Get a postage stamp size (appropriate for photon-shooting) given a
1198 minimum surface brightness in photons/pixel out to which to
1199 extend the stamp region.
1201 Parameters
1202 ----------
1203 obj: galsim.GSObject
1204 The GalSim object for which we will call .drawImage.
1205 keep_sb_level: float
1206 The minimum surface brightness (photons/pixel) out to which to
1207 extend the postage stamp, e.g., a value of
1208 sqrt(sky_bg_per_pixel)/3 would be 1/3 the Poisson noise
1209 per pixel from the sky background.
1210 pixel_scale: float [0.2]
1211 The CCD pixel scale in arcsec.
1213 Returns
1214 -------
1215 int: The length N of the desired NxN postage stamp.
1217 Notes
1218 -----
1219 Use of this function should be avoided with PSF implementations that
1220 are costly to evaluate. A roughly equivalent DoubleGaussian
1221 could be used as a proxy.
1223 This function was originally written by Mike Jarvis.
1224 """
1225 # The factor by which to adjust N in each step.
1226 factor = 1.1
1228 # Start with the normal image size from GalSim
1229 N = obj.getGoodImageSize(pixel_scale)
1230 #print('N = ',N)
1232 if isinstance(obj, galsim.RandomKnots):
1233 # If the galaxy is a RandomWalk, extract the underlying profile for this calculation
1234 # rather than using the knotty version, which will pose problems for the xValue function.
1235 obj = obj._profile
1237 # This can be too small for bright stars, so increase it in steps until the edges are
1238 # all below the requested sb level.
1239 # (Don't go bigger than 4096)
1240 Nmax = 4096
1241 while N < Nmax:
1242 # Check the edges and corners of the current square
1243 h = N / 2 * pixel_scale
1244 xvalues = [ obj.xValue(h,0), obj.xValue(-h,0),
1245 obj.xValue(0,h), obj.xValue(0,-h),
1246 obj.xValue(h,h), obj.xValue(h,-h),
1247 obj.xValue(-h,h), obj.xValue(-h,-h) ]
1248 maxval = np.max(xvalues)
1249 #print(N, maxval)
1250 if maxval < keep_sb_level:
1251 break
1252 N *= factor
1254 N = min(N, Nmax)
1256 # This can be quite huge for Devauc profiles, but we don't actually have much
1257 # surface brightness way out in the wings. So cut it back some.
1258 # (Don't go below 64 though.)
1259 while N >= 64 * factor:
1260 # Check the edges and corners of a square smaller by a factor of N.
1261 h = N / (2 * factor) * pixel_scale
1262 xvalues = [ obj.xValue(h,0), obj.xValue(-h,0),
1263 obj.xValue(0,h), obj.xValue(0,-h),
1264 obj.xValue(h,h), obj.xValue(h,-h),
1265 obj.xValue(-h,h), obj.xValue(-h,-h) ]
1266 maxval = np.max(xvalues)
1267 #print(N, maxval)
1268 if maxval > keep_sb_level:
1269 break
1270 N /= factor
1272 return int(N)
1275class ObjectFlags:
1276 """
1277 Class to keep track of the object rendering bit flags. The bits
1278 will be composed as an int for storing in centroid files.
1279 """
1280 def __init__(self, conditions='skipped simple_sed no_silicon fft_rendered'.split()):
1281 """
1282 Parameters
1283 ----------
1284 conditions: list or tuple
1285 The sequence of strings describing the various conditions
1286 to be tracked. The order will determine how the bits
1287 are assigned, so it should be a well-ordered sequence, i.e.,
1288 specifically a list or a tuple.
1289 """
1290 if type(conditions) not in (list, tuple):
1291 raise TypeError("conditions must be a list or a tuple")
1292 if len(conditions) != len(set(conditions)):
1293 raise ValueError("conditions must contain unique entries")
1294 self.flags = {condition: 1<<shift for shift, condition in
1295 enumerate(conditions)}
1296 self.value = 0
1298 def set_flag(self, condition):
1299 """
1300 Set the bit associated with the specified condition.
1302 Parameters
1303 ----------
1304 condition: str
1305 A condition not in the known set will raise a ValueError.
1306 """
1307 try:
1308 self.value |= self.flags[condition]
1309 except KeyError:
1310 raise ValueError("unknown bit flag: %s" % condition)
1312 def unset_flag(self, condition):
1313 """
1314 Unset the bit associated with the specified condition.
1316 Parameters
1317 ----------
1318 condition: str
1319 A condition not in the known set will raise a ValueError.
1320 """
1321 try:
1322 self.value &= ~self.flags[condition]
1323 except KeyError:
1324 raise ValueError("unknown bit flag: %s" % condition)