Beispiel #1
0
def test_rayGrid():
    dist = 10.0
    length = 10.0
    xcos = 0.1
    ycos = 0.2
    zcos = -np.sqrt(1.0 - xcos**2 - ycos**2)
    nside = 10
    wavelength = 500e-9
    flux = 1.2
    medium = batoid.ConstMedium(1.2)

    rays = batoid.rayGrid(
        dist, length, xcos, ycos, zcos, nside, wavelength, flux, medium,
        lattice=True
    )
    assert rays.monochromatic == True
    # Check that all rays are perpendicular to v
    ray0 = rays[0]
    for ray in rays:
        dr = ray.r - ray0.r
        dp = np.dot(dr, ray0.v)
        np.testing.assert_allclose(dp, 0.0, atol=1e-14, rtol=0.0)
        np.testing.assert_allclose(ray.wavelength, wavelength)
        np.testing.assert_allclose(ray.flux, flux)
        np.testing.assert_allclose(np.linalg.norm(ray.v), 1./1.2)
        np.testing.assert_allclose(ray.v[0]*1.2, xcos)
        np.testing.assert_allclose(ray.v[1]*1.2, ycos)

    # Check that ray that intersects at origin is initially dist away.
    # Need the ray that is in the middle in both dimensions...
    idx = np.ravel_multi_index((nside//2, nside//2), (nside, nside))
    rays.propagateInPlace(dist*1.2)
    np.testing.assert_equal(rays[idx].r, [0,0,0])
    # but mean position won't be the origin, since lattice implies off-center
    assert np.linalg.norm(np.mean(rays.r, axis=0)) > 0.5

    # Now try again with lattice flag set to False
    rays = batoid.rayGrid(
        dist, length, xcos, ycos, zcos, nside, wavelength, flux, medium,
        lattice=False
    )
    # "Central" ray will not intersect origin in this case, but average of
    # all rays should be the origin
    idx = np.ravel_multi_index((nside//2, nside//2), (nside, nside))
    rays.propagateInPlace(dist*1.2)
    assert np.linalg.norm(rays[idx].r) > 0.1
    np.testing.assert_allclose(np.linalg.norm(np.mean(rays.r, axis=0)), 0.0, rtol=0, atol=1e-14)

    # If we use an odd nside, then both the central point and the mean will be the origin.
    nside = 11
    rays = batoid.rayGrid(
        dist, length, xcos, ycos, zcos, nside, wavelength, flux, medium,
        lattice=False
    )
    idx = np.ravel_multi_index((nside//2, nside//2), (nside, nside))
    rays.propagateInPlace(dist*1.2)
    np.testing.assert_allclose(rays[idx].r, [0,0,0], rtol=0, atol=1e-14)
    np.testing.assert_allclose(np.linalg.norm(np.mean(rays.r, axis=0)), 0.0, rtol=0, atol=1e-14)
def test_refraction_chromatic():
    import random
    random.seed(577215664)
    wavelength1 = 500e-9
    wavelength2 = 600e-9
    flux = 1.0

    plane = batoid.Plane()
    filename = os.path.join(batoid.datadir, "media", "silica_dispersion.txt")
    wave, n = np.genfromtxt(filename).T
    wave *= 1e-6  # micron -> meters
    table = batoid.Table(wave, n, batoid.Table.Interpolant.linear)
    silica = batoid.TableMedium(table)
    air = batoid.Air()

    thx, thy = 0.001, 0.0001
    dirCos = batoid.utils.gnomonicToDirCos(thx, thy)
    rv1 = batoid.rayGrid(10.0, 1., dirCos[0], dirCos[1], -dirCos[2], 2,
                         wavelength1, flux, silica)
    rv2 = batoid.rayGrid(10.0, 1., dirCos[0], dirCos[1], -dirCos[2], 2,
                         wavelength2, flux, silica)
    rays = []
    for ray in rv1:
        rays.append(ray)
    for ray in rv2:
        rays.append(ray)
    rvCombined = batoid.RayVector(rays)

    rv1r = plane.refract(rv1, silica, air)
    rv2r = plane.refract(rv2, silica, air)
    assert rv1r != rv2r
    rays = []
    for ray in rv1r:
        rays.append(ray)
    for ray in rv2r:
        rays.append(ray)
    rvrCombined1 = batoid.RayVector(rays)

    rvrCombined2 = plane.refract(rvCombined, silica, air)

    assert rvrCombined1 == rvrCombined2

    # Check in-place
    plane.refractInPlace(rv1, silica, air)
    plane.refractInPlace(rv2, silica, air)
    assert rv1 != rv2
    plane.refractInPlace(rvCombined, silica, air)
    rays = []
    for ray in rv1:
        rays.append(ray)
    for ray in rv2:
        rays.append(ray)
    rvCombined2 = batoid.RayVector(rays)

    assert rvCombined == rvCombined2
Beispiel #3
0
def zernike(optic, theta_x, theta_y, wavelength, nx=32, jmax=22, eps=0.0, sphereRadius=None):
    import galsim

    dirCos = gnomicToDirCos(theta_x, theta_y)
    rays = batoid.rayGrid(
        optic.dist, optic.pupilSize,
        dirCos[0], dirCos[1], -dirCos[2],
        nx, wavelength, 1.0, optic.inMedium
    )

    # Propagate to t=optic.dist where rays are close to being
    # in the entrance pupil.  (TODO: actually intersect EP here?)
    rays.propagateInPlace(optic.dist)

    orig_x = np.array(rays.x).reshape(nx,nx)
    orig_y = np.array(rays.y).reshape(nx,nx)
    wf = wavefront(optic, theta_x, theta_y, wavelength, nx=nx, sphereRadius=sphereRadius)
    wfarr = wf.array
    w = ~wfarr.mask

    basis = galsim.zernike.zernikeBasis(
            jmax, orig_x[w], orig_y[w],
            R_outer=optic.pupilSize/2, R_inner=optic.pupilSize/2*eps
    )
    coefs, _, _, _ = np.linalg.lstsq(basis.T, wfarr[w], rcond=-1)

    return np.array(coefs)
Beispiel #4
0
def zernike(optic, wavelength, theta_x, theta_y, jmax=22, nx=32, eps=0.0):
    import galsim.zernike as zern

    xcos = np.sin(theta_x)
    ycos = np.sin(theta_y)
    zcos = -np.sqrt(1.0 - xcos**2 - ycos**2)

    rays = batoid.rayGrid(optic.dist, optic.pupilSize, xcos, ycos, zcos, nx,
                          wavelength, optic.inMedium)

    orig_x = np.array(rays.x)
    orig_y = np.array(rays.y)

    wf = wavefront(optic, wavelength, rays=rays)

    w = ~wf.mask

    basis = zern.zernikeBasis(jmax,
                              orig_x[w],
                              orig_y[w],
                              R_outer=optic.pupilSize / 2,
                              R_inner=optic.pupilSize / 2 * eps)
    coefs, _, _, _ = np.linalg.lstsq(basis.T, wf[w])

    return coefs
Beispiel #5
0
def test_traceReverse():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    telescope = batoid.Optic.fromYaml("HSC.yaml")

    init_rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9, 1.0, batoid.ConstMedium(1.0))
    forward_rays = telescope.trace(init_rays)

    # Now, turn the result rays around and trace backwards
    forward_rays = forward_rays.propagatedToTime(40.0)
    reverse_rays = batoid.RayVector(
        [batoid.Ray(r.r, -r.v, -r.t, r.wavelength) for r in forward_rays]
    )

    final_rays = telescope.traceReverse(reverse_rays)
    # propagate all the way to t=0
    final_rays = final_rays.propagatedToTime(0.0)
    final_rays.toCoordSysInPlace(batoid.globalCoordSys)

    w = np.where(np.logical_not(final_rays.vignetted))[0]
    for idx in w:
        np.testing.assert_allclose(init_rays[idx].x, final_rays[idx].x)
        np.testing.assert_allclose(init_rays[idx].y, final_rays[idx].y)
        np.testing.assert_allclose(init_rays[idx].z, final_rays[idx].z)
        np.testing.assert_allclose(init_rays[idx].vx, -final_rays[idx].vx)
        np.testing.assert_allclose(init_rays[idx].vy, -final_rays[idx].vy)
        np.testing.assert_allclose(init_rays[idx].vz, -final_rays[idx].vz)
        np.testing.assert_allclose(final_rays[idx].t, 0)
Beispiel #6
0
def test_optic():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9, 1.0, batoid.ConstMedium(1.0))

    nrays = len(rays)
    print("Tracing {} rays.".format(nrays))
    t_fast = 0.0
    t_slow = 0.0

    telescope = batoid.Optic.fromYaml("HSC.yaml")
    do_pickle(telescope)

    t0 = time.time()

    rays_fast = telescope.trace(rays)
    t1 = time.time()
    rays_slow = batoid.RayVector([telescope.trace(r) for r in rays])
    t2 = time.time()

    assert rays_fast == rays_slow
    t_fast = t1 - t0
    t_slow = t2 - t1
    print("Fast trace: {:5.3f} s".format(t_fast))
    print("            {} rays per second".format(int(nrays/t_fast)))
    print("Slow trace: {:5.3f} s".format(t_slow))
    print("            {} rays per second".format(int(nrays/t_slow)))
Beispiel #7
0
def test_optic():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9, 1.0,
                          batoid.ConstMedium(1.0))

    nrays = len(rays)
    print("Tracing {} rays.".format(nrays))
    t_fast = 0.0
    t_slow = 0.0

    fn = os.path.join(batoid.datadir, "HSC", "HSC.yaml")
    config = yaml.load(open(fn))
    telescope = batoid.parse.parse_optic(config['opticalSystem'])
    do_pickle(telescope)

    t0 = time.time()

    rays_fast, _ = telescope.trace(rays)
    t1 = time.time()
    rays_slow = batoid.RayVector([telescope.trace(r)[0] for r in rays])
    t2 = time.time()

    assert rays_fast == rays_slow
    t_fast = t1 - t0
    t_slow = t2 - t1
    print("Fast trace: {:5.3f} s".format(t_fast))
    print("            {} rays per second".format(int(nrays / t_fast)))
    print("Slow trace: {:5.3f} s".format(t_slow))
    print("            {} rays per second".format(int(nrays / t_slow)))
Beispiel #8
0
def test_traceReverse():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    fn = os.path.join(batoid.datadir, "HSC", "HSC.yaml")
    config = yaml.load(open(fn))
    telescope = batoid.parse.parse_optic(config['opticalSystem'])

    init_rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9,
                               1.0, batoid.ConstMedium(1.0))
    forward_rays, _ = telescope.trace(init_rays, outCoordSys=batoid.CoordSys())

    # Now, turn the result rays around and trace backwards
    forward_rays = forward_rays.propagatedToTime(40.0)
    reverse_rays = batoid.RayVector(
        [batoid.Ray(r.r, -r.v, -r.t, r.wavelength) for r in forward_rays])

    final_rays, _ = telescope.traceReverse(reverse_rays,
                                           outCoordSys=batoid.CoordSys())
    # propagate all the way to t=0
    final_rays = final_rays.propagatedToTime(0.0)

    w = np.where(np.logical_not(final_rays.vignetted))[0]
    for idx in w:
        np.testing.assert_allclose(init_rays[idx].x, final_rays[idx].x)
        np.testing.assert_allclose(init_rays[idx].y, final_rays[idx].y)
        np.testing.assert_allclose(init_rays[idx].z, final_rays[idx].z)
        np.testing.assert_allclose(init_rays[idx].vx, -final_rays[idx].vx)
        np.testing.assert_allclose(init_rays[idx].vy, -final_rays[idx].vy)
        np.testing.assert_allclose(init_rays[idx].vz, -final_rays[idx].vz)
        np.testing.assert_allclose(final_rays[idx].t, 0)
Beispiel #9
0
def initialize(ngrid=25, theta_x=1.):
    DESI_fn = os.path.join(batoid.datadir, 'DESI', 'DESI.yaml')
    config = yaml.safe_load(open(DESI_fn))
    telescope = batoid.parse.parse_optic(config['opticalSystem'])
    dirCos = batoid.utils.gnomonicToDirCos(np.deg2rad(theta_x), 0.)
    rays = batoid.rayGrid(telescope.dist, telescope.pupilSize, dirCos[0],
                          dirCos[1], -dirCos[2], ngrid, 500e-9, 1.0,
                          telescope.inMedium)
    return telescope, telescope.traceFull(rays)
Beispiel #10
0
def initialize(ngrid=25, theta_x=1.):
    telescope = batoid.Optic.fromYaml("DESI.yaml")
    dirCos = batoid.utils.gnomonicToDirCos(np.deg2rad(theta_x), 0.)
    rays = batoid.rayGrid(
        telescope.backDist, telescope.pupilSize,
        dirCos[0], dirCos[1], dirCos[2],
        ngrid, 500e-9, 1.0, telescope.inMedium
    )
    return telescope, telescope.traceFull(rays)
Beispiel #11
0
def wavefront(optic, theta_x, theta_y, wavelength, nx=32, sphereRadius=None):
    """Compute wavefront.

    Parameters
    ----------
    optic : batoid.Optic
        Optic for which to compute wavefront.
    theta_x, theta_y : float
        Field of incoming rays (gnomic projection)
    wavelength : float
        Wavelength of incoming rays
    nx : int, optional
        Size of ray grid to generate to compute wavefront.  Default: 32
    sphereRadius : float, optional
        Radius of reference sphere in meters.  If None, then use optic.sphereRadius.

    Returns
    -------
    wavefront : batoid.Lattice
        A batoid.Lattice object containing the wavefront values in waves and
        the primitive lattice vectors of the entrance pupil grid in meters.
    """
    dirCos = gnomicToDirCos(theta_x, theta_y)
    rays = batoid.rayGrid(
        optic.dist, optic.pupilSize,
        dirCos[0], dirCos[1], -dirCos[2],
        nx, wavelength, 1.0, optic.inMedium
    )

    if sphereRadius is None:
        sphereRadius = optic.sphereRadius

    outCoordSys = batoid.CoordSys()
    optic.traceInPlace(rays, outCoordSys=outCoordSys)
    w = np.where(1-rays.vignetted)[0]
    point = np.mean(rays.r[w], axis=0)

    # We want to place the vertex of the reference sphere one radius length away from the
    # intersection point.  So transform our rays into that coordinate system.
    transform = batoid.CoordTransform(
            outCoordSys, batoid.CoordSys(point+np.array([0,0,sphereRadius])))
    transform.applyForwardInPlace(rays)

    sphere = batoid.Sphere(-sphereRadius)
    sphere.intersectInPlace(rays)

    w = np.where(1-rays.vignetted)[0]
    # Should potentially try to make the reference time w.r.t. the chief ray instead of the mean
    # of the good (unvignetted) rays.
    t0 = np.mean(rays.t[w])

    arr = np.ma.masked_array((t0-rays.t)/wavelength, mask=rays.vignetted).reshape(nx, nx)
    primitiveVectors = np.vstack([[optic.pupilSize/nx, 0], [0, optic.pupilSize/nx]])
    return batoid.Lattice(arr, primitiveVectors)
Beispiel #12
0
def fpPosition(optic,
               theta_x,
               theta_y,
               wavelength,
               nx=32,
               projection='postel'):
    dirCos = fieldToDirCos(theta_x, theta_y, projection=projection)
    rays = batoid.rayGrid(optic.dist, optic.pupilSize, dirCos[0], dirCos[1],
                          -dirCos[2], nx, wavelength, 1.0, optic.inMedium)
    optic.traceInPlace(rays)
    rays.trimVignettedInPlace()
    return np.mean(rays.x), np.mean(rays.y)
Beispiel #13
0
def test_inplace():
    rays = batoid.rayGrid(10.0, 10.0, 0.1, 0.1, 1.0, 1024, 500e-9, 1.0, batoid.ConstMedium(1.2))
    print("propagating {} rays.".format(len(rays)))
    rays2 = batoid.RayVector(rays)

    t0 = time.time()
    outrays = rays.propagatedToTime(11.2)
    t1 = time.time()
    rays2.propagateInPlace(11.2)
    t2 = time.time()

    print("immutable propagation took {:6.3f} seconds.".format(t1-t0))
    print("in-place propagation took  {:6.3f} seconds.".format(t2-t1))

    assert outrays == rays2
Beispiel #14
0
def wavefront(optic,
              wavelength,
              theta_x=0,
              theta_y=0,
              nx=32,
              rays=None,
              saveRays=False,
              sphereRadius=None):
    if rays is None:
        xcos = np.sin(theta_x)
        ycos = np.sin(theta_y)
        zcos = -np.sqrt(1.0 - xcos**2 - ycos**2)

        rays = batoid.rayGrid(optic.dist, optic.pupilSize, xcos, ycos, zcos,
                              nx, wavelength, optic.inMedium)
    if saveRays:
        rays = batoid.RayVector(rays)
    if sphereRadius is None:
        sphereRadius = optic.sphereRadius

    outCoordSys = batoid.CoordSys()
    optic.traceInPlace(rays, outCoordSys=outCoordSys)
    goodRays = batoid._batoid.trimVignetted(rays)
    point = np.array(
        [np.mean(goodRays.x),
         np.mean(goodRays.y),
         np.mean(goodRays.z)])

    # We want to place the vertex of the reference sphere one radius length away from the
    # intersection point.  So transform our rays into that coordinate system.
    transform = batoid.CoordTransform(
        outCoordSys, batoid.CoordSys(point + np.array([0, 0, sphereRadius])))
    transform.applyForwardInPlace(rays)

    sphere = batoid.Sphere(-sphereRadius)
    sphere.intersectInPlace(rays)
    goodRays = batoid._batoid.trimVignetted(rays)
    # Should potentially try to make the reference time w.r.t. the chief ray instead of the mean
    # of the good (unvignetted) rays.
    t0 = np.mean(goodRays.t0)

    ts = rays.t0[:]
    isV = rays.isVignetted[:]
    ts -= t0
    ts /= wavelength
    wf = np.ma.masked_array(ts, mask=isV)
    return wf
Beispiel #15
0
def test_traceFull():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9, 1.0, batoid.ConstMedium(1.0))

    nrays = len(rays)
    print("Tracing {} rays.".format(nrays))

    telescope = batoid.Optic.fromYaml("HSC.yaml")

    tf = telescope.traceFull(rays)
    rays = telescope.trace(rays)

    assert rays == tf['D']['out']
Beispiel #16
0
def test_rSplit():
    for i in range(100):
        R = np.random.normal(0.7, 0.8)
        conic = np.random.uniform(-2.0, 1.0)
        ncoef = np.random.randint(0, 4)
        coefs = [np.random.normal(0, 1e-10) for i in range(ncoef)]
        asphere = batoid.Asphere(R, conic, coefs)

        rays = batoid.rayGrid(10, 2*R, 0.0, 0.0, -1.0, 16, 500e-9, 1.0, batoid.Air())
        coating = batoid.SimpleCoating(0.9, 0.1)
        reflectedRays = asphere.reflect(rays, coating)
        m1 = batoid.Air()
        m2 = batoid.ConstMedium(1.1)
        refractedRays = asphere.refract(rays, m1, m2, coating)
        reflectedRays2, refractedRays2 = asphere.rSplit(rays, m1, m2, coating)

        assert reflectedRays == reflectedRays2
        assert refractedRays == refractedRays2
Beispiel #17
0
def dkdu(optic, wavelength, theta_x, theta_y, nx=16):
    """Calculate derivative of outgoing ray k-vector with respect to incoming ray
    pupil coordinate.


    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.

    Returns
    -------
    dkdu : (2, 2), ndarray
        Jacobian transformation matrix for converting between (kx, ky) of rays impacting the focal
        plane and initial field angle.

    Notes
    -----
        This is essentially the plate scale in radians per meter.
    """
    rays = batoid.rayGrid(optic.dist, optic.pupilSize, theta_x, theta_y, -1.0,
                          nx, wavelength, optic.inMedium)

    pupilRays = batoid._batoid.propagatedToTimesMany(rays,
                                                     np.zeros_like(rays.x))
    ux = np.array(pupilRays.x)
    uy = np.array(pupilRays.y)

    optic.traceInPlace(rays)
    w = np.where(1 - rays.isVignetted)[0]
    ux = ux[w]
    uy = uy[w]

    kx = rays.kx[w]
    ky = rays.ky[w]

    soln = bivariate_fit(ux, uy, kx, ky)
    return soln[1:]
Beispiel #18
0
def dkdu(optic, theta_x, theta_y, wavelength, nx=16):
    """Calculate derivative of outgoing ray k-vector with respect to incoming ray
    pupil coordinate.


    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians (gnomic tangent plane projection)
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.

    Returns
    -------
    dkdu : (2, 2), ndarray
        Jacobian transformation matrix for converting between (kx, ky) of rays impacting the focal
        plane and initial field angle (gnomic tangent plane projection).
    """
    dirCos = gnomicToDirCos(theta_x, theta_y)
    rays = batoid.rayGrid(
        optic.dist, optic.pupilSize,
        dirCos[0], dirCos[1], -dirCos[2],
        nx, wavelength, 1.0, optic.inMedium
    )

    pupilRays = rays.propagatedToTime(0.0)
    ux = np.array(pupilRays.x)
    uy = np.array(pupilRays.y)

    optic.traceInPlace(rays)
    w = np.where(1-rays.vignetted)[0]
    ux = ux[w]
    uy = uy[w]

    kx = rays.kx[w]
    ky = rays.ky[w]

    soln = bilinear_fit(ux, uy, kx, ky)
    return soln[1:]
Beispiel #19
0
def newFFTPSF(optic, wavelength, theta_x, theta_y, nx=64, pupilFactor=2):
    xcos = np.sin(theta_x)
    ycos = np.sin(theta_y)
    zcos = -np.sqrt(1.0 - xcos**2 - ycos**2)

    # We could incorporate pupilFactor here, but since we know that all the
    # additional rays should vignette, we'll delay until later.
    rays = batoid.rayGrid(optic.dist, optic.pupilSize, xcos, ycos, zcos, nx,
                          wavelength, optic.inMedium)

    wf = wavefront(optic, wavelength, rays=rays).reshape(nx, nx)

    padSize = nx * pupilFactor
    expwf = np.zeros((padSize, padSize), dtype=np.complex128)
    start = padSize // 2 - nx // 2
    stop = padSize // 2 + nx // 2
    expwf[start:stop, start:stop][~wf.mask] = np.exp(2j * np.pi * wf[~wf.mask])
    arr = np.abs(np.fft.fftshift(np.fft.fft2(np.fft.fftshift(expwf))))**2

    # Now the tricky part, what is the coordinate system for the output array?
    du = optic.pupilSize / nx
    N = nx * pupilFactor
    u1 = np.array([1, 0]) * du
    u2 = np.array([0, 1]) * du

    # Convert pupil coords to Fourier focal plane coords
    _dkdu = dkdu(optic, wavelength, theta_x, theta_y, nx=nx)
    q1 = _dkdu[0, 0] * u1 + _dkdu[0, 1] * u2
    q2 = _dkdu[1, 0] * u1 + _dkdu[1, 1] * u2

    # Use reciprical lattice vectors for focal plane coords...
    norm = 2 * np.pi / (q1[0] * q2[1] - q1[1] * q2[0]) / N
    x1 = norm * np.array([q2[1], -q2[0]])
    x2 = norm * np.array([q1[1], -q1[0]])

    # Also finally convert to sky coords
    dthdx = dthdr(optic, wavelength, theta_x, theta_y, nx=nx)
    th1 = dthdx[0, 0] * x1 + dthdx[0, 1] * x2
    th2 = dthdx[1, 0] * x1 + dthdx[1, 1] * x2

    # th1 and th2 now indicate the field angles corresponding to rows and columns.
    return th1, th2, arr
Beispiel #20
0
def dkdu(optic, theta_x, theta_y, wavelength, nx=16, projection='postel'):
    """Calculate derivative of outgoing ray k-vector with respect to incoming ray
    pupil coordinate.

    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.
    projection : {'postel', 'zemax', 'gnomonic', 'stereographic', 'lambert', 'orthographic'}
        Projection used to convert field angle to direction cosines.

    Returns
    -------
    dkdu : (2, 2), ndarray
        Jacobian transformation matrix for converting between (kx, ky) of rays impacting the focal
        plane and initial field angle (gnomonic tangent plane projection).
    """
    dirCos = fieldToDirCos(theta_x, theta_y, projection=projection)
    rays = batoid.rayGrid(optic.backDist, optic.pupilSize, dirCos[0],
                          dirCos[1], dirCos[2], nx, wavelength, 1.0,
                          optic.inMedium)

    pupilRays = rays.copy().propagate(0.0)
    ux = np.array(pupilRays.x)
    uy = np.array(pupilRays.y)

    optic.trace(rays)
    w = np.where(1 - rays.vignetted)[0]
    ux = ux[w]
    uy = uy[w]

    kx = rays.kx[w]
    ky = rays.ky[w]

    soln = bilinear_fit(ux, uy, kx, ky)
    return soln[1:]
Beispiel #21
0
def test_traceFull():
    if __name__ == '__main__':
        nside = 128
    else:
        nside = 32

    rays = batoid.rayGrid(20, 12.0, 0.005, 0.005, -1.0, nside, 500e-9, 1.0,
                          batoid.ConstMedium(1.0))

    nrays = len(rays)
    print("Tracing {} rays.".format(nrays))

    fn = os.path.join(batoid.datadir, "HSC", "HSC.yaml")
    config = yaml.load(open(fn))
    telescope = batoid.parse.parse_optic(config['opticalSystem'])

    tf = telescope.traceFull(rays)
    rays, _ = telescope.trace(rays)

    assert rays == tf[-1]['out']
Beispiel #22
0
def test_thread():
    telescope = batoid.Optic.fromYaml("HSC.yaml")

    rayGrid = batoid.rayGrid(telescope.backDist, telescope.pupilSize, 0.0, 0.0,
                             -1.0, 32, 750e-9, 1.0, telescope.inMedium)

    batoid._batoid.setNThread(4)
    assert batoid._batoid.getNThread() == 4

    rays4 = telescope.trace(rayGrid)

    batoid._batoid.setNThread(2)
    assert batoid._batoid.getNThread() == 2

    rays2 = telescope.trace(rayGrid)

    batoid._batoid.setNThread(1)
    assert batoid._batoid.getNThread() == 1

    rays1 = telescope.trace(rayGrid)

    assert rays1 == rays2 == rays4
Beispiel #23
0
def test_thread():
    fn = os.path.join(batoid.datadir, "HSC", "HSC.yaml")
    config = yaml.load(open(fn))
    telescope = batoid.parse.parse_optic(config['opticalSystem'])

    rayGrid = batoid.rayGrid(telescope.dist, telescope.pupilSize, 0.0, 0.0,
                             -1.0, 32, 750e-9, 1.0, telescope.inMedium)

    batoid._batoid.setNThread(4)
    assert batoid._batoid.getNThread() == 4

    rays4, _ = telescope.trace(rayGrid)

    batoid._batoid.setNThread(2)
    assert batoid._batoid.getNThread() == 2

    rays2, _ = telescope.trace(rayGrid)

    batoid._batoid.setNThread(1)
    assert batoid._batoid.getNThread() == 1

    rays1, _ = telescope.trace(rayGrid)

    assert rays1 == rays2 == rays4
Beispiel #24
0
def huygensPSF(optic,
               theta_x=None,
               theta_y=None,
               wavelength=None,
               nx=None,
               projection='postel',
               dx=None,
               dy=None,
               nxOut=None):
    r"""Compute a PSF via the Huygens construction.

    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float, optional
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.
    projection : {'postel', 'zemax', 'gnomonic', 'stereographic', 'lambert', 'orthographic'}
        Projection used to convert field angle to direction cosines.
    dx, dy : float, optional
        Lattice scales to use for PSF evaluation locations.  Default, use fftPSF lattice.

    Returns
    -------
    psf : batoid.Lattice
        The PSF.

    Notes
    -----
    The Huygens construction is to evaluate the PSF as

    I(x) \propto \Sum_u exp(i phi(u)) exp(i k(u).r)

    The u are assumed to uniformly sample the entrance pupil, but not include any rays that get
    vignetted before they reach the focal plane.  The phis are the phases of the exit rays evaluated
    at a single arbitrary time.  The k(u) indicates the conversion of the uniform entrance pupil
    samples into nearly (though not exactly) uniform samples in k-space of the output rays.

    The output locations where the PSF is evaluated are governed by dx, dy and nx.  If dx and dy are
    None, then the same lattice as in fftPSF will be used.  If dx and dy are scalars, then a lattice
    with primitive vectors [dx, 0] and [0, dy] will be used.  If dx and dy are 2-vectors, then those
    will be the primitive vectors of the output lattice.
    """
    from numbers import Real

    if dx is None:
        primitiveU = np.array([[optic.pupilSize / nx, 0],
                               [0, optic.pupilSize / nx]])
        primitiveK = dkdu(optic,
                          theta_x,
                          theta_y,
                          wavelength,
                          projection=projection).dot(primitiveU)
        pad_factor = 2
        primitiveX = np.vstack(
            reciprocalLatticeVectors(primitiveK[0], primitiveK[1],
                                     pad_factor * nx))
    elif isinstance(dx, Real):
        if dy is None:
            dy = dx
        primitiveX = np.vstack([[dx, 0], [0, dy]])
        pad_factor = 1
    else:
        primitiveX = np.vstack([dx, dy])
        pad_factor = 1

    if nxOut is None:
        nxOut = nx

    dirCos = fieldToDirCos(theta_x, theta_y, projection=projection)
    rays = batoid.rayGrid(optic.backDist,
                          optic.pupilSize,
                          dirCos[0],
                          dirCos[1],
                          dirCos[2],
                          nx,
                          wavelength=wavelength,
                          flux=1,
                          medium=optic.inMedium,
                          lattice=True)

    amplitudes = np.zeros((nxOut * pad_factor, nxOut * pad_factor),
                          dtype=np.complex128)
    out = batoid.Lattice(
        np.zeros((nxOut * pad_factor, nxOut * pad_factor), dtype=float),
        primitiveX)

    rays = optic.trace(rays)
    rays.trimVignetted()
    # Need transpose to conform to numpy [y,x] ordering convention
    xs = out.coords[..., 0].T + np.mean(rays.x)
    ys = out.coords[..., 1].T + np.mean(rays.y)
    zs = np.zeros_like(xs)

    points = np.concatenate([aux[..., None] for aux in (xs, ys, zs)], axis=-1)
    time = rays[0].t
    for idx in np.ndindex(amplitudes.shape):
        amplitudes[idx] = rays.sumAmplitude(points[idx], time)
    out.array = np.abs(amplitudes)**2
    return out
Beispiel #25
0
def wavefront(optic,
              theta_x,
              theta_y,
              wavelength,
              nx=32,
              projection='postel',
              sphereRadius=None,
              lattice=False,
              _addedWF=None):
    """Compute wavefront.

    Parameters
    ----------
    optic : batoid.Optic
        Optic for which to compute wavefront.
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float, optional
        Wavelength in meters
    nx : int, optional
        Size of ray grid to generate to compute wavefront.  Default: 32
    projection : {'postel', 'zemax', 'gnomonic', 'stereographic', 'lambert', 'orthographic'}
        Projection used to convert field angle to direction cosines.
    sphereRadius : float, optional
        Radius of reference sphere in meters.  If None, then use optic.sphereRadius.
    lattice : bool, optional
        If true, then decenter the grid so it spans (-N/2, N/2+1), as appropriate
        for Fourier transforms.

    Returns
    -------
    wavefront : batoid.Lattice
        A batoid.Lattice object containing the wavefront values in waves and
        the primitive lattice vectors of the entrance pupil grid in meters.
    """
    dirCos = fieldToDirCos(theta_x, theta_y, projection=projection)
    rays = batoid.rayGrid(optic.backDist,
                          optic.pupilSize,
                          dirCos[0],
                          dirCos[1],
                          dirCos[2],
                          nx,
                          wavelength,
                          1.0,
                          optic.inMedium,
                          lattice=lattice)

    if sphereRadius is None:
        sphereRadius = optic.sphereRadius

    optic.trace(rays)
    w = np.where(1 - rays.vignetted)[0]
    point = np.mean(rays.r[w], axis=0)

    # We want to place the vertex of the reference sphere one radius length away from the
    # intersection point.  So transform our rays into that coordinate system.
    targetCoordSys = rays.coordSys.shiftLocal(point +
                                              np.array([0, 0, sphereRadius]))
    rays.toCoordSys(targetCoordSys)

    sphere = batoid.Sphere(-sphereRadius)
    sphere.intersect(rays)

    w = np.where(1 - rays.vignetted)[0]
    # Should potentially try to make the reference time w.r.t. the chief ray instead of the mean
    # of the good (unvignetted) rays.
    t0 = np.mean(rays.t[w])

    arr = np.ma.masked_array((t0 - rays.t) / wavelength,
                             mask=rays.vignetted).reshape(nx, nx)
    if _addedWF is not None:
        arr += _addedWF
    primitiveVectors = np.vstack([[optic.pupilSize / nx, 0],
                                  [0, optic.pupilSize / nx]])
    return batoid.Lattice(arr, primitiveVectors)
Beispiel #26
0
def drdth(optic, theta_x, theta_y, wavelength, nx=16, projection='postel'):
    """Calculate derivative of focal plane coord with respect to field angle.

    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.
    projection : {'postel', 'zemax', 'gnomonic', 'stereographic', 'lambert', 'orthographic'}
        Projection used to convert field angle to direction cosines.

    Returns
    -------
    drdth : (2, 2), ndarray
        Jacobian transformation matrix for converting between (theta_x, theta_y)
        and (x, y) on the focal plane.

    Notes
    -----
        This is the Jacobian of pixels -> tangent plane, (and importantly, not pixels -> ra/dec).
        It should be *close* to the inverse plate scale though, especially near the center of the
        tangent plane projection.
    """
    # We just use a finite difference approach here.
    dth = 1e-5

    # Make direction cosine vectors
    nominalCos = fieldToDirCos(theta_x, theta_y, projection=projection)
    dthxCos = fieldToDirCos(theta_x + dth, theta_y, projection=projection)
    dthyCos = fieldToDirCos(theta_x, theta_y + dth, projection=projection)

    rays = batoid.rayGrid(optic.backDist,
                          optic.pupilSize,
                          nominalCos[0],
                          nominalCos[1],
                          nominalCos[2],
                          nx,
                          wavelength=wavelength,
                          flux=1,
                          medium=optic.inMedium)

    rays_x = batoid.rayGrid(optic.backDist,
                            optic.pupilSize,
                            dthxCos[0],
                            dthxCos[1],
                            dthxCos[2],
                            nx,
                            wavelength=wavelength,
                            flux=1,
                            medium=optic.inMedium)

    rays_y = batoid.rayGrid(optic.backDist,
                            optic.pupilSize,
                            dthyCos[0],
                            dthyCos[1],
                            dthyCos[2],
                            nx,
                            wavelength=wavelength,
                            flux=1,
                            medium=optic.inMedium)

    optic.trace(rays)
    optic.trace(rays_x)
    optic.trace(rays_y)

    rays.trimVignetted()
    rays_x.trimVignetted()
    rays_y.trimVignetted()

    # meters / radian
    drx_dthx = (np.mean(rays_x.x) - np.mean(rays.x)) / dth
    drx_dthy = (np.mean(rays_y.x) - np.mean(rays.x)) / dth
    dry_dthx = (np.mean(rays_x.y) - np.mean(rays.y)) / dth
    dry_dthy = (np.mean(rays_y.y) - np.mean(rays.y)) / dth

    return np.array([[drx_dthx, drx_dthy], [dry_dthx, dry_dthy]])
Beispiel #27
0
def drdth(optic, theta_x, theta_y, wavelength, nx=16):
    """Calculate derivative of focal plane coord with respect to field angle.

    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians (gnomic tangent plane projection)
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.

    Returns
    -------
    drdth : (2, 2), ndarray
        Jacobian transformation matrix for converting between (theta_x, theta_y)
        and (x, y) on the focal plane.

    Notes
    -----
        This is the Jacobian of pixels -> tangent plane, (and importantly, not pixels -> ra/dec).
        It should be *close* to the inverse plate scale though, especially near the center of the
        tangent plane projection.
    """
    # We just use a finite difference approach here.
    dth = 1e-5

    # Make direction cosine vectors
    nominalCos = gnomicToDirCos(theta_x, theta_y)
    dthxCos = gnomicToDirCos(theta_x + dth, theta_y)
    dthyCos = gnomicToDirCos(theta_x, theta_y+ dth)

    # Flip the dirCos z-components so rays are headed downwards
    rays = batoid.rayGrid(optic.dist, optic.pupilSize,
        nominalCos[0], nominalCos[1], -nominalCos[2],
        nx, wavelength=wavelength, flux=1, medium=optic.inMedium)

    rays_x = batoid.rayGrid(optic.dist, optic.pupilSize,
        dthxCos[0], dthxCos[1], -dthxCos[2],
        nx, wavelength=wavelength, flux=1, medium=optic.inMedium)

    rays_y = batoid.rayGrid(optic.dist, optic.pupilSize,
        dthyCos[0], dthyCos[1], -dthyCos[2],
        nx, wavelength=wavelength, flux=1, medium=optic.inMedium)

    optic.traceInPlace(rays)
    optic.traceInPlace(rays_x)
    optic.traceInPlace(rays_y)

    rays.trimVignettedInPlace()
    rays_x.trimVignettedInPlace()
    rays_y.trimVignettedInPlace()

    # meters / radian
    drx_dthx = (np.mean(rays_x.x) - np.mean(rays.x))/dth
    drx_dthy = (np.mean(rays_y.x) - np.mean(rays.x))/dth
    dry_dthx = (np.mean(rays_x.y) - np.mean(rays.y))/dth
    dry_dthy = (np.mean(rays_y.y) - np.mean(rays.y))/dth

    return np.array([[drx_dthx, dry_dthx], [drx_dthy, dry_dthy]])
Beispiel #28
0
def drdth(optic, wavelength, theta_x, theta_y, nx=16):
    """Calculate derivative of focal plane coord with respect to field angle.

    Parameters
    ----------
    optic : batoid.Optic
        Optical system
    theta_x, theta_y : float
        Field angle in radians
    wavelength : float
        Wavelength in meters
    nx : int, optional
        Size of ray grid to use.

    Returns
    -------
    drdth : (2, 2), ndarray
        Jacobian transformation matrix for converting between (theta_x, theta_y)
        and (x, y) on the focal plane.

    Notes
    -----
        This is essentially the inverse plate scale in meters per radian.
    """
    # We just use a finite difference approach here.
    dth = 1e-5

    nominalCos = [
        np.sin(theta_x),
        np.sin(theta_y),
        -np.sqrt(1.0 - np.sin(theta_x)**2 - np.sin(theta_y)**2)
    ]
    dthxCos = [
        np.sin(theta_x + dth),
        np.sin(theta_y),
        -np.sqrt(1.0 - np.sin(theta_x + dth)**2 - np.sin(theta_y)**2)
    ]
    dthyCos = [
        np.sin(theta_x),
        np.sin(theta_y + dth),
        -np.sqrt(1.0 - np.sin(theta_x)**2 - np.sin(theta_y + dth)**2)
    ]

    rays = batoid.rayGrid(optic.dist,
                          optic.pupilSize,
                          nominalCos[0],
                          nominalCos[1],
                          nominalCos[2],
                          nx,
                          wavelength=wavelength,
                          medium=optic.inMedium)

    rays_x = batoid.rayGrid(optic.dist,
                            optic.pupilSize,
                            dthxCos[0],
                            dthxCos[1],
                            dthxCos[2],
                            nx,
                            wavelength=wavelength,
                            medium=optic.inMedium)

    rays_y = batoid.rayGrid(optic.dist,
                            optic.pupilSize,
                            dthyCos[0],
                            dthyCos[1],
                            dthyCos[2],
                            nx,
                            wavelength=wavelength,
                            medium=optic.inMedium)

    optic.traceInPlace(rays)
    optic.traceInPlace(rays_x)
    optic.traceInPlace(rays_y)

    batoid.trimVignettedInPlace(rays)
    batoid.trimVignettedInPlace(rays_x)
    batoid.trimVignettedInPlace(rays_y)

    # meters / radian
    drx_dthx = (np.mean(rays_x.x) - np.mean(rays.x)) / dth
    drx_dthy = (np.mean(rays_y.x) - np.mean(rays.x)) / dth
    dry_dthx = (np.mean(rays_x.y) - np.mean(rays.y)) / dth
    dry_dthy = (np.mean(rays_y.y) - np.mean(rays.y)) / dth

    return np.array([[drx_dthx, dry_dthx], [drx_dthy, dry_dthy]])