You signed in with another tab or window. Reload to refresh your session.You signed out in another tab or window. Reload to refresh your session.You switched accounts on another tab or window. Reload to refresh your session.Dismiss alert
# Affine decompositions with multi-indices and additional parameters
2
+
3
+
In this example, we want to explore the capabilities of the central [`AffineDecomposition`](@ref) type.
4
+
To advance the previous examples where we covered the magnetization — a very simple observable that consists of only one affine term and a parameter-independent coefficient — we now turn to observables where the indices are multi-indices ``r = (r_1, \dots, r_d)`` and the coefficients can depend on additional parameters ``p``, aside from the ``\bm{\mu}`` parameter points that are present in the Hamiltonian:
So let us continue the first example where we have generated an ``L=6`` XXZ surrogate `basis` with a reduced Hamiltonian `h` using an exact diagonalization solver.
67
+
Now the task is to implement the double-sum in ``\mathcal{S}``, as well as the ``k``-dependency in the coefficients.
68
+
The double-sum can be encoded by putting all ``S^z_r S^z_{r'}`` combinations into a ``L \times L`` matrix:
One feature of the structure factor that also shows up in many other affine decompositions with double-sums is that the term indices commute, i.e. ``O_{r,r'} = O_{r',r}``.
83
+
In that case, only the upper triangular matrix has to be computed since ``B^\dagger O_{r,r'} B = B^\dagger O_{r',r} B`` are the same in the compressed affine decomposition.
84
+
So let's create the [`AffineDecomposition`](@ref) and compress, exploiting this symmetry using the `symmetric_terms` keyword argument:
In the online evaluation of the structure factor, we then need to define some wavevector values and compute the structure factor at each of them.
92
+
With the `grid_online` from before, this reads:
93
+
94
+
```@example multi_ad; continued = true
95
+
wavevectors = [0.0, π/4, π/2, π]
96
+
sf = [zeros(size(grid_online)) for _ in 1:length(wavevectors)]
97
+
for (idx, μ) in pairs(grid_online)
98
+
_, φ_rb = solve(h, basis.metric, μ, fulldiag)
99
+
for (i, k) in enumerate(wavevectors)
100
+
sf[i][idx] = sum(eachcol(φ_rb)) do u
101
+
real(dot(u, sfspin(k), u))
102
+
end / size(φ_rb, 2)
103
+
end
104
+
end
105
+
```
106
+
107
+
Here we again see the convenience of measuring observables in the online stage;
108
+
adding more wavevector values does not significantly increase the computational cost, since it corresponds to a mere reevaluation of the coefficient functions and small vector-matrix products.
109
+
Finally, let us see how the structure factor behaves for the different wavevector values:
It can be nicely seen that the spin structure factor indicates the ferromagnetic phase at ``k=0`` and then moves through the magnetization plateaus until it reaches the antiferromagnetic plateau at ``k=\pi``.
Copy file name to clipboardExpand all lines: docs/src/examples/xxz_dmrg.md
+25-25Lines changed: 25 additions & 25 deletions
Display the source diff
Display the rich diff
Original file line number
Diff line number
Diff line change
@@ -1,7 +1,7 @@
1
1
# Greedy basis assembly using DMRG
2
2
3
-
As a follow-up example, we now want to showcase how to compute an `RBasis` by means of the [density matrix renormalization group](https://tensornetwork.org/mps/algorithms/dmrg/) (DMRG).
4
-
To that end, we utilize the `ITensors.jl` library which, among other things, efficiently implements DMRG.
3
+
As a follow-up example, we now want to showcase how to compute a reduced basis by means of the [density matrix renormalization group](https://tensornetwork.org/mps/algorithms/dmrg/) (DMRG).
4
+
To that end, we utilize the ITensors.jl package which, among other things, efficiently implements DMRG.
5
5
We will see that, while we need to adjust the way we set up the model Hamiltonian as well as our solver, most steps stay the same.
6
6
Again, we treat the one-dimensional ``S=1/2`` XXZ model from the previous example.
7
7
@@ -19,7 +19,7 @@ using ReducedBasis
19
19
```
20
20
21
21
now featuring `ITensors`.
22
-
To build the Hamiltonian terms as MPOs, we make use of the [`OpSum()`](https://itensor.github.io/ITensors.jl/stable/tutorials/DMRG.html) object that automatically produces an MPO from a string of operators.
22
+
To build the Hamiltonian terms as MPOs, we make use of the `ITensors.OpSum()` object that automatically produces an MPO from a string of operators.
23
23
The affine MPO terms are then stored in an [`AffineDecomposition`](@ref) as [`ApproxMPO`](@ref)s which also include possible truncation keyword arguments:
24
24
25
25
```@example xxz_dmrg; continued = true
@@ -35,29 +35,27 @@ function xxz_chain(sites::IndexSet; kwargs...)
35
35
end
36
36
magn_term += "Sz", length(sites) # Add last magnetization term
So let us instantiate such an MPO Hamiltonian where we also specify a singular value `cutoff`, which is passed to the [`ApproxMPO`](@ref) objects:
48
46
49
47
```@example xxz_dmrg; continued = true
50
-
L = 12
48
+
L = 18
51
49
sites = siteinds("S=1/2", L)
52
50
H = xxz_chain(sites; cutoff=1e-9)
53
51
```
54
52
55
-
We now chose a bigger system size, since the tensor format allows for efficient low rank approximations (hence the `cutoff`) that buy us a substantial performance advantage when going to larger systems.
53
+
Notice that we can now choose a bigger system size (which is still very small here), since the tensor format allows for efficient low rank approximations (hence the `cutoff`) that buy us a substantial performance advantage when going to larger systems.
56
54
57
55
## Using the [`DMRG`](@ref) solver for basis assembly
58
56
59
57
Having created our Hamiltonian in MPO format, we now need a solver that is able to compute ground states from MPOs.
60
-
The corresponding ground state will also be obtained in a tensor format, namely as *matrix product states* (MPS).
58
+
The corresponding ground state will also be obtained in a tensor format, namely as a *matrix product state* (MPS).
61
59
This is achieved by `ITensors.dmrg` which is wrapped in the [`DMRG`](@ref) solver type:
While the implemented DMRG solver is capable of solving degenerate ground state, we here opt for non-degenerate settings (i.e. `n_states=1` and `tol_degeneracy=0.0`), since one encounters a ``L+1``-fold degeneracy on the parameter domain, where the degenerate DMRG solver can produce instable results for larger ``L``.
67
+
While the implemented DMRG solver is capable of solving degenerate ground states, we here opt for non-degenerate settings, i.e. `n_states=1` and `tol_degeneracy=0.0`.
68
+
(We do this due to a ``L+1``-fold degeneracy on the parameter domain, where the degenerate DMRG solver can produce instable results for larger ``L``.)
70
69
71
70
As discussed in the last example, we need a way to orthogonalize the reduced basis.
72
-
Due to the MPS format that the snapshots will have, we cannot use QR decompositions anymore and resort to a different method, [`EigenDecomposition`](@ref), featuring an eigenvalue decomposition of the snapshot overlap matrix:
71
+
Due to the MPS format that the snapshots will have, we cannot use QR decompositions anymore and resort to a different method, [`EigenDecomposition`](@ref), featuring an eigenvalue decomposition of the snapshot overlap matrix that can drop insignificant snapshots that fall below a cutoff:
Now with different types for `H`, the solver and the orthogonalizer, we call `assemble` using the `greedy` strategy and training grid from the last example:
84
+
Now with different types for the Hamiltonian, the solver and the orthogonalizer, we call `assemble` using the `greedy` strategy and training grid from the last example:
86
85
87
86
```@example xxz_dmrg; continued = true
88
-
basis, h, info = assemble(H, grid_train, greedy, dm, edcomp)
87
+
info = assemble(H, grid_train, greedy, dm, edcomp)
88
+
basis = info.basis; h = info.h_cache.h;
89
89
```
90
90
91
91
The returned `basis` now has snapshot vectors of `ITensors.MPS` type, which we have to keep in mind when we want to compress observables.
92
-
That is to say, the observables have to be constructed as `AffineDecompositions` with `ApproxMPO` terms as for the Hamiltonian.
92
+
That is to say, the observables have to be constructed as [`AffineDecomposition`](@ref)s with [`ApproxMPO`](@ref) terms as we did for the Hamiltonian.
93
93
Again, we want to compute the magnetization so that we can reuse the third term of `H`:
94
94
95
95
```@example xxz_dmrg; continued = true
96
96
M = AffineDecomposition([H.terms[3]], μ -> [2 / L])
And at that point, we can continue as before since we have arrived at the online phase where we only operate in the low-dimensional reduced basis space, agnostic of the previously used solver method.
116
-
We have to make sure, however, to choose matching degeneracy settings for the [`FullDiagonalization`](@ref) solver in the online phase, which we do via the matching constructor:
115
+
And at that point, we continue as before since we have arrived at the online phase where we only operate in the low-dimensional reduced basis space, agnostic of the snapshot solver method.
116
+
We have to make sure, however, to choose matching degeneracy settings for the [`FullDiagonalization`](@ref) solver in the online phase:
0 commit comments