def test_nonempty_matrix_vector(self): n = 10 m = trackcpp.Matrix() for i in range(n): self.ml.append(m) self.assertEqual(len(self.ml), n)
def append(self, value): if isinstance(value, _trackcpp.Matrix): self._ml.append(value) elif isinstance(value, _numpy.ndarray): m = _trackcpp.Matrix() for line in value: line_as_floats = [float(x) for x in line] m.append(line_as_floats) self._ml.append(m) elif self._is_list_of_lists(value): m = _trackcpp.Matrix() for line in value: m.append(line) self._ml.append(m) else: raise TrackingException('can only append matrix-like objects')
def test_append_trackcpp_matrix(self): m = trackcpp.Matrix() for line in self.matrix: m.append(line) self.ml.append(m) comparison = self.ml[0] == numpy.array(self.matrix) self.assertTrue(numpy.all(comparison))
def test_matrix_list_init_arg(self): m_cpp = trackcpp.Matrix() for line in self.m: m_cpp.append(line) ml_cpp = trackcpp.CppMatrixVector() ml_cpp.append(m_cpp) ml = pyaccel.tracking.MatrixList(ml_cpp) self.assertEqual(len(ml), 1) comparison = ml[0] == numpy.array(self.m) self.assertTrue(numpy.all(comparison))
def find_m44(accelerator, indices=None, energy_offset=0.0, closed_orbit=None): """Calculate 4D transfer matrices of elements in an accelerator. Keyword arguments: accelerator indices energy_offset closed_orbit Return values: m44 cumul_trans_matrices -- values at the start of each lattice element """ if indices is None: indices = list(range(len(accelerator))) if closed_orbit is None: # calcs closed orbit if it was not passed. fixed_point_guess = _trackcpp.CppDoublePos() fixed_point_guess.de = energy_offset _closed_orbit = _trackcpp.CppDoublePosVector() r = _trackcpp.track_findorbit4(accelerator._accelerator, _closed_orbit, fixed_point_guess) if r > 0: raise TrackingException(_trackcpp.string_error_messages[r]) else: _closed_orbit = _4Numpy2CppDoublePosVector(closed_orbit, de=energy_offset) _cumul_trans_matrices = _trackcpp.CppMatrixVector() _m44 = _trackcpp.Matrix() _v0 = _trackcpp.CppDoublePos() r = _trackcpp.track_findm66(accelerator._accelerator, _closed_orbit, _cumul_trans_matrices, _m44, _v0) if r > 0: raise TrackingException(_trackcpp.string_error_messages[r]) m44 = _CppMatrix24Numpy(_m44) if indices == 'm44': return m44 cumul_trans_matrices = [] for i in range(len(_cumul_trans_matrices)): cumul_trans_matrices.append(_CppMatrix24Numpy( _cumul_trans_matrices[i])) return m44, cumul_trans_matrices
def find_m66(accelerator, indices=None, closed_orbit=None): """Calculate 6D transfer matrices of elements in an accelerator. Keyword arguments: accelerator indices closed_orbit Return values: m66 cumul_trans_matrices -- values at the start of each lattice element """ if indices is None: indices = list(range(len(accelerator))) if closed_orbit is None: # Closed orbit is calculated by trackcpp _closed_orbit = _trackcpp.CppDoublePosVector() else: _closed_orbit = _Numpy2CppDoublePosVector(closed_orbit) _cumul_trans_matrices = _trackcpp.CppMatrixVector() _m66 = _trackcpp.Matrix() _v0 = _trackcpp.CppDoublePos() r = _trackcpp.track_findm66(accelerator._accelerator, _closed_orbit, _cumul_trans_matrices, _m66, _v0) if r > 0: raise TrackingException(_trackcpp.string_error_messages[r]) m66 = _CppMatrix2Numpy(_m66) if indices == 'm66': return m66 cumul_trans_matrices = MatrixList(_cumul_trans_matrices) return m66, cumul_trans_matrices
def calc_twiss(accelerator=None, init_twiss=None, fixed_point=None, indices='open', energy_offset=None): """Return Twiss parameters of uncoupled dynamics. Keyword arguments: accelerator -- Accelerator object init_twiss -- Twiss parameters at the start of first element fixed_point -- 6D position at the start of first element indices -- Open or closed energy_offset -- float denoting the energy deviation (used only for periodic solutions). Returns: tw -- list of Twiss objects (closed orbit data is in the objects vector) m66 -- one-turn transfer matrix """ if indices == 'open': closed_flag = False elif indices == 'closed': closed_flag = True else: raise OpticsException("invalid value for 'indices' in calc_twiss") _m66 = _trackcpp.Matrix() _twiss = _trackcpp.CppTwissVector() if init_twiss is not None: ''' as a transport line: uses init_twiss ''' _init_twiss = init_twiss._t if fixed_point is None: _fixed_point = _init_twiss.co else: raise OpticsException( 'arguments init_twiss and fixed_point are mutually exclusive') r = _trackcpp.calc_twiss(accelerator._accelerator, _fixed_point, _m66, _twiss, _init_twiss, closed_flag) else: ''' as a periodic system: try to find periodic solution ''' if accelerator.harmonic_number == 0: raise OpticsException( 'Either harmonic number was not set or calc_twiss was' 'invoked for transport line without initial twiss') if fixed_point is None: _closed_orbit = _trackcpp.CppDoublePosVector() _fixed_point_guess = _trackcpp.CppDoublePos() if energy_offset is not None: _fixed_point_guess.de = energy_offset if not accelerator.cavity_on and not accelerator.radiation_on: r = _trackcpp.track_findorbit4(accelerator._accelerator, _closed_orbit, _fixed_point_guess) elif not accelerator.cavity_on and accelerator.radiation_on: raise OpticsException( 'The radiation is on but the cavity is off') else: r = _trackcpp.track_findorbit6(accelerator._accelerator, _closed_orbit, _fixed_point_guess) if r > 0: raise _tracking.TrackingException( _trackcpp.string_error_messages[r]) _fixed_point = _closed_orbit[0] else: _fixed_point = _tracking._Numpy2CppDoublePos(fixed_point) if energy_offset is not None: _fixed_point.de = energy_offset r = _trackcpp.calc_twiss(accelerator._accelerator, _fixed_point, _m66, _twiss) if r > 0: raise OpticsException(_trackcpp.string_error_messages[r]) twiss = TwissList(_twiss) m66 = _tracking._CppMatrix2Numpy(_m66) return twiss, m66