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.597445e+02     1.635877e+02
 * time: 0.0006971359252929688
     1     1.540207e+02     1.529545e+02
 * time: 0.002201080322265625
     2     1.168600e+02     1.706955e+02
 * time: 0.003958940505981445
     3     4.719435e+01     1.074110e+02
 * time: 0.005813121795654297
     4     1.990498e+01     6.233995e+01
 * time: 0.0075130462646484375
     5     1.446872e+01     5.470987e+01
 * time: 0.008764028549194336
     6     1.028559e+01     3.463054e+01
 * time: 0.009814977645874023
     7     9.239520e+00     1.018030e+01
 * time: 0.010827064514160156
     8     6.585533e+00     1.100819e+01
 * time: 0.011834144592285156
     9     3.948380e+00     1.142955e+01
 * time: 0.012946128845214844
    10     2.122278e+00     2.391355e+00
 * time: 0.01401209831237793
    11     1.416845e+00     1.553289e+00
 * time: 0.015066146850585938
    12     1.369622e+00     1.573589e+00
 * time: 0.01585698127746582
    13     1.250320e+00     7.940685e-01
 * time: 0.01666712760925293
    14     1.184392e+00     8.198353e-01
 * time: 0.017460107803344727
    15     1.158106e+00     4.092584e-01
 * time: 0.018252134323120117
    16     1.148319e+00     1.811436e-01
 * time: 0.019036054611206055
    17     1.145951e+00     1.127054e-01
 * time: 0.019865989685058594
    18     1.145023e+00     7.041936e-02
 * time: 0.02069711685180664
    19     1.144439e+00     4.518698e-02
 * time: 0.021502017974853516
    20     1.144231e+00     6.553943e-02
 * time: 0.022083044052124023
    21     1.144195e+00     6.369244e-02
 * time: 0.022686004638671875
    22     1.144130e+00     3.549207e-02
 * time: 0.02433609962463379
    23     1.144091e+00     1.907598e-02
 * time: 0.02518606185913086
    24     1.144058e+00     1.124082e-02
 * time: 0.0260009765625
    25     1.144040e+00     7.059753e-03
 * time: 0.026803970336914062
    26     1.144040e+00     4.578568e-03
 * time: 0.027383089065551758
    27     1.144038e+00     4.604402e-03
 * time: 0.028222084045410156
    28     1.144037e+00     2.416074e-03
 * time: 0.029050111770629883
    29     1.144037e+00     1.055905e-03
 * time: 0.02986311912536621
    30     1.144037e+00     6.336158e-04
 * time: 0.03067493438720703
    31     1.144037e+00     7.362609e-04
 * time: 0.03147411346435547
    32     1.144037e+00     2.860822e-04
 * time: 0.03205394744873047
    33     1.144037e+00     1.964097e-04
 * time: 0.03290200233459473
    34     1.144037e+00     1.304958e-04
 * time: 0.03381705284118652
    35     1.144037e+00     1.078917e-04
 * time: 0.0346379280090332
    36     1.144037e+00     4.127986e-05
 * time: 0.03550004959106445
    37     1.144037e+00     6.999419e-05
 * time: 0.036366939544677734
    38     1.144037e+00     3.934010e-05
 * time: 0.03724193572998047
    39     1.144037e+00     2.635852e-05
 * time: 0.03808093070983887
    40     1.144037e+00     1.961053e-05
 * time: 0.03891801834106445
    41     1.144037e+00     1.166020e-05
 * time: 0.03971505165100098
    42     1.144037e+00     7.174118e-06
 * time: 0.04052400588989258
    43     1.144037e+00     6.365032e-06
 * time: 0.0413510799407959
    44     1.144037e+00     4.842403e-06
 * time: 0.04215693473815918
    45     1.144037e+00     3.860060e-06
 * time: 0.04297900199890137
    46     1.144037e+00     3.154550e-06
 * time: 0.044181108474731445
    47     1.144037e+00     2.319926e-06
 * time: 0.0450289249420166
    48     1.144037e+00     1.045393e-06
 * time: 0.04588508605957031
    49     1.144037e+00     6.560508e-07
 * time: 0.04670095443725586
    50     1.144037e+00     2.975209e-07
 * time: 0.04751706123352051
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]:
8.338746851893258e-16

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]:
3.326664185386434e-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.00022348749912162638