Optimization¶
import some modules¶
[1]:
from dftpy.ions import Ions
from dftpy.field import DirectField
from dftpy.grid import DirectGrid
from dftpy.functional import LocalPseudo, Functional, TotalFunctional
from dftpy.formats import io
from dftpy.math_utils import ecut2nr
from dftpy.time_data import TimeData
from dftpy.optimization import Optimization
from dftpy.mpi import sprint
pseudopotential file¶
[2]:
path_pp='../DATA/'
file1='Al_lda.oe01.recpot'
PP_list = {'Al': path_pp+file1}
build the ions or read from file¶
[3]:
from ase.build import bulk
atoms = bulk('Al', 'fcc', a=4.05, cubic=True)
ions = Ions.from_ase(atoms)
# ions = io.read(posfile)
make a grid¶
[4]:
nr = ecut2nr(ecut=35, lattice=ions.cell)
grid = DirectGrid(lattice=ions.cell, nr=nr)
sprint('The final grid size is ', nr)
The final grid size is [20 20 20]
build local pseudo, and generate guess density¶
[5]:
PSEUDO = LocalPseudo(grid = grid, ions=ions, PP_list=PP_list)
rho_ini = DirectField(grid=grid)
rho_ini[:] = ions.get_ncharges()/ions.cell.volume
setting key: Al -> ../DATA/Al_lda.oe01.recpot
instance KEDF, XC and HARTREE functionals¶
[6]:
KE = Functional(type='KEDF',name='TFvW')
XC = Functional(type='XC',name='LDA')
HARTREE = Functional(type='HARTREE')
instance DFTpy evaluator¶
[7]:
evaluator = TotalFunctional(KE=KE, XC=XC, HARTREE=HARTREE, PSEUDO=PSEUDO)
instance and execute DFTpy density optimizer¶
[8]:
optimization_options = {'econv' : 1e-6*ions.nat}
opt = Optimization(EnergyEvaluator=evaluator, optimization_options = optimization_options,
optimization_method = 'TN')
rho = opt.optimize_rho(guess_rho=rho_ini)
Step Energy(a.u.) dE dP Nd Nls Time(s)
0 -8.090977710718E+00 -8.090978E+00 7.877088E-01 1 1 1.902819E-02
1 -8.273226052119E+00 -1.822483E-01 7.033208E-02 2 1 3.294325E-02
2 -8.280858141000E+00 -7.632089E-03 4.803522E-03 7 1 5.704927E-02
3 -8.281101062862E+00 -2.429219E-04 3.640550E-04 5 1 7.545233E-02
4 -8.281135978745E+00 -3.491588E-05 3.590863E-05 6 1 9.528708E-02
5 -8.281138625688E+00 -2.646943E-06 2.500624E-06 5 1 1.117182E-01
6 -8.281138996216E+00 -3.705278E-07 4.567996E-08 8 1 1.352971E-01
#### Density Optimization Converged ####
Chemical potential (a.u.): 0.30115196532053357
Chemical potential (eV) : 8.194762380333897
evaluate final energy¶
[9]:
energy = evaluator.Energy(rho=rho, ions=ions)
print('Energy (a.u.)', energy)
Energy (a.u.) -8.2811389962163
print the timing¶
[10]:
TimeData.output(lprint=True, sort='cost')
--------------------------------Time information--------------------------------
Label Cost(s) Number Avg. Cost(s)
ewald.Energy_corr 0.0001 1 0.0001
CBspline._calc_PME_Qarray 0.0016 1 0.0016
ewald.Energy_rec_PME 0.0025 1 0.0025
LocalPseudo.local_PP 0.0045 1 0.0045
ewald.Energy_real_fast2 0.0054 1 0.0054
TF 0.0066 41 0.0002
ewald.energy 0.0080 41 0.0002
LDA 0.0125 41 0.0003
FFT 0.0320 84 0.0004
InvFFT 0.0359 83 0.0004
Hartree.compute 0.0364 41 0.0009
vW 0.0458 41 0.0011
Optimize 0.1377 1 0.1377
[11]:
rho.write('rho.xsf', ions=ions)
rho.write('rho.cube', ions=ions)
Visualize with scikit-image and matplotlib¶
!pip install scikit-image matplotlib
[12]:
from dftpy.visualize import view
Visualize with VESTA¶
[15]:
view(ions=ions, data=rho, viewer='vesta')
save ./.dftpy.xsf in darwin platform