def main(argv):
    """
    Make images to be used for characterizing the brighter-fatter effect
      - Each fits file is 5 x 5 postage stamps.
      - Each postage stamp is 40 x 40 pixels.
      - There are 3 sets of 5 images each.  The 5 images are at 5 different flux levels
      - The three sets are (bf_1) B-F off, (bf_2) B-F on, diffusion off, (bf_3) B-F and diffusion on
      - Each image is in output/bf_set/bf_nfile.fits, where set ranges from 1-3 and nfile ranges from 1-5.
    """
    logging.basicConfig(format="%(message)s", level=logging.INFO, stream=sys.stdout)
    logger = logging.getLogger("bf_plots")

    # Add the wavelength info
    bppath = "../../share/bandpasses/"
    sedpath = "../../share/"
    sed = galsim.SED(os.path.join(sedpath, 'CWW_E_ext.sed'), 'nm', 'flambda').thin()

    # Add the directions (seems to work - CL)
    fratio = 1.2
    obscuration = 0.2
    seed = 12345
    assigner = galsim.FRatioAngles(fratio, obscuration, seed)
    bandpass = galsim.Bandpass(os.path.join(bppath, 'LSST_r.dat'), 'nm').thin()
    rng3 = galsim.BaseDeviate(1234)
    sampler = galsim.WavelengthSampler(sed, bandpass, rng3)


    # Define some parameters we'll use below.
    # Normally these would be read in from some parameter file.

    nx_tiles = 10                   #
    ny_tiles = 10                   #
    stamp_xsize = 40                #
    stamp_ysize = 40                #

    random_seed = 6424512           #

    pixel_scale = 0.2               # arcsec / pixel
    sky_level = 0.01                # ADU / arcsec^2

    # Make output directory if not already present.
    if not os.path.isdir('output'):
        os.mkdir('output')

    gal_sigma = 0.2     # arcsec
    psf_sigma = 0.01     # arcsec
    pixel_scale = 0.2  # arcsec / pixel
    noise = 0.01        # standard deviation of the counts in each pixel

    shift_radius = 0.2              # arcsec (=pixels)

    logger.info('Starting bf_plots using:')
    logger.info('    - image with %d x %d postage stamps',nx_tiles,ny_tiles)
    logger.info('    - postage stamps of size %d x %d pixels',stamp_xsize,stamp_ysize)
    logger.info('    - Centroid shifts up to = %.2f pixels',shift_radius)

    rng = galsim.BaseDeviate(5678)    
    sensor1 = galsim.Sensor()
    sensor2 = galsim.SiliconSensor(rng=rng, diffusion_factor=0.0)
    sensor3 = galsim.SiliconSensor(rng=rng)

    for set in range(1,4):
        starttime = time.time()
        exec("sensor = sensor%d"%set)
        for nfile in range(1,6):
            # Make bf_x directory if not already present.
            if not os.path.isdir('output/bf_%d'%set):
                os.mkdir('output/bf_%d'%set)

            gal_file_name = os.path.join('output','bf_%d/bf_%d.fits'%(set,nfile))
            sex_file_name = os.path.join('output','bf_%d/bf_%d_SEX.fits.cat.reg'%(set,nfile))
            sexfile = open(sex_file_name, 'w')
            gal_flux = 2.0e5 * nfile    # total counts on the image
            # Define the galaxy profile
            gal = galsim.Gaussian(flux=gal_flux, sigma=gal_sigma)
            logger.debug('Made galaxy profile')

            # Define the PSF profile
            psf = galsim.Gaussian(flux=1., sigma=psf_sigma) # PSF flux should always = 1
            logger.debug('Made PSF profile')

            # This profile is placed with different orientations and noise realizations
            # at each postage stamp in the gal image.
            gal_image = galsim.ImageF(stamp_xsize * nx_tiles-1 , stamp_ysize * ny_tiles-1,
                                      scale=pixel_scale)
            psf_image = galsim.ImageF(stamp_xsize * nx_tiles-1 , stamp_ysize * ny_tiles-1,
                                      scale=pixel_scale)

            shift_radius_sq = shift_radius**2

            first_in_pair = True  # Make pairs that are rotated by 90 degrees

            k = 0
            for iy in range(ny_tiles):
                for ix in range(nx_tiles):
                    # The normal procedure for setting random numbers in GalSim is to start a new
                    # random number generator for each object using sequential seed values.
                    # This sounds weird at first (especially if you were indoctrinated by Numerical 
                    # Recipes), but for the boost random number generator we use, the "random" 
                    # number sequences produced from sequential initial seeds are highly uncorrelated.
                    # 
                    # The reason for this procedure is that when we use multiple processes to build
                    # our images, we want to make sure that the results are deterministic regardless
                    # of the way the objects get parcelled out to the different processes. 
                    #
                    # Of course, this script isn't using multiple processes, so it isn't required here.
                    # However, we do it nonetheless in order to get the same results as the config
                    # version of this demo script (demo5.yaml).
                    ud = galsim.UniformDeviate(random_seed+k)

                    # Any kind of random number generator can take another RNG as its first 
                    # argument rather than a seed value.  This makes both objects use the same
                    # underlying generator for their pseudo-random values.
                    #gd = galsim.GaussianDeviate(ud, sigma=gal_ellip_rms)

                    # The -1's in the next line are to provide a border of
                    # 1 pixel between postage stamps
                    b = galsim.BoundsI(ix*stamp_xsize+1 , (ix+1)*stamp_xsize-1, 
                                       iy*stamp_ysize+1 , (iy+1)*stamp_ysize-1)
                    sub_gal_image = gal_image[b]
                    sub_psf_image = psf_image[b]

                    # Great08 randomized the locations of the two galaxies in each pair,
                    # but for simplicity, we just do them in sequential postage stamps.

                    if first_in_pair:
                        # Use a random orientation:
                        beta = ud() * 2. * math.pi * galsim.radians

                        # Determine the ellipticity to use for this galaxy.
                        ellip = 0.0
                        first_in_pair = False
                    else:
                        # Use the previous ellip and beta + 90 degrees
                        beta += 90 * galsim.degrees
                        first_in_pair = True

                    # Make a new copy of the galaxy with an applied e1/e2-type distortion 
                    # by specifying the ellipticity and a real-space position angle
                    this_gal = gal#gal.shear(e=ellip, beta=beta)

                    # Apply a random shift_radius:
                    rsq = 2 * shift_radius_sq
                    while (rsq > shift_radius_sq):
                        dx = (2*ud()-1) * shift_radius
                        dy = (2*ud()-1) * shift_radius
                        rsq = dx**2 + dy**2

                    this_gal = this_gal.shift(dx,dy)
                    # Note that the shifted psf that we create here is purely for the purpose of being able
                    # to draw a separate, shifted psf image.  We do not use it when convolving the galaxy
                    # with the psf.
                    this_psf = psf.shift(dx,dy)

                    # Make the final image, convolving with the (unshifted) psf
                    final_gal = galsim.Convolve([psf,this_gal])

                    # Draw the image

                    if ix == 0 and iy == 0:
                        final_gal.drawImage(sub_gal_image, method = 'phot', sensor=sensor, surface_ops=[sampler, assigner], rng = rng, save_photons = True)
                        photon_file = os.path.join('output','bf_%d/bf_%d_nx_%d_ny_%d_photon_file.fits'%(set,nfile,ix,iy))
                        sub_gal_image.photons.write(photon_file)
                    else:
                        final_gal.drawImage(sub_gal_image, method = 'phot', sensor=sensor, surface_ops=[sampler, assigner], rng = rng)
                    
                    # Now add an appropriate amount of noise to get our desired S/N
                    # There are lots of definitions of S/N, but here is the one used by Great08
                    # We use a weighted integral of the flux:
                    #   S = sum W(x,y) I(x,y) / sum W(x,y)
                    #   N^2 = Var(S) = sum W(x,y)^2 Var(I(x,y)) / (sum W(x,y))^2
                    # Now we assume that Var(I(x,y)) is constant so
                    #   Var(I(x,y)) = noise_var
                    # We also assume that we are using a matched filter for W, so W(x,y) = I(x,y).
                    # Then a few things cancel and we find that
                    # S/N = sqrt( sum I(x,y)^2 / noise_var )
                    #
                    # The above procedure is encapsulated in the function image.addNoiseSNR which
                    # sets the flux appropriately given the variance of the noise model.
                    # In our case, noise_var = sky_level_pixel
                    sky_level_pixel = sky_level * pixel_scale**2
                    noise = galsim.PoissonNoise(ud, sky_level=sky_level_pixel)
                    #sub_gal_image.addNoiseSNR(noise, gal_signal_to_noise)

                    # Draw the PSF image
                    # No noise on PSF images.  Just draw it as is.
                    this_psf.drawImage(sub_psf_image)

                    # For first instance, measure moments
                    """
                    if ix==0 and iy==0:
                        psf_shape = sub_psf_image.FindAdaptiveMom()
                        temp_e = psf_shape.observed_shape.e
                        if temp_e > 0.0:
                            g_to_e = psf_shape.observed_shape.g / temp_e
                        else:
                            g_to_e = 0.0
                        logger.info('Measured best-fit elliptical Gaussian for first PSF image: ')
                        logger.info('  g1, g2, sigma = %7.4f, %7.4f, %7.4f (pixels)',
                                    g_to_e*psf_shape.observed_shape.e1,
                                    g_to_e*psf_shape.observed_shape.e2, psf_shape.moments_sigma)
                    """
                    x = b.center().x
                    y = b.center().y
                    logger.info('Galaxy (%d,%d): center = (%.0f,%0.f)  (e,beta) = (%.4f,%.3f)',
                                ix,iy,x,y,ellip,beta/galsim.radians)
                    k = k+1
                    sexline = 'circle %f %f %f\n'%(x+dx/pixel_scale,y+dy/pixel_scale,gal_sigma/pixel_scale)
                    sexfile.write(sexline)

            sexfile.close()
            logger.info('Done making images of postage stamps')

            # Now write the images to disk.
            #psf_image.write(psf_file_name)
            #logger.info('Wrote PSF file %s',psf_file_name)

            gal_image.write(gal_file_name)
            logger.info('Wrote image to %r',gal_file_name)  # using %r adds quotes around filename for us

        finishtime = time.time()
        print("Time to complete set %d = %.2f seconds\n"%(set, finishtime-starttime))
Exemple #2
0
def test_simple():
    """Test the default Sensor class that acts basically like not passing any sensor object.
    """

    # Start with photon shooting, since that's the most typical way that sensors are used.
    obj = galsim.Gaussian(flux=10000, sigma=1.3)

    # We'll draw the same object using SiliconSensor, Sensor, and the default (sensor=None)
    im1 = galsim.ImageD(64, 64, scale=0.3)  # Refefence image with sensor=None
    im2 = galsim.ImageD(64, 64, scale=0.3)  # Use sensor=simple

    rng1 = galsim.BaseDeviate(5678)
    rng2 = galsim.BaseDeviate(5678)

    simple = galsim.Sensor()

    # Start with photon shooting, since that's more straightforward.
    obj.drawImage(im1, method='phot', poisson_flux=False, rng=rng1)
    obj.drawImage(im2, method='phot', poisson_flux=False, sensor=simple, rng=rng2)

    # Should be exactly equal
    np.testing.assert_array_equal(im2.array, im1.array)

    # Fluxes should all equal obj.flux
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=6)

    # Now test fft drawing, which is more complicated with possible temporaries and subsampling.

    im1 = galsim.ImageD(64, 64, scale=0.3)  # Reference image with sensor=None
    im2 = galsim.ImageD(64, 64, scale=0.3)  # Use sensor=simple
    im3 = galsim.ImageD(64, 64, scale=0.3)  # Use sensor=simple, no subsampling
    im4 = galsim.ImageCD(64, 64, scale=0.3) # Equivalent to image2, but checks using a temporary.
                                            # Also check add_to_image=True with im5.
    im5 = galsim.ImageD(64, 64, scale=0.3)  # Check manually convolving by the pixel.
    im6 = galsim.ImageD(64, 64, scale=0.3)  # Check manually convolving by the pixel, n_subsample=1

    # The rng shouldn't matter anymore for these, so just use the default rng=None
    obj.drawImage(im1, method='fft')
    obj.drawImage(im2, method='fft', sensor=simple)
    obj.drawImage(im3, method='fft', sensor=simple, n_subsample=1)
    obj.drawImage(im4, method='fft', sensor=simple, add_to_image=True)

    obj_with_pixel = galsim.Convolve(obj, galsim.Pixel(0.3))
    obj_with_pixel.drawImage(im5, method='no_pixel', sensor=simple)
    obj_with_pixel.drawImage(im6, method='no_pixel', sensor=simple, n_subsample=1)

    # Fluxes should all equal obj.flux
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im4.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im5.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im6.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im4.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im5.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im6.added_flux, obj.flux, decimal=3)

    # im1 and im2 are not precisely equal, since im2 was made with subsampling and then binning,
    # but with a largish object relative to the pixel, it's very close. (cf. similar test below
    # in test_silicon_fft, where the agreement is not so good.)
    print('max diff between im1, im2 with fft = ',np.max(np.abs(im2.array-im1.array)))
    np.testing.assert_almost_equal(im2.array, im1.array, decimal=10)

    # With no subsampling it should be nearly perfect (although this would be expected to be worse
    # when done with a real Sensor model).
    print('max diff without subsampling = ',np.max(np.abs(im3.array-im1.array)))
    np.testing.assert_almost_equal(im3.array, im1.array, decimal=12)

    # Using a temporary (and add_to_image) shouldn't affect anything for the D -> CD case.
    print('max diff with temporary = ',np.max(np.abs(im4.array-im2.array)))
    np.testing.assert_almost_equal(im4.array.real, im2.array, decimal=12)

    # Manual convolution should be identical to what 'fft' does automatically.
    print('max diff with manual pixel conv = ',np.max(np.abs(im5.array-im2.array)))
    #np.testing.assert_almost_equal(im5.array, im2.array, decimal=12)
    print('max diff with manual pixel conv, no subsampling = ',np.max(np.abs(im6.array-im3.array)))
    np.testing.assert_almost_equal(im6.array, im3.array, decimal=12)

    do_pickle(simple)
Exemple #3
0
def test_bf_slopes():
    """Test the brighter-fatter slopes
    with both the B-F effect and diffusion turned on and off.
    """
    from scipy import stats
    simple = galsim.Sensor()

    init_flux = 200000
    obj = galsim.Gaussian(flux=init_flux, sigma=0.3)

    num_fluxes = 5
    x_moments = np.zeros([num_fluxes, 3])
    y_moments = np.zeros([num_fluxes, 3])
    fluxes = np.zeros([num_fluxes])

    for fluxmult in range(num_fluxes):
        rng1 = galsim.BaseDeviate(5678)
        rng2 = galsim.BaseDeviate(5678)
        rng3 = galsim.BaseDeviate(5678)
        # silicon1 has diffusion turned off, silicon2 has it turned on.
        silicon1 = galsim.SiliconSensor(rng=rng1, diffusion_factor=0.0)
        silicon2 = galsim.SiliconSensor(rng=rng2)

        # We'll draw the same object using SiliconSensor, Sensor, and the default (sensor=None)
        im1 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=silicon1 (diffusion off)
        im2 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=silicon2 (diffusion on)
        im3 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=simple

        obj1 = obj * (fluxmult + 1)

        obj1.drawImage(im1, method='phot', poisson_flux=False, sensor=silicon1, rng=rng1)
        obj1.drawImage(im2, method='phot', poisson_flux=False, sensor=silicon2, rng=rng2)
        obj1.drawImage(im3, method='phot', poisson_flux=False, sensor=simple, rng=rng3)

        print('Moments Mx, My, Mxx, Myy, Mxy for im1, flux = %f:'%obj1.flux)
        mom = galsim.utilities.unweighted_moments(im1)
        x_moments[fluxmult,0] = mom['Mxx']
        y_moments[fluxmult,0] = mom['Myy']
        print('Moments Mx, My, Mxx, Myy, Mxy for im2, flux = %f:'%obj1.flux)
        mom = galsim.utilities.unweighted_moments(im2)
        x_moments[fluxmult,1] = mom['Mxx']
        y_moments[fluxmult,1] = mom['Myy']
        print('Moments Mx, My, Mxx, Myy, Mxy for im3, flux = %f:'%obj1.flux)
        mom = galsim.utilities.unweighted_moments(im3)
        x_moments[fluxmult,2] = mom['Mxx']
        y_moments[fluxmult,2] = mom['Myy']
        fluxes[fluxmult] = im1.array.max()
    print('fluxes = ',fluxes)
    print('x_moments = ',x_moments[:,0])
    print('y_moments = ',y_moments[:,0])
    x_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,x_moments[:,0])
    y_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,y_moments[:,0])
    x_slope *= 50000.0 * 100.0
    y_slope *= 50000.0 * 100.0
    print('With BF turned on, diffusion off, x_slope = %.3f, y_slope = %.3f %% per 50K e-'%(
            x_slope, y_slope ))
    assert x_slope > 0.5
    assert y_slope > 0.5
    x_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,x_moments[:,1])
    y_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,y_moments[:,1])
    x_slope *= 50000.0 * 100.0
    y_slope *= 50000.0 * 100.0
    print('With BF turned on, diffusion on, x_slope = %.3f, y_slope = %.3f %% per 50K e-'%(
            x_slope, y_slope ))
    assert x_slope > 0.5
    assert y_slope > 0.5
    x_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,x_moments[:,2])
    y_slope, intercept, r_value, p_value, std_err = stats.linregress(fluxes,y_moments[:,2])
    x_slope *= 50000.0 * 100.0
    y_slope *= 50000.0 * 100.0
    print('With BF turned off, x_slope = %.3f, y_slope = %.3f %% per 50K e-'%(x_slope, y_slope ))
    assert abs(x_slope) < 0.5
    assert abs(y_slope) < 0.5
Exemple #4
0
def test_silicon_fft():
    """Test that drawing with method='fft' also works for SiliconSensor.
    """
    # Lower this somewhat so we get more accurate fluxes from the FFT.
    # (Still only accurate to 3 d.p. though.)
    gsparams = galsim.GSParams(maxk_threshold=1.e-5)
    obj = galsim.Gaussian(flux=3539, sigma=0.3, gsparams=gsparams)

    # We'll draw the same object using SiliconSensor, Sensor, and the default (sensor=None)
    im1 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=silicon
    im2 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=simple
    im3 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=None

    rng = galsim.BaseDeviate(5678)
    silicon = galsim.SiliconSensor(rng=rng, diffusion_factor=0.0)
    simple = galsim.Sensor()

    obj.drawImage(im1, method='fft', sensor=silicon, rng=rng)
    obj.drawImage(im2, method='fft', sensor=simple, rng=rng)
    obj.drawImage(im3, method='fft')

    printval(im1,im2)

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak         radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    # First, im2 and im3 should be almost exactly equal.  Not precisely, since im2 was made with
    # subsampling and then binning, so the FFT ringing is different (im3 is worse in this regard,
    # since it used convolution with a larger pixel).  So 3 digits is all we manage to get here.
    np.testing.assert_almost_equal(im2.array, im3.array, decimal=3)

    # im1 should be similar, but not equal
    np.testing.assert_almost_equal(im1.array/obj.flux, im2.array/obj.flux, decimal=2)

    # Fluxes should all equal obj.flux
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=3)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=3)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=3)

    # Sizes are all about equal since flux is not large enough for B/F to be significant
    sigma_r = 1./np.sqrt(obj.flux) * im1.scale
    np.testing.assert_allclose(r1, r2, atol=2.*sigma_r)
    np.testing.assert_allclose(r2, r3, atol=2.*sigma_r)

    # Repeat with 20X more photons where the brighter-fatter effect should kick in more.
    obj *= 200
    obj.drawImage(im1, method='fft', sensor=silicon, rng=rng)
    obj.drawImage(im2, method='fft', sensor=simple, rng=rng)
    obj.drawImage(im3, method='fft')

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak          radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    # Fluxes should still be fine.
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=1)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=1)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=1)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=1)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=1)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=1)

    # Sizes for 2,3 should be about equal, but 1 should be larger.
    sigma_r = 1./np.sqrt(obj.flux) * im1.scale
    print('check |r2-r3| = %f <? %f'%(np.abs(r2-r3), 2.*sigma_r))
    np.testing.assert_allclose(r2, r3, atol=2.*sigma_r)
    print('check |r1-r3| = %f >? %f'%(np.abs(r1-r3), 2.*sigma_r))
    assert r1-r3 > 2.*sigma_r
Exemple #5
0
def test_silicon():
    """Test the basic construction and use of the SiliconSensor class.
    """

    # Note: Use something quite small in terms of npixels so the B/F effect kicks in without
    # requiring a ridiculous number of photons
    obj = galsim.Gaussian(flux=10000, sigma=0.3)

    # We'll draw the same object using SiliconSensor, Sensor, and the default (sensor=None)
    im1 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=silicon
    im2 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=simple
    im3 = galsim.ImageD(64, 64, scale=0.3)  # Will use sensor=None

    rng1 = galsim.BaseDeviate(5678)
    rng2 = galsim.BaseDeviate(5678)
    rng3 = galsim.BaseDeviate(5678)

    silicon = galsim.SiliconSensor(rng=rng1, diffusion_factor=0.0)
    simple = galsim.Sensor()

    # Start with photon shooting, since that's more straightforward.
    obj.drawImage(im1, method='phot', poisson_flux=False, sensor=silicon, rng=rng1)
    obj.drawImage(im2, method='phot', poisson_flux=False, sensor=simple, rng=rng2)
    obj.drawImage(im3, method='phot', poisson_flux=False, rng=rng3)

    # First, im2 and im3 should be exactly equal.
    np.testing.assert_array_equal(im2.array, im3.array)

    # im1 should be similar, but not equal
    np.testing.assert_almost_equal(im1.array/obj.flux, im2.array/obj.flux, decimal=2)

    # Now use a different seed for 3 to see how much of the variation is just from randomness.
    rng3.seed(234241)
    obj.drawImage(im3, method='phot', poisson_flux=False, rng=rng3)

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak          radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    # Fluxes should all equal obj.flux
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=6)

    # Sizes are all about equal since flux is not large enough for B/F to be significant
    # Variance of Irr for Gaussian with Poisson noise is
    # Var(Irr) = Sum(I r^4) = 4Irr [using Gaussian kurtosis = 8sigma^2, Irr = 2sigma^2]
    # r = sqrt(Irr/flux), so sigma(r) = 1/2 r sqrt(Var(Irr))/Irr = 1/sqrt(flux)
    # Use 2sigma for below checks.
    sigma_r = 1. / np.sqrt(obj.flux) * im1.scale
    np.testing.assert_allclose(r1, r2, atol=2.*sigma_r)
    np.testing.assert_allclose(r2, r3, atol=2.*sigma_r)

    # Repeat with 100X more photons where the brighter-fatter effect should kick in more.
    obj *= 100
    rng1 = galsim.BaseDeviate(5678)
    rng2 = galsim.BaseDeviate(5678)
    rng3 = galsim.BaseDeviate(5678)

    obj.drawImage(im1, method='phot', poisson_flux=False, sensor=silicon, rng=rng1)
    obj.drawImage(im2, method='phot', poisson_flux=False, sensor=simple, rng=rng2)
    obj.drawImage(im3, method='phot', poisson_flux=False, rng=rng3)

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak          radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    # Fluxes should still be fine.
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=6)

    # Sizes for 2,3 should be about equal, but 1 should be larger.
    sigma_r = 1. / np.sqrt(obj.flux) * im1.scale
    print('check |r2-r3| = %f <? %f'%(np.abs(r2-r3), 2.*sigma_r))
    np.testing.assert_allclose(r2, r3, atol=2.*sigma_r)
    print('check r1 - r3 = %f > %f due to brighter-fatter'%(r1-r2,sigma_r))
    assert r1 - r3 > 2*sigma_r

    # Check that it is really responding to flux, not number of photons.
    # Using fewer shot photons will mean each one encapsulates several electrons at once.
    obj.drawImage(im1, method='phot', n_photons=1000, poisson_flux=False, sensor=silicon,
                  rng=rng1)
    obj.drawImage(im2, method='phot', n_photons=1000, poisson_flux=False, sensor=simple,
                  rng=rng2)
    obj.drawImage(im3, method='phot', n_photons=1000, poisson_flux=False, rng=rng3)

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak          radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=6)

    print('check r1 - r3 = %f > %f due to brighter-fatter'%(r1-r2,sigma_r))
    assert r1 - r3 > 2*sigma_r

    # Can also get the stronger BF effect with the strength parameter.
    obj /= 100  # Back to what it originally was.
    rng1 = galsim.BaseDeviate(5678)
    rng2 = galsim.BaseDeviate(5678)
    rng3 = galsim.BaseDeviate(5678)

    silicon = galsim.SiliconSensor(dir='lsst_itl', strength=100., rng=rng1, diffusion_factor=0.0)
    obj.drawImage(im1, method='phot', poisson_flux=False, sensor=silicon, rng=rng1)
    obj.drawImage(im2, method='phot', poisson_flux=False, sensor=simple, rng=rng2)
    obj.drawImage(im3, method='phot', poisson_flux=False, rng=rng3)

    r1 = im1.calculateMomentRadius(flux=obj.flux)
    r2 = im2.calculateMomentRadius(flux=obj.flux)
    r3 = im3.calculateMomentRadius(flux=obj.flux)
    print('Flux = %.0f:  sum        peak          radius'%obj.flux)
    print('im1:         %.1f     %.2f       %f'%(im1.array.sum(),im1.array.max(), r1))
    print('im2:         %.1f     %.2f       %f'%(im2.array.sum(),im2.array.max(), r2))
    print('im3:         %.1f     %.2f       %f'%(im3.array.sum(),im3.array.max(), r3))

    # Fluxes should still be fine.
    np.testing.assert_almost_equal(im1.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.array.sum(), obj.flux, decimal=6)
    np.testing.assert_almost_equal(im1.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im2.added_flux, obj.flux, decimal=6)
    np.testing.assert_almost_equal(im3.added_flux, obj.flux, decimal=6)

    # Sizes for 2,3 should be about equal, but 1 should be larger.
    sigma_r = 1. / np.sqrt(obj.flux) * im1.scale
    print('check |r2-r3| = %f <? %f'%(np.abs(r2-r3), 2.*sigma_r))
    np.testing.assert_allclose(r2, r3, atol=2.*sigma_r)
    print('check r1 - r3 = %f > %f due to brighter-fatter'%(r1-r2,sigma_r))
    assert r1 - r3 > 2*sigma_r / 100

    # Check the construction with an explicit directory.
    s0 = galsim.SiliconSensor(rng=rng1)
    dir = os.path.join(galsim.meta_data.share_dir, 'sensors', 'lsst_itl')
    s1 = galsim.SiliconSensor(dir=dir, strength=1.0, rng=rng1, diffusion_factor=1.0, qdist=3,
                              nrecalc=10000)
    assert s0 == s1
    s1 = galsim.SiliconSensor(dir, 1.0, rng1, 1.0, 3, 10000)
    assert s0 == s1
    s2 = galsim.SiliconSensor(rng=rng1, dir='lsst_itl')
    assert s0 == s2
    s3 = galsim.SiliconSensor(rng=rng1, strength=10.)
    s4 = galsim.SiliconSensor(rng=rng1, diffusion_factor=2.0)
    s5 = galsim.SiliconSensor(rng=rng1, qdist=4)
    s6 = galsim.SiliconSensor(rng=rng1, nrecalc=12345)
    s7 = galsim.SiliconSensor(dir=dir, strength=1.5, rng=rng1, diffusion_factor=1.3, qdist=4,
                              nrecalc=12345)
    for s in [ s3, s4, s5, s6, s7 ]:
        assert silicon != s
        assert s != s0

    do_pickle(s0)
    do_pickle(s1)
    do_pickle(s7)

    try:
        np.testing.assert_raises(IOError, galsim.SiliconSensor, dir='junk')
        np.testing.assert_raises(IOError, galsim.SiliconSensor, dir='output')
        np.testing.assert_raises(RuntimeError, galsim.SiliconSensor, rng=3.4)
        np.testing.assert_raises(TypeError, galsim.SiliconSensor, 'lsst_itl', 'hello')
    except ImportError:
        print('The assert_raises tests require nose')
Exemple #6
0
def test_treerings():
    """Test the addition of tree rings.
    compare image positions with the simple sensor to
    a SiliconSensor with no tree rings and six
    different additions of tree rings.
    """
    # Set up the different sensors.
    treering_amplitude = 0.5
    rng1 = galsim.BaseDeviate(5678)
    rng2 = galsim.BaseDeviate(5678)
    rng3 = galsim.BaseDeviate(5678)
    rng4 = galsim.BaseDeviate(5678)
    rng5 = galsim.BaseDeviate(5678)
    rng6 = galsim.BaseDeviate(5678)
    rng7 = galsim.BaseDeviate(5678)
    sensor0 = galsim.Sensor()
    sensor1 = galsim.SiliconSensor(rng=rng1)
    tr2 = galsim.SiliconSensor.simple_treerings(treering_amplitude, 250.)
    sensor2 = galsim.SiliconSensor(rng=rng2, treering_func=tr2,
                                   treering_center=galsim.PositionD(-1000.0,0.0))
    sensor3 = galsim.SiliconSensor(rng=rng3, treering_func=tr2,
                                   treering_center=galsim.PositionD(1000.0,0.0))
    sensor4 = galsim.SiliconSensor(rng=rng4, treering_func=tr2,
                                   treering_center=galsim.PositionD(0.0,-1000.0))
    tr5 = galsim.SiliconSensor.simple_treerings(treering_amplitude, 250., r_max=2000, dr=1.)
    sensor5 = galsim.SiliconSensor(rng=rng5, treering_func=tr5,
                                   treering_center=galsim.PositionD(0.0,1000.0))

    # Now test the ability to read in a user-specified function
    tr6 = galsim.LookupTable.from_func(treering_function, x_min=0.0, x_max=2000)
    sensor6 = galsim.SiliconSensor(rng=rng6, treering_func=tr6,
                                   treering_center=galsim.PositionD(-1000.0,0.0))

    # Now test the ability to read in a lookup table
    tr7 = galsim.LookupTable.from_file('tree_ring_lookup.dat', amplitude=treering_amplitude)
    sensor7 = galsim.SiliconSensor(rng=rng7, treering_func=tr7,
                                   treering_center=galsim.PositionD(-1000.0,0.0))

    sensors = [sensor0, sensor1, sensor2, sensor3, sensor4, sensor5, sensor6, sensor7]
    names = ['Sensor()',
             'SiliconSensor, no TreeRings',
             'SiliconSensor, TreeRingCenter= (-1000,0)',
             'SiliconSensor, TreeRingCenter= (1000,0)',
             'SiliconSensor, TreeRingCenter= (0,-1000)',
             'SiliconSensor, TreeRingCenter= (0,1000)',
             'SiliconSensor with specified function, TreeRingCenter= (-1000,0)',
             'SiliconSensor with lookup table, TreeRingCenter= (-1000,0)']
    centers = [None, None,
               (-1000,0),
               (1000,0),
               (0,-1000),
               (0,1000),
               (-1000,0),
               (-1000,0)]

    init_flux = 200000
    obj = galsim.Gaussian(flux=init_flux, sigma=0.3)

    im = galsim.ImageD(10,10, scale=0.3)
    obj.drawImage(im)
    ref_mom = galsim.utilities.unweighted_moments(im)
    print('Reference Moments Mx, My = (%f, %f):'%(ref_mom['Mx'], ref_mom['My']))

    for sensor, name, center in zip(sensors, names, centers):
        im = galsim.ImageD(10,10, scale=0.3)
        obj.drawImage(im, method='phot', sensor=sensor)
        mom = galsim.utilities.unweighted_moments(im)
        print('%s, Moments Mx, My = (%f, %f):'%(name, mom['Mx'], mom['My']))
        if center is None:
            np.testing.assert_almost_equal(mom['Mx'], ref_mom['Mx'], decimal = 1)
            np.testing.assert_almost_equal(mom['My'], ref_mom['My'], decimal = 1)
        else:
            np.testing.assert_almost_equal(ref_mom['Mx'] + treering_amplitude * center[0] / 1000,
                                           mom['Mx'], decimal=1)
            np.testing.assert_almost_equal(ref_mom['My'] + treering_amplitude * center[1] / 1000,
                                           mom['My'], decimal=1)

    assert_raises(TypeError, galsim.SiliconSensor, treering_func=lambda x:np.cos(x))
    assert_raises(TypeError, galsim.SiliconSensor, treering_func=tr7, treering_center=(3,4))
yesbf_filename = 'spectrum_sim_gausshalf_bftrue.fits'

rng = galsim.BaseDeviate(5678)

# multiply the total flux by a scalar
scalar = 0.5
# transform the spectrum image into a galsim object
spectrum_image = galsim.Image(smeared_spectrum2d * scalar,
                              scale=1.0)  # scale is pixel/pixel
# interpolate the image so GalSim can manipulate it
spectrum_interpolated = galsim.InterpolatedImage(spectrum_image)
spectrum_interpolated.drawImage(
    image=spectrum_image,
    method='phot',
    # center=(15, 57),
    sensor=galsim.Sensor())

print('image center:', spectrum_image.center)
print('image true center:', spectrum_image.true_center)
spectrum_image.write(nobf_filename)
if display:
    galsim_sensor_image = spectrum_image.array.copy()

# now do it again, but with the BF effect
# spectrum_image = galsim.Image(smeared_spectrum2d, scale=.25)  # scale is angstroms/pixel
# interpolate the image so GalSim can manipulate it
# spectrum_interpolated = galsim.InterpolatedImage(spectrum_image)

spectrum_interpolated.drawImage(
    image=spectrum_image,
    method='phot',
Exemple #8
0
def main(argv):
    """
    About as simple as it gets:
      - Use a circular Gaussian profile for the galaxy.
      - Convolve it by a circular Gaussian PSF.
      - Add Gaussian noise to the image.
    """
    # In non-script code, use getLogger(__name__) at module scope instead.
    logging.basicConfig(format="%(message)s", level=logging.INFO, stream=sys.stdout)
    logger = logging.getLogger("demo1")

    # MOD crank up the flux
    data_index = 1
    flux_factor = 1
    gal_flux = 2.e6 * flux_factor   # total counts on the image
    gal_sigma = 2.     # arcsec
    psf_sigma = 1.     # arcsec
    # MOD decrease the resolution
    pixel_scale = 2.0  # arcsec / pixel
    noise = 30.        # standard deviation of the counts in each pixel

    logger.info('Starting demo script 1 using:')
    logger.info('    - circular Gaussian galaxy (flux = %.1e, sigma = %.1f),',gal_flux,gal_sigma)
    logger.info('    - circular Gaussian PSF (sigma = %.1f),',psf_sigma)
    logger.info('    - pixel scale = %.2f,',pixel_scale)
    logger.info('    - Gaussian noise (sigma = %.2f).',noise)

    # Define the galaxy profile
    gal = galsim.Gaussian(flux=gal_flux, sigma=gal_sigma)
    logger.debug('Made galaxy profile')

    # Define the PSF profile
    psf = galsim.Gaussian(flux=1., sigma=psf_sigma) # PSF flux should always = 1
    logger.debug('Made PSF profile')

    # Final profile is the convolution of these
    # Can include any number of things in the list, all of which are convolved
    # together to make the final flux profile.
    final = galsim.Convolve([gal, psf])
    logger.debug('Convolved components into final profile')

    # Draw the image with a particular pixel scale, given in arcsec/pixel.
    # The returned image has a member, added_flux, which is gives the total flux actually added to
    # the image.  One could use this value to check if the image is large enough for some desired
    # accuracy level.  Here, we just ignore it.
    # MOD use the 'phot' method with an ideal sensor
    image = final.drawImage(scale=pixel_scale,
                            method='phot',
                            sensor=galsim.Sensor())
    logger.debug('Made image of the profile: flux = %f, added_flux = %f', gal_flux, image.added_flux)
    file_name = os.path.join('output', 'spot_nobf_' + str(data_index) + '.fits')

    # Write the image to a file
    if not os.path.isdir('output'):
        os.mkdir('output')

    # Note: if the file already exists, this will overwrite it.
    image.write(file_name)
    logger.info('Wrote Ideal image to %r' % file_name)  # using %r adds quotes around filename for us

    results = image.FindAdaptiveMom()

    logger.info('HSM reports that the image has observed shape and size:')
    logger.info('    e1 = %.3f, e2 = %.3f, sigma = %.3f (pixels)', results.observed_shape.e1,
                results.observed_shape.e2, results.moments_sigma)
    logger.info('Expected values in the limit that pixel response and noise are negligible:')
    logger.info('    e1 = %.3f, e2 = %.3f, sigma = %.3f', 0.0, 0.0,
                math.sqrt(gal_sigma ** 2 + psf_sigma ** 2) / pixel_scale)

    sensor_name = 'lsst_e2v_50_32'
    # sensor_name = 'lsst_itl_50_32'
    # MOD use the 'phot' method, with a e2v sensor
    rng = galsim.BaseDeviate(5678)
    image = final.drawImage(scale=pixel_scale,
                            method='phot',
                            sensor=galsim.SiliconSensor(name=sensor_name,
                                                        transpose=False,
                                                        rng=rng,
                                                        diffusion_factor=0.0))
    logger.debug('Made image of the profile: flux = %f, added_flux = %f', gal_flux, image.added_flux)
    file_name = os.path.join('output', 'spot_' + sensor_name + '_bf_' + str(data_index) + '.fits')

    # Add Gaussian noise to the image with specified sigma
    # MOD actually, don't add any noise
    # image.addNoise(galsim.GaussianNoise(sigma=noise))
    # logger.debug('Added Gaussian noise')
    logger.debug('No noise added')

    # Write the image to a file
    if not os.path.isdir('output'):
        os.mkdir('output')

    # Note: if the file already exists, this will overwrite it.
    image.write(file_name)
    logger.info('Wrote bfed image to %r' % file_name)  # using %r adds quotes around filename for us

    results = image.FindAdaptiveMom()

    logger.info('HSM reports that the image has observed shape and size:')
    logger.info('    e1 = %.3f, e2 = %.3f, sigma = %.3f (pixels)', results.observed_shape.e1,
                results.observed_shape.e2, results.moments_sigma)
    logger.info('Expected values in the limit that pixel response and noise are negligible:')
    logger.info('    e1 = %.3f, e2 = %.3f, sigma = %.3f', 0.0, 0.0,
                math.sqrt(gal_sigma**2 + psf_sigma**2)/pixel_scale)