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