COMPASS  5.0.0
End-to-end AO simulation tool using GPU acceleration
shesha.supervisor.optimizers.modalBasis.ModalBasis Class Reference

This optimizer class handles all the modal basis and DM Influence functions related operations. More...

Inheritance diagram for shesha.supervisor.optimizers.modalBasis.ModalBasis:
Collaboration diagram for shesha.supervisor.optimizers.modalBasis.ModalBasis:

Public Member Functions

def __init__ (self, config, dms, target)
 
csr_matrix compute_influ_basis (self, int dm_index)
 Computes and return the influence function phase basis of the specified DM as a sparse matrix. More...
 
np.ndarray compute_modes_to_volts_basis (self, str modal_basis_type, *bool merged=False, int nbpairs=None, bool return_delta=False)
 Computes a given modal basis ("KL2V", "Btt", "Btt_petal") and return the 2 transfer matrices. More...
 
np.ndarray compute_btt_basis (self, *bool merged=False, int nbpairs=None, bool return_delta=False)
 Computes the so-called Btt modal basis. More...
 
np.ndarray compute_merged_influ (self, int dm_index, *int nbpairs=None)
 Used to compute merged IF from each side of the spider for an ELT case (Petalling Effect) More...
 
np.ndarray compute_btt_petal (self)
 Computes a Btt modal basis with Pistons filtered. More...
 
np.ndarray compute_phase_to_modes (self, np.ndarray modal_basis)
 Return the phase to modes matrix by using the given modal basis. More...
 

Public Attributes

 config
 
 dms
 
 target
 
 slaved_actus
 
 selected_actus
 
 couples_actus
 
 index_under_spiders
 
 modal_basis
 
 projection_matrix
 Instantiate a ModalBasis object. More...
 
 btt
 

Detailed Description

This optimizer class handles all the modal basis and DM Influence functions related operations.

Definition at line 50 of file modalBasis.py.

Constructor & Destructor Documentation

◆ __init__()

def shesha.supervisor.optimizers.modalBasis.ModalBasis.__init__ (   self,
  config,
  dms,
  target 
)

Definition at line 103 of file modalBasis.py.

103  nbpairs: int = None, return_delta: bool = False) -> np.ndarray:
104  """ Computes a given modal basis ("KL2V", "Btt", "Btt_petal") and return the 2 transfer matrices
105 
106  Parameters:
107  modal_basis_type : (str) : modal basis to compute ("KL2V", "Btt", "Btt_petal")
108 
109  merged : (bool, optional) :
110 
111  nbpairs : (int, optional) :
112 
113  Return:

Member Function Documentation

◆ compute_btt_basis()

np.ndarray shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_btt_basis (   self,
*bool   merged = False,
int   nbpairs = None,
bool   return_delta = False 
)

Computes the so-called Btt modal basis.

The <merged> flag allows merto merge 2x2 the actuators influence functions for actuators on each side of the spider (ELT case)

Parameters
merged (bool, optional) : If True, merge 2x2 the actuators influence functions for actuators on each side of the spider (ELT case). Default is False
nbpairs (int, optional) : Default is None. TODO : description
return_delta (bool, optional) : If False (default), the function returns Btt (modes to volts matrix), and P (volts to mode matrix). If True, returns delta = IF.T.dot(IF) / N instead of P
Return
Btt (np.ndarray) : Btt modes to volts matrix
projection_matrix (np.ndarray) : volts to Btt modes matrix

Definition at line 199 of file modalBasis.py.

199 
200  if (merged):
201  btt2 = np.zeros((len(boolarray) + 2, self.btt.shape[1]))
202  btt2[np.r_[~boolarray, True, True], :] = self.btt
203  btt2[couples_actus[:, 1], :] = btt2[couples_actus[:, 0], :]
204 
205  P2 = np.zeros((self.btt.shape[1], len(boolarray) + 2))
206  P2[:, np.r_[~boolarray, True, True]] = self.projection_matrix
207  P2[:, couples_actus[:, 1]] = P2[:, couples_actus[:, 0]]
208  self.btt = btt2
209  self.projection_matrix = P2
210 
211  return self.btt, self.projection_matrix
212 
213  def compute_merged_influ(self, dm_index : int, *, nbpairs: int = None) -> np.ndarray:
214  """ Used to compute merged IF from each side of the spider
215  for an ELT case (Petalling Effect)
216 
217  Parameters:
218  dm_index : (int) : DM index
219 
220  nbpairs : (int, optional) : Default is None. TODO : description
221 
222  Return:
223  pairs : (np.ndarray) : TODO description
224 
225  discard : (list) : TODO description
226  """
227  p_geom = self.config.p_geom
228 
229 
230  cent = p_geom.pupdiam / 2. + 0.5
231  p_tel = self.config.p_tel
232  p_tel.t_spiders = 0.51
233  spup = mkP.make_pupil(p_geom.pupdiam, p_geom.pupdiam, p_tel, cent,
234  cent).astype(np.float32).T
235 
236  p_tel.t_spiders = 0.
237  spup2 = mkP.make_pupil(p_geom.pupdiam, p_geom.pupdiam, p_tel, cent,
238  cent).astype(np.float32).T
Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_btt_petal()

np.ndarray shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_btt_petal (   self)

Computes a Btt modal basis with Pistons filtered.

Return
Btt (np.ndarray) : Btt modes to volts matrix
P (np.ndarray) : volts to Btt modes matrix

Definition at line 359 of file modalBasis.py.

359  for i in range(nbmode):
360  self.dms.set_command((modal_basis[:, i]).copy())
361  # self.next(see_atmos=False)
362  self.target.raytrace(0, dms=self.dms, ncpa=False, reset=True)
363  phase = self.target.get_tar_phase(0, pupil=True)
364  # Normalisation pour les unites rms en microns !!!
365  norm = np.sqrt(np.sum((phase)**2) / S)
366  phase_to_modes[i] = phase / norm
367  return phase_to_modes
Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_influ_basis()

csr_matrix shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_influ_basis (   self,
int  dm_index 
)

Computes and return the influence function phase basis of the specified DM as a sparse matrix.

Parameters
dm_index (int) : Index of the DM
Return
influ_sparse (csr_matrix) : influence function phases

Definition at line 123 of file modalBasis.py.

123  self.config.p_atmos, self.config.p_tel)
124  fnz = util.first_non_zero(self.modal_basis, axis=0)
125  # Computing the sign of the first non zero element
126  #sig = np.sign(modal_basis[[fnz, np.arange(modal_basis.shape[1])]])
127  sig = np.sign(self.modal_basis[tuple([
Here is the caller graph for this function:

◆ compute_merged_influ()

np.ndarray shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_merged_influ (   self,
int  dm_index,
*int   nbpairs = None 
)

Used to compute merged IF from each side of the spider for an ELT case (Petalling Effect)

Parameters
dm_index (int) : DM index
nbpairs (int, optional) : Default is None. TODO : description
Return
pairs (np.ndarray) : TODO description
discard (list) : TODO description

Definition at line 252 of file modalBasis.py.

252  discard = np.zeros(len(dm_posx), dtype=np.bool)
253  pairs = []
254 
255  # For each of the k pieces of the spider
256  for k, px_list in enumerate(px_list_spider):
257  pts = np.c_[px_list[1],
258  px_list[0]] # x,y coord of pixels of the spider piece
259  # line_eq = [a, b]
260  # Which minimizes leqst squares of aa*x + bb*y = 1
261  line_eq = np.linalg.pinv(pts).dot(np.ones(pts.shape[0]))
262  aa, bb = line_eq[0], line_eq[1]
263 
264  # Find any point of the fitted line.
265  # For simplicity, the intercept with one of the axes x = 0 / y = 0
266  if np.abs(bb) < np.abs(aa): # near vertical
267  one_point = np.array([1 / aa, 0.])
268  else: # otherwise
269  one_point = np.array([0., 1 / bb])
270 
271  # Rotation that aligns the spider piece to the horizontal
272  rotation = np.array([[-bb, aa], [-aa, -bb]]) / (aa**2 + bb**2)**.5
273 
274  # Rotated the spider mask
275  rotated_px = rotation.dot(pts.T - one_point[:, None])
276  # Min and max coordinates along the spider length - to filter actuators that are on
277  # 'This' side of the pupil and not the other side
278  min_u, max_u = rotated_px[0].min() - 5. * pitch, rotated_px[0].max(
279  ) + 5. * pitch
280 
281  # Rotate the actuators
282  rotated_actus = rotation.dot(dm_pos_mat - one_point[:, None])
283  sel_good_side = (rotated_actus[0] > min_u) & (rotated_actus[0] < max_u)
284  threshold = 0.05
285  # Actuators below this piece of spider
286  sel_discard = (np.abs(rotated_actus[1]) < threshold * pitch) & sel_good_side
287  discard |= sel_discard
288 
289  # Actuator 'near' this piece of spider
290  sel_pairable = (np.abs(rotated_actus[1]) > threshold * pitch) & \
291  (np.abs(rotated_actus[1]) < 1. * pitch) & \
292  sel_good_side
293 
294  pairable_index = np.where(sel_pairable)[0] # Indices of these actuators
295  u_coord = rotated_actus[
296  0, sel_pairable] # Their linear coord along the spider major axis
297 
298  order = np.sort(u_coord) # Sort by linear coordinate
299  order_index = pairable_index[np.argsort(
300  u_coord)] # And keep track of original indexes
301 
302  # i = 0
303  # while i < len(order) - 1:
304  if (nbpairs is None):
305  i = 0
306  ii = len(order) - 1
307  else:
308  i = len(order) // 2 - nbpairs
309  ii = len(order) // 2 + nbpairs
310  while (i < ii):
311  # Check if next actu in sorted order is very close
312  # Some lonely actuators may be hanging in this list
313  if np.abs(order[i] - order[i + 1]) < .2 * pitch:
314  pairs += [(order_index[i], order_indx[i + 1])]
315  i += 2
316  else:
317  i += 1
318  print('To discard: %u actu' % np.sum(discard))
319  print('%u pairs to slave' % len(pairs))
320  if np.sum(discard) == 0:
321  discard = []
322  else:
323  list(np.where(discard)[0])
324  return np.asarray(pairs), list(np.where(discard)[0])
325 
326  def compute_btt_petal(self) -> np.ndarray:
327  """ Computes a Btt modal basis with Pistons filtered
328 
329  Return:
330  Btt : (np.ndarray) : Btt modes to volts matrix
331 
332  P : (np.ndarray) : volts to Btt modes matrix
333  """
334  pzt_index = np.where([d.type is scons.DmType.PZT for d in self.config.p_dms])[0][0]
335  influ_pzt = self.compute_influ_basis(pzt_index)
336  petal_dm_index = np.where([
337  d.influ_type is scons.InfluType.PETAL for d in self.config.p_dms
338  ])[0][0]
339  influ_petal = self.compute_influ_basis(petal_dm_index)
340  tt_index = np.where([d.type is scons.DmType.TT for d in self.config.p_dms])[0][0]
341  influ_tt = self.compute_influ_basis(tt_index).toarray()
342 
343  self.modal_basis, self.projection_matrix = basis.compute_btt(influ_pzt.T, influ_tt.T, influ_petal=influ_petal)
344  return self.modal_basis, self.projection_matrix
345 
346  def compute_phase_to_modes(self, modal_basis: np.ndarray) -> np.ndarray:
347  """ Return the phase to modes matrix by using the given modal basis
348 
349  Parameters:
350  modal_basis : (np.ndarray) : Modal basis matrix
351 
Here is the caller graph for this function:

◆ compute_modes_to_volts_basis()

np.ndarray shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_modes_to_volts_basis (   self,
str  modal_basis_type,
*bool   merged = False,
int   nbpairs = None,
bool   return_delta = False 
)

Computes a given modal basis ("KL2V", "Btt", "Btt_petal") and return the 2 transfer matrices.

Parameters
modal_basis_type (str) : modal basis to compute ("KL2V", "Btt", "Btt_petal")
merged (bool, optional) :
nbpairs (int, optional) :
Return
modal_basis (np.ndarray) : modes to volts matrix
projection_matrix (np.ndarray) : volts to modes matrix (None if "KL")

Definition at line 142 of file modalBasis.py.

142  ])]) # pour remove le future warning!
143  self.modal_basis *= sig[None, :]
144  elif (modal_basis_type == "Btt_petal"):
145  print("Computing Btt with a Petal basis...")
146  self.modal_basis, self.projection_matrix = self.compute_btt_petal()
147  else:
148  raise ArgumentError("Unsupported modal basis")
149 
150  return self.modal_basis, self.projection_matrix
151 
152  def compute_btt_basis(self, *, merged: bool = False, nbpairs: int = None,
153  return_delta: bool = False) -> np.ndarray:
154  """ Computes the so-called Btt modal basis. The <merged> flag allows merto merge
155  2x2 the actuators influence functions for actuators on each side of the spider (ELT case)
156 
157  Parameters:
158  merged : (bool, optional) : If True, merge 2x2 the actuators influence functions for
159  actuators on each side of the spider (ELT case). Default
160  is False
161 
162  nbpairs : (int, optional) : Default is None. TODO : description
163 
164  return_delta : (bool, optional) : If False (default), the function returns
165  Btt (modes to volts matrix),
166  and P (volts to mode matrix).
167  If True, returns delta = IF.T.dot(IF) / N
168  instead of P
169 
170  Return:
171  Btt : (np.ndarray) : Btt modes to volts matrix
172 
173  projection_matrix : (np.ndarray) : volts to Btt modes matrix
174  """
175  dms_basis = basis.compute_IFsparse(self.dms.dms, self.config.p_dms, self.config.p_geom)
176  influ_basis = dms_basis[:-2,:]
177  tt_basis = dms_basis[-2:,:].toarray()
Here is the call graph for this function:

◆ compute_phase_to_modes()

np.ndarray shesha.supervisor.optimizers.modalBasis.ModalBasis.compute_phase_to_modes (   self,
np.ndarray  modal_basis 
)

Return the phase to modes matrix by using the given modal basis.

Parameters
modal_basis (np.ndarray) : Modal basis matrix
Return
phase_to_modes (np.ndarray) : phase to modes matrix

Definition at line 380 of file modalBasis.py.

Member Data Documentation

◆ btt

shesha.supervisor.optimizers.modalBasis.ModalBasis.btt

Definition at line 233 of file modalBasis.py.

◆ config

shesha.supervisor.optimizers.modalBasis.ModalBasis.config

(config) : Configuration parameters module

Definition at line 104 of file modalBasis.py.

◆ couples_actus

shesha.supervisor.optimizers.modalBasis.ModalBasis.couples_actus

TODO : docstring

Definition at line 109 of file modalBasis.py.

◆ dms

shesha.supervisor.optimizers.modalBasis.ModalBasis.dms

(DmCompass) : DmCompass instance

Definition at line 105 of file modalBasis.py.

◆ index_under_spiders

shesha.supervisor.optimizers.modalBasis.ModalBasis.index_under_spiders

TODO : docstring

Definition at line 110 of file modalBasis.py.

◆ modal_basis

shesha.supervisor.optimizers.modalBasis.ModalBasis.modal_basis

(np.ndarray) : Last modal basis computed

Definition at line 111 of file modalBasis.py.

◆ projection_matrix

shesha.supervisor.optimizers.modalBasis.ModalBasis.projection_matrix

Instantiate a ModalBasis object.

(np.ndarray) : Last projection_matrix computed

Parameters
config (config) : Configuration parameters module
dms (DmCompass) : DmCompass instance
target (TargetCompass) : TargetCompass instance

Definition at line 112 of file modalBasis.py.

◆ selected_actus

shesha.supervisor.optimizers.modalBasis.ModalBasis.selected_actus

TODO : docstring

Definition at line 108 of file modalBasis.py.

◆ slaved_actus

shesha.supervisor.optimizers.modalBasis.ModalBasis.slaved_actus

TODO : docstring

Definition at line 107 of file modalBasis.py.

◆ target

shesha.supervisor.optimizers.modalBasis.ModalBasis.target

(TargetCompass) : TargetCompass instance

Definition at line 106 of file modalBasis.py.


The documentation for this class was generated from the following file: