Circuit Class
Circuit
- class SQcircuit.Circuit(elements, flux_dist='junctions')[source]
- Bases: - object- Class that contains circuit properties and builds the Hamiltonian using the theory discussed in the original SQcircuit paper. Provides methods to calculate: - Eigenvalues and eigenvectors 
- Phase coordinate representation of eigenvectors 
- Coupling operators 
- Matrix elements 
- Decoherence rates 
- Gradients of Hamiltonian, eigenvalues/vectors , and Decoherence 
 - Parameters:
- elements ( - Dict[- Tuple[- int,- int],- List[- Element]]) – A dictionary that contains the circuit’s elements at each edge of the circuit.
- flux_dist ( - str) – Provide the method of distributing the external fluxes. If- flux_distis- "all", SQcircuit assign the external fluxes based on the capacitor of each inductive element (This option is necessary for time-dependent external fluxes). If- flux_distis- "inductor"SQcircuit finds the external flux distribution by assuming the capacitor of the inductors are much smaller than the junction capacitors, If- flux_distis- "junction"it is the other way around.
 
 - update()[source]
- Update the circuit Hamiltonian to reflect in-place changes made to the scalar values used for circuit elements (ex. C, L, J…). - Return type:
- None
 
 - property efreqs: ndarray | Tensor
- Eigenfrequencies in the chosen frequency unit for SQcircuit. If the SQcircuit engine is - PyTorch, the efreqs will be in- Tensorformat; otherwise, they will be in- ndarrayformat.
 - property evecs: List[Qobj] | Tensor
- List of circuit eigenvectors. If the SQcircuit engine is - PyTorch, each eigenvector will be in- Tensorformat; otherwise, they will be in- Qutip.Qobjformat.
 - property trunc_nums: List[int]
- List of truncation numbers of the circuit. For harmonic modes, these are N where the Hilbert space is 0, 1, …, (N-1) and for charge modes these are N where the Hilbert space is -(N-1), …, 0, …, (N-1). 
 - property parameters: List[Tensor]
- The values of the elements in the circuit which require gradient. The parameters can be set by a list of new values for each of the elements. Only valid if - get_optim_mode() = Trueor we are using- PyTorchengine of SQcircuit.
 - property parameters_grad: List[Tensor | None] | Tensor
- Return the gradients of the tensors in - .parameters. If all values are not- None, it is returned as a stacked- Tensor, otherwise as a list of individual values.
 - property parameters_dict: OrderedDict[Tuple[Element | Loop, Tensor]]
- The dictionary of (element, value) pairs for the elements in the circuit which require gradient. 
 - property parameters_elems: List[Element | Loop]
- The elements in the circuit which require gradient. 
 - get_params_type()[source]
- List of the types for each element in the circuit’s parameters. - Return type:
- List[- Union[- Type[- Element],- Type[- Loop]]]
 
 - zero_parameters_grad()[source]
- Set the gradient of all values in self.parameters to None. - Return type:
- None
 
 - safecopy(save_eigs=False)[source]
- Return a copy of - self, explicitly detaching and cloning all tensor values in the circuit (which are element and loop values). Eigenvalues/vectors are either discarded or detached and cloned based on the value of- save_eigs.- Parameters:
- save_eigs – Whether to retain the eigenvalues/vectors in the copied version of the circuit. 
- Return type:
- Self
- Returns:
- Deepcopy of self. 
 
 - picklecopy()[source]
- Helper function which returns a shallow copy of - selfwith- ._toggle_fullcopy = False. Use for pickling circuit to save memory.- Return type:
- Self
- Returns:
- Copy of self with - ._toggle_fullcopy = False.
 
 - description(tp=None, _test=False)[source]
- Print out Hamiltonian and a listing of the modes (whether they are harmonic or charge modes with the frequency for each harmonic mode), Hamiltonian parameters, and external flux values. Constructs a symbolic Hamiltonian, which is stored in .descrip_vars[‘H’]. - Parameters:
- tp ( - Optional[- str]) – If- Noneprints out the output as Latex if SQcircuit is running in a Jupyter notebook and as text if SQcircuit is running in Python terminal. If- tpis- "ltx", the output is in Latex format, and if- tpis- "txt"the output is in text format.
- _test ( - bool) – if True, return the entire description as string text (use only for testing the function).
 
- Return type:
- Optional[- str]
- Returns:
- The text of the description as a string, if - _testis- True.
 
 - loop_description(_test=False)[source]
- Print out the external flux distribution over inductive elements. - Parameters:
- _test ( - bool) – if True, return the entire description as string text. (use only for testing the function)
- Return type:
- Optional[- str]
- Returns:
- The text of the external flux distribution, if - _testis- True.
 
 - set_trunc_nums(nums)[source]
- Set the truncation numbers for each mode. - Parameters:
- nums ( - List[- int]) – A list that contains the truncation numbers for each mode. Harmonic modes with truncation number N are 0, 1 , …, (N-1), and charge modes with truncation number N are -(N-1), …, 0, …, (N-1).
- Return type:
- None
 
 - set_charge_offset(mode, ng)[source]
- Set the charge offset for each charge mode. - Parameters:
- mode ( - int) – An integer that specifies the charge mode. To see, which mode is a charge mode, one can use- description()method.
- ng ( - float) – The charge offset.
 
- Return type:
- None
 
 - set_charge_noise(mode, A)[source]
- Set the charge noise for each charge mode. - Parameters:
- mode ( - int) – An integer that specifies the charge mode. To see which mode is a charge mode, we can use- description()method.
- A ( - float) – The charge noise.
 
- Return type:
- None
 
 - charge_op(mode, basis='FC')[source]
- Return charge operator for specific mode in the Fock/Charge basis or the eigenbasis. - Parameters:
- mode ( - int) – Integer that specifies the mode number.
- basis ( - str) – String that specifies the basis. It can be either- 'FC'for original Fock/Charge basis or- 'eig'for eigenbasis.
 
- Return type:
- Qobj
- Returns:
- The charge operato for the - i``th mode in the basis specified by ``basis.
 
 - op(typ, keywords)[source]
- Get a saved circuit operator of type - typ, specified by keywords given in the- keywordsdict, as a backpropagatable- Tensorobject when- .get_optim_mode()is- True. Currently supports the following operators:- 'sin_half'
 - Parameters:
- typ ( - str) – Type of saved operator.
- keywords ( - Dict) – Dictionary specifying which operator of type- typto return.
 
- Return type:
- Union[- Qobj,- Tensor]
- Returns:
- The typ operator of the circuit specified by - keywords.
 
 - diag(n_eig)[source]
- Diagonalize the Hamiltonian of the circuit and return the eigenfrequencies and eigenvectors of the circuit up to specified number of eigenvalues. - Parameters:
- n_eig ( - int) – Number of eigenvalues to output. The lower- n_eig, the faster- SQcircuitfinds the eigenvalues.
- Return type:
- Tuple[- Union[- ndarray,- Tensor],- List[- Union[- Qobj,- Tensor]]]
- Returns:
- efreqs:
- ndarray of eigenfrequencies in frequency unit of SQcircuit (gigahertz by default). 
- evecs:
- List of eigenvectors in qutip.Qobj or Tensor format, depending on optimization mode. 
 
 
 - truncate_circuit(K, heuristic=False)[source]
- Set truncation numbers of circuit to - k=ceil(K^{1/n})for all modes, where- nis the number of modes in the circuit. If- heuristicis true, then the truncation number for each harmonic mode is weighted by setting- k_i = k * prod(omega_j^(1/n))/omega_iAll charge modes are left with truncation number- Kas above.- Parameters:
- K ( - int) – Total truncation number
- heuristic – Whether to use a heurstic to set harmonic mode truncations 
 
- Return type:
- List[- int]
- Returns:
- trunc_nums:
- List of truncation numbers for each mode of circuit 
 
 
 - check_convergence(eig_vec_idx=1, t=10, threshold=1e-05)[source]
- Check whether the diagonalization of the circuit has converged. - Parameters:
- eig_vec_idx – Index of eigenvector to use to test convergence. 
- t – Number of entries of eigenvector to use to test convergence. 
- threshold – Cutoff for convergence. 
 
- Returns:
- convergence_succeeded:
- Truthy value of whether the circuit converged 
- epsilon:
- Calculated value for convergence test 
 
 
 - coord_transform(var_type)[source]
- Return the transformation of the coordinates as ndarray for each type of variables, either charge or flux. - Parameters:
- var_type ( - str) – The type of the variables that can be either- "charge"or- "flux".
- Return type:
- ndarray
- Returns:
- Matrix giving coordinate transformation for - var_typecoordinates.
 
 - hamiltonian()[source]
- Returns the transformed hamiltonian of the circuit as - qutip.Qobjformat.- Return type:
- Qobj
- Returns:
- Circuit Hamiltonian. 
 
 - eig_phase_coord(k, grid)[source]
- Return the phase coordinate representations of the eigenvectors as ndarray. - Parameters:
- k ( - int) – The eigenvector index. For example, we set it to 0 for the ground state and 1 for the first excited state.
- grid ( - Sequence[- ndarray]) – A list that contains the range of values of phase φ for which we want to evaluate the wavefunction.
 
- Return type:
- ndarray
- Returns:
- Phase coordinate representation of the - k``th eigenvector over the values of φ provided in ``grid.
 
 - coupling_op(ctype, nodes)[source]
- Return the capacitive or inductive coupling operator related to the specified nodes. The output has the - qutip.Qobjformat.- Parameters:
- ctype ( - str) – Coupling type which is either- "capacitive"or- "inductive".
- nodes ( - Tuple[- int,- int]) – A tuple of circuit nodes to which we want to couple.
 
- Return type:
- Qobj
- Returns:
- Coupling operator of type - ctypebetween nodes in- nodes.
 
 - matrix_elements(ctype, nodes, states)[source]
- Return the matrix element of two eigenstates for either capacitive or inductive coupling. - Parameters:
- ctype ( - str) – Coupling type which is either- "capacitive"or- "inductive".
- nodes ( - Tuple[- int,- int]) – A tuple of circuit nodes to which we want to couple.
- states ( - Tuple[- int,- int]) – A tuple of indices of eigenstates for which we want to calculate the matrix element.
 
- Return type:
- float
- Returns:
- Matrix element between eigenstates in - statesfor coupling operator of type- ctypebetween nodes in- nodes.
 
 - dec_rate(dec_type, states, total=True)[source]
- Return the decoherence rate in [1/s] between each two eigenstates for different types of depolarization and dephasing. - Parameters:
- dec_type ( - str) – decoherence type that can be:- "capacitive"for capacitive loss;- "inductive"for inductive loss; “quasiparticle” for quasiparticle loss;- "charge"for charge noise,- "flux"for flux noise; and- "cc"for critical current noise.
- states ( - Tuple[- int,- int]) – A tuple of eigenstate indices, for which we want to calculate the decoherence rate. For example, for- states=(0, 1), we calculate the decoherence rate between the ground state and the first excited state.
- total ( - bool) – if False return a decoherence rate associated with a transition from state m to state n for- states=(m, n). if True return a decoherence rate associated with both m to n and n to m transitions.
 
- Return type:
- float
- Returns:
- Decoherence/dephasing rate between - statesspecified by- dec_type.
 
 - get_partial_omega(el, m, subtract_ground=True, _B_idx=None)[source]
- Return the gradient of the eigen angular frequency with respect to elements or loop as - qutip.Qobjformat.- Parameters:
- el ( - Union[- Capacitor,- Inductor,- Junction,- Loop]) – Element of a circuit that can be either- Capacitor,- Inductor,- Junction, or- Loop.
- m ( - int) – Integer specifies the eigenvalue. for example- m=0specifies the ground state and- m=1specifies the first excited state.
- subtract_ground ( - bool) – If- True, it subtracts the ground state frequency from the desired frequency.
- _B_idx ( - Optional[- int]) – Optional integer to indicate which row of the B matrix (per-element external flux distribution) to use. This specifies which JJ of the circuit to consider specifically (ex. for critical current noise calculation).
 
- Return type:
- float
- Returns:
- Partial derivative of eigenfrequency - mwith respect to- el, in units of angular frequency.
 
 - get_partial_vec(el, m, epsilon=1e-12)[source]
- Return the gradient of the eigenvectors with respect to elements or loop as - qutip.Qobjformat.- Parameters:
- el ( - Union[- Element,- Loop]) – Element of a circuit that can be either- Capacitor,- Inductor,- Junction, or- Loop.
- m ( - int) – Integer specifies the eigenvalue. for example- m=0specifies the ground state and- m=1specifies the first excited state.
 
- Return type:
- Qobj
- Returns:
- Partial derivative of the - m``th eigenvector, with respect to ``el.