Questaal Home
Navigation

Nb/Fe/Nb Metallic Trilayer, Electronic Structure and Landauer-Buttiker Transport

This tutorial studies the Landauer-Buttiker transport through a Nb/Fe/Nb trilayer. The setup of the trilayer is complicated, and as a first step the Nb/Fe superlattice (with periodic boundary conditions) is constructed, as explained in this tutorial.

Nb and Fe are both bcc but are strongly lattice mismatched, and it necessitated rotating them to axes different from the cubic axes. See the Nb/Fe superlattice tutorial for its construction and the relaxation of atoms at the interface.

This tutorial begins where the superlattice tutorial leaves off. It uses the Atomic Spheres Approximation, and the principal layer code designed to calculate Landauer-Buttiker transmission through a multilayer.

There is an additional complication when the ASA is used — namely sphere packing for overlapping spheres becomes a challenge for a complex interface such as this. The early parts of this tutorial modify the original SL by adding empty spheres, before proceeding to the (nonperiodic) Nb/Fe/Nb trilayer.


Table of Contents

Energy bands of Nb (3×1) reconstructed cell in the ASA

This section requires files init.nb and site3.nb.

cp init.nb init.nbasa
blm --nl=3 --nfile --nk=12,12,12 --asa nbasa
cat actrl.nbasa | sed 's/LMX=2/LMX=2  GRP2=1/'  | awk '{ if ($1 == "OPTIONS") {print "HAM   PMIN=-1 QASA=1"}; print}' > ctrl.nbasa
cp site3.nb site3.nbasa
lmstr -vfile=3 ctrl.nbasa
lm -vnit=0 -vfile=3 ctrl.nbasa
lm -vnit=30 -vfile=3 ctrl.nbasa --rs=0,1
lmchk ctrl.nbasa -vfile=3 --syml~n=31~q=0,1,0,0,0,0,1/2,0,-1/2
lm -vnit=30 -vfile=3 ctrl.nbasa --band~fn=syml
echo -6,6,5,10 | plbnds -fplot~sh~scl=.5 -ef=0 -scl=13.6 --nocol -dat=red nbasa
fplot -f plot.plbnds
open fplot.ps

Energy bands of Fe (4×1) reconstructed cell in the ASA This section requires files init.fe and site4.fe.

cp init.fe init.feasa
blm --nfile --nk=10,5,13 --mag --asa feasa
cat actrl.feasa | sed 's/\(MMOM.*\)/\1 GRP2=1/'  | awk '{ if ($1 == "OPTIONS") {print "HAM   PMIN=-1 QASA=1"}; print}' > ctrl.feasa
cp site4.fe site4.feasa
lmstr -vfile=4 ctrl.feasa
lm -vnit=0 -vfile=4 ctrl.feasa
lm -vnit=30 -vfile=4 ctrl.feasa --rs=0,1
cp syml.nbasa syml.feasa
lm -vnit=30 -vfile=4 ctrl.feasa --band~fn=syml
echo -6,6,5,10 | plbnds -fplot~sh~scl=.5 -ef=0 -scl=13.6 --nocol -dat=red -spin1 feasa
fplot -f plot.plbnds
open fplot.ps

Preliminaries

This tutorial makes extensive use of many of the Questaal packages. including lmscell, lmplan lm, lmgf, lmpg, and miscellaneous utilities such as vextract. They are assumed to be in your path, or in the directory where you cloned the Questaal package (assumed here to be ~/lm).

It invokes mpix to execute MPI jobs. mpix is a proxy for mpirun, the default command that launches MPI jobs. We assume your system has 16 processors on a core, but you can use any number, or run in serial.

For drawing postscript files, this tutorial assumes you will use the Apple open. Substitute your postscript viewer of choice for open.

Input files

This tutorial uses the following files from Nb/Fe superlattice tutorial init.nb, site3.nb, init.fe, site3.fe all found in this section, and site4.nbfe, created in this section. If you want to skip that tutorial, you can

Copy file init.nb to init.nbasa:

# Init file for Nb.  Lattice constant at 0K
LATTICE
    ALAT=6.226 PLAT= -.5 .5 .5  .5 -.5 .5  .5 .5 -.5
SPEC ATOM=Nb  PZ=0,0,0
SITE
    ATOM=Nb   X=0 0 0

File site3.nb

% site-data vn=3.0 fast io=31 nbas=3 alat=6.226 plat= -0.5 0.5 0.5 1.5 1.5 -1.5 0.0 0.0 1.0
#                        pos
 Nb        0.0000000   0.0000000   0.0000000
 Nb        0.5000000   0.5000000  -0.5000000
 Nb        1.0000000   1.0000000  -1.0000000

File init.fe

# Init file for Fe.  Lattice constant at 0K
LATTICE
    ALAT=5.404 PLAT= -.5 .5 .5  .5 -.5 .5  .5 .5 -.5
SPEC
    ATOM=Fe  MMOM=0,0,2.2
SITE
    ATOM=Fe   X=0 0 0

File site4.fe

% site-data vn=3.0 fast io=15 nbas=4 alat=6.226 plat= -0.5 0.5 0.5 1.5 1.5 -1.5 0.18425022 -0.50338096 0.68763117
#                        pos
 Fe        0.0000000   0.0000000   0.0000000
 Fe        0.0000000   0.7500000   0.0000000
 Fe        0.5000000   1.0000000  -0.5000000
 Fe        1.0000000   1.2500000  -1.0000000

Introduction and Basic Strategy

This tutorial shows how to make a self-consistent density starting from the superlattice developed in the Nb/Fe superlattice tutorial, within the Atomic Spheres approximation.

Initially the ASA-LDA code lm code is then applied to the relaxed superlattice to confirm the reliability of the ASA.

Next the crystal Green’s function package lmgf as setup to the layer Green’s function package lmpg.

Finally a structure is set up for trilayer geometry (… Nb / Fe / Nb …), which is related to but distinct from the Nb/Fe superlattice, as it is not periodic along the stacking axis. Transmission calculations through the trilayer are performed using lmpg.

In the Nb/Fe superlattice tutorial a superlattice of Nb and Fe was constructed; structural information is contained in site4.nbfe. Starting from it this tutorial will do:

  • Make a ctrl and site file suitable for the ASA, with periodic boundary conditions. The interface has an additional complication in that open spaces appear, which the ASA must deal with by packing with additional empty spheres. Care must be taken to fill space, but not allow the sphere overlaps become too large.
  • Repeat the self-consistent calculation of the interface code with the ASA-LDA code lm. lm is more approximate than the full LDA code lmf, but it is very efficient, and provides a stepping stone to final objective, transport through a trilayer
  • Repeat the calculation once again with the crystal Green’s function code, lmgf, and demonstrate that the spectral functions are similar to the bands created by lm.
  • Setup the trilayer geometry.
  • Calculate transmission with the transport code, lmpg.

1. Self-consistent ASA-LDA energy bands for reconstructed Nb

It is instructive to perform a band calculation for Nb in the (3×1) reconstructed superlattice in the ASA approximation, to compare it to the full LDA result generated here.

This step isn’t essential for the superlattice, but it is useful as a sanity check we will need the ASA for transport and need to confirm its reliability for bulk Nb.

cp init.nb init.nbasa
blm --nl=3 --nfile --nk=12,12,12 --asa nbasa
cat actrl.nbasa | sed 's/LMX=2/LMX=2  GRP2=1/'  | awk '{ if ($1 == "OPTIONS") {print "HAM   PMIN=-1 QASA=1"}; print}' > ctrl.nbasa

blm makes an acceptable ctrl file, but we make a few tweaks in anticipation of the fact that ultimately we will use lmpg for transport. This code doesn’t allow down-folding, and it uses KKR style accumulation of weights. Compare the changes and refer to the input file documentation for the meaning of each new tag.

Check the sphere overlaps

cp site3.nb site3.nbasa
lmchk -vfile=3 ctrl.nbasa --shell

Because of the reconstruction, the three Nb split into two classes. Each class has 8 nearest neighbors, 6 second neighbors, 12 third neighbors, etc. characteristic of the bcc latttice. The MT sphere has been enlarged relative to the full-potential result so that the sphere volume matches the unit cell volume, as the ASA requires (look for  Sum of sphere volumes).

Make the ASA density self-consistent in the reconstructed structure

lmstr -vfile=3 ctrl.nbasa
lm -vnit=0 -vfile=3 ctrl.nbasa
lm -vnit=30 -vfile=3 ctrl.nbasa
lmchk ctrl.nbasa -vfile=3 --syml~n=31~q=0,1,0,0,0,0,1/2,0,-1/2
lm -vnit=30 -vfile=3 ctrl.nbasa --band~fn=syml

Draw the ASA bands:

echo -6,6,5,10 | plbnds -fplot~sh~scl=.5 -ef=0 -scl=13.6 --nocol -dat=red -spin1 nbasa
open fplot.ps

The ASA energy bands are quite similar to the full LDA ones, and sufficiently close for the present purposes. Note the bands cross the Fermi surface at similar points with similar slope.

Blue: LDA bands generated by lmf. Red: ASA-LDA generated by lm.

2. Self-consistent ASA-LDA energy bands for reconstructed Fe

As in the Nb case, it is instructive to perform a band calculation for ferromagnetic Fe in the (4×1) reconstructed superlattice in the ASA approximation, to compare it to the full LDA result generated here.

cp init.fe init.feasa
blm --nfile --nk=10,5,13 --mag --asa feasa
cat actrl.feasa | sed 's/\(MMOM.*\)/\1 GRP2=1/'  | awk '{ if ($1 == "OPTIONS") {print "HAM   PMIN=-1 QASA=1"}; print}' > ctrl.feasa

Check the sphere overlaps

cp site4.fe site4.feasa
lmchk -vfile=4 ctrl.feasa

The MT sphere has been enlarged so that the sphere volume matches the unit cell volume, as the ASA requires.

Make the ASA density self-consistent in the reconstructed structure

lmstr -vfile=4 ctrl.feasa
lm -vnit=0 -vfile=4 ctrl.feasa
lm -vnit=30 -vfile=4 ctrl.feasa
cp syml.fe syml.feasa
lm -vnit=30 -vfile=4 ctrl.feasa --band~fn=syml

Draw the ASA bands:

echo -6,6,5,10 | plbnds -fplot~sh~scl=.5 -ef=0 -scl=13.6 --nocol -dat=red feasa
open fplot.ps

As with Nb, the ASA energy bands are quite similar to the full LDA ones, and indeed significantly closer to each other than to bands computed from a higher level theory such as GW.

Majority (left) and minority (right) energy bands of Fe. blue: LDA generated by lmf; red: ASA-LDA generated by lm.

FP and ASA energy bands of Fe compared

3. Self-consistent ASA Calculation of the Nb/Fe interface

The ASA is more difficult to set up because of the constraint that the spheres fill space: the sum-of-sphere volumes must equal the volume of the unit cell. Also In the end regions we want to keep the spheres the same as the bulk Nb. At the interface there is considerable reconstruction which makes good sphere packing difficult.

Create a directory asa below your current working directory to keep the files distinct.

We require the structure generated by the Nb/Fe superlattice tutorial. Either copy the structure file of the relaxed interface, site4.nbfe, to file site9.asa in your working directory, or paste the contexts of the box below into it.

% site-data vn=3.0 fast io=15 nbas=24 alat=6.226 plat= -0.5 0.5 0.5 1.5 1.5 -1.5 0.55275066 -1.51014288 6.06289351
#                        pos
 Nb        0.0000000   0.0000000   0.0000000
 Nb        0.5000000   0.5000000  -0.5000000
 Nb        1.0000000   1.0000000  -1.0000000
 Nb        0.0038386  -0.0166966   0.9910505
 Nb        0.5161225   0.4990571   0.5150855
 Nb        0.9928133   1.0176382   0.0048461
 Fe       -0.0259839   0.7809822   1.9780095
 Fe        0.5160204   0.9996464   1.5211879
 Fe        0.9792297   1.2189900   0.9795957
 Fe       -0.0399873  -0.0010873   1.9603949
 Fe        0.1509317   0.2442758   2.6493636
 Fe        0.6537761   0.4947430   2.1670347
 Fe        1.1725712   0.7480391   1.6834229
 Fe        0.1626577  -0.5019115   2.6567494
 Fe        0.3308379  -0.2533275   3.3370495
 Fe        0.8421375  -0.0133975   2.8192889
 Fe        1.3371460   0.2434396   2.3440138
 Fe        0.3278900  -0.9987992   3.3586058
 Nb        0.5490728  -1.5111615   4.0608585
 Nb        1.0623082  -1.0379833   3.5524052
 Nb        1.5400241  -0.4828374   3.0627345
 Nb        0.5527507  -1.5101429   5.0628935
 Nb        1.0527507  -1.0101429   4.5628935
 Nb        1.5527507  -0.5101429   4.0628935

A working ctrl file is needed to use Questaal’s tools. The ctrl file in the Nb/Fe superlattice tutorial was designed for lmf, and it is a bit different from the ASA.

Here we start with minimum information from that tutorial and build ctrl from scratch. Make a preliminary file with

blm --nfile --rdsite=site9 --nosort --shorten=no --ctrl=ctrl asa

To manage sphere packing in the ASA we must distinguish sites immediately at the interface from those in the bulk. Recall the superlattice geometry:

  Layer    Nb 1      Nb 2      Fe 1       Fe 2     Fe 3      Nb 3      Nb 4
  sites    1-3       4-6       7-10       11-14    15-18     19-21     22-24

The frontier sites are 4-6 and 19-21 (Nb) and 7-10 and 15-18 (Fe). Use the lmscell, utility to assign species Nbx and Fex to these files, or alternatively just use a text editor and do it by hand. The following creates site8.asa

lmscell ctrl.asa --stack~setlbl@spid=Nbx@targ=4:6,19:21~setlbl@spid=Fex@targ=7:10,15:18~wsite@short@fn=site8

site8.asa differs from site9.asa only in that some labels changed.

To obtain a site file with adequate sphere packing, empty spheres will need to be added. This is not trivial. The automatic empty sphere finder in the blm utility greatly facilitates this step, and even with the blm it takes some tinkering. The set of switches was found to work:

blm --pr35 --nl=3 --gf --nfile --rdsite=site8 --nosort --shorten=no --nk=9,4,2 --mag --asa~pfree~rmaxs=7.2 --wsite@fn=site7 asa -vfile=4 --findes~rmin=.9~shorten=0,0,1~omax=0~sclwsr=20~spec~nspmx=4~nesmx=36 --omax=.18,.26,.26
cp actrl.asa ctrl.asa

Section needs to be completed …

In construction of the ctrl and sites files, attention needs to be paid to the following:

  1. The total sphere volume equals the unit cell volume. Confirm that the autogenerated files meet this requirement:
    lmchk ctrl.asa -vfile=7 | grep 'Sum of sphere'
    
  2. No sphere overlap is too large. –omax=.18,.26,.26 imposes constraints on the atom-atom, atom-ES, and ES-ES overlaps, respectively. To see what these overlaps actually are, do:
    lmchk ctrl.asa -vfile=7 | grep 'max ovlp'
    
  3. The Nb sphere radius should be the same as the bulk value. This isn’t strictly necessary, but if this condition is not met, the density inside the Nb spheres will not be exactly the same as the bulk, and we should impose that radius as a constraint. lmchk can do this; by telling it to resize the spheres (tag SPEC_SCLWSR=11). We must also constrain the Nb sphere radius to that of bulk nb, R=3.065511 (see ctrl.nbasa in the Nb/Fe superlattice tutorial). Temporarily the adjustments ctrl.asa, e.g.

    cp ctrl.asa ctrl.bk
    sed -i -e 's/R= 3.[0-9]*/R= 3.065511 CSTRMX=1/' -e 's/# SCLWSR/  SCLWSR/' ctrl.asa
    lmchk ctrl.asa -vfile=7 | grep -A 10 'spec  name'
    

    You should see the following table:

    spec  name        old rmax    new rmax     ratio
      1   Nb          3.065511    3.065511    1.000000
      2   Nbx         2.730129    2.730129    1.000000
      3   Fex         2.496030    2.496030    1.000000
      4   Fe          2.719829    2.719829    1.000000
      5   E           1.542848    1.715533    1.111926
      6   E1          1.483502    1.710466    1.152992
      7   E2          1.432208    1.651325    1.152992
      8   E3          1.081334    1.202364    1.111926
    

    The E spheres were resized, but the others were not. Restore ctrl.asa from ctrl.bk and use your text editor to substitute sphere radii in the table above for the original ones. Check once more the total sphere volume and the overlaps

    lmchk ctrl.asa -vfile=7 | grep 'Sum of sphere'
    lmchk ctrl.asa -vfile=7 | grep 'max ovlp'
    

    You can confirm that the maximum overlap between real atoms remains at 10%, while the overlaps with E sites increased. This the best we can do without adding more very small E spheres and it is not necessary.

  4. This step is not essential, but makes matters simpler, if the sites are ordered by increasing projection along the line normal to the interface. lmscell can do this automatically:
    lmscell -vfile=7 ctrl.asa --stack~sort@h3~wsite@short@fn=site6
    

    Note that h3 evaluates to the projection of a given site normal to the plane given by PLAT1×PLAT2. The rearranged sites in site6.asa clarifies where the E spheres actually go.

ctrl.asa:

# Autogenerated from init.nbfe using:
# blm --pr35 --gf --nfile --rdsite=site4 --nosort --shorten=no --nk=9,4,2 --mag --asa~pfree~rmaxs=7 --wsite@fn=null nbfe -vfile=4 --findes~rmin=.9~shorten=0,0,1~omax=0~sclwsr=20~spec~nspmx=4~nesmx=36 --omax=.18,.26,.26

# Variables entering into expressions parsed by input
% const nit=10 asa=t les=1
% const met=(asa?1:5)
% const nsp=2 so=0
% const lxcf=2 lxcf1=0 lxcf2=0     # for PBE use: lxcf=0 lxcf1=101 lxcf2=130
% const ccor=T sx=0 gamma=sx scr=4 # ASA-specific
% const rmaxs=7                    # Range for structure constants
% const gfmode=1 nz=16 ef=0 c3=t   # lmgf-specific variables
% const nk1=9 nk2=4 nk3=2 beta=.3

VERS  LM:7 ASA:7 # FP:7
IO    SHOW=f HELP=f IACTIV=f VERBOS=35,35  OUTPUT=*
EXPRESS
# Lattice vectors and site positions
%ifdef file>=1
  file=   site{file}
%endif
  file=   site
# file= essite

# Self-consistency
  nit=    {abs(nit)}               # Maximum number of iterations
#  mix=    A2,b={beta},n=5,k=5,0
#   mix=    A7,b={beta},n=3,k=3;A2,b={beta/10},n=2,k=2,elind=0.001
#   mix=    A7,b={beta},n=7,k=7;A2,b={beta/10},n=2,k=2,elind=0
  mix=    B2,b={beta},k=7          # Charge density mixing parameters
  conv=   1e-5                     # Convergence tolerance (energy)
  convc=  3e-5/3                   # tolerance in RMS (output-input) density

# Brillouin zone
  nkabc=  {nk1},{nk2},{nk3}        # 1 to 3 values
  metal=  {met}                    # Management of k-point integration weights in metals

# Potential
  nspin=  {nsp}                    # 2 for spin polarized calculations
  so=     {so}                     # 1 turns on spin-orbit coupling
  xcfun=  {lxcf},{lxcf1},{lxcf2}   # set lxcf=0 for libxc functionals

SYMGRP  find
HAM   PMIN=-1                      # Constrain Pnu > P_free
      ELIND=1

OPTIONS
      ASA[ CCOR={ccor} TWOC=F ADNF=F GAMMA={gamma}] SCR={scr} SX={sx}
STR   RMAXS={rmaxs}                # Range of strux

GF    MODE={gfmode} GFOPTS={?~c3~p3;~p2;}
BZ    EMESH={nz},10,-1,{ef},.5,.3  # nz-pts;contour mode;emin;emax;ecc;bunching

.SPEC
# SCLWSR=11 WSRMAX=3.3 OMAX1=0.18 0.26 0.26
  ATOM=Nb       Z= 41  R= 3.065511  LMX=2  IDXDN=0,0,0,2 # IDMOD=0,1
  ATOM=Nbx      Z= 41  R= 2.730129  LMX=2  IDXDN=0,0,0,2 # IDMOD=0,1
  ATOM=Fex      Z= 26  R= 2.496030  LMX=2
  ATOM=Fe       Z= 26  R= 2.719829  LMX=2
  ATOM=E        Z=  0  R= 1.715533  LMX=0
  ATOM=E1       Z=  0  R= 1.710466  LMX=0
  ATOM=E2       Z=  0  R= 1.651325  LMX=0
  ATOM=E3       Z=  0  R= 1.202364  LMX=0

SPEC
# SCLWSR=11 WSRMAX=3.3 OMAX1=0.18 0.26 0.26
  ATOM=Nb       Z= 41  R= 3.065511  LMX=2  IDXDN=0,0,0,2  IDMOD=0,1
  ATOM=Nbx      Z= 41  R= 2.683855  LMX=2  IDXDN=0,0,0,2  IDMOD=0,1
  ATOM=Fex      Z= 26  R= 2.453724  LMX=2  MMOM=0,0,2
  ATOM=Fe       Z= 26  R= 2.673730  LMX=2  MMOM=0,0,2
  ATOM=E        Z=  0  R= 1.715533  LMX=0
  ATOM=E1       Z=  0  R= 1.796268  LMX=0
  ATOM=E2       Z=  0  R= 1.737788  LMX=0
  ATOM=E3       Z=  0  R= 1.202364  LMX=0

START CNTROL={nit<=0} BEGMOM={nit<=0}

site.asa:

% site-data vn=3.0 fast io=15 nbas=60 alat=6.226 plat= -0.5 0.5 0.5 1.5 1.5 -1.5 0.55275066 -1.51014288 6.06289351
#                        pos
 Nb        0.0000000   0.0000000   0.0000000
 Nb        0.5000000   0.5000000  -0.5000000
 Nb        1.0000000   1.0000000  -1.0000000
 E3       -0.4737501  -0.7394575   1.0371136
 E3        0.5357464   0.7459879   0.0428640
 E3        0.0444594   0.2488781   0.5344094
 E3       -0.4454942  -0.2743963   1.0326326
 E3        0.5547866   0.2622121   0.0406214
 E3        0.0405665  -0.5050013   0.7628725
 E3        0.5391316  -0.0233702   0.2705096
 E3       -0.4670669  -0.9928333   1.2831687
 E3        0.2968868  -0.0115608   0.5300687
 E3        0.8001271   0.5039323   0.0363901
 Nbx       0.0038386  -0.0166966   0.9910505
 Nbx       0.9928133   1.0176382   0.0048461
 Nbx       0.5161225   0.4990571   0.5150855
 E3       -0.0164014  -0.4956354   1.2589020
 E3       -0.3739799  -0.9880550   1.6531767
 E3        0.7909603   0.0066238   0.5321685
 E3        0.4664294  -0.3020557   0.9890124
 E3       -0.0259534  -0.7485791   1.4813952
 E3        0.5611756   0.1357335   0.9604225
 E3        0.0229530  -0.2263954   1.4986452
 E2       -0.5404845  -0.6328016   2.0620827
 E1       -0.4216957  -1.1627362   1.9898101
 Fex      -0.0399873  -0.0010873   1.9603949
 Fex      -0.0259839   0.7809822   1.9780095
 Fex       0.9792297   1.2189900   0.9795957
 Fex       0.5160204   0.9996464   1.5211879
 Fe        0.1509317   0.2442758   2.6493636
 Fe        0.1626577  -0.5019115   2.6567494
 Fe        0.6537761   0.4947430   2.1670347
 Fe        1.1725712   0.7480391   1.6834229
 Fex       0.8421375  -0.0133975   2.8192889
 Fex       0.3308379  -0.2533275   3.3370495
 Fex       1.3371460   0.2434396   2.3440138
 Fex       0.3278900  -0.9987992   3.3586058
 E3       -0.3619496  -1.1309320   4.2693144
 E         0.8837140  -0.5118108   3.3420287
 E3       -0.0787140  -1.4964656   4.3085914
 E3        0.6610712  -0.7050789   3.5688063
 E3        0.2767879  -1.0210851   4.0233808
 E3        0.5795804  -1.0270209   3.8087107
 Nbx       1.5400241  -0.4828374   3.0627345
 Nbx       0.5490728  -1.5111615   4.0608585
 Nbx       1.0623082  -1.0379833   3.5524052
 E3        1.0240526  -0.5127517   3.7389527
 E3        0.5165514  -1.0225847   4.2903860
 E3        0.0116892  -1.4962605   4.7978325
 E3        0.3010275  -1.4799238   4.5283927
 E3        1.2996213  -0.5585176   3.5297990
 E3        1.0168605  -0.2808990   4.0027595
 E3       -0.0005736  -1.2373961   5.0284631
 E3       -0.0037828  -1.7497071   5.0358071
 E3        0.5024477  -0.7478516   4.5295767
 E3        0.5170890  -1.2513967   4.5213959
 E3        1.0304754  -0.7795029   4.0211892
 Nb        0.5527507  -1.5101429   5.0628935
 Nb        1.0527507  -1.0101429   4.5628935
 Nb        1.5527507  -0.5101429   4.0628935

Remove any atom files that might be present and Make the system self-consistent:

touch nb.asa
rm {nb,fe,e}*.asa
lm asa -vnit=0 --pr30,20 --rs=0,1
lm asa --rsedit~rs~rs@ib=1@src=1:3@fn=rsta2~rs@ib=14@src=1:3@fn=rsta2~~rs@ib=44@src=1:3@fn=rsta2~rs@ib=58@src=1:3@fn=rsta2~rs@ib=1@src=1:3@fn=rsta2~rs@ib=26@src=1:4@fn=rsta3~rs@ib=30@src=1:4@fn=rsta3~rs@ib=34@src=1:4@fn=rsta3
c=1:4@fn=rsta4
cp rsta4.asa rsta.asa
lm asa -vnit=0 --pr30,20 --rs=1,0
lmstr asa
rm -f save.asa mixm.asa sv.asa log.asa
mpix -n 16 lm -vmet=3 ctrl.asa -vccor=f -vnit=400 --iactiv -vscr=4 --zerq~all~qout --pr30,21 > out  &

Convergence to self-consistency is slow, as a consequence of the screening algorithm. This is necessary to protect against low-energy, long-wave charge oscillations along the long axis. That this is an important issue can readily be seen from Poisson’s equation of a periodic cell in Fourier space:

Away from self-consistency there will be a deviation . This generates a change in the potential . Long cells have are vectors of low G, and changes to the potential from get greatly amplified. The solution is to screen out these low Fourier components using a model dielectric function. But the ASA has no plane-wave representation of the density, so constructing a model dielectric function is not so simple.

A model construction is nevertheless implemented, and it protects against these fluctuations, but at the same time it conceals true physical changes in them. This limits the rate of convergence. It can be dealt with by some extra tricks, but for the present we simply iterate a very large number of times to circumvent the problem.

The self-consistent magnetic moment is close to 24 μB, or about 2 μB per Fe atom (inspect the last line of save.asa). This is close to the full LDA result, 23.6 μB, and is a bit less than the moment of 12 Fe atoms in bulk, 12×2.2 = 26.4 μB. It shows that the Fe moment isn’t quite rigid, and gets extinguished somewhat at the Fe/Nb interface.

As a last step run lm one final iteration, to save in file rsta.asa essential information that can reconstruct the self-consistent density:

mpix -n 16 lm -vmet=3 ctrl.asa -vccor=f -vnit=1 --iactiv -vscr=4 --pr40,21 --rs=0,1 --quit=rho > out.asa

That rsta.asa contains sufficient information to reconstruct the density, remove the atom files containing potentials and do

rm {nb,fe,e}*.asa
lm asa -vnit=0 --pr30,20 --rs=1,0
mpix -n 16 lm -vmet=3 ctrl.asa -vccor=f -vnit=1 --iactiv -vscr=4 --pr40,21 --quit=rho > out

Files out and out.asa are essentially identical.

ASA Energy bands

Make the energy bands in a manner similar to the full LDA. As a first step, extract the a list of orbitals in the Hamiltonian associated with Nb, and ones with Fe, so as to color the bands according to the Nb or Fe character.

lm asa --pr55 --quit=ham | grep -A 70 'Orbital positions'

Inspect the the table. We can associate orbitals 1:64 and 188:244 with Nb, and orbitals 74:172 with Fe.

Generate the energy bands

cp ../syml.nbfe syml.asa
lm -vmet=3 ctrl.asa -vccor=f -vnit=1 --iactiv -vscr=4 --pr40,21 --quit=rho --band~col=1:64,188:244~col2=74:172~fn=syml

and create and view postscript figure for the majority and minority bands

echo -5,5,5,10 | plbnds -spin1 -fplot~sh~scl=.7~shft=-.28 -ef=0 -scl=13.6 -dat=up -lt=1,col=0,0,0,colw=1,0,0,colw2=0,1,0 bnds.asa
open fplot.ps
echo -5,5,5,10 | plbnds -spin2 -fplot~sh~scl=.7~shft=.45 -ef=0 -scl=13.6 -dat=dn -lt=1,col=0,0,0,colw=1,0,0,colw2=0,1,0 bnds.asa
open fplot.ps

Majority (left) and minority (right) energy bands of the Nb/Fe superlattice.   Red: Nb character; green: Fe character

ASA energy bands of the Fe/Nb interface

Compare to the full LDA band structure. You can see that they are very similar.

Crystal Green’s function calculation of the Nb/Fe superlattice

The ASA suite also has a crystal Green’s function code, lmgf. This page has a basic tutorial describing lmgf and its operation. lmgf is of particular interest here mainly as a segue into the Principal layer Green’s function code, lmpg.

lmgf operates in a manner approximately similar to lm but it generates Green’s functions instead of wave functions. lmgf reads the same potential parameters as lm but constructs the Green’s function instead of a hamiltonian. The program flow is similar to lm, only the “crystal” branch is replaced by a branch that obtains multipole moments through Green’s functions. Multipole moments provide a compact way to represent the charge density in the ASA: the ability to generate them enables lmgf to perform self-consistent calculations.

As noted in the Introductory tutorial, lmgf requires some additional input, namely a GF category and an EMESH tag in the BZ category. It is strongly advised that you run through that tutorial before doing this section.

The elliptical energy contour, normally used for Brillouin zone integrations needed, e.g. to find the Fermi level and the output moments, requires information about the Fermi level. It is easiest to get an estimate from lm, since the two methods are closely linked. From the output of lm, the Fermi level is displayed in a table with a line

 BZINTS: Fermi energy:      0.154488; 156.000000 electrons;  D(Ef):  323.085
 ...

The Fermi level appropriate to lmgf should be close to 0.154 Ry, but not identically so because lmgf uses a different integration scheme and the third order Green’s function is close to, but identical with the 3rd order ASA hamiltonian.

Following the Introductory tutorial run lmgf to get the Fermi level satisfying charge neutrality

lmgf ctrl.asa -vccor=f -vnit=1 -vscr=4 --pr40,21 --quit=rho -vgamma=t -vef=.154

Inspect file vshft.asa. It found a constant potential shift, vconst=-0.0004575 must be subtracted from the estimate for the Fermi level. The Fermi level then should be ef-vconst=0.1544 Ry, perhaps fortuitously close to the Fermi level generated by lm. The close agreement should give use some confidence that the Green’s function is well described.

Spectral Function from lmgf

A Green’s functions does not supply quasiparticle levels, but it does contain the spectral function . The analog of the band structure can be got by drawing as a “heat map” on a grid. Since in this case there are no alloy or many-body scattering processes, the spectral function is fully coherent with unit Z factor and sharp peaks at the QP levels. They show up as bright spots on a heat map.

lmgf will make data for a heat map in lieu of the band structure, using the --band switch. By passing the Fermi level within the --band switch, lmgf will shift the spectral function to make EF coincide with 0. This is useful for comparison to the band structures drawn previously.

mpix -n 16 lmgf -vmet=3 ctrl.asa -vccor=f -vnit=1 -vscr=4 --pr40,21 --quit=rho -vgamma=t  -vef=.1544 --band~ef=.1544~fn=syml2

Draw a figure with

~/lm/utils/SpectralFunction.sh  -w -5:5
open asa_UP.png
open asa_DN.png

Majority (left) and minority (right) spectral functions bands of the Nb/Fe superlattice

NbFe SF up figure.

NbFe SF dn figure.

Spectral functions correspond well to the energy bands generated by lm, Note in particular the narrow Fe d bands are bunched together and show greater intensity (red). They correspond to the Fe d bands, seen as nearly flat bands in green near the Fermi level.

Nb/Fe/Nb Trilayer: lmpg

… to be completed

Footnotes and references