Gross-Pitaevskii equation in one dimension

In this example we will use DFTK to solve the Gross-Pitaevskii equation, and use this opportunity to explore a few internals.

The model

The Gross-Pitaevskii equation (GPE) is a simple non-linear equation used to model bosonic systems in a mean-field approach. Denoting by ψ the effective one-particle bosonic wave function, the time-independent GPE reads in atomic units: $$ H ψ = \left(-\frac12 Δ + V + 2 C |ψ|^2\right) ψ = μ ψ \qquad \|ψ\|_{L^2} = 1 $$ where C provides the strength of the boson-boson coupling. It's in particular a favorite model of applied mathematicians because it has a structure simpler than but similar to that of DFT, and displays interesting behavior (especially in higher dimensions with magnetic fields, see Gross-Pitaevskii equation with magnetism).

We wish to model this equation in 1D using DFTK. First we set up the lattice. For a 1D case we supply two zero lattice vectors,

In [1]:
a = 10
lattice = a .* [[1 0 0.]; [0 0 0]; [0 0 0]];

which is special cased in DFTK to support 1D models.

For the potential term V we just pick a harmonic potential. The real-space grid is in [0,1) in fractional coordinates( see Lattices and lattice vectors), therefore:

In [2]:
pot(x) = (x - a/2)^2;

We setup each energy term in sequence: kinetic, potential and nonlinear term. For the non-linearity we use the PowerNonlinearity(C, α) term of DFTK. This object introduces an energy term C ∫ ρ(r)^α dr to the total energy functional, thus a potential term α C ρ^{α-1}. In our case we thus need the parameters

In [3]:
C = 1.0
α = 2;

... and with this build the model

In [4]:
using DFTK
using LinearAlgebra

n_electrons = 1  # Increase this for fun
terms = [Kinetic(),
         ExternalFromReal(r -> pot(r[1])),
         PowerNonlinearity(C, α),
]
model = Model(lattice; n_electrons=n_electrons, terms=terms,
              spin_polarization=:spinless);  # use "spinless electrons"

We discretize using a moderate Ecut (For 1D values up to 5000 are completely fine) and run a direct minimization algorithm:

In [5]:
basis = PlaneWaveBasis(model, Ecut=500, kgrid=(1, 1, 1))
scfres = direct_minimization(basis, tol=1e-8) # This is a constrained preconditioned LBFGS
scfres.energies
Iter     Function value   Gradient norm 
     0     1.803271e+02     1.687193e+02
 * time: 0.0005929470062255859
     1     1.601613e+02     1.408108e+02
 * time: 0.0019040107727050781
     2     1.260739e+02     1.559715e+02
 * time: 0.0030488967895507812
     3     5.005600e+01     1.034580e+02
 * time: 0.004617929458618164
     4     2.562358e+01     7.426721e+01
 * time: 0.005960941314697266
     5     7.775187e+00     4.265292e+00
 * time: 0.007016897201538086
     6     5.938181e+00     9.590303e+00
 * time: 0.007807016372680664
     7     5.303142e+00     2.436536e+01
 * time: 0.00855398178100586
     8     3.146115e+00     1.014774e+01
 * time: 0.009292840957641602
     9     2.180215e+00     2.588234e+00
 * time: 0.0101318359375
    10     1.714000e+00     2.342090e+00
 * time: 0.010967016220092773
    11     1.385824e+00     1.603570e+00
 * time: 0.011739969253540039
    12     1.262433e+00     1.010181e+00
 * time: 0.012336015701293945
    13     1.193454e+00     8.628967e-01
 * time: 0.012928962707519531
    14     1.168538e+00     9.568073e-01
 * time: 0.013535022735595703
    15     1.151669e+00     6.331508e-01
 * time: 0.014189958572387695
    16     1.148250e+00     5.083292e-01
 * time: 0.014854907989501953
    17     1.146282e+00     1.382512e-01
 * time: 0.015578031539916992
    18     1.146192e+00     1.878342e-01
 * time: 0.01607799530029297
    19     1.145136e+00     1.445251e-01
 * time: 0.016679048538208008
    20     1.144409e+00     1.156383e-01
 * time: 0.01729106903076172
    21     1.144244e+00     3.846202e-02
 * time: 0.017920970916748047
    22     1.144143e+00     2.496072e-02
 * time: 0.01853489875793457
    23     1.144075e+00     1.052256e-02
 * time: 0.01916193962097168
    24     1.144048e+00     9.322276e-03
 * time: 0.019819974899291992
    25     1.144041e+00     5.145441e-03
 * time: 0.020519018173217773
    26     1.144039e+00     2.955804e-03
 * time: 0.021149873733520508
    27     1.144038e+00     2.301091e-03
 * time: 0.02176189422607422
    28     1.144037e+00     4.981745e-03
 * time: 0.022357940673828125
    29     1.144037e+00     3.642429e-03
 * time: 0.022958993911743164
    30     1.144037e+00     1.534142e-03
 * time: 0.0235898494720459
    31     1.144037e+00     6.220448e-04
 * time: 0.024198055267333984
    32     1.144037e+00     3.882335e-04
 * time: 0.024888992309570312
    33     1.144037e+00     3.317794e-04
 * time: 0.025586843490600586
    34     1.144037e+00     1.604869e-04
 * time: 0.026239871978759766
    35     1.144037e+00     1.243124e-04
 * time: 0.0268399715423584
    36     1.144037e+00     1.524377e-04
 * time: 0.02746295928955078
    37     1.144037e+00     6.338290e-05
 * time: 0.02789592742919922
    38     1.144037e+00     3.728940e-05
 * time: 0.0285189151763916
    39     1.144037e+00     3.866383e-05
 * time: 0.029130935668945312
    40     1.144037e+00     1.915725e-05
 * time: 0.029855966567993164
    41     1.144037e+00     1.303086e-05
 * time: 0.030536890029907227
    42     1.144037e+00     4.742708e-06
 * time: 0.03122091293334961
    43     1.144037e+00     3.254390e-06
 * time: 0.03182697296142578
    44     1.144037e+00     2.711997e-06
 * time: 0.03243088722229004
    45     1.144037e+00     2.139196e-06
 * time: 0.03302288055419922
    46     1.144037e+00     1.146952e-06
 * time: 0.03362083435058594
    47     1.144037e+00     1.153891e-06
 * time: 0.034214019775390625
    48     1.144037e+00     8.805688e-07
 * time: 0.03488302230834961
    49     1.144037e+00     7.859604e-07
 * time: 0.03557300567626953
    50     1.144037e+00     6.630402e-07
 * time: 0.03619790077209473
    51     1.144037e+00     2.651445e-07
 * time: 0.036778926849365234
    52     1.144037e+00     1.375474e-07
 * time: 0.03736686706542969
Out[5]:
Energy breakdown:
    Kinetic             0.2682057 
    ExternalFromReal    0.4707475 
    PowerNonlinearity   0.4050836 

    total               1.144036852755 

Internals

We use the opportunity to explore some of DFTK internals.

Extract the converged density and the obtained wave function:

In [6]:
ρ = real(scfres.ρ)[:, 1, 1, 1]  # converged density, first spin component
ψ_fourier = scfres.ψ[1][:, 1];    # first kpoint, all G components, first eigenvector

Transform the wave function to real space and fix the phase:

In [7]:
ψ = G_to_r(basis, basis.kpoints[1], ψ_fourier)[:, 1, 1]
ψ /= (ψ[div(end, 2)] / abs(ψ[div(end, 2)]));

Check whether ψ is normalised:

In [8]:
x = a * vec(first.(DFTK.r_vectors(basis)))
N = length(x)
dx = a / N  # real-space grid spacing
@assert sum(abs2.(ψ)) * dx  1.0

The density is simply built from ψ:

In [9]:
norm(scfres.ρ - abs2.(ψ))
Out[9]:
1.1045183068376976e-15

We summarize the ground state in a nice plot:

In [10]:
using Plots

p = plot(x, real.(ψ), label="real(ψ)")
plot!(p, x, imag.(ψ), label="imag(ψ)")
plot!(p, x, ρ, label="ρ")
Out[10]:

The energy_hamiltonian function can be used to get the energy and effective Hamiltonian (derivative of the energy with respect to the density matrix) of a particular state (ψ, occupation). The density ρ associated to this state is precomputed and passed to the routine as an optimization.

In [11]:
E, ham = energy_hamiltonian(basis, scfres.ψ, scfres.occupation; ρ=scfres.ρ)
@assert E.total == scfres.energies.total

Now the Hamiltonian contains all the blocks corresponding to kpoints. Here, we just have one kpoint:

In [12]:
H = ham.blocks[1];

H can be used as a linear operator (efficiently using FFTs), or converted to a dense matrix:

In [13]:
ψ11 = scfres.ψ[1][:, 1] # first kpoint, first eigenvector
Hmat = Array(H) # This is now just a plain Julia matrix,
#                which we can compute and store in this simple 1D example
@assert norm(Hmat * ψ11 - H * ψ11) < 1e-10

Let's check that ψ11 is indeed an eigenstate:

In [14]:
norm(H * ψ11 - dot(ψ11, H * ψ11) * ψ11)
Out[14]:
2.0076092799647823e-7

Build a finite-differences version of the GPE operator H, as a sanity check:

In [15]:
A = Array(Tridiagonal(-ones(N - 1), 2ones(N), -ones(N - 1)))
A[1, end] = A[end, 1] = -1
K = A / dx^2 / 2
V = Diagonal(pot.(x) + C .* α .* (ρ.^(α-1)))
H_findiff = K + V;
maximum(abs.(H_findiff*ψ - (dot(ψ, H_findiff*ψ) / dot(ψ, ψ)) * ψ))
Out[15]:
0.00022341788149503677