6. Splitting the valence/conduction of silicon

In the previous tutorial, we have seen how to use the parallel_transport function to automate the construction of WFs for isolated manifold. In this tutorial, we will again use this technique, combined with disentanglement, to Wannierize the valence and conduction manifolds of silicon, separately.

Usually, Wannierization of valence or conduction manifold is difficult because there is no prior knowledge of the bonding or anti-bonding orbitals, except for simple system.

For the valence manifold, since it is isolated, we can

  • use random Gaussian as initial projections, then run many iterations of maxiaml localization,
  • use SCDM for the isolated manifold, then run maximal localization.

However, for the conduction manifold, it is entangled with higher energy bands, so the random Gaussian is not applicable. In such case, we can

  • use SCDM with specified $\mu$ and $\sigma$ for the conduction manifold, then run maximal localization.

However, the quality of Wannierization sensitively depends on the $\mu$ and $\sigma$, and also SCDM works in real space which is more memory consuming.

Now, our solution is

  • run the normal disentanglement on valence + conduction manifold, the initial guess can be the analytic spdf projectors, or using numerical atomic orbitals (NAOs) from pseudopotentials (for plane-wave code like QE, this automates the initial projections for the disentanglement).
  • separate the valence + conduction manifold into two, by diagonalizing the Wannier gauge Hamiltonian, the eigenvectors are unitary matrices, in which the upper (lower) part maps to the valence (conduction) manifold. Thus we have two separated submanifolds.
  • run the parallel_transport on the separated submanifolds, to construct the parallel transport gauge (PTG) for valence and conduction, respectively. We can further run optimal rotation and maximal localization to smoothen the gauge.

The disentanglement generates a smooth "total" manifold, helping us getting rid of higher energy states, and the PTG automates the gauge construction for the two separated submanifolds. In final, we have an automated approach to Wannierize the valence and/or conduction manifolds without human intervention.

Outline

  1. construct a Model, by reading the win, amn, mmn, and eig files
  2. split the Model by diagonalizing the Wannier gauge Hamiltonian
  3. run parallel_transport on the two Models to construct PTG for valence and conduction, respectively
  4. maxiaml localize to smoothen the PTG
  5. interpolate band structures
Tip

This is a HTML version of the tutorial, you can download corresponding

Preparation

Load the package

using Wannier
using WannierPlots

Model generation

We will use the read_w90 function to read the win, amn, mmn, and eig files, and construct a Model that abstracts the calculation

model = read_w90("si2")
lattice: Å
  a1:  0.00000  2.71527  2.71527
  a2:  2.71527  0.00000  2.71527
  a3:  2.71527  2.71527  0.00000

atoms: fractional
  Si:  0.00000  0.00000  0.00000
  Si:  0.25000  0.25000  0.25000

n_bands: 16
n_wann : 8
kgrid  : 4 4 4
n_kpts : 64
n_bvecs: 8

b-vectors:
         [bx, by, bz] / Å⁻¹                weight
  1       0.28925   -0.28925    0.28925    1.49401
  2      -0.28925   -0.28925   -0.28925    1.49401
  3      -0.28925    0.28925   -0.28925    1.49401
  4      -0.28925   -0.28925    0.28925    1.49401
  5       0.28925    0.28925   -0.28925    1.49401
  6      -0.28925    0.28925    0.28925    1.49401
  7       0.28925   -0.28925   -0.28925    1.49401
  8       0.28925    0.28925    0.28925    1.49401

Disentanglement

First let's disentangle the valence + conduction manifold,

U = disentangle(model);
[ Info: Initial spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1    -0.00000     0.00000     0.00000     1.33119
   2     0.00001     0.00011    -0.00000     1.92293
   3    -0.00000     0.00011     0.00011     1.92297
   4     0.00001    -0.00000     0.00011     1.92303
   5     1.35763     1.35763     1.35763     1.33119
   6     1.35762     1.35752     1.35763     1.92293
   7     1.35763     1.35752     1.35752     1.92297
   8     1.35762     1.35763     1.35752     1.92303
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =     9.69628
   Ω̃   =     4.50396
   ΩOD =     4.49427
   ΩD  =     0.00968
   Ω   =    14.20024


[ Info: Initial spread (with states freezed)
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.00000     0.00000     0.00000     1.79953
   2     0.00000    -0.00000    -0.00000     2.44042
   3     0.00000    -0.00000     0.00000     2.44042
   4    -0.00000    -0.00000    -0.00000     2.44042
   5     1.35763     1.35763     1.35763     1.79953
   6     1.35763     1.35763     1.35763     2.44042
   7     1.35763     1.35763     1.35763     2.44042
   8     1.35763     1.35763     1.35763     2.44042
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    12.82012
   Ω̃   =     5.42148
   ΩOD =     5.34421
   ΩD  =     0.07726
   Ω   =    18.24160


Iter     Function value   Gradient norm
     0     1.824160e+01     6.135873e-01
 * time: 0.0042569637298583984
     1     1.609863e+01     8.063143e-02
 * time: 0.03941798210144043
     2     1.571861e+01     4.479970e-02
 * time: 0.06448197364807129
     3     1.559753e+01     3.171032e-02
 * time: 0.10131692886352539
     4     1.551187e+01     1.252183e-02
 * time: 0.21619391441345215
     5     1.548086e+01     8.310836e-03
 * time: 0.25153398513793945
     6     1.547552e+01     4.120118e-03
 * time: 0.27658891677856445
     7     1.547423e+01     1.505832e-03
 * time: 0.31236886978149414
     8     1.547405e+01     5.299607e-04
 * time: 0.34804487228393555
     9     1.547403e+01     4.244924e-04
 * time: 0.38411903381347656
    10     1.547402e+01     2.528288e-04
 * time: 0.409639835357666
    11     1.547402e+01     8.185063e-05
 * time: 0.4459500312805176
    12     1.547402e+01     4.605264e-05
 * time: 0.5533130168914795
[ Info: Final spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.00000     0.00000     0.00000     1.44670
   2     0.00000     0.00000    -0.00000     2.09677
   3    -0.00000    -0.00000     0.00000     2.09677
   4    -0.00000    -0.00000    -0.00000     2.09677
   5     1.35763     1.35763     1.35763     1.44670
   6     1.35763     1.35763     1.35763     2.09677
   7     1.35763     1.35763     1.35763     2.09677
   8     1.35763     1.35763     1.35763     2.09677
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    10.05117
   Ω̃   =     5.42285
   ΩOD =     5.36573
   ΩD  =     0.05711
   Ω   =    15.47402

Splitting the model

Then we split the model into two, by diagonalizing the Wannier gauge Hamiltonian, we need to pass it the number of valence bands, n_val,

n_val = 4
(model_v, Uv), (model_c, Uc) = split_model(model, n_val)
2-element Vector{Tuple{Any, Any}}:
 (lattice: Å
  a1:  0.00000  2.71527  2.71527
  a2:  2.71527  0.00000  2.71527
  a3:  2.71527  2.71527  0.00000

atoms: fractional
  Si:  0.00000  0.00000  0.00000
  Si:  0.25000  0.25000  0.25000

n_bands: 4
n_wann : 4
kgrid  : 4 4 4
n_kpts : 64
n_bvecs: 8

b-vectors:
         [bx, by, bz] / Å⁻¹                weight
  1       0.28925   -0.28925    0.28925    1.49401
  2      -0.28925   -0.28925   -0.28925    1.49401
  3      -0.28925    0.28925   -0.28925    1.49401
  4      -0.28925   -0.28925    0.28925    1.49401
  5       0.28925    0.28925   -0.28925    1.49401
  6      -0.28925    0.28925    0.28925    1.49401
  7       0.28925   -0.28925   -0.28925    1.49401
  8       0.28925    0.28925    0.28925    1.49401, [-0.894002702921691 - 0.4414971779865484im 1.4083478477605758e-10 - 1.6136884027789426e-10im 1.8248865681011806e-10 - 2.696676689073669e-10im 8.528847239264198e-10 + 5.206305234474865e-10im; 3.8801508659345084e-10 + 1.312138130001942e-9im -0.11452048143542194 - 0.4631901009837881im 0.5257864281268069 + 0.19518429428353512im 0.08288800652856984 + 0.6087726607646692im; … ; 6.49639648027564e-10 + 1.7840748763106313e-10im -4.812986820969236e-9 + 3.1054570351505215e-9im 5.637282779951185e-9 + 2.6778505413831644e-9im 1.5626138820610927e-9 + 1.2668022879702624e-9im; 3.58173017012043e-9 - 2.003785583966899e-9im -4.813279497195566e-9 + 1.1642188766855521e-9im -4.244262009553412e-9 - 4.674619996647422e-10im -0.12953778624257484 + 0.06323972261193567im;;; 0.2664630951158803 + 0.9598314815603557im 0.007118931011431268 - 0.00696619829597428im 4.803887797950273e-10 - 5.658801272781643e-10im 4.098273778485805e-10 + 2.322998244712203e-10im; -0.003559608912873241 - 0.015102706282721479im 0.7286597971823794 - 0.6589841681703532im 1.8002962593622772e-8 - 5.073511264020493e-9im 7.949762199775097e-9 + 1.0302985160679763e-8im; … ; 9.677549515153478e-10 - 6.105207993104395e-10im -2.096908867071439e-9 - 5.580308556672876e-10im -0.0176236867562885 + 0.030871930529712086im 0.05421342280326157 + 0.001858352828750496im; 4.137972750896516e-9 + 1.799911585995144e-9im -3.590244730223988e-8 - 3.225774395982599e-8im -0.01951821527406029 - 0.050612106482135086im 0.023393690676321143 - 0.02676577891925633im;;; -0.27891927734293515 - 0.9530923218615222im -2.582700524240568e-9 + 8.419681693831469e-9im 3.64645139120758e-8 + 2.8521438924890404e-8im -7.481358038374672e-8 + 6.696313233499476e-8im; 5.183390013699773e-9 - 8.221401759653682e-9im 0.8959335097338064 - 0.4078965769631317im -1.6366019132477102e-7 - 4.306946428578446e-8im -6.541540084209365e-8 - 9.138867484179964e-8im; … ; -5.8821661551505825e-9 - 7.092986007886968e-9im 0.09669698419862097 + 0.08555688558026782im -5.2069525254318876e-8 - 3.1989871374079277e-7im 1.095486477719682e-7 - 1.8736321587967148e-7im; 0.023820357736336024 + 0.05718437156548336im 9.524340057905919e-9 + 4.187990421865672e-9im -2.1706854874796937e-7 - 1.376443845656832e-7im 3.7043416472921864e-7 - 4.035562622798368e-7im;;; … ;;; 0.7872046548624827 - 0.6088856030329838im -0.0037274385455807885 - 0.004388773891115048im -0.001605326609046026 - 0.012928823969105496im 2.7094233041528377e-9 + 1.3071910554654092e-9im; -0.0024865055892907273 - 0.00855334338708921im 0.9370818232330018 - 0.31955669536805814im 0.06722540616127048 + 0.0173028969106163im 1.981134379877069e-8 + 1.674503767170175e-8im; … ; -3.5844421603490506e-9 - 6.604466391860573e-9im 2.3299533735206817e-8 - 3.920371795758541e-8im -5.242666090960984e-9 + 2.722013383322293e-9im 0.11261983498552722 - 0.0002571298450818882im; 0.003448348545767397 - 0.008348384549399802im 0.015802349765101073 + 0.0056960288726943115im -0.027803099787172856 - 0.037012104791959516im 2.636313155893825e-8 - 2.196584238725692e-8im;;; -0.45505627978927965 + 0.8847818421435386im -1.4739961751164543e-9 + 4.247110663158082e-9im -0.005061137331093856 + 0.009840511700358625im 3.665394629436784e-9 + 1.3792780215943135e-8im; -4.3383215908576744e-9 - 1.956012435355703e-9im 0.9651249269868651 + 0.22246608832038522im 2.6986965828701415e-8 + 2.676924015787563e-8im -4.811213514586297e-9 - 1.3840720046143066e-8im; … ; 5.864705878161457e-9 - 3.954445898930521e-8im 1.6985415670943016e-8 + 9.858012062841267e-10im 1.7747601392353033e-8 - 2.6419174218226782e-8im -0.08451645229153253 - 0.13402670748752574im; 1.5875834024139945e-9 - 3.7177818074107805e-8im -0.0005767705903412851 - 0.0115073855477851im 1.0049234406856594e-8 - 2.1703693689478776e-8im -7.1891785524704995e-9 + 6.65361537521141e-9im;;; 0.7913216689933823 - 0.6050529482873804im 0.009117820340494659 + 0.00400904593380233im 2.2093412528400954e-9 + 1.740594891318179e-9im 1.205816202460965e-9 - 7.05180772778712e-10im; -0.014310240563215254 - 0.005998296992461932im 0.10477018352389192 + 0.9768461714329542im -2.438766063383388e-9 - 2.1149402022937462e-8im -1.0549048269118859e-8 - 2.9280723656301924e-9im; … ; 7.185891478423279e-11 - 2.7901009763160896e-10im 4.385038911941548e-10 - 8.442083953291385e-10im -0.016083473099792536 + 0.025034123750787327im -0.03434942484387839 - 0.04627044110737001im; 2.049550205919011e-8 - 1.0535206343346919e-8im -3.348533837787378e-8 - 3.9481457864336995e-8im 0.03740438797672109 + 0.04383771996274882im 0.027999771996675804 - 0.010069684111009905im])
 (lattice: Å
  a1:  0.00000  2.71527  2.71527
  a2:  2.71527  0.00000  2.71527
  a3:  2.71527  2.71527  0.00000

atoms: fractional
  Si:  0.00000  0.00000  0.00000
  Si:  0.25000  0.25000  0.25000

n_bands: 4
n_wann : 4
kgrid  : 4 4 4
n_kpts : 64
n_bvecs: 8

b-vectors:
         [bx, by, bz] / Å⁻¹                weight
  1       0.28925   -0.28925    0.28925    1.49401
  2      -0.28925   -0.28925   -0.28925    1.49401
  3      -0.28925    0.28925   -0.28925    1.49401
  4      -0.28925   -0.28925    0.28925    1.49401
  5       0.28925    0.28925   -0.28925    1.49401
  6      -0.28925    0.28925    0.28925    1.49401
  7       0.28925   -0.28925   -0.28925    1.49401
  8       0.28925    0.28925    0.28925    1.49401, [3.678613811636494e-11 + 1.0802962389017304e-10im 5.666452985396249e-10 + 4.6259072897313395e-10im -7.698032076490095e-11 + 1.8310541141933174e-10im -3.5235761912412135e-11 - 1.1062176772105295e-10im; 3.1878919077500037e-9 + 3.876777164427694e-9im 1.1605601294844389e-8 + 8.539662074849056e-9im -2.744863366708112e-7 + 5.96510080964425e-7im -1.4509510772362457e-10 + 1.120494130960083e-9im; … ; 8.932012724425882e-9 + 6.630172463811699e-9im 9.818638451588797e-11 + 5.70590705718868e-9im -1.8097570246109727e-8 - 1.3639506233349171e-8im 0.12597113637302335 + 0.1770601508288769im; -1.2656604835639797e-9 - 7.766553698178936e-11im 1.1589847947078431e-9 - 2.137975110073251e-9im -5.943352790341144e-7 - 6.024747051067262e-8im 4.782345356444402e-9 - 1.3675644877927692e-9im;;; -0.0027287510091086493 - 0.014270478283865323im 1.550626981234116e-9 + 3.82406307516645e-10im -2.1611333185249367e-10 + 8.053079143962366e-10im -0.005272112548761865 - 0.0021048253249410682im; -0.003361212015712718 - 0.02229867912153427im 1.8233552839114902e-8 + 1.875218504249511e-8im -7.2806519653241634e-9 + 1.3568069378516166e-8im -0.022714204870479815 - 0.010121075838996477im; … ; -3.603207769635486e-10 - 2.9281992938619986e-10im -0.13914568605970293 - 0.07665615104940901im -0.032724667822814635 + 0.26099380802079875im -1.0566780686743492e-9 + 1.0836750603848065e-9im; 1.3533052327692941e-8 + 2.4556687722543473e-9im 0.17503414009765833 - 0.19634566258331385im 0.14927199697515903 + 0.05436539825525473im 6.040758351301037e-8 - 1.1279015359259533e-7im;;; 6.933853284490038e-10 - 1.3958888596573377e-9im -1.1926510321089915e-9 - 2.9304781241326567e-9im -1.9467946867641512e-8 - 7.31454015284357e-9im -0.0006728675825196033 - 0.0022992645210359913im; 0.018104616153934158 - 0.008242586135421025im 1.4895666413999783e-8 - 1.7620064509256878e-9im 1.5463366101127355e-8 - 1.5389409974536279e-9im -9.260652852253116e-10 - 1.3125862458432664e-9im; … ; -0.11489276739863655 - 0.10165632795526464im 7.491946620341431e-8 + 1.4420943610604708e-7im 8.84573562643497e-8 + 1.5344027934737044e-7im 2.5038245561682696e-8 - 8.143718934474755e-9im; -1.0158457594390044e-8 - 3.451053162543223e-9im 1.89539003377207e-8 + 4.4480072751118703e-8im 2.8453830733572895e-7 + 9.269005279491711e-8im -0.11158297656069838 - 0.26787072201799894im;;; … ;;; 0.011535898739932687 - 0.006091427133098737im 0.006431059691891568 + 0.006076076030305239im -4.839900371383571e-9 + 1.2977795440461205e-8im 0.004532648521396114 - 0.0018039052121527274im; 0.002643969671168811 + 0.023857599238924866im -0.02107648237318979 + 0.009879056003432474im -4.142649242855431e-8 - 2.3090042195072876e-8im -1.6861334172835196e-5 - 0.005085808527842378im; … ; -5.657583708925154e-8 - 6.114304555827502e-9im -2.766799724332062e-7 - 1.873408517771635e-7im 0.0004203061135273495 + 0.184033521847266im 1.299622334182897e-7 - 3.654869910886304e-9im; -0.04356792679975066 + 0.0688415667931971im -0.32232448647029743 - 0.07761364880801433im -1.7014952171789717e-7 - 6.685463339374295e-7im -0.2488566413879828 + 0.3132730412941224im;;; -2.727015341243578e-9 - 5.1156419910616345e-9im 7.068822305278388e-9 - 2.2910474628545785e-8im -0.004510160355871182 + 0.00876922636856443im 2.789952137371987e-9 - 5.362865873603795e-9im; -2.073630798006839e-8 - 1.1483744557626101e-7im 0.0005393964137048694 + 0.00012436635137648677im -1.491475903891981e-8 - 7.785568089530015e-9im 0.024127911790728984 + 0.005561545833676872im; … ; -1.4019912924143063e-6 + 5.244577149232654e-7im -2.8220944395188258e-9 + 2.7566927040594105e-8im -1.0039793542534511e-7 + 2.909260780252526e-7im -2.0234480273036903e-8 + 1.4325239151271704e-8im; 1.8952692871498668e-6 - 9.271214677843525e-7im 0.015031661345467491 + 0.30012219915242977im -9.287518911056669e-8 + 3.637479870088889e-7im 0.024084829983951478 + 0.48048934377056385im;;; 0.012224727433499125 - 0.007851666827181858im 1.2238235025332585e-9 - 1.7834435132995895e-10im 4.551118069214883e-9 - 1.5796184728303966e-9im 1.2160668863628321e-5 + 0.005676734311796952im; 0.020015004882040943 + 0.010388860202373701im -5.795354481193922e-9 - 3.801210909590518e-9im -2.5545081218142884e-8 - 1.208966141175239e-8im -0.02154000939683698 + 0.0124257503932447im; … ; 7.23726430363849e-10 - 1.3432574576337086e-9im 0.2858582294121772 - 0.08561822397106772im -0.06608239212100842 + 0.031844847730647885im -3.2168699007172713e-9 + 2.004261058717154e-9im; -4.2417959877857786e-8 + 8.251831332121284e-9im 0.07332021794776017 - 0.00226320112203956im 0.29309615381272863 - 0.05603628135798971im 6.089151343631455e-8 - 8.99527090536911e-8im])

Returned are two separated Models, and the two corresponding gauge transformation matrices from the total manifold to the two separated submanifolds. This is useful if you want to further rotate the gauge of some operators, or if you want to rotate the UNK files for plotting WFs, in that case you need to pass these two matrices to the function split_unk.

Next, we construct PTG for the two Models,

Uv, _ = parallel_transport(model_v)
model_v.U .= Uv;
parallel transport
[ Info: Filling (kx,0,0)
[ Info: Filling (kx,ky,0)
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = -0.004577663482537961
[ Info: Filling (k1,k2,k3)
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = 0.04924394292523244
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = -0.0222702466235844
initial error = 0.5998
final error   = 0.1618

and

Uc, _ = parallel_transport(model_c)
model_c.U .= Uc;
parallel transport
[ Info: Filling (kx,0,0)
[ Info: Filling (kx,ky,0)
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = -0.06880416726126168
[ Info: Filling (k1,k2,k3)
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = -0.0650573287936393
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
[ Info: Pole chosen by Barycenter of path
┌ Info: Chern number
└   m = -0.016374288184603365
initial error = 0.6084
final error   = 0.189

and take a look at the spread

omega(model_v)
omega(model_c)
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.12266     0.08267     0.38530     4.67622
   2    -0.00147    -0.16371    -0.17590     4.96181
   3    -0.00357     0.34725     0.07612     4.66569
   4    -0.06925    -0.26337    -0.28859     5.16882
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    11.25427
   Ω̃   =     8.21827
   ΩOD =     6.81834
   ΩD  =     1.39993
   Ω   =    19.47254
Tip

There is a split_wannierize function that combines the above steps, however, for pedagogical purpose we run them manually.

Further maximal localization

We can further run maximal localization to smoothen the gauge,

Uv = max_localize(model_v)
model_v.U .= Uv;
[ Info: Initial spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1    -0.00414    -0.06009    -0.09507     2.28300
   2     0.03873     0.23203    -0.02785     3.09473
   3     0.06834     0.03812     0.38726     3.00592
   4    -0.07576    -0.19926    -0.25918     3.08638
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =     5.38374
   Ω̃   =     6.08629
   ΩOD =     5.85950
   ΩD  =     0.22680
   Ω   =    11.47003


Iter     Function value   Gradient norm
     0     1.147003e+01     6.488197e-01
 * time: 0.0014870166778564453
     1     1.020453e+01     3.767061e-02
 * time: 0.011614084243774414
     2     9.975831e+00     2.070231e-02
 * time: 0.021627187728881836
     3     9.939701e+00     1.357925e-02
 * time: 0.031783103942871094
     4     9.907616e+00     2.561955e-02
 * time: 0.04177999496459961
     5     9.767324e+00     4.729370e-02
 * time: 0.05472111701965332
     6     9.514631e+00     5.809008e-02
 * time: 0.06488609313964844
     7     9.013009e+00     4.985240e-02
 * time: 0.07481908798217773
     8     8.482475e+00     6.823858e-02
 * time: 0.08490419387817383
     9     7.969027e+00     5.785297e-02
 * time: 0.09505820274353027
    10     6.977420e+00     8.881997e-02
 * time: 0.105010986328125
    11     6.360951e+00     5.122475e-02
 * time: 0.1121680736541748
    12     6.050457e+00     3.714282e-02
 * time: 0.12244820594787598
    13     5.907667e+00     2.509691e-02
 * time: 0.13277101516723633
    14     5.831826e+00     2.250381e-02
 * time: 0.1439990997314453
    15     5.780070e+00     1.405661e-02
 * time: 0.1550910472869873
    16     5.765670e+00     7.922847e-03
 * time: 0.16260409355163574
    17     5.761047e+00     3.477379e-03
 * time: 0.1735830307006836
    18     5.759176e+00     3.414329e-03
 * time: 0.18464899063110352
    19     5.757903e+00     2.157822e-03
 * time: 0.19580698013305664
    20     5.757275e+00     1.757323e-03
 * time: 0.20682716369628906
    21     5.757036e+00     7.337036e-04
 * time: 0.21785902976989746
    22     5.756973e+00     5.143370e-04
 * time: 0.28983020782470703
    23     5.756948e+00     2.764406e-04
 * time: 0.2996561527252197
    24     5.756935e+00     2.153013e-04
 * time: 0.30963921546936035
    25     5.756929e+00     1.709295e-04
 * time: 0.3197669982910156
    26     5.756925e+00     1.354181e-04
 * time: 0.3298771381378174
    27     5.756924e+00     6.193683e-05
 * time: 0.3400561809539795
    28     5.756924e+00     2.520679e-05
 * time: 0.35010814666748047
    29     5.756924e+00     2.572938e-05
 * time: 0.36008310317993164
[ Info: Final spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1    -0.67882    -0.67882     0.67881     1.43925
   2    -0.67882     0.67881    -0.67881     1.43935
   3     0.67881     0.67883     0.67882     1.43771
   4     0.67882    -0.67881    -0.67882     1.44061
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =     5.38374
   Ω̃   =     0.37319
   ΩOD =     0.37319
   ΩD  =     0.00000
   Ω   =     5.75692

and

Uc = max_localize(model_c)
model_c.U .= Uc;
[ Info: Initial spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.12266     0.08267     0.38530     4.67622
   2    -0.00147    -0.16371    -0.17590     4.96181
   3    -0.00357     0.34725     0.07612     4.66569
   4    -0.06925    -0.26337    -0.28859     5.16882
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    11.25427
   Ω̃   =     8.21827
   ΩOD =     6.81834
   ΩD  =     1.39993
   Ω   =    19.47254


Iter     Function value   Gradient norm
     0     1.947254e+01     5.380885e-01
 * time: 0.0014748573303222656
     1     1.694777e+01     1.100288e-01
 * time: 0.011771917343139648
     2     1.620877e+01     5.806310e-02
 * time: 0.021925926208496094
     3     1.600008e+01     3.685524e-02
 * time: 0.03201603889465332
     4     1.591740e+01     2.023350e-02
 * time: 0.042122840881347656
     5     1.586866e+01     2.721928e-02
 * time: 0.05202889442443848
     6     1.573774e+01     2.502966e-02
 * time: 0.06218981742858887
     7     1.553735e+01     2.791503e-02
 * time: 0.07219886779785156
     8     1.532713e+01     4.696656e-02
 * time: 0.08213186264038086
     9     1.505181e+01     4.372920e-02
 * time: 0.09226393699645996
    10     1.460699e+01     6.126373e-02
 * time: 0.10251998901367188
    11     1.398381e+01     6.352369e-02
 * time: 0.11300396919250488
    12     1.339260e+01     6.758271e-02
 * time: 0.1245729923248291
    13     1.308261e+01     4.371015e-02
 * time: 0.1355748176574707
    14     1.282837e+01     4.676811e-02
 * time: 0.14701581001281738
    15     1.267524e+01     3.244126e-02
 * time: 0.15802383422851562
    16     1.258424e+01     2.141517e-02
 * time: 0.16915488243103027
    17     1.253841e+01     1.967834e-02
 * time: 0.18019890785217285
    18     1.251652e+01     8.469666e-03
 * time: 0.19142603874206543
    19     1.250812e+01     7.368176e-03
 * time: 0.20259785652160645
    20     1.250232e+01     5.107745e-03
 * time: 0.21372199058532715
    21     1.249965e+01     3.967977e-03
 * time: 0.22469592094421387
    22     1.249841e+01     2.316879e-03
 * time: 0.2357029914855957
    23     1.249790e+01     1.689382e-03
 * time: 0.2467188835144043
    24     1.249759e+01     1.269381e-03
 * time: 0.25778985023498535
    25     1.249743e+01     1.058153e-03
 * time: 0.26897192001342773
    26     1.249734e+01     6.166992e-04
 * time: 0.2801778316497803
    27     1.249730e+01     3.514423e-04
 * time: 0.29138684272766113
    28     1.249728e+01     2.925587e-04
 * time: 0.3677358627319336
    29     1.249728e+01     1.929891e-04
 * time: 0.37770795822143555
    30     1.249727e+01     1.497478e-04
 * time: 0.3877270221710205
    31     1.249727e+01     1.363214e-04
 * time: 0.39771389961242676
    32     1.249727e+01     6.374996e-05
 * time: 0.40775084495544434
    33     1.249727e+01     5.605883e-05
 * time: 0.4177539348602295
[ Info: Final spread
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.67887     0.67884     0.67878     3.12428
   2     0.67877    -0.67893    -0.67880     3.12427
   3    -0.67883     0.67879    -0.67883     3.12439
   4    -0.67881    -0.67870     0.67885     3.12433
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    11.25427
   Ω̃   =     1.24300
   ΩOD =     1.24300
   ΩD  =     0.00000
   Ω   =    12.49727
Tip

For such simple case (few bands, very little kpoints) a direct maximal localization is sufficient. But for more complex scenarios, running an opt_rotate might be helpful.

Comparing spreads

The valence + conduction,

model.U .= U;
omega(model)
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.00000     0.00000     0.00000     1.44670
   2     0.00000     0.00000    -0.00000     2.09677
   3    -0.00000    -0.00000     0.00000     2.09677
   4    -0.00000    -0.00000    -0.00000     2.09677
   5     1.35763     1.35763     1.35763     1.44670
   6     1.35763     1.35763     1.35763     2.09677
   7     1.35763     1.35763     1.35763     2.09677
   8     1.35763     1.35763     1.35763     2.09677
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    10.05117
   Ω̃   =     5.42285
   ΩOD =     5.36573
   ΩD  =     0.05711
   Ω   =    15.47402

the valence,

omega(model_v)
  WF     center [rx, ry, rz]/Š             spread/Ų
   1    -0.67882    -0.67882     0.67881     1.43925
   2    -0.67882     0.67881    -0.67881     1.43935
   3     0.67881     0.67883     0.67882     1.43771
   4     0.67882    -0.67881    -0.67882     1.44061
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =     5.38374
   Ω̃   =     0.37319
   ΩOD =     0.37319
   ΩD  =     0.00000
   Ω   =     5.75692

the conduction,

omega(model_c)
  WF     center [rx, ry, rz]/Š             spread/Ų
   1     0.67887     0.67884     0.67878     3.12428
   2     0.67877    -0.67893    -0.67880     3.12427
   3    -0.67883     0.67879    -0.67883     3.12439
   4    -0.67881    -0.67870     0.67885     3.12433
Sum spread: Ω = ΩI + Ω̃, Ω̃ = ΩOD + ΩD
   ΩI  =    11.25427
   Ω̃   =     1.24300
   ΩOD =     1.24300
   ΩD  =     0.00000
   Ω   =    12.49727

Look, we have well-localized WFs for all the three cases, the initial 8 s,p orbitals are separated into 4 bonding and 4 anti-bonding orbitals! 🚀

Band structure

Finally, let's compare band interpolation,

using PlotlyJS

the valence + conduction,

interp_model = Wannier.InterpModel(model)
kpi, E = interpolate(interp_model)
(StaticArraysCore.SVector{3, Float64}[[0.0, 0.0, 0.0], [0.005, 0.0, 0.005], [0.01, 0.0, 0.01], [0.015, 0.0, 0.015], [0.02, 0.0, 0.02], [0.024999999999999998, 0.0, 0.024999999999999998], [0.03, 0.0, 0.03], [0.034999999999999996, 0.0, 0.034999999999999996], [0.04, 0.0, 0.04], [0.04499999999999999, 0.0, 0.04499999999999999]  …  [0.5, 0.045000000000000005, 0.545], [0.5, 0.04, 0.5399999999999999], [0.5, 0.034999999999999996, 0.535], [0.5, 0.029999999999999985, 0.53], [0.5, 0.02499999999999998, 0.525], [0.5, 0.020000000000000025, 0.52], [0.5, 0.015000000000000017, 0.515], [0.5, 0.010000000000000012, 0.51], [0.5, 0.005000000000000006, 0.505], [0.5, 0.0, 0.5]], [-5.798515544177007 -5.797983061606554 … -1.6351007785473197 -1.6353082316090068; 6.202325333825 6.1984987174204775 … -1.6351007730836997 -1.635308231608975; … ; 8.703475344557578 8.705889267969335 … 16.246495461035114 16.247897361941277; 9.514549688154366 9.518709092702471 … 16.24649551718928 16.247897361941423])

the valence,

interp_model_v = Wannier.InterpModel(model_v)
_, Ev = interpolate(interp_model_v)
(StaticArraysCore.SVector{3, Float64}[[0.0, 0.0, 0.0], [0.005, 0.0, 0.005], [0.01, 0.0, 0.01], [0.015, 0.0, 0.015], [0.02, 0.0, 0.02], [0.024999999999999998, 0.0, 0.024999999999999998], [0.03, 0.0, 0.03], [0.034999999999999996, 0.0, 0.034999999999999996], [0.04, 0.0, 0.04], [0.04499999999999999, 0.0, 0.04499999999999999]  …  [0.5, 0.045000000000000005, 0.545], [0.5, 0.04, 0.5399999999999999], [0.5, 0.034999999999999996, 0.535], [0.5, 0.029999999999999985, 0.53], [0.5, 0.02499999999999998, 0.525], [0.5, 0.020000000000000025, 0.52], [0.5, 0.015000000000000017, 0.515], [0.5, 0.010000000000000012, 0.51], [0.5, 0.005000000000000006, 0.505], [0.5, 0.0, 0.5]], [-5.683474999400311 -5.683160312458814 … -1.6189465249916333 -1.6185917276750166; 7.003492063601115 7.000297044503034 … -1.6174972033692596 -1.6185916252845773; 7.003492125235773 7.002570226301684 … 3.3584415294668695 3.3596391220216146; 7.462665898770761 7.461472339853736 … 3.359881377050272 3.3596391239841967])

the conduction,

interp_model_c = Wannier.InterpModel(model_c)
_, Ec = interpolate(interp_model_c)
(StaticArraysCore.SVector{3, Float64}[[0.0, 0.0, 0.0], [0.005, 0.0, 0.005], [0.01, 0.0, 0.01], [0.015, 0.0, 0.015], [0.02, 0.0, 0.02], [0.024999999999999998, 0.0, 0.024999999999999998], [0.03, 0.0, 0.03], [0.034999999999999996, 0.0, 0.034999999999999996], [0.04, 0.0, 0.04], [0.04499999999999999, 0.0, 0.04499999999999999]  …  [0.5, 0.045000000000000005, 0.545], [0.5, 0.04, 0.5399999999999999], [0.5, 0.034999999999999996, 0.535], [0.5, 0.029999999999999985, 0.53], [0.5, 0.02499999999999998, 0.525], [0.5, 0.020000000000000025, 0.52], [0.5, 0.015000000000000017, 0.515], [0.5, 0.010000000000000012, 0.51], [0.5, 0.005000000000000006, 0.505], [0.5, 0.0, 0.5]], [8.703475344555013 8.697558978534058 … 13.506724376369926 13.506720002703762; 8.703475344555052 8.70641428175449 … 13.506730395212003 13.506720096470826; 8.703475344560456 8.712052268359766 … 16.247428090393655 16.24789736194152; 10.069264200454386 10.070263551329775 … 16.2474304574677 16.247897361960934])

and compare valence,

P = plot_band_diff(kpi, E, Ev)

and conduction,

P = plot_band_diff(kpi, E, Ec)
Note

Again, the interpolation quality is not good (especially for the conduction which has much larger spread WFs), bacause of the very coarse kgrid, increase the kgrid and compare the resulting band structures.

Bravo, we have successfully Wannierized the valence and conduction manifolds, by splitting the gauge of a valence + conduction calculation!


This page was generated using Literate.jl.