[0b990d] | 1 | //
|
---|
| 2 | // gaussbas.cc
|
---|
| 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 | #ifdef __GNUC__
|
---|
| 29 | #pragma implementation
|
---|
| 30 | #endif
|
---|
| 31 |
|
---|
| 32 | #include <stdio.h>
|
---|
| 33 | #include <stdexcept>
|
---|
| 34 |
|
---|
| 35 | #include <scconfig.h>
|
---|
| 36 | #ifdef HAVE_SSTREAM
|
---|
| 37 | # include <sstream>
|
---|
| 38 | #else
|
---|
| 39 | # include <strstream.h>
|
---|
| 40 | #endif
|
---|
| 41 |
|
---|
| 42 | #include <util/keyval/keyval.h>
|
---|
| 43 | #include <util/misc/formio.h>
|
---|
| 44 | #include <util/misc/newstring.h>
|
---|
| 45 | #include <util/state/stateio.h>
|
---|
| 46 | #include <math/scmat/blocked.h>
|
---|
| 47 | #include <chemistry/molecule/molecule.h>
|
---|
| 48 | #include <chemistry/qc/basis/gaussshell.h>
|
---|
| 49 | #include <chemistry/qc/basis/gaussbas.h>
|
---|
| 50 | #include <chemistry/qc/basis/files.h>
|
---|
| 51 | #include <chemistry/qc/basis/cartiter.h>
|
---|
| 52 | #include <chemistry/qc/basis/transform.h>
|
---|
| 53 | #include <chemistry/qc/basis/integral.h>
|
---|
| 54 |
|
---|
| 55 | using namespace std;
|
---|
| 56 | using namespace sc;
|
---|
| 57 |
|
---|
| 58 | static ClassDesc GaussianBasisSet_cd(
|
---|
| 59 | typeid(GaussianBasisSet),"GaussianBasisSet",3,"public SavableState",
|
---|
| 60 | 0, create<GaussianBasisSet>, create<GaussianBasisSet>);
|
---|
| 61 |
|
---|
| 62 | static bool
|
---|
| 63 | skip_atom(bool skip_ghosts, bool include_q,
|
---|
| 64 | const Ref<Molecule> &mol, int iatom)
|
---|
| 65 | {
|
---|
| 66 | if (skip_ghosts && mol->charge(iatom) == 0.0) return true;
|
---|
| 67 | // charges do not have basis functions
|
---|
| 68 | if (!include_q && mol->atom_symbol(iatom) == "Q") return true;
|
---|
| 69 | return false;
|
---|
| 70 | }
|
---|
| 71 |
|
---|
| 72 | GaussianBasisSet::GaussianBasisSet(const Ref<KeyVal>&topkeyval)
|
---|
| 73 | {
|
---|
| 74 | molecule_ << topkeyval->describedclassvalue("molecule");
|
---|
| 75 | if (molecule_.null()) {
|
---|
| 76 | ExEnv::err0() << indent << "GaussianBasisSet: no \"molecule\"\n";
|
---|
| 77 | abort();
|
---|
| 78 | }
|
---|
| 79 |
|
---|
| 80 | // see if the user requests pure am or cartesian functions
|
---|
| 81 | int pure;
|
---|
| 82 | pure = topkeyval->booleanvalue("puream");
|
---|
| 83 | if (topkeyval->error() != KeyVal::OK) pure = -1;
|
---|
| 84 |
|
---|
| 85 | // construct a keyval that contains the basis library
|
---|
| 86 | Ref<KeyVal> keyval;
|
---|
| 87 |
|
---|
| 88 | if (topkeyval->exists("basisfiles")) {
|
---|
| 89 | Ref<MessageGrp> grp = MessageGrp::get_default_messagegrp();
|
---|
| 90 | Ref<ParsedKeyVal> parsedkv = new ParsedKeyVal();
|
---|
| 91 | char *in_char_array;
|
---|
| 92 | if (grp->me() == 0) {
|
---|
| 93 | #ifdef HAVE_SSTREAM
|
---|
| 94 | ostringstream ostrs;
|
---|
| 95 | #else
|
---|
| 96 | ostrstream ostrs;
|
---|
| 97 | #endif
|
---|
| 98 | // Look at the basisdir and basisfiles variables to find out what
|
---|
| 99 | // basis set files are to be read in. The files are read on node
|
---|
| 100 | // 0 only.
|
---|
| 101 | ParsedKeyVal::cat_files("basis",topkeyval,ostrs);
|
---|
| 102 | #ifdef HAVE_SSTREAM
|
---|
| 103 | int n = 1 + strlen(ostrs.str().c_str());
|
---|
| 104 | in_char_array = strcpy(new char[n],ostrs.str().c_str());
|
---|
| 105 | #else
|
---|
| 106 | ostrs << ends;
|
---|
| 107 | in_char_array = ostrs.str();
|
---|
| 108 | int n = ostrs.pcount();
|
---|
| 109 | #endif
|
---|
| 110 | grp->bcast(n);
|
---|
| 111 | grp->bcast(in_char_array, n);
|
---|
| 112 | }
|
---|
| 113 | else {
|
---|
| 114 | int n;
|
---|
| 115 | grp->bcast(n);
|
---|
| 116 | in_char_array = new char[n];
|
---|
| 117 | grp->bcast(in_char_array, n);
|
---|
| 118 | }
|
---|
| 119 | parsedkv->parse_string(in_char_array);
|
---|
| 120 | delete[] in_char_array;
|
---|
| 121 | Ref<KeyVal> libkeyval = parsedkv.pointer();
|
---|
| 122 | keyval = new AggregateKeyVal(topkeyval,libkeyval);
|
---|
| 123 | }
|
---|
| 124 | else {
|
---|
| 125 | keyval = topkeyval;
|
---|
| 126 | }
|
---|
| 127 |
|
---|
| 128 | // if there isn't a matrixkit in the input, init2() will assign the
|
---|
| 129 | // default matrixkit
|
---|
| 130 | matrixkit_ << keyval->describedclassvalue("matrixkit");
|
---|
| 131 |
|
---|
| 132 | // Bases keeps track of what basis set data bases have already
|
---|
| 133 | // been read in. It also handles the conversion of basis
|
---|
| 134 | // names to file names.
|
---|
| 135 | BasisFileSet bases(keyval);
|
---|
| 136 | init(molecule_,keyval,bases,1,pure);
|
---|
| 137 | }
|
---|
| 138 |
|
---|
| 139 | GaussianBasisSet::GaussianBasisSet(UnitType)
|
---|
| 140 | {
|
---|
| 141 | molecule_ = new Molecule;
|
---|
| 142 | molecule_->add_atom(molecule()->atominfo()->string_to_Z("Q"),
|
---|
| 143 | 0.0, 0.0, 0.0, // xyz
|
---|
| 144 | "dummy", // label
|
---|
| 145 | 0.0, // no mass
|
---|
| 146 | 1, 0.0 // no charge
|
---|
| 147 | );
|
---|
| 148 | name_ = new_string("Unit");
|
---|
| 149 | label_ = new_string(name_);
|
---|
| 150 | shell_ = new GaussianShell*[1];
|
---|
| 151 |
|
---|
| 152 | double *exp = new double[1];
|
---|
| 153 | int *am = new int[1];
|
---|
| 154 | int *pure = new int[1];
|
---|
| 155 | double **c = new double*[1];
|
---|
| 156 | *c = new double[1];
|
---|
| 157 | exp[0] = 0.0;
|
---|
| 158 | am[0] = 0;
|
---|
| 159 | pure[0] = 0;
|
---|
| 160 | c[0][0] = 1.0;
|
---|
| 161 | shell_[0] = new GaussianShell(1,1,exp,am,pure,c,
|
---|
| 162 | GaussianShell::Unnormalized,
|
---|
| 163 | false);
|
---|
| 164 |
|
---|
| 165 | ncenter_ = 1;
|
---|
| 166 | nshell_ = 1;
|
---|
| 167 | center_to_nshell_.push_back(1);
|
---|
| 168 | init2(0,1);
|
---|
| 169 | }
|
---|
| 170 |
|
---|
| 171 | GaussianBasisSet::GaussianBasisSet(const GaussianBasisSet& gbs) :
|
---|
| 172 | molecule_(gbs.molecule_),
|
---|
| 173 | matrixkit_(gbs.matrixkit_),
|
---|
| 174 | basisdim_(gbs.basisdim_),
|
---|
| 175 | ncenter_(gbs.ncenter_),
|
---|
| 176 | nshell_(gbs.nshell_)
|
---|
| 177 | {
|
---|
| 178 | int i,j;
|
---|
| 179 |
|
---|
| 180 | name_ = new_string(gbs.name_);
|
---|
| 181 | label_ = new_string(gbs.label_);
|
---|
| 182 |
|
---|
| 183 | center_to_nshell_.resize(ncenter_);
|
---|
| 184 | for (i=0; i < ncenter_; i++) {
|
---|
| 185 | center_to_nshell_[i] = gbs.center_to_nshell_[i];
|
---|
| 186 | }
|
---|
| 187 |
|
---|
| 188 | shell_ = new GaussianShell*[nshell_];
|
---|
| 189 | for (i=0; i<nshell_; i++) {
|
---|
| 190 | const GaussianShell& gsi = gbs(i);
|
---|
| 191 |
|
---|
| 192 | int nc=gsi.ncontraction();
|
---|
| 193 | int np=gsi.nprimitive();
|
---|
| 194 |
|
---|
| 195 | int *ams = new int[nc];
|
---|
| 196 | int *pure = new int[nc];
|
---|
| 197 | double *exps = new double[np];
|
---|
| 198 | double **coefs = new double*[nc];
|
---|
| 199 |
|
---|
| 200 | for (j=0; j < nc; j++) {
|
---|
| 201 | ams[j] = gsi.am(j);
|
---|
| 202 | pure[j] = gsi.is_pure(j);
|
---|
| 203 | coefs[j] = new double[np];
|
---|
| 204 | for (int k=0; k < np; k++)
|
---|
| 205 | coefs[j][k] = gsi.coefficient_unnorm(j,k);
|
---|
| 206 | }
|
---|
| 207 |
|
---|
| 208 | for (j=0; j < np; j++)
|
---|
| 209 | exps[j] = gsi.exponent(j);
|
---|
| 210 |
|
---|
| 211 | shell_[i] = new GaussianShell(nc, np, exps, ams, pure, coefs,
|
---|
| 212 | GaussianShell::Unnormalized);
|
---|
| 213 | }
|
---|
| 214 |
|
---|
| 215 | init2();
|
---|
| 216 | }
|
---|
| 217 |
|
---|
| 218 | GaussianBasisSet::GaussianBasisSet(const char* name,
|
---|
| 219 | const char* label,
|
---|
| 220 | const Ref<Molecule>& molecule,
|
---|
| 221 | const Ref<SCMatrixKit>& matrixkit,
|
---|
| 222 | const RefSCDimension& basisdim,
|
---|
| 223 | const int ncenter,
|
---|
| 224 | const int nshell,
|
---|
| 225 | GaussianShell** shell,
|
---|
| 226 | const std::vector<int>& center_to_nshell) :
|
---|
| 227 | molecule_(molecule),
|
---|
| 228 | matrixkit_(matrixkit),
|
---|
| 229 | basisdim_(basisdim),
|
---|
| 230 | ncenter_(ncenter),
|
---|
| 231 | nshell_(nshell),
|
---|
| 232 | shell_(shell),
|
---|
| 233 | center_to_nshell_(center_to_nshell)
|
---|
| 234 | {
|
---|
| 235 | name_ = new_string(name);
|
---|
| 236 | label_ = new_string(label);
|
---|
| 237 |
|
---|
| 238 | init2();
|
---|
| 239 | }
|
---|
| 240 |
|
---|
| 241 | Ref<GaussianBasisSet>
|
---|
| 242 | GaussianBasisSet::operator+(const Ref<GaussianBasisSet>& B)
|
---|
| 243 | {
|
---|
| 244 | GaussianBasisSet* b = B.pointer();
|
---|
| 245 | if (molecule_.pointer() != b->molecule_.pointer())
|
---|
| 246 | throw std::runtime_error("GaussianBasisSet::concatenate -- cannot concatenate basis sets, molecules are different");
|
---|
| 247 |
|
---|
| 248 | Ref<Molecule> molecule = molecule_;
|
---|
| 249 | Ref<SCMatrixKit> matrixkit = matrixkit_;
|
---|
| 250 | const int ncenter = ncenter_;
|
---|
| 251 | const int nshell = nshell_ + b->nshell_;
|
---|
| 252 | std::vector<int> center_to_nshell(ncenter);
|
---|
| 253 |
|
---|
| 254 | GaussianShell** shell = new GaussianShell*[nshell];
|
---|
| 255 | int* func_per_shell = new int[nshell];
|
---|
| 256 |
|
---|
| 257 | for(int c=0; c<ncenter; c++) {
|
---|
| 258 |
|
---|
| 259 | int ns1 = center_to_nshell_[c];
|
---|
| 260 | int ns2 = b->center_to_nshell_[c];
|
---|
| 261 | int ns = ns1+ns2;
|
---|
| 262 | int s1off = center_to_shell_[c];
|
---|
| 263 | int s2off = b->center_to_shell_[c];
|
---|
| 264 | int soff = s1off + s2off;
|
---|
| 265 | center_to_nshell[c] = ns;
|
---|
| 266 |
|
---|
| 267 | for (int i=0; i<ns; i++) {
|
---|
| 268 | const GaussianShell* gsi;
|
---|
| 269 | if (i < ns1)
|
---|
| 270 | gsi = shell_[s1off + i];
|
---|
| 271 | else
|
---|
| 272 | gsi = b->shell_[s2off + i - ns1];
|
---|
| 273 |
|
---|
| 274 | int nc=gsi->ncontraction();
|
---|
| 275 | int np=gsi->nprimitive();
|
---|
| 276 | func_per_shell[soff + i] = gsi->nfunction();
|
---|
| 277 |
|
---|
| 278 | int *ams = new int[nc];
|
---|
| 279 | int *pure = new int[nc];
|
---|
| 280 | double *exps = new double[np];
|
---|
| 281 | double **coefs = new double*[nc];
|
---|
| 282 |
|
---|
| 283 | for (int j=0; j < nc; j++) {
|
---|
| 284 | ams[j] = gsi->am(j);
|
---|
| 285 | pure[j] = gsi->is_pure(j);
|
---|
| 286 | coefs[j] = new double[np];
|
---|
| 287 | for (int k=0; k < np; k++)
|
---|
| 288 | coefs[j][k] = gsi->coefficient_unnorm(j,k);
|
---|
| 289 | }
|
---|
| 290 |
|
---|
| 291 | for (int j=0; j < np; j++)
|
---|
| 292 | exps[j] = gsi->exponent(j);
|
---|
| 293 |
|
---|
| 294 | shell[soff + i] = new GaussianShell(nc, np, exps, ams, pure, coefs,
|
---|
| 295 | GaussianShell::Unnormalized);
|
---|
| 296 | }
|
---|
| 297 | }
|
---|
| 298 |
|
---|
| 299 | int nbas = nbasis() + b->nbasis();
|
---|
| 300 | RefSCDimension basisdim
|
---|
| 301 | = new SCDimension(nbas, nshell, func_per_shell, "basis set dimension");
|
---|
| 302 |
|
---|
| 303 | const char* A_name = name();
|
---|
| 304 | const char* B_name = B->name();
|
---|
| 305 | const char* AplusB_name = 0;
|
---|
| 306 | if (!A_name && !B_name) {
|
---|
| 307 | ostringstream oss;
|
---|
| 308 | oss << "[" << A_name << "]+[" << B_name << "]";
|
---|
| 309 | std::string tmpname = oss.str();
|
---|
| 310 | AplusB_name = strcpy(new char[tmpname.size()+1],tmpname.c_str());
|
---|
| 311 | }
|
---|
| 312 | const char* AplusB_label = 0;
|
---|
| 313 | if (AplusB_name) {
|
---|
| 314 | AplusB_label = AplusB_name;
|
---|
| 315 | }
|
---|
| 316 | else {
|
---|
| 317 | ostringstream oss;
|
---|
| 318 | const char* A_label = label();
|
---|
| 319 | const char* B_label = B->label();
|
---|
| 320 | oss << "[" << A_label << "]+[" << B_label << "]";
|
---|
| 321 | std::string tmpname = oss.str();
|
---|
| 322 | AplusB_label = strcpy(new char[tmpname.size()+1],tmpname.c_str());
|
---|
| 323 | }
|
---|
| 324 | Ref<GaussianBasisSet> AplusB
|
---|
| 325 | = new GaussianBasisSet(AplusB_name, AplusB_label, molecule,
|
---|
| 326 | matrixkit, basisdim, ncenter,
|
---|
| 327 | nshell, shell, center_to_nshell);
|
---|
| 328 |
|
---|
| 329 | delete[] func_per_shell;
|
---|
| 330 | delete[] AplusB_name;
|
---|
| 331 | if (AplusB_name != AplusB_label)
|
---|
| 332 | delete[] AplusB_label;
|
---|
| 333 |
|
---|
| 334 | return AplusB;
|
---|
| 335 | }
|
---|
| 336 |
|
---|
| 337 | GaussianBasisSet::GaussianBasisSet(StateIn&s):
|
---|
| 338 | SavableState(s)
|
---|
| 339 | {
|
---|
| 340 | matrixkit_ = SCMatrixKit::default_matrixkit();
|
---|
| 341 |
|
---|
| 342 | if (s.version(::class_desc<GaussianBasisSet>()) < 3) {
|
---|
| 343 | // read the duplicate length saved in older versions
|
---|
| 344 | int junk;
|
---|
| 345 | s.get(junk);
|
---|
| 346 | }
|
---|
| 347 | s.get(center_to_nshell_);
|
---|
| 348 |
|
---|
| 349 | molecule_ << SavableState::restore_state(s);
|
---|
| 350 | basisdim_ << SavableState::restore_state(s);
|
---|
| 351 |
|
---|
| 352 |
|
---|
| 353 | ncenter_ = center_to_nshell_.size();
|
---|
| 354 | s.getstring(name_);
|
---|
| 355 | s.getstring(label_);
|
---|
| 356 |
|
---|
| 357 | nshell_ = 0;
|
---|
| 358 | int i;
|
---|
| 359 | for (i=0; i<ncenter_; i++) {
|
---|
| 360 | nshell_ += center_to_nshell_[i];
|
---|
| 361 | }
|
---|
| 362 |
|
---|
| 363 | shell_ = new GaussianShell*[nshell_];
|
---|
| 364 | for (i=0; i<nshell_; i++) {
|
---|
| 365 | shell_[i] = new GaussianShell(s);
|
---|
| 366 | }
|
---|
| 367 |
|
---|
| 368 | init2();
|
---|
| 369 | }
|
---|
| 370 |
|
---|
| 371 | void
|
---|
| 372 | GaussianBasisSet::save_data_state(StateOut&s)
|
---|
| 373 | {
|
---|
| 374 | s.put(center_to_nshell_);
|
---|
| 375 |
|
---|
| 376 | SavableState::save_state(molecule_.pointer(),s);
|
---|
| 377 | SavableState::save_state(basisdim_.pointer(),s);
|
---|
| 378 |
|
---|
| 379 | s.putstring(name_);
|
---|
| 380 | s.putstring(label_);
|
---|
| 381 | for (int i=0; i<nshell_; i++) {
|
---|
| 382 | shell_[i]->save_object_state(s);
|
---|
| 383 | }
|
---|
| 384 | }
|
---|
| 385 |
|
---|
| 386 | void
|
---|
| 387 | GaussianBasisSet::init(Ref<Molecule>&molecule,
|
---|
| 388 | Ref<KeyVal>&keyval,
|
---|
| 389 | BasisFileSet& bases,
|
---|
| 390 | int have_userkeyval,
|
---|
| 391 | int pur)
|
---|
| 392 | {
|
---|
| 393 | int pure, havepure;
|
---|
| 394 | pure = pur;
|
---|
| 395 | if (pur == -1) {
|
---|
| 396 | havepure = 0;
|
---|
| 397 | }
|
---|
| 398 | else {
|
---|
| 399 | havepure = 1;
|
---|
| 400 | }
|
---|
| 401 |
|
---|
| 402 | int skip_ghosts = keyval->booleanvalue("skip_ghosts");
|
---|
| 403 | bool missing_ok = keyval->booleanvalue("missing_ok");
|
---|
| 404 | bool include_q = keyval->booleanvalue("include_q");
|
---|
| 405 |
|
---|
| 406 | name_ = keyval->pcharvalue("name");
|
---|
| 407 | int have_custom, nelement;
|
---|
| 408 |
|
---|
| 409 | if (keyval->exists("basis")) {
|
---|
| 410 | have_custom = 1;
|
---|
| 411 | nelement = keyval->count("element");
|
---|
| 412 | }
|
---|
| 413 | else {
|
---|
| 414 | have_custom = 0;
|
---|
| 415 | nelement = 0;
|
---|
| 416 | if (!name_) {
|
---|
| 417 | ExEnv::err0() << indent
|
---|
| 418 | << "GaussianBasisSet: No name given for basis set\n";
|
---|
| 419 | abort();
|
---|
| 420 | }
|
---|
| 421 | }
|
---|
| 422 |
|
---|
| 423 | // Construct label_
|
---|
| 424 | if (name_)
|
---|
| 425 | label_ = new_string(name_);
|
---|
| 426 | else {
|
---|
| 427 | if (have_custom) {
|
---|
| 428 | ostringstream oss;
|
---|
| 429 | Ref<AtomInfo> atominfo = molecule->atominfo();
|
---|
| 430 | // If element is given then construct label_ using element symbol and basis name
|
---|
| 431 | // combinations, e.g. "{ [Fe S1] [Ni S2] [C aug-cc-pVDZ] }"
|
---|
| 432 | if (nelement) {
|
---|
| 433 | oss << "{ ";
|
---|
| 434 | for(int e=0; e<nelement; e++) {
|
---|
| 435 | char* tmpelementname = keyval->pcharvalue("element", e);
|
---|
| 436 | int Z = atominfo->string_to_Z(tmpelementname);
|
---|
| 437 | std::string elemsymbol = atominfo->symbol(Z);
|
---|
| 438 | char* basisname = keyval->pcharvalue("basis", e);
|
---|
| 439 | oss << "[" << elemsymbol << " " << basisname << "] ";
|
---|
| 440 | }
|
---|
| 441 | oss << "}";
|
---|
| 442 | }
|
---|
| 443 | // If element is not given then construct label_ using basis names for each atom
|
---|
| 444 | // e.g. "[ aug-cc-pVDZ cc-pVDZ cc-pVDZ ]"
|
---|
| 445 | else {
|
---|
| 446 | int natom = molecule->natom();
|
---|
| 447 | oss << "[ ";
|
---|
| 448 | for(int a=0; a<natom; a++) {
|
---|
| 449 | char* basisname = keyval->pcharvalue("basis", a);
|
---|
| 450 | oss << basisname << " ";
|
---|
| 451 | }
|
---|
| 452 | oss << "]";
|
---|
| 453 | }
|
---|
| 454 | std::string label = oss.str();
|
---|
| 455 | label_ = new char[label.size() + 1];
|
---|
| 456 | strcpy(label_,label.c_str());
|
---|
| 457 | }
|
---|
| 458 | }
|
---|
| 459 |
|
---|
| 460 |
|
---|
| 461 | // construct prefixes for each atom: :basis:<atom>:<basisname>:<shell#>
|
---|
| 462 | // and read in the shell
|
---|
| 463 | nbasis_ = 0;
|
---|
| 464 | int ishell = 0;
|
---|
| 465 | ncenter_ = molecule->natom();
|
---|
| 466 | int iatom;
|
---|
| 467 | for (iatom=0; iatom < ncenter_; iatom++) {
|
---|
| 468 | if (skip_atom(skip_ghosts,include_q,molecule,iatom)) continue;
|
---|
| 469 | int Z = molecule->Z(iatom);
|
---|
| 470 | // see if there is a specific basisname for this atom
|
---|
| 471 | char* sbasisname = 0;
|
---|
| 472 | if (have_custom && !nelement) {
|
---|
| 473 | sbasisname = keyval->pcharvalue("basis",iatom);
|
---|
| 474 | }
|
---|
| 475 | else if (nelement) {
|
---|
| 476 | int i;
|
---|
| 477 | for (i=0; i<nelement; i++) {
|
---|
| 478 | char *tmpelementname = keyval->pcharvalue("element", i);
|
---|
| 479 | int tmpZ = molecule->atominfo()->string_to_Z(tmpelementname);
|
---|
| 480 | if (tmpZ == Z) {
|
---|
| 481 | sbasisname = keyval->pcharvalue("basis", i);
|
---|
| 482 | break;
|
---|
| 483 | }
|
---|
| 484 | delete[] tmpelementname;
|
---|
| 485 | }
|
---|
| 486 | }
|
---|
| 487 | if (!sbasisname) {
|
---|
| 488 | if (!name_) {
|
---|
| 489 | ExEnv::err0()
|
---|
| 490 | << indent << "GaussianBasisSet: no basis name for atom "
|
---|
| 491 | << iatom
|
---|
| 492 | << " (Z=" <<molecule->atominfo()->name(Z) << ")"
|
---|
| 493 | << std::endl;
|
---|
| 494 | abort();
|
---|
| 495 | }
|
---|
| 496 | sbasisname = strcpy(new char[strlen(name_)+1],name_);
|
---|
| 497 | }
|
---|
| 498 | std::string name(molecule->atominfo()->name(Z));
|
---|
| 499 | ishell += count_shells_(keyval, name.c_str(),
|
---|
| 500 | sbasisname, bases, havepure, pure, missing_ok);
|
---|
| 501 | delete[] sbasisname;
|
---|
| 502 | }
|
---|
| 503 | nshell_ = ishell;
|
---|
| 504 | shell_ = new GaussianShell*[nshell_];
|
---|
| 505 | ishell = 0;
|
---|
| 506 | center_to_nshell_.resize(ncenter_);
|
---|
| 507 | for (iatom=0; iatom<ncenter_; iatom++) {
|
---|
| 508 | if (skip_atom(skip_ghosts,include_q,molecule,iatom)) {
|
---|
| 509 | center_to_nshell_[iatom] = 0;
|
---|
| 510 | continue;
|
---|
| 511 | }
|
---|
| 512 | int Z = molecule->Z(iatom);
|
---|
| 513 | // see if there is a specific basisname for this atom
|
---|
| 514 | char* sbasisname = 0;
|
---|
| 515 | if (have_custom && !nelement) {
|
---|
| 516 | sbasisname = keyval->pcharvalue("basis",iatom);
|
---|
| 517 | }
|
---|
| 518 | else if (nelement) {
|
---|
| 519 | int i;
|
---|
| 520 | for (i=0; i<nelement; i++) {
|
---|
| 521 | char *tmpelementname = keyval->pcharvalue("element", i);
|
---|
| 522 | int tmpZ = molecule->atominfo()->string_to_Z(tmpelementname);
|
---|
| 523 | if (tmpZ == Z) {
|
---|
| 524 | sbasisname = keyval->pcharvalue("basis", i);
|
---|
| 525 | break;
|
---|
| 526 | }
|
---|
| 527 | delete[] tmpelementname;
|
---|
| 528 | }
|
---|
| 529 | }
|
---|
| 530 | if (!sbasisname) {
|
---|
| 531 | if (!name_) {
|
---|
| 532 | ExEnv::err0()
|
---|
| 533 | << indent << "GaussianBasisSet: no basis name for atom "
|
---|
| 534 | << iatom
|
---|
| 535 | << " (Z=" <<molecule->atominfo()->name(Z) << ")"
|
---|
| 536 | << std::endl;
|
---|
| 537 | abort();
|
---|
| 538 | }
|
---|
| 539 | sbasisname = strcpy(new char[strlen(name_)+1],name_);
|
---|
| 540 | }
|
---|
| 541 |
|
---|
| 542 | int ishell_old = ishell;
|
---|
| 543 | std::string name(molecule->atominfo()->name(Z));
|
---|
| 544 | get_shells_(ishell, keyval, name.c_str(),
|
---|
| 545 | sbasisname, bases, havepure, pure, missing_ok);
|
---|
| 546 |
|
---|
| 547 | center_to_nshell_[iatom] = ishell - ishell_old;
|
---|
| 548 |
|
---|
| 549 | delete[] sbasisname;
|
---|
| 550 | }
|
---|
| 551 |
|
---|
| 552 | // delete the name_ if the basis set is customized
|
---|
| 553 | if (have_custom) {
|
---|
| 554 | delete[] name_;
|
---|
| 555 | name_ = 0;
|
---|
| 556 | }
|
---|
| 557 |
|
---|
| 558 | // finish with the initialization steps that don't require any
|
---|
| 559 | // external information
|
---|
| 560 | init2(skip_ghosts,include_q);
|
---|
| 561 | }
|
---|
| 562 |
|
---|
| 563 | double
|
---|
| 564 | GaussianBasisSet::r(int icenter, int xyz) const
|
---|
| 565 | {
|
---|
| 566 | return molecule_->r(icenter,xyz);
|
---|
| 567 | }
|
---|
| 568 |
|
---|
| 569 | void
|
---|
| 570 | GaussianBasisSet::init2(int skip_ghosts,bool include_q)
|
---|
| 571 | {
|
---|
| 572 | // center_to_shell_ and shell_to_center_
|
---|
| 573 | shell_to_center_.resize(nshell_);
|
---|
| 574 | center_to_shell_.resize(ncenter_);
|
---|
| 575 | center_to_nbasis_.resize(ncenter_);
|
---|
| 576 | int ishell = 0;
|
---|
| 577 | for (int icenter=0; icenter<ncenter_; icenter++) {
|
---|
| 578 | if (skip_atom(skip_ghosts,include_q,molecule(),icenter)) {
|
---|
| 579 | center_to_shell_[icenter] = -1;
|
---|
| 580 | center_to_nbasis_[icenter] = 0;
|
---|
| 581 | continue;
|
---|
| 582 | }
|
---|
| 583 | int j;
|
---|
| 584 | center_to_shell_[icenter] = ishell;
|
---|
| 585 | center_to_nbasis_[icenter] = 0;
|
---|
| 586 | for (j = 0; j<center_to_nshell_[icenter]; j++) {
|
---|
| 587 | center_to_nbasis_[icenter] += shell_[ishell+j]->nfunction();
|
---|
| 588 | }
|
---|
| 589 | ishell += center_to_nshell_[icenter];
|
---|
| 590 | for (j = center_to_shell_[icenter]; j<ishell; j++) {
|
---|
| 591 | shell_to_center_[j] = icenter;
|
---|
| 592 | }
|
---|
| 593 | }
|
---|
| 594 |
|
---|
| 595 | // compute nbasis_ and shell_to_function_[]
|
---|
| 596 | shell_to_function_.resize(nshell_);
|
---|
| 597 | shell_to_primitive_.resize(nshell_);
|
---|
| 598 | nbasis_ = 0;
|
---|
| 599 | nprim_ = 0;
|
---|
| 600 | for (ishell=0; ishell<nshell_; ishell++) {
|
---|
| 601 | shell_to_function_[ishell] = nbasis_;
|
---|
| 602 | shell_to_primitive_[ishell] = nprim_;
|
---|
| 603 | nbasis_ += shell_[ishell]->nfunction();
|
---|
| 604 | nprim_ += shell_[ishell]->nprimitive();
|
---|
| 605 | }
|
---|
| 606 |
|
---|
| 607 | // would like to do this in function_to_shell(), but it is const
|
---|
| 608 | int n = nbasis();
|
---|
| 609 | int nsh = nshell();
|
---|
| 610 | function_to_shell_.resize(n);
|
---|
| 611 | int ifunc = 0;
|
---|
| 612 | for (int i=0; i<nsh; i++) {
|
---|
| 613 | int nfun = operator()(i).nfunction();
|
---|
| 614 | for (int j=0; j<nfun; j++) {
|
---|
| 615 | function_to_shell_[ifunc] = i;
|
---|
| 616 | ifunc++;
|
---|
| 617 | }
|
---|
| 618 | }
|
---|
| 619 |
|
---|
| 620 | // figure out if any shells are spherical harmonics
|
---|
| 621 | has_pure_ = false;
|
---|
| 622 | for(int i=0; i<nsh; i++)
|
---|
| 623 | has_pure_ = has_pure_ || shell_[i]->has_pure();
|
---|
| 624 |
|
---|
| 625 | if (matrixkit_.null())
|
---|
| 626 | matrixkit_ = SCMatrixKit::default_matrixkit();
|
---|
| 627 |
|
---|
| 628 | so_matrixkit_ = new BlockedSCMatrixKit(matrixkit_);
|
---|
| 629 |
|
---|
| 630 | if (basisdim_.null()) {
|
---|
| 631 | int nb = nshell();
|
---|
| 632 | int *bs = new int[nb];
|
---|
| 633 | for (int s=0; s < nb; s++)
|
---|
| 634 | bs[s] = shell(s).nfunction();
|
---|
| 635 | basisdim_ = new SCDimension(nbasis(), nb, bs, "basis set dimension");
|
---|
| 636 | delete[] bs;
|
---|
| 637 | }
|
---|
| 638 | }
|
---|
| 639 |
|
---|
| 640 | void
|
---|
| 641 | GaussianBasisSet::set_matrixkit(const Ref<SCMatrixKit>& mk)
|
---|
| 642 | {
|
---|
| 643 | matrixkit_ = mk;
|
---|
| 644 | so_matrixkit_ = new BlockedSCMatrixKit(matrixkit_);
|
---|
| 645 | }
|
---|
| 646 |
|
---|
| 647 |
|
---|
| 648 | int
|
---|
| 649 | GaussianBasisSet::count_shells_(Ref<KeyVal>& keyval, const char* element, const char* basisname, BasisFileSet& bases,
|
---|
| 650 | int havepure, int pure, bool missing_ok)
|
---|
| 651 | {
|
---|
| 652 | int nshell = 0;
|
---|
| 653 | char keyword[KeyVal::MaxKeywordLength];
|
---|
| 654 |
|
---|
| 655 | sprintf(keyword,":basis:%s:%s",element,basisname);
|
---|
| 656 | bool exists = keyval->exists(keyword);
|
---|
| 657 | if (!exists) {
|
---|
| 658 | keyval = bases.keyval(keyval, basisname);
|
---|
| 659 | exists = keyval->exists(keyword);
|
---|
| 660 | if (!exists) {
|
---|
| 661 | if (missing_ok) return 0;
|
---|
| 662 | ExEnv::err0() << indent
|
---|
| 663 | << scprintf("GaussianBasisSet::count_shells_ couldn't find \"%s\":\n", keyword);
|
---|
| 664 | keyval->errortrace(ExEnv::err0());
|
---|
| 665 | throw std::runtime_error("GaussianBasisSet::count_shells_ -- couldn't find the basis set");
|
---|
| 666 | }
|
---|
| 667 | }
|
---|
| 668 |
|
---|
| 669 | // Check if the basis set is an array of shells
|
---|
| 670 | keyval->count(keyword);
|
---|
| 671 | if (keyval->error() != KeyVal::OK) {
|
---|
| 672 | nshell = count_even_temp_shells_(keyval, element, basisname, havepure, pure);
|
---|
| 673 | }
|
---|
| 674 | else {
|
---|
| 675 | recursively_get_shell(nshell, keyval, element, basisname,
|
---|
| 676 | bases, havepure, pure, 0, false);
|
---|
| 677 | }
|
---|
| 678 |
|
---|
| 679 | return nshell;
|
---|
| 680 | }
|
---|
| 681 |
|
---|
| 682 | void
|
---|
| 683 | GaussianBasisSet::get_shells_(int& ishell, Ref<KeyVal>& keyval, const char* element, const char* basisname, BasisFileSet& bases,
|
---|
| 684 | int havepure, int pure, bool missing_ok)
|
---|
| 685 | {
|
---|
| 686 | char keyword[KeyVal::MaxKeywordLength];
|
---|
| 687 |
|
---|
| 688 | sprintf(keyword,":basis:%s:%s:am",
|
---|
| 689 | element,basisname);
|
---|
| 690 | if (keyval->exists(keyword)) {
|
---|
| 691 | get_even_temp_shells_(ishell, keyval, element, basisname, havepure, pure);
|
---|
| 692 | }
|
---|
| 693 | else {
|
---|
| 694 | recursively_get_shell(ishell, keyval, element, basisname,
|
---|
| 695 | bases, havepure, pure, 1, missing_ok);
|
---|
| 696 | }
|
---|
| 697 | }
|
---|
| 698 |
|
---|
| 699 |
|
---|
| 700 | int
|
---|
| 701 | GaussianBasisSet::count_even_temp_shells_(Ref<KeyVal>& keyval, const char* element, const char* basisname,
|
---|
| 702 | int havepure, int pure)
|
---|
| 703 | {
|
---|
| 704 | int nshell = 0;
|
---|
| 705 | char keyword[KeyVal::MaxKeywordLength];
|
---|
| 706 |
|
---|
| 707 | sprintf(keyword,":basis:%s:%s:am",element,basisname);
|
---|
| 708 | if (!keyval->exists(keyword)) {
|
---|
| 709 | sprintf(keyword,":basis:%s:%s",element,basisname);
|
---|
| 710 | ExEnv::err0() << indent
|
---|
| 711 | << scprintf("GaussianBasisSet::count_even_temp_shells_ -- couldn't read \"%s\":\n", keyword);
|
---|
| 712 | throw std::runtime_error("GaussianBasisSet::count_even_temp_shells_ -- basis set specification is invalid");
|
---|
| 713 | }
|
---|
| 714 |
|
---|
| 715 | // count the number of even-tempered primitive blocks
|
---|
| 716 | int nblocks = keyval->count(keyword) - 1;
|
---|
| 717 | if (keyval->error() != KeyVal::OK) {
|
---|
| 718 | ExEnv::err0() << indent
|
---|
| 719 | << scprintf("GaussianBasisSet::count_even_temp_shells_ -- couldn't read \"%s\":\n", keyword);
|
---|
| 720 | throw std::runtime_error("GaussianBasisSet::count_even_temp_shells_ -- failed to read am");
|
---|
| 721 | }
|
---|
| 722 | if (nblocks == -1)
|
---|
| 723 | return 0;
|
---|
| 724 |
|
---|
| 725 | sprintf(keyword,":basis:%s:%s:nprim", element, basisname);
|
---|
| 726 | int j = keyval->count(keyword) - 1;
|
---|
| 727 | if (nblocks != j) {
|
---|
| 728 | ExEnv::err0() << indent
|
---|
| 729 | << scprintf("GaussianBasisSet::count_even_temp_shells_ -- problem reading \"%s\":\n", keyword);
|
---|
| 730 | throw std::runtime_error("GaussianBasisSet::count_even_temp_shells_ -- am and nprim have different dimensions");
|
---|
| 731 | }
|
---|
| 732 |
|
---|
| 733 | for(int b=0; b<=nblocks; b++) {
|
---|
| 734 | sprintf(keyword,":basis:%s:%s:nprim:%d", element, basisname, b);
|
---|
| 735 | int nprim = keyval->intvalue(keyword);
|
---|
| 736 | if (nprim <= 0) {
|
---|
| 737 | ExEnv::err0() << indent
|
---|
| 738 | << scprintf("GaussianBasisSet::count_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 739 | throw std::runtime_error("GaussianBasisSet::count_shells_ -- the number of primitives has to be positive");
|
---|
| 740 | }
|
---|
| 741 | nshell += nprim;
|
---|
| 742 | }
|
---|
| 743 |
|
---|
| 744 | return nshell;
|
---|
| 745 | }
|
---|
| 746 |
|
---|
| 747 |
|
---|
| 748 | void
|
---|
| 749 | GaussianBasisSet::get_even_temp_shells_(int& ishell, Ref<KeyVal>& keyval, const char* element, const char* basisname,
|
---|
| 750 | int havepure, int pure)
|
---|
| 751 | {
|
---|
| 752 | char keyword[KeyVal::MaxKeywordLength];
|
---|
| 753 |
|
---|
| 754 | // count the number of even-tempered primitive blocks
|
---|
| 755 | sprintf(keyword,":basis:%s:%s:am",
|
---|
| 756 | element,basisname);
|
---|
| 757 | int nblocks = keyval->count(keyword) - 1;
|
---|
| 758 | if (keyval->error() != KeyVal::OK) {
|
---|
| 759 | ExEnv::err0() << indent
|
---|
| 760 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- couldn't read \"%s\":\n", keyword);
|
---|
| 761 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- failed to read am");
|
---|
| 762 | }
|
---|
| 763 | if (nblocks == -1)
|
---|
| 764 | return;
|
---|
| 765 |
|
---|
| 766 | sprintf(keyword,":basis:%s:%s:nprim", element, basisname);
|
---|
| 767 | int j = keyval->count(keyword) - 1;
|
---|
| 768 | if (nblocks != j) {
|
---|
| 769 | ExEnv::err0() << indent
|
---|
| 770 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem reading \"%s\":\n", keyword);
|
---|
| 771 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- am and nprim have different dimensions");
|
---|
| 772 | }
|
---|
| 773 |
|
---|
| 774 | sprintf(keyword,":basis:%s:%s:last_exp", element, basisname);
|
---|
| 775 | bool have_last_exp = keyval->exists(keyword);
|
---|
| 776 |
|
---|
| 777 | sprintf(keyword,":basis:%s:%s:first_exp", element, basisname);
|
---|
| 778 | bool have_first_exp = keyval->exists(keyword);
|
---|
| 779 |
|
---|
| 780 | sprintf(keyword,":basis:%s:%s:exp_ratio", element, basisname);
|
---|
| 781 | bool have_exp_ratio = keyval->exists(keyword);
|
---|
| 782 |
|
---|
| 783 | if ( !have_first_exp && !have_last_exp )
|
---|
| 784 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- neither last_exp nor first_exp has been specified");
|
---|
| 785 |
|
---|
| 786 | if ( have_first_exp && have_last_exp && have_exp_ratio)
|
---|
| 787 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- only two of (last_exp,first_exp,exp_ratio) can be specified");
|
---|
| 788 |
|
---|
| 789 | if ( !have_first_exp && !have_last_exp && have_exp_ratio)
|
---|
| 790 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- any two of (last_exp,first_exp,exp_ratio) must be specified");
|
---|
| 791 |
|
---|
| 792 | for(int b=0; b<=nblocks; b++) {
|
---|
| 793 |
|
---|
| 794 | sprintf(keyword,":basis:%s:%s:nprim:%d", element, basisname, b);
|
---|
| 795 | int nprim = keyval->intvalue(keyword);
|
---|
| 796 | if (nprim <= 0) {
|
---|
| 797 | ExEnv::err0() << indent
|
---|
| 798 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 799 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- the number of primitives has to be positive");
|
---|
| 800 | }
|
---|
| 801 |
|
---|
| 802 | sprintf(keyword,":basis:%s:%s:am:%d", element, basisname, b);
|
---|
| 803 | int l = keyval->intvalue(keyword);
|
---|
| 804 | if (l < 0) {
|
---|
| 805 | ExEnv::err0() << indent
|
---|
| 806 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 807 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- angular momentum has to be non-negative");
|
---|
| 808 | }
|
---|
| 809 |
|
---|
| 810 | double alpha0, alphaN, beta;
|
---|
| 811 | if (have_first_exp) {
|
---|
| 812 | sprintf(keyword,":basis:%s:%s:first_exp:%d", element, basisname, b);
|
---|
| 813 | alpha0 = keyval->doublevalue(keyword);
|
---|
| 814 | if (alpha0 <= 0.0) {
|
---|
| 815 | ExEnv::err0() << indent
|
---|
| 816 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 817 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- orbital exponents have to be positive");
|
---|
| 818 | }
|
---|
| 819 | }
|
---|
| 820 |
|
---|
| 821 | if (have_last_exp) {
|
---|
| 822 | sprintf(keyword,":basis:%s:%s:last_exp:%d", element, basisname, b);
|
---|
| 823 | alphaN = keyval->doublevalue(keyword);
|
---|
| 824 | if (alphaN <= 0.0) {
|
---|
| 825 | ExEnv::err0() << indent
|
---|
| 826 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 827 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- orbital exponents have to be positive");
|
---|
| 828 | }
|
---|
| 829 | }
|
---|
| 830 |
|
---|
| 831 | if (have_last_exp && have_first_exp) {
|
---|
| 832 | if (alphaN > alpha0) {
|
---|
| 833 | ExEnv::err0() << indent
|
---|
| 834 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 835 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- last_exps[i] must be smaller than first_exp[i]");
|
---|
| 836 | }
|
---|
| 837 | if (nprim > 1)
|
---|
| 838 | beta = pow(alpha0/alphaN,1.0/(nprim-1));
|
---|
| 839 | else
|
---|
| 840 | beta = 1.0;
|
---|
| 841 | }
|
---|
| 842 | else {
|
---|
| 843 | sprintf(keyword,":basis:%s:%s:exp_ratio:%d", element, basisname, b);
|
---|
| 844 | beta = keyval->doublevalue(keyword);
|
---|
| 845 | if (beta <= 1.0) {
|
---|
| 846 | ExEnv::err0() << indent
|
---|
| 847 | << scprintf("GaussianBasisSet::get_even_temp_shells_ -- problem with \"%s\":\n", keyword);
|
---|
| 848 | throw std::runtime_error("GaussianBasisSet::get_even_temp_shells_ -- exponent ratio has to be greater than 1.0");
|
---|
| 849 | }
|
---|
| 850 | if (have_last_exp)
|
---|
| 851 | alpha0 = alphaN * pow(beta,nprim-1);
|
---|
| 852 | }
|
---|
| 853 |
|
---|
| 854 | double alpha = alpha0;
|
---|
| 855 | for(int p=0; p<nprim; p++, alpha /= beta ) {
|
---|
| 856 | int* am = new int[1];
|
---|
| 857 | double* exps = new double[1];
|
---|
| 858 | double** coeffs = new double*[1];
|
---|
| 859 | coeffs[0] = new double[1];
|
---|
| 860 |
|
---|
| 861 | exps[0] = alpha;
|
---|
| 862 | am[0] = l;
|
---|
| 863 | coeffs[0][0] = 1.0;
|
---|
| 864 |
|
---|
| 865 | if (l <= 1)
|
---|
| 866 | shell_[ishell] = new GaussianShell(1,1,exps,am,GaussianShell::Cartesian,coeffs,GaussianShell::Normalized);
|
---|
| 867 | else if (havepure)
|
---|
| 868 | shell_[ishell] = new GaussianShell(1,1,exps,am,GaussianShell::Pure,coeffs,GaussianShell::Normalized);
|
---|
| 869 | else
|
---|
| 870 | shell_[ishell] = new GaussianShell(1,1,exps,am,GaussianShell::Cartesian,coeffs,GaussianShell::Normalized);
|
---|
| 871 | ishell++;
|
---|
| 872 |
|
---|
| 873 | // Recompute beta at each step to maintain accuracy
|
---|
| 874 | if (have_last_exp && have_first_exp) {
|
---|
| 875 | int nprim_left = nprim - p;
|
---|
| 876 | if (nprim_left > 1)
|
---|
| 877 | beta = pow(alpha/alphaN,1.0/(nprim_left-1));
|
---|
| 878 | else
|
---|
| 879 | beta = 1.0;
|
---|
| 880 | }
|
---|
| 881 | }
|
---|
| 882 | }
|
---|
| 883 | }
|
---|
| 884 |
|
---|
| 885 | void
|
---|
| 886 | GaussianBasisSet::
|
---|
| 887 | recursively_get_shell(int&ishell,Ref<KeyVal>&keyval,
|
---|
| 888 | const char*element,
|
---|
| 889 | const char*basisname,
|
---|
| 890 | BasisFileSet &bases,
|
---|
| 891 | int havepure,int pure,
|
---|
| 892 | int get, bool missing_ok)
|
---|
| 893 | {
|
---|
| 894 | char keyword[KeyVal::MaxKeywordLength],prefix[KeyVal::MaxKeywordLength];
|
---|
| 895 |
|
---|
| 896 | sprintf(keyword,":basis:%s:%s",
|
---|
| 897 | element,basisname);
|
---|
| 898 | int count = keyval->count(keyword);
|
---|
| 899 | if (keyval->error() != KeyVal::OK) {
|
---|
| 900 | keyval = bases.keyval(keyval, basisname);
|
---|
| 901 | }
|
---|
| 902 | count = keyval->count(keyword);
|
---|
| 903 | if (keyval->error() != KeyVal::OK) {
|
---|
| 904 | if (missing_ok) return;
|
---|
| 905 | ExEnv::err0() << indent
|
---|
| 906 | << scprintf("GaussianBasisSet:: couldn't find \"%s\":\n", keyword);
|
---|
| 907 | keyval->errortrace(ExEnv::err0());
|
---|
| 908 | throw std::runtime_error("GaussianBasisSet::recursively_get_shell -- couldn't find the basis set");
|
---|
| 909 | }
|
---|
| 910 | if (!count) return;
|
---|
| 911 | for (int j=0; j<count; j++) {
|
---|
| 912 | sprintf(prefix,":basis:%s:%s",
|
---|
| 913 | element,basisname);
|
---|
| 914 | Ref<KeyVal> prefixkeyval = new PrefixKeyVal(keyval,prefix,j);
|
---|
| 915 | if (prefixkeyval->exists("get")) {
|
---|
| 916 | char* newbasis = prefixkeyval->pcharvalue("get");
|
---|
| 917 | if (!newbasis) {
|
---|
| 918 | ExEnv::err0() << indent << "GaussianBasisSet: "
|
---|
| 919 | << scprintf("error processing get for \"%s\"\n", prefix);
|
---|
| 920 | keyval->errortrace(ExEnv::err0());
|
---|
| 921 | exit(1);
|
---|
| 922 | }
|
---|
| 923 | recursively_get_shell(ishell,keyval,element,newbasis,bases,
|
---|
| 924 | havepure,pure,get,missing_ok);
|
---|
| 925 | delete[] newbasis;
|
---|
| 926 | }
|
---|
| 927 | else {
|
---|
| 928 | if (get) {
|
---|
| 929 | if (havepure) shell_[ishell] = new GaussianShell(prefixkeyval,pure);
|
---|
| 930 | else shell_[ishell] = new GaussianShell(prefixkeyval);
|
---|
| 931 | }
|
---|
| 932 | ishell++;
|
---|
| 933 | }
|
---|
| 934 | }
|
---|
| 935 | }
|
---|
| 936 |
|
---|
| 937 | GaussianBasisSet::~GaussianBasisSet()
|
---|
| 938 | {
|
---|
| 939 | delete[] name_;
|
---|
| 940 | delete[] label_;
|
---|
| 941 |
|
---|
| 942 | int ii;
|
---|
| 943 | for (ii=0; ii<nshell_; ii++) {
|
---|
| 944 | delete shell_[ii];
|
---|
| 945 | }
|
---|
| 946 | delete[] shell_;
|
---|
| 947 | }
|
---|
| 948 |
|
---|
| 949 | int
|
---|
| 950 | GaussianBasisSet::max_nfunction_in_shell() const
|
---|
| 951 | {
|
---|
| 952 | int i;
|
---|
| 953 | int max = 0;
|
---|
| 954 | for (i=0; i<nshell_; i++) {
|
---|
| 955 | if (max < shell_[i]->nfunction()) max = shell_[i]->nfunction();
|
---|
| 956 | }
|
---|
| 957 | return max;
|
---|
| 958 | }
|
---|
| 959 |
|
---|
| 960 | int
|
---|
| 961 | GaussianBasisSet::max_ncontraction() const
|
---|
| 962 | {
|
---|
| 963 | int i;
|
---|
| 964 | int max = 0;
|
---|
| 965 | for (i=0; i<nshell_; i++) {
|
---|
| 966 | if (max < shell_[i]->ncontraction()) max = shell_[i]->ncontraction();
|
---|
| 967 | }
|
---|
| 968 | return max;
|
---|
| 969 | }
|
---|
| 970 |
|
---|
| 971 | int
|
---|
| 972 | GaussianBasisSet::max_angular_momentum() const
|
---|
| 973 | {
|
---|
| 974 | int i;
|
---|
| 975 | int max = 0;
|
---|
| 976 | for (i=0; i<nshell_; i++) {
|
---|
| 977 | int maxshi = shell_[i]->max_angular_momentum();
|
---|
| 978 | if (max < maxshi) max = maxshi;
|
---|
| 979 | }
|
---|
| 980 | return max;
|
---|
| 981 | }
|
---|
| 982 |
|
---|
| 983 | int
|
---|
| 984 | GaussianBasisSet::max_cartesian() const
|
---|
| 985 | {
|
---|
| 986 | int i;
|
---|
| 987 | int max = 0;
|
---|
| 988 | for (i=0; i<nshell_; i++) {
|
---|
| 989 | int maxshi = shell_[i]->max_cartesian();
|
---|
| 990 | if (max < maxshi) max = maxshi;
|
---|
| 991 | }
|
---|
| 992 | return max;
|
---|
| 993 | }
|
---|
| 994 |
|
---|
| 995 | int
|
---|
| 996 | GaussianBasisSet::max_ncartesian_in_shell(int aminc) const
|
---|
| 997 | {
|
---|
| 998 | int i;
|
---|
| 999 | int max = 0;
|
---|
| 1000 | for (i=0; i<nshell_; i++) {
|
---|
| 1001 | int maxshi = shell_[i]->ncartesian_with_aminc(aminc);
|
---|
| 1002 | if (max < maxshi) max = maxshi;
|
---|
| 1003 | }
|
---|
| 1004 | return max;
|
---|
| 1005 | }
|
---|
| 1006 |
|
---|
| 1007 | int
|
---|
| 1008 | GaussianBasisSet::max_nprimitive_in_shell() const
|
---|
| 1009 | {
|
---|
| 1010 | int i;
|
---|
| 1011 | int max = 0;
|
---|
| 1012 | for (i=0; i<nshell_; i++) {
|
---|
| 1013 | if (max < shell_[i]->nprimitive()) max = shell_[i]->nprimitive();
|
---|
| 1014 | }
|
---|
| 1015 | return max;
|
---|
| 1016 | }
|
---|
| 1017 |
|
---|
| 1018 | int
|
---|
| 1019 | GaussianBasisSet::max_am_for_contraction(int con) const
|
---|
| 1020 | {
|
---|
| 1021 | int i;
|
---|
| 1022 | int max = 0;
|
---|
| 1023 | for (i=0; i<nshell_; i++) {
|
---|
| 1024 | if (shell_[i]->ncontraction() <= con) continue;
|
---|
| 1025 | int maxshi = shell_[i]->am(con);
|
---|
| 1026 | if (max < maxshi) max = maxshi;
|
---|
| 1027 | }
|
---|
| 1028 | return max;
|
---|
| 1029 | }
|
---|
| 1030 |
|
---|
| 1031 | int
|
---|
| 1032 | GaussianBasisSet::function_to_shell(int func) const
|
---|
| 1033 | {
|
---|
| 1034 | return function_to_shell_[func];
|
---|
| 1035 | }
|
---|
| 1036 |
|
---|
| 1037 | int
|
---|
| 1038 | GaussianBasisSet::ncenter() const
|
---|
| 1039 | {
|
---|
| 1040 | return ncenter_;
|
---|
| 1041 | }
|
---|
| 1042 |
|
---|
| 1043 | int
|
---|
| 1044 | GaussianBasisSet::nshell_on_center(int icenter) const
|
---|
| 1045 | {
|
---|
| 1046 | return center_to_nshell_[icenter];
|
---|
| 1047 | }
|
---|
| 1048 |
|
---|
| 1049 | int
|
---|
| 1050 | GaussianBasisSet::nbasis_on_center(int icenter) const
|
---|
| 1051 | {
|
---|
| 1052 | return center_to_nbasis_[icenter];
|
---|
| 1053 | }
|
---|
| 1054 |
|
---|
| 1055 | int
|
---|
| 1056 | GaussianBasisSet::shell_on_center(int icenter, int ishell) const
|
---|
| 1057 | {
|
---|
| 1058 | return center_to_shell_[icenter] + ishell;
|
---|
| 1059 | }
|
---|
| 1060 |
|
---|
| 1061 | const GaussianShell&
|
---|
| 1062 | GaussianBasisSet::operator()(int icenter,int ishell) const
|
---|
| 1063 | {
|
---|
| 1064 | return *shell_[center_to_shell_[icenter] + ishell];
|
---|
| 1065 | }
|
---|
| 1066 |
|
---|
| 1067 | GaussianShell&
|
---|
| 1068 | GaussianBasisSet::operator()(int icenter,int ishell)
|
---|
| 1069 | {
|
---|
| 1070 | return *shell_[center_to_shell_[icenter] + ishell];
|
---|
| 1071 | }
|
---|
| 1072 |
|
---|
| 1073 | int
|
---|
| 1074 | GaussianBasisSet::equiv(const Ref<GaussianBasisSet> &b)
|
---|
| 1075 | {
|
---|
| 1076 | if (nshell() != b->nshell()) return 0;
|
---|
| 1077 | for (int i=0; i<nshell(); i++) {
|
---|
| 1078 | if (!shell_[i]->equiv(b->shell_[i])) return 0;
|
---|
| 1079 | }
|
---|
| 1080 | return 1;
|
---|
| 1081 | }
|
---|
| 1082 |
|
---|
| 1083 | void
|
---|
| 1084 | GaussianBasisSet::print_brief(ostream& os) const
|
---|
| 1085 | {
|
---|
| 1086 | os << indent
|
---|
| 1087 | << "GaussianBasisSet:" << endl << incindent
|
---|
| 1088 | << indent << "nbasis = " << nbasis_ << endl
|
---|
| 1089 | << indent << "nshell = " << nshell_ << endl
|
---|
| 1090 | << indent << "nprim = " << nprim_ << endl;
|
---|
| 1091 | if (name_) {
|
---|
| 1092 | os << indent
|
---|
| 1093 | << "name = \"" << name_ << "\"" << endl;
|
---|
| 1094 | }
|
---|
| 1095 | else {
|
---|
| 1096 | os << indent
|
---|
| 1097 | << "label = \"" << label_ << "\"" << endl;
|
---|
| 1098 | }
|
---|
| 1099 | os << decindent;
|
---|
| 1100 | }
|
---|
| 1101 |
|
---|
| 1102 | void
|
---|
| 1103 | GaussianBasisSet::print(ostream& os) const
|
---|
| 1104 | {
|
---|
| 1105 | print_brief(os);
|
---|
| 1106 | if (!SCFormIO::getverbose(os)) return;
|
---|
| 1107 |
|
---|
| 1108 | os << incindent;
|
---|
| 1109 |
|
---|
| 1110 | // Loop over centers
|
---|
| 1111 | int icenter = 0;
|
---|
| 1112 | int ioshell = 0;
|
---|
| 1113 | for (icenter=0; icenter < ncenter_; icenter++) {
|
---|
| 1114 | os << endl << indent
|
---|
| 1115 | << scprintf(
|
---|
| 1116 | "center %d: %12.8f %12.8f %12.8f, nshell = %d, shellnum = %d\n",
|
---|
| 1117 | icenter,
|
---|
| 1118 | r(icenter,0),
|
---|
| 1119 | r(icenter,1),
|
---|
| 1120 | r(icenter,2),
|
---|
| 1121 | center_to_nshell_[icenter],
|
---|
| 1122 | center_to_shell_[icenter]);
|
---|
| 1123 | for (int ishell=0; ishell < center_to_nshell_[icenter]; ishell++) {
|
---|
| 1124 | os << indent
|
---|
| 1125 | << scprintf("Shell %d: functionnum = %d, primnum = %d\n",
|
---|
| 1126 | ishell,shell_to_function_[ioshell],shell_to_primitive_[ioshell]);
|
---|
| 1127 | os << incindent;
|
---|
| 1128 | operator()(icenter,ishell).print(os);
|
---|
| 1129 | os << decindent;
|
---|
| 1130 | ioshell++;
|
---|
| 1131 | }
|
---|
| 1132 | }
|
---|
| 1133 |
|
---|
| 1134 | os << decindent;
|
---|
| 1135 | }
|
---|
| 1136 |
|
---|
| 1137 | /////////////////////////////////////////////////////////////////////////////
|
---|
| 1138 | // GaussianBasisSet::ValueData
|
---|
| 1139 |
|
---|
| 1140 | GaussianBasisSet::ValueData::ValueData(
|
---|
| 1141 | const Ref<GaussianBasisSet> &basis,
|
---|
| 1142 | const Ref<Integral> &integral)
|
---|
| 1143 | {
|
---|
| 1144 | maxam_ = basis->max_angular_momentum();
|
---|
| 1145 |
|
---|
| 1146 | civec_ = new CartesianIter *[maxam_+1];
|
---|
| 1147 | sivec_ = new SphericalTransformIter *[maxam_+1];
|
---|
| 1148 | for (int i=0; i<=maxam_; i++) {
|
---|
| 1149 | civec_[i] = integral->new_cartesian_iter(i);
|
---|
| 1150 | if (i>1) sivec_[i] = integral->new_spherical_transform_iter(i);
|
---|
| 1151 | else sivec_[i] = 0;
|
---|
| 1152 | }
|
---|
| 1153 | }
|
---|
| 1154 |
|
---|
| 1155 | GaussianBasisSet::ValueData::~ValueData()
|
---|
| 1156 | {
|
---|
| 1157 | for (int i=0; i<=maxam_; i++) {
|
---|
| 1158 | delete civec_[i];
|
---|
| 1159 | delete sivec_[i];
|
---|
| 1160 | }
|
---|
| 1161 | delete[] civec_;
|
---|
| 1162 | delete[] sivec_;
|
---|
| 1163 | }
|
---|
| 1164 |
|
---|
| 1165 | /////////////////////////////////////////////////////////////////////////////
|
---|
| 1166 |
|
---|
| 1167 | // Local Variables:
|
---|
| 1168 | // mode: c++
|
---|
| 1169 | // c-file-style: "CLJ"
|
---|
| 1170 | // End:
|
---|