ContinuousDeformableMirror
- class poppy.ContinuousDeformableMirror(dm_shape=(10, 10), actuator_spacing=None, influence_func=None, name='DM', include_actuator_print_through=False, actuator_print_through_file=None, actuator_mask_file=None, radius=<Quantity 1. m>, flip_x=False, flip_y=False, include_factor_of_two=False, **kwargs)[source]
Bases:
AnalyticOpticalElement
Generic deformable mirror, of the continuous face sheet variety.
- Parameters:
- dm_shapetuple with 2 elements
Number of actuators across the clear aperture in each dimension
- actuator_spacingfloat or astropy Quantity with dimension length
Spacing between adjacent actuators as seen in that plane
- influence_funcInfluence function filename. Optional. If not supplied,
a Gaussian model approximately representative of an influence function for a Boston MEMS DMs will be used. This parameter can let you provide a more detailed model for your particular mirror.
- radiusfloat or Quantity with dimension length
radius of the clear aperture of the DM. Note, the reflective clear aperture may potentially be larger than the controllable active aperture, depending on the details of your particular hardware. This parameter does not affect any of the actuator models, only the reflective area for wavefront amplitude.
- flip_x, flip_yBool
Flip the orientation of the X or Y axes of the DM actuators. Useful if your device is not oriented such that the origin is at lower left as seen in the pupil.
- include_factor_of_twoBool
include the factor of two due to reflection in the OPD function (optional, default False). If this is set False (default), actuator commands are interpreted as being in units of desired wavefront error directly; the returned WFE will be directly proportional to the requested values (convolved with the actuator response function etc). If this is set to True, then the actuator commands are interpreted as being in physical surface units, and the WFE is therefore a factor of two larger. The returned WFE will be twice the amplitude of the requested values (convolved with the actuator response function etc.)
- Additionally, the standard parameters for shift_x and shift_y can be accepted and
- will be handled by the **kwargs mechanism. Note, rotation is not yet supported for DMs,
- and shifts are currently rounded to integer pixels in the sampled wavefront, rather
- than being applied as floating point values. Shifts are specified in physical units
- of meters transverse motion in the DM plane.
- Note: Keeping track of actuator locations and spacing is subtle. This note
- assumes you are familiar with ‘counting fenceposts’ and off-by-one errors.
- This class follows the convention adopted by Boston Micromachines:
- If there are N actuators across a given distance, there are N-1 spaces between
- actuators across that distance (plus half an actuator space border on the outside,
- which is still controlled by the outermost actuators.). The controlled portion of
- the pupil thus has diameter = (N-1)*actuator_spacing.
- However, for display purposes etc it is often useful if the displayed pupil extends
- over the full N actuators, so we set the `pupil_diam` attribute to N*actuator_spacing,
- or the reflective radius, whichever is larger.
Attributes Summary
DM actuator geometry - i.e. how many actuators per axis .
The commanded surface shape of the deformable mirror, in meters.
Methods Summary
annotate
([marker])Overplot actuator coordinates on some already-existing pupil display
annotate_grid
([linestyle, color])display
([annotate, grid, what, crosshairs])Display an Analytic optic by first computing it onto a grid. Parameters ---------- wavelength : float Wavelength to evaluate this optic's properties at npix : int Number of pixels to use when sampling the optical element. what : str What to display: 'intensity', 'surface' or 'phase', or 'both' ax : matplotlib.Axes instance Axes to display into nrows, row : integers # of rows and row index for subplot display crosshairs : bool Display crosshairs indicating the center? colorbar : bool Show colorbar? colorbar_orientation : bool Desired orientation, horizontal or vertical? Default is horizontal if only 1 row of plots, else vertical opd_vmax : float Max value for OPD image display, in meters. title : string Plot label annotate : bool Draw annotations on plot grid : bool Show the grid for the DM actuator spacing.
display_actuators
([annotate, grid, what, ...])Display the optical surface, viewed as discrete actuators
flatten
()Flatten the DM by setting all actuators to zero piston
get_act_coordinates
([one_d, ...])Y and X coordinates for the actuators
get_opd
(wave)Return the surface shape OPD for the optic.
get_transmission
(wave)Note that this is the amplitude transmission, not the total intensity transmission.
set_actuator
(actx, acty, new_value)Set an individual actuator of the DM. Parameters ------------- actx, acty : integers Coordinates of the actuator you wish to control new_value : float Desired surface height for that actuator, in meters by default or use astropy Units to specify another unit if desired. Example ----------- dm.set_actuator(12, 22, 123.4*u.nm).
set_surface
(new_surface)Set the entire surface shape of the DM.
Attributes Documentation
- dm_shape
DM actuator geometry - i.e. how many actuators per axis
- surface
The commanded surface shape of the deformable mirror, in meters. This is the input to the DM. See the .opd property for the output.
Methods Documentation
- annotate(marker='+', **kwargs)[source]
Overplot actuator coordinates on some already-existing pupil display
- display(annotate=False, grid=False, what='opd', crosshairs=False, *args, **kwargs)[source]
Display an Analytic optic by first computing it onto a grid. Parameters ———- wavelength : float
Wavelength to evaluate this optic’s properties at
- npixint
Number of pixels to use when sampling the optical element.
- whatstr
What to display: ‘intensity’, ‘surface’ or ‘phase’, or ‘both’
- axmatplotlib.Axes instance
Axes to display into
- nrows, rowintegers
# of rows and row index for subplot display
- crosshairsbool
Display crosshairs indicating the center?
- colorbarbool
Show colorbar?
- colorbar_orientationbool
Desired orientation, horizontal or vertical? Default is horizontal if only 1 row of plots, else vertical
- opd_vmaxfloat
Max value for OPD image display, in meters.
- titlestring
Plot label
- annotatebool
Draw annotations on plot
- gridbool
Show the grid for the DM actuator spacing
- display_actuators(annotate=False, grid=True, what='opd', crosshairs=False, *args, **kwargs)[source]
Display the optical surface, viewed as discrete actuators
- Parameters:
- annotatebool
Annotate coordinates and types of actuators on the display? Default false.
- gridbool
Annotate grid of actuators on the display? Default true.
- whatstring
What to display: ‘intensity’ transmission, ‘opd’, or ‘both’
- crosshairsbool
Draw crosshairs on plot to indicate the origin.
- get_act_coordinates(one_d=False, include_transformations=False)[source]
Y and X coordinates for the actuators
- Parameters:
- one_dbool
Return 1-dimensional arrays of coordinates per axis? Default is to return 2D arrays with same shape as full array.
- include_transformationsbool
Return apparent coordinates after rotations, etc of the DM.
- Returns:
- y_act, x_actfloat ndarrays
actuator coordinates, in units of meters
- get_opd(wave)[source]
Return the surface shape OPD for the optic. Interpolates from the current optic surface state onto the desired coordinates for the wave.
- get_transmission(wave)[source]
Note that this is the amplitude transmission, not the total intensity transmission.
- set_actuator(actx, acty, new_value)[source]
Set an individual actuator of the DM. Parameters ————- actx, acty : integers
Coordinates of the actuator you wish to control
- new_valuefloat
Desired surface height for that actuator, in meters by default or use astropy Units to specify another unit if desired.
Example
dm.set_actuator(12, 22, 123.4*u.nm)