Diagonalization
The diag subpackage provides eigensolvers, diversity selection, and
projected Hamiltonian construction utilities.
Eigensolvers
eigensolver — Backward-compatible re-exports
This module re-exports symbols that were split into
eigenvalue and
davidson for backward compatibility.
- qvartools.diag.eigen.eigensolver.solve_generalized_eigenvalue(H, S, k=1, which='SA', use_gpu=False, davidson_threshold=500)[source]
Solve the generalized eigenvalue problem Hv = ESv.
Dispatches to the appropriate backend depending on matrix format, size, and whether GPU acceleration is requested.
Solver selection (from highest to lowest priority):
GPU — if use_gpu is
Trueand CuPy is available.Sparse — if H or S is a SciPy sparse matrix.
Davidson — if the matrix dimension
nexceeds davidson_threshold (iterative, much faster than dense for large systems when only a few eigenvalues are needed).Dense — full eigendecomposition via
scipy.linalg.eigh.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).S (
np.ndarrayorscipy.sparse.spmatrix) – Overlap matrix, shape(n, n). Must be positive-definite; useregularize_overlap_matrix()if necessary.k (
int, optional) – Number of lowest eigenvalues / eigenvectors to compute (default1).which (
str, optional) – Eigenvalue selection criterion:"SA"for smallest algebraic (default).use_gpu (
bool, optional) – IfTrue, attempt to use CuPy for GPU-accelerated computation. Falls back to CPU if CuPy is unavailable.davidson_threshold (
int, optional) – Minimum matrix dimension above which the Davidson iterative solver is preferred over full dense decomposition (default500). Ignored when the matrix is sparse.
- Returns:
eigenvalues (
np.ndarray) – The lowestkeigenvalues, shape(k,), sorted ascending.eigenvectors (
np.ndarray) – Corresponding eigenvectors, shape(n, k).
- Raises:
ValueError – If
HandShave incompatible shapes ork < 1.RuntimeError – If the eigensolver fails to converge.
- Return type:
Examples
>>> H = np.diag([1.0, 2.0, 3.0]) >>> S = np.eye(3) >>> vals, vecs = solve_generalized_eigenvalue(H, S, k=2) >>> vals array([1., 2.])
- qvartools.diag.eigen.eigensolver.compute_ground_state_energy(H, use_gpu=False)[source]
Compute the ground-state energy (lowest eigenvalue) of a Hamiltonian.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration via CuPy.
- Returns:
The ground-state energy.
- Return type:
Examples
>>> H = np.diag([3.0, 1.0, 2.0]) >>> compute_ground_state_energy(H) 1.0
- qvartools.diag.eigen.eigensolver.analyze_spectrum(H, k=6, use_gpu=False)[source]
Compute spectral statistics of a Hamiltonian.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).k (
int, optional) – Number of lowest eigenvalues to compute (default6).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration.
- Returns:
Dictionary with keys:
"eigenvalues": np.ndarray – the lowestkeigenvalues."gaps": np.ndarray – consecutive eigenvalue gaps, shape(k-1,)."first_excited_gap": float – gap between ground and first excited state."ground_state_energy": float – the lowest eigenvalue.
- Return type:
Examples
>>> H = np.diag([1.0, 3.0, 7.0]) >>> info = analyze_spectrum(H, k=3) >>> info["first_excited_gap"] 2.0
- qvartools.diag.eigen.eigensolver.regularize_overlap_matrix(S, threshold=1e-06, use_gpu=False)[source]
Regularize an overlap matrix to ensure positive-definiteness.
Diagonalizes
S, clamps eigenvalues belowthresholdtothreshold, and reconstructs a well-conditioned overlap matrix.- Parameters:
S (
np.ndarrayorscipy.sparse.spmatrix) – Overlap matrix, shape(n, n).threshold (
float, optional) – Minimum allowed eigenvalue (default1e-6).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration via CuPy.
- Returns:
Regularized overlap matrix in CSR format.
- Return type:
Examples
>>> S = np.array([[1.0, 0.99], [0.99, 1.0]]) >>> S_reg = regularize_overlap_matrix(S, threshold=0.1)
- class qvartools.diag.eigen.eigensolver.DavidsonSolver(max_iterations=100, tolerance=1e-08, max_subspace_size=20)[source]
Bases:
objectIterative Davidson eigensolver for large sparse Hermitian matrices.
The Davidson algorithm extends the power method by building a subspace of approximate eigenvectors and using a preconditioner (the diagonal of the matrix) to accelerate convergence.
- Parameters:
Examples
>>> from scipy.sparse import random as sparse_random >>> H = sparse_random(100, 100, density=0.1, format="csr") >>> H = H + H.T # symmetrize >>> solver = DavidsonSolver(max_iterations=200) >>> eigenvalues, eigenvectors = solver.solve(H, k=2)
- solve(H, k=1)[source]
Compute the lowest
keigenvalues and eigenvectors ofH.- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hermitian matrix, shape(n, n).k (
int, optional) – Number of lowest eigenpairs to compute (default1).
- Returns:
eigenvalues (
np.ndarray) – Lowestkeigenvalues, shape(k,), sorted ascending.eigenvectors (
np.ndarray) – Corresponding eigenvectors, shape(n, k).
- Raises:
ValueError – If
kexceeds the matrix dimension.RuntimeError – If the solver does not converge within
max_iterations.
- Return type:
- class qvartools.diag.eigen.davidson.DavidsonSolver(max_iterations=100, tolerance=1e-08, max_subspace_size=20)[source]
Bases:
objectIterative Davidson eigensolver for large sparse Hermitian matrices.
The Davidson algorithm extends the power method by building a subspace of approximate eigenvectors and using a preconditioner (the diagonal of the matrix) to accelerate convergence.
- Parameters:
Examples
>>> from scipy.sparse import random as sparse_random >>> H = sparse_random(100, 100, density=0.1, format="csr") >>> H = H + H.T # symmetrize >>> solver = DavidsonSolver(max_iterations=200) >>> eigenvalues, eigenvectors = solver.solve(H, k=2)
- solve(H, k=1)[source]
Compute the lowest
keigenvalues and eigenvectors ofH.- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hermitian matrix, shape(n, n).k (
int, optional) – Number of lowest eigenpairs to compute (default1).
- Returns:
eigenvalues (
np.ndarray) – Lowestkeigenvalues, shape(k,), sorted ascending.eigenvectors (
np.ndarray) – Corresponding eigenvectors, shape(n, k).
- Raises:
ValueError – If
kexceeds the matrix dimension.RuntimeError – If the solver does not converge within
max_iterations.
- Return type:
eigenvalue — Eigenvalue problem solvers and overlap regularization
Provides functions for solving standard and generalized eigenvalue problems arising from projected Hamiltonian construction. Supports both CPU (SciPy) and GPU (CuPy) backends.
Functions
- solve_generalized_eigenvalue
Solve the generalized eigenvalue problem Hv = ESv.
- compute_ground_state_energy
Extract the ground-state energy from a Hamiltonian matrix.
- analyze_spectrum
Compute eigenvalues, gaps, and spectral statistics.
- regularize_overlap_matrix
Regularize an overlap matrix to ensure positive-definiteness.
- qvartools.diag.eigen.eigenvalue.solve_generalized_eigenvalue(H, S, k=1, which='SA', use_gpu=False, davidson_threshold=500)[source]
Solve the generalized eigenvalue problem Hv = ESv.
Dispatches to the appropriate backend depending on matrix format, size, and whether GPU acceleration is requested.
Solver selection (from highest to lowest priority):
GPU — if use_gpu is
Trueand CuPy is available.Sparse — if H or S is a SciPy sparse matrix.
Davidson — if the matrix dimension
nexceeds davidson_threshold (iterative, much faster than dense for large systems when only a few eigenvalues are needed).Dense — full eigendecomposition via
scipy.linalg.eigh.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).S (
np.ndarrayorscipy.sparse.spmatrix) – Overlap matrix, shape(n, n). Must be positive-definite; useregularize_overlap_matrix()if necessary.k (
int, optional) – Number of lowest eigenvalues / eigenvectors to compute (default1).which (
str, optional) – Eigenvalue selection criterion:"SA"for smallest algebraic (default).use_gpu (
bool, optional) – IfTrue, attempt to use CuPy for GPU-accelerated computation. Falls back to CPU if CuPy is unavailable.davidson_threshold (
int, optional) – Minimum matrix dimension above which the Davidson iterative solver is preferred over full dense decomposition (default500). Ignored when the matrix is sparse.
- Returns:
eigenvalues (
np.ndarray) – The lowestkeigenvalues, shape(k,), sorted ascending.eigenvectors (
np.ndarray) – Corresponding eigenvectors, shape(n, k).
- Raises:
ValueError – If
HandShave incompatible shapes ork < 1.RuntimeError – If the eigensolver fails to converge.
- Return type:
Examples
>>> H = np.diag([1.0, 2.0, 3.0]) >>> S = np.eye(3) >>> vals, vecs = solve_generalized_eigenvalue(H, S, k=2) >>> vals array([1., 2.])
- qvartools.diag.eigen.eigenvalue.compute_ground_state_energy(H, use_gpu=False)[source]
Compute the ground-state energy (lowest eigenvalue) of a Hamiltonian.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration via CuPy.
- Returns:
The ground-state energy.
- Return type:
Examples
>>> H = np.diag([3.0, 1.0, 2.0]) >>> compute_ground_state_energy(H) 1.0
- qvartools.diag.eigen.eigenvalue.analyze_spectrum(H, k=6, use_gpu=False)[source]
Compute spectral statistics of a Hamiltonian.
- Parameters:
H (
np.ndarrayorscipy.sparse.spmatrix) – Hamiltonian matrix, shape(n, n).k (
int, optional) – Number of lowest eigenvalues to compute (default6).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration.
- Returns:
Dictionary with keys:
"eigenvalues": np.ndarray – the lowestkeigenvalues."gaps": np.ndarray – consecutive eigenvalue gaps, shape(k-1,)."first_excited_gap": float – gap between ground and first excited state."ground_state_energy": float – the lowest eigenvalue.
- Return type:
Examples
>>> H = np.diag([1.0, 3.0, 7.0]) >>> info = analyze_spectrum(H, k=3) >>> info["first_excited_gap"] 2.0
- qvartools.diag.eigen.eigenvalue.regularize_overlap_matrix(S, threshold=1e-06, use_gpu=False)[source]
Regularize an overlap matrix to ensure positive-definiteness.
Diagonalizes
S, clamps eigenvalues belowthresholdtothreshold, and reconstructs a well-conditioned overlap matrix.- Parameters:
S (
np.ndarrayorscipy.sparse.spmatrix) – Overlap matrix, shape(n, n).threshold (
float, optional) – Minimum allowed eigenvalue (default1e-6).use_gpu (
bool, optional) – IfTrue, attempt GPU acceleration via CuPy.
- Returns:
Regularized overlap matrix in CSR format.
- Return type:
Examples
>>> S = np.array([[1.0, 0.99], [0.99, 1.0]]) >>> S_reg = regularize_overlap_matrix(S, threshold=0.1)
Projected Hamiltonian
- class qvartools.diag.eigen.projected_hamiltonian.ProjectedHamiltonianBuilder(hamiltonian, config=None)[source]
Bases:
objectBuild the projected Hamiltonian H_ij = <x_i|H|x_j> in a sampled basis.
Uses the Hamiltonian’s
get_connections()method for efficient sparse construction: for each basis state |x_j>, the connected states and their matrix elements are retrieved, and a hash-based lookup determines which connections lie within the sampled basis.- Parameters:
hamiltonian (
Hamiltonian) – The Hamiltonian operator (must implementdiagonal_elementandget_connections).config (
ProjectedHamiltonianConfigorNone, optional) – Construction hyperparameters. IfNone, defaults are used.
Examples
>>> from qvartools.hamiltonians.spin import HeisenbergHamiltonian >>> ham = HeisenbergHamiltonian(num_sites=6, J=1.0) >>> builder = ProjectedHamiltonianBuilder(ham) >>> basis = torch.tensor([[1,0,1,0,1,0], [0,1,0,1,0,1]]) >>> H_proj = builder.build(basis) >>> H_proj.shape (2, 2)
- build(basis_states)[source]
Construct the projected Hamiltonian matrix.
Uses vectorised batch diagonal computation and hash-based connection matching for efficiency.
- Parameters:
basis_states (
torch.Tensor) – Sampled basis configurations, shape(n_basis, n_sites). Each row is a computational-basis state.- Returns:
Projected Hamiltonian of shape
(n_basis, n_basis).- Return type:
- Raises:
ValueError – If
basis_stateshas fewer than 1 row or the wrong number of columns.
Selection
- class qvartools.diag.selection.diversity_selection.DiversitySelector(config, reference, n_orbitals)[source]
Bases:
objectSelect a diverse, representative subset of configurations.
Configurations are bucketed by excitation rank relative to a reference state (e.g. Hartree-Fock). Within each bucket, the top-weighted configurations are greedily selected while enforcing a minimum pairwise Hamming distance. An optional DPP kernel can be applied for maximal coverage.
- Parameters:
config (
DiversityConfig) – Selection hyperparameters.reference (
torch.Tensor) – Reference configuration (e.g. Hartree-Fock ground state), shape(n_orbitals,)with binary entries.n_orbitals (
int) – Number of orbitals (must matchreference.shape[0]).
- Raises:
ValueError – If
referencelength does not matchn_orbitals.
Examples
>>> cfg = DiversityConfig(max_configs=100) >>> ref = torch.tensor([1, 1, 1, 0, 0, 0]) >>> selector = DiversitySelector(cfg, ref, n_orbitals=6) >>> selected, stats = selector.select(all_configs, weights)
- select(configs, weights=None)[source]
Select a diverse subset of configurations.
- Parameters:
configs (
torch.Tensor) – Pool of candidate configurations, shape(n_pool, n_orbitals)with binary entries.weights (
torch.TensororNone, optional) – NQS importance weights for each configuration, shape(n_pool,). IfNone, uniform weights are used.
- Returns:
selected (
torch.Tensor) – Selected configurations, shape(n_selected, n_orbitals).stats (
dict) – Selection statistics with keys:"total_pool": int – size of the input pool."n_selected": int – number of selected configurations."rank_counts": dict – count per excitation rank in pool."rank_selected": dict – count per rank in selection."rank_quotas": dict – allocated quota per rank.
- Return type:
- class qvartools.diag.selection.diversity_selection.DiversityConfig(max_configs=1000, rank_0_fraction=0.05, rank_1_fraction=0.15, rank_2_fraction=0.4, rank_3_fraction=0.25, rank_4_plus_fraction=0.15, min_hamming_distance=2, use_dpp_selection=False, dpp_kernel_scale=0.1)[source]
Bases:
objectHyperparameters for diversity-aware configuration selection.
The selector partitions configurations by excitation rank relative to a reference state and allocates per-rank quotas according to the given fractions. Within each rank bucket, configurations are selected by NQS importance weight and filtered to enforce a minimum pairwise Hamming distance.
- Parameters:
max_configs (
int) – Maximum number of configurations to select.rank_0_fraction (
float) – Fraction ofmax_configsallocated to rank-0 (reference) configs.rank_1_fraction (
float) – Fraction allocated to single-excitation configurations.rank_2_fraction (
float) – Fraction allocated to double-excitation configurations.rank_3_fraction (
float) – Fraction allocated to triple-excitation configurations.rank_4_plus_fraction (
float) – Fraction allocated to quadruple-and-higher excitations.min_hamming_distance (
int) – Minimum pairwise Hamming distance between any two selected configs.use_dpp_selection (
bool) – IfTrue, use a determinantal point process kernel for maximal diversity within each rank bucket (more expensive).dpp_kernel_scale (
float) – Length-scale parameter for the DPP Gaussian kernel.
Examples
>>> cfg = DiversityConfig(max_configs=500, min_hamming_distance=3) >>> cfg.max_configs 500
excitation_rank — Excitation-rank and Hamming-distance utilities
Provides functions for computing excitation ranks, Hamming distances, and bit-packing of binary configuration tensors.
Functions
- compute_excitation_rank
Count the number of bit flips between a configuration and a reference.
- compute_hamming_distance
Hamming distance between two binary configuration tensors.
- bitpack_configs
Pack binary configuration tensors into int64 words for fast Hamming.
- bitpacked_hamming
Vectorized Hamming distance via popcount on bit-packed tensors.
- qvartools.diag.selection.excitation_rank.compute_excitation_rank(config, reference)[source]
Count the number of bit flips from a reference configuration.
The excitation rank is the Hamming distance between the configuration and the reference, which for a fermionic system corresponds to the number of single-particle excitations (orbital substitutions).
- Parameters:
config (
torch.Tensor) – Binary configuration vector, shape(n_orbitals,).reference (
torch.Tensor) – Reference configuration (e.g. Hartree-Fock), shape(n_orbitals,).
- Returns:
Number of differing bits (excitation rank).
- Return type:
- Raises:
ValueError – If
configandreferencehave different shapes.
Examples
>>> ref = torch.tensor([1, 1, 0, 0]) >>> cfg = torch.tensor([1, 0, 1, 0]) >>> compute_excitation_rank(cfg, ref) 2
- qvartools.diag.selection.excitation_rank.compute_hamming_distance(config1, config2)[source]
Compute the Hamming distance between two binary configurations.
- Parameters:
config1 (
torch.Tensor) – First binary configuration, shape(n_orbitals,).config2 (
torch.Tensor) – Second binary configuration, shape(n_orbitals,).
- Returns:
Number of positions where the two configurations differ.
- Return type:
- Raises:
ValueError – If
config1andconfig2have different shapes.
Examples
>>> a = torch.tensor([1, 0, 1, 0]) >>> b = torch.tensor([1, 1, 0, 0]) >>> compute_hamming_distance(a, b) 2
- qvartools.diag.selection.excitation_rank.bitpack_configs(configs)[source]
Pack binary configurations into int64 words for O(1) Hamming distance.
Each configuration of
n_orbitalsbinary values is packed intoceil(n_orbitals / 63)int64 words (63 bits per word to avoid sign bit issues).- Parameters:
configs (
torch.Tensor) – Binary configuration matrix, shape(n_configs, n_orbitals)with entries in{0, 1}.- Returns:
Bit-packed tensor, shape
(n_configs, n_words), dtypeint64.- Return type:
Examples
>>> cfgs = torch.tensor([[1, 0, 1, 1], [0, 1, 0, 0]]) >>> packed = bitpack_configs(cfgs) >>> packed.shape torch.Size([2, 1])
- qvartools.diag.selection.excitation_rank.bitpacked_hamming(packed, idx_a, idx_b)[source]
Compute vectorized Hamming distances via popcount on bit-packed configs.
- Parameters:
packed (
torch.Tensor) – Bit-packed configurations, shape(n_configs, n_words), dtypeint64, as returned bybitpack_configs().idx_a (
torch.Tensor) – Index tensor for the first set of configurations, shape(n_pairs,).idx_b (
torch.Tensor) – Index tensor for the second set of configurations, shape(n_pairs,).
- Returns:
Hamming distances, shape
(n_pairs,), dtypeint64.- Return type:
Examples
>>> cfgs = torch.tensor([[1, 0, 1, 1], [0, 1, 0, 0], [1, 1, 0, 0]]) >>> packed = bitpack_configs(cfgs) >>> idx_a = torch.tensor([0, 0]) >>> idx_b = torch.tensor([1, 2]) >>> bitpacked_hamming(packed, idx_a, idx_b) tensor([4, 2])
Bitstring and basis-set utility functions for SKQD postprocessing.
Provides helpers for converting between bitstring and integer representations, accumulating measurement results across Krylov steps, filtering low-probability states, and computing basis-set overlap metrics.
- qvartools.diag.selection.bitstring.bitstring_to_int(bitstring)[source]
Convert a binary bitstring to its integer representation.
- Parameters:
bitstring (
str) – String of'0'and'1'characters (e.g."0110").- Returns:
Integer value of the bitstring.
- Return type:
Examples
>>> bitstring_to_int("0110") 6
- qvartools.diag.selection.bitstring.int_to_bitstring(value, num_bits)[source]
Convert an integer to a zero-padded binary bitstring.
- Parameters:
- Returns:
Binary string of length
num_bits.- Return type:
Examples
>>> int_to_bitstring(6, 4) '0110'
- qvartools.diag.selection.bitstring.get_basis_states_as_array(measurement_results, num_qubits)[source]
Convert measurement results to an array of unique basis-state integers.
- Parameters:
- Returns:
Sorted array of unique basis-state integers, dtype
int64.- Return type:
np.ndarray
Examples
>>> results = {"01": 5, "10": 3, "01": 2} >>> get_basis_states_as_array(results, num_qubits=2) array([1, 2])
- qvartools.diag.selection.bitstring.calculate_cumulative_results(all_measurement_results)[source]
Calculate cumulative measurement results across Krylov steps.
For step k, the cumulative results include all unique bitstrings from steps 0, 1, …, k with their total counts.
- Parameters:
all_measurement_results (
listofdict) – One measurement dictionary per Krylov step, each mapping bitstring to count.- Returns:
Cumulative measurement dictionaries. Entry k contains the union of all bitstrings observed in steps 0 through k.
- Return type:
Examples
>>> step0 = {"00": 3, "01": 2} >>> step1 = {"01": 1, "10": 4} >>> cumulative = calculate_cumulative_results([step0, step1]) >>> cumulative[1] {'00': 3, '01': 3, '10': 4}
- qvartools.diag.selection.bitstring.filter_high_probability_states(measurement_results, threshold=0.0, max_states=None)[source]
Filter measurement results to keep only high-probability states.
- Parameters:
measurement_results (
dictofstrtoint) – Mapping from bitstring to occurrence count.threshold (
float, optional) – Minimum empirical probability for a state to be retained (default0.0, i.e. keep all).max_states (
intorNone, optional) – If notNone, keep at most this many states (the highest-count states are preferred).
- Returns:
Filtered measurement dictionary.
- Return type:
Examples
>>> counts = {"00": 90, "01": 5, "10": 3, "11": 2} >>> filter_high_probability_states(counts, threshold=0.05) {'00': 90, '01': 5}
- qvartools.diag.selection.bitstring.compute_basis_overlap(basis1, basis2)[source]
Compute the overlap fraction between two basis sets.
Returns the fraction of states in
basis1that are also present inbasis2.- Parameters:
basis1 (
np.ndarray) – First basis as an array of state integers.basis2 (
np.ndarray) – Second basis as an array of state integers.
- Returns:
Overlap fraction in the range
[0, 1]. Returns0.0whenbasis1is empty.- Return type:
Examples
>>> b1 = np.array([0, 1, 2, 3]) >>> b2 = np.array([2, 3, 4, 5]) >>> compute_basis_overlap(b1, b2) 0.5
- qvartools.diag.selection.bitstring.estimate_ground_state_sparsity(ground_state, threshold=1e-06)[source]
Estimate sparsity metrics of a ground-state wavefunction.
- Parameters:
ground_state (
np.ndarray) – Ground-state wavefunction vector, shape(dim,).threshold (
float, optional) – Minimum probability|c_i|^2for a component to be counted as significant (default1e-6).
- Returns:
Sparsity metrics with keys:
"n_significant": int – number of components abovethreshold."sparsity_ratio": float – fraction of Hilbert space with significant weight."concentration": float – total probability weight in the top 10 %% of components."total_dimension": int – total Hilbert-space dimension.
- Return type:
Examples
>>> psi = np.array([0.9, 0.1, 0.0, 0.0]) >>> stats = estimate_ground_state_sparsity(psi, threshold=1e-4) >>> stats["n_significant"] 2
- qvartools.diag.selection.bitstring.merge_basis_sets(*bases)[source]
Merge multiple basis sets into one sorted, unique set.
- Parameters:
*bases (
np.ndarray) – Variable number of 1-D arrays of basis-state integers.- Returns:
Sorted array of unique basis states, dtype
int64.- Return type:
np.ndarray
Examples
>>> b1 = np.array([0, 1, 2]) >>> b2 = np.array([2, 3]) >>> merge_basis_sets(b1, b2) array([0, 1, 2, 3])