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

QMData.h

Go to the documentation of this file.
00001 /***************************************************************************
00002  *cr                                                                       
00003  *cr            (C) Copyright 1995-2019 The Board of Trustees of the           
00004  *cr                        University of Illinois                       
00005  *cr                         All Rights Reserved                        
00006  *cr                                                                   
00007  ***************************************************************************/
00008 
00009 /***************************************************************************
00010  * RCS INFORMATION:
00011  *
00012  *      $RCSfile: QMData.h,v $
00013  *      $Author: johns $        $Locker:  $             $State: Exp $
00014  *      $Revision: 1.90 $       $Date: 2020/02/24 18:23:52 $
00015  *
00016  ***************************************************************************
00017  * DESCRIPTION:
00018  *
00019  * The QMData class, which stores all QM calculation data that
00020  * are not dependent on the timestep (the latter are handled by
00021  * the QMTimestep class).
00022  * These data are, for instance, basis set and calculation input
00023  * parameters such as method and multiplicity.
00024  * Note that the wavefunctions are stored in QMTimestep but that
00025  * a signature of all wavefunctions that occur in the trajectory
00026  * is kept here in QMData. The signature is needed for setting up
00027  * the Orbital representation GUI. 
00028  *
00029  * Basis set storage
00030  * -----------------
00031  * The basis set data are stored in hierarchical data structures
00032  * for convenient access and better readability of the non 
00033  * performance critical code. However, the actual orbital
00034  * computation (performed by the Orbital class) needs simple
00035  * contiguous arrays which is why we keep those too. Since we
00036  * have only one basis set per molecule this redundant storage
00037  * should not be too memory demanding.
00038  *
00039  * Basis set normalization
00040  * -----------------------
00041  * The general functional form of a normalized Gaussian-type
00042  * orbital basis function in cartesian coordinates is:
00043  *
00044  * Phi_GTO(x,y,z;a,i,j,k) = (2a/pi)^(3/4) * sqrt[(8a)^L]
00045  *          * sqrt[i!j!k!/((2i)!(2j)!(2k)!)]
00046  *          * x^i * y^j * z^k * exp[-a(x^2+y^2+z^2)]
00047  *
00048  * where x,y,z are the coordinates, a is the primitive exponent
00049  * and i,j,k are non-negative integers controlling the cartesian
00050  * shape of the GTO. For instance:
00051  * S  -> (0,0,0)
00052  * Px -> (1,0,0)
00053  * Py -> (0,1,0)
00054  * Dxx-> (2,0,0)
00055  * L=i+j+k denotes the shell type with L={0,1,2,3,...} for
00056  * {S,P,D,F,...} shells.
00057  * The contracted GTO is then defined by a linear combination of
00058  * GTO primitives
00059  *
00060  * Phi_CGTO = sum[c_p * Phi_GTO(x,y,z;a_p,i,j,k)]
00061  *
00062  * where c_p and a_p are the contraction coefficient and the
00063  * exponent for primitive p.
00064  *
00065  * The 3rd line in the expression for Phi_GTO is the actual
00066  * basis function while the first line shell-type dependent part
00067  * of the normalization factor and the second line comprises the 
00068  * angular momentum dependent part of the normalization.
00069  *
00070  * Finally, the wavefunction Psi for an orbital is a linear
00071  * combination of contracted basis functions Phi_CGTO.
00072  *
00073  * Psi = sum[c_k*Phi_CGTO_k]
00074  *
00075  * I.e. each wavefunction coefficient c_k must be multiplied with
00076  * the corresponding Phi_CGTO_k.
00077  *
00078  * Conceptionally it would be straightforward to just generate
00079  * an array with all Phi_CGTO_k. But for a high performance
00080  * implementation where fast memory is precious (e.g. constant
00081  * memory on GPUs) this is not efficient, especially for basis
00082  * sets with a large number of primitives per shell.
00083  * We have only one pair of contraction coefficient and exponent
00084  * per primitive but several cartesian components such as
00085  * Px, Py, Pz. Hence, we can multiply the shell dependent part
00086  * of the normalization factor into the contraction coefficient.
00087  * This is done in normalize_basis() which is called at the end
00088  * of init_basis(). For the angular momentum dependent part we
00089  * create a lookup table in init_angular_norm_factors() which 
00090  * is factored into the wavefunction coefficients when an Orbital
00091  * object is created for rendering purpose.
00092  *
00093  ***************************************************************************/
00094 #ifndef QMDATA_H
00095 #define QMDATA_H
00096 
00097 
00098 #define GUI_WAVEF_TYPE_CANON    0
00099 #define GUI_WAVEF_TYPE_GEMINAL  1
00100 #define GUI_WAVEF_TYPE_MCSCFNAT 2
00101 #define GUI_WAVEF_TYPE_MCSCFOPT 3
00102 #define GUI_WAVEF_TYPE_CINAT    4
00103 #define GUI_WAVEF_TYPE_LOCAL    5
00104 #define GUI_WAVEF_TYPE_OTHER    6
00105 
00106 #define GUI_WAVEF_SPIN_ALPHA    0
00107 #define GUI_WAVEF_SPIN_BETA     1
00108 
00109 #define GUI_WAVEF_EXCI_GROUND   0
00110 
00111 // NOTE: this has to be kept in sync with the corresponding table
00112 //       in the QM molfile plugin reader. 
00113 // XXX   These should all be added to the molfile_plugin header
00114 //       When copied from the plugin to the internal data structure,
00115 //       we can either arrange that the molfile plugin constants and these
00116 //       are always identical, or we can use a switch() block to translate 
00117 //       between them if they ever diverge.
00118 #define SPD_D_SHELL -5
00119 #define SPD_P_SHELL -4
00120 #define SPD_S_SHELL -3
00121 #define SP_S_SHELL  -2
00122 #define SP_P_SHELL  -1
00123 #define S_SHELL      0
00124 #define P_SHELL      1
00125 #define D_SHELL      2
00126 #define F_SHELL      3
00127 #define G_SHELL      4
00128 #define H_SHELL      5
00129 #define I_SHELL      6
00130 
00131 #define QMDATA_BUFSIZ              81 
00132 #define QMDATA_BUFSTRFMT      "%.80s" 
00133 #define QMDATA_BIGBUFSIZ         8192 
00134 #define QMDATA_BIGBUFSTRFMT "%.8191s" 
00135 
00136 // forward declarations.
00137 class Orbital;
00138 class QMTimestep;
00139 class Timestep;
00140 class Molecule;
00141 
00143 enum {
00144   QMSTATUS_UNKNOWN        = -1, // don't know
00145   QMSTATUS_OPT_CONV       =  0, // optimization converged
00146   QMSTATUS_SCF_NOT_CONV   =  1, // SCF convergence failed
00147   QMSTATUS_OPT_NOT_CONV   =  2, // optimization not converged
00148   QMSTATUS_FILE_TRUNCATED =  3  // file was truncated
00149 };
00150 
00152 enum {
00153   SCFTYPE_UNKNOWN = -1, // we have no info about the method used
00154   SCFTYPE_NONE    =  0, // calculation didn't make use of SCF at all
00155   SCFTYPE_RHF     =  1, // Restricted Hartree-Fock
00156   SCFTYPE_UHF     =  2, // unrestricted Hartree-Fock
00157   SCFTYPE_ROHF    =  3, // restricted open-shell Hartree-Fock
00158   SCFTYPE_GVB     =  4, // generalized valence bond orbitals
00159   SCFTYPE_MCSCF   =  5, // multi-configuration SCF
00160   SCFTYPE_FF      =  6  // classical force-field based sim.
00161 };
00162 
00164 enum {
00165   RUNTYPE_UNKNOWN    = 0,
00166   RUNTYPE_ENERGY     = 1, // single point energy evaluation
00167   RUNTYPE_OPTIMIZE   = 2, // geometry optimization
00168   RUNTYPE_SADPOINT   = 3, // saddle point search
00169   RUNTYPE_HESSIAN    = 4, // Hessian/frequency calculation
00170   RUNTYPE_SURFACE    = 5, // potential surface scan
00171   RUNTYPE_GRADIENT   = 6, // energy gradient calculation
00172   RUNTYPE_MEX        = 7, // minimum energy crossing
00173   RUNTYPE_DYNAMICS   = 8, // Any type of molecular dynamics
00174                           // e.g. Born-Oppenheimer, Car-Parinello,
00175                           // or classical MD
00176   RUNTYPE_PROPERTIES = 9  // Properties were calculated from a
00177                           // wavefunction that was read from file
00178 };
00179 
00181 enum {
00182   QMCHARGE_UNKNOWN  = 0, 
00183   QMCHARGE_MULLIKEN = 1, // Mulliken charges
00184   QMCHARGE_LOWDIN   = 2, // Lowdin charges
00185   QMCHARGE_ESP      = 3, // elstat. potential fitting
00186   QMCHARGE_NPA      = 4  // natural population analysis
00187 };
00188 
00189 
00217 enum {
00218   WAVE_CANON,    WAVE_GEMINAL, WAVE_MCSCFNAT,
00219   WAVE_MCSCFOPT, WAVE_CINATUR,
00220   WAVE_PIPEK,    WAVE_BOYS,    WAVE_RUEDEN,
00221   WAVE_NAO,      WAVE_PNAO,    WAVE_NHO, 
00222   WAVE_PNHO,     WAVE_NBO,     WAVE_PNBO, 
00223   WAVE_PNLMO,    WAVE_NLMO,    WAVE_MOAO, 
00224   WAVE_NATO,     WAVE_UNKNOWN
00225 };
00226 
00227 
00228 // Signature of a wavefunction, used to identify
00229 // a wavefunction thoughout the trajectory.
00230 // The orbital metadata are needed by the GUI to
00231 // determine what is the union of orbitals showing up
00232 // over the trajetory. It is displayed as the available
00233 // orbital list in the representations window.
00234 typedef struct {
00235   int type;
00236   int spin;
00237   int exci;
00238   char info[QMDATA_BUFSIZ];
00239 
00240   int max_avail_orbs; // maximum number of available orbitals
00241                       // over all frames.
00242   int *orbids;        // list of existing orbital IDs
00243   float *orbocc;      // occupancies of all available orbitals
00244 } wavef_signa_t;
00245 
00246 
00247 // Basis set data structures
00248 // -------------------------
00249 
00251 typedef struct {
00252   float expon;       
00253   float coeff;       
00254 } prim_t;
00255 
00257 typedef struct {
00258   int numprims;      
00259   int type;          
00260   int combo;         
00261 
00262 
00263   int num_cart_func; 
00264 
00265 
00266 
00267   float *norm_fac;   
00268   float *basis;      
00269   prim_t *prim;      
00270 } shell_t;
00271 
00273 typedef struct {
00274   int atomicnum;     
00275   int numshells;     
00276   shell_t *shell;    
00277 } basis_atom_t;
00278 
00281 class QMData {
00282  public:
00283   // Orbital data
00284   int num_wave_f;      
00285 
00286 
00287 
00288 
00289 
00290 
00291 
00292   int num_basis;       
00293 
00294 
00295 
00296   int num_wavef_signa; 
00297 
00298 
00299  private:
00300   int num_types;            
00301   int num_atoms;            
00302 
00303   int *atom_types;          
00304   int *atom_sort;           
00305   int *num_shells_per_atom; 
00306 
00307   int *wave_offset;         
00308 
00309 
00310 
00311 
00312 
00313   int num_shells;           
00314   int *num_prim_per_shell;  
00315 
00316   int *angular_momentum;    
00317 
00318 
00319 
00320 
00321   int highest_shell;        
00322 
00323 
00324   float **norm_factors;     
00325 
00326 
00327 
00328 
00329 
00330 
00331   float *basis_array;       
00332 
00333 
00334 
00335 
00336 
00337 
00338 
00339 
00340 
00341 
00342 
00343   int *atom_basis;          
00344 
00345   int *shell_types;         
00346 
00347 
00348 
00349 
00350 
00351 
00352   basis_atom_t *basis_set;  
00353 
00354 
00355 
00356 
00357 
00358 
00359   
00360 
00361   wavef_signa_t *wavef_signa; 
00362 
00363 
00364 
00365   // Hessian and frequency data.
00366   // These are actually timestep dependent but in almost all practical
00367   // cases they are just printed for the last step at the end of a calculation. 
00368   int nimag;                
00369   int nintcoords;           
00370 
00371   double *carthessian;     
00372 
00373 
00374   double *inthessian;      
00375 
00376 
00377   float *wavenumbers;      
00378   float *intensities;      
00379   float *normalmodes;      
00380   int   *imagmodes;        
00381 
00382 
00383 public:
00384   // QM run info
00385   int runtype;             
00386   int scftype;             
00387   int status;              
00388   int nproc;               
00389   int memory;              
00390   int num_electrons;       
00391   int totalcharge;         
00392   //int num_orbitals_A;      ///< number of alpha orbitals
00393   //int num_orbitals_B;      ///< number of beta orbitals
00394 
00395   char *basis_string;                  
00396   char runtitle[QMDATA_BIGBUFSIZ];     
00397   char geometry[QMDATA_BUFSIZ];        
00398 
00399   char version_string[QMDATA_BUFSIZ];  
00400 
00401   QMData(int natoms, int nbasis, int shells, int nwave); 
00402   ~QMData(void);                       
00403 
00406   void init_electrons(Molecule *mol, int totcharge);
00407 
00408   // ====================================
00409   // Functions for basis set initializion
00410   // ====================================
00411 
00413   void delete_basis_set();
00414 
00417   int init_basis(Molecule *mol, int num_basis_atoms, const char *string,
00418                  const float *basis, const int *atomic_numbers,
00419                  const int *nshells, const int *nprim,
00420                  const int *shelltypes);
00421 
00423   int create_unique_basis(Molecule *mol, int num_basis_atoms);
00424 
00427   void normalize_basis();
00428 
00429 #if 0
00430   // XXX: unused
00432   void sort_atoms_by_type();
00433 #endif
00434 
00437   void init_angular_norm_factors();
00438 
00439 
00440   // =================
00441   // Basis set acccess
00442   // =================
00443 
00445   const basis_atom_t* get_basis(int atom=0) const;
00446 
00448   const shell_t*      get_basis(int atom, int shell) const;
00449 
00450   const int* get_atom_types()          const { return atom_types; }
00451   const int* get_num_shells_per_atom() const { return num_shells_per_atom; }
00452   const int* get_num_prim_per_shell()  const { return num_prim_per_shell; }
00453 
00454   int get_num_atoms()                { return num_atoms; }
00455   int get_num_types()                { return num_types; }
00456   int get_num_shells()               { return num_shells; }
00457   int get_num_wavecoeff_per_atom(int atom) const;
00458 
00460   int get_wave_offset(int atom, int shell=0) const;
00461 
00463   const char* get_shell_type_str(const shell_t *shell);
00464 
00465   char* get_angular_momentum_str(int atom, int shell, int mom) const;
00466   void  set_angular_momentum_str(int atom, int shell, int mom, const char *tag);
00467   void  set_angular_momentum(int atom, int shell, int mom, int *pow);
00468   int   get_angular_momentum(int atom, int shell, int mom, int comp);
00469   int   set_angular_momenta(const int *angmom);
00470 
00471 
00472   // ========================
00473   // Hessian and normal modes
00474   // ========================
00475 
00476   int get_num_imag()       { return nimag; }
00477   int get_num_intcoords()  { return nintcoords; }
00478   void set_carthessian(int numcartcoords, double *array);
00479   void set_inthessian(int numintcoords, double *array);
00480   void set_normalmodes(int numcart, float *array);
00481   void set_wavenumbers(int numcart, float *array);
00482   void set_intensities(int numcart, float *array);
00483   void set_imagmodes(int numimag, int *array);
00484   const double* get_carthessian() const { return carthessian; }
00485   const double* get_inthessian()  const { return inthessian;  }
00486   const float*  get_normalmodes() const { return normalmodes; }
00487   const float*  get_wavenumbers() const { return wavenumbers; }
00488   const float*  get_intensities() const { return intensities; }
00489   const int*    get_imagmodes()   const { return imagmodes;   }
00490 
00491 
00492   // =====================
00493   // Calculation meta data
00494   // =====================
00495 
00497   const char *get_scftype_string(void) const;
00498 
00500   const char *get_runtype_string(void) const;
00501   
00503   const char* get_status_string();
00504 
00505 
00506   // =======================
00507   // Wavefunction signatures
00508   // =======================
00509 
00512   int assign_wavef_id(int type, int spin, int excit, char *info,
00513                       wavef_signa_t *(&signa_ts), int &numsig);
00514 
00518   int find_wavef_id_from_gui_specs(int type, int spin, int exci);
00519 
00521   int compare_wavef_guitype_to_type(int guitype, int type);
00522 
00524   int has_wavef_guitype(int guitype);
00525 
00527   int has_wavef_spin(int spin);
00528 
00530   int has_wavef_exci(int exci);
00531 
00534   int has_wavef_signa(int type, int spin, int exci);
00535 
00538   int get_highest_excitation(int guitype);
00539 
00540 
00541   // ========================================================
00542   // Functions dealing with the list of orbitals available 
00543   // for a given wavefunction. Needed for the GUI.
00544   // ========================================================
00545 
00547   void update_avail_orbs(int iwavesig, int norbitals,
00548                          const int *orbids, const float *orbocc);
00549 
00555   int get_max_avail_orbitals(int iwavesig);
00556 
00560   int get_avail_orbitals(int iwavesig, int *(&orbids));
00561 
00564   int get_avail_occupancies(int iwavesig, float *(&orbocc));
00565 
00569   int get_orbital_label_from_gui_index(int iwavesig, int iorb);
00570 
00573   int has_orbital(int iwavesig, int orbid);
00574 
00575 #define ANGS_TO_BOHR 1.88972612478289694072f
00576 
00577   int expand_atompos(const float *atompos,
00578                      float *(&expandedpos));
00579   int expand_basis_array(float *(&expandedbasis), int *(&numprims));
00580 
00581   void compute_overlap_integrals(Timestep *ts,
00582                         const float *expandedbasis,
00583                         const int *numprims, float *(&overlap_matrix));
00584 
00585   int mullikenpop(Timestep *ts, int iwavesig, 
00586                   const float *expandedbasis,
00587                   const int *numprims);
00588 
00589 
00590   double localization_orbital_change(int orbid1, int orbid2,
00591                                         const float *C,
00592                                         const float *S);
00593 
00594   float pair_localization_measure(int num_orbitals,
00595                                   int orbid1, int orbid2,
00596                                   const float *C, const float *S);
00597 
00598   double mean_localization_measure(int num_orbitals, const float *C,
00599                                   const float *S);
00600 
00603   void rotate_2x2_orbitals(float *C, int orbid1, int orbid2,
00604                            double gamma);
00605 
00608   double localization_rotation_angle(const float *C, const float *S,
00609                                     int orbid1, int orbid2);
00610 
00612   double gross_atom_mulliken_pop(const float *C, const float *S,
00613                                 int atom, int orbid) const;
00614 
00615   double orb_pair_gross_atom_mulliken_pop(const float *C, const float *S,
00616                                      int atom, int orbid1, int orbid2);
00617 
00620   int orblocalize(Timestep *ts, int waveid, const float *expandedbasis, 
00621                   const int *numprims);
00622 
00625   Orbital* create_orbital(int iwave, int orbid, float *pos,
00626                           QMTimestep *qmts);
00627 
00628 };
00629 
00630 
00631 
00632 #endif
00633 

Generated on Fri Oct 24 02:45:22 2025 for VMD (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002