Atmospheric Turbulence

Finite Phase Screens

Creation of phase screens of a defined size with Von Karmen Statistics.

aotools.turbulence.phasescreen.ft_phase_screen(r0, N, delta, L0, l0, FFT=None, seed=None)[source]

Creates a random phase screen with Von Karmen statistics. (Schmidt 2010)

Parameters:
  • r0 (float) – r0 parameter of scrn in metres
  • N (int) – Size of phase scrn in pxls
  • delta (float) – size in Metres of each pxl
  • L0 (float) – Size of outer-scale in metres
  • l0 (float) – inner scale in metres
Returns:

numpy array representing phase screen

Return type:

ndarray

aotools.turbulence.phasescreen.ft_sh_phase_screen(r0, N, delta, L0, l0, FFT=None, seed=None)[source]

Creates a random phase screen with Von Karmen statistics with added sub-harmonics to augment tip-tilt modes. (Schmidt 2010)

Parameters:
  • r0 (float) – r0 parameter of scrn in metres
  • N (int) – Size of phase scrn in pxls
  • delta (float) – size in Metres of each pxl
  • L0 (float) – Size of outer-scale in metres
  • l0 (float) – inner scale in metres
Returns:

numpy array representing phase screen

Return type:

ndarray

Infinite Phase Screens

An implementation of the “infinite phase screen”, as deduced by Francois Assemat and Richard W. Wilson, 2006.

class aotools.turbulence.infinitephasescreen.PhaseScreenVonKarman(nx_size, pixel_scale, r0, L0, random_seed=None, n_columns=2)[source]

Bases: aotools.turbulence.infinitephasescreen.PhaseScreen

A “Phase Screen” for use in AO simulation with Von Karmon statistics.

This represents the phase addition light experiences when passing through atmospheric turbulence. Unlike other phase screen generation techniques that translate a large static screen, this method keeps a small section of phase, and extends it as necessary for as many steps as required. This can significantly reduce memory consumption at the expense of more processing power required.

The technique is described in a paper by Assemat and Wilson, 2006. It essentially assumes that there are two matrices, “A” and “B”, that can be used to extend an existing phase screen. A single row or column of new phase can be represented by

X = A.Z + B.b

where X is the new phase vector, Z is some number of columns of the existing screen, and b is a random vector with gaussian statistics.

This object calculates the A and B matrices using an expression of the phase covariance when it is initialised. Calculating A is straightforward through the relationship:

A = Cov_xz . (Cov_zz)^(-1).

B is less trivial.

BB^t = Cov_xx - A.Cov_zx

(where B^t is the transpose of B) is a symmetric matrix, hence B can be expressed as

B = UL,

where U and L are obtained from the svd for BB^t

U, w, U^t = svd(BB^t)

L is a diagonal matrix where the diagonal elements are w^(1/2).

On initialisation an initial phase screen is calculated using an FFT based method. When add_row is called, a new vector of phase is added to the phase screen using nCols columns of previous phase. Assemat & Wilson claim that two columns are adequate for good atmospheric statistics. The phase in the screen data is always accessed as <phasescreen>.scrn.

Parameters:
  • nx_size (int) – Size of phase screen (NxN)
  • pixel_scale (float) – Size of each phase pixel in metres
  • r0 (float) – fried parameter (metres)
  • L0 (float) – Outer scale (metres)
  • n_columns (int, optional) – Number of columns to use to continue screen, default is 2
add_row()

Adds a new row to the phase screen and removes old ones.

scrn

The current phase map held in the PhaseScreen object.

class aotools.turbulence.infinitephasescreen.PhaseScreenKolmogorov(nx_size, pixel_scale, r0, L0, random_seed=None, stencil_length_factor=4)[source]

Bases: aotools.turbulence.infinitephasescreen.PhaseScreen

A “Phase Screen” for use in AO simulation using the Fried method for Kolmogorov turbulence.

This represents the phase addition light experiences when passing through atmospheric turbulence. Unlike other phase screen generation techniques that translate a large static screen, this method keeps a small section of phase, and extends it as neccessary for as many steps as required. This can significantly reduce memory consuption at the expense of more processing power required.

The technique is described in a paper by Assemat and Wilson, 2006 and expanded upon by Fried, 2008. It essentially assumes that there are two matrices, “A” and “B”, that can be used to extend an existing phase screen. A single row or column of new phase can be represented by

X = A.Z + B.b

where X is the new phase vector, Z is some data from the existing screen, and b is a random vector with gaussian statistics.

This object calculates the A and B matrices using an expression of the phase covariance when it is initialised. Calculating A is straightforward through the relationship:

A = Cov_xz . (Cov_zz)^(-1).

B is less trivial.

BB^t = Cov_xx - A.Cov_zx

(where B^t is the transpose of B) is a symmetric matrix, hence B can be expressed as

B = UL,

where U and L are obtained from the svd for BB^t

U, w, U^t = svd(BB^t)

L is a diagonal matrix where the diagonal elements are w^(1/2).

The Z data is taken from points in a “stencil” defined by Fried that samples the entire screen. An additional “reference point” is also considered, that is picked from a point separate from teh stencil and applied on each iteration such that the new phase equation becomes:

On initialisation an initial phase screen is calculated using an FFT based method. When add_row is called, a new vector of phase is added to the phase screen. The phase in the screen data is always accessed as <phasescreen>.scrn.

Parameters:
  • nx_size (int) – Size of phase screen (NxN)
  • pixel_scale (float) – Size of each phase pixel in metres
  • r0 (float) – fried parameter (metres)
  • L0 (float) – Outer scale (metres)
  • random_seed (int, optional) – seed for the random number generator
  • stencil_length_factor (int, optional) – How much longer is the stencil than the desired phase? default is 4
add_row()

Adds a new row to the phase screen and removes old ones.

scrn

The current phase map held in the PhaseScreen object.

Slope Covariance Matrix Generation

Slope covariance matrix routines for AO systems observing through Von Karmon turbulence. Such matrices have a variety of uses, though they are especially useful for creating ‘tomographic reconstructors’ that can reconstruct some ‘psuedo’ WFS measurements in a required direction (where there might be an interesting science target but no guide stars), given some actual measurements in other directions (where the some suitable guide stars are).

Warning

This code has been tested qualitatively and seems OK, but needs more rigorous testing.

class aotools.turbulence.slopecovariance.CovarianceMatrix(n_wfs, pupil_masks, telescope_diameter, subap_diameters, gs_altitudes, gs_positions, wfs_wavelengths, n_layers, layer_altitudes, layer_r0s, layer_L0s, threads=1)[source]

Bases: object

A creator of slope covariance matrices in Von Karmon turbulence, based on the paper by Martin et al, SPIE, 2012.

Given a list of paramters describing an AO WFS system and the atmosphere above the telescope, this class can compute the covariance matrix between all the WFS measurements. This can support LGS sources that exist at a finite altitude. When computing the covariance matrix, Python’s multiprocessing module is used to spread the work between different processes and processing cores.

On initialisation, this class performs some initial calculations and parameter sorting. To create the covariance matrix run the make_covariace_matrix method. This may take some time depending on your paramters…

Parameters:
  • n_wfs (int) – Number of wavefront sensors present.
  • pupil_masks (ndarray) – A map of the pupil for each WFS which is nx_subaps by ny_subaps. 1 if subap active, 0 if not.
  • telescope_diameter (float) – Diameter of the telescope
  • subap_diameters (ndarray) – The diameter of the sub-apertures for each WFS in metres
  • gs_altitudes (ndarray) – Reciprocal (1/metres) of the Guide star alitude for each WFS
  • gs_positions (ndarray) – X,Y position of each WFS in arcsecs. Array shape (Wfs, 2)
  • wfs_wavelengths (ndarray) – Wavelength each WFS observes
  • n_layers (int) – The number of atmospheric turbulence layers
  • layer_altitudes (ndarray) – The altitude of each turbulence layer in meters
  • layer_r0s (ndarray) – The Fried parameter of each turbulence layer
  • layer_L0s (ndarray) – The outer-scale of each layer in metres
  • threads (int, optional) – Number of processes to use for calculation. default is 1
make_covariance_matrix()[source]

Calculate and build the covariance matrix

Returns:Covariance Matrix
Return type:ndarray
make_tomographic_reconstructor(svd_conditioning=0)[source]

Creats a tomohraphic reconstructor from the covariance matrix as in Vidal, 2010. See the documentation for the function create_tomographic_covariance_reconstructor in this module. Assumes the 1st WFS given is the one for which reconstruction is required to.

Parameters:svd_conditioning (float) – Conditioning for the SVD used in inversion.
Returns:A tomohraphic reconstructor.
Return type:ndarray
aotools.turbulence.slopecovariance.calculate_structure_function(phase, nbOfPoint=None, step=None)[source]

Compute the structure function of an 2D array, along the first dimension. SF defined as sf[j]= < (phase - phase_shifted_by_j)^2 > Translated from a YAO function.

Parameters:
  • phase (ndarray, 2d) – 2d-array
  • nbOfPoint (int) – final size of the structure function vector. Default is phase.shape[1] / 4
  • step (int) – step in pixel when computing the sf. (step * sampling_phase) gives the sf sampling in meters. Default is 1
Returns:

values for the structure function of the data.

Return type:

ndarray, float

aotools.turbulence.slopecovariance.calculate_wfs_seperations(n_subaps1, n_subaps2, wfs1_positions, wfs2_positions)[source]

Calculates the seperation between all sub-apertures in two WFSs

Parameters:
  • n_subaps1 (int) – Number of sub-apertures in WFS 1
  • n_subaps2 (int) – Number of sub-apertures in WFS 2
  • wfs1_positions (ndarray) – Array of the X, Y positions of the centre of each sub-aperture with respect to the centre of the telescope pupil
  • wfs2_positions (ndarray) – Array of the X, Y positions of the centre of each sub-aperture with respect to the centre of the telescope pupil
Returns:

2-D Array of sub-aperture seperations

Return type:

ndarray

aotools.turbulence.slopecovariance.create_tomographic_covariance_reconstructor(covariance_matrix, n_onaxis_subaps, svd_conditioning=0)[source]

Calculates a tomographic reconstructor using the method of Vidal, JOSA A, 2010.

Uses a slope covariance matrix to generate a reconstruction matrix that will convert a collection of measurements from WFSs observing in a variety of directions (the “off-axis” directions) to the measurements that would have been observed by a WFS observing in a different direction (the “on-axis” direction. The given covariance matrix must include the covariance between all WFS measurements, including the “psuedo” WFS pointing in the on-axis direction. It is assumed that the covariance of measurements from this on-axis direction are first in the covariance matrix.

To create the tomographic reconstructor it is neccessary to invert the covariance matrix between off-axis WFS measurements. An SVD based psueod inversion is used for this (numpy.linalg.pinv). A conditioning to this SVD may be required to filter potentially unwanted modes.

Parameters:
  • covariance_matrix (ndarray) – A covariance matrix between WFSs
  • n_onaxis_subaps (int) – Number of sub-aperture in on-axis WFS
  • svd_conditioning (float, optional) – Conditioning for SVD inversion (default is 0)

Returns:

aotools.turbulence.slopecovariance.mirror_covariance_matrix(cov_mat)[source]

Mirrors a covariance matrix around the axis of the diagonal.

Parameters:
  • cov_mat (ndarray) – The covariance matrix to mirror
  • n_subaps (ndarray) – Number of sub-aperture in each WFS
aotools.turbulence.slopecovariance.structure_function_kolmogorov(separation, r0)[source]

Compute the Kolmogorov phase structure function

Parameters:
  • separation (ndarray, float) – float or array of data representing separations between points
  • r0 (float) – Fried parameter for atmosphere
Returns:

Structure function for separation(s)

Return type:

ndarray, float

aotools.turbulence.slopecovariance.structure_function_vk(seperation, r0, L0)[source]

Computes the Von Karmon structure function of atmospheric turbulence

Parameters:
  • seperation (ndarray, float) – float or array of data representing seperations between points
  • r0 (float) – Fried parameter for atmosphere
  • L0 (float) – Outer scale of turbulence
Returns:

Structure function for seperation(s)

Return type:

ndarray, float

aotools.turbulence.slopecovariance.wfs_covariance(n_subaps1, n_subaps2, wfs1_positions, wfs2_positions, wfs1_diam, wfs2_diam, r0, L0)[source]

Calculates the covariance between 2 WFSs

Parameters:
  • n_subaps1 (int) – number of sub-apertures in WFS 1
  • n_subaps2 (int) – number of sub-apertures in WFS 1
  • wfs1_positions (ndarray) – Central position of each sub-apeture from telescope centre for WFS 1
  • wfs2_positions (ndarray) – Central position of each sub-apeture from telescope centre for WFS 2
  • wfs1_diam – Diameter of WFS 1 sub-apertures
  • wfs2_diam – Diameter of WFS 2 sub-apertures
  • r0 – Fried parameter of turbulence
  • L0 – Outer scale of turbulence
Returns:

slope covariance of X with X , slope covariance of Y with Y, slope covariance of X with Y

aotools.turbulence.slopecovariance.wfs_covariance_mpwrap(args)[source]

Temporal Power Spectra

Turbulence gradient temporal power spectra calculation and plotting

author:Andrew Reeves
date:September 2016
aotools.turbulence.temporal_ps.calc_slope_temporalps(slope_data)[source]

Calculates the temporal power spectra of the loaded centroid data.

Calculates the Fourier transform over the number frames of the centroid data, then returns the square of the mean of all sub-apertures, for x and y. This is a temporal power spectra of the slopes, and should adhere to a -11/3 power law for the slopes in the wind direction, and -14/3 in the direction tranverse to the wind direction. See Conan, 1995 for more.

The phase slope data should be split into X and Y components, with all X data first, then Y (or visa-versa).

Parameters:slope_data (ndarray) – 2-d array of shape (…, nFrames, nCentroids)
Returns:The temporal power spectra of the data for X and Y components.
Return type:ndarray
aotools.turbulence.temporal_ps.fit_tps(tps, t_axis, D, V_guess=10, f_noise_guess=20, A_guess=9, tps_err=None, plot=False)[source]

Runs minimization routines to get t0.

Parameters:
  • tps (ndarray) – The temporal power spectrum to fit
  • axis (str) – fit parallel (‘par’) or perpendicular (‘per’) to wind direction
  • D (float) – (Sub-)Aperture diameter

Returns:

aotools.turbulence.temporal_ps.get_tps_time_axis(frame_rate, n_frames)[source]
Parameters:
  • frame_rate (float) – Frame rate of detector observing slope gradients (Hz)
  • n_frames – (int): Number of frames used for temporal power spectrum
Returns:

Time values for temporal power spectra plots

Return type:

ndarray

aotools.turbulence.temporal_ps.plot_tps(slope_data, frame_rate)[source]

Generates a plot of the temporal power spectrum/a for a data set of phase gradients

Parameters:
  • slope_data (ndarray) – 2-d array of shape (…, nFrames, nCentroids)
  • frame_rate (float) – Frame rate of detector observing slope gradients (Hz)
Returns:

The computed temporal power spectrum/a, and the time axis data

Return type:

tuple

Atmospheric Parameter Conversions

Functions for converting between different atmospheric parameters,

aotools.turbulence.atmos_conversions.cn2_to_r0(cn2, lamda=5e-07)[source]

Calculates r0 from the integrated Cn2 value

Parameters:
  • cn2 (float) – integrated Cn2 value in m^2/3
  • lamda – wavelength
Returns:

r0 in cm

aotools.turbulence.atmos_conversions.cn2_to_seeing(cn2, lamda=5e-07)[source]

Calculates the seeing angle from the integrated Cn2 value

Parameters:
  • cn2 (float) – integrated Cn2 value in m^2/3
  • lamda – wavelength
Returns:

seeing angle in arcseconds

aotools.turbulence.atmos_conversions.coherenceTime(cn2, v, lamda=5e-07)[source]

Calculates the coherence time from profiles of the Cn2 and wind velocity

Parameters:
  • cn2 (array) – Cn2 profile in m^2/3
  • v (array) – profile of wind velocity, same altitude scale as cn2
  • lamda – wavelength
Returns:

coherence time in seconds

aotools.turbulence.atmos_conversions.isoplanaticAngle(cn2, h, lamda=5e-07)[source]

Calculates the isoplanatic angle from the Cn2 profile

Parameters:
  • cn2 (array) – Cn2 profile in m^2/3
  • h (Array) – Altitude levels of cn2 profile in m
  • lamda – wavelength
Returns:

isoplanatic angle in arcseconds

aotools.turbulence.atmos_conversions.r0_from_slopes(slopes, wavelength, subapDiam)[source]

Measures the value of R0 from a set of WFS slopes.

Uses the equation in Saint Jaques, 1998, PhD Thesis, Appendix A to calculate the value of atmospheric seeing parameter, r0, that would result in the variance of the given slopes.

Parameters:
  • slopes (ndarray) – A 3-d set of slopes in radians, of shape (dimension, nSubaps, nFrames)
  • wavelength (float) – The wavelegnth of the light observed
  • subapDiam (float) –
Returns:

An estimate of r0 for that dataset.

Return type:

float

aotools.turbulence.atmos_conversions.r0_to_cn2(r0, lamda=5e-07)[source]

Calculates integrated Cn2 value from r0

Parameters:
  • r0 (float) – r0 in cm
  • lamda – wavelength
Returns:

integrated Cn2 value in m^2/3

Return type:

cn2 (float)

aotools.turbulence.atmos_conversions.r0_to_seeing(r0, lamda=5e-07)[source]

Calculates the seeing angle from r0

Parameters:
  • r0 (float) – Freid’s parameter in cm
  • lamda – wavelength
Returns:

seeing angle in arcseconds

aotools.turbulence.atmos_conversions.seeing_to_cn2(seeing, lamda=5e-07)[source]

Calculates the integrated Cn2 value from the seeing

Parameters:
  • seeing (float) – seeing in arcseconds
  • lamda – wavelength
Returns:

integrated Cn2 value in m^2/3

aotools.turbulence.atmos_conversions.seeing_to_r0(seeing, lamda=5e-07)[source]

Calculates r0 from seeing

Parameters:
  • seeing (float) – seeing angle in arcseconds
  • lamda – wavelength
Returns:

Freid’s parameter in cm

Return type:

r0 (float)

aotools.turbulence.atmos_conversions.slope_variance_from_r0(r0, wavelength, subapDiam)[source]

Uses the equation in Saint Jaques, 1998, PhD Thesis, Appendix A to calculate the slope variance resulting from a value of r0.

Parameters:
  • r0 (float) – Fried papamerter of turubulence in metres
  • wavelength (float) – Wavelength of light in metres (where 1e-9 is 1nm)
  • subapDiam (float) – Diameter of the aperture in metres
Returns:

The expected slope variance for a given r0 ValueError