def Subaru_frontend(empty_lamda, grid_size, PASSVALUE): """ propagates instantaneous complex E-field thru Subaru from the primary through the AO188 AO system in loop over wavelength range this function is called a 'prescription' by proper uses PyPROPER3 to generate the complex E-field at the source, then propagates it through atmosphere, then telescope, to the focal plane the AO simulator happens here this does not include the observation of the wavefront by the detector :returns spectral cube at instantaneous time in the focal_plane() """ # print("Propagating Broadband Wavefront Through Subaru") # Initialize the Wavefront in Proper wfo = opx.Wavefronts() wfo.initialize_proper() # Atmosphere # atmos has only effect on phase delay, not intensity wfo.loop_collection(atmos.add_atmos, PASSVALUE['iter'], plane_name='atmosphere') # Defines aperture (baffle-before primary) # Obscurations (Secondary and Spiders) wfo.loop_collection(opx.add_obscurations, d_primary=tp.d_nsmyth, d_secondary=tp.d_secondary, legs_frac=0.05) wfo.loop_collection(proper.prop_circular_aperture, **{'radius': tp.entrance_d / 2}) # clear inside, dark outside wfo.loop_collection( proper.prop_define_entrance, plane_name='entrance_pupil') # normalizes abs intensity if ap.companion: # Must do this after all calls to prop_define_entrance wfo.loop_collection(opx.offset_companion) wfo.loop_collection(proper.prop_circular_aperture, **{'radius': tp.entrance_d / 2 }) # clear inside, dark outside # Test Sampling if sp.verbose: opx.check_sampling(PASSVALUE['iter'], wfo, "initial", getframeinfo(stack()[0][0]), units='mm') # Testing Primary Focus (instead of propagating to focal plane) # wfo.loop_collection(opx.prop_pass_lens, tp.flen_nsmyth, tp.flen_nsmyth) # test only going to prime focus ######################################## # Subaru Propagation ####################################### # Effective Primary # CPA from Effective Primary wfo.loop_collection(aber.add_aber, step=PASSVALUE['iter'], lens_name='ao188-OAP1') # Zernike Aberrations- Low Order # wfo.loop_collection(aber.add_zern_ab, tp.zernike_orders, aber.randomize_zern_values(tp.zernike_orders)) wfo.loop_collection(opx.prop_pass_lens, tp.flen_nsmyth, tp.dist_nsmyth_ao1) ######################################## # AO188 Propagation ######################################## # # AO188-OAP1 wfo.loop_collection(aber.add_aber, step=PASSVALUE['iter'], lens_name='ao188-OAP1') wfo.loop_collection(opx.prop_pass_lens, tp.fl_ao1, tp.dist_ao1_dm) # AO System if tp.use_ao: WFS_map = ao.open_loop_wfs(wfo) wfo.loop_collection( ao.deformable_mirror, WFS_map, PASSVALUE['iter'], plane_name='woofer', debug=sp.debug ) # don't use PASSVALUE['WFS_map'] here because open loop # ------------------------------------------------ wfo.loop_collection(proper.prop_propagate, tp.dist_dm_ao2) # AO188-OAP2 wfo.loop_collection(aber.add_aber, step=PASSVALUE['iter'], lens_name='ao188-OAP2') # wfo.loop_collection(aber.add_zern_ab, tp.zernike_orders, aber.randomize_zern_values(tp.zernike_orders)/2) wfo.loop_collection(opx.prop_pass_lens, tp.fl_ao2, tp.dist_oap2_focus) ######################################## # Focal Plane # ####################################### # Check Sampling in focal plane if sp.verbose: wfo.loop_collection(opx.check_sampling, PASSVALUE['iter'], "focal plane", getframeinfo(stack()[0][0]), units='nm') # wfo.focal_plane fft-shifts wfo from Fourier Space (origin==lower left corner) to object space (origin==center) cpx_planes, sampling = wfo.focal_plane() print(f"Finished datacube at timestep = {PASSVALUE['iter']}") return cpx_planes, sampling
def deformable_mirror(wf, WFS_map, iter, previous_output=None, apodize=False, plane_name='', debug=False): """ combine different DM actuator commands into single map to send to prop_dm prop_dm needs an input map of n_actuators x n_actuators in units of actuator command height. quick_ao will handle the conversion to actuator command height, and the CDI probe must be scaled in cdi.probe_amp in params in units of m. Each subroutine is also responsible for creating a map of n_actuators x n_actuators spacing. prop_dm handles the resampling of this map onto the wavefront, including the influence function. Its some wizardry that happens in c, and presumably it is taken care of so you don't have to worry about it. In the call to proper.prop_dm, we apply the flag tp.fit_dm, which switches between two 'modes' of proper's DM surface fitting. If FALSE, the DM is driven to the heights specified by dm_map, and the influence function will act on these heights to define the final surface shape applied to the DM, which may differ substantially from the initial heights specified by dm_map. If TRUE, proper will iterate applying the influence function to the input heights, and adjust the heights until the difference between the influenced-map and input map meets some proper-defined convergence criterea. Setting tp.fit_dm=TRUE will obviously slow down the code, but will (likely) more accurately represent a well-calibrated DM response function. much of this code copied over from example from Proper manual on pg 94 :param wf: single wavefront :param WFS_map: wavefront sensor map, should be in units of phase delay :param previous_output: :param iter: the current index of iteration (which timestep this is) :param plane_name: name of plane (should be 'woofer' or 'tweeter' for best functionality) :return: nothing is returned, but the probe map has been applied to the DM via proper.prop_dm. DM plane post DM application can be saved via the sp.save_list functionality """ assert np.logical_xor(WFS_map is None, previous_output is None) # AO Actuator Count from DM Type if plane_name == 'tweeter' and hasattr(tp, 'act_tweeter'): nact = tp.act_tweeter elif plane_name == 'woofer' and hasattr(tp, 'act_woofer'): nact = tp.act_woofer else: nact = tp.ao_act # DM Coordinates nact_across_pupil = nact - 2 # number of full DM actuators across pupil (oversizing DM extent) dm_xc = ( nact / 2 ) # The location of the optical axis (center of the wavefront) on the DM in dm_yc = ( nact / 2 ) # actuator units. First actuator is centered on (0.0, 0.0). The 0.5 is a # parameter introduced/tuned by Rupert to remove weird errors (address this). # KD verified this needs to be here or else suffer weird errors 9/19 # TODO address/remove the 0.5 in DM x,y coordinates ############################ # Creating DM Surface Map ############################ d_beam = 2 * proper.prop_get_beamradius(wf) # beam diameter act_spacing = d_beam / nact_across_pupil # actuator spacing [m] ####### # AO ####### if previous_output is not None and WFS_map is None: dm_map = update_dm(previous_output) else: dm_map = quick_ao(wf, nact, WFS_map[wf.iw]) ######### # Waffle ######### if tp.satelite_speck['apply'] and plane_name is not 'woofer': waffle = make_speckle_kxy(tp.satelite_speck['xloc'], tp.satelite_speck['yloc'], tp.satelite_speck['amp'], tp.satelite_speck['phase']) waffle += make_speckle_kxy(tp.satelite_speck['xloc'], -tp.satelite_speck['yloc'], tp.satelite_speck['amp'], tp.satelite_speck['phase']) dm_map += waffle ####### # CDI ###### if cdi.use_cdi and plane_name == cdi.which_DM: theta = cdi.phase_series[iter] if not np.isnan(theta): # dprint(f"Applying CDI probe, lambda = {wfo.wsamples[iw]*1e9:.2f} nm") cdi.save_tseries(iter, datetime.datetime.now()) probe = config_probe(theta, nact, iw=wf.iw, ib=wf.ib, tstep=iter) dm_map = dm_map + probe # Add Probe to DM map ######################### # Applying Piston Error ######################### if tp.piston_error: mean_dm_map = np.mean(np.abs(dm_map)) var = 1e-4 # 1e-11 dm_map = dm_map + np.random.normal(0, var, (dm_map.shape[0], dm_map.shape[1])) ######################### # proper.prop_dm ######################### dmap = proper.prop_dm(wf, dm_map, dm_xc, dm_yc, act_spacing, FIT=tp.fit_dm) # if debug and wf.iw == 0 and wf.ib == 0 and iter == 0: dprint(plane_name) check_sampling(wf, iter, plane_name + ' DM pupil plane', getframeinfo(stack()[0][0]), units='mm') quick2D(WFS_map[wf.iw], title=f"WFS map after masking", zlabel='unwrapped phase (rad)', vlim=[-3 * np.pi, 3 * np.pi]) fig, ax = plt.subplots(1, 1) cax = ax.imshow(dm_map * 1e9, interpolation='none', origin='lower') plt.title(f'{plane_name} dm_map (actuator coordinates)') cb = plt.colorbar(cax) cb.set_label('nm') plt.show() post_ao = unwrap_phase( proper.prop_get_phase(wf)) * wf.lamda / (2 * np.pi) # quick2D(pre_ao_dist*1e9, title='unwrapped wavefront before DM', zlabel='nm', show=False) # , vlim=(-0.5e-7,0.5e-7)) # quick2D(np.abs(pre_ao_amp)**2, title='Pre-AO Intensity', show=False)#, vlim=(-0.5e-7,0.5e-7)) # quick2D(dmap, title='the phase map prop_dm is applying', zlabel='distance (m)', show=False)#, vlim=(-0.5e-7,0.5e-7)) # plt.figure() # plt.plot(pre_ao_dist[len(pre_ao_dist)//2], label=f'pre_ao 1D cut, row {len(pre_ao_dist)//2}') # plt.plot(2*dmap[len(dmap)//2], label=f'dmap 1D cut (x2), row {len(dmap)//2}') # plt.plot((pre_ao_dist + (2*dmap))[len(dmap)//2], label='difference') # plt.legend() # plt.xlim(sp.grid_size//2*np.array([1-sp.beam_ratio*1.1, 1+sp.beam_ratio*1.1])) # quick2D(pre_ao + (2*dmap), title='diff', zlabel='m', show=False, vlim=(-0.5e-7,0.5e-7)) # quick2D(post_ao, title='unwrapped wavefront after DM', zlabel='m', show=True, vlim=(-0.5e-7,0.5e-7)) # quick2D(np.abs(proper.prop_get_amplitude(wf))**2, title='wavefront after DM intensity', show=False) # quick2D(proper.prop_get_phase(wf), title='wavefront after DM in phase units', zlabel='Phase', # show=True) # colormap='sunlight', if apodize: hardmask_pupil(wf) return dmap