[0b990d] | 1 | //
|
---|
| 2 | // gaussbas.h
|
---|
| 3 | //
|
---|
| 4 | // Copyright (C) 1996 Limit Point Systems, Inc.
|
---|
| 5 | //
|
---|
| 6 | // Author: Curtis Janssen <cljanss@limitpt.com>
|
---|
| 7 | // Maintainer: LPS
|
---|
| 8 | //
|
---|
| 9 | // This file is part of the SC Toolkit.
|
---|
| 10 | //
|
---|
| 11 | // The SC Toolkit is free software; you can redistribute it and/or modify
|
---|
| 12 | // it under the terms of the GNU Library General Public License as published by
|
---|
| 13 | // the Free Software Foundation; either version 2, or (at your option)
|
---|
| 14 | // any later version.
|
---|
| 15 | //
|
---|
| 16 | // The SC Toolkit is distributed in the hope that it will be useful,
|
---|
| 17 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 18 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 19 | // GNU Library General Public License for more details.
|
---|
| 20 | //
|
---|
| 21 | // You should have received a copy of the GNU Library General Public License
|
---|
| 22 | // along with the SC Toolkit; see the file COPYING.LIB. If not, write to
|
---|
| 23 | // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
|
---|
| 24 | //
|
---|
| 25 | // The U.S. Government is granted a limited license as per AL 91-7.
|
---|
| 26 | //
|
---|
| 27 |
|
---|
| 28 | #ifndef _chemistry_qc_basis_gaussbas_h
|
---|
| 29 | #define _chemistry_qc_basis_gaussbas_h
|
---|
| 30 |
|
---|
| 31 | #ifdef __GNUC__
|
---|
| 32 | #pragma interface
|
---|
| 33 | #endif
|
---|
| 34 |
|
---|
| 35 | #include <vector>
|
---|
| 36 | #include <iostream>
|
---|
| 37 |
|
---|
| 38 | #include <util/state/state.h>
|
---|
| 39 | #include <util/keyval/keyval.h>
|
---|
| 40 | #include <math/scmat/matrix.h>
|
---|
| 41 | #include <math/scmat/vector3.h>
|
---|
| 42 | #include <chemistry/molecule/molecule.h>
|
---|
| 43 |
|
---|
| 44 | namespace sc {
|
---|
| 45 |
|
---|
| 46 | class GaussianShell;
|
---|
| 47 | class BasisFileSet;
|
---|
| 48 | class Integral;
|
---|
| 49 |
|
---|
| 50 | class CartesianIter;
|
---|
| 51 | class SphericalTransformIter;
|
---|
| 52 |
|
---|
| 53 | /** The GaussianBasisSet class is used describe a basis set composed of
|
---|
| 54 | atomic gaussian orbitals. Inputs for common basis sets are included in the
|
---|
| 55 | MPQC distribution. They have been obtained from the EMSL Basis Set
|
---|
| 56 | Database and translated into the MPQC format. The citation for this
|
---|
| 57 | database is below. The technical citation for each basis set is listed in
|
---|
| 58 | the individual basis set data files, in MPQC's <tt>lib/basis</tt>
|
---|
| 59 | directory.
|
---|
| 60 |
|
---|
| 61 | Following is a table with available basis sets listing the supported
|
---|
| 62 | elements for each basis and the number of basis functions for H, \f$n_0\f$,
|
---|
| 63 | first row, \f$n_1\f$, and second row, \f$n_2\f$, atoms. Basis sets with
|
---|
| 64 | non-alpha-numerical characters in their name must be given in quotes.
|
---|
| 65 |
|
---|
| 66 | <table>
|
---|
| 67 | <tr><td>Basis Set<td>Elements<td>\f$n_0\f$<td>\f$n_1\f$<td>\f$n_2\f$
|
---|
| 68 | <tr><td><tt>STO-2G</tt><td>H-Ca<td>1<td>5<td>9
|
---|
| 69 | <tr><td><tt>STO-3G</tt><td>H-Kr<td>1<td>5<td>9
|
---|
| 70 | <tr><td><tt>STO-3G*</tt><td>H-Ar<td>1<td>5<td>14
|
---|
| 71 | <tr><td><tt>STO-6G</tt><td>H-Kr<td>1<td>5<td>9
|
---|
| 72 | <tr><td><tt>MINI (Huzinaga)</tt><td>H-Ca<td>1<td>5<td>9
|
---|
| 73 | <tr><td><tt>MINI (Scaled)</tt><td>H-Ca<td>1<td>5<td>9
|
---|
| 74 | <tr><td><tt>MIDI (Huzinaga)</tt><td>H-Na, Al-K<td>2<td>9<td>13
|
---|
| 75 | <tr><td><tt>DZ (Dunning)</tt><td>H, Li, B-Ne, Al-Cl<td>2<td>10<td>18
|
---|
| 76 | <tr><td><tt>DZP (Dunning)</tt><td>H, Li, B-Ne, Al-Cl<td>5<td>16<td>24
|
---|
| 77 | <tr><td><tt>DZP + Diffuse (Dunning)</tt><td>H, B-Ne<td>6<td>19<td>
|
---|
| 78 | <tr><td><tt>3-21G</tt><td>H-Kr<td>2<td>9<td>13
|
---|
| 79 | <tr><td><tt>3-21G*</tt><td>H-Ar<td>2<td>9<td>19
|
---|
| 80 | <tr><td><tt>3-21++G</tt><td>H-Ar<td>3<td>13<td>17
|
---|
| 81 | <tr><td><tt>3-21++G*</tt><td>H-Ar<td>3<td>13<td>23
|
---|
| 82 | <tr><td><tt>4-31G</tt><td>H-Ne, P-Cl<td>2<td>9<td>13
|
---|
| 83 | <tr><td><tt>6-31G</tt><td>H-Zn<td>2<td>9<td>13
|
---|
| 84 | <tr><td><tt>6-31G*</tt><td>H-Zn<td>2<td>15<td>19
|
---|
| 85 | <tr><td><tt>6-31G**</tt><td>H-Zn<td>5<td>15<td>19
|
---|
| 86 | <tr><td><tt>6-31+G*</tt><td>H-Ar<td>2<td>19<td>23
|
---|
| 87 | <tr><td><tt>6-31++G</tt><td>H-Ca<td>3<td>13<td>17
|
---|
| 88 | <tr><td><tt>6-31++G*</tt><td>H-Ar<td>3<td>19<td>23
|
---|
| 89 | <tr><td><tt>6-31++G**</tt><td>H-Ar<td>6<td>19<td>23
|
---|
| 90 | <tr><td><tt>6-311G</tt><td>H-Ca, Ga-Kr<td>3<td>13<td>21
|
---|
| 91 | <tr><td><tt>6-311G*</tt><td>H-Ca, Ga-Kr<td>3<td>18<td>26
|
---|
| 92 | <tr><td><tt>6-311G**</tt><td>H-Ca, Ga-Kr<td>6<td>18<td>26
|
---|
| 93 | <tr><td><tt>6-311G(2df,2pd)</tt><td>H-Ne, K, Ca<td>14<td>30<td>
|
---|
| 94 | <tr><td><tt>6-311++G**</tt><td>H-Ne<td>7<td>22<td>
|
---|
| 95 | <tr><td><tt>6-311++G(2d,2p)</tt><td>H-Ca<td>10<td>27<td>35
|
---|
| 96 | <tr><td><tt>6-311++G(3df,3pd)</tt><td>H-Ar<td>18<td>39<td>47
|
---|
| 97 | <tr><td><tt>cc-pVDZ</tt><td>H-Ar, Ca, Ga-Kr<td>5<td>14<td>18
|
---|
| 98 | <tr><td><tt>cc-pVTZ</tt><td>H-Ar, Ca, Ga-Kr<td>14<td>30<td>34
|
---|
| 99 | <tr><td><tt>cc-pVQZ</tt><td>H-Ar, Ca, Ga-Kr<td>30<td>55<td>59
|
---|
| 100 | <tr><td><tt>cc-pV5Z</tt><td>H-Ar, Ca, Ga-Kr<td>55<td>91<td>95
|
---|
| 101 | <tr><td><tt>cc-pV6Z</tt><td>H, He, B-Ne, Al-Ar<td>91<td>140<td>144
|
---|
| 102 | <tr><td><tt>aug-cc-pVDZ</tt><td>H, He, B-Ne, Al-Ar, Ga-Kr<td>9<td>23<td>27
|
---|
| 103 | <tr><td><tt>aug-cc-pVTZ</tt><td>H, He, B-Ne, Al-Ar, Ga-Kr<td>23<td>46<td>50
|
---|
| 104 | <tr><td><tt>aug-cc-pVQZ</tt><td>H, He, B-Ne, Al-Ar, Ga-Kr<td>46<td>80<td>84
|
---|
| 105 | <tr><td><tt>aug-cc-pV5Z</tt><td>H, He, B-Ne, Al-Ar, Ga-Kr<td>80<td>127<td>131
|
---|
| 106 | <tr><td><tt>aug-cc-pV6Z</tt><td>H, He, B-Ne, Al-Ar<td>127<td>189<td>193
|
---|
| 107 | <tr><td><tt>cc-pCVDZ</tt><td>Li, B-Ar<td><td>18<td>27
|
---|
| 108 | <tr><td><tt>cc-pCVTZ</tt><td>Li, B-Ar<td><td>43<td>59
|
---|
| 109 | <tr><td><tt>cc-pCVQZ</tt><td>Li, B-Ar<td><td>84<td>109
|
---|
| 110 | <tr><td><tt>cc-pCV5Z</tt><td>B-Ne<td><td>145<td>
|
---|
| 111 | <tr><td><tt>aug-cc-pCVDZ</tt><td>B-F, Al-Ar<td><td>27<td>36
|
---|
| 112 | <tr><td><tt>aug-cc-pCVTZ</tt><td>B-Ne, Al-Ar<td><td>59<td>75
|
---|
| 113 | <tr><td><tt>aug-cc-pCVQZ</tt><td>B-Ne, Al-Ar<td><td>109<td>134
|
---|
| 114 | <tr><td><tt>aug-cc-pCV5Z</tt><td>B-F<td><td>181<td>
|
---|
| 115 | <tr><td><tt>NASA Ames ANO</tt><td>H, B-Ne, Al, P, Ti, Fe, Ni<td>30<td>55<td>59
|
---|
| 116 | <tr><td><tt>pc-0</tt><td>H, C-F, Si-Cl<td>2<td>9<td>13
|
---|
| 117 | <tr><td><tt>pc-1</tt><td>H, C-F, Si-Cl<td>5<td>14<td>18
|
---|
| 118 | <tr><td><tt>pc-2</tt><td>H, C-F, Si-Cl<td>14<td>30<td>34
|
---|
| 119 | <tr><td><tt>pc-3</tt><td>H, C-F, Si-Cl<td>34<td>64<td>64
|
---|
| 120 | <tr><td><tt>pc-4</tt><td>H, C-F, Si-Cl<td>63<td>109<td>105
|
---|
| 121 | <tr><td><tt>pc-0-aug</tt><td>H, C-F, Si-Cl<td>3<td>13<td>17
|
---|
| 122 | <tr><td><tt>pc-1-aug</tt><td>H, C-F, Si-Cl<td>9<td>23<td>27
|
---|
| 123 | <tr><td><tt>pc-2-aug</tt><td>H, C-F, Si-Cl<td>23<td>46<td>50
|
---|
| 124 | <tr><td><tt>pc-3-aug</tt><td>H, C-F, Si-Cl<td>50<td>89<td>89
|
---|
| 125 | <tr><td><tt>pc-4-aug</tt><td>H, C-F, Si-Cl<td>88<td>145<td>141
|
---|
| 126 | </table>
|
---|
| 127 |
|
---|
| 128 | All basis sets except for the pc-<i>n</i> and pc-<i>n</i>-aug basis sets
|
---|
| 129 | were obtained from the Extensible Computational Chemistry
|
---|
| 130 | Environment Basis Set Database, Version 12/03/03, as developed and
|
---|
| 131 | distributed by the Molecular Science Computing Facility, Environmental and
|
---|
| 132 | Molecular Sciences Laboratory which is part of the Pacific Northwest
|
---|
| 133 | Laboratory, P.O. Box 999, Richland, Washington 99352, USA, and funded by
|
---|
| 134 | the U.S. Department of Energy. The Pacific Northwest Laboratory is a
|
---|
| 135 | multi-program laboratory operated by Battelle Memorial Institute for the
|
---|
| 136 | U.S. Department of Energy under contract DE-AC06-76RLO 1830. Contact David
|
---|
| 137 | Feller or Karen Schuchardt for further information.
|
---|
| 138 |
|
---|
| 139 | The pc-<i>n</i> and pc-<i>n</i>-aug basis sets are the polarization
|
---|
| 140 | consistent basis sets of Frank Jensen. See J. Chem. Phys. 115 (2001) 9113;
|
---|
| 141 | J. Chem. Phys. 116 (2002) 7372; J. Chem. Phys. 117 (2002) 9234; and
|
---|
| 142 | J. Chem. Phys. 121 (2004) 3463.
|
---|
| 143 |
|
---|
| 144 | */
|
---|
| 145 | class GaussianBasisSet: public SavableState
|
---|
| 146 | {
|
---|
| 147 | private:
|
---|
| 148 | // nonnull if keyword "name" was provided
|
---|
| 149 | char* name_;
|
---|
| 150 | // same as name_ if name_!=0, else something else
|
---|
| 151 | char* label_;
|
---|
| 152 | GaussianShell** shell_;
|
---|
| 153 | std::vector<int> shell_to_function_;
|
---|
| 154 | std::vector<int> function_to_shell_;
|
---|
| 155 |
|
---|
| 156 | Ref<Molecule> molecule_;
|
---|
| 157 |
|
---|
| 158 | Ref<SCMatrixKit> matrixkit_;
|
---|
| 159 | Ref<SCMatrixKit> so_matrixkit_;
|
---|
| 160 | RefSCDimension basisdim_;
|
---|
| 161 |
|
---|
| 162 | int ncenter_;
|
---|
| 163 |
|
---|
| 164 | std::vector<int> shell_to_center_;
|
---|
| 165 | std::vector<int> shell_to_primitive_;
|
---|
| 166 | std::vector<int> center_to_shell_;
|
---|
| 167 | std::vector<int> center_to_nshell_;
|
---|
| 168 | std::vector<int> center_to_nbasis_;
|
---|
| 169 |
|
---|
| 170 | int nshell_;
|
---|
| 171 | int nbasis_;
|
---|
| 172 | int nprim_;
|
---|
| 173 | bool has_pure_;
|
---|
| 174 |
|
---|
| 175 | GaussianBasisSet(const char* name, const char* label, const Ref<Molecule>& molecule,
|
---|
| 176 | const Ref<SCMatrixKit>& matrixkit,
|
---|
| 177 | const RefSCDimension& basisdim,
|
---|
| 178 | const int ncenter, const int nshell,
|
---|
| 179 | GaussianShell** shell,
|
---|
| 180 | const std::vector<int>& center_to_nshell);
|
---|
| 181 |
|
---|
| 182 | // Counts shells in this basis for this chemical element
|
---|
| 183 | int count_shells_(Ref<KeyVal>& keyval, const char* elemname, const char* sbasisname, BasisFileSet& bases,
|
---|
| 184 | int havepure, int pure, bool missing_ok);
|
---|
| 185 | // Constructs this basis
|
---|
| 186 | void get_shells_(int& ishell, Ref<KeyVal>& keyval, const char* elemname, const char* sbasisname, BasisFileSet& bases,
|
---|
| 187 | int havepure, int pure, bool missing_ok);
|
---|
| 188 | // Counts shells in an even-tempered primitive basis
|
---|
| 189 | int count_even_temp_shells_(Ref<KeyVal>& keyval, const char* elemname, const char* sbasisname,
|
---|
| 190 | int havepure, int pure);
|
---|
| 191 | // Constructs an even-tempered primitive basis
|
---|
| 192 | void get_even_temp_shells_(int& ishell, Ref<KeyVal>& keyval, const char* elemname, const char* sbasisname,
|
---|
| 193 | int havepure, int pure);
|
---|
| 194 | // Constructs basis set specified as an array of shells
|
---|
| 195 | void recursively_get_shell(int&,Ref<KeyVal>&,
|
---|
| 196 | const char*,const char*,BasisFileSet&,
|
---|
| 197 | int,int,int,bool missing_ok);
|
---|
| 198 |
|
---|
| 199 | void init(Ref<Molecule>&,Ref<KeyVal>&,
|
---|
| 200 | BasisFileSet&,
|
---|
| 201 | int have_userkeyval,
|
---|
| 202 | int pure);
|
---|
| 203 | void init2(int skip_ghosts=0,bool include_q=0);
|
---|
| 204 |
|
---|
| 205 | protected:
|
---|
| 206 | GaussianBasisSet(const GaussianBasisSet&);
|
---|
| 207 | virtual void set_matrixkit(const Ref<SCMatrixKit>&);
|
---|
| 208 |
|
---|
| 209 | public:
|
---|
| 210 | /** This holds scratch data needed to compute basis function values. */
|
---|
| 211 | class ValueData {
|
---|
| 212 | protected:
|
---|
| 213 | CartesianIter **civec_;
|
---|
| 214 | SphericalTransformIter **sivec_;
|
---|
| 215 | int maxam_;
|
---|
| 216 | public:
|
---|
| 217 | ValueData(const Ref<GaussianBasisSet> &, const Ref<Integral> &);
|
---|
| 218 | ~ValueData();
|
---|
| 219 | CartesianIter **civec() { return civec_; }
|
---|
| 220 | SphericalTransformIter **sivec() { return sivec_; }
|
---|
| 221 | };
|
---|
| 222 |
|
---|
| 223 | /// This can be given to a CTOR to construct a unit basis function.
|
---|
| 224 | enum UnitType {Unit};
|
---|
| 225 |
|
---|
| 226 | /** The KeyVal constructor.
|
---|
| 227 |
|
---|
| 228 | <dl>
|
---|
| 229 |
|
---|
| 230 | <dt><tt>molecule</tt><dd> The gives a Molecule object. The is no
|
---|
| 231 | default.
|
---|
| 232 |
|
---|
| 233 | <dt><tt>puream</tt><dd> If this boolean parameter is true then 5D,
|
---|
| 234 | 7F, etc. will be used. Otherwise all cartesian functions will be
|
---|
| 235 | used. The default depends on the particular basis set.
|
---|
| 236 |
|
---|
| 237 | <dt><tt>name</tt><dd> This is a string giving the name of the basis
|
---|
| 238 | set. The above table of basis sets gives some of the recognized
|
---|
| 239 | basis set names. It may be necessary to put the name in double
|
---|
| 240 | quotes. There is no default.
|
---|
| 241 |
|
---|
| 242 | <dt><tt>basis</tt><dd> This is a vector of basis set names that can
|
---|
| 243 | give a different basis set to each atom in the molecule. If the
|
---|
| 244 | element vector is given, then it gives different basis sets to
|
---|
| 245 | different elements. The default is to give every atom the basis
|
---|
| 246 | set specified in name.
|
---|
| 247 |
|
---|
| 248 | <dt><tt>element</tt><dd> This is a vector of elements. If it is
|
---|
| 249 | given then it must have the same number of entries as the basis
|
---|
| 250 | vector.
|
---|
| 251 |
|
---|
| 252 | <dt><tt>basisdir</tt><dd> A string giving a directory where basis
|
---|
| 253 | set data files are to be sought. See the text below for a complete
|
---|
| 254 | description of what directories are consulted.
|
---|
| 255 |
|
---|
| 256 | <dt><tt>basisfiles</tt><dd> Each keyword in this vector of files is
|
---|
| 257 | appended to the directory specified with basisdir and basis set
|
---|
| 258 | data is read from them.
|
---|
| 259 |
|
---|
| 260 | <dt><tt>matrixkit</tt><dd> Specifies a SCMatrixKit object. It is
|
---|
| 261 | usually not necessary to give this keyword, as the default action
|
---|
| 262 | should get the correct SCMatrixKit.
|
---|
| 263 |
|
---|
| 264 | </dl>
|
---|
| 265 |
|
---|
| 266 | Several files in various directories are checked for basis set
|
---|
| 267 | data. First, basis sets can be given by the user in the basis
|
---|
| 268 | section at the top level of the main input file. Next, if a path
|
---|
| 269 | is given with the basisdir keyword, then all of the files given
|
---|
| 270 | with the basisfiles keyword are read in after appending their names
|
---|
| 271 | to the value of basisdir. Basis sets can be given in these files
|
---|
| 272 | in the basis section at the top level as well. If the named basis
|
---|
| 273 | set still cannot be found, then GaussianBasisSet will try convert
|
---|
| 274 | the basis set name to a file name and check first in the directory
|
---|
| 275 | given by basisdir. Next it checks for the environment variable
|
---|
| 276 | SCLIBDIR. If it is set it will look for the basis file in
|
---|
| 277 | $SCLIBDIR/basis. Otherwise it will look in the source code
|
---|
| 278 | distribution in the directory SC/lib/basis. If the executable has
|
---|
| 279 | changed machines or the source code has be moved, then it may be
|
---|
| 280 | necessary to copy the library files to your machine and set the
|
---|
| 281 | SCLIBDIR environmental variable.
|
---|
| 282 |
|
---|
| 283 | The basis set itself is also given in the ParsedKeyVal format. There are two
|
---|
| 284 | recognized formats for basis sets:
|
---|
| 285 | <dl>
|
---|
| 286 |
|
---|
| 287 | <dt>array of shells<dd> One must specify the keyword :basis: followed by the
|
---|
| 288 | lowercase atomic name followed by : followed by the basis set name
|
---|
| 289 | (which may need to be placed inside double quotes). The value for the keyword
|
---|
| 290 | is an array of shells. Each shell
|
---|
| 291 | reads the following keywords:
|
---|
| 292 |
|
---|
| 293 | <dl>
|
---|
| 294 |
|
---|
| 295 | <dt><tt>type</tt><dd> This is a vector that describes each
|
---|
| 296 | component of this shell. For each element the following two
|
---|
| 297 | keywords are read:
|
---|
| 298 |
|
---|
| 299 | <dl>
|
---|
| 300 |
|
---|
| 301 | <dt><tt>am</tt><dd> The angular momentum of the component. This
|
---|
| 302 | can be given as the letter designation, s, p, d, etc. There is
|
---|
| 303 | no default.
|
---|
| 304 |
|
---|
| 305 | <dt><tt>puream</tt><dd> If this boolean parameter is true then
|
---|
| 306 | 5D, 7F, etc. shells are used. The default is false. This
|
---|
| 307 | parameter can be overridden in the GaussianBasisSet
|
---|
| 308 | specification.
|
---|
| 309 |
|
---|
| 310 | </dl>
|
---|
| 311 |
|
---|
| 312 | <dt><tt>exp</tt><dd> This is a vector giving the exponents of the
|
---|
| 313 | primitive Gaussian functions.
|
---|
| 314 |
|
---|
| 315 | <dt><tt>coef</tt><dd> This is a matrix giving the coeffients of the
|
---|
| 316 | primitive Gaussian functions. The first index gives the component
|
---|
| 317 | number of the shell and the second gives the primitive number.
|
---|
| 318 |
|
---|
| 319 | </dl>
|
---|
| 320 |
|
---|
| 321 | <dt><dd>An example might be easier to understand. This is a basis set
|
---|
| 322 | specificition for STO-2G carbon:
|
---|
| 323 |
|
---|
| 324 | <pre>
|
---|
| 325 | basis: (
|
---|
| 326 | carbon: "STO-2G": [
|
---|
| 327 | (type: [(am = s)]
|
---|
| 328 | { exp coef:0 } = {
|
---|
| 329 | 27.38503303 0.43012850
|
---|
| 330 | 4.87452205 0.67891353
|
---|
| 331 | })
|
---|
| 332 | (type: [(am = p) (am = s)]
|
---|
| 333 | { exp coef:1 coef:0 } = {
|
---|
| 334 | 1.13674819 0.04947177 0.51154071
|
---|
| 335 | 0.28830936 0.96378241 0.61281990
|
---|
| 336 | })
|
---|
| 337 | ]
|
---|
| 338 | )
|
---|
| 339 | </pre>
|
---|
| 340 |
|
---|
| 341 | <dt>basis set of even-tempered primitive Gaussians<dd>
|
---|
| 342 | Such basis set format is given as a group of keywords. The name of the group is :basis: followed by the
|
---|
| 343 | lowercase atomic name followed by : followed by the basis set name
|
---|
| 344 | (which may need to be placed inside double quotes).
|
---|
| 345 | The group of keywords must contain vectors <tt>am</tt> and <tt>nprim</tt>,
|
---|
| 346 | which specify the angular momentum and the number of shells in each
|
---|
| 347 | block of even-tempered primitives. In addition, one must provide any
|
---|
| 348 | two of the following vectors:
|
---|
| 349 |
|
---|
| 350 | <dl>
|
---|
| 351 | <dt><tt>first_exp</tt><dd> The exponent of the "tightest" primitive Gaussian in the block.
|
---|
| 352 |
|
---|
| 353 | <dt><tt>last_exp</tt><dd> The exponent of the most "diffuse" primitive Gaussian in the block.
|
---|
| 354 |
|
---|
| 355 | <dt><tt>exp_ratio</tt><dd> The ratio of exponents of consecutive primitive Gaussians
|
---|
| 356 | in the block.
|
---|
| 357 |
|
---|
| 358 | </dl>
|
---|
| 359 |
|
---|
| 360 | <dt><dd>Note that the dimensions of each vector must be the same.
|
---|
| 361 |
|
---|
| 362 | Here's an example of a basis set composed of 2 blocks of even-tempered s-functions
|
---|
| 363 | and 1 block of even-tempered p-functions.
|
---|
| 364 |
|
---|
| 365 | <pre>
|
---|
| 366 | basis: (
|
---|
| 367 | neon: "20s5s13p":(
|
---|
| 368 |
|
---|
| 369 | am = [ 0 0 1 ]
|
---|
| 370 | nprim = [ 20 5 13 ]
|
---|
| 371 | first_exp = [ 1000.0 0.1 70.0 ]
|
---|
| 372 | last_exp = [ 1.0 0.01 0.1 ]
|
---|
| 373 |
|
---|
| 374 | )
|
---|
| 375 | )
|
---|
| 376 | </pre>
|
---|
| 377 |
|
---|
| 378 | </dl>
|
---|
| 379 |
|
---|
| 380 | */
|
---|
| 381 | GaussianBasisSet(const Ref<KeyVal>&);
|
---|
| 382 | /** This can be given GaussianBasisSet::Unit to construct a basis
|
---|
| 383 | set with a single basis function that is one everywhere. This
|
---|
| 384 | can be used with integral evaluators to compute certain classes
|
---|
| 385 | of integrals, with limitations. */
|
---|
| 386 | GaussianBasisSet(UnitType);
|
---|
| 387 | GaussianBasisSet(StateIn&);
|
---|
| 388 | virtual ~GaussianBasisSet();
|
---|
| 389 |
|
---|
| 390 | /** Returns a GaussianBasisSet object that consists of the basis
|
---|
| 391 | functions for each atom in this followed by the basis functions in
|
---|
| 392 | B for the corresponding atom. The Molecule object for the two
|
---|
| 393 | basis sets must be identical. */
|
---|
| 394 | Ref<GaussianBasisSet> operator+(const Ref<GaussianBasisSet>& B);
|
---|
| 395 |
|
---|
| 396 | void save_data_state(StateOut&);
|
---|
| 397 |
|
---|
| 398 | /// Return the name of the basis set (is nonnull only if keyword "name" was provided)
|
---|
| 399 | const char* name() const { return name_; }
|
---|
| 400 | /** Return the label of the basis set. label() return the same string as name() if
|
---|
| 401 | keyword "name" was provided, otherwise it is a unique descriptive string which
|
---|
| 402 | can be arbitrarily long. */
|
---|
| 403 | const char* label() const { if (name()) { return name(); } else { return label_; } }
|
---|
| 404 |
|
---|
| 405 | /// Return the Molecule object.
|
---|
| 406 | Ref<Molecule> molecule() const { return molecule_; }
|
---|
| 407 | /// Returns the SCMatrixKit that is to be used for AO bases.
|
---|
| 408 | Ref<SCMatrixKit> matrixkit() { return matrixkit_; }
|
---|
| 409 | /// Returns the SCMatrixKit that is to be used for SO bases.
|
---|
| 410 | Ref<SCMatrixKit> so_matrixkit() { return so_matrixkit_; }
|
---|
| 411 | /// Returns the SCDimension object for the dimension.
|
---|
| 412 | RefSCDimension basisdim() { return basisdim_; }
|
---|
| 413 |
|
---|
| 414 | /// Return the number of centers.
|
---|
| 415 | int ncenter() const;
|
---|
| 416 | /// Return the number of shells.
|
---|
| 417 | int nshell() const { return nshell_; }
|
---|
| 418 | /// Return the number of shells on the given center.
|
---|
| 419 | int nshell_on_center(int icenter) const;
|
---|
| 420 | /** Return an overall shell number, given a center and the shell number
|
---|
| 421 | on that center. */
|
---|
| 422 | int shell_on_center(int icenter, int shell) const;
|
---|
| 423 | /// Return the center on which the given shell is located.
|
---|
| 424 | int shell_to_center(int ishell) const { return shell_to_center_[ishell]; }
|
---|
| 425 | /// Return the overall index of the first primitive from the given shell
|
---|
| 426 | int shell_to_primitive(int ishell) const {return shell_to_primitive_[ishell]; }
|
---|
| 427 | /// Return the number of basis functions.
|
---|
| 428 | int nbasis() const { return nbasis_; }
|
---|
| 429 | /// Return the number of basis functions on the given center.
|
---|
| 430 | int nbasis_on_center(int icenter) const;
|
---|
| 431 | /// Return the number of primitive Gaussians.
|
---|
| 432 | int nprimitive() const { return nprim_; }
|
---|
| 433 | /// Return true if basis contains solid harmonics Gaussians
|
---|
| 434 | int has_pure() const { return has_pure_; }
|
---|
| 435 |
|
---|
| 436 | /// Return the maximum number of functions that any shell has.
|
---|
| 437 | int max_nfunction_in_shell() const;
|
---|
| 438 | /** Return the maximum number of Cartesian functions that any shell has.
|
---|
| 439 | The optional argument is an angular momentum increment. */
|
---|
| 440 | int max_ncartesian_in_shell(int aminc=0) const;
|
---|
| 441 | /// Return the maximum number of primitive Gaussian that any shell has.
|
---|
| 442 | int max_nprimitive_in_shell() const;
|
---|
| 443 | /// Return the highest angular momentum in any shell.
|
---|
| 444 | int max_angular_momentum() const;
|
---|
| 445 | /// Return the maximum number of Gaussians in a contraction in any shell.
|
---|
| 446 | int max_ncontraction() const;
|
---|
| 447 | /** Return the maximum angular momentum found in the given contraction
|
---|
| 448 | number for any shell. */
|
---|
| 449 | int max_am_for_contraction(int con) const;
|
---|
| 450 | /// Return the maximum number of Cartesian functions in any shell.
|
---|
| 451 | int max_cartesian() const;
|
---|
| 452 |
|
---|
| 453 | /// Return the number of the first function in the given shell.
|
---|
| 454 | int shell_to_function(int i) const { return shell_to_function_[i]; }
|
---|
| 455 | /// Return the shell to which the given function belongs.
|
---|
| 456 | int function_to_shell(int i) const;
|
---|
| 457 |
|
---|
| 458 | /// Return a reference to GaussianShell number i.
|
---|
| 459 | const GaussianShell& operator()(int i) const { return *shell_[i]; }
|
---|
| 460 | /// Return a reference to GaussianShell number i.
|
---|
| 461 | GaussianShell& operator()(int i) { return *shell_[i]; }
|
---|
| 462 | /// Return a reference to GaussianShell number i.
|
---|
| 463 | const GaussianShell& operator[](int i) const { return *shell_[i]; }
|
---|
| 464 | /// Return a reference to GaussianShell number i.
|
---|
| 465 | GaussianShell& operator[](int i) { return *shell_[i]; }
|
---|
| 466 | /// Return a reference to GaussianShell number i.
|
---|
| 467 | const GaussianShell& shell(int i) const { return *shell_[i]; }
|
---|
| 468 | /// Return a reference to GaussianShell number i.
|
---|
| 469 | GaussianShell& shell(int i) { return *shell_[i]; }
|
---|
| 470 |
|
---|
| 471 | /// Return a reference to GaussianShell number ishell on center icenter.
|
---|
| 472 | const GaussianShell& operator()(int icenter,int ishell) const;
|
---|
| 473 | /// Return a reference to GaussianShell number ishell on center icenter.
|
---|
| 474 | GaussianShell& operator()(int icenter,int ishell);
|
---|
| 475 | /// Return a reference to GaussianShell number j on center i.
|
---|
| 476 | const GaussianShell& shell(int i,int j) const { return operator()(i,j); }
|
---|
| 477 | /// Return a reference to GaussianShell number j on center i.
|
---|
| 478 | GaussianShell& shell(int i,int j) { return operator()(i,j); }
|
---|
| 479 |
|
---|
| 480 | /** The location of center icenter. The xyz argument is 0 for x, 1 for
|
---|
| 481 | y, and 2 for z. */
|
---|
| 482 | double r(int icenter,int xyz) const;
|
---|
| 483 |
|
---|
| 484 | /** Compute the values for this basis set at position r. The
|
---|
| 485 | basis_values argument must be vector of length nbasis. */
|
---|
| 486 | int values(const SCVector3& r, ValueData *, double* basis_values) const;
|
---|
| 487 | /** Like values(...), but computes gradients of the basis function
|
---|
| 488 | values, too. The g_values argument must be vector of length
|
---|
| 489 | 3*nbasis. The data will be written in the order bf1_x, bf1_y,
|
---|
| 490 | bf1_z, ... */
|
---|
| 491 | int grad_values(const SCVector3& r, ValueData *,
|
---|
| 492 | double*g_values,double* basis_values=0) const;
|
---|
| 493 | /** Like values(...), but computes first and second derivatives of the
|
---|
| 494 | basis function values, too. h_values must be vector of length
|
---|
| 495 | 6*nbasis. The data will be written in the order bf1_xx, bf1_yx,
|
---|
| 496 | bf1_yy, bf1_zx, bf1_zy, bf1_zz, ... */
|
---|
| 497 | int hessian_values(const SCVector3& r, ValueData *, double *h_values,
|
---|
| 498 | double*g_values=0,double* basis_values=0) const;
|
---|
| 499 | /** Compute the values for the given shell functions at position r.
|
---|
| 500 | See the other values(...) members for more information. */
|
---|
| 501 | int shell_values(const SCVector3& r, int sh,
|
---|
| 502 | ValueData *, double* basis_values) const;
|
---|
| 503 | /** Like values(...), but computes gradients of the shell function
|
---|
| 504 | values, too. See the other grad_values(...)
|
---|
| 505 | members for more information. */
|
---|
| 506 | int grad_shell_values(const SCVector3& r, int sh,
|
---|
| 507 | ValueData *,
|
---|
| 508 | double*g_values, double* basis_values=0) const;
|
---|
| 509 | /** Like values(...), but computes first and second derivatives of the
|
---|
| 510 | shell function values, too. See the other hessian_values(...)
|
---|
| 511 | members for more information. */
|
---|
| 512 | int hessian_shell_values(const SCVector3& r, int sh,
|
---|
| 513 | ValueData *, double *h_values,
|
---|
| 514 | double*g_values=0,double* basis_values=0) const;
|
---|
| 515 |
|
---|
| 516 | /// Returns true if this and the argument are equivalent.
|
---|
| 517 | int equiv(const Ref<GaussianBasisSet> &b);
|
---|
| 518 |
|
---|
| 519 | /// Print a brief description of the basis set.
|
---|
| 520 | void print_brief(std::ostream& =ExEnv::out0()) const;
|
---|
| 521 | /// Print a detailed description of the basis set.
|
---|
| 522 | void print(std::ostream& =ExEnv::out0()) const;
|
---|
| 523 | };
|
---|
| 524 |
|
---|
| 525 | }
|
---|
| 526 |
|
---|
| 527 | #endif
|
---|
| 528 |
|
---|
| 529 | // Local Variables:
|
---|
| 530 | // mode: c++
|
---|
| 531 | // c-file-style: "CLJ"
|
---|
| 532 | // End:
|
---|