Modular Framework

The types in CtrlVQE are designed to be as flexible as possible, while remaining performant. We did this by designing the DeviceType to represent the complete Hamiltonian,

\[\hat H(\vec x,t) = \hat H_0(\vec x) + \sum_i \hat V_i (\vec x,t)\]

There is an interface to follow drift and drive terms are handled correctly and efficiently in time evolution, but it is completely and entirely your choice how to implement them. But prioritizing flexibility and performance has a tendency to sacrifice readability and reusability.

This package attempts to find a middle ground, by splitting each bit of the Hamiltonian equation into chunks. That is, there is one type (AlgebraType) which defines the space on which each operator $\hat H$, $\hat H_0$, $\hat V_i$ acts. There is another (DriftType) defining $\hat H_0$. There is another (DriveType) defining $\hat V_i$.

And to help keep the drives as modular as possible, there is a MapperType to map "device" parameters to individual drift and drive. A more complete equation - which is itself less readable, but it makes the code so much easier to work with - would be

\[\hat H(\vec x,t) = \hat H_0(h(\vec x)) + \sum_i \hat V_i(f_i(\vec x), t)\]

For example, let's say you're curious how a single three-qubit interaction term in the static Hamiltonian of a transmon device would alter the minimal evolution time. This is a tiny change, but going from Basics, it would seem to require implementing a whole new device type, and to do that correctly, you'd have to be sure to implement the dozen or so methods required for new device types. You'll presumably just copy and paste from an existing implementation, which is generally considered bad practice and has a tendency to produce unexpected bugs (hence the desire for code reusability). And you'll generate a file with hundreds of lines of code, a real pain to organize, document, and debug (hence the desire for code readability).

Going from here, you'd simply implement a new drift type, with two methods, and be done.

Core Interface

CtrlVQE.ModularFrameworkModule
ModularFramework

Provides a framework for designing devices and energy functions in a less monolithic way. Different behaviors are delegated to different types.

Behavior for devices are divided into:

  • AlgebraType
  • DriftType
  • DriveType
  • ParameterMap

These can be implemented by the user and then combined into a pre-existing concrete type, usually LocalDevice.

Behavior for energy functions are divided into:

  • ReferenceType
  • MeasurementType

These can be implemented by the user and then combined with an integration, device, and evolution using the pre-existing concrete type Energy.

Each type has a validate method defined, with the following signatures. (See the function's docstring for details.)

validate(::AlgebraType)
validate(::DriftType; algebra::AlgebraType)
validate(::DriveType; algebra::AlgebraType, grid::IntegrationType, t)
validate(::ParameterMap; device::DeviceType)

validate(::ReferenceType; device::DeviceType)
validate(::MeasurementType; grid::IntegrationType, device::DeviceType, t)

Types may optionally implement the Prototypes interface. Since these are often inter-dependent (e.g. a prototypical DipoleDrive may depend on the DriftType), prototypes defined here are implemented in a special file, Prototypes.jl, and users defining their own types are advised to do something similar for the sake of writing concise doctests.

source
CtrlVQE.ModularFramework.AlgebraTypeType
AlgebraType{m,n}

Defines a Hilbert space and an operator basis.

Delegated the following methods:

  • Devices.nlevels
  • Devices.nqubits
  • Devices.noperators
  • Devices.localalgebra

Type Parameters

  • m: the number of levels in each qubit
  • n: the number of qubits

Implementation

Subtypes A must implement the following methods:

  • Devices.noperators(::Type{A}): how many unique operators define an algebra on one qubit.
  • Devices.localalgebra(::A; result): a 4d array ā where ā[:,:,σ,q] is the σ'th algebraic operator on qubit q, in the bare basis.
source
CtrlVQE.ModularFramework.DriftTypeType
DriftType{A}

Defines the static Hamiltonian $\hat H_0$.

Delegated the following methods:

  • Devices.qubithamiltonian
  • Devices.staticcoupling

Implementation

Subtypes H must implement the following methods:

  • Devices.qubithamiltonian(::H, ā, q::Int; result): the static components of the device Hamiltonian local to qubit q.
  • Devices.staticcoupling(::H, ā; result): the static components of the device Hamiltonian nonlocal to any one qubit.

The arg is a 4d array with all algebra operators, like the one returned by Devices.localalgebra (but not necessarily in a local basis).

source
CtrlVQE.ModularFramework.DriveTypeType
DriveType{A}

Defines a drive term $\hat V_i$.

Delegated the following methods:

  • Devices.ngrades
  • Devices.driveoperator
  • Devices.gradeoperator
  • Devices.gradient

Implementation

Subtypes V must implement the Parameters interface, and the following methods:

  • ngrades(::Type{<:V}): the number of distinct gradient operators for this drive type.
  • driveoperator(::V, ā, t::Real; result): the distinct drive operator for this channel at time t.
  • gradeoperator(::V, ā, j::Int, t::Real; result): the distinct gradient operator indexed by j at time t.
  • gradient(::V, grid::Integrations.IntegrationType, ϕ; result): the gradient vector for each variational parameter in the device.

Note that the drive index i is omitted from the interface for driveoperator, and the grade index j is with respect to just this channel. That means j will only ever take values between 1 and ngrades(V). This is true of both the j in the signature of gradeoperator and the column indices of ϕ in the signature of gradient.

source
CtrlVQE.ModularFramework.LocalDriveType
LocalDrive{F,A}

Extension of DriveType for $\hat V_i$ which act on a single qubit.

Delegated drivequbit.

Implementation

Subtypes D must implement the DriveType interface, and the following method:

  • Devices.drivequbit(::D): index of the qubit on which this drive is applied.

Note that the drive index i is omitted from the interface for drivequbit.

source
CtrlVQE.ModularFramework.MeasurementTypeType
MeasurementType{B,O}

A protocol to measure scalars from a statevector.

Type Parameters

  • B<:Bases.BasisType: the basis that a measurement takes place in.
  • O<:Operators.OperatorType: the frame rotation that a measurement takes place in. Specifically, identifies the static operator in the interacture picture, so STATIC identifies the dressed frame and IDENTITY identifies the lab frame.

Implementation

TODO: Rather different now. Be sure to mention overriding initbasis/initframe if needed.

Subtypes M must implement the following methods:

  • measure(::M, device, ψ): measures a state given in the initbasis.
  • observables(::M, device; result): constructs the observables in the initbasis.
  • CostFunctions.nobservables(::Type{M}): the number of distinct observables involved in a measurement.
  • Devices.gradient(::M, device, grid, ϕ, ψ; result): calculates the gradient of device parameters, given gradient signals.

In Devices.gradient, ϕ[:,j,k] contains the jth gradient signal $ϕ_j(t)$ evaluated at each point in grid for observable k, while ψ contains the wavefunction itself, evolved to the end of the grid and rotated into the measurement basis. For the most part, implementing types will simply delegate to device, whose gradient method expects the 2d array ϕ[:,k]. When M consists of more than one observable, it will need to decide how to combine the resulting gradients into one.

While initbasis and initframe have default implementations for singleton basis and operator types, these methods will need to be overridden if your measurement is ever needed with more complex types.

source
CtrlVQE.ModularFramework.ParameterMapType
ParameterMap

Enumerates the protocol for mapping device parameters to each drive.

Delegated Parameters.names.

Implementation

Subtypes P must implement the following methods:

  • Parameters.names(::P, device): constructs a list of human-readable names for each parameter in the device, which is a device implementing the Parameters interface. Certain implementations of ParameterMap may have additional requirements.
  • sync!(::P, device): mutate the internal parameters of a device to match those of its drives.
  • map_values(::P, device, i::Int; result): computes drive parameters from device parameters.
  • map_gradients(::P, device, i::Int; result): computes drive gradients with respect to device parameters.
source
CtrlVQE.ModularFramework.ReferenceTypeType
ReferenceType{B}

A protocol to prepare an initial statevector.

Type Parameters

  • B<:Bases.BasisType: the basis that a reference was defined in

Implementation

Subtypes R must implement the following methods:

  • prepare(::R, device; result): constructs the initial statevector in the initbasis.

While initbasis has a default implementation for singleton basis types, this method will need to be overridden if your reference is ever needed with more complex types.

source
CtrlVQE.ModularFramework.initframeMethod
initframe(measurement)

Fetch the frame that this measurement takes place in.

Specifically, this function fetches the static operator in the interacture picture, so STATIC identifies the dressed frame and IDENTITY identifies the lab frame.

source
CtrlVQE.ModularFramework.map_gradientsFunction
map_gradients(pmap::ParameterMap, device::ModularDevice, i::Int; result)

Compute gradients for parameters in a drive term, with respect to each device parameter.

Parameters

  • pmap: the ParameterMap defining the family of functions to map parameters.
  • device: the device, giving the device parameters (via Parameters.values). Certain implementations of ParameterMap may have additional requirements.
  • i: identifies which function in the family (i.e. indexes the drive).

Returns

A matrix g, such that g[k,j] is $∂_{x_k} y_j$.

source
CtrlVQE.ModularFramework.map_valuesFunction
map_values(pmap::ParameterMap, device, i::Int; result)

Compute the parameters for a drive term, as a function of all device parameters.

Parameters

  • pmap: the ParameterMap defining the family of functions to map parameters.
  • device: the device, giving the device parameters (via Parameters.values). Certain implementations of ParameterMap may have additional requirements.
  • i: identifies which function in the family (i.e. indexes the drive).
source
CtrlVQE.ModularFramework.measureFunction
measure(measurement, device, ψ)
measure(measurement, device, ψ, t)
measure(measurement, device, basis, ψ)
measure(measurement, device, basis, ψ, t)

Measure the state ψ, provided in the given basis, at time t.

Parameters

  • measurement::MeasurementType: the measurement protocol.
  • device::DeviceType: the device being measured.
  • basis::BasisType: the basis that ψ is represented in. Defaults to initbasis when omitted.
  • ψ::AbstractVector: the state to measure.
  • t::Real: the time at which ψ is being measured, i.e. the time over which to evolve the frame operator. When t is omitted, the frame rotation is skipped.
source
CtrlVQE.ModularFramework.observablesFunction
observables(measurement, device; result=nothing)
observables(measurement, device, t; result=nothing)
observables(measurement, device, basis; result=nothing)
observables(measurement, device, basis, t; result=nothing)

Constructs the Hermitian observables involved in this measurement.

Parameters

  • measurement::MeasurementType: the measurement protocol.
  • device::DeviceType: the device being measured.
  • basis::BasisType: the basis the observables are represented in. Defaults to initbasis when omitted.
  • t::Real: the time at which the measurement takes place, i.e. the time over which to evolve the frame operator. When t is omitted, the frame rotation is skipped.

The returned result is also written to the result kwarg when provided.

Returns

A 3darray indexed such that result[:,:,k] is the kth matrix, the same format required by the parameter Ō in Evolutions.gradientsignal. The size of the third dimension must be equal to nobservables(measurement).

source
CtrlVQE.ModularFramework.prepareFunction
prepare(reference, device; result=nothing)
prepare(reference, device, basis; result=nothing)

Prepare the state and represent it as a statevector in the given basis.

Parameters

  • reference::ReferenceType: the state preparation protocol.
  • device::DeviceType: the device for which the state is being prepared.
  • basis::BasisType: the basis to represent the reference state in (defaults to initbasis when omitted).

The returned result is also written to the result kwarg when provided.

source
CtrlVQE.ModularFramework.sync!Function
sync!(device)
sync!(pmap::ParameterMap, device)

Mutate the internal parameters of a device to match those of its drives.

The one-parameter signature requires device to have a property pmap, which is the parameter map that will be used for dispatch. Both signatures require device to have the property x which is a vector of all the parameters of device. Certain implementations of ParameterMap may have additional requirements.

This function may resize x.

source
CtrlVQE.ModularFramework.LocalDevices.LocalDeviceType
LocalDevice(algebra, drift, drives, pmap, x)

A totally modular device, with some restrictions (see below).

Parameters

  • algebra::AlgebraType: the algebra, defining the Hilbert space.

  • drift::DriftType: the drift term, describing the static Hamiltonian.

  • drives::Vector{<:DriveType} a vector of drive terms, describing the time-dependent Hamiltonian.

  • pmap::ParameterMap: the parameter map, defining the relationship between device and drive parameters.

  • x::Vector{<:AbstractFloat}: the initial device parameters

    LocalDevice(F, algebra, drift, drives, pmap)

The typical constructor for a LocalDevice.

Parameters

  • F: the float type for this device (typically Float64).
  • algebra::AlgebraType: the algebra, defining the Hilbert space.
  • drift::DriftType: the drift term, describing the static Hamiltonian.
  • drives::Vector{<:DriveType} a vector of drive terms, describing the time-dependent Hamiltonian.
  • pmap::ParameterMap: the parameter map, defining the relationship between device and drive parameters.

The device parameters are initialized to be consistent with drives (though the degree of consistency is defined by pmap).

Restrictions

  • Drift can't be parametric

    • I could program around this, but it would add an unnecessary amount of code for a feature I may never use.
    • In principle, another type could relax this restriction if:
      • compatible drift types implement the Parameters interface.
      • compatible parameter map types implement suitable versions of map_values and map_gradients
      • the bind! method clears relevant methods from the Memoization.jl cache.
      • all loops in this file with parameter manipulation include a step for the drift.
  • Only works with LocalDrives

    • This restriction is necessary only due to lack of multiple inheritance in Julia. In other words, I have to subtype LocallyDrivenDevice in order to get the fast evolution, and that means I have to restrict the drives to those compatible with it. But none of the code here depends on it, other than the parts implementing the interface for LocallyDrivenDevice.
    • In principle, another type could relax this restriction if:
      • you copy and paste the whole file.
      • switch out LocallDrivenDevice for LocalDevice
      • switch out LocalDrive with DriveType.
      • delete the parts with drivequbit and gradequbit.
source
CtrlVQE.ModularFramework.Energies.EnergyType
Energy(evolution, device, grid, reference, measurement)

A totally modular energy function.

Parameters

  • evolution::EvolutionType: the time integration algorithm
  • device::DeviceType: the device that gets evolved in time
  • grid::IntegrationType: the time grid to integrate over
  • reference::ReferenceType: the reference state to initialize the device to
  • measurement::MeasurementType: the protocol to extract energies from evolved states
source
CtrlVQE.ModularFramework.DrivePenalties.DrivePenaltyType
DrivePenalty(device, penalties::AbstractVector{<:CostFunctionType})

The sum of several penalty functions acting on a device's drives.

In order for this type to be usable, each penalty in penalties must act on parameters corresponding to the drives contained in device.drives, and device parameters must be mapped onto drives according to device.pmap. This interface is suitable for the Device types in ModularFramework.

This struct records the last values computed (value or gradient) for each cost function in its values and gradients fields.

Beware that the components field is a vector with abstract eltype. This is a relatively mild form of type instability; certain compiled code could hypothetically be sub-optimal. But it does not appear to cause any noticeable disadvantage; see Examples/CompositeCostFunctions.

source

Device Modules

Algebras

CtrlVQE.ModularFramework.TruncatedBosonicAlgebras.TruncatedBosonicAlgebraType
TruncatedBosonicAlgebra{m,n}()

An algebra of n distinguishable bosonic modes, each represented with m levels.

This algebra is useful for transmon architectures.

julia> using CtrlVQE.ModularFramework;

julia> A = TruncatedBosonicAlgebra{3,2};

julia> validate(A());

julia> nlevels(A)
3
julia> nqubits(A)
2
julia> noperators(A)
1
source
CtrlVQE.ModularFramework.PauliAlgebras.PauliAlgebraType
PauliAlgebra{n}()

An algebra of n spinors, represented with the Pauli spin matrices.

julia> using CtrlVQE.ModularFramework;

julia> A = PauliAlgebra{4};

julia> validate(A());

julia> nlevels(A)
2
julia> nqubits(A)
4
julia> noperators(A)
3
source

Drifts

CtrlVQE.ModularFramework.TransmonDrifts.TransmonDriftType
TransmonDrift(ω, δ, g, quples)

A static Hamiltonian for architectures of n transmons with fixed-couplings.

``\hat H0 = \sumq ωq \hat a^\daggerq \hat a_q

  • \sumq δq \hat a^\daggerq \hat a^\daggerq \hat aq \hat aq
  • \sum{⟨pq⟩} g{pq} (\hat a^\daggerp \hat aq + \hat a^\daggerq \hat ap)``

Parameters

  • ω: a vector of n resonance frequencies
  • δ: a vector of n anharmonicities
  • g: a vector of coupling strengths
  • quples: a vector of Quple identifying the qubit pairs associated with each g[i]
julia> using CtrlVQE.ModularFramework;

julia> A = TruncatedBosonicAlgebra{3,2};

julia> drift = TransmonDrift{A}([4.82, 4.84], [0.30, 0.30], [0.02], [Quple(1,2)]);

julia> validate(drift; algebra=A());

julia> ā0 = Devices.localalgebra(A());

julia> Devices.qubithamiltonian(drift, ā0, 1)
3×3 Matrix{ComplexF64}:
 0.0+0.0im   0.0+0.0im   0.0+0.0im
 0.0+0.0im  4.82+0.0im   0.0+0.0im
 0.0+0.0im   0.0+0.0im  9.34+0.0im

julia> Devices.qubithamiltonian(drift, ā0, 2)
3×3 Matrix{ComplexF64}:
 0.0+0.0im   0.0+0.0im   0.0+0.0im
 0.0+0.0im  4.84+0.0im   0.0+0.0im
 0.0+0.0im   0.0+0.0im  9.38+0.0im
source

Drives

CtrlVQE.ModularFramework.DipoleDrives.DipoleDriveType
DipoleDrive(q, ω, Ω, Δ)

A drive term representing interaction of a bosonic mode with an electric dipole, in the rotating wave approximation.

$\hat V = Ω(t) e^{i [Δ(t)+ω] t} \hat a_q + {\rm h.t.}$

Parameters

  • q: the integer index identifying which qubit this drive applies to
  • ω: the resonance frequency of the qubit
  • Ω: a SignalType representing the complex drive strength over time
  • Δ: a SignalType representing the detuning over time
julia> grid = TemporalLattice(20.0, 400);

julia> Ω = Constant(2.0+1.0im);

julia> Δ = Constant(0.0);

julia> using CtrlVQE.ModularFramework;

julia> A = TruncatedBosonicAlgebra{3,2};

julia> drive = DipoleDrive{A}(1, 4.82, Ω, Δ);

julia> validate(drive; algebra=A(), grid=grid, t=10.0);

julia> ā0 = Devices.localalgebra(A());

julia> Devices.driveoperator(drive, ā0, 10.0)
3×3 Matrix{ComplexF64}:
        0.0+0.0im      -0.0693931-2.23499im         0.0+0.0im
 -0.0693931+2.23499im         0.0+0.0im      -0.0981367-3.16075im
        0.0+0.0im      -0.0981367+3.16075im         0.0+0.0im
source

Parameter Mappers

CtrlVQE.ModularFramework.DisjointMappers.DisjointMapperType
DISJOINT

Useful for when parameters for each drive are completely independent of one another.

In this case, all the parameters for drive 1 are followed by all the parameters for drive 2, and so on.

julia> using CtrlVQE.ModularFramework;

julia> pmap = DISJOINT;

julia> device = Prototype(LocalDevice{Float64}; n=3, pmap=pmap);

julia> validate(pmap; device=device);

julia> map_values(pmap, device, 1)
2-element Vector{Float64}:
 0.0
 0.0

julia> Parameters.count(device)
6
source
CtrlVQE.ModularFramework.LinearMappers.LinearMapperType
LinearMapper(encoding::Vector{F}, size::Vector{Int})
LinearMapper(A::Array{F,3})

Associate each device parameter with the weight of a basis vector for each drive.

The implemetation of this type is designed around the assumption that each drive has the same number of parameters.

The linear map A is a 3d array such that y[i] = A[:,:,i] * x, where y[i] is the parameter vector for drive i and x is the parameter vector for the whole device.

A is represented internally as a data vector encoding and shape vector size, so that it can be resized in-place.

julia> using CtrlVQE.ModularFramework;

julia> pmap = LinearMapper(ones(2,4,3));

julia> device = Prototype(LocalDevice{Float64}; n=3, pmap=pmap);

julia> validate(pmap; device=device);

julia> map_values(pmap, device, 1)
2-element Vector{Float64}:
 0.0
 0.0

julia> Parameters.count(device)
4
source
CtrlVQE.ModularFramework.LinearMappers.addbasisvector!Method
addbasisvector!(pmap, a)

Update a LinearMapper to add in a new device parameter.

You'll pretty much always need to sync! a device after calling this.

Parameters

  • pmap: the LinearMapper object to mutate.
  • a::AbstractMatrix: a[:,i] contains the mapping from the new device parameter to the parameters of drive i.
source

Energy Modules

References

CtrlVQE.ModularFramework.KetReferences.KetReferenceType
KetReference(basis::B, ket)

Represents a basis state of the given basis B.

Parameters

  • basis::Bases.BasisType: the basis this ket is defined in.
  • ket: a vector of integers representing the ket. For example, [0,1,1] represents state |011⟩, the fourth basis state. Note that index 1 of the ket is the most significant bit.
julia> using CtrlVQE.ModularFramework;

julia> reference = KetReference(BARE, [0,1]);

julia> device = Prototype(LocalDevice{Float64}; n=2);

julia> validate(reference; device=device);

julia> prepare(reference, device)
4-element Vector{ComplexF64}:
 0.0 + 0.0im
 1.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
source
CtrlVQE.ModularFramework.DenseReferences.DenseReferenceType
DenseReference(basis::B, statevector; m, n)

Represents an arbitrary statevector of the given basis B.

Parameters

  • basis::Bases.BasisType: the basis where statevector represents the reference.
  • statevector::AbstractVector: a dense statevector.

The kwargs m and n specify the number of levels and number of qubits for which the statevector is defined. The reference state may be prepared on a device with the same number of qubits, and at least as many levels. Typically, only one of m and n are provided, and the other is inferred from the length of statevector. If neither is provided, m defaults to 2. If both are provided, an error will be thrown if they are not consistent with the length of statevector.

julia> using CtrlVQE.ModularFramework;

julia> reference = DenseReference(BARE, [0,1,0,0]);

julia> device = Prototype(LocalDevice{Float64}; n=2);

julia> validate(reference; device=device);

julia> prepare(reference, device)
4-element Vector{ComplexF64}:
 0.0 + 0.0im
 1.0 + 0.0im
 0.0 + 0.0im
 0.0 + 0.0im
source

Measurements

CtrlVQE.ModularFramework.PauliMeasurements.PauliMeasurementType
PauliMeasurement(basis, frame, X, Z, c)

Represents a bare measurement of a linear combination of Paulis, without any logical projection prior to frame rotation.

Paulis are encoded symplectically, as two length-n bitstrings labeled x and z. 1s in x indicate a Pauli X appears on that bit; 1s in z indicate the same for Z. When a bit has both an X and a Z, it is interpreted as having a Y. (You could think of it as Z*X = iY, except the implementation here ignores phase.)

Parameters

  • basis: the BasisType identifying the basis observable is written in.
  • frame: the OperatorType identifying the frame where measurements are conducted.
  • X: a vector of integers whose bitstrings give each x.
  • Z: a vector of integers whose bitstrings give each z.
  • c: a vector of floats giving the coefficients for each pauli word.
source
CtrlVQE.ModularFramework.PauliMeasurements.PauliMeasurementMethod
PauliMeasurement(basis, frame, observable; eps=1e-10)

A constructor accepting a dense matrix observable.

The Pauli coefficients are computed once from (a not very efficient implementation of) the Hilbert-Schmidt inner product.

Parameters

  • basis: the BasisType identifying the basis observable is written in.
  • frame: the OperatorType identifying the frame where measurements are conducted.
  • observable: the dense matrix observable.
  • eps: smallest coefficient to keep before discarding as negligible.
source
CtrlVQE.ModularFramework.PauliMeasurements.PauliMeasurementMethod
PauliMeasurement(basis, frame; paulis...)

A constructor accepting coefficients for each Pauli word.

Parameters

  • basis: the BasisType identifying the basis observable is written in.
  • frame: the OperatorType identifying the frame where measurements are conducted.

Keyword Arguments

  • Each key is a consistent-length string consisting only of "I", "X", "Y", "Z".
  • Each value is the float coefficient for that Pauli word.
julia> using CtrlVQE.ModularFramework;

julia> I = [1 0; 0 1]; X = [0 1; 1 0]; Z = [1 0; 0 -1];

julia> observable = (0.5 .* kron(X,Z)) .+ (1.0 .* kron(I, Z))
4×4 Matrix{Float64}:
 1.0   0.0  0.5   0.0
 0.0  -1.0  0.0  -0.5
 0.5   0.0  1.0   0.0
 0.0  -0.5  0.0  -1.0

julia> measurement = PauliMeasurement(BARE, STATIC; XZ=0.5, IZ=1.0);

julia> device = Prototype(LocalDevice{Float64}; n=2);

julia> validate(measurement; device=device);

julia> Ō = observables(measurement, device);

julia> size(Ō)
(4, 4, 1)

julia> observable ≈ reshape(Ō, (4,4))
true
source
CtrlVQE.ModularFramework.DenseMeasurements.DenseMeasurementType
DenseMeasurement(basis, frame, observable)

Represents a bare measurement of a matrix observable in a given basis and frame, without any logical projection prior to frame rotation.

Parameters

  • observable: the dense matrix observable.
  • basis: the BasisType identifying the basis observable is written in.
  • frame: the OperatorType identifying the frame where measurements are conducted.

Specifically, the operator identified by frame is the one we rotate by for duration t to move from the "lab" frame to the interaction picture. Use IDENTITY for the lab frame itslef, and STATIC for the dressed frame.

julia> using CtrlVQE.ModularFramework;

julia> I = [1 0; 0 1]; X = [0 1; 1 0]; Z = [1 0; 0 -1];

julia> observable = (0.5 .* kron(X,Z)) .+ (1.0 .* kron(I, Z))
4×4 Matrix{Float64}:
 1.0   0.0  0.5   0.0
 0.0  -1.0  0.0  -0.5
 0.5   0.0  1.0   0.0
 0.0  -0.5  0.0  -1.0

julia> measurement = DenseMeasurement(BARE, STATIC, observable);

julia> device = Prototype(LocalDevice{Float64}; n=2);

julia> validate(measurement; device=device);

julia> Ō = observables(measurement, device);

julia> size(Ō)
(4, 4, 1)

julia> observable ≈ reshape(Ō, (4,4))
true
source