Class QuantumCircuit

Inheritance Relationships

Derived Type

Class Documentation

class QuantumCircuit

量子回路のクラス

量子回路を管理するクラス。QuantumGateクラスをリストとして持ち、種々の操作を行う。 管理する量子ゲートは量子回路の解放時にすべて解放される。

Subclassed by ParametricQuantumCircuit

Public Functions

explicit QuantumCircuit(UINT qubit_count)

空の量子回路を作成する

パラメータ

qubit_count -- [in] 量子ビットの数

戻り値

生成された量子回路のインスタンス

QuantumCircuit *copy() const

量子回路のディープコピーを生成する

戻り値

量子回路のディープコピー

virtual ~QuantumCircuit()

デストラクタ

virtual void add_gate(QuantumGateBase *gate)

量子ゲートを回路の末尾に追加する

量子ゲートを回路に追加する。 追加した量子ゲートは量子回路の解放時に開放される。

パラメータ

gate -- [in] 追加する量子ゲート

virtual void add_gate(QuantumGateBase *gate, UINT index)

量子ゲートを回路の指定位置に追加する。

量子ゲートを回路の指定した位置に追加する。 追加した量子ゲートは量子回路の解放時に開放される。

パラメータ
  • gate -- [in] 追加する量子ゲート

  • index -- [in] 追加する位置

virtual void add_gate_copy(const QuantumGateBase *gate)

量子ゲートを回路の末尾に追加する

与えられた量子ゲートのコピーを回路に追加する。 add_gateに比べコピーが発生する分低速な一方、引数で与えたゲートを再利用できる。

パラメータ

gate -- [in] 追加する量子ゲート

virtual void add_gate_copy(const QuantumGateBase *gate, UINT index)

量子ゲートを回路の指定位置に追加する。

与えらた量子ゲートを回路の指定した位置に追加する。

パラメータ
  • gate -- [in] 追加する量子ゲート

  • index -- [in] 追加する位置

virtual void add_noise_gate(QuantumGateBase *gate, std::string noise_type, double noise_prob)

ノイズ付き量子ゲートを回路の末尾に追加する

ノイズ付き量子ゲートを回路に追加する。 追加したノイズ付き量子ゲートは量子回路の解放時に開放される。

パラメータ
  • gate -- [in] 追加するノイズ付き量子ゲート

  • noise_type -- [in] 追加するノイズの種類

  • noise_prob -- [in] ノイズが発生する確率

void add_noise_gate_copy(QuantumGateBase *gate, std::string noise_type, double noise_prob)

ノイズ付き量子ゲートのコピーを回路の末尾に追加する

ノイズ付き量子ゲートのコピーを回路に追加する。 add_noise_gateに比べコピーが発生する分低速な一方、引数で与えたゲートを再利用できる。

パラメータ
  • gate -- [in] 追加するノイズ付き量子ゲート

  • noise_type -- [in] 追加するノイズの種類

  • noise_prob -- [in] ノイズが発生する確率

virtual void remove_gate(UINT index)

量子回路からゲートを削除する。

削除した量子ゲートは解放される。

パラメータ

index -- [in] 削除するゲートの位置

inline void merge_circuit(const QuantumCircuit *circuit)

量子回路をマージする。

引数で与えた量子回路のゲートを後ろに追加していく。 マージされた側の量子回路に変更を加えてもマージした側の量子回路には変更は加わらないことに注意する。 circuit1.add_circuit(circuit2) circuit2.add_gate(gate) # これをしても、circuit1にgateは追加されない

パラメータ

circuit -- [in] マージする量子回路

void update_quantum_state(QuantumStateBase *state)

量子状態を更新する

順番にすべての量子ゲートを作用する。量子状態の初期化などは行わない。

パラメータ

state -- [inout] 作用する量子状態

void update_quantum_state(QuantumStateBase *state, UINT start_index, UINT end_index)

量子回路の指定範囲のみを用いて量子状態をを更新する

添え字がstart_indexからend_index-1までの量子ゲートを順番に量子ゲートを作用する。量子状態の初期化などは行わない。

パラメータ
  • state -- [inout] 作用する量子状態

  • start_index -- [in] 開始位置

  • end_index -- [in] 修了位置

bool is_Clifford() const

量子回路がCliffordかどうかを判定する。

全ての量子ゲートがCliffordである場合にtrueと判定される。 Non-Cliffordゲートが複数あり積がCliffordとなっている場合もfalseとして判定される点に注意。

戻り値
  • true -- Clifford

  • false -- Non-Clifford

bool is_Gaussian() const

量子回路がFermionic Gaussianかどうかを判定する。

全ての量子ゲートがFermionic Gaussianである場合にtrueと判定される。 Non-Gaussianゲートが複数あり、結果としてGaussianとなっている場合もNon-Gaussianとして判定される点に注意。

戻り値
  • true -- Fermionic Gaussian

  • false -- Non-fermionic Gaussian

UINT calculate_depth() const

量子回路のdepthを計算する。

ここでいうdepthとは、可能な限り量子ゲートを並列実行した時に掛かるステップ数を指す。

戻り値

量子回路のdepth

virtual std::string to_string() const

量子回路のデバッグ情報の文字列を生成する

戻り値

生成した文字列

virtual void add_X_gate(UINT target_index)

\(X\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_Y_gate(UINT target_index)

\(Y\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_Z_gate(UINT target_index)

\(Z\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_H_gate(UINT target_index)

Hadamard gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_S_gate(UINT target_index)

\(S\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_Sdag_gate(UINT target_index)

\(S^{\dagger}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_T_gate(UINT target_index)

\(T\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_Tdag_gate(UINT target_index)

\(T^{\dagger}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_sqrtX_gate(UINT target_index)

\(\sqrt{X}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_sqrtXdag_gate(UINT target_index)

\(\sqrt{X}^{\dagger}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_sqrtY_gate(UINT target_index)

\(\sqrt{Y}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_sqrtYdag_gate(UINT target_index)

\(\sqrt{Y}^{\dagger}\) gateを追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_P0_gate(UINT target_index)

0状態への射影演算を追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_P1_gate(UINT target_index)

1状態への射影演算を追加する。

パラメータ

target_index -- [in] 作用する量子ビットの添え字

virtual void add_CNOT_gate(UINT control_index, UINT target_index)

CNOT gateを追加する。

パラメータ
  • control_index -- [in] 作用するcontrol qubitの添え字

  • target_index -- [in] 作用するtarget qubitの添え字

virtual void add_CZ_gate(UINT control_index, UINT target_index)

Control-Z gateを追加する。

パラメータ
  • control_index -- [in] 作用するcontrol qubitの添え字

  • target_index -- [in] 作用するtarget qubitの添え字

virtual void add_SWAP_gate(UINT target_index1, UINT target_index2)

SWAP gateを追加する。

パラメータ
  • target_index1 -- [in] 作用するtarget qubitの添え字

  • target_index2 -- [in] 作用するもう一方のtarget qubitの添え字

virtual void add_RX_gate(UINT target_index, double angle)

X-rotation gateを追加する。

ゲートの表記は \( R_X(\theta) = \exp(i\theta X) \)になっている。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RY_gate(UINT target_index, double angle)

Y-rotation gateを追加する。

ゲートの表記は \( R_Y(\theta) = \exp(i\theta Y) \)になっている。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RZ_gate(UINT target_index, double angle)

Z-rotation gateを追加する。

ゲートの表記は \( R_Z(\theta) = \exp(i\theta Z) \)になっている。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotInvX_gate(UINT target_index, double angle)

X-rotation gateを追加する。

ゲートの表記は \( R_X(\theta) = \exp(i\theta X) \)になっている。 一般的な表記に対して逆向き,qulacsのRXと同じ向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotInvY_gate(UINT target_index, double angle)

Y-rotation gateを追加する。

ゲートの表記は \( R_Y(\theta) = \exp(i\theta Y) \)になっている。 一般的な表記に対して逆向き,qulacsのRYと同じ向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotInvZ_gate(UINT target_index, double angle)

Z-rotation gateを追加する。

ゲートの表記は \( R_Z(\theta) = \exp(i\theta Z) \)になっている。 一般的な表記に対して逆向き,qulacsのRZと同じ向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotX_gate(UINT target_index, double angle)

X-rotation gateを追加する。

ゲートの表記は \( R_X(\theta) = \exp(-i\theta X) \)になっている。 一般的な表記と同じ向き,qulacsのRXに対して逆向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotY_gate(UINT target_index, double angle)

Y-rotation gateを追加する。

ゲートの表記は \( R_Y(\theta) = \exp(-i\theta Y) \)になっている。 一般的な表記と同じ向き,qulacsのRYに対して逆向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_RotZ_gate(UINT target_index, double angle)

Z-rotation gateを追加する。

ゲートの表記は \( R_Z(\theta) = \exp(-i\theta Z) \)になっている。 一般的な表記と同じ向き,qulacsのRZに対して逆向きである。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • angle -- [in] 回転角 \(\theta\)

virtual void add_U1_gate(UINT target_index, double phi)

OpenQASMのu1 gateを追加する。

ゲートの表記はIBMQのページを参照。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • phi -- [in] 回転角 \(\phi\)

virtual void add_U2_gate(UINT target_index, double phi, double psi)

OpenQASMのu2 gateを追加する。

ゲートの表記はIBMQのページを参照。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • phi -- [in] 回転角 \(\phi\)

  • psi -- [in] 回転角 \(\psi\)

virtual void add_U3_gate(UINT target_index, double phi, double psi, double lambda)

OpenQASMのu3 gateを追加する。

ゲートの表記はIBMQのページを参照。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • phi -- [in] 回転角 \(\phi\)

  • psi -- [in] 回転角 \(\psi\)

  • lambda -- [in] 回転角 \(\lambda\)

virtual void add_multi_Pauli_gate(std::vector<UINT> target_index_list, std::vector<UINT> pauli_id_list)

n-qubitパウリゲートを追加する。

n-qubitパウリゲートを作用する。 パウリ演算子は \((I,X,Y,Z) \mapsto (0,1,2,3)\)と対応している。 例えば、 \(X_3 Y_2 Z_4\)であれば、target_index_list = {3,2,4}, pauli_id_list = {1,2,3}である。 1-qubitパウリゲートとn-qubitパウリゲートの計算コストはほぼ同じため、パウリゲートのテンソル積を作用する場合はパウリゲートとして作用した方が処理が高速になる。

パラメータ
  • target_index_list -- [in] 作用するtarget qubitの添え字のリスト

  • pauli_id_list -- [in] target_index_listに対応したパウリ演算子のid

virtual void add_multi_Pauli_gate(const PauliOperator &pauli_operator)

n-qubitパウリゲートを追加する。

n-qubitパウリゲートを作用する。

パラメータ

pauli_operator -- [in] 追加するパウリ演算子

virtual void add_multi_Pauli_rotation_gate(std::vector<UINT> target_index_list, std::vector<UINT> pauli_id_list, double angle)

n-qubitパウリ回転ゲートを追加する。

n-qubitパウリ回転ゲートを作用する。 パウリ演算子は{I,X,Y,Z} = {0,1,2,3}と対応している。 例えば、 \(\exp(i\theta X_3 Y_2 Z_4)\)であれば、target_index_list = {3,2,4}, pauli_id_list = {1,2,3}, angle = \(\theta\)とする。 1-qubitパウリゲートとn-qubitパウリゲートの計算コストはほぼ同じため、パウリゲートのテンソル積を作用する場合はパウリゲートとして作用した方が処理が高速になる。

パラメータ
  • target_index_list -- [in] 作用するtarget qubitの添え字のリスト

  • pauli_id_list -- [in] target_index_listに対応したパウリ演算子のid

  • angle -- [in] 回転角

virtual void add_multi_Pauli_rotation_gate(const PauliOperator &pauli_operator)

n-qubitパウリ回転ゲートを追加する。 n-qubitパウリ回転ゲートを作用する。 回転角はPauliOperatorの係数を用いる。

パラメータ

pauli_operator -- [in] 追加するパウリ演算子

virtual void add_diagonal_observable_rotation_gate(const Observable &observable, double angle)

n-qubitオブザーバブル回転ゲートを追加する。(対角のみ)

n-qubitオブザーバブル回転ゲートを作用する。ここで用いるオブザーバブルは、対角である必要がある。

パラメータ
  • observable -- [in] 追加するオブザーバブル

  • angle -- [in] 回転角

virtual void add_observable_rotation_gate(const Observable &observable, double angle, UINT num_repeats = 0)

n-qubitオブザーバブル回転ゲートを追加する。

Suzuki-Trotter展開によりn-qubitオブザーバブル回転ゲートを作用する。ここで用いるオブザーバブルは、対角でなくてもよい。

パラメータ
  • observable -- [in] 追加するオブザーバブル

  • angle -- [in] 回転角 \( \theta \)

  • num_repeats -- [in] Trotter展開をする際の分割数 \(N\)。指定しない場合は、関数内部で \( \#qubit \cdot \theta / N = 0.01\) となるように設定される。

virtual void add_dense_matrix_gate(UINT target_index, const ComplexMatrix &matrix)

1-qubitのdenseな行列のゲートを追加する。

denseな行列はユニタリである必要はなく、射影演算子やクラウス演算子でも良い。

パラメータ
  • target_index -- [in] 作用するtarget qubitの添え字

  • matrix -- [in] 作用する2*2の行列

virtual void add_dense_matrix_gate(std::vector<UINT> target_index_list, const ComplexMatrix &matrix)

multi qubitのdenseな行列のゲートを追加する。

denseな行列はユニタリである必要はなく、射影演算子やクラウス演算子でも良い。 matrixの次元はtarget_index_listのサイズを \(m\)としたとき \(2^m\)でなければいけない。

パラメータ
  • target_index_list -- [in] 作用するtarget qubitの添え字のリスト

  • matrix -- [in] 作用する行列

virtual void add_random_unitary_gate(std::vector<UINT> target_index_list)

multi qubitのランダムユニタリゲートを追加する。

パラメータ
  • target_index_list -- [in] 作用するtarget qubitの添え字のリスト

  • seed -- [in] 乱数のseed値

virtual void add_random_unitary_gate(std::vector<UINT> target_index_list, UINT seed)
virtual QuantumCircuit *get_inverse(void)

Public Members

const UINT &qubit_count

量子ビットの数

const std::vector<QuantumGateBase*> &gate_list

量子ゲートのリスト

Protected Functions

QuantumCircuit(const QuantumCircuit &obj)
QuantumCircuit &operator=(const QuantumCircuit&) = delete

Protected Attributes

std::vector<QuantumGateBase*> _gate_list
UINT _qubit_count

Friends

friend std::ostream &operator<<(std::ostream &os, const QuantumCircuit&)

量子回路のデバッグ情報を出力する。

戻り値

受け取ったストリーム

friend std::ostream &operator<<(std::ostream &os, const QuantumCircuit *gate)

量子回路のデバッグ情報を出力する。

戻り値

受け取ったストリーム