LSSTApplications  19.0.0-14-gb0260a2+599893a4c6,20.0.0+126303c00d,20.0.0+2f3d0e5c40,20.0.0+36ef800059,20.0.0+5ac7adcc0c,20.0.0+693a64958a,20.0.0+bebc1f60e8,20.0.0+cad136aba6,20.0.0+e2e26847c2,20.0.0+e69b5d60e7,20.0.0-1-g10df615+11b215b765,20.0.0-1-g253301a+36ef800059,20.0.0-1-g2b7511a+bebc1f60e8,20.0.0-1-g4d801e7+aeeb640673,20.0.0-1-g5b95a8c+f111d5f02f,20.0.0-1-g660595b+f45b7d88f4,20.0.0-1-gc96f8cb+e3b38461e6,20.0.0-1-gd1c87d7+85c46248f3,20.0.0-1-gedffbd8+d0b27f8bcb,20.0.0-16-g111fe95+e3b38461e6,20.0.0-16-g18096c8+d1a4df0137,20.0.0-16-g233ea98+a4df35922d,20.0.0-17-ga9337b4+41f27cfd54,20.0.0-2-g4dae9ad+e3b38461e6,20.0.0-2-g7818986+85c46248f3,20.0.0-2-gec03fae+ff10c6d78d,20.0.0-28-g282f9e7e+feda6aebd8,20.0.0-3-g4cc78c6+63636aeed8,20.0.0-3-g6a8623c+d1a4df0137,20.0.0-3-g750bffe+f5427621ce,20.0.0-4-gfea843c+f45b7d88f4,20.0.0-5-g357b56b+f45b7d88f4,20.0.0-5-gfcebe35+e2b15ed341,20.0.0-52-g73d9071+9bf1eb8e0a,20.0.0-7-gcda7bf1+773ba852cb,20.0.0-8-g4540fe2a+952f6d3c43,20.0.0-9-g61a2a9a3d+14f89e4eca,w.2020.40
LSSTDataManagementBasePackage
makeFakeGalaxy.py
Go to the documentation of this file.
1 import warnings
2 
3 import numpy as np
4 from astropy.io import fits
5 
6 import galsim
7 
8 
9 def makeGalaxy(flux, gal, psfImage,
10  galType='sersic', cosmosCat=None,
11  drawMethod='no_pixel', trunc=10.0,
12  transform=None, addShear=False,
13  calib=None,
14  sersic_prec=0.01, addPoisson=False):
15  """
16  Function called by task to make galaxy images
17 
18  INPUTS:
19  flux = calibrated total flux
20  gal = galaxy parameters in record (np.recarray)
21  psfImage = np.ndarray of psfImage
22  galType = type of galaxy we want to make, this has to agree
23  with what's in the record array, options now are:
24  'sersic' (single sersic),
25  'dsersic' (double sersic), and
26  'real' (for making galaxies from real HST images)
27  'cosmos' (Using GalSim.COSMOSCatalog)
28 
29  All the necessary keywords need to be in the fits catalog,
30  including maybe drawMethod and trunc...
31  """
32  if galType == 'sersic':
33  return galSimFakeSersic(flux, gal, psfImage=psfImage,
34  trunc=trunc,
35  drawMethod=drawMethod,
36  returnObj=False,
37  transform=transform,
38  addShear=addShear,
39  addPoisson=addPoisson)
40 
41  if galType == 'dsersic':
42  # TODO: addShear option is not available for double Sersic yet
43  (comp1, comp2) = parseDoubleSersic(flux, gal)
44  return galSimFakeDoubleSersic(comp1, comp2, psfImage=psfImage,
45  trunc=trunc,
46  drawMethod=drawMethod,
47  returnObj=False,
48  transform=transform,
49  addShear=addShear,
50  addPoisson=addPoisson)
51 
52  if galType == 'real':
53  # TODO: For real galaxies, we need to decide which to use: index in the
54  # catalog or Object ID. Now, I just use Index
55  (real_galaxy_catalog, index) = parseRealGalaxy(gal)
56  if index >= 0:
57  index = index
58  random = False
59  else:
60  index = None
61  random = True
62  return galSimRealGalaxy(flux, real_galaxy_catalog,
63  index=index,
64  psfImage=psfImage,
65  random=random,
66  returnObj=False,
67  drawMethod=drawMethod,
68  transform=transform,
69  addPoisson=addPoisson)
70 
71  if galType == 'cosmos':
72  if cosmosCat is None or calib is None:
73  raise Exception("# No COSMOSCatalog() provided!")
74  return galSimFakeCosmos(cosmosCat, calib, gal,
75  psfImage=psfImage,
76  returnObj=False,
77  sersic_prec=sersic_prec,
78  drawMethod=drawMethod,
79  transform=transform,
80  addShear=addShear,
81  addPoisson=addPoisson)
82 
83 
84 def parseRealGalaxy(gal):
85  """Place holder for real galaxi."""
86  try:
87  index = gal['index']
88  except KeyError:
89  index = -1
90 
91  # Get the catalog name and the directory
92  try:
93  cat_name = gal['cat_name']
94  except KeyError:
95  raise KeyError('Can not find the name of the catlog')
96  try:
97  cat_dir = gal['cat_dir']
98  except KeyError:
99  cat_dir = None
100 
101  real_galaxy_catalog = galsim.RealGalaxyCatalog(cat_name,
102  dir=cat_dir)
103 
104  return real_galaxy_catalog, index
105 
106 
107 def parseDoubleSersic(tflux, gal):
108  """
109  Parse the input total flux [tflux] and parameter record array
110  [gal] into two parameter records for each component [comp1, comp2]
111  """
112 
113  # Check if this is a real 2-Sersic record
114  try:
115  frac1 = float(gal['b2t'])
116  except KeyError:
117  raise KeyError("No b2t parameter is found in the record!!")
118 
119  # Make sure the flux fraction of the first component is reasonable
120  if (frac1 <= 0) or (frac1 >= 1):
121  raise Exception("b2t should be > 0 and <1 !!")
122  flux1, flux2 = (tflux * frac1), (tflux * (1.0 - frac1))
123 
124  # Check, then read in other parameters
125  galID = gal["ID"]
126  # Effective radius
127  try:
128  reff1, reff2 = float(gal["reff1"]), float(gal["reff2"])
129  except KeyError:
130  raise KeyError("reff1 or reff2 is found in the record!!")
131  # Sersic index
132  try:
133  nser1, nser2 = float(gal["sersic_n1"]), float(gal["sersic_n2"])
134  except KeyError:
135  raise KeyError("sersic_n1 or sersic_n2 is found in the record!!")
136  # Axis ratio
137  try:
138  ba1, ba2 = float(gal["b_a1"]), float(gal["b_a2"])
139  except KeyError:
140  raise KeyError("b_a1 or b_a2 is found in the record!!")
141  # Position angle
142  try:
143  pa1, pa2 = float(gal["theta1"]), float(gal["theta2"])
144  except KeyError:
145  raise KeyError("theta1 or theta2 is found in the record!!")
146 
147  comp1 = np.array((galID, flux1, nser1, reff1, ba1, pa1),
148  dtype=[('ID', 'int'),
149  ('mag', 'float'),
150  ('sersic_n', 'float'),
151  ('reff', 'float'),
152  ('b_a', 'float'),
153  ('theta', 'float')])
154  comp2 = np.array((galID, flux2, nser2, reff2, ba2, pa2),
155  dtype=[('ID', 'int'),
156  ('mag', 'float'),
157  ('sersic_n', 'float'),
158  ('reff', 'float'),
159  ('b_a', 'float'),
160  ('theta', 'float')])
161 
162  return comp1, comp2
163 
164 
165 def arrayToGSObj(imgArr, scale=1.0, norm=False):
166  """
167  Convert an input 2-D array into a GalSim Image object
168 
169  TODO : Check the scale here
170  According to the GalSim Doxygen
171  If provided, use this as the pixel scale for the Image;
172  this will override the pixel scale stored by the provided Image.
173  If scale is None, then take the provided image's pixel scale.
174  [default: None]
175  """
176  if norm:
177  return galsim.InterpolatedImage(galsim.image.Image(imgArr),
178  scale=scale,
179  normalization="flux")
180  else:
181  return galsim.InterpolatedImage(galsim.image.Image(imgArr),
182  scale=scale)
183 
184 
185 def galSimDrawImage(galObj, size=0, scale=1.0, method="no_pixel",
186  addPoisson=False):
187  """
188  "Draw" a GalSim Object into an GalSim Image using certain method, and with
189  certain size
190 
191  TODO : Think about the scale here:
192  By default scale=None
193  According to GalSim Doxygen :
194  If provided, use this as the pixel scale for the image. If scale is None
195  and image != None, then take the provided image's pixel scale.
196  If scale is None and image == None, then use the Nyquist scale.
197  If scale <= 0 (regardless of image), then use the Nyquist scale.
198  """
199 
200  # Generate an "Image" object for the model
201  if size > 0:
202  imgTemp = galsim.image.Image(size, size)
203  galImg = galObj.drawImage(imgTemp, scale=scale, method=method)
204  else:
205  galImg = galObj.drawImage(scale=scale, method=method)
206 
207  # Just an option for test
208  if addPoisson:
209  galImg.addNoise(galsim.PoissonNoise())
210 
211  # Return the Numpy array version of the image
212  return galImg.array
213 
214 
215 def galSimConvolve(galObj, psfObj, size=0, scale=1.0, method="no_pixel",
216  returnObj=False):
217  """
218  Just do convolution using GalSim
219 
220  Make sure the inputs are both GalSim GSObj
221  The galaxy model should be the first one, and the PSF object is the second
222  one; Returns a imgArr or GSObj
223  """
224  outObj = galsim.Convolve([galObj, psfObj])
225 
226  if returnObj:
227  return outObj
228  else:
229  outArr = galSimDrawImage(galObj, size=size,
230  scale=scale, method=method)
231  return outArr
232 
233 
234 def galSimAdd(galObjList, size=0, scale=1.0, method="no_pixel",
235  returnArr=False):
236  """
237  Just add a list of GSObjs together using GalSim
238  Make sure all elements in the input list are GSObjs
239  """
240  if len(galObjList) < 2:
241  raise Exception("Should be more than one GSObjs to add !")
242 
243  outObj = galsim.Add(galObjList)
244 
245  if returnArr:
246  outArr = galSimDrawImage(outObj, size=size, scale=scale,
247  method=method)
248  return outArr
249  else:
250  return outObj
251 
252 
253 def plotFakeGalaxy(galObj, galID=None, suffix=None,
254  size=0, addPoisson=False):
255  """
256  Generate a PNG image of the model
257  By default, the image will be named as 'fake_galaxy.png'
258  """
259 
260  import matplotlib.pyplot as plt
261 
262  if galID is None:
263  outPNG = 'fake_galaxy'
264  else:
265  outPNG = 'fake_galaxy_%i' % galID
266  if suffix is not None:
267  outPNG = outPNG + '_' + suffix.strip() + '.png'
268 
269  plt.figure(1, figsize=(8, 8))
270 
271  # Use "fft" just to be fast
272  plt.imshow(np.arcsinh(galSimDrawImage(galObj, size=size,
273  method="no_pixel",
274  addPoisson=addPoisson,
275  scale=1.0)))
276  plt.savefig(outPNG)
277 
278 
279 def galSimFakeCosmos(cosmosCat, calib, gal,
280  psfImage=None, plotFake=False,
281  returnObj=True, sersic_prec=0.02,
282  drawMethod='no_pixel', scale=1.0,
283  transform=None, addShear=False,
284  addPoisson=False):
285  """
286  Generate fake galaxy using galSim.COSMOSCatalog objects.
287  """
288  # Select the useful parameter catalog
289  indexUse = cosmosCat.orig_index
290  paramCat = cosmosCat.param_cat
291  objectUse = paramCat[indexUse]
292 
293  galFound = np.where(objectUse['IDENT'] == gal['COSMOS_ID'])[0]
294  if len(galFound) == 0:
295  warnings.warn("# Find no match for %d" % gal['COSMOS_ID'])
296  elif len(galFound) > 1:
297  warnings.warn("# Multiple match for %d" % gal['COSMOS_ID'])
298 
299  galIndex = galFound[0]
300  cosObj = cosmosCat.makeGalaxy(index=galIndex,
301  sersic_prec=sersic_prec,
302  gal_type='parametric')
303  hstFlux = cosObj.flux
304 
305  # If necessary, apply addtion shear (e.g. for weak lensing test)
306  if addShear:
307  try:
308  g1 = float(gal['g1'])
309  g2 = float(gal['g2'])
310  cosObj = cosObj.shear(g1=g1, g2=g2)
311  except ValueError:
312  warnings.warn("Can not find g1 or g2 in the input!\n",
313  " No shear has been added!")
314 
315  # Do the transformation from sky to pixel coordinates, if given
316  if transform is not None:
317  cosObj = cosObj.transform(*tuple(transform.ravel()))
318 
319  # TODO : Should check the flux
320  """
321  The flux of the galaxy corresponds to a 1 second exposure time with the
322  Hubble Space Telescope. Users who wish to simulate F814W images with a
323  different telescope and an exposure time longer than 1 second should
324  multiply by that exposure time, and by the square of the ratio of the
325  effective diameter of their telescope compared to that of HST.
326  (Effective diameter may differ from the actual diameter if there is
327  significant obscuration.)
328 
329  The catalog returns objects that are appropriate for HST in 1 second
330  exposures. So for our telescope we scale up by the relative area and
331  exposure time. Note that what is important is the *effective* area after
332  taking into account obscuration. For HST, the telescope diameter is 2.4
333  but there is obscuration (a linear factor of 0.33). Here, we assume that
334  the telescope we're simulating effectively has no obscuration factor.
335  We're also ignoring the pi/4 factor since it appears in the numerator and
336  denominator, so we use area = diam^2.
337  """
338  """
339  hstEffArea = (2.4 ** 2) * (1.0 - 0.33 ** 2)
340  subaruEffArea = (8.2 ** 2) * (1.0 - 0.1 ** 2)
341  fluxScaling = (subaruEffArea / hstEffArea)
342  """
343  hstMag = -2.5 * np.log10(hstFlux) + 25.94
344  # hscFlux = 10.0 ** ((27.0 - hstMag) / 2.5)
345  hscFlux = calib.getFlux(hstMag)
346  # cosObj *= (hscFlux / hstFlux)
347  cosObj = cosObj.withFlux(float(hscFlux))
348 
349  # Convolve the Sersic model using the provided PSF image
350  if psfImage is not None:
351  # Convert the PSF Image Array into a GalSim Object
352  # Norm=True by default
353  psfObj = arrayToGSObj(psfImage, norm=True)
354  cosFinal = galsim.Convolve([cosObj, psfObj])
355  else:
356  cosFinal = cosObj
357 
358  # Make a PNG figure of the fake galaxy to check if everything is Ok
359  if plotFake:
360  plotFakeGalaxy(cosFinal, galID=gal['ID'])
361 
362  # Now, by default, the function will just return the GSObj
363  if returnObj:
364  return cosFinal
365  else:
366  return galSimDrawImage(cosFinal,
367  method=drawMethod, scale=scale,
368  addPoisson=addPoisson)
369 
370 
371 def galSimFakeSersic(flux, gal, psfImage=None, scaleRad=False, returnObj=True,
372  expAll=False, devAll=False, plotFake=False, trunc=0,
373  drawMethod="no_pixel", addPoisson=False, scale=1.0,
374  transform=None, addShear=False):
375  """
376  Make a fake single Sersic galaxy using the galSim.Sersic function
377 
378  Inputs: total flux of the galaxy, and a record array that stores the
379  necessary parameters [reffPix, nSersic, axisRatio, posAng]
380 
381  Output: a 2-D image array of the galaxy model OR
382  a GalSim object of the model
383 
384  Options:
385  psfImage: PSF image for convolution
386  trunc: Flux of Sersic models will truncate at trunc * reffPix
387  radius; trunc=0 means no truncation
388  drawMethod: The method for drawImage: ['auto', 'fft', 'real_space']
389  addPoisson: Add Poisson noise
390  plotFake: Generate a PNG figure of the model
391  expAll: Input model will be seen as nSersic=1
392  devAll: Input model will be seen as nSersic=4
393  returnObj: If TRUE, will return the GSObj
394  """
395  reff = float(gal["reff"])
396  posAng = float(gal["theta"])
397  axisRatio = float(gal["b_a"])
398  nSersic = float(gal["sersic_n"])
399 
400  # Truncate the flux at trunc x reff
401  if trunc > 0:
402  trunc = trunc * reff
403 
404  # Make sure Sersic index is not too large
405  if nSersic > 6.0:
406  raise ValueError("Sersic index is too large! Should be <= 6.0")
407  # Check the axisRatio value
408  if axisRatio <= 0.05:
409  raise ValueError("Axis Ratio is too small! Should be >= 0.05")
410 
411  # Make the Sersic model based on flux, re, and Sersic index
412  if nSersic == 1.0 or expAll:
413  if scaleRad:
414  serObj = galsim.Exponential(scale_radius=reff)
415  else:
416  serObj = galsim.Exponential(half_light_radius=reff)
417  if expAll:
418  print(" * Treated as a n=1 Exponential disk : %d" % (gal["ID"]))
419  elif nSersic == 4.0 or devAll:
420  serObj = galsim.DeVaucouleurs(half_light_radius=reff, trunc=trunc)
421  if devAll:
422  print(" * Treated as a n=4 De Vaucouleurs model: %d" % (gal["ID"]))
423  elif nSersic <= 0.9:
424  serObj = galsim.Sersic(nSersic, half_light_radius=reff)
425  else:
426  serObj = galsim.Sersic(nSersic, half_light_radius=reff,
427  trunc=trunc)
428 
429  # If necessary, apply the Axis Ratio (q=b/a) using the Shear method
430  if axisRatio < 1.0:
431  serObj = serObj.shear(q=axisRatio, beta=(0.0 * galsim.degrees))
432 
433  # If necessary, apply the Position Angle (theta) using the Rotate method
434  # if posAng != 0.0 or posAng != 180.0:
435  serObj = serObj.rotate((90.0 - posAng) * galsim.degrees)
436 
437  # If necessary, apply addtion shear (e.g. for weak lensing test)
438  if addShear:
439  try:
440  g1 = float(gal['g1'])
441  g2 = float(gal['g2'])
442  serObj = serObj.shear(g1=g1, g2=g2)
443  except ValueError:
444  warnings.warn("Can not find g1 or g2 in the input!\n",
445  " No shear has been added!")
446 
447  # Do the transformation from sky to pixel coordinates, if given
448  if transform is not None:
449  serObj = serObj.transform(*tuple(transform.ravel()))
450 
451  # Pass the flux to the object
452  serObj = serObj.withFlux(float(flux))
453 
454  # Convolve the Sersic model using the provided PSF image
455  if psfImage is not None:
456  # Convert the PSF Image Array into a GalSim Object
457  # Norm=True by default
458  psfObj = arrayToGSObj(psfImage, norm=True)
459  serFinal = galsim.Convolve([serObj, psfObj])
460  else:
461  serFinal = serObj
462 
463  # Make a PNG figure of the fake galaxy to check if everything is Ok
464  if plotFake:
465  plotFakeGalaxy(serFinal, galID=gal['ID'])
466 
467  # Now, by default, the function will just return the GSObj
468  if returnObj:
469  return serFinal
470  else:
471  return galSimDrawImage(serFinal, method=drawMethod, scale=scale,
472  addPoisson=addPoisson)
473 
474 
475 def galSimFakeDoubleSersic(comp1, comp2, psfImage=None, trunc=0,
476  returnObj=True, devExp=False,
477  plotFake=False, drawMethod='auto',
478  addPoisson=False, scale=1.0, transform=None,
479  addShear=False):
480  """
481  Make a fake double Sersic galaxy using the galSim.Sersic function
482 
483  Inputs: total flux of the galaxy, and a record array that stores the
484  necessary parameters [reffPix, nSersic, axisRatio, posAng]
485 
486  Output: a 2-D image array of the galaxy model OR
487  a GalSim object of the model
488 
489  Options:
490  psfImage: PSF image for convolution
491  trunc: Flux of Sersic models will truncate at trunc * reffPix
492  radius; trunc=0 means no truncation
493  drawMethod: The method for drawImage: ['auto', 'fft', 'real_space']
494  addPoisson: Add Poisson noise
495  plotFake: Generate a PNG figure of the model
496  devexp: The first component will be seen as a nSersic=4 bulge;
497  And, the second one will be seen as a nSersic=1 disk
498  returnObj: If TRUE, will return the GSObj
499  """
500 
501  # Get the flux of both components
502  flux1 = float(comp1['mag'])
503  flux2 = float(comp2['mag'])
504  # tflux = flux1 + flux2
505 
506  # If devExp = True : Treat the first component as an n=4 DeVaucouleurs
507  # and, the second component as an n=1 Exponential disk
508  if devExp:
509  serModel1 = galSimFakeSersic(flux1, comp1, returnObj=True, devAll=True,
510  trunc=trunc, addShear=addShear)
511  serModel2 = galSimFakeSersic(flux2, comp2, returnObj=True, expAll=True,
512  trunc=trunc, addShear=addShear)
513  else:
514  serModel1 = galSimFakeSersic(flux1, comp1, returnObj=True, trunc=trunc,
515  addShear=addShear)
516  serModel2 = galSimFakeSersic(flux2, comp2, returnObj=True, trunc=trunc,
517  addShear=addShear)
518 
519  # Combine these two components
520  doubleSersic = galSimAdd([serModel1, serModel2])
521 
522  # Do the transformation from sky to pixel coordinates, if given
523  if transform is not None:
524  doubleSersic = doubleSersic.transform(*tuple(transform.ravel()))
525 
526  # Convolve the Sersic model using the provided PSF image
527  if psfImage is not None:
528  # Convert the PSF Image Array into a GalSim Object
529  # Norm=True by default
530  psfObj = arrayToGSObj(psfImage, norm=True)
531  dserFinal = galsim.Convolve([doubleSersic, psfObj])
532  else:
533  dserFinal = doubleSersic
534 
535  # Make a PNG figure of the fake galaxy to check if everything is Ok
536  if plotFake:
537  if devExp:
538  plotFakeGalaxy(dserFinal, galID=comp1['ID'], suffix='devexp')
539  else:
540  plotFakeGalaxy(dserFinal, galID=comp1['ID'], suffix='double')
541 
542  # Now, by default, the function will just return the GSObj
543  if returnObj:
544  return dserFinal
545  else:
546  return galSimDrawImage(dserFinal, method=drawMethod, scale=scale,
547  addPoisson=addPoisson)
548 
549 
550 def galSimRealGalaxy(flux, real_galaxy_catalog, index=None, psfImage=None,
551  random=False, returnObj=True, plotFake=False,
552  drawMethod='auto', addPoisson=False, scale=1.0,
553  transform=None):
554  """
555  Real galaxy.
556  """
557 
558  if index is None:
559  random = True
560  realObj = galsim.RealGalaxy(real_galaxy_catalog, index=index,
561  random=random)
562  index = realObj.index
563 
564  # Pass the flux to the object
565  realObj = realObj.withFlux(flux)
566 
567  # Do the transformation from sky to pixel coordinates, if given
568  if transform is not None:
569  realObj = realObj.transform(*tuple(transform.ravel()))
570 
571  # Convolve the Sersic model using the provided PSF image
572  if psfImage is not None:
573  # Convert the PSF Image Array into a GalSim Object
574  # Norm=True by default
575  psfObj = arrayToGSObj(psfImage, norm=True)
576  realFinal = galsim.Convolve([realObj, psfObj])
577  else:
578  realFinal = realFinal
579 
580  # Make a PNG figure of the fake galaxy to check if everything is Ok
581  if plotFake:
582  plotFakeGalaxy(realFinal, galID=index, suffix='realga')
583 
584  # Now, by default, the function will just return the GSObj
585  if returnObj:
586  return realFinal
587  else:
588  return galSimDrawImage(realFinal, method=drawMethod, scale=scale,
589  addPoisson=addPoisson)
590 
591 
592 def testMakeFake(galList, asciiTab=False, single=True, double=True, real=True):
593  """Test the makeFake functions."""
594  # Make a fake Gaussian PSF
595  psfGaussian = galsim.Gaussian(fwhm=2.0)
596  psfImage = psfGaussian.drawImage().array
597 
598  # Test SingleSersic
599  if single:
600  if asciiTab:
601  galData = np.loadtxt(galList, dtype=[('ID', 'int'),
602  ('mag', 'float'),
603  ('sersic_n', 'float'),
604  ('reff', 'float'),
605  ('b_a', 'float'),
606  ('theta', 'float')])
607  else:
608  galData = fits.open(galList)[1].data
609 
610  for igal, gal in enumerate(galData):
611 
612  flux = 10.0 ** ((27.0 - gal['mag']) / 2.5)
613 
614  print('\n---------------------------------')
615  print(" Input Flux : ", flux)
616  print(" Input Parameters : ", gal["sersic_n"], gal["reff"])
617  print(" ", gal["b_a"], gal["theta"])
618 
619  galArray = galSimFakeSersic(flux, gal, psfImage=psfImage,
620  plotFake=True, returnObj=False,
621  trunc=12.0, drawMethod="no_pixel")
622 
623  print(" Output Flux : ", np.sum(galArray))
624  print(" Shape of the Output Array : ", galArray.shape)
625  print('---------------------------------')
626 
627  # Test DoubleSersic
628  if double:
629  if asciiTab:
630  raise Exception("For now, only FITS input is allowed !!")
631  else:
632  galData = fits.open(galList)[1].data
633 
634  for igal, gal in enumerate(galData):
635 
636  flux = 10.0 ** ((27.0 - gal['mag']) / 2.5)
637 
638  print('\n---------------------------------')
639  print(" Input Flux : ", flux)
640 
641  (comp1, comp2) = parseDoubleSersic(flux, gal)
642 
643  # TODO: Get error when the axis ratio is small: 0.2?
644  # RuntimeError: Solve error: Too many iterations in
645  # bracketLowerWithLimit()
646  # It seems like that GalSim has some issues with highly elliptical
647  # objects. Although different, see this one:
648  # https://github.com/GalSim-developers/GalSim/issues/384
649  # It seems that b/a = 0.25 is fine, so right now, just change the
650  # lower limit of b/a to 0.25
651 
652  print(" Flux for Component 1 : ", comp1['mag'])
653  print(" Flux for Component 2 : ", comp2['mag'])
654  print(" Comp 1 Parameters : %5.2f %8.2f" % (comp1["sersic_n"],
655  comp1["reff"]))
656  print(" %5.2f %8.2f" % (comp1["b_a"],
657  comp1["theta"]))
658  print(" Comp 2 Parameters : %5.2f %8.2f" % (comp2["sersic_n"],
659  comp2["reff"]))
660  print(" %5.2f %8.2f" % (comp2["b_a"],
661  comp2["theta"]))
662 
663  doubleArray = galSimFakeDoubleSersic(comp1, comp2,
664  psfImage=psfImage,
665  trunc=12, returnObj=False,
666  devExp=True, plotFake=True,
667  drawMethod='no_pixel')
668 
669  print(" Output Flux : ", np.sum(doubleArray))
670  print(" Shape of the Output Array : ", doubleArray.shape)
671  print('---------------------------------')
672 
673  # Test RealGalaxy
674  if real:
675  # Make a special PSF for real galaxy
676  psfReal = galsim.Gaussian(fwhm=0.2)
677  psfRealImage = psfReal.drawImage().array
678  # TODO : Scale seems to be a problem, should check again !!
679 
680  if asciiTab:
681  raise Exception("For now, only FITS input is allowed !!")
682  else:
683  galData = fits.open(galList)[1].data
684 
685  for igal, gal in enumerate(galData):
686 
687  (real_galaxy_catalog, index) = parseRealGalaxy(gal)
688  if index >= 0:
689  index = index
690  random = False
691  else:
692  index = None
693  random = True
694 
695  flux = 10.0 ** ((27.0 - gal['mag']) / 2.5)
696 
697  realArray = galSimRealGalaxy(flux, real_galaxy_catalog,
698  index=index,
699  psfImage=psfRealImage,
700  random=random,
701  plotFake=True,
702  returnObj=False,
703  drawMethod='no_pixel')
704 
705  print('\n---------------------------------')
706  print(" Input Flux : ", flux)
707 
708  print(" Output Flux : ", np.sum(realArray))
709  print(" Shape of the Output Array : ", realArray.shape)
710  print('---------------------------------')
lsst.synpipe.makeFakeGalaxy.testMakeFake
def testMakeFake(galList, asciiTab=False, single=True, double=True, real=True)
Definition: makeFakeGalaxy.py:592
lsst.synpipe.makeFakeGalaxy.galSimDrawImage
def galSimDrawImage(galObj, size=0, scale=1.0, method="no_pixel", addPoisson=False)
Definition: makeFakeGalaxy.py:185
lsst.synpipe.makeFakeGalaxy.arrayToGSObj
def arrayToGSObj(imgArr, scale=1.0, norm=False)
Definition: makeFakeGalaxy.py:165
lsst.synpipe.makeFakeGalaxy.galSimAdd
def galSimAdd(galObjList, size=0, scale=1.0, method="no_pixel", returnArr=False)
Definition: makeFakeGalaxy.py:234
lsst.synpipe.makeFakeGalaxy.galSimFakeCosmos
def galSimFakeCosmos(cosmosCat, calib, gal, psfImage=None, plotFake=False, returnObj=True, sersic_prec=0.02, drawMethod='no_pixel', scale=1.0, transform=None, addShear=False, addPoisson=False)
Definition: makeFakeGalaxy.py:279
lsst.synpipe.makeFakeGalaxy.plotFakeGalaxy
def plotFakeGalaxy(galObj, galID=None, suffix=None, size=0, addPoisson=False)
Definition: makeFakeGalaxy.py:253
lsst.synpipe.makeFakeGalaxy.parseDoubleSersic
def parseDoubleSersic(tflux, gal)
Definition: makeFakeGalaxy.py:107
lsst.synpipe.makeFakeGalaxy.galSimFakeSersic
def galSimFakeSersic(flux, gal, psfImage=None, scaleRad=False, returnObj=True, expAll=False, devAll=False, plotFake=False, trunc=0, drawMethod="no_pixel", addPoisson=False, scale=1.0, transform=None, addShear=False)
Definition: makeFakeGalaxy.py:371
lsst.synpipe.makeFakeGalaxy.parseRealGalaxy
def parseRealGalaxy(gal)
Definition: makeFakeGalaxy.py:84
lsst.synpipe.makeFakeGalaxy.makeGalaxy
def makeGalaxy(flux, gal, psfImage, galType='sersic', cosmosCat=None, drawMethod='no_pixel', trunc=10.0, transform=None, addShear=False, calib=None, sersic_prec=0.01, addPoisson=False)
Definition: makeFakeGalaxy.py:9
lsst.synpipe.makeFakeGalaxy.galSimConvolve
def galSimConvolve(galObj, psfObj, size=0, scale=1.0, method="no_pixel", returnObj=False)
Definition: makeFakeGalaxy.py:215
lsst.synpipe.makeFakeGalaxy.galSimRealGalaxy
def galSimRealGalaxy(flux, real_galaxy_catalog, index=None, psfImage=None, random=False, returnObj=True, plotFake=False, drawMethod='auto', addPoisson=False, scale=1.0, transform=None)
Definition: makeFakeGalaxy.py:550
lsst.synpipe.makeFakeGalaxy.galSimFakeDoubleSersic
def galSimFakeDoubleSersic(comp1, comp2, psfImage=None, trunc=0, returnObj=True, devExp=False, plotFake=False, drawMethod='auto', addPoisson=False, scale=1.0, transform=None, addShear=False)
Definition: makeFakeGalaxy.py:475