1
+ /* ******************************************************************************
2
+ * Copyright (c) 2022 - 2025 NVIDIA Corporation & Affiliates. *
3
+ * All rights reserved. *
4
+ * *
5
+ * This source code and the accompanying materials are made available under *
6
+ * the terms of the Apache License 2.0 which accompanies this distribution. *
7
+ ******************************************************************************/
8
+
9
+ // Compile and run with:
10
+ // ```
11
+ // nvq++ --target dynamics cavity_qed.cpp -o a.out && ./a.out
12
+ // ```
13
+
14
+ #include " cudaq/algorithms/evolve.h"
15
+ #include " cudaq/dynamics_integrators.h"
16
+ #include " cudaq/operators.h"
17
+ #include " export_csv_helper.h"
18
+ #include < cudaq.h>
19
+
20
+ int main () {
21
+
22
+ // `delta` represents the detuning between the two qubits.
23
+ // In physical terms, detuning is the energy difference (or frequency offset)
24
+ // between qubit levels. Detuning term (in angular frequency units).
25
+ double delta = 100 * 2 * M_PI;
26
+ // `J` is the static coupling strength between the two qubits.
27
+ // This terms facilitates energy exchange between the qubits, effectively
28
+ // coupling their dynamics.
29
+ double J = 7 * 2 * M_PI;
30
+ // `m_12` models spurious electromagnetic `crosstalk`.
31
+ // `Crosstalk` is an unwanted interaction , here represented as a fraction of
32
+ // the drive strength applied to the second qubit.
33
+ double m_12 = 0.2 ;
34
+ // `Omega` is the drive strength applied to the qubits.
35
+ // A driving field can induce transitions between qubit states.
36
+ double Omega = 20 * 2 * M_PI;
37
+
38
+ // For a spin-1/2 system, the raising operator S^+ and lowering operator S^-
39
+ // are defined as: S^+ = 0.5 * (X + `iY`) and S^- = 0.5 * (X - `iY`) These
40
+ // operators allow transitions between the spin states (|0> and |1>).
41
+ auto spin_plus = [](int degree) {
42
+ return 0.5 *
43
+ (cudaq::spin_operator::x (degree) +
44
+ std::complex<double >(0.0 , 1.0 ) * cudaq::spin_operator::y (degree));
45
+ };
46
+
47
+ auto spin_minus = [](int degree) {
48
+ return 0.5 *
49
+ (cudaq::spin_operator::x (degree) -
50
+ std::complex<double >(0.0 , 1.0 ) * cudaq::spin_operator::y (degree));
51
+ };
52
+
53
+ // The Hamiltonian describes the energy and dynamics of our 2-qubit system.
54
+ // It consist of several parts:
55
+ // 1. Detuning term for qubit 0: (delta / 2) * Z. This sets the energy
56
+ // splitting for qubit 0.
57
+ // 2. Exchange interaction: J * (S^-_1 * S^+_0 + S^+_1 * S^-_0). This couples
58
+ // the two qubits, enabling excitation transfer.
59
+ // 3. Drive on qubit 0: Omega * X. A control field that drives transition in
60
+ // qubit 0.
61
+ // 4. `Crosstalk` drive on qubit 1: m_12 * Omega * X. A reduces drive on qubit
62
+ // 1 due to electromagnetic `crosstalk`.
63
+ auto hamiltonian =
64
+ (delta / 2.0 ) * cudaq::spin_operator::z (0 ) +
65
+ J * (spin_minus (1 ) * spin_plus (0 ) + spin_plus (1 ) * spin_minus (0 )) +
66
+ Omega * cudaq::spin_operator::x (0 ) +
67
+ m_12 * Omega * cudaq::spin_operator::x (1 );
68
+
69
+ // Each qubit is a 2-level system (dimension 2).
70
+ // The composite system (two qubits) has a total Hilbert space dimension of 2
71
+ // * 2 = 4.
72
+ cudaq::dimension_map dimensions{{0 , 2 }, {1 , 2 }};
73
+
74
+ // Build the initial state
75
+ // psi_00 represents the state |00> (both qubits in the ground state).
76
+ // psi_10 represents the state |10> (first qubit excited, second qubit in the
77
+ // ground state).
78
+ std::vector<std::complex<double >> psi00_data = {
79
+ {1.0 , 0.0 }, {0.0 , 0.0 }, {0.0 , 0.0 }, {0.0 , 0.0 }};
80
+ std::vector<std::complex<double >> psi10_data = {
81
+ {0.0 , 0.0 }, {1.0 , 0.0 }, {0.0 , 0.0 }, {0.0 , 0.0 }};
82
+
83
+ // Two initial state vectors for the 2-qubit system (dimension 4)
84
+ auto psi_00 = cudaq::state::from_data (psi00_data);
85
+ auto psi_10 = cudaq::state::from_data (psi10_data);
86
+
87
+ // Create a list of time steps for the simulation.
88
+ // Here we use 1001 points linearly spaced between time 0 and 1.
89
+ // This schedule will be used to integrate the time evolution of the system.
90
+ const int num_steps = 1001 ;
91
+ std::vector<double > steps = cudaq::linspace (0.0 , 1.0 , num_steps);
92
+ cudaq::schedule schedule (steps);
93
+
94
+ // Use Runge-`Kutta` integrator (4`th` order) to solve the time-dependent
95
+ // evolution. `dt` is the integration time step, and `order` sets the accuracy
96
+ // of the integrator method.
97
+ cudaq::integrators::runge_kutta integrator (4 , 0.0001 );
98
+
99
+ // The observables are the spin components along the x, y, and z directions
100
+ // for both qubits. These observables will be measured during the evolution.
101
+ auto observables = {cudaq::spin_operator::x (0 ), cudaq::spin_operator::y (0 ),
102
+ cudaq::spin_operator::z (0 ), cudaq::spin_operator::x (1 ),
103
+ cudaq::spin_operator::y (1 ), cudaq::spin_operator::z (1 )};
104
+
105
+ // Evolution with 2 initial states
106
+ // We evolve the system under the defined Hamiltonian for both initial states
107
+ // simultaneously. No collapsed operators are provided (closed system
108
+ // evolution). The evolution returns expectation values for all defined
109
+ // observables at each time step.
110
+ auto evolution_results =
111
+ cudaq::evolve (hamiltonian, dimensions, schedule, {psi_00, psi_10},
112
+ integrator, {}, observables, true );
113
+
114
+ // Retrieve the evolution result corresponding to each initial state.
115
+ auto &evolution_result_00 = evolution_results[0 ];
116
+ auto &evolution_result_10 = evolution_results[1 ];
117
+
118
+ // Lambda to extract expectation values for a given observable index
119
+ auto get_expectation = [](int idx, auto &result) -> std::vector<double > {
120
+ std::vector<double > expectations;
121
+
122
+ auto all_exps = result.get_expectation_values ().value ();
123
+ for (auto exp_vals : all_exps) {
124
+ expectations.push_back ((double )exp_vals[idx]);
125
+ }
126
+ return expectations;
127
+ };
128
+
129
+ // For the two `evolutions`, extract the six observable trajectories.
130
+ // For the |00> initial state, we extract the expectation trajectories for
131
+ // each observable.
132
+ auto result_00_0 = get_expectation (0 , evolution_result_00);
133
+ auto result_00_1 = get_expectation (1 , evolution_result_00);
134
+ auto result_00_2 = get_expectation (2 , evolution_result_00);
135
+ auto result_00_3 = get_expectation (3 , evolution_result_00);
136
+ auto result_00_4 = get_expectation (4 , evolution_result_00);
137
+ auto result_00_5 = get_expectation (5 , evolution_result_00);
138
+
139
+ // Similarly, for the |10> initial state:
140
+ auto result_10_0 = get_expectation (0 , evolution_result_10);
141
+ auto result_10_1 = get_expectation (1 , evolution_result_10);
142
+ auto result_10_2 = get_expectation (2 , evolution_result_10);
143
+ auto result_10_3 = get_expectation (3 , evolution_result_10);
144
+ auto result_10_4 = get_expectation (4 , evolution_result_10);
145
+ auto result_10_5 = get_expectation (5 , evolution_result_10);
146
+
147
+ // Export the results to a `CSV` file
148
+ // Export the Z-component of qubit 1's expectation values for both initial
149
+ // states. The `CSV` file "cross_resonance_z.`csv`" will have time versus (Z1)
150
+ // data for both |00> and |10> initial conditions.
151
+ export_csv (" cross_resonance_z.csv" , " time" , steps, " <Z1>_00" , result_00_5,
152
+ " <Z1>_10" , result_10_5);
153
+ // Export the Y-component of qubit 1's expectation values for both initial
154
+ // states. The `CSV` file "cross_resonance_y.`csv`" will have time versus (Y1)
155
+ // data.
156
+ export_csv (" cross_resonance_y.csv" , " time" , steps, " <Y1>_00" , result_00_4,
157
+ " <Y1>_10" , result_10_4);
158
+
159
+ std::cout
160
+ << " Simulation complete. The results are saved in cross_resonance_z.csv "
161
+ " and cross_resonance_y.csv files."
162
+ << std::endl;
163
+ return 0 ;
164
+ }
0 commit comments