Skip to content

Commit 7a0806e

Browse files
Improve SciMLOperator handling of sesolveProblem
1 parent 54ded38 commit 7a0806e

File tree

3 files changed

+8
-4
lines changed

3 files changed

+8
-4
lines changed

src/time_evolution/sesolve.jl

Lines changed: 5 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -36,6 +36,10 @@ function _generate_sesolve_kwargs(e_ops, progress_bar::Val{false}, tlist, kwargs
3636
return _generate_sesolve_kwargs_with_callback(tlist, kwargs)
3737
end
3838

39+
_sesolve_make_U_QobjEvo(H::QuantumObjectEvolution{<:MatrixOperator}) =
40+
QobjEvo(MatrixOperator(-1im * H.data.A), dims = H.dims, type = Operator)
41+
_sesolve_make_U_QobjEvo(H) = QobjEvo(H, -1im)
42+
3943
@doc raw"""
4044
sesolveProblem(
4145
H::Union{AbstractQuantumObject{DT1,OperatorQuantumObject},Tuple},
@@ -88,7 +92,7 @@ function sesolveProblem(
8892

8993
tlist = convert(Vector{_FType(ψ0)}, tlist) # Convert it to support GPUs and avoid type instabilities for OrdinaryDiffEq.jl
9094

91-
H_evo = QobjEvo(H, -1im) # pre-multiply by -i
95+
H_evo = _sesolve_make_U_QobjEvo(H) # Multiply by -i
9296
isoper(H_evo) || throw(ArgumentError("The Hamiltonian must be an Operator."))
9397
check_dims(H_evo, ψ0)
9498

src/time_evolution/time_evolution_dynamical.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -564,8 +564,8 @@ function _DSF_mcsolve_Affect!(integrator)
564564
@. e_ops0 = get_data(e_ops2)
565565
@. c_ops0 = get_data(c_ops2)
566566
H_nh = lmul!(convert(eltype(ψt), 0.5im), mapreduce(op -> op' * op, +, c_ops0))
567-
# By doing this, we are assuming that the system is time-independent and f is a ScaledOperator
568-
copyto!(integrator.f.f.L.A, H(op_l2, dsf_params).data - H_nh)
567+
# By doing this, we are assuming that the system is time-independent and f is a MatrixOperator
568+
copyto!(integrator.f.f.A, lmul!(-1im, H(op_l2, dsf_params).data - H_nh))
569569
return u_modified!(integrator, true)
570570
end
571571

test/core-test/time_evolution.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -98,7 +98,7 @@
9898
sol_mc_string = sprint((t, s) -> show(t, "text/plain", s), sol_mc)
9999
sol_sse_string = sprint((t, s) -> show(t, "text/plain", s), sol_sse)
100100
@test prob_me.f.f isa MatrixOperator
101-
@test prob_mc.f.f isa ScaledOperator # TODO: can be optimized as MatrixOperator
101+
@test prob_mc.f.f isa MatrixOperator
102102
@test sum(abs.(sol_mc.expect .- sol_me.expect)) / length(t_l) < 0.1
103103
@test sum(abs.(sol_mc2.expect .- sol_me.expect)) / length(t_l) < 0.1
104104
@test sum(abs.(vec(expect_mc_states_mean) .- vec(sol_me.expect))) / length(t_l) < 0.1

0 commit comments

Comments
 (0)