step2_many_body

libra_py.workflows.nbra.step2_many_body.compute_cube_ks_overlaps(cubefiles_prev, params)[source]

This function computes overlaps between cube files of two time steps. In order to not read the cube files twice it returns the cube files of the current step, and gets the cube files of the previous step.

Parameters
  • cubefiles_prev (list) – The list containing th cube files of the previous step.

  • params (dict) –

    curr_step (int): The current time step.

    isUKS (int): This parameter is set for spin restricted and unrestricted calculations. When it is

    set to 1 it means that unrestricted calculations were set in the input file otherwise it is restricted.

    nprocs (int): The number of processors used to read the cube files and perform the integration.

Returns

The list of the read current step cube files.

S_ks_prev (2D numpy array): The overlap matrix of the wavefunctions for the previous time step.

S_ks_curr (2D numpy array): The overlap matrix of the wavefunctions for the current time step.

St_ks (2D numpy array): The overlap matrix between wavefunctions of the two time step.

Return type

cubefiles_curr (list)

libra_py.workflows.nbra.step2_many_body.curr_and_final_step_job(istep, fstep, njobs, njob)[source]

This function is used to determine the initial and final step of a job when distributing the molecular dynamics trajectory over a number of jobs.

Parameters
  • istep (int) – The initial time step for molecular dynamics trajectory which the user wants to start the calculations with. This paramter is set in the input file for submitting the jobs.

  • fstep (int) – The final time step for molecular dynamics trajectory which the user wants to start the calculations with. This paramter is set in the input file for submitting the jobs.

  • njobs (int) – The number of jobs specified by user.

  • njob (int) – The job number.

Returns

The job initial time step.

job_final_step (int): The job final time step.

Return type

job_init_step (int)

libra_py.workflows.nbra.step2_many_body.form_Hvib_real(params)[source]

This function forms the real part of the vibronic Hamiltonian by inserting the energies on the diagonal of a zero matrix. The Hvib is in two-spinor format containing the alpha and beta states energies.

Parameters

params (distionary) –

logfile_directory (str): The directory of the output files by CP2K.

time (int): The time step of the molecular dynamics. For single point calculations

it should be set to 0.

min_band (int): The minimum state number.

max_band (int): The maximum state number.

isUKS (int): This flag is used whenever the UKS calculations were called in CP2K input. if it is set to

1 then UKS calculations were called and if not it will consider only alpha energies.

Returns

The diagonal matrix containing the energies of alpha and beta spins if

UKS is set to True and the energies of alpha spin in a block form matrix.

Return type

Hvib_ks_re (2D numpy array)

libra_py.workflows.nbra.step2_many_body.get_excitation_analysis_output(params)[source]

This function reads the information of the excited states from the log file of the single point calculations.

Parameters

params (dict) –

logfile_directory (str): The log files directory.

es_software (str): The name of the software used to calculate the energy calculations.

curr_step (int): The current time step of the calculations.

isUKS (int): This parameter is set for spin restricted and unrestricted calculations. When it is

set to 1 it means that unrestricted calculations were set in the input file otherwise it is restricted.

Returns

The excitation energies of the curr_step.

ci_basis_raw (numpy array): The excited states which contains the occupied and virtual orbitals.

ci_coefficients_raw_unnorm (numpy array): The excited states CI coefficients.

spin_components (numpy array): Contains the excited states spin components (‘alp’ for alpha spin and ‘bet’ for beta spin)

Return type

excitation_energies (numpy array)

libra_py.workflows.nbra.step2_many_body.integrate_cube_set(cubefiles_set_1, cubefiles_set_2, dv)[source]

This function comutes the overlap matrix between two set of cube files.

Parameters
  • cubefiles_set_1 (list) – The list of cube files of the curr_step.

  • cubefiles_set_2 (list) – The list of cube files of the previous step.

  • dv (float) – The integration element obtained from the cube files.

Returns

The overlap between the two set of cube files.

Return type

overlap_matrix (numpy 2D array)

libra_py.workflows.nbra.step2_many_body.normalize_ci_coefficients(ci_coefficients_raw_unnorm)[source]

This funciton normalizes the list of configuration interaction (CI) coefficients.

Parameters

ci_coefficients_raw_unnorm (list) – The list containing the lists of unnormalized CI coefficients.

Returns

The list containing the lists of normalized CI coefficients.

Return type

ci_coefficients_raw_norm (list)

libra_py.workflows.nbra.step2_many_body.reindex_cp2k_sd_states(ks_orbital_homo_index, ks_orbital_indicies, sd_basis_states, sd_format=2)[source]

ks_orbital_homo_index: Index of the homo ks orbital, from 1 ks_orbital_indicies: Range of the considered ks orbtials. Ex) [8,9,10,11], where 9 is homo orbtial index (from 1) sd_basis_states( list of lists of lists ): A list of Slater determinants, where each slater determinant is a excitation in the Kohn-Sham

basis. This function assumes that all Kohn-Sham excitations are for alpha electrons. To differentiate between alpha and beta excitations, elements of sd_basis_states contain spin information.

Ex) sd_basis_states[0] = [ [9,10], “alp” ]

sd_basis_states[1] = [ [9,10], “bet” ]

Returns

Reindexed Slater determinant basis states in terms of the ks_orbital_homo_index. The SD states returned can be either in one of two ways

1. All alpha orbitals before beta orbtials. This format is still a fixed slot format. It is used when used both alpha and beta spin channels Ex.) ks_orbital_homo_index = 68, ks_orbital_indicies = [67, 68, 69, 70, 71, 72]

The ground state will be: [1, 2, -7, -8] [ [68,69], “alp” ] –> [1, 3, -7, -8] [ [68,69], “bet” ] –> [1, 2, -7, -9]

2. Fixed slot format. this format used only the alpha spin channel Ex.) ks_orbital_homo_index = 68, ks_orbital_indicies = [67, 68, 69, 70, 71, 72]

The ground state will be: [1, -1, 2, -2] [ [68,69], “alp” ] –> [1, -1, 3, -2] [ [68,69], “bet” ] –> [1, -1, 2, -3]

libra_py.workflows.nbra.step2_many_body.run_step2_many_body(params)[source]

This function is the main function which runs the following calculations:

  • Executes calls to perform CP2K, DFTB+, or Gaussian electronic structure calculations

  • Computing the overlap, time-overlap, and energy eigenvalue matrices for the KS basis

  • Molecular orbital visualization using VMD

Parameters

params (dictionary) –

min_band (integer): The minimum state number.

max_band (integer): The maximum state number.

ks_orbital_homo_index (integer): The index of the Kohn-Sham HOMO orbital (State numbers start from 1).

nsteps_this_job (integer): The number of steps in this job.

es_software (string): The electronic structure calculation software.

  • Note: The current code works with CP2K but it can be generalized to other electronic structure

    calculation software packages which are capable of producing cube files and performing TD-DFT calculations such as Gaussian. So this parameter is for future implementations.

es_software_exe (string): The CP2K executable file or path.

es_software_input_template (str): The full path to CP2K input template.

project_name (string): The project name for CP2K input.

isUKS (integer): This is the flag for unrestricted spin calculations. When it is set to

1 it will perform the nrestricted calculations otherwise it will consider the restricted case.

trajectory_xyz_filename (string): The full path to trajectory xyz file.

njob (integer): The current job number.

nprocs (integer): The number of processors to be used.

res_dir (string): The full path to res directory where the output data are stored.

logfile_directory (string): The path to where the log files are stored.

Returns

None