DFTB_methods¶
-
libra_py.DFTB_methods.
cube_generator_dftbplus
(project_name, time_step, min_band, max_band, waveplot_exe, isUKS)[source]¶ This function generates the cube files by first forming the ‘fchk’ file from ‘chk’ file. Then it will generate the cube files from min_band to max_band the same as CP2K naming. This will helps us to use the read_cube and integrate_cube functions easier.
- Parameters
project_name (string) – The project_name.
time_step (integer) – The time step.
min_band (integer) – The minimum state number.
max_band (integer) – The maximum state number.
waveplot_exe (str) – Location of the executable for the waveplot program of the dftb+ software package
isUKS (integer) – The unrestricted spin calculation flag. 1 is for spin unrestricted calculations. Other numbers are for spin restricted calculations
- Returns
None
-
libra_py.DFTB_methods.
dftb_distribute
(istep, fstep, nsteps_this_job, trajectory_xyz_file, dftb_input, waveplot_input, curr_job_number)[source]¶ Distributes dftb jobs for trivial parallelization
Make sure that your dftb input file has absolute paths to the following input parameters:
- SlaterKosterFiles = Type2FileNames {
Prefix = “/panasas/scratch/grp-alexeyak/brendan/dftbp_development/” Separator = “-” Suffix = “.skf”
}
- Parameters
istep (integer) – The initial time step in the trajectory xyz file.
fstep (integer) – The final time step in the trajectory xyz file.
nsteps_this_job (integer) – The number of steps for this job.
trajectory_xyz_file (string) – The full path to trajectory xyz file.
dftb_input (string) – the dftb_input_template.hsd file
waveplot_input (string) – the input file for the waveplot subpackage of dftbplus for generating cubefiles
curr_job_number (integer) – The current job number.
- Returns
None
-
libra_py.DFTB_methods.
get_dftb_ks_energies
(params)[source]¶ This function reads the band.out file generated from a dftb+ computations. The band.out file stores the ks energies and their occupation.
Args:
- Returns
The vector consisting of the KS energies from min_band to max_band.
total_energy (float): The total energy obtained from the log file.
- Return type
E (1D numpy array)
-
libra_py.DFTB_methods.
get_dftb_matrices
(filename, act_sp1=None, act_sp2=None)[source]¶ Get the matrices printed out by the DFTB+
- Parameters
filename (string) – the name of the file to read, usually these are any of the files: “hamsqr1.dat”, “hamsqr2.dat”, etc. or “oversqrt.dat”. Produced with the input containing option “WriteHS = Yes”
act_sp1 (list of ints or None) – the row active space to extract from the original files Indices here start from 0. If set to None - the number of AOs will be determined automatically from the file. [default: None]
act_sp2 (list of ints or None) – the cols active space to extract from the original files Indices here start from 0. If set to None - the number of AOs will be determined automatically from the file. [default: None]
- Returns
- X: where N = len(act_sp1) and M = len(act_sp2)
These are the corresponding property matrices (converted to the complex type) for each k-point, such that X[ikpt] is a CMATRIX(N, M) containing the overlaps/Hamiltonians for the k-point
`ikpt`
.
- Return type
list of CMATRIX(N, M)
Warning
So far tested only for a single k-point!
-
libra_py.DFTB_methods.
get_energy_forces
(filename, nat)[source]¶ Get forces from the input file
- Parameters
filename (string) – the name of the file to read, usually this is a the file “detailed.out” produced with the input containing option “CalculateForces = Yes “
nat (int) – the number of atoms in the system
- Returns
( E_ex, F, Flst ), where:
E_ex ( double ): excitation energy for the given state [ units: a.u. ]
F ( MATRIX(ndof, 1) ): the forces acting on all atoms [ units: a.u. of force ]
Flst ( list of [x, y, z] ): the forces acting on all atoms [ units: a.u. of force ]
- Return type
Warning
it is likely not gonna work for files other than “detailed.out”
-
libra_py.DFTB_methods.
make_dftb_input
(dftb_input_template_filename, istate, timestep=- 1)[source]¶ This function creates an input file for DFTB+ package from a template file, it changes several placeholder lines to ensure the calculations are done for required electronic states
- Parameters
dftb_input_template_filename (strings) – the name of the generic input file (template) for DFTB+ calculations
istate (int) – the index of the state, for which the calculations will be carried out
timestep (int) – the index of a given timestep along a MD trajectory. If not provided (-1), then the defaul name of “tmp.gen” will be used as the coordinate file name. If a time index is given, the name will be “coord-<timestep>.gen”.
- Returns
just creates the files
- Return type
-
libra_py.DFTB_methods.
read_dftb_output
(natoms, istate)[source]¶ This file reads in the total energy (essentially the reference, ground state energy), the excitation energies, and the corresponding forces on all atoms and return them in a digital format for further processing.
-
libra_py.DFTB_methods.
read_dftbplus_TRA_file
(params)[source]¶ This function reads TRA.dat files generated from TD-DFTB calculations using DFTB+ and returns the TD-DFTB excitation energies and CI-like coefficients
- Parameters
params (dictionary) –
parameters dictionary
- logfile_name ( string ): this is actualyl the TRA.dat file, but we keep it as logfile_name for ease, as many other similar
type functions (for cp2k, gaussian) rely on this format internally.
number_of_states ( int ): how many ci states to consider tolerance ( float ): the tolerance for accepting SDs are members of the CI wavefunctions isUKS ( Boolean ): flag for whether or not a spin-polarized Kohn-Sham basis is being used. TRUE means
that a spin-polarized Kohn-Sham basis is being used.
- Returns
The excitation energies in the log file. ci_basis ( list ): The CI-basis configuration. ci_coefficients ( list ): The coefficients of the CI-states. spin_components (list): The spin component of the excited states.
- Return type
excitation_energies ( list )
-
libra_py.DFTB_methods.
run_dftb_adi
(q, params_, full_id)[source]¶ This function executes the DFTB+ quantum chemistry calculations and returns the key properties needed for dynamical calculations.
- Parameters
q (MATRIX(ndof, ntraj)) – coordinates of the particle [ units: Bohr ]
params (dictionary) –
model parameters
- params[“labels”] ( list of strings ): the labels of atomic symbolc - for all atoms,
and in a order that is consistent with the coordinates (in triples) stored in q. The number of this labels is natoms, such that ndof = 3 * natoms. [ Required ]
- params[“nstates”] ( int ): the total number of electronic states
in this model [ default: 1 - just the ground state]
- params[“dftb_input_template_filename”] ( string ): the name of the input file template
[ default: “dftb_input_template.hsd” ]
- params[“dftbp_exe”] ( string ): the full path to the DFTB+ executable
[ defaut: “dftb+” ]
- params[“xyz2gen_exe”] ( string ): the full path to the xyz2gen executable
(part of the DFTB+ package) [ defaut: “xyz2gen” ] Note: sometimes, especially if you are using conda-installed Python, you may need to edit the “xyz2gen” file in your DFTB+ installation to change the topmost line to point to the correct python executable (e.g. #! /home/alexey/miniconda2/envs/py37/bin/python)
- Returns
obj, with the members:
obj.ham_adi ( CMATRIX(nstates,nstates) ): adiabatic Hamiltonian
- obj.d1ham_adi ( list of ndof CMATRIX(nstates, nstates) objects ):
derivatives of the adiabatic Hamiltonian w.r.t. the nuclear coordinate
- Return type
PyObject
-
libra_py.DFTB_methods.
xyz_traj2gen_ovlp
(infile, outfile, md_iter1, md_iter2, sys_type)[source]¶ This file converts the xyz trajectory file in the format produced by DFTB+ to a gen file containing the superimposed md_iter1-th and md_iter2-th steps geometries
- Parameters
infile (string) – the name of the input xyz trajectory file
outfile (string) – the name of a output gen single point file
md_iter1 (int) – index of the first timeframe to extract
md_iter2 (int) – index of the second timeframe to extract
sys_type (string) – “C” or “P” - the system type - cluster or periodic
- Returns
but creates a file
- Return type
none
-
libra_py.DFTB_methods.
xyz_traj2gen_sp
(infile, outfile, md_iter, sys_type)[source]¶ This file converts the xyz trajectory file in the format produced by DFTB+ to a gen file containing the md_iter-th step geometry
- Parameters
infile (string) – the name of the input xyz trajectory file
outfile (string) – the name of a output gen single point file
md_iter (int) – index of the timeframe to extract
sys_type (string) – “C” or “P” - the system type - cluster or periodic
- Returns
but creates a file
- Return type
none