Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages  

QMData Class Reference

QM data management class. Concerns QM related data that are not timestep dependent. More...

#include <QMData.h>

List of all members.

Public Methods

 QMData (int natoms, int nbasis, int shells, int nwave)
 constructor:. More...

 ~QMData (void)
 destructor. More...

void init_electrons (Molecule *mol, int totcharge)
 Initialize total charge, multiplicity and everything that can be derived from that. More...

void delete_basis_set ()
 Free memory of the basis set. More...

int init_basis (Molecule *mol, int num_basis_atoms, const char *string, const float *basis, const int *atomic_numbers, const int *nshells, const int *nprim, const int *shelltypes)
 Populate basis set data and organize them into hierarcical data structures. More...

int create_unique_basis (Molecule *mol, int num_basis_atoms)
 Collapse basis set so that we have one atomic basis set per atom type. More...

void normalize_basis ()
 Multiply contraction coefficients with the shell dependent part of the normalization factor. More...

void init_angular_norm_factors ()
 Initialize table of normalization factors containing different factors for each shell and its cartesian functions. More...

const basis_atom_tget_basis (int atom=0) const
 Get basis set for an atom. More...

const shell_tget_basis (int atom, int shell) const
 Get basis set for a shell. More...

const int * get_atom_types () const
const int * get_num_shells_per_atom () const
const int * get_num_prim_per_shell () const
int get_num_atoms ()
int get_num_types ()
int get_num_shells ()
int get_num_wavecoeff_per_atom (int atom) const
int get_wave_offset (int atom, int shell=0) const
 Get the offset in the wavefunction array for given shell. More...

const char * get_shell_type_str (const shell_t *shell)
 Get shell type letter (S, P, D, F, ...) followed by '\0'. More...

char * get_angular_momentum_str (int atom, int shell, int mom) const
void set_angular_momentum_str (int atom, int shell, int mom, const char *tag)
void set_angular_momentum (int atom, int shell, int mom, int *pow)
int get_angular_momentum (int atom, int shell, int mom, int comp)
int set_angular_momenta (const int *angmom)
int get_num_imag ()
int get_num_intcoords ()
void set_carthessian (int numcartcoords, double *array)
void set_inthessian (int numintcoords, double *array)
void set_normalmodes (int numcart, float *array)
void set_wavenumbers (int numcart, float *array)
void set_intensities (int numcart, float *array)
void set_imagmodes (int numimag, int *array)
const double * get_carthessian () const
const double * get_inthessian () const
const float * get_normalmodes () const
const float * get_wavenumbers () const
const float * get_intensities () const
const int * get_imagmodes () const
const char * get_scftype_string (void) const
 Translate the scftype constant into a string. More...

const char * get_runtype_string (void) const
 Translate the runtype constant into a string. More...

const char * get_status_string ()
 Get status of SCF and optimization convergence. More...

int assign_wavef_id (int type, int spin, int excit, char *info, wavef_signa_t *(&signa_ts), int &numsig)
 Determines and returns a unique ID for each wavefuntion based on it's signature (type, spin, excitation, info). More...

int find_wavef_id_from_gui_specs (int type, int spin, int exci)
 Matches the wavefunction selected from the GUI buttons (type, spin, excitation) with a wavefunction IDtag. Returns -1 if no such wavefunction exists. More...

int compare_wavef_guitype_to_type (int guitype, int type)
 Return 1 if the wavefunction type from the GUI matches the given type. More...

int has_wavef_guitype (int guitype)
 Return 1 if we have a wavefunction with the given GUI_WAVEF_TYPE_*. More...

int has_wavef_spin (int spin)
 Return 1 if we have any wavefunction with the given spin. More...

int has_wavef_exci (int exci)
 Return 1 if we have any wavefunction with the given excitation. More...

int has_wavef_signa (int type, int spin, int exci)
 Return 1 if we have any wavefunction with the given signature (type, spin, and excitation). More...

int get_highest_excitation (int guitype)
 Get the highest excitation for any wavefunction with the given type. More...

void update_avail_orbs (int iwavesig, int norbitals, const int *orbids, const float *orbocc)
 Update the list of available orbitals. More...

int get_max_avail_orbitals (int iwavesig)
 Return the maximum number of available orbitals for the given wavefunction over all frames Can be used to determine the number of orbitals to be displayed in the GUI. Returns -1 if requested wavefunction doesn't exist. More...

int get_avail_orbitals (int iwavesig, int *(&orbids))
 Get IDs of all orbitals available for the given wavefunction. iwavesig is the index of the wavefunction signature. Returns 1 upon success, 0 otherwise. More...

int get_avail_occupancies (int iwavesig, float *(&orbocc))
 Get the occupancies of all available orbitals for the given wavefunction. More...

int get_orbital_label_from_gui_index (int iwavesig, int iorb)
 For the given wavefunction signature return the iorb-th orbital label. Used to translate from the GUI list of available orbitals the unique orbital label. More...

int has_orbital (int iwavesig, int orbid)
 Return 1 if the given wavefunction has an orbital with ID orbid in any frame. More...

int expand_atompos (const float *atompos, float *(&expandedpos))
int expand_basis_array (float *(&expandedbasis), int *(&numprims))
void compute_overlap_integrals (Timestep *ts, const float *expandedbasis, const int *numprims, float *(&overlap_matrix))
 Compute electronic overlap integrals Sij. Memory for overlap_matrix will be allocated. More...

int mullikenpop (Timestep *ts, int iwavesig, const float *expandedbasis, const int *numprims)
double localization_orbital_change (int orbid1, int orbid2, const float *C, const float *S)
float pair_localization_measure (int num_orbitals, int orbid1, int orbid2, const float *C, const float *S)
double mean_localization_measure (int num_orbitals, const float *C, const float *S)
void rotate_2x2_orbitals (float *C, int orbid1, int orbid2, double gamma)
 Perform 2x2 rotation of orbital pair. Equations (10a) and (10b) in Pipek & Mezey (1989). More...

double localization_rotation_angle (const float *C, const float *S, int orbid1, int orbid2)
 Rotation angle gamma for the 2 by 2 orbital rotations Equation (9) in Boughton & Pulay (1993). More...

double gross_atom_mulliken_pop (const float *C, const float *S, int atom, int orbid) const
 Return the gross Mulliken population of orbital i on atom A. More...

double orb_pair_gross_atom_mulliken_pop (const float *C, const float *S, int atom, int orbid1, int orbid2)
int orblocalize (Timestep *ts, int waveid, const float *expandedbasis, const int *numprims)
 Create a new wavefunction object based on existing wavefunction "waveid" with orbitals localized using the Pipek-Mezey algorithm. More...

Orbitalcreate_orbital (int iwave, int orbid, float *pos, QMTimestep *qmts)
 Create a new Orbital and return the pointer. User is responsible for deleting. More...


Public Attributes

int num_wave_f
 max. size of the wave_function array (stored in QMTimestep). The actual number of MOs present can be different for each frame but this is the maximum number of possible occupied and virtual orbitals, i.e. this is the number of contracted cartesian gaussian basis functions or the size of the secular equation. More...

int num_basis
 Number of the {exp, coeff) pairs in basis array. This is NOT the same as above since each (P, D, F, ...) shell consists of (3, 6, 10, ...) cartesian functions! More...

int num_wavef_signa
 total number of wavefunctions with different signatures (over all frames). More...

int runtype
 runtype for internal use. More...

int scftype
 scftype for internal use. More...

int status
 SCF and optimization status. More...

int nproc
 number of processors used. More...

int memory
 amount of memory used in MByte. More...

int num_electrons
 number of electrons. More...

int totalcharge
 total charge of system. More...

char * basis_string
 basis name as "nice" string. More...

char runtitle [QMDATA_BIGBUFSIZ]
 title of run. (large for Gaussian). More...

char geometry [QMDATA_BUFSIZ]
 type of provided geometry, e.g. UNIQUE, ZMT, CART, ... More...

char version_string [QMDATA_BUFSIZ]
 QM code version information. More...


Detailed Description

QM data management class. Concerns QM related data that are not timestep dependent.

Definition at line 281 of file QMData.h.


Constructor & Destructor Documentation

QMData::QMData int    natoms,
int    nbasis,
int    shells,
int    nwave
 

constructor:.

Definition at line 49 of file QMData.C.

References NULL, num_wavef_signa, QMSTATUS_UNKNOWN, runtype, RUNTYPE_UNKNOWN, scftype, SCFTYPE_UNKNOWN, and status.

QMData::~QMData void   
 

destructor.

Definition at line 79 of file QMData.C.

References basis_string, delete_basis_set, and num_wavef_signa.


Member Function Documentation

int QMData::assign_wavef_id int    type,
int    spin,
int    excit,
char *    info,
wavef_signa_t *&    signa_ts,
int &    numsig
 

Determines and returns a unique ID for each wavefuntion based on it's signature (type, spin, excitation, info).

Definition at line 864 of file QMData.C.

References wavef_signa_t::exci, wavef_signa_t::info, NULL, num_wavef_signa, QMDATA_BUFSIZ, wavef_signa_t::spin, and wavef_signa_t::type.

Referenced by QMTimestep::add_wavefunction.

int QMData::compare_wavef_guitype_to_type int    guitype,
int    type
 

Return 1 if the wavefunction type from the GUI matches the given type.

Definition at line 974 of file QMData.C.

References GUI_WAVEF_TYPE_CANON, GUI_WAVEF_TYPE_CINAT, GUI_WAVEF_TYPE_GEMINAL, GUI_WAVEF_TYPE_LOCAL, GUI_WAVEF_TYPE_MCSCFNAT, GUI_WAVEF_TYPE_MCSCFOPT, GUI_WAVEF_TYPE_OTHER, WAVE_BOYS, WAVE_CANON, WAVE_CINATUR, WAVE_GEMINAL, WAVE_MCSCFNAT, WAVE_MCSCFOPT, WAVE_PIPEK, WAVE_RUEDEN, and WAVE_UNKNOWN.

Referenced by find_wavef_id_from_gui_specs, get_highest_excitation, and has_wavef_guitype.

void QMData::compute_overlap_integrals Timestep   ts,
const float *    expandedbasis,
const int *    numprims,
float *&    overlap_matrix
 

Compute electronic overlap integrals Sij. Memory for overlap_matrix will be allocated.

Definition at line 1463 of file QMData.C.

References expand_atompos, get_overlap_matrix, NULL, num_wave_f, and Timestep::pos.

Referenced by molinfo_get, mullikenpop, and orblocalize.

Orbital * QMData::create_orbital int    iwave,
int    orbid,
float *    pos,
QMTimestep   qmts
 

Create a new Orbital and return the pointer. User is responsible for deleting.

Definition at line 1911 of file QMData.C.

References QMTimestep::get_wavecoeffs, num_basis, and num_wave_f.

int QMData::create_unique_basis Molecule   mol,
int    num_basis_atoms
 

Collapse basis set so that we have one atomic basis set per atom type.

Definition at line 384 of file QMData.C.

References BaseMolecule::atom, basis_atom_t::atomicnum, MolAtom::atomicnumber, shell_t::basis, prim_t::coeff, compare_atomic_basis, copy_atomic_basis, delete_basis_set, prim_t::expon, get_shell_type_str, num_basis, shell_t::num_cart_func, shell_t::numprims, basis_atom_t::numshells, shell_t::prim, basis_atom_t::shell, and shell_t::type.

Referenced by init_basis.

void QMData::delete_basis_set  
 

Free memory of the basis set.

Definition at line 112 of file QMData.C.

References NULL, basis_atom_t::numshells, shell_t::prim, and basis_atom_t::shell.

Referenced by create_unique_basis, and ~QMData.

int QMData::expand_atompos const float *    atompos,
float *&    expandedpos
 

Definition at line 1148 of file QMData.C.

References n, num_wave_f, basis_atom_t::numshells, basis_atom_t::shell, shell_t::type, and z.

Referenced by compute_overlap_integrals.

int QMData::expand_basis_array float *&    expandedbasis,
int *&    numprims
 

Definition at line 1183 of file QMData.C.

References shell_t::basis, prim_t::expon, n, shell_t::norm_fac, shell_t::num_cart_func, num_wave_f, shell_t::numprims, basis_atom_t::numshells, shell_t::prim, basis_atom_t::shell, and shell_t::type.

Referenced by VMDApp::molecule_orblocalize, and molinfo_get.

int QMData::find_wavef_id_from_gui_specs int    type,
int    spin,
int    exci
 

Matches the wavefunction selected from the GUI buttons (type, spin, excitation) with a wavefunction IDtag. Returns -1 if no such wavefunction exists.

Definition at line 947 of file QMData.C.

References compare_wavef_guitype_to_type, wavef_signa_t::exci, num_wavef_signa, and wavef_signa_t::spin.

Referenced by GraphicsFltkRepOrbital::regen_excitationlist, and GraphicsFltkRepOrbital::regen_orbitallist.

int QMData::get_angular_momentum int    atom,
int    shell,
int    mom,
int    comp
 

Definition at line 688 of file QMData.C.

References get_basis, and get_wave_offset.

char * QMData::get_angular_momentum_str int    atom,
int    shell,
int    mom
const
 

Definition at line 732 of file QMData.C.

References get_wave_offset, and NULL.

Referenced by molinfo_get.

const int* QMData::get_atom_types   const [inline]
 

Definition at line 450 of file QMData.h.

int QMData::get_avail_occupancies int    iwavesig,
float *&    orbocc
 

Get the occupancies of all available orbitals for the given wavefunction.

Definition at line 1112 of file QMData.C.

References wavef_signa_t::max_avail_orbs, num_wavef_signa, and wavef_signa_t::orbocc.

int QMData::get_avail_orbitals int    iwavesig,
int *&    orbids
 

Get IDs of all orbitals available for the given wavefunction. iwavesig is the index of the wavefunction signature. Returns 1 upon success, 0 otherwise.

Definition at line 1098 of file QMData.C.

References wavef_signa_t::max_avail_orbs, num_wavef_signa, and wavef_signa_t::orbids.

const shell_t * QMData::get_basis int    atom,
int    shell
const
 

Get basis set for a shell.

Definition at line 611 of file QMData.C.

References NULL.

const basis_atom_t * QMData::get_basis int    atom = 0 const
 

Get basis set for an atom.

Definition at line 603 of file QMData.C.

References NULL.

Referenced by get_angular_momentum, molinfo_get, and Wavefunction::sort_wave_coefficients.

const double* QMData::get_carthessian   const [inline]
 

Definition at line 484 of file QMData.h.

int QMData::get_highest_excitation int    guitype
 

Get the highest excitation for any wavefunction with the given type.

Definition at line 1024 of file QMData.C.

References compare_wavef_guitype_to_type, wavef_signa_t::exci, and num_wavef_signa.

Referenced by GraphicsFltkRepOrbital::regen_excitationlist.

const int* QMData::get_imagmodes   const [inline]
 

Definition at line 489 of file QMData.h.

const float* QMData::get_intensities   const [inline]
 

Definition at line 488 of file QMData.h.

const double* QMData::get_inthessian   const [inline]
 

Definition at line 485 of file QMData.h.

int QMData::get_max_avail_orbitals int    iwavesig
 

Return the maximum number of available orbitals for the given wavefunction over all frames Can be used to determine the number of orbitals to be displayed in the GUI. Returns -1 if requested wavefunction doesn't exist.

Definition at line 1089 of file QMData.C.

References wavef_signa_t::max_avail_orbs, and num_wavef_signa.

Referenced by molinfo_get, and GraphicsFltkRepOrbital::regen_orbitallist.

const float* QMData::get_normalmodes   const [inline]
 

Definition at line 486 of file QMData.h.

int QMData::get_num_atoms   [inline]
 

Definition at line 454 of file QMData.h.

Referenced by Wavefunction::sort_wave_coefficients.

int QMData::get_num_imag   [inline]
 

Definition at line 476 of file QMData.h.

Referenced by molinfo_get.

int QMData::get_num_intcoords   [inline]
 

Definition at line 477 of file QMData.h.

Referenced by molinfo_get.

const int* QMData::get_num_prim_per_shell   const [inline]
 

Definition at line 452 of file QMData.h.

int QMData::get_num_shells   [inline]
 

Definition at line 456 of file QMData.h.

Referenced by molinfo_get.

const int* QMData::get_num_shells_per_atom   const [inline]
 

Definition at line 451 of file QMData.h.

Referenced by Wavefunction::sort_wave_coefficients.

int QMData::get_num_types   [inline]
 

Definition at line 455 of file QMData.h.

int QMData::get_num_wavecoeff_per_atom int    atom const
 

Definition at line 619 of file QMData.C.

References shell_t::num_cart_func, basis_atom_t::numshells, and basis_atom_t::shell.

Referenced by Wavefunction::density_matrix, gross_atom_mulliken_pop, mullikenpop, orb_pair_gross_atom_mulliken_pop, and Wavefunction::population_matrix.

int QMData::get_orbital_label_from_gui_index int    iwavesig,
int    iorb
 

For the given wavefunction signature return the iorb-th orbital label. Used to translate from the GUI list of available orbitals the unique orbital label.

Definition at line 1128 of file QMData.C.

References wavef_signa_t::max_avail_orbs, num_wavef_signa, and wavef_signa_t::orbids.

const char * QMData::get_runtype_string void    const
 

Translate the runtype constant into a string.

Definition at line 800 of file QMData.C.

References RUNTYPE_DYNAMICS, RUNTYPE_ENERGY, RUNTYPE_GRADIENT, RUNTYPE_HESSIAN, RUNTYPE_MEX, RUNTYPE_OPTIMIZE, RUNTYPE_PROPERTIES, RUNTYPE_SADPOINT, and RUNTYPE_SURFACE.

const char * QMData::get_scftype_string void    const
 

Translate the scftype constant into a string.

Definition at line 818 of file QMData.C.

References SCFTYPE_FF, SCFTYPE_GVB, SCFTYPE_MCSCF, SCFTYPE_NONE, SCFTYPE_RHF, SCFTYPE_ROHF, and SCFTYPE_UHF.

const char * QMData::get_shell_type_str const shell_t   shell
 

Get shell type letter (S, P, D, F, ...) followed by '\0'.

Definition at line 655 of file QMData.C.

References shell_t::type.

Referenced by create_unique_basis, and molinfo_get.

const char * QMData::get_status_string  
 

Get status of SCF and optimization convergence.

Definition at line 834 of file QMData.C.

References QMSTATUS_FILE_TRUNCATED, QMSTATUS_OPT_CONV, QMSTATUS_OPT_NOT_CONV, QMSTATUS_SCF_NOT_CONV, and status.

int QMData::get_wave_offset int    atom,
int    shell = 0
const
 

Get the offset in the wavefunction array for given shell.

Definition at line 635 of file QMData.C.

Referenced by Wavefunction::density_matrix, get_angular_momentum, get_angular_momentum_str, gross_atom_mulliken_pop, molinfo_get, orb_pair_gross_atom_mulliken_pop, Wavefunction::population_matrix, set_angular_momentum, and set_angular_momentum_str.

const float* QMData::get_wavenumbers   const [inline]
 

Definition at line 487 of file QMData.h.

double QMData::gross_atom_mulliken_pop const float *    C,
const float *    S,
int    atom,
int    orbid
const
 

Return the gross Mulliken population of orbital i on atom A.

Definition at line 1846 of file QMData.C.

References get_num_wavecoeff_per_atom, get_wave_offset, and num_wave_f.

Referenced by localization_orbital_change, localization_rotation_angle, mean_localization_measure, and pair_localization_measure.

int QMData::has_orbital int    iwavesig,
int    orbid
 

Return 1 if the given wavefunction has an orbital with ID orbid in any frame.

Definition at line 1137 of file QMData.C.

References wavef_signa_t::max_avail_orbs, num_wavef_signa, and wavef_signa_t::orbids.

Referenced by GraphicsFltkRepOrbital::regen_orbitallist.

int QMData::has_wavef_exci int    exci
 

Return 1 if we have any wavefunction with the given excitation.

Definition at line 1000 of file QMData.C.

References wavef_signa_t::exci, and num_wavef_signa.

int QMData::has_wavef_guitype int    guitype
 

Return 1 if we have a wavefunction with the given GUI_WAVEF_TYPE_*.

Definition at line 964 of file QMData.C.

References compare_wavef_guitype_to_type, and num_wavef_signa.

Referenced by GraphicsFltkRepOrbital::regen_excitationlist, and GraphicsFltkRepOrbital::regen_wavefunctypes.

int QMData::has_wavef_signa int    type,
int    spin,
int    exci
 

Return 1 if we have any wavefunction with the given signature (type, spin, and excitation).

Definition at line 1011 of file QMData.C.

References wavef_signa_t::exci, num_wavef_signa, wavef_signa_t::spin, and wavef_signa_t::type.

int QMData::has_wavef_spin int    spin
 

Return 1 if we have any wavefunction with the given spin.

Definition at line 990 of file QMData.C.

References num_wavef_signa, and wavef_signa_t::spin.

Referenced by GraphicsFltkRepOrbital::regen_wavefunctypes.

void QMData::init_angular_norm_factors  
 

Initialize table of normalization factors containing different factors for each shell and its cartesian functions.

Definition at line 570 of file QMData.C.

References fac.

Referenced by init_basis.

int QMData::init_basis Molecule   mol,
int    num_basis_atoms,
const char *    string,
const float *    basis,
const int *    atomic_numbers,
const int *    nshells,
const int *    nprim,
const int *    shelltypes
 

Populate basis set data and organize them into hierarcical data structures.

Definition at line 171 of file QMData.C.

References basis_atom_t::atomicnum, shell_t::basis, basis_string, prim_t::coeff, shell_t::combo, create_unique_basis, prim_t::expon, init_angular_norm_factors, shell_t::norm_fac, normalize_basis, num_basis, shell_t::num_cart_func, shell_t::numprims, basis_atom_t::numshells, shell_t::prim, basis_atom_t::shell, SP_P_SHELL, SP_S_SHELL, SPD_D_SHELL, SPD_P_SHELL, SPD_S_SHELL, and shell_t::type.

Referenced by MolFilePlugin::read_qm_data.

void QMData::init_electrons Molecule   mol,
int    totcharge
 

Initialize total charge, multiplicity and everything that can be derived from that.

Definition at line 130 of file QMData.C.

References BaseMolecule::atom, MolAtom::atomicnumber, num_electrons, scftype, SCFTYPE_RHF, SCFTYPE_ROHF, SCFTYPE_UHF, and totalcharge.

Referenced by MolFilePlugin::read_qm_data.

double QMData::localization_orbital_change int    orbid1,
int    orbid2,
const float *    C,
const float *    S
 

Definition at line 1713 of file QMData.C.

References gross_atom_mulliken_pop, and orb_pair_gross_atom_mulliken_pop.

Referenced by orblocalize.

double QMData::localization_rotation_angle const float *    C,
const float *    S,
int    orbid1,
int    orbid2
 

Rotation angle gamma for the 2 by 2 orbital rotations Equation (9) in Boughton & Pulay (1993).

Definition at line 1803 of file QMData.C.

References gross_atom_mulliken_pop, orb_pair_gross_atom_mulliken_pop, and VMD_PI.

Referenced by orblocalize.

double QMData::mean_localization_measure int    num_orbitals,
const float *    C,
const float *    S
 

Definition at line 1756 of file QMData.C.

References gross_atom_mulliken_pop.

Referenced by orblocalize.

int QMData::mullikenpop Timestep   ts,
int    iwavesig,
const float *    expandedbasis,
const int *    numprims
 

Definition at line 1500 of file QMData.C.

References compute_overlap_integrals, Wavefunction::density_matrix, Wavefunction::get_num_coeffs, get_num_wavecoeff_per_atom, QMTimestep::get_wavef_index, QMTimestep::get_wavefunction, mm_mul, num_wavef_signa, Wavefunction::population_matrix, and Timestep::qm_timestep.

void QMData::normalize_basis  
 

Multiply contraction coefficients with the shell dependent part of the normalization factor.

Definition at line 533 of file QMData.C.

References shell_t::basis, prim_t::coeff, prim_t::expon, shell_t::numprims, basis_atom_t::numshells, shell_t::prim, basis_atom_t::shell, shell_t::type, and VMD_PI.

Referenced by init_basis.

double QMData::orb_pair_gross_atom_mulliken_pop const float *    C,
const float *    S,
int    atom,
int    orbid1,
int    orbid2
 

Definition at line 1870 of file QMData.C.

References get_num_wavecoeff_per_atom, get_wave_offset, and num_wave_f.

Referenced by localization_orbital_change, and localization_rotation_angle.

int QMData::orblocalize Timestep   ts,
int    waveid,
const float *    expandedbasis,
const int *    numprims
 

Create a new wavefunction object based on existing wavefunction "waveid" with orbitals localized using the Pipek-Mezey algorithm.

Definition at line 1601 of file QMData.C.

References QMTimestep::add_wavefunction, compute_overlap_integrals, Wavefunction::get_coeffs, Wavefunction::get_excitation, Wavefunction::get_multiplicity, Wavefunction::get_num_occupied_double, Wavefunction::get_spin, QMTimestep::get_wavef_index, QMTimestep::get_wavefunction, localization_orbital_change, localization_rotation_angle, mean_localization_measure, NULL, num_wave_f, num_wavef_signa, pair_localization_measure, Timestep::qm_timestep, rotate_2x2_orbitals, and WAVE_PIPEK.

Referenced by VMDApp::molecule_orblocalize.

float QMData::pair_localization_measure int    num_orbitals,
int    orbid1,
int    orbid2,
const float *    C,
const float *    S
 

Definition at line 1732 of file QMData.C.

References gross_atom_mulliken_pop.

Referenced by orblocalize.

void QMData::rotate_2x2_orbitals float *    C,
int    orbid1,
int    orbid2,
double    gamma
 

Perform 2x2 rotation of orbital pair. Equations (10a) and (10b) in Pipek & Mezey (1989).

Definition at line 1782 of file QMData.C.

References num_wave_f.

Referenced by orblocalize.

int QMData::set_angular_momenta const int *    angmom
 

Definition at line 664 of file QMData.C.

References num_wave_f.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_angular_momentum int    atom,
int    shell,
int    mom,
int *    pow
 

Definition at line 671 of file QMData.C.

References get_wave_offset.

void QMData::set_angular_momentum_str int    atom,
int    shell,
int    mom,
const char *    tag
 

Definition at line 699 of file QMData.C.

References get_wave_offset.

void QMData::set_carthessian int    numcartcoords,
double *    array
 

Definition at line 754 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_imagmodes int    numimag,
int *    array
 

Definition at line 785 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_intensities int    numcart,
float *    array
 

Definition at line 779 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_inthessian int    numintcoords,
double *    array
 

Definition at line 760 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_normalmodes int    numcart,
float *    array
 

Definition at line 767 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::set_wavenumbers int    numcart,
float *    array
 

Definition at line 773 of file QMData.C.

Referenced by MolFilePlugin::read_qm_data.

void QMData::update_avail_orbs int    iwavesig,
int    norbitals,
const int *    orbids,
const float *    orbocc
 

Update the list of available orbitals.

Definition at line 1046 of file QMData.C.

References wavef_signa_t::max_avail_orbs, wavef_signa_t::orbids, and wavef_signa_t::orbocc.

Referenced by QMTimestep::add_wavefunction.


Member Data Documentation

char* QMData::basis_string
 

basis name as "nice" string.

Definition at line 395 of file QMData.h.

Referenced by init_basis, and ~QMData.

char QMData::geometry[QMDATA_BUFSIZ]
 

type of provided geometry, e.g. UNIQUE, ZMT, CART, ...

Definition at line 397 of file QMData.h.

Referenced by MolFilePlugin::read_qm_data.

int QMData::memory
 

amount of memory used in MByte.

Definition at line 389 of file QMData.h.

Referenced by MolFilePlugin::read_qm_data.

int QMData::nproc
 

number of processors used.

Definition at line 388 of file QMData.h.

Referenced by MolFilePlugin::read_qm_data.

int QMData::num_basis
 

Number of the {exp, coeff) pairs in basis array. This is NOT the same as above since each (P, D, F, ...) shell consists of (3, 6, 10, ...) cartesian functions!

Definition at line 292 of file QMData.h.

Referenced by create_orbital, create_unique_basis, init_basis, and Wavefunction::sort_wave_coefficients.

int QMData::num_electrons
 

number of electrons.

Definition at line 390 of file QMData.h.

Referenced by QMTimestep::add_wavefunction, and init_electrons.

int QMData::num_wave_f
 

max. size of the wave_function array (stored in QMTimestep). The actual number of MOs present can be different for each frame but this is the maximum number of possible occupied and virtual orbitals, i.e. this is the number of contracted cartesian gaussian basis functions or the size of the secular equation.

Definition at line 284 of file QMData.h.

Referenced by compute_overlap_integrals, create_orbital, expand_atompos, expand_basis_array, gross_atom_mulliken_pop, molinfo_get, orb_pair_gross_atom_mulliken_pop, orblocalize, rotate_2x2_orbitals, set_angular_momenta, and Wavefunction::sort_wave_coefficients.

int QMData::num_wavef_signa
 

total number of wavefunctions with different signatures (over all frames).

Definition at line 296 of file QMData.h.

Referenced by assign_wavef_id, find_wavef_id_from_gui_specs, get_avail_occupancies, get_avail_orbitals, get_highest_excitation, get_max_avail_orbitals, get_orbital_label_from_gui_index, has_orbital, has_wavef_exci, has_wavef_guitype, has_wavef_signa, has_wavef_spin, molinfo_get, mullikenpop, orblocalize, QMData, and ~QMData.

char QMData::runtitle[QMDATA_BIGBUFSIZ]
 

title of run. (large for Gaussian).

Definition at line 396 of file QMData.h.

Referenced by MolFilePlugin::read_qm_data.

int QMData::runtype
 

runtype for internal use.

Definition at line 385 of file QMData.h.

Referenced by QMData, and MolFilePlugin::read_qm_data.

int QMData::scftype
 

scftype for internal use.

Definition at line 386 of file QMData.h.

Referenced by QMTimestep::add_wavefunction, init_electrons, QMData, and MolFilePlugin::read_qm_data.

int QMData::status
 

SCF and optimization status.

Definition at line 387 of file QMData.h.

Referenced by get_status_string, QMData, and MolFilePlugin::read_qm_data.

int QMData::totalcharge
 

total charge of system.

Definition at line 391 of file QMData.h.

Referenced by init_electrons.

char QMData::version_string[QMDATA_BUFSIZ]
 

QM code version information.

Definition at line 399 of file QMData.h.

Referenced by MolFilePlugin::read_qm_data.


The documentation for this class was generated from the following files:
Generated on Sat Apr 20 02:45:28 2024 for VMD (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002