// // gaussbas.h // // Copyright (C) 1996 Limit Point Systems, Inc. // // Author: Curtis Janssen // Maintainer: LPS // // This file is part of the SC Toolkit. // // The SC Toolkit is free software; you can redistribute it and/or modify // it under the terms of the GNU Library General Public License as published by // the Free Software Foundation; either version 2, or (at your option) // any later version. // // The SC Toolkit is distributed in the hope that it will be useful, // but WITHOUT ANY WARRANTY; without even the implied warranty of // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the // GNU Library General Public License for more details. // // You should have received a copy of the GNU Library General Public License // along with the SC Toolkit; see the file COPYING.LIB. If not, write to // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA. // // The U.S. Government is granted a limited license as per AL 91-7. // #ifndef _chemistry_qc_basis_gaussbas_h #define _chemistry_qc_basis_gaussbas_h #ifdef __GNUC__ #pragma interface #endif #include #include #include #include #include #include #include namespace sc { class GaussianShell; class BasisFileSet; class Integral; class CartesianIter; class SphericalTransformIter; /** The GaussianBasisSet class is used describe a basis set composed of atomic gaussian orbitals. Inputs for common basis sets are included in the MPQC distribution. They have been obtained from the EMSL Basis Set Database and translated into the MPQC format. The citation for this database is below. The technical citation for each basis set is listed in the individual basis set data files, in MPQC's lib/basis directory. Following is a table with available basis sets listing the supported elements for each basis and the number of basis functions for H, \f$n_0\f$, first row, \f$n_1\f$, and second row, \f$n_2\f$, atoms. Basis sets with non-alpha-numerical characters in their name must be given in quotes.
Basis SetElements\f$n_0\f$\f$n_1\f$\f$n_2\f$
STO-2GH-Ca159
STO-3GH-Kr159
STO-3G*H-Ar1514
STO-6GH-Kr159
MINI (Huzinaga)H-Ca159
MINI (Scaled)H-Ca159
MIDI (Huzinaga)H-Na, Al-K2913
DZ (Dunning)H, Li, B-Ne, Al-Cl21018
DZP (Dunning)H, Li, B-Ne, Al-Cl51624
DZP + Diffuse (Dunning)H, B-Ne619
3-21GH-Kr2913
3-21G*H-Ar2919
3-21++GH-Ar31317
3-21++G*H-Ar31323
4-31GH-Ne, P-Cl2913
6-31GH-Zn2913
6-31G*H-Zn21519
6-31G**H-Zn51519
6-31+G*H-Ar21923
6-31++GH-Ca31317
6-31++G*H-Ar31923
6-31++G**H-Ar61923
6-311GH-Ca, Ga-Kr31321
6-311G*H-Ca, Ga-Kr31826
6-311G**H-Ca, Ga-Kr61826
6-311G(2df,2pd)H-Ne, K, Ca1430
6-311++G**H-Ne722
6-311++G(2d,2p)H-Ca102735
6-311++G(3df,3pd)H-Ar183947
cc-pVDZH-Ar, Ca, Ga-Kr51418
cc-pVTZH-Ar, Ca, Ga-Kr143034
cc-pVQZH-Ar, Ca, Ga-Kr305559
cc-pV5ZH-Ar, Ca, Ga-Kr559195
cc-pV6ZH, He, B-Ne, Al-Ar91140144
aug-cc-pVDZH, He, B-Ne, Al-Ar, Ga-Kr92327
aug-cc-pVTZH, He, B-Ne, Al-Ar, Ga-Kr234650
aug-cc-pVQZH, He, B-Ne, Al-Ar, Ga-Kr468084
aug-cc-pV5ZH, He, B-Ne, Al-Ar, Ga-Kr80127131
aug-cc-pV6ZH, He, B-Ne, Al-Ar127189193
cc-pCVDZLi, B-Ar1827
cc-pCVTZLi, B-Ar4359
cc-pCVQZLi, B-Ar84109
cc-pCV5ZB-Ne145
aug-cc-pCVDZB-F, Al-Ar2736
aug-cc-pCVTZB-Ne, Al-Ar5975
aug-cc-pCVQZB-Ne, Al-Ar109134
aug-cc-pCV5ZB-F181
NASA Ames ANOH, B-Ne, Al, P, Ti, Fe, Ni305559
pc-0H, C-F, Si-Cl2913
pc-1H, C-F, Si-Cl51418
pc-2H, C-F, Si-Cl143034
pc-3H, C-F, Si-Cl346464
pc-4H, C-F, Si-Cl63109105
pc-0-augH, C-F, Si-Cl31317
pc-1-augH, C-F, Si-Cl92327
pc-2-augH, C-F, Si-Cl234650
pc-3-augH, C-F, Si-Cl508989
pc-4-augH, C-F, Si-Cl88145141
All basis sets except for the pc-n and pc-n-aug basis sets were obtained from the Extensible Computational Chemistry Environment Basis Set Database, Version 12/03/03, as developed and distributed by the Molecular Science Computing Facility, Environmental and Molecular Sciences Laboratory which is part of the Pacific Northwest Laboratory, P.O. Box 999, Richland, Washington 99352, USA, and funded by the U.S. Department of Energy. The Pacific Northwest Laboratory is a multi-program laboratory operated by Battelle Memorial Institute for the U.S. Department of Energy under contract DE-AC06-76RLO 1830. Contact David Feller or Karen Schuchardt for further information. The pc-n and pc-n-aug basis sets are the polarization consistent basis sets of Frank Jensen. See J. Chem. Phys. 115 (2001) 9113; J. Chem. Phys. 116 (2002) 7372; J. Chem. Phys. 117 (2002) 9234; and J. Chem. Phys. 121 (2004) 3463. */ class GaussianBasisSet: public SavableState { private: // nonnull if keyword "name" was provided char* name_; // same as name_ if name_!=0, else something else char* label_; GaussianShell** shell_; std::vector shell_to_function_; std::vector function_to_shell_; Ref molecule_; Ref matrixkit_; Ref so_matrixkit_; RefSCDimension basisdim_; int ncenter_; std::vector shell_to_center_; std::vector shell_to_primitive_; std::vector center_to_shell_; std::vector center_to_nshell_; std::vector center_to_nbasis_; int nshell_; int nbasis_; int nprim_; bool has_pure_; GaussianBasisSet(const char* name, const char* label, const Ref& molecule, const Ref& matrixkit, const RefSCDimension& basisdim, const int ncenter, const int nshell, GaussianShell** shell, const std::vector& center_to_nshell); // Counts shells in this basis for this chemical element int count_shells_(Ref& keyval, const char* elemname, const char* sbasisname, BasisFileSet& bases, int havepure, int pure, bool missing_ok); // Constructs this basis void get_shells_(int& ishell, Ref& keyval, const char* elemname, const char* sbasisname, BasisFileSet& bases, int havepure, int pure, bool missing_ok); // Counts shells in an even-tempered primitive basis int count_even_temp_shells_(Ref& keyval, const char* elemname, const char* sbasisname, int havepure, int pure); // Constructs an even-tempered primitive basis void get_even_temp_shells_(int& ishell, Ref& keyval, const char* elemname, const char* sbasisname, int havepure, int pure); // Constructs basis set specified as an array of shells void recursively_get_shell(int&,Ref&, const char*,const char*,BasisFileSet&, int,int,int,bool missing_ok); void init(Ref&,Ref&, BasisFileSet&, int have_userkeyval, int pure); void init2(int skip_ghosts=0,bool include_q=0); protected: GaussianBasisSet(const GaussianBasisSet&); virtual void set_matrixkit(const Ref&); public: /** This holds scratch data needed to compute basis function values. */ class ValueData { protected: CartesianIter **civec_; SphericalTransformIter **sivec_; int maxam_; public: ValueData(const Ref &, const Ref &); ~ValueData(); CartesianIter **civec() { return civec_; } SphericalTransformIter **sivec() { return sivec_; } }; /// This can be given to a CTOR to construct a unit basis function. enum UnitType {Unit}; /** The KeyVal constructor.
molecule
The gives a Molecule object. The is no default.
puream
If this boolean parameter is true then 5D, 7F, etc. will be used. Otherwise all cartesian functions will be used. The default depends on the particular basis set.
name
This is a string giving the name of the basis set. The above table of basis sets gives some of the recognized basis set names. It may be necessary to put the name in double quotes. There is no default.
basis
This is a vector of basis set names that can give a different basis set to each atom in the molecule. If the element vector is given, then it gives different basis sets to different elements. The default is to give every atom the basis set specified in name.
element
This is a vector of elements. If it is given then it must have the same number of entries as the basis vector.
basisdir
A string giving a directory where basis set data files are to be sought. See the text below for a complete description of what directories are consulted.
basisfiles
Each keyword in this vector of files is appended to the directory specified with basisdir and basis set data is read from them.
matrixkit
Specifies a SCMatrixKit object. It is usually not necessary to give this keyword, as the default action should get the correct SCMatrixKit.
Several files in various directories are checked for basis set data. First, basis sets can be given by the user in the basis section at the top level of the main input file. Next, if a path is given with the basisdir keyword, then all of the files given with the basisfiles keyword are read in after appending their names to the value of basisdir. Basis sets can be given in these files in the basis section at the top level as well. If the named basis set still cannot be found, then GaussianBasisSet will try convert the basis set name to a file name and check first in the directory given by basisdir. Next it checks for the environment variable SCLIBDIR. If it is set it will look for the basis file in $SCLIBDIR/basis. Otherwise it will look in the source code distribution in the directory SC/lib/basis. If the executable has changed machines or the source code has be moved, then it may be necessary to copy the library files to your machine and set the SCLIBDIR environmental variable. The basis set itself is also given in the ParsedKeyVal format. There are two recognized formats for basis sets:
array of shells
One must specify the keyword :basis: followed by the lowercase atomic name followed by : followed by the basis set name (which may need to be placed inside double quotes). The value for the keyword is an array of shells. Each shell reads the following keywords:
type
This is a vector that describes each component of this shell. For each element the following two keywords are read:
am
The angular momentum of the component. This can be given as the letter designation, s, p, d, etc. There is no default.
puream
If this boolean parameter is true then 5D, 7F, etc. shells are used. The default is false. This parameter can be overridden in the GaussianBasisSet specification.
exp
This is a vector giving the exponents of the primitive Gaussian functions.
coef
This is a matrix giving the coeffients of the primitive Gaussian functions. The first index gives the component number of the shell and the second gives the primitive number.
An example might be easier to understand. This is a basis set specificition for STO-2G carbon:
        basis: (
         carbon: "STO-2G": [
          (type: [(am = s)]
           {      exp      coef:0 } = {
              27.38503303 0.43012850
               4.87452205 0.67891353
           })
          (type: [(am = p) (am = s)]
           {     exp      coef:1     coef:0 } = {
               1.13674819 0.04947177 0.51154071
               0.28830936 0.96378241 0.61281990
           })
         ]
        )
        
basis set of even-tempered primitive Gaussians
Such basis set format is given as a group of keywords. The name of the group is :basis: followed by the lowercase atomic name followed by : followed by the basis set name (which may need to be placed inside double quotes). The group of keywords must contain vectors am and nprim, which specify the angular momentum and the number of shells in each block of even-tempered primitives. In addition, one must provide any two of the following vectors:
first_exp
The exponent of the "tightest" primitive Gaussian in the block.
last_exp
The exponent of the most "diffuse" primitive Gaussian in the block.
exp_ratio
The ratio of exponents of consecutive primitive Gaussians in the block.
Note that the dimensions of each vector must be the same. Here's an example of a basis set composed of 2 blocks of even-tempered s-functions and 1 block of even-tempered p-functions.
        basis: (
         neon: "20s5s13p":(

           am = [ 0 0 1 ]
           nprim = [ 20 5 13 ]
           first_exp = [ 1000.0 0.1  70.0 ]
           last_exp =  [    1.0 0.01  0.1 ]

         )
        )
        
*/ GaussianBasisSet(const Ref&); /** This can be given GaussianBasisSet::Unit to construct a basis set with a single basis function that is one everywhere. This can be used with integral evaluators to compute certain classes of integrals, with limitations. */ GaussianBasisSet(UnitType); GaussianBasisSet(StateIn&); virtual ~GaussianBasisSet(); /** Returns a GaussianBasisSet object that consists of the basis functions for each atom in this followed by the basis functions in B for the corresponding atom. The Molecule object for the two basis sets must be identical. */ Ref operator+(const Ref& B); void save_data_state(StateOut&); /// Return the name of the basis set (is nonnull only if keyword "name" was provided) const char* name() const { return name_; } /** Return the label of the basis set. label() return the same string as name() if keyword "name" was provided, otherwise it is a unique descriptive string which can be arbitrarily long. */ const char* label() const { if (name()) { return name(); } else { return label_; } } /// Return the Molecule object. Ref molecule() const { return molecule_; } /// Returns the SCMatrixKit that is to be used for AO bases. Ref matrixkit() { return matrixkit_; } /// Returns the SCMatrixKit that is to be used for SO bases. Ref so_matrixkit() { return so_matrixkit_; } /// Returns the SCDimension object for the dimension. RefSCDimension basisdim() { return basisdim_; } /// Return the number of centers. int ncenter() const; /// Return the number of shells. int nshell() const { return nshell_; } /// Return the number of shells on the given center. int nshell_on_center(int icenter) const; /** Return an overall shell number, given a center and the shell number on that center. */ int shell_on_center(int icenter, int shell) const; /// Return the center on which the given shell is located. int shell_to_center(int ishell) const { return shell_to_center_[ishell]; } /// Return the overall index of the first primitive from the given shell int shell_to_primitive(int ishell) const {return shell_to_primitive_[ishell]; } /// Return the number of basis functions. int nbasis() const { return nbasis_; } /// Return the number of basis functions on the given center. int nbasis_on_center(int icenter) const; /// Return the number of primitive Gaussians. int nprimitive() const { return nprim_; } /// Return true if basis contains solid harmonics Gaussians int has_pure() const { return has_pure_; } /// Return the maximum number of functions that any shell has. int max_nfunction_in_shell() const; /** Return the maximum number of Cartesian functions that any shell has. The optional argument is an angular momentum increment. */ int max_ncartesian_in_shell(int aminc=0) const; /// Return the maximum number of primitive Gaussian that any shell has. int max_nprimitive_in_shell() const; /// Return the highest angular momentum in any shell. int max_angular_momentum() const; /// Return the maximum number of Gaussians in a contraction in any shell. int max_ncontraction() const; /** Return the maximum angular momentum found in the given contraction number for any shell. */ int max_am_for_contraction(int con) const; /// Return the maximum number of Cartesian functions in any shell. int max_cartesian() const; /// Return the number of the first function in the given shell. int shell_to_function(int i) const { return shell_to_function_[i]; } /// Return the shell to which the given function belongs. int function_to_shell(int i) const; /// Return a reference to GaussianShell number i. const GaussianShell& operator()(int i) const { return *shell_[i]; } /// Return a reference to GaussianShell number i. GaussianShell& operator()(int i) { return *shell_[i]; } /// Return a reference to GaussianShell number i. const GaussianShell& operator[](int i) const { return *shell_[i]; } /// Return a reference to GaussianShell number i. GaussianShell& operator[](int i) { return *shell_[i]; } /// Return a reference to GaussianShell number i. const GaussianShell& shell(int i) const { return *shell_[i]; } /// Return a reference to GaussianShell number i. GaussianShell& shell(int i) { return *shell_[i]; } /// Return a reference to GaussianShell number ishell on center icenter. const GaussianShell& operator()(int icenter,int ishell) const; /// Return a reference to GaussianShell number ishell on center icenter. GaussianShell& operator()(int icenter,int ishell); /// Return a reference to GaussianShell number j on center i. const GaussianShell& shell(int i,int j) const { return operator()(i,j); } /// Return a reference to GaussianShell number j on center i. GaussianShell& shell(int i,int j) { return operator()(i,j); } /** The location of center icenter. The xyz argument is 0 for x, 1 for y, and 2 for z. */ double r(int icenter,int xyz) const; /** Compute the values for this basis set at position r. The basis_values argument must be vector of length nbasis. */ int values(const SCVector3& r, ValueData *, double* basis_values) const; /** Like values(...), but computes gradients of the basis function values, too. The g_values argument must be vector of length 3*nbasis. The data will be written in the order bf1_x, bf1_y, bf1_z, ... */ int grad_values(const SCVector3& r, ValueData *, double*g_values,double* basis_values=0) const; /** Like values(...), but computes first and second derivatives of the basis function values, too. h_values must be vector of length 6*nbasis. The data will be written in the order bf1_xx, bf1_yx, bf1_yy, bf1_zx, bf1_zy, bf1_zz, ... */ int hessian_values(const SCVector3& r, ValueData *, double *h_values, double*g_values=0,double* basis_values=0) const; /** Compute the values for the given shell functions at position r. See the other values(...) members for more information. */ int shell_values(const SCVector3& r, int sh, ValueData *, double* basis_values) const; /** Like values(...), but computes gradients of the shell function values, too. See the other grad_values(...) members for more information. */ int grad_shell_values(const SCVector3& r, int sh, ValueData *, double*g_values, double* basis_values=0) const; /** Like values(...), but computes first and second derivatives of the shell function values, too. See the other hessian_values(...) members for more information. */ int hessian_shell_values(const SCVector3& r, int sh, ValueData *, double *h_values, double*g_values=0,double* basis_values=0) const; /// Returns true if this and the argument are equivalent. int equiv(const Ref &b); /// Print a brief description of the basis set. void print_brief(std::ostream& =ExEnv::out0()) const; /// Print a detailed description of the basis set. void print(std::ostream& =ExEnv::out0()) const; }; } #endif // Local Variables: // mode: c++ // c-file-style: "CLJ" // End: