Understanding AIA temperature response with aiapy

Understanding AIA temperature response with aiapy#

In this example, we will use the aiapy package to convolve the Fe XVIII contribution function with the instrument response function to understand the temperature sensitivity of the 94 angstrom channel.

import astropy.constants as const
import astropy.units as u
import matplotlib.pyplot as plt
import numpy as np

from aiapy.response import Channel
from astropy.visualization import quantity_support
from scipy.interpolate import interp1d

import fiasco

# sphinx_gallery_thumbnail_number = 2

quantity_support()
<astropy.visualization.units.quantity_support.<locals>.MplQuantityConverter object at 0x7f5c91dc3f90>

First, create the Channel object for the 94 angstrom channel and compute the wavelength response.

Files Downloaded:   0%|          | 0/1 [00:00<?, ?file/s]

aiapy.aia_V8_all_fullinst.genx:   0%|          | 0.00/5.43M [00:00<?, ?B/s]

aiapy.aia_V8_all_fullinst.genx:   6%|▌         | 306k/5.43M [00:00<00:01, 3.05MB/s]

aiapy.aia_V8_all_fullinst.genx:  60%|██████    | 3.27M/5.43M [00:00<00:00, 18.7MB/s]


Files Downloaded: 100%|██████████| 1/1 [00:00<00:00,  2.90file/s]
Files Downloaded: 100%|██████████| 1/1 [00:00<00:00,  2.89file/s]

Plot the wavelength response

plt.plot(ch.wavelength, response)
plt.xlim(ch.channel + [-4, 4]*u.angstrom)
plt.axvline(x=ch.channel, ls='--', color='k')
plt.title(f'{ch.channel.to_string(format="latex")} Wavelength Response')
plt.show()
$94 \; \mathrm{\mathring{A}}$ Wavelength Response

Next, we construct for the Ion object for Fe 18. Note that by default we use the coronal abundances of Feldman et al. [FMS+92].

temperature = 10.**(np.arange(4.5, 8, 0.05)) * u.K
fe18 = fiasco.Ion('Fe 18', temperature)

Compute the contribution function,

\[G(n,T,\lambda) = 0.83\mathrm{Ab}\frac{hc}{\lambda}N_{\lambda}A_{\lambda}f\frac{1}{n}\]

for each transition of Fe 18 at constant pressure of \(10^{15}\) K \(\mathrm{cm}^{-3}\). Note that we use the couple_density_to_temperature keyword such that temperature and density vary along the same axis.

constant_pressure = 1e15 * u.K * u.cm**(-3)
density = constant_pressure / fe18.temperature
g = fe18.contribution_function(density, include_protons=False,
                               couple_density_to_temperature=True)

Get the corresponding transition wavelengths

transitions = fe18.transitions.wavelength[~fe18.transitions.is_twophoton]
energy = const.h * const.c / transitions / u.photon

Interpolate the response function to the transition wavelengths

In order to compute the temperature response function \(K_c\), we integrate the wavelength response function \(R_c\), weighted by the contribution function \(G(\lambda,T)\), over wavelength,

\[K_c(T) = \int_0^{\infty}\mathrm{d}\lambda\,G(\lambda,T)R_c(\lambda)\]

We divide by \(hc/\lambda\) in order to convert from units of energy to photons. The factor of \(0.83\) is a relative scaling factor for the abundance of H and is not included in the temperature responses computed by aia_get_response.pro. For more more information on the AIA wavelength response calculation, see Boerner et al. [BEL+12].

K = (g / energy * response_transitions).sum(axis=2) / (4*np.pi*u.steradian) / 0.83
K = K.squeeze().to('cm5 ct pix-1 s-1')

Plot the effective temperature response function for the 94 channel.

plt.plot(temperature, K)
plt.xscale('log')
plt.yscale('log')
plt.ylim(1e-29, 1e-26)
plt.xlim(1e5, 3e7)
plt.title(f'{ch.channel.to_string(format="latex")} Response for {fe18.ion_name_roman}')
plt.show()
$94 \; \mathrm{\mathring{A}}$ Response for Fe XVIII

Note that this represents only the contribution of Fe 18 to the 94 angstrom bandpass. In order to construct the full response functions, one would need to repeat this procedure for every ion listed in the CHIANTI atomic database and for all AIA channels.

Total running time of the script: (0 minutes 11.244 seconds)

Gallery generated by Sphinx-Gallery