multilayer_fit#

esis.flights.f1.optics.primaries.materials.multilayer_fit()[source]#

Modify the result of multilayer_witness_fit() to have the correct substrate.

The witness sample has a silicon substrate and the primary has a glass substrate. So this function changes the substrate from silicon to glass and modifies the roughness to be consistent with the actual primary mirror.

Examples

Plot the theoretical reflectivity of this multilayer stack vs. the theoretical reflectivity of multilayer_witness_fit().

import numpy as np
import matplotlib.pyplot as plt
import astropy.units as u
import astropy.visualization
import named_arrays as na
import optika
from esis.flights.f1.optics import primaries

# Define a grid of wavelength samples
wavelength = na.geomspace(100, 1000, axis="wavelength", num=1001) * u.AA

# Define a grid of incidence angles
angle = 4 * u.deg

# Define the light rays incident on the multilayer stack
rays = optika.rays.RayVectorArray(
    wavelength=wavelength,
    direction=na.Cartesian3dVectorArray(
        x=np.sin(angle),
        y=0,
        z=np.cos(angle),
    ),
)

# Initialize the multilayer stacks
multilayer_witness_fit = primaries.materials.multilayer_witness_fit()
multilayer_fit = primaries.materials.multilayer_fit()

# Define the vector normal to the multilayer stack
normal = na.Cartesian3dVectorArray(0, 0, -1)

# Compute the reflectivity of the multilayer for the given incident rays
reflectivity_witness = multilayer_witness_fit.efficiency(rays, normal)
reflectivity_fit = multilayer_fit.efficiency(rays, normal)

# Plot the reflectivities as a function of wavelength
with astropy.visualization.quantity_support():
    fig, ax = plt.subplots(constrained_layout=True)
    na.plt.plot(
        wavelength,
        reflectivity_witness,
        ax=ax,
        axis="wavelength",
        label="witness fit",
    );
    na.plt.plot(
        wavelength,
        reflectivity_fit,
        ax=ax,
        axis="wavelength",
        label="primary fit",
    );
    ax.set_xlabel(f"wavelength ({ax.get_xlabel()})");
    ax.set_ylabel("reflectivity");
    ax.legend();

# Print the fitted multilayer stack
multilayer_fit
MultilayerMirror(
    layers=[
        Layer(
            chemical='SiO2',
            thickness=1. nm,
            interface=ErfInterfaceProfile(
                width=3.44068978 nm,
            ),
            kwargs_plot={'color': 'tab:blue', 'alpha': 0.3},
            x_label=None,
        ),
        Layer(
            chemical='SiC',
            thickness=21.27827884 nm,
            interface=ErfInterfaceProfile(
                width=3.44068978 nm,
            ),
            kwargs_plot={'color': 'tab:blue', 'alpha': 0.5},
            x_label=None,
        ),
        Layer(
            chemical='Cr',
            thickness=5. nm,
            interface=ErfInterfaceProfile(
                width=1.77840007 nm,
            ),
            kwargs_plot={'color': 'tab:orange', 'alpha': 0.5},
            x_label=None,
        ),
    ],
    substrate=Layer(
        chemical='SiO2',
        thickness=30. mm,
        interface=ErfInterfaceProfile(
            width=5.53 Angstrom,
        ),
        kwargs_plot={'color': 'gray', 'alpha': 0.5},
        x_label=None,
    ),
)
../_images/esis.flights.f1.optics.primaries.materials.multilayer_fit_0_2.png

Plot the reflectivity of this multilayer stack as a function of incidence angle.

import numpy as np
import matplotlib.pyplot as plt
import astropy.units as u
import astropy.visualization
import named_arrays as na
import optika
from esis.flights.f1.optics import primaries

# Define a grid of wavelength samples
wavelength = na.geomspace(100, 1000, axis="wavelength", num=1001) * u.AA

# Define a grid of incidence angles
angle = na.linspace(0, 20, axis="angle", num=3) * u.deg

# Define the light rays incident on the multilayer stack
rays = optika.rays.RayVectorArray(
    wavelength=wavelength,
    direction=na.Cartesian3dVectorArray(
        x=np.sin(angle),
        y=0,
        z=np.cos(angle),
    ),
)

# Initialize the multilayer stack
multilayer = primaries.materials.multilayer_fit()

# Compute the reflectivity of the multilayer for the given incident rays
reflectivity = multilayer.efficiency(
    rays=rays,
    normal=na.Cartesian3dVectorArray(0, 0, -1),
)

# Plot the reflectivity of the multilayer as a function of wavelength
with astropy.visualization.quantity_support():
    fig, ax = plt.subplots(constrained_layout=True)
    na.plt.plot(
        wavelength,
        reflectivity,
        ax=ax,
        axis="wavelength",
        label=angle,
    );
    ax.set_xlabel(f"wavelength ({ax.get_xlabel()})");
    ax.set_ylabel("reflectivity");
    ax.legend();
../_images/esis.flights.f1.optics.primaries.materials.multilayer_fit_1_0.png

Return type:

MultilayerMirror