SST Lab Dokuwiki Header header picture

ユーザ用ツール

サイト用ツール


info:dacapo

Jacapo

Python interface to the Fortran DACAPO code

Jacapo Jacapo

functions

read(ncfile)

read(ncfile) return atoms and calculator from ncfile

>>> atoms, calc = read('co.nc')

class

class Jacapo
  Python interface to the Fortran DACAPO code

Constractor

Jacapo

__init__(self, nc='out.nc', outnc=None, debug=30, stay_alive=False, **kwargs)

Parameters

Name Type Description
nc string output netcdf file, or input file if nc already exists
outnc string output file. by default equal to nc
debug integer logging debug level

Valid kwargs

Name Type Description
atoms ASE.Atoms instance atoms is an ase.Atoms object that will be attached to this calculator.
pw integer sets planewave cutoff
dw integer sets density cutoff
kpts iterable set chadi-cohen, monkhorst-pack kpt grid,
e.g. kpts = (2,2,1) or explicit list of kpts
spinpol Boolean sets whether spin-polarization is used or not.
fixmagmom float set the magnetic moment of the unit cell. only used in spin polarize calculations
ftfloatset the Fermi temperature used in occupation smearing
xcstringset the exchange-correlation functional. one of ['PZ','VWN','PW91','PBE','RPBE','revPBE'],
dipolebooleanturn the dipole correction on (True) or off (False)
dictionary dictionary of parameters to fine-tune behavior
  {'status':False,
   'mixpar':0.2,
   'initval':0.0,
   'adddipfield':0.0,
   'position':None}
nbandsintegerset the number of bands
symmetryBooleanTurn symmetry reduction on (True) or off (False)
stressBooleanTurn stress calculation on (True) or off (False)
debug level for logging could be something like logging.DEBUG or an integer 0-50. The higher the integer, the less information you see set debug level (0 = off, 10 = extreme)

Initialize the Jacapo calculator Modification of the nc file only occurs at calculate time if needed

>>> calc = Jacapo('CO.nc')

reads the calculator from CO.nc if it exists or minimally initializes CO.nc with dimensions if it does not exist.

>>> calc = Jacapo('CO.nc', pw=300)

reads the calculator from CO.nc or initializes it if it does not exist and changes the planewave cutoff energy to 300eV

 >>> atoms = Jacapo.read_atoms('CO.nc')

returns the atoms in the netcdffile CO.nc, with the calculator attached to it.

>>> atoms, calc = read('CO.nc')

__str__

str(self) pretty-print the calculator and atoms.

we read everything directly from the ncfile to prevent triggering any calculations

atoms_are_equal

atoms_are_equal(self, atoms) comparison of atoms to self.atoms using tolerances to account for float/double differences and float math.

attach_child

attach_child(self, child)

calculate

calculate(self) run a calculation.

you have to be a little careful with code in here. Use the calculation_required function to tell if a calculation is required. It is assumed here that if you call this, you mean it.

calculation_required

calculation_required(self, atoms=None, quantities=None) determines if a calculation is needed.

return True if a calculation is needed to get up to date data. return False if no calculation is needed.

quantities is here because of the ase interface. ==== delete_ncattdimvar —- delete_ncattdimvar(self, ncf, ncattrs=None, ncdims=None, ncvars=None) helper function to delete attributes, dimensions and variables in a netcdffile

this functionality is not implemented for some reason in netcdf, so the only way to do this is to copy all the attributes, dimensions, and variables to a new file, excluding the ones you want to delete and then rename the new file.

if you delete a dimension, all variables with that dimension are also deleted.

execute_external_dynamics

execute_external_dynamics(self, nc=None, txt=None, stoppfile='stop', stopprogram=None) Implementation of the stay alive functionality with socket communication between dacapo and python. Known limitations: It is not possible to start 2 independent Dacapo calculators from the same python process, since the python PID is used as identifier for the script[PID].py file.

execute_parent_calculation

execute_parent_calculation(self) Implementation of an extra level of parallelization, where one jacapo calculator spawns several dacapo.run processes. This is used for NEBs parallelized over images.

get_ados

get_ados(self, **kwargs)

attempt at maintaining backward compatibility with get_ados returning data

Now when we call calc.get_ados() it will return settings,

and calc.get_ados(atoms=[],…) should return data

get_ados_data

get_ados_data(self, atoms, orbitals, cutoff, spin) get atom projected data

:Parameters:

atoms

    list of atom indices (integers)

orbitals

    list of strings
    ['s','p','d'],
    ['px','py','pz']
    ['d_zz', 'dxx-yy', 'd_xy', 'd_xz', 'd_yz']

cutoff : string

  cutoff radius you want the results for 'short' or 'infinite'

spin

  : list of integers
  spin you want the results for
  [0] or [1] or [0,1] for both

returns (egrid, ados) egrid has the fermi level at 0 eV

get_all_eigenvalues

get_all_eigenvalues(self, spin=0) return all the eigenvalues at all the kpoints for a spin.

:Parameters:

spin : integer

  which spin the eigenvalues are for

get_ascii_debug

get_ascii_debug(self) Return the debug settings in Dacapo

get_atoms

get_atoms(self) return the atoms attached to a calculator()

get_bz_k_points

get_bz_k_points(self) return list of kpoints in the Brillouin zone

get_calculate_stress

get_calculate_stress(self) return whether stress is calculated or not

get_cd

get_cd = get_charge_density(self, spin=0)

get_charge_density

get_charge_density(self, spin=0) return x,y,z,charge density data

x,y,z are grids sampling the unit cell cd is the charge density data

netcdf documentation::

ChargeDensity(number_of_spin,

              hardgrid_dim3,
              hardgrid_dim2,
              hardgrid_dim1)
ChargeDensity:Description = "realspace charge density" ;
        ChargeDensity:unit = "-e/A^3" ;

get_charge_mixing

get_charge_mixing(self) return charge mixing parameters

get_convergence

get_convergence(self) return convergence settings for Dacapo

get_debug

get_debug(self) Return the python logging level

get_decoupling

get_decoupling(self) return the electrostatic decoupling parameters

get_dipole

get_dipole(self) return dictionary of parameters if the DipoleCorrection was used

get_dipole_moment

get_dipole_moment(self, atoms=None) return dipole moment of unit cell

Defined by the vector connecting the center of electron charge density to the center of nuclear charge density.

Units = eV*angstrom

1 Debye = 0.208194 eV*angstrom get_dw(self) return the density wave cutoff get_ef = get_fermi_level(self) get_effective_potential(self, spin=1) returns the realspace local effective potential for the spin. the units of the potential are eV

:Parameters:

spin : integer

   specify which spin you want, 0 or 1

get_eigenvalues(self, kpt=0, spin=0) return the eigenvalues for a kpt and spin

:Parameters:

kpt : integer

  index of the IBZ kpoint

spin : integer

  which spin the eigenvalues are for

get_electronic_minimization(self) get method and diagonalizations per band for electronic minimization algorithms get_electronic_temperature = get_ft(self)

get_electrostatic_potential

get_electrostatic_potential(self, spin=0) get electrostatic potential

Netcdf documentation::

double ElectrostaticPotential(number_of_spin,

                              hardgrid_dim3,
                              hardgrid_dim2,
                              hardgrid_dim1) ;
       ElectrostaticPotential:
           Description = "realspace local effective potential" ;
           unit = "eV" ;

get_ensemble_coefficients

get_ensemble_coefficients(self) returns exchange correlation ensemble coefficients

get_esp

get_esp = get_electrostatic_potential(self, spin=0) get_external_dipole(self) return the External dipole settings

get_extpot

get_extpot(self) return the external potential set in teh calculator

get_extracharge

get_extracharge(self) Return the extra charge set in teh calculator

get_fermi_level

get_fermi_level(self) return Fermi level

get_fftgrid

get_fftgrid(self) return soft and hard fft grids

get_fixmagmom

get_fixmagmom(self) returns the value of FixedMagneticMoment

get_forces

get_forces(self, atoms=None) Calculate atomic forces

get_ft

get_ft(self) return the FermiTemperature used in the calculation

get_ibz_k_points

get_ibz_k_points = get_ibz_kpoints(self)

get_ibz_kpoints

get_ibz_kpoints(self) return list of kpoints in the irreducible brillouin zone

get_k_point_weights

get_k_point_weights(self) return the weights on the IBZ kpoints

get_kpts

get_kpts(self) return the BZ kpts

get_kpts_type

get_kpts_type(self) return the kpt grid type

get_magnetic_moment

get_magnetic_moment(self, atoms=None) calculates the magnetic moment (Bohr-magnetons) of the supercell

get_magnetic_moments

get_magnetic_moments(self, atoms=None) return magnetic moments on each atom after the calculation is run

get_mdos

get_mdos(self) return multicentered projected dos parameters

get_mdos_data

get_mdos_data(self, spin=0, cutoffradius='infinite') returns data from multicentered projection

returns (mdos, rotmat)

the rotation matrices are retrieved from the text file. I am not sure what you would do with these, but there was a note about them in the old documentation so I put the code to retrieve them here. the syntax for the return value is: rotmat[atom#][label] returns the rotation matrix for the center on the atom# for label. I do not not know what the label refers to.

get_nbands

get_nbands(self) return the number of bands used in the calculation get_nc(self) return the ncfile used for output get_ncoutput(self) returns the control variables for the ncfile get_number_of_bands = get_nbands(self) get_number_of_electrons = get_valence(self, atoms=None) get_number_of_grid_points(self) return soft fft grid get_number_of_iterations(self) get_number_of_spins(self) if spin-polarized returns 2, if not returns 1 get_occ = get_occupation_numbers(self, kpt=0, spin=0) get_occupation_numbers(self, kpt=0, spin=0) return occupancies of eigenstates for a kpt and spin

:Parameters:

kpt : integer

  index of the IBZ kpoint you want the occupation of

spin : integer

  0 or 1

get_occupationstatistics(self) return occupation statistics method get_potential_energy(self, atoms=None, force_consistent=False) return the potential energy. get_pseudo_wave_function(self, band=0, kpt=0, spin=0, pad=True) return the pseudo wavefunction get_pseudopotentials(self) get pseudopotentials set for atoms attached to calculator get_psp(self, sym=None, z=None) get the pseudopotential filename from the psp database

:Parameters:

sym : string

  the chemical symbol of the species

z : integer

  the atomic number of the species

you can only specify sym or z. Returns the pseudopotential filename, not the full path. get_psp_nuclear_charge(self, psp) get the nuclear charge of the atom from teh psp-file.

This is not the same as the atomic number, nor is it necessarily the negative of the number of valence electrons, since a psp may be an ion. this function is needed to compute centers of ion charge for the dipole moment calculation.

We read in the valence ion configuration from the psp file and add up the charges in each shell. get_psp_valence(self, psp) get the psp valence charge on an atom from the pspfile. get_pw(self) return the planewave cutoff used get_reciprocal_bloch_function(self, band=0, kpt=0, spin=0) return the reciprocal bloch function. Need for Jacapo Wannier class. get_reciprocal_fft_index(self, kpt=0) return the Wave Function FFT Index get_scratch(self) finds an appropriate scratch directory for the calculation get_spin_polarized(self) Return True if calculate is spin-polarized or False if not get_spinpol(self) Returns the spin polarization setting, either True or False get_status(self) get status of calculation from ncfile. usually one of: 'new', 'aborted' 'running' 'finished' None get_stay_alive(self) return the stay alive settings get_stress(self, atoms=None) get stress on the atoms.

you should have set up the calculation to calculate stress first.

returns [sxx, syy, szz, syz, sxz, sxy] get_symmetry(self) return the type of symmetry used get_txt(self) return the txt file used for output get_ucgrid(self, dims) Return X,Y,Z grids for uniform sampling of the unit cell

dims = (n0,n1,n2)

n0 points along unitcell vector 0 n1 points along unitcell vector 1 n2 points along unitcell vector 2 get_valence(self, atoms=None) return the total number of valence electrons for the atoms. valence electrons are read directly from the pseudopotentials.

the psp filenames are stored in the ncfile. They may be just the name of the file, in which case the psp may exist in the same directory as the ncfile, or in $DACAPOPATH, or the psp may be defined by an absolute or relative path. This function deals with all these possibilities. get_wannier_localization_matrix(self, nbands, dirG, kpoint, nextkpoint, G_I, spin) return wannier localization matrix get_wave_function(self, band=0, kpt=0, spin=0) return the wave function. This is the pseudo wave function divided by volume. get_wf = get_wave_function(self, band=0, kpt=0, spin=0) get_xc(self) return the self-consistent exchange-correlation functional used

returns a string get_xc_energies(self, *functional) Get energies for different functionals self-consistent and non-self-consistent.

:Parameters:

functional : strings

  some set of 'PZ','VWN','PW91','PBE','revPBE', 'RPBE'

This function returns the self-consistent energy and/or energies associated with various functionals. The functionals are currently PZ,VWN,PW91,PBE,revPBE, RPBE. The different energies may be useful for calculating improved adsorption energies as in B. Hammer, L.B. Hansen and J.K. Norskov, Phys. Rev. B 59,7413. Examples: get_xcenergies() #returns all the energies get_xcenergies('PBE') # returns the PBE total energy get_xcenergies('PW91','PBE','revPBE') # returns a # list of energies in the order asked for initial_wannier(self, initialwannier, kpointgrid, fixedstates, edf, spin) return initial wannier initnc(self, ncfile=None) create an ncfile with minimal dimensions in it

this makes sure the dimensions needed for other set functions exist when needed. read_only_atoms(self, ncfile) read only the atoms from an existing netcdf file. Used to initialize a calculator from a ncfilename.

:Parameters:

ncfile : string

  name of file to read from.

return ASE.Atoms with no calculator attached or None if no atoms found restart(self) Restart the calculator by deleting nc dimensions that will be rewritten on the next calculation. This is sometimes required when certain dimensions change related to unitcell size changes planewave/densitywave cutoffs and kpt changes. These can cause fortran netcdf errors if the data does not match the pre-defined dimension sizes.

also delete all the output from previous calculation. set(self, **kwargs) set a parameter

parameter is stored in dictionary that is processed later if a calculation is need. set_ados(self, energywindow=(-15, 5), energywidth=0.2, npoints=250, cutoff=1.0) setup calculation of atom-projected density of states

:Parameters:

energywindow : (float, float)

  sets (emin,emax) in eV referenced to the Fermi level

energywidth : float

  the gaussian used in smearing

npoints : integer

  the number of points to sample the DOS at

cutoff : float

  the cutoff radius in angstroms for the integration.

set_ascii_debug(self, level) set the debug level for Dacapo

:Parameters:

level : string

  one of 'Off', 'MediumLevel', 'HighLevel'

set_atoms(self, atoms) attach an atoms to the calculator and update the ncfile

:Parameters:

atoms

  ASE.Atoms instance

set_calculate_stress

set_calculate_stress(self, stress=True) Turn on stress calculation

:Parameters:

stress : boolean

  set_calculate_stress(True) calculates stress
  set_calculate_stress(False) do not calculate stress

set_charge_mixing

set_charge_mixing(self, method='Pulay', mixinghistory=10, mixingcoeff=0.1, precondition='No', updatecharge='Yes') set density mixing method and parameters

:Parameters:

method : string

  'Pulay' for Pulay mixing. only one supported now

mixinghistory : integer

  number of iterations to mix
  Number of charge residual vectors stored for generating
  the Pulay estimate on the self-consistent charge density,
  see Sec. 4.2 in Kresse/Furthmuller:
  Comp. Mat. Sci. 6 (1996) p34ff

mixingcoeff : float

  Mixing coefficient for Pulay charge mixing, corresponding
  to A in G$^1$ in Sec. 4.2 in Kresse/Furthmuller:
  Comp. Mat. Sci. 6 (1996) p34ff
              
precondition : string
  'Yes' or 'No'
  
  * "Yes" : Kerker preconditiong is used,
     i.e. q$_0$ is different from zero, see eq. 82
     in Kresse/Furthmuller: Comp. Mat. Sci. 6 (1996).
     The value of q$_0$ is fix to give a damping of 20
     of the lowest q vector.
  
  * "No" : q$_0$ is zero and mixing is linear (default).

updatecharge : string

  'Yes' or 'No'
  
  * "Yes" : Perform charge mixing according to
     ChargeMixing:Method setting
    
  * "No" : Freeze charge to initial value.
     This setting is useful when evaluating the Harris-Foulkes
     density functional

set_convergence

set_convergence(self, energy=1e-05, density=0.0001, occupation=0.001, maxsteps=None, maxtime=None) set convergence criteria for stopping the dacapo calculator.

:Parameters:

energy : float

  set total energy change (eV) required for stopping

density : float

  set density change required for stopping

occupation : float

  set occupation change required for stopping

maxsteps : integer

  specify maximum number of steps to take

maxtime : integer

  specify maximum number of hours to run.

Autopilot not supported here.

set_debug

set_debug(self, debug) set debug level for python logging

debug should be an integer from 0-100 or one of the logging constants like logging.DEBUG, logging.WARN, etc…

set_decoupling

set_decoupling(self, ngaussians=3, ecutoff=100, gausswidth=0.35) Decoupling activates the three dimensional electrostatic decoupling. Based on paper by Peter E. Bloechl: JCP 103 page7422 (1995).

:Parameters:

ngaussians : int

  The number of gaussian functions per atom
  used for constructing the model charge of the system

ecutoff : int

  The cut off energy (eV) of system charge density in
  g-space used when mapping constructing the model change of
  the system, i.e. only charge density components below
  ECutoff enters when constructing the model change.

gausswidth : float

  The width of the Gaussians defined by
  $widthofgaussian*1.5^(n-1)$  $n$=(1 to numberofgaussians)

set_dipole

set_dipole(self, status=True, mixpar=0.2, initval=0.0, adddipfield=0.0, position=None) turn on and set dipole correction scheme

:Parameters:

status : Boolean

  True turns dipole on. False turns Dipole off

mixpar : float

  Mixing Parameter for the the dipole correction field
  during the electronic minimization process. If instabilities
  occur during electronic minimization, this value may be
  decreased.

initval : float

  initial value to start at

adddipfield : float

  additional dipole field to add
  units : V/ang
  External additive, constant electrostatic field along
  third unit cell vector, corresponding to an external
  dipole layer. The field discontinuity follows the position
  of the dynamical dipole correction, i.e. if
  DipoleCorrection:DipoleLayerPosition is set, the field
  discontinuity is at this value, otherwise it is at the
  vacuum position farthest from any other atoms on both
  sides of the slab.

position : float

  scaled coordinates along third unit cell direction.
  If this attribute is set, DACAPO will position the
  compensation dipole layer plane in at the provided value.
  If this attribute is not set, DACAPO will put the compensation
  dipole layer plane in the vacuum position farthest from any
  other atoms on both sides of the slab. Do not set this to
  0.0

calling set_dipole() sets all default values.

set_dw

set_dw(self, dw) set the density wave cutoff energy.

:Parameters:

dw : integer

  the density wave cutoff

The function checks to make sure it is not less than the planewave cutoff.

Density_WaveCutoff describes the kinetic energy neccesary to represent a wavefunction associated with the total density, i.e. G-vectors for which $ ert G ert^2$ $<$ 4*Density_WaveCutoff will be used to describe the total density (including augmentation charge and partial core density). If Density_WaveCutoff is equal to PlaneWaveCutoff this implies that the total density is as soft as the wavefunctions described by the kinetic energy cutoff PlaneWaveCutoff. If a value of Density_WaveCutoff is specified (must be larger than or equal to PlaneWaveCutoff) the program will run using two grids, one for representing the wavefunction density (softgrid_dim) and one representing the total density (hardgrid_dim). If the density can be reprensented on the same grid as the wavefunction density Density_WaveCutoff can be chosen equal to PlaneWaveCutoff (default).

set_electronic_minimization

set_electronic_minimization(self, method='eigsolve', diagsperband=2) set the eigensolver method

Selector for which subroutine to use for electronic minimization

Recognized options : “resmin”, “eigsolve” and “rmm-diis”.

* “resmin” : Power method (Lennart Bengtson), can only handle

 k-point parallization.

* “eigsolve : Block Davidson algorithm

 (Claus Bendtsen et al).

* “rmm-diis : Residual minimization

 method (RMM), using DIIS (direct inversion in the iterate
 subspace) The implementaion follows closely the algorithm
 outlined in Kresse and Furthmuller, Comp. Mat. Sci, III.G/III.H

:Parameters:

method : string

  should be 'resmin', 'eigsolve' or 'rmm-diis'

diagsperband : int

  The number of diagonalizations per band for
  electronic minimization algorithms (maps onto internal
  variable ndiapb). Applies for both
  ElectronicMinimization:Method = "resmin" and "eigsolve".
  default value = 2

set_external_dipole

set_external_dipole(self, value, position=None) Externally imposed dipole potential. This option overwrites DipoleCorrection if set.

:Parameters:

value : float

  units of volts

position : float

  scaled coordinates along third unit cell direction.
  if None, the compensation dipole layer plane in the
  vacuum position farthest from any other atoms on both
  sides of the slab. Do not set to 0.0.

set_extpot

set_extpot(self, potgrid) add external potential of value

see this link before using this https://listserv.fysik.dtu.dk/pipermail/campos/2003-August/000657.html

:Parameters:

potgrid : np.array with shape (nx,ny,nz)

  the shape must be the same as the fft soft grid
  the value of the potential to add

you have to know both of the fft grid dimensions ahead of time! if you know what you are doing, you can set the fft_grid you want before hand with: calc.set_fftgrid1)

set_extracharge

set_extracharge(self, val) add extra charge to unit cell

:Parameters:

val : float

  extra electrons to add or subtract from the unit cell

Fixed extra charge in the unit cell (i.e. deviation from charge neutrality). This assumes a compensating, positive constant backgound charge (jellium) to forge overall charge neutrality.

set_fftgrid

set_fftgrid(self, soft=None, hard=None) sets the dimensions of the FFT grid to be used

:Parameters:

soft : (n1,n2,n3) integers

  make a n1 x n2 x n3 grid

hard : (n1,n2,n3) integers

  make a n1 x n2 x n3 grid

»> calc.set_fftgrid(soft=[42,44,46]) sets the soft and hard grid dimensions to 42,44,46

»> calc.set_fftgrid(soft=[42,44,46],hard=[80,84,88]) sets the soft grid dimensions to 42,44,46 and the hard grid dimensions to 80,84,88

These are the fast FFt grid numbers listed in fftdimensions.F

data list_of_fft /2, 4, 6, 8, 10, 12, 14, 16, 18, 20, & 22,24, 28, 30,32, 36, 40, 42, 44, 48, & 56,60, 64, 66, 70, 72, 80, 84, 88, 90, & 96,108,110,112,120,126,128,132,140,144,154, & 160,168,176,180,192,198,200, & 216,240,264,270,280,288,324,352,360,378,384,400,432, & 450,480,540,576,640/

otherwise you will get some errors from mis-dimensioned variables.

this is usually automatically set by Dacapo.

set_fixmagmom

set_fixmagmom(self, fixmagmom=None) set a fixed magnetic moment for a spin polarized calculation

:Parameters:

fixmagmom : float

  the magnetic moment of the cell in Bohr magnetons

set_ft

set_ft(self, ft) set the Fermi temperature for occupation smearing

:Parameters:

ft : float

  Fermi temperature in kT (eV)

Electronic temperature, corresponding to gaussian occupation statistics. Device used to stabilize the convergence towards the electronic ground state. Higher values stabilizes the convergence. Values in the range 0.1-1.0 eV are recommended, depending on the complexity of the Fermi surface (low values for d-metals and narrow gap semiconducters, higher for free electron-like metals).

set_kpts

set_kpts(self, kpts) set the kpt grid.

Parameters:

kpts: (n1,n2,n3) or [k1,k2,k3,…] or one of these chadi-cohen sets:

* cc6_1x1 * cc12_2x3 * cc18_sq3xsq3 * cc18_1x1 * cc54_sq3xsq3 * cc54_1x1 * cc162_sq3xsq3 * cc162_1x1

(n1,n2,n3) creates an n1 x n2 x n3 monkhorst-pack grid, [k1,k2,k3,…] creates a kpt-grid based on the kpoints defined in k1,k2,k3,…

There is also a possibility to have Dacapo (fortran) create the Kpoints in chadi-cohen or monkhorst-pack form. To do this you need to set the KpointSetup.gridtype attribute, and KpointSetup.

KpointSetup = [3,0,0] KpointSetup.gridtype = 'ChadiCohen'

KpointSetup(1) Chadi-Cohen k-point set 1 6 k-points 1×1 2 18-kpoints sqrt(3)*sqrt(3) 3 18-kpoints 1×1 4 54-kpoints sqrt(3)*sqrt(3) 5 54-kpoints 1×1 6 162-kpoints 1×1 7 12-kpoints 2×3 8 162-kpoints 3xsqrt 3

or KpointSetup = [4,4,4] KpointSetup.gridtype = 'MonkhorstPack' we do not use this functionality. set_mdos(self, mcenters=None, energywindow=(-15, 5), energywidth=0.2, numberenergypoints=250, cutoffradius=1.0) Setup multicentered projected DOS.

mcenters

 a list of tuples containing (atom#,l,m,weight)
 (0,0,0,1.0) specifies (atom 0, l=0, m=0, weight=1.0) an s orbital
 on atom 0
 
 (1,1,1,1.0) specifies (atom 1, l=1, m=1, weight=1.0) a p orbital
 with m = +1 on atom 0

l=0 s-orbital

 l=1 p-orbital
 l=2 d-orbital

m in range of ( -l … 0 … +l )

The direction cosines for which the spherical harmonics are

 set up are using the next different atom in the list
 (cyclic) as direction pointer, so the z-direction is chosen
 along the direction to this next atom. At the moment the
 rotation matrices is only given in the text file, you can
 use grep'MUL: Rmatrix' out_o2.txt to get this information.
 

adapated from old MultiCenterProjectedDOS.py set_nbands(self, nbands=None) Set the number of bands. a few unoccupied bands are recommended.

:Parameters:

nbands : integer

  the number of bands.
  

if nbands = None the function returns with nothing done. At calculate time, if there are still no bands, they will be set by:

the number of bands is calculated as $nbands=nvalence*0.65 + 4$

set_nc

set_nc(self, nc='out.nc')

set filename for the netcdf and text output for this calculation

:Parameters:

nc : string

  filename for netcdf file
      

if the ncfile attached to the calculator is changing, the old file will be copied to the new file if it doesn not exist so that all the calculator details are preserved. Otherwise, the

if the ncfile does not exist, it will get initialized.

the text file will have the same basename as the ncfile, but with a .txt extension.

set_ncoutput

set_ncoutput(self, wf=None, cd=None, efp=None, esp=None)

set the output of large variables in the netcdf output file

:Parameters:

wf : string

  controls output of wavefunction. values can
  be 'Yes' or 'No'

cd : string

  controls output of charge density. values can
  be 'Yes' or 'No'

efp : string

  controls output of effective potential. values can
  be 'Yes' or 'No'

esp : string

  controls output of electrostatic potential. values can
  be 'Yes' or 'No'

set_occupationstatistics(self, method) set the method used for smearing the occupations.

:Parameters:

method : string

  one of 'FermiDirac' or 'MethfesselPaxton'
  Currently, the Methfessel-Paxton scheme (PRB 40, 3616 (1989).)
  is implemented to 1th order (which is recommemded by most authors).
  'FermiDirac' is the default

set_parent(self, parent) set_pseudopotentials(self, pspdict) Set all the pseudopotentials from a dictionary.

The dictionary should have this form::

{symbol1: path1,

   symbol2: path2}

set_psp

set_psp(self, sym=None, z=None, psp=None)

set the pseudopotential file for a species or an atomic number.

:Parameters:

sym : string

 chemical symbol of the species

z : integer

  the atomic number of the species

psp : string

  filename of the pseudopotential

you can only set sym or z.

examples::

set_psp('N',psp='pspfile')

set_psp(z=6,psp='pspfile')

set_psp_database(self, xc=None) get the xc-dependent psp database

:Parameters:

xc : string

 one of 'PW91', 'PBE', 'revPBE', 'RPBE', 'PZ'

not all the databases are complete, and that means some psp do not exist.

note: this function is not supported fully. only pw91 is imported now. Changing the xc at this point results in loading a nearly empty database, and I have not thought about how to resolve that

set_pw(self, pw)

set the planewave cutoff.

:Parameters:

pw : integer

 the planewave cutoff in eV
 

this function checks to make sure the density wave cutoff is greater than or equal to the planewave cutoff.

set_spinpol(self, spinpol=False)

set Spin polarization.

:Parameters:

spinpol : Boolean

 set_spinpol(True)  spin-polarized.
 set_spinpol(False) no spin polarization, default

Specify whether to perform a spin polarized or unpolarized calculation.

set_status(self, status)

set the status flag in the netcdf file

:Parameters:

status : string

  status flag, e.g. 'new', 'finished'

set_stay_alive(self, value)

set the stay alive setting

set_symmetry(self, val=False)

set how symmetry is used to reduce k-points

:Parameters:

val : Boolean

 set_sym(True) Maximum symmetry is used
 set_sym(False) No symmetry is used

This variable controls the if and how DACAPO should attempt using symmetry in the calculation. Imposing symmetry generally speeds up the calculation and reduces numerical noise to some extent. Symmetry should always be applied to the maximum extent, when ions are not moved. When relaxing ions, however, the symmetry of the equilibrium state may be lower than the initial state. Such an equilibrium state with lower symmetry is missed, if symmetry is imposed. Molecular dynamics-like algorithms for ionic propagation will generally not break the symmetry of the initial state, but some algorithms, like the BFGS may break the symmetry of the initial state. Recognized options:

“Off”: No symmetry will be imposed, apart from time inversion symmetry in recipical space. This is utilized to reduce the k-point sampling set for Brillouin zone integration and has no influence on the ionic forces/motion.

“Maximum”: DACAPO will look for symmetry in the supplied atomic structure and extract the highest possible symmetry group. During the calculation, DACAPO will impose the found spatial symmetry on ionic forces and electronic structure, i.e. the symmetry will be conserved during the calculation.

set_xc(self, xc)

Set the self-consistent exchange-correlation functional

:Parameters:

xc : string

 Must be one of 'PZ', 'VWN', 'PW91', 'PBE', 'revPBE', 'RPBE'

Selects which density functional to use for exchange-correlation when performing electronic minimization (the electronic energy is minimized with respect to this selected functional) Notice that the electronic energy is also evaluated non-selfconsistently by DACAPO for other exchange-correlation functionals Recognized options :

* “PZ” (Perdew Zunger LDA-parametrization) * “VWN” (Vosko Wilk Nusair LDA-parametrization) * “PW91” (Perdew Wang 91 GGA-parametrization) * “PBE” (Perdew Burke Ernzerhof GGA-parametrization) * “revPBE” (revised PBE/1 GGA-parametrization) * “RPBE” (revised PBE/2 GGA-parametrization)

option “PZ” is not allowed for spin polarized calculation; use “VWN” instead.

strip(self)

remove all large memory nc variables not needed for anything I use very often. update_input_parameters(self) read in all the input parameters from the netcdfile

write(self, new=False)

write out everything to the ncfile : get_nc()

new determines whether to delete any existing ncfile, and rewrite it.

write_input(self)

write out input parameters as needed

you must define a self._set_keyword function that does all the actual writing.

write_nc

write_nc(self, nc=None, atoms=None)

write out atoms to a netcdffile.

This does not write out the calculation parameters!

Parameters

NameTypeDescription
ncstring ncfilename to write to. this file will get clobbered if it already exists.
atomsASE.Atoms atoms to write. if None use the attached atoms if no atoms are attached only the calculator is written out.

the ncfile is always opened in 'a' mode.

note

it is good practice to use the atoms argument to make sure that the geometry you mean gets written! Otherwise, the atoms in the calculator is used, which may be different than the external copy of the atoms.

Static methods

read_atoms

read_atoms(filename)

read atoms and calculator from an existing netcdf file.

Parameters
NameTypeDescription
filenamestringname of file to read from
example
  >>> atoms = Jacapo.read_atoms(ncfile)
  >>> calc = atoms.get_calculator()

This method is here for legacy purposes. I used to use it alot.

1)
n1,n2,n3
info/dacapo.txt · 最終更新: 2022/08/23 13:34 by 127.0.0.1

Donate Powered by PHP Valid HTML5 Valid CSS Driven by DokuWiki