[gpaw-users] Trouble with spin-polarized calculation convergence NO molecule

Michael Walter Michael.Walter at fmf.uni-freiburg.de
Fri May 18 09:39:13 CEST 2012


Dear Andre,

it works quite nicely like this:

------------------------------------------------------
import sys
from gpaw.cluster import Cluster

# get initial staructure and move to the box center
import os
import sys

from ase import Atom, optimize

from gpaw import GPAW
from gpaw.cluster import Cluster

h = 0.25
R0 = 1.

ext = ''
fname='PBE' + ext + '.gpw'
ftraj='relax' + ext + '.traj'

try:
   s = Cluster(ftraj)
except:
   s = Cluster([Atom('N'), Atom('O',[0,0,R0])])
s.minimal_box(3., h=h)

c = GPAW(xc='PBE', h=h, nbands=-8, width=0.1)
s.set_calculator(c)

# Make a trajectory file:
dyn = optimize.FIRE(s)
dyn.run(fmax=0.05)
------------------------------------------------------

It might be that this is not what you want, though. The occupations are
(last step):

Fermi Level: -4.06658
 Band   Eigenvalues  Occupancy
   0    -31.02655     2.00000
   1    -16.55047     2.00000
   2    -12.48119     2.00000
   3    -12.48119     2.00000
   4    -11.44573     2.00000
   5     -3.95671     0.50000
   6     -3.95671     0.50000
   7      2.58220     0.00000
   8      3.84096     0.00000
   9      4.49919     0.00000
  10      4.49919     0.00000
  11      6.15199     0.00000
  12      6.38558     0.00000

In case you want to occupy state 5 with one electron, you have to break the
symmetry somehow. It's a radical, so spin polarization can lower the energy.

Best,
Michael

2012/5/17 Clayborne, Penee A. <clayborne at anl.gov>

> Hello.  I am trying a spin polarized calculation for the NO molecule.
> However, after numerous attempts, I am having difficulty with the
> convergence.  Any recommendations?
>
> from ase import *
> from ase.all import *
> from ase.optimize import QuasiNewton
> from gpaw import GPAW, FermiDirac
> from ase import Atoms
> from ase.io import read, write, PickleTrajectory
> from gpaw import *
> import numpy
> import numpy as np
>
> atoms = read('NO.xyz')
> atoms.center(vacuum=6.0)
>
> calc = GPAW(h=0.18,
>            maxiter=450,
>            occupations=FermiDirac(0.05),
>            spinpol=True,
>            usesymm=False,
>            charge=00.0,
>            nbands=-20,
>            xc='PBE',
>            mixer=Mixer(0.2, 5),
>            convergence={'energy':0.05, 'density':1.0e-5,
> 'eigenstates':1.0e-7},
>           txt='NO.out')
>
> atoms.set_calculator(calc)
>
> traj = PickleTrajectory('NO.traj', 'w', atoms)
>
> calc.set(convergence={'bands':-20})
> e = atoms.get_potential_energy()
>
> relax = QuasiNewton(atoms, trajectory='NO.traj')
> relax.attach(traj.write)
> relax.run(fmax=0.05)
>
> write('relax_NO.xyz',atoms)
> calc.write('NO.gpw', mode='all')
>
> Thanks in advance!
>
> P. Andre Clayborne, Ph.D.
> Postdoctoral Appointee
> Center for Nanoscale Materials
> Argonne National Laboratory, Bldg. 440
> Argonne, IL 60439
>
>
> _______________________________________________
> gpaw-users mailing list
> gpaw-users at listserv.fysik.dtu.dk
> https://listserv.fysik.dtu.dk/mailman/listinfo/gpaw-users
>



-- 
------------------------------------------
PD Dr Michael Walter
Address: Freiburger Materialforschungszentrum
         Stefan-Meier-Straße 21
         D-79104 Freiburg i. Br.
         Germany
Tel.: +49 761 203 4758 and +49 761 203 7695
Fax: +49 761 203 4701
email: Michael.Walter at fmf.uni-freiburg.de
www: http://omnibus.uni-freiburg.de/~mw767
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://listserv.fysik.dtu.dk/pipermail/gpaw-users/attachments/20120518/39e2a056/attachment-0001.html 


More information about the gpaw-users mailing list