radiation_pattern¶
-
pyoof.aperture.radiation_pattern(K_coeff, I_coeff, d_z, wavel, illum_func, telgeo, resolution, box_factor)[source] [edit on github]¶ Spectrum or (field) radiation pattern, \(F(u, v)\), it is the FFT2 computation of the aperture distribution, \(\underline{E_\mathrm{a}} (x, y)\), in a rectangular grid. Passing the majority of arguments to the aperture distribution except the FFT2 resolution.
Parameters: K_coeff :
ndarrayConstants coefficients, \(K_{n\ell}\), for each of them there is only one Zernike circle polynomial, \(U^\ell_n(\varrho, \varphi)\). The coefficients are between \([-2, 2]\).
I_coeff :
ndarrayList which contains 4 parameters, the illumination amplitude, \(A_{E_\mathrm{a}}\), the illumination taper, \(c_\mathrm{dB}\) and the two coordinate offset, \((x_0, y_0)\). The illumination coefficients must be listed as follows,
I_coeff = [i_amp, c_dB, x0, y0].d_z :
floatRadial offset, \(d_z\), added to the sub-reflector in meters. This characteristic measurement adds the classical interference pattern to the beam maps, normalized squared (field) radiation pattern, which is an out-of-focus property. It is usually of the order of centimeters.
wavel :
floatWavelength, \(\lambda\), of the observation in meters.
illum_func :
functionIllumination function, \(E_\mathrm{a}(x, y)\), to be evaluated with the key I_coeff. The illumination functions available are
illum_pedestalandillum_gauss.telgeo :
listList that contains the blockage distribution, optical path difference (OPD) function, and the primary radius (
float) in meters. The list must have the following order,telego = [block_dist, opd_func, pr].resolution :
intFast Fourier Transform resolution for a rectangular grid. The input value has to be greater or equal to the telescope resolution and with power of 2 for faster FFT processing. It is recommended a value higher than
resolution = 2 ** 8.box_factor :
intRelated to the FFT resolution (resolution key), defines the image pixel size level. It depends on the primary radius,
pr, of the telescope, e.g. abox_factor = 5returnsx = np.linspace(-5 * pr, 5 * pr, resolution), an array to be used in the FFT2 (fft2).Returns: u_shift :
ndarray\(u\) wave-vector in 1 / m units. It belongs to the \(x\) coordinate in meters from the aperture distribution, \(\underline{E_\mathrm{a}}(x, y)\).
v_shift :
ndarray\(v\) wave-vector in 1 / m units. It belongs to the y coordinate in meters from the aperture distribution, \(\underline{E_\mathrm{a}}(x, y)\).
F_shift :
ndarrayOutput from the FFT2 (
fft2), \(F(u, v)\), unnormalized solution in a grid, defined by resolution and box_factor keys.Notes
The (field) radiation pattern is the direct Fourier Transform in two dimensions of the aperture distribution, hence,
\[F(u, v) = \mathcal{F} \left[ \underline{E_\mathrm{a}}(x, y) \right].\]