Public API

QcInput

class pyqchem.qc_input.QchemInput(molecule, jobtype='sp', method=None, exchange=None, correlation=None, unrestricted=None, basis='6-31G', basis2=None, thresh=14, scf_convergence=8, max_scf_cycles=50, scf_algorithm='diis', purecart=None, ras_roots=None, ras_do_hole=True, ras_do_part=True, ras_act=None, ras_act_orb=None, ras_elec=None, ras_elec_alpha=None, ras_elec_beta=None, ras_occ=None, ras_spin_mult=1, ras_sts_tm=False, ras_fod=False, ras_natorb=False, ras_natorb_state=None, ras_print=1, ras_diabatization_scheme=None, ras_diabatization_states=None, ras_guess=None, use_reduced_ras_guess=False, ras_omega=400, ras_srdft=None, ras_srdft_damp=0.5, ras_srdft_exc=None, ras_srdft_cor=None, ras_srdft_spinpol=0, calc_soc=False, state_analysis=False, ee_singlets=False, ee_triplets=False, cc_trans_prop=False, cc_symmetry=True, cc_e_conv=None, cc_t_conv=None, eom_davidson_conv=5, cis_convergence=6, cis_n_roots=None, cis_singlets=False, cis_triplets=False, cis_ampl_anal=False, loc_cis_ov_separate=False, er_cis_numstate=0, boys_cis_numstate=0, cis_diabath_decompose=False, max_cis_cycles=30, localized_diabatization=None, sts_multi_nroots=None, cc_state_to_opt=None, cis_state_deriv=None, RPA=False, set_iter=30, gui=0, geom_opt_dmax=300, geom_opt_update=-1, geom_opt_linear_angle=165, geom_opt_coords=-1, geom_opt_tol_gradient=300, geom_opt_tol_displacement=1200, geom_opt_tol_energy=100, geom_opt_max_cycles=50, geom_opt_constrains=None, solvent_method=None, solvent_params=None, pcm_params=None, rpath_coords=0, rpath_direction=1, rpath_max_cycles=20, rpath_max_stepsize=150, rpath_tol_displacement=5000, symmetry=True, sym_ignore=False, nto_pairs=None, n_frozen_core=None, n_frozen_virt=None, n_frozen_virtual=0, mom_start=False, reorder_orbitals=None, namd_nsurfaces=None, scf_print=None, scf_guess=None, scf_energies=None, scf_density=None, scf_guess_mix=False, hessian=None, sym_tol=5, mem_total=2000, mem_static=64, skip_scfman=False, extra_rem_keywords=None, extra_sections=None)[source]

Handles the Q-Chem input info

get_copy()[source]

Get a copy of the input

Returns:
get_txt()[source]

get qchem input in plain text

Return string:qchem input in plain text
update_input(dictionary)[source]

Update the input from data in a dictionary Note: already existing parameters will be overwritten

Parameters:dictionary – parameters to add
pyqchem.qc_input.normalize_values(value)[source]

Set all string values (including keys and values of dictionaries) to lower case

Parameters:value – the values
Returns:normalized values

Structure

class pyqchem.structure.Structure(coordinates=None, symbols=None, atomic_numbers=None, connectivity=None, charge=0, multiplicity=1, name=None)[source]

Structure object containing all the geometric data of the molecule

alpha_electrons

returns the alpha electrons

Returns:number of alpha electrons
beta_electrons

returns the number of beta electrons

Returns:number of beta electrons
charge

returns the charge :return: the charge

get_atomic_masses()[source]

get the atomic masses of the atoms of the molecule

Returns:list of atomic masses
get_atomic_numbers()[source]

get the atomic numbers of the atoms of the molecule

Returns:list with the atomic numbers
get_connectivity(thresh=1.2)[source]

get the connectivity as a list of pairs of indices of atoms from atomic radii

Parameters:thresh – radii threshold used to determine the connectivity
Returns:
get_coordinates(fragment=None)[source]

gets the cartesian coordinates

Parameters:fragment – list of atoms that are part of the fragment
Returns:coordinates list
get_number_of_atoms()[source]

get the number of atoms

Returns:number of atoms
get_point_symmetry()[source]

Returns the point group of the molecule using pymatgen

Returns:point symmetry label
get_symbols()[source]

get the atomic element symbols of the atoms of the molecule

Returns:list of symbols
get_valence_electrons()[source]

get number of valence electrons

Returns:number of valence electrons
get_xyz(title='')[source]

generates a XYZ formatted file

Parameters:title – title of the molecule
Returns:string with the formatted XYZ file
multiplicity

returns the multiplicity

Returns:the multiplicity
name

returns the name :return: structure name

number_of_electrons

returns the total number of electrons

Returns:number of total electrons
set_coordinates(coordinates)[source]

sets the cartessian coordinates

Parameters:coordinates – cartesian coordinates matrix

QcCore

pyqchem.qchem_core.generate_additional_files(input_qchem, work_dir)[source]

Generate additional files on scratch (work dir) for special calculations

Parameters:
  • input_qchem – QChem input object
  • work_dir – scratch directory
pyqchem.qchem_core.get_output_from_qchem(input_qchem, processors=1, use_mpi=False, scratch=None, read_fchk=False, return_electronic_structure=False, parser=None, parser_parameters=None, force_recalculation=False, fchk_only=False, store_full_output=False, delete_scratch=True, remote=None, scratch_read_level=0)[source]

Runs qchem and returns the output in the following format:

  1. If return_electronic_structure is requested:
    [output, parsed_fchk]
  2. If return_electronic_structure is not requested:
    [output]
Note: if parser is set then output contains a dictionary with the parsed info
else output contains the q-chem output in plain text
Parameters:
  • input_qchem – QcInput object containing the Q-Chem input
  • processors – number of threads/processors to use in the calculation
  • use_mpi – If False use OpenMP (threads) else use MPI (processors)
  • scratch – Full Q-Chem scratch directory path. If None read from $QCSCRATCH
  • return_electronic_structure – if True, returns the parsed FCHK file containing the electronic structure
  • read_fchk – same as return_electronic_structure (to be deprecated)
  • parser – function to use to parse the Q-Chem output
  • parser_parameters – additional parameters that parser function may have
  • force_recalculation – Force to recalculate even identical calculation has already performed
  • fchk_only – If true, returns electronic structure data from cache ignoring output (to be deprecated)
  • remote – dictionary containing the data for remote calculation (beta)
  • store_full_output – store full output in plain text in pkl file
  • delete_scratch – delete all scratch files when calculation is finished
Returns:

output [, electronic_structure]

pyqchem.qchem_core.local_run(input_file_name, work_dir, fchk_file, use_mpi=False, processors=1)[source]

Run Q-Chem locally

Parameters:
  • input_file_name – Q-Chem input file in plain text format
  • work_dir – Scratch directory where calculation run
  • fchk_file – filename of fchk
  • use_mpi – use mpi instead of openmp
Returns:

output, err: Q-Chem standard output and standard error

pyqchem.qchem_core.local_run_stream(input_file_name, work_dir, fchk_file, use_mpi=False, processors=1, print_stream=True)[source]

Run Q-Chem locally

Parameters:
  • input_file_name – Q-Chem input file in plain text format
  • work_dir – Scratch directory where calculation run
  • fchk_file – filename of fchk
  • use_mpi – use mpi instead of openmp
  • print_stream – set True to print output stream during execution
Returns:

output, err: Q-Chem standard output and standard error

pyqchem.qchem_core.parse_output(get_output_function)[source]

to be deprecated

Parameters:get_output_function
Returns:parsed output
pyqchem.qchem_core.remote_run(input_file_name, work_dir, fchk_file, remote_params, use_mpi=False, processors=1)[source]

Run Q-Chem remotely

Parameters:
  • input_file – Q-Chem input file in plain text format
  • work_dir – Scratch directory where calculation run
  • fchk_file – filename of fchk
  • remote_params – connection parameters for paramiko
  • use_mpi – use mpi instead of openmp
Returns:

output, err: Q-Chem standard output and standard error

pyqchem.qchem_core.retrieve_additional_files(input_qchem, data_fchk, work_dir, scratch_read_level=0)[source]

retrieve data from files in scratch data (on development, currently for test only)

Parameters:
  • input_qchem – QChem input object
  • data_fchk – FCHK parsed dictionary
  • work_dir – scratch directory
  • scratch_read_level – defines what data to retrieve
Returns:

dictionary with additional data

Utils

pyqchem.utils.get_inertia(structure)[source]

returns the inertia moments and main axis of inertia (in rows)

Parameters:structure – Structure object containg the molecule
Returns:eigenvalues, eigenvectors
pyqchem.utils.get_plane(coords, direction=None)[source]

Returns the center and normal vector of tha plane formed by a list of atomic coordinates

Parameters:coords – List of atomic coordinates
Returns:center, normal_vector
pyqchem.utils.get_sdm(matrix_1, matrix_2)[source]

get differences square matrix between to matrices :param matrix_1: the matrix 1 :param matrix_2: the matrix 2 :return: difference square matrix

pyqchem.utils.is_rasci_transition(configuration, reference, n_electron=1, max_jump=10)[source]

Determine if a configuration corresponds to a transition of n_electron

Parameters:
  • configuration – dictionary containing the configuration to be analyzed
  • reference – reference configuration (in general lowest energy Slater determinant)
  • n_electron
  • max_jump – Restrict to transitions with jumps less or equal to max_jump orbitals
Returns:

True if conditions are met, otherwise False

pyqchem.utils.is_transition(configuration, reference, n_electron=1, max_jump=10)[source]

Determine if a configuration corresponds to a transition of n_electron

Parameters:
  • configuration – dictionary containing the configuration to be analyzed
  • reference – reference configuration (in general lowest energy Slater determinant)
  • n_electron – number of electrons in the transition
  • max_jump – Restrict to transitions with jumps less or equal to max_jump orbitals
Returns:

True if conditions are met, otherwise False

pyqchem.utils.reorder_coefficients(occupations, coefficients)[source]

Reorder the coefficients according to occupations. Occupated orbitals will be grouped at the beginning non occupied will be attached at the end

Parameters:
  • occupations – list on integers (0 or 1) or list of Boolean
  • coefficients – dictionary containing the molecular orbitals coefficients {‘alpha’: coeff, ‘beta:’ coeff}. coeff should be a list of lists (Norb x NBas)
Returns:

pyqchem.tools.get_geometry_from_pubchem(entry, type='name')[source]

Get structure form PubChem database

Parameters:
  • entry – entry data
  • type – data type: ‘name’, ‘cid’
Returns:

Structure

pyqchem.tools.plot_rasci_state_configurations(states)[source]

Prints :param states: parsed data (excited states) dictionary entry from RASCI calculation :return: None

pyqchem.tools.print_excited_states(parsed_data, include_conf_rasci=False, include_mulliken_rasci=False)[source]

Prints excited states in nice format. It works for CIS/TDDFT/RASCI methods

Parameters:
  • parsed_data – parsed data (excited states) dictionary entry from CIS/RASCI/TDDFT calculation
  • include_conf_rasci – print also configuration data (only RASCI method)
  • include_mulliken_rasci – print also mulliken analysis (only RASCI method)
Returns:

None

pyqchem.tools.rotate_coordinates(coordinates, angle, axis, atoms_list=None, center=(0, 0, 0))[source]

Rotate the coordinates (or range of coordinates) with respect a given axis

Parameters:
  • coordinates – coordinates to rotate
  • angle – rotation angle in radians
  • axis – rotation axis
  • atoms_list – list of atoms to rotate (if None then rotate all)
Returns:

rotated coordinates

pyqchem.tools.submit_notice(message, service='pushbullet', pb_token=None, sp_url=None, gc_key=None, gc_token=None, gc_thread=None, slack_token=None, slack_channel=None)[source]

Submit a notification using webhooks

Parameters:
  • message – The message to send
  • service – pushbullet, samepage, google_chat
  • pb_token – pushbullet token
  • sp_url – samepage url
  • gc_key – google chat key
  • gc_token – google chat token
  • gc_thread – google chat thread
  • slack_token – slack bot token (xoxb-xxx.xxx.xxx),
  • slack_channel – slack channel
Returns:

server response

pyqchem.tools.geometry.get_angle(coordinates, atoms)[source]

Compute the angle between 3 atoms

Parameters:
  • coordinates – list of coordinates of the molecule
  • atoms – list of 3 atom indices to use to calculate the angle (from 1 to N)
Returns:

the angle

pyqchem.tools.geometry.get_dihedral(coordinates, atoms)[source]

Compute the dihedral angle between 4 atoms

Parameters:
  • coordinates – list of coordinates of the molecule
  • atoms – list of 4 atom indices to use to calculate the dihedral angle (from 1 to N)
Returns:

the dihedral angle

pyqchem.tools.geometry.get_distance(coordinates, atoms)[source]

Compute the distance between 2 atoms

Parameters:
  • coordinates – list of coordinates of the molecule
  • atoms – list of 2 atom indices to use to calculate the distance (from 1 to N)
Returns:

the distance

pyqchem.tools.geometry.unit_vector(vector)[source]

Compute unit vector from general vector

Parameters:vector – the vector
Returns:the unit vector