├── README └── ramsey_experiment_suite.py /README: -------------------------------------------------------------------------------- 1 | A list of functions for master equation quantum optics simulations based on python QuTiP. 2 | 3 | Simulations include: carrier, red, and blue sideband Rabi flops, free evolution, and Ramsey experiments. 4 | 5 | Can initiallize the ion state to pure states, as well as mixed states, e.g. thernmal, or displaced thermal states. All operations on the ion are a result of a numerical solution of the master equation for the Hamiltonian 6 | evolution + the Lindblad master equation operators. 7 | 8 | This module uses the QuTiP library of functions. See QuTip installation instructions for dependencies. 9 | -------------------------------------------------------------------------------- /ramsey_experiment_suite.py: -------------------------------------------------------------------------------- 1 | from qutip import tensor, qeye, basis, destroy, thermal_dm, mesolve, tomography, hinton, displace, squeez 2 | from pylab import * 3 | import numpy as np 4 | from numpy import pi, sin, cos, exp, sqrt, abs, sum, shape, arange as rang 5 | from scipy.special.orthogonal import eval_genlaguerre as laguerre 6 | from scipy.sparse import csr_matrix, dia_matrix 7 | from scipy.optimize import minimize_scalar # needs scipy 0.9.0 or newer 8 | from scipy.optimize import fminbound 9 | 10 | def collapse_operators(N, n_th_a, gamma_motion, gamma_motion_phi, gamma_atom): 11 | '''Collapse operators for the master equation of a single atom and a harmonic oscillator 12 | @ var N: size of the harmonic oscillator Hilbert space 13 | @ var n_th: temperature of the noise bath in quanta 14 | @ var gamma_motion: heating rate of the motion 15 | @ var gamma_motion_phi: dephasing rate of the motion 16 | @ var gamma_atom: decay rate of the atom 17 | 18 | returns: list of collapse operators for master equation solution of atom + harmonic oscillator 19 | ''' 20 | a = tensor(destroy(N), qeye(2)) 21 | sm = tensor(qeye(N), destroy(2)) 22 | c_op_list = [] 23 | 24 | rate = gamma_motion * (1 + n_th_a) 25 | if rate > 0.0: 26 | c_op_list.append(sqrt(rate) * a) 27 | 28 | rate = gamma_motion * n_th_a 29 | if rate > 0.0: 30 | c_op_list.append(sqrt(rate) * a.dag()) 31 | 32 | rate = gamma_motion_phi 33 | if rate > 0.0: 34 | c_op_list.append(sqrt(rate) * a.dag() * a) 35 | 36 | rate = gamma_atom 37 | if rate > 0.0: 38 | c_op_list.append(sqrt(rate) * sm) 39 | return c_op_list 40 | 41 | def population_operators(N): 42 | '''Population operators for the master equation 43 | @ var N: size of the oscillator Hilbert space 44 | 45 | returns: list of population operators for the harmonic oscillator and the atom 46 | ''' 47 | a = tensor(destroy(N), qeye(2)) 48 | sm = tensor( qeye(N), destroy(2)) 49 | return [a.dag()*a, sm.dag()*sm] 50 | 51 | 52 | def lamb_dicke(wavelength, frequency, lasertheta): 53 | '''computes the Lamb Dicke parameter 54 | @ var theta: laser projection angle in degrees 55 | @ var wavelength: laser wavelength in meters 56 | @ var frequency: trap frequency in Hz 57 | 58 | returns: Lamb-Dicke parameter 59 | ''' 60 | k = 2.*pi/wavelength 61 | return k*sqrt(hbar/(2*mass*2*pi*frequency))*np.abs(cos(pi*lasertheta/180)) 62 | 63 | def rabi_coupling(n, m, eta): 64 | '''return the Rabi matrix elements for n = 0...n-1 65 | @ var n: size of the Hilbert space 66 | @ var m: order m of the sideband 67 | @ var eta: Lamd-Dicke parameter 68 | 69 | returns: Wnm Rabi matrix elements for n = 0,...n-1 70 | ''' 71 | L = np.array([laguerre(ii, abs(m), eta**2) for ii in range(n)]) 72 | fctr = np.array([factorial_simplified(ii, abs(m))**(0.5) for ii in range(n)]) 73 | Wmn = exp(-1./2.*eta**2) * eta**abs(m) * fctr * L 74 | if m<0: 75 | Wmn = np.concatenate((np.zeros(-m), Wmn[:m])) 76 | return abs(Wmn) 77 | 78 | def factorial_simplified(n, m): 79 | '''return the factor n!/(n+m)! with n>0, m = 0,\pm 1, ...''' 80 | if not isinstance(m , (long, int)): 81 | print 'factorial_simplified only accepts integer arguments' 82 | return NaN 83 | if n+m<0: 84 | print 'n+m cannot be negative. Returning 0.' 85 | return 0. 86 | if m<0: 87 | return (n+m+1)*factorial_simplified(n,m+1) 88 | elif m>0: 89 | return (1./(n+m))*factorial_simplified(n,m-1) 90 | else: 91 | return 1. 92 | 93 | def carrier_flop(rho0, W, eta, delta, theta, phi, c_op_list = [], return_op_list = []): 94 | ''' Return values of atom and motion populations during carrier Rabi flop 95 | for rotation angles theta. Calls numerical solution of master equation. 96 | @ var rho0: initial density matrix 97 | @ var W: bare Rabi frequency 98 | @ var eta: Lamb-Dicke parameter 99 | @ var delta: detuning between atom and motion 100 | @ var theta: list of Rabi rotation angles (i.e. theta, or g*time) 101 | @ var phi: phase of the input laser pulse 102 | @ var c_op_list: list of collapse operators for the master equation treatment 103 | @ var return_op_list: list of population operators the values of which will be returned 104 | 105 | returns: time, populations of motional mode and atom 106 | ''' 107 | N = shape(rho0.data)[0]/2 # assume N Fock states and two atom states 108 | a = tensor(destroy(N), qeye(2)) 109 | Wc = qeye(N) 110 | Wc.data = csr_matrix( qeye(N).data.dot( np.diag(rabi_coupling(N,0,eta) ) ) ) 111 | sm = tensor( Wc, destroy(2)) 112 | 113 | # use the rotating wave approxiation 114 | H = delta * a.dag() * a + \ 115 | (1./2.)* W * (sm.dag()*exp(1j*phi) + sm*exp(-1j*phi)) 116 | 117 | if hasattr(theta, '__len__'): 118 | if len(theta)>1: # I need to be able to pass a list of length zero and not get an error 119 | time = theta/W 120 | else: 121 | time = theta/W 122 | 123 | output = mesolve(H, rho0, time, c_op_list, return_op_list) 124 | return time, output 125 | 126 | def bsb_flop(rho0, W, eta, delta, theta, phi, c_op_list = [], return_op_list = []): 127 | ''' Return values of atom and motion populations during blue sideband Rabi flop 128 | for rotation angles theta. Calls numerical solution of master equation for the 129 | anti Jaynes-Cummings Hamiltonian. 130 | @ var rho0: initial density matrix 131 | @ var W: bare Rabi frequency 132 | @ var delta: detuning between atom and motion 133 | @ var theta: list of Rabi rotation angles (i.e. theta, or g*time) 134 | @ var phi: phase of the input laser pulse 135 | @ var c_op_list: list of collapse operators for the master equation treatment 136 | @ var return_op_list: list of population operators the values of which will be returned 137 | 138 | returns: time, populations of motional mode and atom 139 | ''' 140 | N = shape(rho0.data)[0]/2 # assume N Fock states and two atom states 141 | a = tensor(destroy(N), qeye(2)) 142 | sm = tensor( qeye(N), destroy(2)) 143 | Wbsb = destroy(N).dag() 144 | Wbsb.data = csr_matrix( destroy(N).dag().data.dot( np.diag( rabi_coupling(N,1,eta) / np.sqrt(np.linspace(1,N,N)) ) ) ) 145 | Absb = tensor(Wbsb.dag(), qeye(2)) 146 | # use the rotating wave approxiation 147 | # Note that the regular a, a.dag() is used for the time evolution of the oscillator 148 | # Absb is the destruction operator including the state dependent coupling strength 149 | H = delta * a.dag() * a + \ 150 | (1./2.)* W * (Absb.dag() * sm.dag()*exp(1j*phi) + Absb * sm*exp(-1j*phi)) 151 | #print Absb.dag() * sm.dag() + Absb * sm 152 | 153 | if hasattr(theta, '__len__'): 154 | if len(theta)>1: # I need to be able to pass a list of length one and not get an error (needed for Ramsey) 155 | time = theta/(eta*W) 156 | else: 157 | time = theta/(eta*W) 158 | 159 | output = mesolve(H, rho0, time, c_op_list, return_op_list) 160 | return time, output 161 | 162 | def rsb_flop(rho0, W, eta, delta, theta, phi, c_op_list = [], return_op_list = []): 163 | ''' Return values of atom and motion populations during red sideband Rabi flop 164 | for rotation angles theta. Calls numerical solution of master equation for the 165 | Jaynes-Cummings Hamiltonian. 166 | @ var rho0: initial density matrix 167 | @ var W: bare Rabi frequency 168 | @ var delta: detuning between atom and motion 169 | @ var theta: list of Rabi rotation angle (i.e. theta, or g*time) 170 | @ var phi: phase of the input laser pulse 171 | @ var c_op_list: list of collapse operators for the master equation treatment 172 | @ var return_op_list: list of population operators the values of which will be returned 173 | 174 | returns: time, populations of motional mode and atom 175 | ''' 176 | N = shape(rho0.data)[0]/2 # assume N Fock states and two atom states 177 | a = tensor(destroy(N), qeye(2)) 178 | sm = tensor( qeye(N), destroy(2)) 179 | Wrsb = destroy(N) 180 | one_then_zero = ([float(x<1) for x in range(N)]) 181 | Wrsb.data = csr_matrix( destroy(N).data.dot( np.diag( rabi_coupling(N,-1,eta) / np.sqrt(one_then_zero+np.linspace(0,N-1,N)) ) ) ) 182 | Arsb = tensor(Wrsb, qeye(2)) 183 | # use the rotating wave approxiation 184 | # Note that the regular a, a.dag() is used for the time evolution of the oscillator 185 | # Arsb is the destruction operator including the state dependent coupling strength 186 | H = delta * a.dag() * a + \ 187 | (1./2.) * W * (Arsb.dag() * sm * exp(1j*phi) + Arsb * sm.dag() * exp(-1j*phi)) 188 | 189 | if hasattr(theta, '__len__'): 190 | if len(theta)>1: # I need to be able to pass a list of length zero and not get an error 191 | time = theta/(eta*W) 192 | else: 193 | time = theta/(eta*W) 194 | 195 | output = mesolve(H, rho0, time, c_op_list, return_op_list) 196 | return time, output 197 | 198 | def flop_phase(flop_func, is_sideband, rho0, W, eta, delta, c_op_list = [], return_op_list = []): 199 | '''find the rotation angle for maximum atom excitation using the Rabi flopping function defined in flop_func 200 | @ var flop_func: the Rabi flopping function name 201 | @ var is_sideband: flag for whether you are looking at a motional sideband flop (True means sideband) 202 | @ var rho0: density matrix of the input state 203 | @ var W: bare Rabi frequency 204 | @ var delta: detuning between atom and motion 205 | @ var c_op_list: list of collapse operators for the master equation treatment 206 | @ var return_op_list: list of population operators the values of which will be returned 207 | 208 | returns: the 'pi' phase 209 | 210 | Attention: when using with carrier_flop on a pure state without diddipation, you might get 2*n*pi time 211 | ''' 212 | C = int(is_sideband) 213 | f = lambda x: -abs(flop_func(rho0, W, eta, delta, [0, x], 0., c_op_list, return_op_list)[1].expect[1][-1]) 214 | res = minimize_scalar(f, bounds = (0., 2.*pi), method = 'bounded') 215 | return res.x 216 | 217 | def free_evolution(rho0, delta, time, c_op_list = [], return_op_list = []): 218 | '''Free evolution of an atom + oscillator system under a master equation 219 | @ var rho0: initial density matrix 220 | @ var delta: detunning 221 | @ var time: list of free evolution times 222 | @var c_op_list: collapse operators for the master equation treatment 223 | @ var return_op_list: population operators the values of which will be returned 224 | 225 | returns: time, populations of motional mode and atom 226 | ''' 227 | N = shape(rho0.data)[0]/2 # assume N Fock states and two atom states 228 | a = tensor(destroy(N), qeye(2)) 229 | H = delta * a.dag() * a 230 | output = mesolve(H, rho0, time, c_op_list, return_op_list) 231 | return time, output 232 | 233 | def ramsey_experiment(rho0, W, eta, delta, timeslist, c_op_list = [], return_op_list = []): 234 | '''Ramsey experiment where the laser is detuned and I scan the Ramsey time to get phase depencdence 235 | @ var: rho0: initial density matrix 236 | @ var W: bare Rabi frequency 237 | @ var eta: lamb Dicke parameter 238 | @ var delta: detuning 239 | @ var timeslist: list of Ramsey times (not including the pi/2 rotations) 240 | @ var c_op_list: collapse operators for the master equation treatment 241 | @ var return_op_list: population operators the values of which will be returned 242 | 243 | returns: populations of motional mode and atom for different Ramsey times 244 | ''' 245 | Theta = flop_phase(bsb_flop, True, rho0, W, eta, delta, c_op_list, return_op_list)/2 246 | thetas = np.linspace(0,Theta,2) 247 | phi = 0 248 | on_equator = bsb_flop(rho0, W, eta, delta, thetas, phi, c_op_list, []) 249 | rho = on_equator[1].states[-1] 250 | after_wait = free_evolution(rho, delta, timeslist, c_op_list, []) 251 | end_populations = [] 252 | for ii in range(len(timeslist)): 253 | rho = after_wait[1].states[ii] 254 | on_northpole = bsb_flop(rho, W, eta, delta, thetas, phi, c_op_list, return_op_list) 255 | print on_northpole[1].expect[1][1] 256 | end_populations.append(on_northpole[1].expect[1][1]) 257 | return end_populations 258 | 259 | def visualize_dm(rho, atomIncluded = True): 260 | ''' make a hinton diagram for density matrix of an atom + motion (atomIncluded = True), or only motion (atomIncluded = False) 261 | @ var rho: the density matrix 262 | @ var atomIncluded: flag to determine whether the input density matrix has an atom state included or not 263 | ''' 264 | if atomIncluded: 265 | N = shape(rho0.data)[0]/2 # assume N Fock states and two atom states 266 | lbls_list = [[str(d) for d in range(N)], ["u", "d"]] 267 | else: 268 | N = shape(rho.data)[0] # assume N Fock states 269 | lbls_list = [[str(d) for d in range(N)]] 270 | xlabels = [] 271 | for inds in tomography._index_permutations([len(lbls) for lbls in lbls_list]): 272 | xlabels.append("".join([lbls_list[k][inds[k]] for k in range(len(lbls_list))])) 273 | ax = hinton(rho, xlabels=xlabels, ylabels=xlabels) 274 | show() 275 | 276 | def displaced_dm(rho,alpha): 277 | '''apply a displacement operator to a density matrix 278 | @ var rho: the input density matrix for the motion only 279 | @ var alpha: the displacement parameter 280 | 281 | returns rho_displaced: a density matrix Qobj 282 | ''' 283 | N = shape(rho.data)[0] # assume N Fock states and two atom states 284 | D = displace(N, alpha) 285 | return D * rho * D.dag() 286 | 287 | def squeezed_dm(rho,sp): 288 | '''apply a squeeze operator to a density matrix 289 | @ var rho: the input density matrix for the motion only 290 | @ var sp: the squeeze parameter 291 | 292 | returns rho_sqz: a density matrix Qobj 293 | ''' 294 | N = shape(rho.data)[0] # assume N Fock states and two atom states 295 | S = squeez(N, sp) 296 | return S * rho * S.dag() 297 | 298 | if __name__ == '__main__': 299 | # Configure parameters 300 | pi = np.pi 301 | hbar = 1.05457173e-34 # planck constant 302 | mass = 40*1.67262178e-27 # ion mass 303 | frequency = 8.0e5 # trap frequency 304 | wavelength = 729.e-9 # laser wavelength 305 | lasertheta = 7. # projection of the laser to the motion axis in degrees 306 | bare_rabi_freq = 2*pi*200.e3 # bare Rabi frequency 307 | eta = lamb_dicke(wavelength, frequency, lasertheta) # Lamb-Dicke parameter 308 | 309 | if 1: # test carrier, bsb, rsb flops 310 | detuning = 2*pi*0.e3 # detuning from the driven transition 311 | N_noise = 1.e5 # noise equilibration temperature in quanta (make sure it is >> 1) 312 | gamma_motion = 1e-3 /N_noise # motion heating rate (1./N_noise means 1 quantum/sec ) 313 | gamma_motion_phi = 1e-3 # motion dephasing rate 314 | gamma_atom = 1e0 # atom dissipation rate 315 | N_hilbert = 200 # number of motional Fock states in Hilbert space 316 | 317 | # Atom state 318 | atom_up = basis(2,1) 319 | atom_down = basis(2,0) 320 | 321 | # Fock motional states 322 | #N_fock = 0 323 | #motionfock = basis(N_hilbert,N_fock) 324 | #rho_up = tensor( motionfock*motionfock.dag(), atom_up*atom_up.dag()) 325 | #rho_down = tensor( motionfock*motionfock.dag(), atom_down*atom_down.dag()) 326 | # Thermal motional state 327 | N_motion = 20 # Thermal component of the initial state 328 | motion_dm = thermal_dm(N_hilbert, N_motion) 329 | alpha = 10. # Coherent displacement amplitude 330 | motion_dm = displaced_dm(thermal_dm(N_hilbert, N_motion), alpha) 331 | rho_up = tensor(motion_dm, atom_up*atom_up.dag()) 332 | rho_down = tensor(motion_dm, atom_down*atom_down.dag()) 333 | 334 | # collapse operators 335 | c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 336 | # return operators 337 | r_ops = population_operators(N_hilbert) 338 | 339 | 340 | thetalist = linspace(0, 2*pi, 100) 341 | phi = 0. 342 | mesolution_c = carrier_flop(rho_down, bare_rabi_freq, eta, detuning, thetalist, phi, c_ops, r_ops) 343 | mesolution_b = bsb_flop(rho_down, bare_rabi_freq, eta, detuning, thetalist, phi, c_ops, r_ops) 344 | mesolution_r = rsb_flop(rho_down, bare_rabi_freq, eta, detuning, thetalist, phi, c_ops, r_ops) 345 | #print mesolution[1].states[-1] 346 | tlist_c = mesolution_c[0] 347 | populations_c = mesolution_c[1] 348 | tlist_b = mesolution_b[0] 349 | populations_b = mesolution_b[1] 350 | tlist_r = mesolution_r[0] 351 | populations_r = mesolution_r[1] 352 | #plot the results 353 | #plot(tlist, populations.expect[0]) 354 | plot(tlist_c, populations_c.expect[1], label = 'carrier') 355 | plot(tlist_b, populations_b.expect[1], label = 'bsb') 356 | plot(tlist_r, populations_r.expect[1], label = 'rsb') 357 | 358 | legend(loc = 1) 359 | xlabel('Time (s)') 360 | ylabel('Prob(D)') 361 | title('Rabi oscillations') 362 | show() 363 | 364 | if 0: # Ramsey experiment with thermal state and high dephasing rate 365 | detuning = 2*pi*3.e3 # detuning from the driven transition 366 | N_noise = 1.e5 # noise equilibration temperature in quanta (make sure it is >> 1) 367 | gamma_motion = 1e0 /N_noise # motion heating rate (1./N_noise means 1 quantum/sec ) 368 | gamma_motion_phi = 1e0 # motion dephasing rate 369 | gamma_atom = 1e0 # atom dissipation rate 370 | N_hilbert = 10 # number of motional Fock states in Hilbert space 371 | 372 | # Atom state 373 | atom_up = basis(2,1) 374 | atom_down = basis(2,0) 375 | 376 | # Motion state 377 | N_motion = 0. # Thermal component of the initial state 378 | motion_dm = thermal_dm(N_hilbert, N_motion) 379 | rho0 = tensor(motion_dm, atom_down*atom_down.dag()) # start with a thermal motional state and atom in ground state 380 | 381 | # collapse operators 382 | #c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 383 | # return operators 384 | r_ops = population_operators(N_hilbert) 385 | 386 | times = linspace(0., 8e-4, 100) 387 | 388 | gamma_motion_phi = 1e1 # motion dephasing rate 389 | c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 390 | ramsey_thermal = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 391 | plot(times,ramsey_thermal, color = 'blue', label = 'nbar = 0, Tphi = 100 ms') 392 | 393 | gamma_motion_phi = 1e2 394 | c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 395 | ramsey_thermal = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 396 | plot(times,ramsey_thermal, color = 'red', label = 'nbar = 0, Tphi = 10 ms') 397 | 398 | gamma_motion_phi = 1e3 399 | c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 400 | ramsey_thermal = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 401 | plot(times,ramsey_thermal, color = 'green', label = 'nbar = 0, Tphi = 1 ms') 402 | 403 | gamma_motion_phi = 1e4 404 | c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 405 | ramsey_thermal = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 406 | plot(times,ramsey_thermal, color = 'magenta', label = 'nbar = 0, Tphi = 100 us') 407 | 408 | #c_ops = collapse_operators(N_hilbert, N_noise, gamma_motion, gamma_motion_phi, gamma_atom) 409 | #ramsey_thermal = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 410 | #plot(times,ramsey_thermal, color = 'red', label = 'nbar = 15, Tphi = 500 us') 411 | 412 | legend(loc=1) 413 | xlabel('Time') 414 | ylabel('Prob(D)') 415 | title('Ramsey fringes') 416 | show() 417 | 418 | 419 | if 0: # Ramsey experiment with thermal, displaced thermal, squeezed-displaced thermal 420 | times = linspace(0., 1e-3, 300) 421 | 422 | N_motion = 5. # Thermal component of the initial state 423 | alpha = 0. # Coherent displacement amplitude 424 | motion_dm = displaced_dm(thermal_dm(N_hilbert, N_motion), alpha) 425 | rho0 = tensor(motion_dm, atom_down*atom_down.dag()) # start with a thermal motional state and atom in ground state 426 | ramsey_undisplaced = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 427 | plot(times,ramsey_undisplaced, color = 'blue', label = 'undisplaced state') 428 | N_motion = 5. # Thermal component of the initial state 429 | alpha = 8. # Coherent displacement amplitude 430 | motion_dm = displaced_dm(thermal_dm(N_hilbert, N_motion), alpha) 431 | rho0 = tensor(motion_dm, atom_down*atom_down.dag()) # start with a thermal motional state and atom in ground state 432 | ramsey_displaced = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 433 | plot(times,ramsey_displaced, color = 'red', label = 'displaced state') 434 | N_motion = 5. # Thermal component of the initial state 435 | alpha = 8. # Coherent displacement amplitude 436 | SP = 1. # Squeezing parameter 437 | motion_dm = squeezed_dm( displaced_dm( thermal_dm(N_hilbert, N_motion), alpha), SP) 438 | rho0 = tensor(motion_dm, atom_down*atom_down.dag()) # start with a thermal motional state and atom in ground state 439 | ramsey_displaced = ramsey_experiment(rho0, bare_rabi_freq, eta, detuning, times, c_ops, r_ops) 440 | plot(times,ramsey_displaced, color = 'green', label = 'squeezed-displaced state') 441 | 442 | legend(loc=1) 443 | xlabel('Time') 444 | ylabel('Prob(D)') 445 | title('Ramsey fringes') 446 | show() 447 | 448 | if 0: # Generate carrier operators 449 | N = 3 450 | a = tensor(destroy(N), qeye(2)) 451 | Wc = qeye(N) 452 | Wc.data = csr_matrix( qeye(N).data.dot( np.diag(rabi_coupling(N,0,eta) ) ) ) 453 | sm = tensor( Wc, destroy(2)) 454 | print Wc 455 | print sm + sm.dag() 456 | if 0: # Generate blue sideband operators 457 | N = 3 458 | a = tensor(destroy(N), qeye(2)) 459 | sm = tensor( qeye(N), destroy(2)) 460 | Wbsb = destroy(N).dag() # I will overwrite Wsb.data in the next line, but for pedagogical reasons... 461 | Wbsb.data = csr_matrix( destroy(N).dag().data.dot( np.diag( rabi_coupling(N,1,eta) / np.sqrt(np.linspace(1,N,N)) ) ) ) 462 | #Wbsb.data = csr_matrix( destroy(N).dag().data.dot( np.diag( rabi_coupling(N,1,eta) ) ) ) 463 | Absb = tensor(Wbsb.dag(), qeye(2)) 464 | print 'Wbsb', Wbsb 465 | print Absb.dag() * sm.dag() + Absb * sm 466 | if 0: # Generate red sideband operators 467 | N = 3 468 | a = tensor(destroy(N), qeye(2)) 469 | sm = tensor( qeye(N), destroy(2)) 470 | Wrsb = destroy(N) 471 | one_then_zero = ([float(x<1) for x in range(N)]) 472 | Wrsb.data = csr_matrix( destroy(N).data.dot( np.diag( rabi_coupling(N,-1,eta) / np.sqrt(one_then_zero+np.linspace(0,N-1,N)) ) ) ) 473 | #Wrsb.data = csr_matrix( destroy(N).data.dot( np.diag( rabi_coupling(N,-1,eta) ) ) ) 474 | Arsb = tensor(Wrsb, qeye(2)) 475 | print 'Wrsb', Wrsb 476 | print Arsb.dag() * sm + Arsb * sm.dag() 477 | 478 | --------------------------------------------------------------------------------