#include <circuit.hpp>

Public Member Functions

 QuantumCircuit (UINT qubit_count)
 
 QuantumCircuit (std::string qasm_path, std::string qasm_loader_script_path="qasmloader.py")
 
QuantumCircuitcopy () const
 
virtual ~QuantumCircuit ()
 
virtual void add_gate (QuantumGateBase *gate)
 
virtual void add_gate (QuantumGateBase *gate, UINT index)
 
virtual void add_gate_copy (const QuantumGateBase &gate)
 
virtual void add_gate_copy (const QuantumGateBase &gate, UINT index)
 
virtual void remove_gate (UINT index)
 
void update_quantum_state (QuantumStateBase *state)
 
void update_quantum_state (QuantumStateBase *state, UINT start_index, UINT end_index)
 
bool is_Clifford () const
 
bool is_Gaussian () const
 
UINT calculate_depth () const
 
virtual std::string to_string () const
 
virtual void add_X_gate (UINT target_index)
 
virtual void add_Y_gate (UINT target_index)
 
virtual void add_Z_gate (UINT target_index)
 
virtual void add_H_gate (UINT target_index)
 
virtual void add_S_gate (UINT target_index)
 
virtual void add_Sdag_gate (UINT target_index)
 
virtual void add_T_gate (UINT target_index)
 
virtual void add_Tdag_gate (UINT target_index)
 
virtual void add_sqrtX_gate (UINT target_index)
 
virtual void add_sqrtXdag_gate (UINT target_index)
 
virtual void add_sqrtY_gate (UINT target_index)
 
virtual void add_sqrtYdag_gate (UINT target_index)
 
virtual void add_P0_gate (UINT target_index)
 
virtual void add_P1_gate (UINT target_index)
 
virtual void add_CNOT_gate (UINT control_index, UINT target_index)
 
virtual void add_CZ_gate (UINT control_index, UINT target_index)
 
virtual void add_SWAP_gate (UINT target_index1, UINT target_index2)
 
virtual void add_RX_gate (UINT target_index, double angle)
 
virtual void add_RY_gate (UINT target_index, double angle)
 
virtual void add_RZ_gate (UINT target_index, double angle)
 
virtual void add_U1_gate (UINT target_index, double phi)
 
virtual void add_U2_gate (UINT target_index, double phi, double psi)
 
virtual void add_U3_gate (UINT target_index, double phi, double psi, double lambda)
 
virtual void add_multi_Pauli_gate (std::vector< UINT > target_index_list, std::vector< UINT > pauli_id_list)
 
virtual void add_multi_Pauli_gate (const PauliOperator &pauli_operator)
 
virtual void add_multi_Pauli_rotation_gate (std::vector< UINT > target_index_list, std::vector< UINT > pauli_id_list, double angle)
 
virtual void add_multi_Pauli_rotation_gate (const PauliOperator &pauli_operator)
 
virtual void add_diagonal_observable_rotation_gate (const Observable &observable, double angle)
 
virtual void add_observable_rotation_gate (const Observable &observable, double angle, UINT num_repeats=0)
 
virtual void add_dense_matrix_gate (UINT target_index, const ComplexMatrix &matrix)
 
virtual void add_dense_matrix_gate (std::vector< UINT > target_index_list, const ComplexMatrix &matrix)
 

Public Attributes

const UINT & qubit_count = _qubit_count
 
const std::vector< QuantumGateBase * > & gate_list = _gate_list
 

Friends

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

Detailed Description

量子回路のクラス

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

Constructor & Destructor Documentation

◆ QuantumCircuit() [1/2]

QuantumCircuit::QuantumCircuit ( UINT  qubit_count)

空の量子回路を作成する

Parameters
[in]qubit_count量子ビットの数
Returns
生成された量子回路のインスタンス

◆ QuantumCircuit() [2/2]

QuantumCircuit::QuantumCircuit ( std::string  qasm_path,
std::string  qasm_loader_script_path = "qasmloader.py" 
)

量子回路をQASMから生成する

QASMのファイルをパスを指定し、QASMに記載されている量子回路を作成する。 QASMのパースのためにqiskitが指定されている必要がある。

Parameters
[in]qasm_pathQASMファイルのパス
[in]qasm_loader_script_pathQASMを読むためのpythonのパス
Returns
生成されたインスタンス。生成でエラーが生じた場合はNULLを返す。

◆ ~QuantumCircuit()

QuantumCircuit::~QuantumCircuit ( )
virtual

デストラクタ

Member Function Documentation

◆ add_CNOT_gate()

void QuantumCircuit::add_CNOT_gate ( UINT  control_index,
UINT  target_index 
)
virtual

CNOT gateを追加する。

Parameters
[in]control_index作用するcontrol qubitの添え字
[in]target_index作用するtarget qubitの添え字

◆ add_CZ_gate()

void QuantumCircuit::add_CZ_gate ( UINT  control_index,
UINT  target_index 
)
virtual

Control-Z gateを追加する。

Parameters
[in]control_index作用するcontrol qubitの添え字
[in]target_index作用するtarget qubitの添え字

◆ add_dense_matrix_gate() [1/2]

void QuantumCircuit::add_dense_matrix_gate ( UINT  target_index,
const ComplexMatrix matrix 
)
virtual

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

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]matrix作用する2*2の行列

◆ add_dense_matrix_gate() [2/2]

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

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

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

Parameters
[in]target_index_list作用するtarget qubitの添え字のリスト
[in]matrix作用する行列

◆ add_diagonal_observable_rotation_gate()

void QuantumCircuit::add_diagonal_observable_rotation_gate ( const Observable observable,
double  angle 
)
virtual

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

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

Parameters
[in]observable追加するオブザーバブル
[in]angle回転角

◆ add_gate() [1/2]

void QuantumCircuit::add_gate ( QuantumGateBase gate)
virtual

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

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

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

◆ add_gate() [2/2]

void QuantumCircuit::add_gate ( QuantumGateBase gate,
UINT  index 
)
virtual

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

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

Parameters
[in]gate追加する量子ゲート
[in]index追加する位置

◆ add_gate_copy() [1/2]

void QuantumCircuit::add_gate_copy ( const QuantumGateBase gate)
virtual

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

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

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

◆ add_gate_copy() [2/2]

void QuantumCircuit::add_gate_copy ( const QuantumGateBase gate,
UINT  index 
)
virtual

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

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

Parameters
[in]gate追加する量子ゲート
[in]index追加する位置

◆ add_H_gate()

void QuantumCircuit::add_H_gate ( UINT  target_index)
virtual

Hadamard gateを追加する。

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

◆ add_multi_Pauli_gate() [1/2]

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

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パウリゲートの計算コストはほぼ同じため、パウリゲートのテンソル積を作用する場合はパウリゲートとして作用した方が処理が高速になる。

Parameters
[in]target_index_list作用するtarget qubitの添え字のリスト
[in]pauli_id_listtarget_index_listに対応したパウリ演算子のid

◆ add_multi_Pauli_gate() [2/2]

void QuantumCircuit::add_multi_Pauli_gate ( const PauliOperator pauli_operator)
virtual

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

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

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

◆ add_multi_Pauli_rotation_gate() [1/2]

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

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パウリゲートの計算コストはほぼ同じため、パウリゲートのテンソル積を作用する場合はパウリゲートとして作用した方が処理が高速になる。

Parameters
[in]target_index_list作用するtarget qubitの添え字のリスト
[in]pauli_id_listtarget_index_listに対応したパウリ演算子のid
[in]angle回転角

◆ add_multi_Pauli_rotation_gate() [2/2]

void QuantumCircuit::add_multi_Pauli_rotation_gate ( const PauliOperator pauli_operator)
virtual

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

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

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

◆ add_observable_rotation_gate()

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

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

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

Parameters
[in]observable追加するオブザーバブル
[in]angle回転角 \( \theta \)
[in]num_repeatsTrotter展開をする際の分割数 \(N\)。指定しない場合は、関数内部で \( \#qubit \cdot \theta / N = 0.01\) となるように設定される。

◆ add_P0_gate()

void QuantumCircuit::add_P0_gate ( UINT  target_index)
virtual

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

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

◆ add_P1_gate()

void QuantumCircuit::add_P1_gate ( UINT  target_index)
virtual

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

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

◆ add_RX_gate()

void QuantumCircuit::add_RX_gate ( UINT  target_index,
double  angle 
)
virtual

X-rotation gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]angle回転角 \(\theta\)

◆ add_RY_gate()

void QuantumCircuit::add_RY_gate ( UINT  target_index,
double  angle 
)
virtual

Y-rotation gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]angle回転角 \(\theta\)

◆ add_RZ_gate()

void QuantumCircuit::add_RZ_gate ( UINT  target_index,
double  angle 
)
virtual

Z-rotation gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]angle回転角 \(\theta\)

◆ add_S_gate()

void QuantumCircuit::add_S_gate ( UINT  target_index)
virtual

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

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

◆ add_Sdag_gate()

void QuantumCircuit::add_Sdag_gate ( UINT  target_index)
virtual

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

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

◆ add_sqrtX_gate()

void QuantumCircuit::add_sqrtX_gate ( UINT  target_index)
virtual

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

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

◆ add_sqrtXdag_gate()

void QuantumCircuit::add_sqrtXdag_gate ( UINT  target_index)
virtual

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

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

◆ add_sqrtY_gate()

void QuantumCircuit::add_sqrtY_gate ( UINT  target_index)
virtual

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

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

◆ add_sqrtYdag_gate()

void QuantumCircuit::add_sqrtYdag_gate ( UINT  target_index)
virtual

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

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

◆ add_SWAP_gate()

void QuantumCircuit::add_SWAP_gate ( UINT  target_index1,
UINT  target_index2 
)
virtual

SWAP gateを追加する。

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

◆ add_T_gate()

void QuantumCircuit::add_T_gate ( UINT  target_index)
virtual

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

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

◆ add_Tdag_gate()

void QuantumCircuit::add_Tdag_gate ( UINT  target_index)
virtual

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

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

◆ add_U1_gate()

void QuantumCircuit::add_U1_gate ( UINT  target_index,
double  phi 
)
virtual

OpenQASMのu1 gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]phi回転角 \(\phi\)

◆ add_U2_gate()

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

OpenQASMのu2 gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]phi回転角 \(\phi\)
[in]psi回転角 \(\psi\)

◆ add_U3_gate()

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

OpenQASMのu3 gateを追加する。

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

Parameters
[in]target_index作用するtarget qubitの添え字
[in]phi回転角 \(\phi\)
[in]psi回転角 \(\psi\)
[in]lambda回転角 \(\lambda\)

◆ add_X_gate()

void QuantumCircuit::add_X_gate ( UINT  target_index)
virtual

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

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

◆ add_Y_gate()

void QuantumCircuit::add_Y_gate ( UINT  target_index)
virtual

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

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

◆ add_Z_gate()

void QuantumCircuit::add_Z_gate ( UINT  target_index)
virtual

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

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

◆ calculate_depth()

UINT QuantumCircuit::calculate_depth ( ) const

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

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

Returns
量子回路のdepth

◆ copy()

QuantumCircuit * QuantumCircuit::copy ( ) const

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

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

◆ is_Clifford()

bool QuantumCircuit::is_Clifford ( ) const

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

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

Return values
trueClifford
falseNon-Clifford

◆ is_Gaussian()

bool QuantumCircuit::is_Gaussian ( ) const

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

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

Return values
trueFermionic Gaussian
falseNon-fermionic Gaussian

◆ remove_gate()

void QuantumCircuit::remove_gate ( UINT  index)
virtual

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

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

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

◆ to_string()

std::string QuantumCircuit::to_string ( ) const
virtual

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

Returns
生成した文字列

◆ update_quantum_state() [1/2]

void QuantumCircuit::update_quantum_state ( QuantumStateBase state)

量子状態を更新する

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

Parameters
[in,out]state作用する量子状態

◆ update_quantum_state() [2/2]

void QuantumCircuit::update_quantum_state ( QuantumStateBase state,
UINT  start_index,
UINT  end_index 
)

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

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

Parameters
[in,out]state作用する量子状態
[in]start_index開始位置
[in]end_index修了位置

Friends And Related Function Documentation

◆ operator<< [1/2]

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

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

Returns
受け取ったストリーム

◆ operator<< [2/2]

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

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

Returns
受け取ったストリーム

Member Data Documentation

◆ gate_list

const std::vector<QuantumGateBase*>& QuantumCircuit::gate_list = _gate_list

量子ゲートのリスト

◆ qubit_count

const UINT& QuantumCircuit::qubit_count = _qubit_count

量子ビットの数


The documentation for this class was generated from the following files:
  • /Users/tenninyan/Soft/qunasys/qulacs/src/cppsim/circuit.hpp
  • /Users/tenninyan/Soft/qunasys/qulacs/src/cppsim/circuit.cpp