[gpaw-users] Funny electronic band structure

Jens Jørgen Mortensen jensj at fysik.dtu.dk
Wed Feb 4 14:12:09 CET 2015


On 02/03/2015 10:59 AM, Darío Fdez-Pello Lois wrote:
> Dear all,
>
> Thanks for your replies, indeed it should be a colon after 'size', I 
> mistyped it somehow.
>
> The requested figures are attached, please let me say that Mo_band.png 
> comes from the QE and ELK calculations while MoBandstructure.png comes 
> from GPAW.

It looks like you have a cubic cell with two atoms.  Try with the 
primitive unit cell:

 >>> from ase.lattice import bulk
 >>> mo = bulk('Mo')

See here:

https://wiki.fysik.dtu.dk/ase/ase/structure.html#ase.lattice.bulk

Jens Jørgen

>
> In addition, here are the scripts:
> -------SCF
> import numpy as NP
> from ase.lattice.cubic import BodyCenteredCubic
> from gpaw import GPAW,PW,MethfesselPaxton,Mixer
> from ase.io <http://ase.io> import read
>
> k = 10
> ecut = 600
> N = 1
> configs = read('MoStrain.traj@:')
> aa = NP.array([config.cell[0, 0] for config in configs])
> a = aa[-1]
>
> mo = BodyCenteredCubic(size=(N,N,N),symbol='Mo',latticeconstant=a)
> mo.calc = GPAW(mode=PW(ecut), xc='PBE', kpts={'size': (k,k,k),'gamma': 
> True}, setups={'Mo': '6'}, occupations=MethfesselPaxton(0.01), 
> mixer=Mixer(0.1, 1, 100), txt='bulk.txt')
> mo.get_potential_energy()
> mo.calc.write('bulk.gpw')
>
> --------BAND
> import numpy as NP
> from ase.io <http://ase.io> import read, write
> from ase.lattice.cubic import BodyCenteredCubic
> from ase.dft.kpoints import ibz_points, get_bandpath
> from gpaw import GPAW, PW, MethfesselPaxton
> import matplotlib
> matplotlib.use('Agg')
> import matplotlib.pyplot as plt
>
> num_bands=12
>
> calc = GPAW('bulk.gpw', fixdensity=True, usesymm=None, 
> convergence={'bands': num_bands}, maxiter=600, txt='band.txt')
> points = ibz_points['bcc']
> G = points['Gamma']
> H = points['H']
> N = points['N']
> P = points['P']
> kpts, x, X = get_bandpath([G, H, N, G, P], calc.atoms.cell, 200)
> calc.set(kpts=kpts)
> calc.get_potential_energy()
> e_kn = NP.array([calc.get_eigenvalues(k) for k in range(len(kpts))])
>
> # Plot the band structure
>
> e_f = calc.get_fermi_level()
> e_kn -= e_f
> emin = -10 #e_kn.min() - 1.0
> emax = 10 #e_kn[:, num_bands].max() + 1.0
>
> plt.figure(figsize=(5, 6))
> for n in range(num_bands):
>     plt.plot(x, e_kn[:, n])
> for p in X:
>     plt.plot([p, p], [emin, emax], 'k-')
> plt.plot([0, X[-1]], [0, 0], 'k-')
> plt.xticks(X, ['$%s$' % n for n in ['\Gamma', 'H', 'N', '\Gamma', 'P']])
> plt.axis(xmin=0, xmax=X[-1], ymin=emin, ymax=emax)
> plt.xlabel('k-vector')
> plt.ylabel('E - E$_F$ (eV)')
> plt.title('Bandstructure of Mo')
> plt.savefig('MoBandstructure.png')
>
> --------------------
> ​END of scripts​
>
>
> ​Best regards,
>
>    Darío​
>
>
> 2015-02-03 8:19 GMT+01:00 Jussi Enkovaara <jussi.enkovaara at csc.fi 
> <mailto:jussi.enkovaara at csc.fi>>:
>
>     Hi Dario,
>     could you be a bit more specific what's wrong with the band
>     stucture (preferrably with a picture of correct / incorrect one),
>     e.g are some bands fine but others not?
>
>     In principle, your script below looks fine, but could you also
>     provide the full scripts for SCF and band structure calculations?
>
>     Best regards,
>     Jussi Enkovaara
>
>
>     On 2015-02-03 00:38, Darío Fdez-Pello Lois wrote:
>
>         Dear all,
>
>         My question is related to the calculation of the electronic
>         band structure
>         of molybdenum, which I was successfully able to reproduce both
>         using a full
>         potential code (ELK) and a plane wave one (Quantum Espesso).
>         Unfortunately,
>         the same doesn't happen when I use GPAW (with PW mode)
>         resulting a quite
>         different structure.
>         Here are the relevant (I guess) lines, of course the
>         obligatory parameter's
>         convergence were peviouslly done:
>
>         [ ... ]
>         mo=BodyCenteredCubic(size(N,N,N),symbol='Mo',latticeconstant=a)
>         mo.calc = GPAW(mode=PW(ecut),
>                xc='PBE',
>                kpts={'size'(k,k,k),'gamma': True},
>                setups={'Mo': '6'},
>                occupations=MethfesselPaxton(0.01),
>                mixer=Mixer(0.1, 1, 100),
>                txt='bulk.txt')
>         [ ... ]
>
>         thereafter I run the band structure script similar to the one
>         given in the
>         tutorial (with obvious changes: bcc instead of fcc, etc.)
>
>         I am very confused at this point and I hope some of you would
>         kindly give
>         me some insight.
>
>         Best,
>
>            Darío.
>
>
>
>
>         _______________________________________________
>         gpaw-users mailing list
>         gpaw-users at listserv.fysik.dtu.dk
>         <mailto:gpaw-users at listserv.fysik.dtu.dk>
>         https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
>
>
>
>     _______________________________________________
>     gpaw-users mailing list
>     gpaw-users at listserv.fysik.dtu.dk
>     <mailto:gpaw-users at listserv.fysik.dtu.dk>
>     https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
>
>
>
>
> _______________________________________________
> gpaw-users mailing list
> gpaw-users at listserv.fysik.dtu.dk
> https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users



More information about the gpaw-users mailing list