[0b990d] | 1 | //
|
---|
| 2 | // molecule.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 <math.h>
|
---|
| 33 | #include <string.h>
|
---|
| 34 | #include <stdio.h>
|
---|
| 35 |
|
---|
| 36 | #include <util/class/scexception.h>
|
---|
| 37 | #include <util/misc/formio.h>
|
---|
| 38 | #include <util/state/stateio.h>
|
---|
| 39 | #include <chemistry/molecule/molecule.h>
|
---|
| 40 | #include <chemistry/molecule/formula.h>
|
---|
| 41 | #include <chemistry/molecule/localdef.h>
|
---|
| 42 | #include <math/scmat/cmatrix.h>
|
---|
| 43 |
|
---|
| 44 | using namespace std;
|
---|
| 45 | using namespace sc;
|
---|
| 46 |
|
---|
| 47 | //////////////////////////////////////////////////////////////////////////
|
---|
| 48 | // Molecule
|
---|
| 49 |
|
---|
| 50 | static ClassDesc Molecule_cd(
|
---|
| 51 | typeid(Molecule),"Molecule",6,"public SavableState",
|
---|
| 52 | create<Molecule>, create<Molecule>, create<Molecule>);
|
---|
| 53 |
|
---|
| 54 | Molecule::Molecule():
|
---|
| 55 | natoms_(0), r_(0), Z_(0), charges_(0), mass_(0), labels_(0)
|
---|
| 56 | {
|
---|
| 57 | pg_ = new PointGroup;
|
---|
| 58 | atominfo_ = new AtomInfo();
|
---|
| 59 | geometry_units_ = new Units("bohr");
|
---|
| 60 | nuniq_ = 0;
|
---|
| 61 | equiv_ = 0;
|
---|
| 62 | nequiv_ = 0;
|
---|
| 63 | atom_to_uniq_ = 0;
|
---|
| 64 | q_Z_ = atominfo_->string_to_Z("Q");
|
---|
| 65 | include_q_ = false;
|
---|
| 66 | include_qq_ = false;
|
---|
| 67 | init_symmetry_info();
|
---|
| 68 | }
|
---|
| 69 |
|
---|
| 70 | Molecule::Molecule(const Molecule& mol):
|
---|
| 71 | natoms_(0), r_(0), Z_(0), charges_(0), mass_(0), labels_(0)
|
---|
| 72 | {
|
---|
| 73 | nuniq_ = 0;
|
---|
| 74 | equiv_ = 0;
|
---|
| 75 | nequiv_ = 0;
|
---|
| 76 | atom_to_uniq_ = 0;
|
---|
| 77 | *this=mol;
|
---|
| 78 | }
|
---|
| 79 |
|
---|
| 80 | Molecule::~Molecule()
|
---|
| 81 | {
|
---|
| 82 | clear();
|
---|
| 83 | }
|
---|
| 84 |
|
---|
| 85 | void
|
---|
| 86 | Molecule::clear()
|
---|
| 87 | {
|
---|
| 88 | if (r_) {
|
---|
| 89 | delete[] r_[0];
|
---|
| 90 | delete[] r_;
|
---|
| 91 | r_ = 0;
|
---|
| 92 | }
|
---|
| 93 | if (labels_) {
|
---|
| 94 | for (int i=0; i<natoms_; i++) {
|
---|
| 95 | delete[] labels_[i];
|
---|
| 96 | }
|
---|
| 97 | delete[] labels_;
|
---|
| 98 | labels_ = 0;
|
---|
| 99 | }
|
---|
| 100 | delete[] charges_;
|
---|
| 101 | charges_ = 0;
|
---|
| 102 | delete[] mass_;
|
---|
| 103 | mass_ = 0;
|
---|
| 104 | delete[] Z_;
|
---|
| 105 | Z_ = 0;
|
---|
| 106 |
|
---|
| 107 | clear_symmetry_info();
|
---|
| 108 | }
|
---|
| 109 |
|
---|
| 110 | void
|
---|
| 111 | Molecule::throw_if_atom_duplicated(int begin, double tol)
|
---|
| 112 | {
|
---|
| 113 | for (int i=begin; i<natoms_; i++) {
|
---|
| 114 | SCVector3 ri(r_[i]);
|
---|
| 115 | for (int j=0; j<i; j++) {
|
---|
| 116 | SCVector3 rj(r_[j]);
|
---|
| 117 | if (ri.dist(rj) < tol) {
|
---|
| 118 | throw InputError("duplicated atom coordinate",
|
---|
| 119 | __FILE__, __LINE__, 0, 0, class_desc());
|
---|
| 120 | }
|
---|
| 121 | }
|
---|
| 122 | }
|
---|
| 123 | }
|
---|
| 124 |
|
---|
| 125 | Molecule::Molecule(const Ref<KeyVal>&input):
|
---|
| 126 | natoms_(0), r_(0), Z_(0), charges_(0), mass_(0), labels_(0)
|
---|
| 127 | {
|
---|
| 128 | nuniq_ = 0;
|
---|
| 129 | equiv_ = 0;
|
---|
| 130 | nequiv_ = 0;
|
---|
| 131 | atom_to_uniq_ = 0;
|
---|
| 132 |
|
---|
| 133 | KeyValValueboolean kvfalse(0);
|
---|
| 134 | include_q_ = input->booleanvalue("include_q",kvfalse);
|
---|
| 135 | include_qq_ = input->booleanvalue("include_qq",kvfalse);
|
---|
| 136 |
|
---|
| 137 | atominfo_ << input->describedclassvalue("atominfo");
|
---|
| 138 | if (atominfo_.null()) atominfo_ = new AtomInfo;
|
---|
| 139 | q_Z_ = atominfo_->string_to_Z("Q");
|
---|
| 140 | if (input->exists("pdb_file")) {
|
---|
| 141 | geometry_units_ = new Units("angstrom");
|
---|
| 142 | char* filename = input->pcharvalue("pdb_file");
|
---|
| 143 | read_pdb(filename);
|
---|
| 144 | delete[] filename;
|
---|
| 145 | }
|
---|
| 146 | else {
|
---|
| 147 | // check for old style units input first
|
---|
| 148 | if (input->booleanvalue("angstrom")
|
---|
| 149 | ||input->booleanvalue("angstroms")
|
---|
| 150 | ||input->booleanvalue("aangstrom")
|
---|
| 151 | ||input->booleanvalue("aangstroms")) {
|
---|
| 152 | geometry_units_ = new Units("angstrom");
|
---|
| 153 | }
|
---|
| 154 | // check for new style units input
|
---|
| 155 | else {
|
---|
| 156 | geometry_units_ = new Units(input->pcharvalue("unit"),
|
---|
| 157 | Units::Steal);
|
---|
| 158 | }
|
---|
| 159 |
|
---|
| 160 | double conv = geometry_units_->to_atomic_units();
|
---|
| 161 |
|
---|
| 162 | // get the number of atoms and make sure that the geometry and the
|
---|
| 163 | // atoms array have the same number of atoms.
|
---|
| 164 | // right now we read in the unique atoms...then we will symmetrize.
|
---|
| 165 | // the length of atoms must still equal the length of geometry, but
|
---|
| 166 | // we'll try to set up atom_labels such that different lengths are
|
---|
| 167 | // possible
|
---|
| 168 | int natom = input->count("geometry");
|
---|
| 169 | if (natom != input->count("atoms")) {
|
---|
| 170 | throw InputError("size of \"geometry\" != size of \"atoms\"",
|
---|
| 171 | __FILE__, __LINE__, 0, 0, class_desc());
|
---|
| 172 | }
|
---|
| 173 |
|
---|
| 174 | int i;
|
---|
| 175 | for (i=0; i<natom; i++) {
|
---|
| 176 | char *name, *label;
|
---|
| 177 | int ghost = input->booleanvalue("ghost",i);
|
---|
| 178 | double charge = input->doublevalue("charge",i);
|
---|
| 179 | int have_charge = input->error() == KeyVal::OK;
|
---|
| 180 | if (ghost) {
|
---|
| 181 | have_charge = 1;
|
---|
| 182 | charge = 0.0;
|
---|
| 183 | }
|
---|
| 184 | add_atom(atominfo_->string_to_Z(name = input->pcharvalue("atoms",i)),
|
---|
| 185 | input->doublevalue("geometry",i,0)*conv,
|
---|
| 186 | input->doublevalue("geometry",i,1)*conv,
|
---|
| 187 | input->doublevalue("geometry",i,2)*conv,
|
---|
| 188 | label = input->pcharvalue("atom_labels",i),
|
---|
| 189 | input->doublevalue("mass",i),
|
---|
| 190 | have_charge, charge
|
---|
| 191 | );
|
---|
| 192 | delete[] name;
|
---|
| 193 | delete[] label;
|
---|
| 194 | }
|
---|
| 195 | }
|
---|
| 196 |
|
---|
| 197 | char *symmetry = input->pcharvalue("symmetry");
|
---|
| 198 | double symtol = input->doublevalue("symmetry_tolerance",
|
---|
| 199 | KeyValValuedouble(1.0e-4));
|
---|
| 200 | nuniq_ = 0;
|
---|
| 201 | equiv_ = 0;
|
---|
| 202 | nequiv_ = 0;
|
---|
| 203 | atom_to_uniq_ = 0;
|
---|
| 204 | if (symmetry && !strcmp(symmetry,"auto")) {
|
---|
| 205 | set_point_group(highest_point_group(symtol), symtol*10.0);
|
---|
| 206 | }
|
---|
| 207 | else {
|
---|
| 208 | pg_ = new PointGroup(input);
|
---|
| 209 |
|
---|
| 210 | // translate to the origin of the symmetry frame
|
---|
| 211 | double r[3];
|
---|
| 212 | for (int i=0; i<3; i++) {
|
---|
| 213 | r[i] = -pg_->origin()[i];
|
---|
| 214 | pg_->origin()[i] = 0;
|
---|
| 215 | }
|
---|
| 216 | translate(r);
|
---|
| 217 |
|
---|
| 218 | if (input->booleanvalue("redundant_atoms")) {
|
---|
| 219 | init_symmetry_info();
|
---|
| 220 | cleanup_molecule(symtol);
|
---|
| 221 | }
|
---|
| 222 | else {
|
---|
| 223 | symmetrize();
|
---|
| 224 | // In case we were given redundant atoms, clean up
|
---|
| 225 | // the geometry so the symmetry is exact.
|
---|
| 226 | cleanup_molecule(symtol);
|
---|
| 227 | }
|
---|
| 228 | }
|
---|
| 229 | delete[] symmetry;
|
---|
| 230 | }
|
---|
| 231 |
|
---|
| 232 | Molecule&
|
---|
| 233 | Molecule::operator=(const Molecule& mol)
|
---|
| 234 | {
|
---|
| 235 | clear();
|
---|
| 236 |
|
---|
| 237 | pg_ = new PointGroup(*(mol.pg_.pointer()));
|
---|
| 238 | atominfo_ = mol.atominfo_;
|
---|
| 239 | geometry_units_ = new Units(mol.geometry_units_->string_rep());
|
---|
| 240 |
|
---|
| 241 | q_Z_ = mol.q_Z_;
|
---|
| 242 | include_q_ = mol.include_q_;
|
---|
| 243 | include_qq_ = mol.include_qq_;
|
---|
| 244 | q_atoms_ = mol.q_atoms_;
|
---|
| 245 | non_q_atoms_ = mol.non_q_atoms_;
|
---|
| 246 |
|
---|
| 247 | natoms_ = mol.natoms_;
|
---|
| 248 |
|
---|
| 249 | if (natoms_) {
|
---|
| 250 | if (mol.mass_) {
|
---|
| 251 | mass_ = new double[natoms_];
|
---|
| 252 | memcpy(mass_,mol.mass_,natoms_*sizeof(double));
|
---|
| 253 | }
|
---|
| 254 | if (mol.charges_) {
|
---|
| 255 | charges_ = new double[natoms_];
|
---|
| 256 | memcpy(charges_,mol.charges_,natoms_*sizeof(double));
|
---|
| 257 | }
|
---|
| 258 | if (mol.labels_) {
|
---|
| 259 | labels_ = new char *[natoms_];
|
---|
| 260 | for (int i=0; i<natoms_; i++) {
|
---|
| 261 | if (mol.labels_[i]) {
|
---|
| 262 | labels_[i] = strcpy(new char[strlen(mol.labels_[i])+1],
|
---|
| 263 | mol.labels_[i]);
|
---|
| 264 | }
|
---|
| 265 | else labels_[i] = 0;
|
---|
| 266 | }
|
---|
| 267 | }
|
---|
| 268 | r_ = new double*[natoms_];
|
---|
| 269 | r_[0] = new double[natoms_*3];
|
---|
| 270 | for (int i=0; i<natoms_; i++) {
|
---|
| 271 | r_[i] = &(r_[0][i*3]);
|
---|
| 272 | }
|
---|
| 273 | memcpy(r_[0], mol.r_[0], natoms_*3*sizeof(double));
|
---|
| 274 | Z_ = new int[natoms_];
|
---|
| 275 | memcpy(Z_, mol.Z_, natoms_*sizeof(int));
|
---|
| 276 | }
|
---|
| 277 |
|
---|
| 278 | init_symmetry_info();
|
---|
| 279 |
|
---|
| 280 | return *this;
|
---|
| 281 | }
|
---|
| 282 |
|
---|
| 283 | void
|
---|
| 284 | Molecule::add_atom(int Z,double x,double y,double z,
|
---|
| 285 | const char *label,double mass,
|
---|
| 286 | int have_charge, double charge)
|
---|
| 287 | {
|
---|
| 288 | int i;
|
---|
| 289 |
|
---|
| 290 | // allocate new arrays
|
---|
| 291 | int *newZ = new int[natoms_+1];
|
---|
| 292 | double **newr = new double*[natoms_+1];
|
---|
| 293 | double *newr0 = new double[(natoms_+1)*3];
|
---|
| 294 | char **newlabels = 0;
|
---|
| 295 | if (label || labels_) {
|
---|
| 296 | newlabels = new char*[natoms_+1];
|
---|
| 297 | }
|
---|
| 298 | double *newcharges = 0;
|
---|
| 299 | if (have_charge || charges_) {
|
---|
| 300 | newcharges = new double[natoms_+1];
|
---|
| 301 | }
|
---|
| 302 | double *newmass = 0;
|
---|
| 303 | if (mass_ || mass != 0.0) {
|
---|
| 304 | newmass = new double[natoms_+1];
|
---|
| 305 | }
|
---|
| 306 |
|
---|
| 307 | // setup the r_ pointers
|
---|
| 308 | for (i=0; i<=natoms_; i++) {
|
---|
| 309 | newr[i] = &(newr0[i*3]);
|
---|
| 310 | }
|
---|
| 311 |
|
---|
| 312 | // copy old data to new arrays
|
---|
| 313 | if (natoms_) {
|
---|
| 314 | memcpy(newZ,Z_,sizeof(int)*natoms_);
|
---|
| 315 | memcpy(newr0,r_[0],sizeof(double)*natoms_*3);
|
---|
| 316 | if (labels_) {
|
---|
| 317 | memcpy(newlabels,labels_,sizeof(char*)*natoms_);
|
---|
| 318 | }
|
---|
| 319 | else if (newlabels) {
|
---|
| 320 | memset(newlabels,0,sizeof(char*)*natoms_);
|
---|
| 321 | }
|
---|
| 322 | if (charges_) {
|
---|
| 323 | memcpy(newcharges,charges_,sizeof(double)*natoms_);
|
---|
| 324 | }
|
---|
| 325 | else if (newcharges) {
|
---|
| 326 | for (i=0; i<natoms_; i++) newcharges[i] = Z_[i];
|
---|
| 327 | }
|
---|
| 328 | if (mass_) {
|
---|
| 329 | memcpy(newmass,mass_,sizeof(double)*natoms_);
|
---|
| 330 | }
|
---|
| 331 | else if (newmass) {
|
---|
| 332 | memset(newmass,0,sizeof(double)*natoms_);
|
---|
| 333 | }
|
---|
| 334 | }
|
---|
| 335 |
|
---|
| 336 | // delete old data
|
---|
| 337 | delete[] Z_;
|
---|
| 338 | if (r_) {
|
---|
| 339 | delete[] r_[0];
|
---|
| 340 | delete[] r_;
|
---|
| 341 | }
|
---|
| 342 | delete[] labels_;
|
---|
| 343 | delete[] charges_;
|
---|
| 344 | delete[] mass_;
|
---|
| 345 |
|
---|
| 346 | // setup new pointers
|
---|
| 347 | Z_ = newZ;
|
---|
| 348 | r_ = newr;
|
---|
| 349 | labels_ = newlabels;
|
---|
| 350 | charges_ = newcharges;
|
---|
| 351 | mass_ = newmass;
|
---|
| 352 |
|
---|
| 353 | // copy info for this atom into arrays
|
---|
| 354 | Z_[natoms_] = Z;
|
---|
| 355 | r_[natoms_][0] = x;
|
---|
| 356 | r_[natoms_][1] = y;
|
---|
| 357 | r_[natoms_][2] = z;
|
---|
| 358 | if (mass_) mass_[natoms_] = mass;
|
---|
| 359 | if (label) {
|
---|
| 360 | labels_[natoms_] = strcpy(new char[strlen(label)+1],label);
|
---|
| 361 | }
|
---|
| 362 | else if (labels_) {
|
---|
| 363 | labels_[natoms_] = 0;
|
---|
| 364 | }
|
---|
| 365 | if (have_charge) {
|
---|
| 366 | charges_[natoms_] = charge;
|
---|
| 367 | }
|
---|
| 368 | else if (charges_) {
|
---|
| 369 | charges_[natoms_] = Z;
|
---|
| 370 | }
|
---|
| 371 |
|
---|
| 372 | if (Z == q_Z_) {
|
---|
| 373 | q_atoms_.push_back(natoms_);
|
---|
| 374 | }
|
---|
| 375 | else {
|
---|
| 376 | non_q_atoms_.push_back(natoms_);
|
---|
| 377 | }
|
---|
| 378 |
|
---|
| 379 | natoms_++;
|
---|
| 380 |
|
---|
| 381 | throw_if_atom_duplicated(natoms_-1);
|
---|
| 382 | }
|
---|
| 383 |
|
---|
| 384 | void
|
---|
| 385 | Molecule::print_parsedkeyval(ostream& os,
|
---|
| 386 | int print_pg,
|
---|
| 387 | int print_unit,
|
---|
| 388 | int number_atoms) const
|
---|
| 389 | {
|
---|
| 390 | int i;
|
---|
| 391 |
|
---|
| 392 | double conv = geometry_units_->from_atomic_units();
|
---|
| 393 |
|
---|
| 394 | if (print_pg) pg_->print(os);
|
---|
| 395 | if (print_unit && geometry_units_->string_rep()) {
|
---|
| 396 | os << indent
|
---|
| 397 | << "unit = \"" << geometry_units_->string_rep() << "\""
|
---|
| 398 | << endl;
|
---|
| 399 | }
|
---|
| 400 | os << indent << "{";
|
---|
| 401 | if (number_atoms) os << scprintf("%3s", "n");
|
---|
| 402 | os << scprintf(" %5s", "atoms");
|
---|
| 403 | if (labels_) os << scprintf(" %11s", "atom_labels");
|
---|
| 404 | int int_charges = 1;
|
---|
| 405 | if (charges_) {
|
---|
| 406 | for (i=0;i<natom();i++) if (charges_[i]!=(int)charges_[i]) int_charges=0;
|
---|
| 407 | if (int_charges) {
|
---|
| 408 | os << scprintf(" %7s", "charge");
|
---|
| 409 | }
|
---|
| 410 | else {
|
---|
| 411 | os << scprintf(" %17s", "charge");
|
---|
| 412 | }
|
---|
| 413 | }
|
---|
| 414 | os << scprintf(" %16s", "")
|
---|
| 415 | << scprintf(" %16s", "geometry ")
|
---|
| 416 | << scprintf(" %16s ", "");
|
---|
| 417 | os << "}={" << endl;
|
---|
| 418 | for (i=0; i<natom(); i++) {
|
---|
| 419 | os << indent;
|
---|
| 420 | if (number_atoms) os << scprintf(" %3d", i+1);
|
---|
| 421 | std::string symbol(atom_symbol(i));
|
---|
| 422 | os << scprintf(" %5s", symbol.c_str());
|
---|
| 423 | if (labels_) {
|
---|
| 424 | const char *lab = labels_[i];
|
---|
| 425 | if (lab == 0) lab = "";
|
---|
| 426 | char *qlab = new char[strlen(lab)+3];
|
---|
| 427 | strcpy(qlab,"\"");
|
---|
| 428 | strcat(qlab,lab);
|
---|
| 429 | strcat(qlab,"\"");
|
---|
| 430 | os << scprintf(" %11s",qlab);
|
---|
| 431 | delete[] qlab;
|
---|
| 432 | }
|
---|
| 433 | if (charges_) {
|
---|
| 434 | if (int_charges) os << scprintf(" %7.4f", charges_[i]);
|
---|
| 435 | else os << scprintf(" %17.15f", charges_[i]);
|
---|
| 436 | }
|
---|
| 437 | os << scprintf(" [% 16.10f", conv * r(i,0))
|
---|
| 438 | << scprintf(" % 16.10f", conv * r(i,1))
|
---|
| 439 | << scprintf(" % 16.10f]", conv * r(i,2))
|
---|
| 440 | << endl;
|
---|
| 441 | }
|
---|
| 442 | os << indent << "}" << endl;
|
---|
| 443 | }
|
---|
| 444 |
|
---|
| 445 | void
|
---|
| 446 | Molecule::print(ostream& os) const
|
---|
| 447 | {
|
---|
| 448 | int i;
|
---|
| 449 |
|
---|
| 450 | MolecularFormula *mf = new MolecularFormula(this);
|
---|
| 451 | os << indent
|
---|
| 452 | << "Molecular formula: " << mf->formula() << endl;
|
---|
| 453 | delete mf;
|
---|
| 454 |
|
---|
| 455 | os << indent << "molecule<Molecule>: (" << endl;
|
---|
| 456 | os << incindent;
|
---|
| 457 | print_parsedkeyval(os);
|
---|
| 458 | os << decindent;
|
---|
| 459 | os << indent << ")" << endl;
|
---|
| 460 |
|
---|
| 461 | os << indent << "Atomic Masses:" << endl;
|
---|
| 462 | for (i=0; i<natom(); i+=5) {
|
---|
| 463 | os << indent;
|
---|
| 464 | for (int j=i; j<i+5 && j<natom(); j++) {
|
---|
| 465 | os << scprintf(" %10.5f", mass(j));
|
---|
| 466 | }
|
---|
| 467 | os << endl;
|
---|
| 468 | }
|
---|
| 469 | }
|
---|
| 470 |
|
---|
| 471 | int
|
---|
| 472 | Molecule::atom_label_to_index(const char *l) const
|
---|
| 473 | {
|
---|
| 474 | int i;
|
---|
| 475 | for (i=0; i<natom(); i++) {
|
---|
| 476 | if (label(i) && !strcmp(l,label(i))) return i;
|
---|
| 477 | }
|
---|
| 478 | return -1;
|
---|
| 479 | }
|
---|
| 480 |
|
---|
| 481 | double*
|
---|
| 482 | Molecule::charges() const
|
---|
| 483 | {
|
---|
| 484 | double*result = new double[natoms_];
|
---|
| 485 | if (charges_) {
|
---|
| 486 | memcpy(result, charges_, sizeof(double)*natom());
|
---|
| 487 | }
|
---|
| 488 | else {
|
---|
| 489 | for (int i=0; i<natom(); i++) result[i] = Z_[i];
|
---|
| 490 | }
|
---|
| 491 | return result;
|
---|
| 492 | }
|
---|
| 493 |
|
---|
| 494 | double
|
---|
| 495 | Molecule::charge(int iatom) const
|
---|
| 496 | {
|
---|
| 497 | if (charges_) return charges_[iatom];
|
---|
| 498 | return Z_[iatom];
|
---|
| 499 | }
|
---|
| 500 |
|
---|
| 501 | double
|
---|
| 502 | Molecule::nuclear_charge() const
|
---|
| 503 | {
|
---|
| 504 | double c = 0.0;
|
---|
| 505 | if (include_q_) {
|
---|
| 506 | for (int i=0; i<natom(); i++) {
|
---|
| 507 | c += charge(i);
|
---|
| 508 | }
|
---|
| 509 | }
|
---|
| 510 | else {
|
---|
| 511 | for (int ii=0; ii<non_q_atoms_.size(); ii++) {
|
---|
| 512 | c += charge(non_q_atoms_[ii]);
|
---|
| 513 | }
|
---|
| 514 | }
|
---|
| 515 | return c;
|
---|
| 516 | }
|
---|
| 517 |
|
---|
| 518 | void Molecule::save_data_state(StateOut& so)
|
---|
| 519 | {
|
---|
| 520 | so.put(include_q_);
|
---|
| 521 | so.put(include_qq_);
|
---|
| 522 | so.put(natoms_);
|
---|
| 523 | SavableState::save_state(pg_.pointer(),so);
|
---|
| 524 | SavableState::save_state(geometry_units_.pointer(),so);
|
---|
| 525 | SavableState::save_state(atominfo_.pointer(),so);
|
---|
| 526 | if (natoms_) {
|
---|
| 527 | so.put(Z_, natoms_);
|
---|
| 528 | so.put_array_double(r_[0], natoms_*3);
|
---|
| 529 | so.put(charges_,natoms_);
|
---|
| 530 | }
|
---|
| 531 | if (mass_) {
|
---|
| 532 | so.put(1);
|
---|
| 533 | so.put_array_double(mass_, natoms_);
|
---|
| 534 | }
|
---|
| 535 | else {
|
---|
| 536 | so.put(0);
|
---|
| 537 | }
|
---|
| 538 | if (labels_){
|
---|
| 539 | so.put(1);
|
---|
| 540 | for (int i=0; i<natoms_; i++) {
|
---|
| 541 | so.putstring(labels_[i]);
|
---|
| 542 | }
|
---|
| 543 | }
|
---|
| 544 | else {
|
---|
| 545 | so.put(0);
|
---|
| 546 | }
|
---|
| 547 | }
|
---|
| 548 |
|
---|
| 549 | Molecule::Molecule(StateIn& si):
|
---|
| 550 | SavableState(si),
|
---|
| 551 | natoms_(0), r_(0), Z_(0), mass_(0), labels_(0)
|
---|
| 552 | {
|
---|
| 553 | if (si.version(::class_desc<Molecule>()) < 4) {
|
---|
| 554 | throw FileOperationFailed("cannot restore from old molecules",
|
---|
| 555 | __FILE__, __LINE__, 0,
|
---|
| 556 | FileOperationFailed::Corrupt,
|
---|
| 557 | class_desc());
|
---|
| 558 | }
|
---|
| 559 | if (si.version(::class_desc<Molecule>()) < 6) {
|
---|
| 560 | include_q_ = false;
|
---|
| 561 | include_qq_ = false;
|
---|
| 562 | }
|
---|
| 563 | else {
|
---|
| 564 | si.get(include_q_);
|
---|
| 565 | si.get(include_qq_);
|
---|
| 566 | }
|
---|
| 567 | si.get(natoms_);
|
---|
| 568 | pg_ << SavableState::restore_state(si);
|
---|
| 569 | geometry_units_ << SavableState::restore_state(si);
|
---|
| 570 | atominfo_ << SavableState::restore_state(si);
|
---|
| 571 | if (natoms_) {
|
---|
| 572 | si.get(Z_);
|
---|
| 573 | r_ = new double*[natoms_];
|
---|
| 574 | r_[0] = new double[natoms_*3];
|
---|
| 575 | si.get_array_double(r_[0],natoms_*3);
|
---|
| 576 | for (int i=1; i<natoms_; i++) {
|
---|
| 577 | r_[i] = &(r_[0][i*3]);
|
---|
| 578 | }
|
---|
| 579 | if (si.version(::class_desc<Molecule>()) > 4) {
|
---|
| 580 | si.get(charges_);
|
---|
| 581 | }
|
---|
| 582 | else {
|
---|
| 583 | charges_ = 0;
|
---|
| 584 | }
|
---|
| 585 | }
|
---|
| 586 | int test;
|
---|
| 587 | si.get(test);
|
---|
| 588 | if (test) {
|
---|
| 589 | mass_ = new double[natoms_];
|
---|
| 590 | si.get_array_double(mass_, natoms_);
|
---|
| 591 | }
|
---|
| 592 | si.get(test);
|
---|
| 593 | if (test){
|
---|
| 594 | labels_ = new char*[natoms_];
|
---|
| 595 | for (int i=0; i<natoms_; i++) {
|
---|
| 596 | si.getstring(labels_[i]);
|
---|
| 597 | }
|
---|
| 598 | }
|
---|
| 599 |
|
---|
| 600 | for (int i=0; i<natoms_; i++) {
|
---|
| 601 | if (Z_[i] == q_Z_) {
|
---|
| 602 | q_atoms_.push_back(i);
|
---|
| 603 | }
|
---|
| 604 | else {
|
---|
| 605 | non_q_atoms_.push_back(i);
|
---|
| 606 | }
|
---|
| 607 | }
|
---|
| 608 |
|
---|
| 609 | nuniq_ = 0;
|
---|
| 610 | equiv_ = 0;
|
---|
| 611 | nequiv_ = 0;
|
---|
| 612 | atom_to_uniq_ = 0;
|
---|
| 613 | init_symmetry_info();
|
---|
| 614 | }
|
---|
| 615 |
|
---|
| 616 | int
|
---|
| 617 | Molecule::atom_to_unique_offset(int iatom) const
|
---|
| 618 | {
|
---|
| 619 | int iuniq = atom_to_uniq_[iatom];
|
---|
| 620 | int nequiv = nequiv_[iuniq];
|
---|
| 621 | for (int i=0; i<nequiv; i++) {
|
---|
| 622 | if (equiv_[iuniq][i] == iatom) return i;
|
---|
| 623 | }
|
---|
| 624 | ExEnv::errn() << "Molecule::atom_to_unique_offset: internal error"
|
---|
| 625 | << endl;
|
---|
| 626 | return -1;
|
---|
| 627 | }
|
---|
| 628 |
|
---|
| 629 | void
|
---|
| 630 | Molecule::set_point_group(const Ref<PointGroup>&ppg, double tol)
|
---|
| 631 | {
|
---|
| 632 | ExEnv::out0() << indent
|
---|
| 633 | << "Molecule: setting point group to " << ppg->symbol()
|
---|
| 634 | << endl;
|
---|
| 635 | pg_ = new PointGroup(*ppg.pointer());
|
---|
| 636 |
|
---|
| 637 | double r[3];
|
---|
| 638 | for (int i=0; i<3; i++) {
|
---|
| 639 | r[i] = -pg_->origin()[i];
|
---|
| 640 | pg_->origin()[i] = 0;
|
---|
| 641 | }
|
---|
| 642 | translate(r);
|
---|
| 643 |
|
---|
| 644 | clear_symmetry_info();
|
---|
| 645 | init_symmetry_info();
|
---|
| 646 |
|
---|
| 647 | cleanup_molecule(tol);
|
---|
| 648 | }
|
---|
| 649 |
|
---|
| 650 | Ref<PointGroup>
|
---|
| 651 | Molecule::point_group() const
|
---|
| 652 | {
|
---|
| 653 | return pg_;
|
---|
| 654 | }
|
---|
| 655 |
|
---|
| 656 | SCVector3
|
---|
| 657 | Molecule::center_of_mass() const
|
---|
| 658 | {
|
---|
| 659 | SCVector3 ret;
|
---|
| 660 | double M;
|
---|
| 661 |
|
---|
| 662 | ret = 0.0;
|
---|
| 663 | M = 0.0;
|
---|
| 664 |
|
---|
| 665 | for (int i=0; i < natom(); i++) {
|
---|
| 666 | double m = mass(i);
|
---|
| 667 | ret += m * SCVector3(r(i));
|
---|
| 668 | M += m;
|
---|
| 669 | }
|
---|
| 670 |
|
---|
| 671 | ret *= 1.0/M;
|
---|
| 672 |
|
---|
| 673 | return ret;
|
---|
| 674 | }
|
---|
| 675 |
|
---|
| 676 | double
|
---|
| 677 | Molecule::nuclear_repulsion_energy()
|
---|
| 678 | {
|
---|
| 679 | double e=0.0;
|
---|
| 680 |
|
---|
| 681 | // non_q non_q terms
|
---|
| 682 | for (int ii=1; ii < non_q_atoms_.size(); ii++) {
|
---|
| 683 | int i = non_q_atoms_[ii];
|
---|
| 684 | SCVector3 ai(r(i));
|
---|
| 685 | double Zi = charge(i);
|
---|
| 686 |
|
---|
| 687 | for (int jj=0; jj < ii; jj++) {
|
---|
| 688 | int j = non_q_atoms_[jj];
|
---|
| 689 | SCVector3 aj(r(j));
|
---|
| 690 | e += Zi * charge(j) / ai.dist(aj);
|
---|
| 691 | }
|
---|
| 692 | }
|
---|
| 693 |
|
---|
| 694 | // non_q q terms
|
---|
| 695 | for (int ii=0; ii < q_atoms_.size(); ii++) {
|
---|
| 696 | int i = q_atoms_[ii];
|
---|
| 697 | SCVector3 ai(r(i));
|
---|
| 698 | double Zi = charge(i);
|
---|
| 699 |
|
---|
| 700 | for (int jj=0; jj < non_q_atoms_.size(); jj++) {
|
---|
| 701 | int j = non_q_atoms_[jj];
|
---|
| 702 | SCVector3 aj(r(j));
|
---|
| 703 | e += Zi * charge(j) / ai.dist(aj);
|
---|
| 704 | }
|
---|
| 705 | }
|
---|
| 706 |
|
---|
| 707 | if (include_qq_) {
|
---|
| 708 | // q q terms
|
---|
| 709 | for (int ii=1; ii < q_atoms_.size(); ii++) {
|
---|
| 710 | int i = q_atoms_[ii];
|
---|
| 711 | SCVector3 ai(r(i));
|
---|
| 712 | double Zi = charge(i);
|
---|
| 713 |
|
---|
| 714 | for (int jj=0; jj < ii; jj++) {
|
---|
| 715 | int j = q_atoms_[jj];
|
---|
| 716 | SCVector3 aj(r(j));
|
---|
| 717 | e += Zi * charge(j) / ai.dist(aj);
|
---|
| 718 | }
|
---|
| 719 | }
|
---|
| 720 | }
|
---|
| 721 |
|
---|
| 722 | return e;
|
---|
| 723 | }
|
---|
| 724 |
|
---|
| 725 | void
|
---|
| 726 | Molecule::nuclear_repulsion_1der(int center, double xyz[3])
|
---|
| 727 | {
|
---|
| 728 | int i,j,k;
|
---|
| 729 | double rd[3],r2;
|
---|
| 730 | double factor;
|
---|
| 731 |
|
---|
| 732 | xyz[0] = 0.0;
|
---|
| 733 | xyz[1] = 0.0;
|
---|
| 734 | xyz[2] = 0.0;
|
---|
| 735 |
|
---|
| 736 | SCVector3 r_center(r(center));
|
---|
| 737 | double Z_center = charge(center);
|
---|
| 738 | bool center_is_Q = (atom_symbol(center) == "Q");
|
---|
| 739 |
|
---|
| 740 | // this handles center = Q or non_Q and atom = non_Q
|
---|
| 741 | for (int ii=0; ii < non_q_atoms_.size(); ii++) {
|
---|
| 742 | int i = non_q_atoms_[ii];
|
---|
| 743 | if (i == center) continue;
|
---|
| 744 | SCVector3 r_i(r(i));
|
---|
| 745 |
|
---|
| 746 | r2 = 0.0;
|
---|
| 747 | for (k=0; k < 3; k++) {
|
---|
| 748 | rd[k] = r_center[k] - r_i[k];
|
---|
| 749 | r2 += rd[k]*rd[k];
|
---|
| 750 | }
|
---|
| 751 | factor = - Z_center * charge(i) * pow(r2,-1.5);
|
---|
| 752 | for (k=0; k<3; k++) {
|
---|
| 753 | xyz[k] += factor * rd[k];
|
---|
| 754 | }
|
---|
| 755 | }
|
---|
| 756 |
|
---|
| 757 | // this handles center = Q or non_Q and atom = Q
|
---|
| 758 | for (int ii=0; ii < q_atoms_.size(); ii++) {
|
---|
| 759 | int i = q_atoms_[ii];
|
---|
| 760 | if (i == center || (!include_qq_ && center_is_Q)) continue;
|
---|
| 761 | SCVector3 r_i(r(i));
|
---|
| 762 |
|
---|
| 763 | r2 = 0.0;
|
---|
| 764 | for (k=0; k < 3; k++) {
|
---|
| 765 | rd[k] = r_center[k] - r_i[k];
|
---|
| 766 | r2 += rd[k]*rd[k];
|
---|
| 767 | }
|
---|
| 768 | factor = - Z_center * charge(i) * pow(r2,-1.5);
|
---|
| 769 | for (k=0; k<3; k++) {
|
---|
| 770 | xyz[k] += factor * rd[k];
|
---|
| 771 | }
|
---|
| 772 | }
|
---|
| 773 | }
|
---|
| 774 |
|
---|
| 775 | void
|
---|
| 776 | Molecule::nuclear_charge_efield(const double *charges,
|
---|
| 777 | const double *position, double *efield)
|
---|
| 778 | {
|
---|
| 779 | double tmp;
|
---|
| 780 | double rd[3];
|
---|
| 781 |
|
---|
| 782 | for (int i=0; i<3; i++) efield[i] = 0.0;
|
---|
| 783 |
|
---|
| 784 | if (include_q_) {
|
---|
| 785 | for (int i=0; i<natoms_; i++) {
|
---|
| 786 | SCVector3 a(r(i));
|
---|
| 787 | tmp = 0.0;
|
---|
| 788 | for (int j=0; j<3; j++) {
|
---|
| 789 | rd[j] = position[j] - a[j];
|
---|
| 790 | tmp += rd[j]*rd[j];
|
---|
| 791 | }
|
---|
| 792 | tmp = charges[i]/(tmp*sqrt(tmp));
|
---|
| 793 | for (int j=0; j<3; j++) {
|
---|
| 794 | efield[j] += rd[j] * tmp;
|
---|
| 795 | }
|
---|
| 796 | }
|
---|
| 797 | }
|
---|
| 798 | else {
|
---|
| 799 | for (int ii=0; ii<non_q_atoms_.size(); ii++) {
|
---|
| 800 | int i = non_q_atoms_[ii];
|
---|
| 801 | SCVector3 a(r(i));
|
---|
| 802 | tmp = 0.0;
|
---|
| 803 | for (int j=0; j<3; j++) {
|
---|
| 804 | rd[j] = position[j] - a[j];
|
---|
| 805 | tmp += rd[j]*rd[j];
|
---|
| 806 | }
|
---|
| 807 | tmp = charges[i]/(tmp*sqrt(tmp));
|
---|
| 808 | for (int j=0; j<3; j++) {
|
---|
| 809 | efield[j] += rd[j] * tmp;
|
---|
| 810 | }
|
---|
| 811 | }
|
---|
| 812 | }
|
---|
| 813 | }
|
---|
| 814 |
|
---|
| 815 | void
|
---|
| 816 | Molecule::nuclear_efield(const double *position, double *efield)
|
---|
| 817 | {
|
---|
| 818 | double tmp;
|
---|
| 819 | double rd[3];
|
---|
| 820 |
|
---|
| 821 | for (int i=0; i<3; i++) efield[i] = 0.0;
|
---|
| 822 |
|
---|
| 823 | if (include_q_) {
|
---|
| 824 | for (int i=0; i<natoms_; i++) {
|
---|
| 825 | SCVector3 a(r(i));
|
---|
| 826 | tmp = 0.0;
|
---|
| 827 | for (int j=0; j<3; j++) {
|
---|
| 828 | rd[j] = position[j] - a[j];
|
---|
| 829 | tmp += rd[j]*rd[j];
|
---|
| 830 | }
|
---|
| 831 | tmp = charge(i)/(tmp*sqrt(tmp));
|
---|
| 832 | for (int j=0; j<3; j++) {
|
---|
| 833 | efield[j] += rd[j] * tmp;
|
---|
| 834 | }
|
---|
| 835 | }
|
---|
| 836 | }
|
---|
| 837 | else {
|
---|
| 838 | for (int ii=0; ii<non_q_atoms_.size(); ii++) {
|
---|
| 839 | int i = non_q_atoms_[ii];
|
---|
| 840 | SCVector3 a(r(i));
|
---|
| 841 | tmp = 0.0;
|
---|
| 842 | for (int j=0; j<3; j++) {
|
---|
| 843 | rd[j] = position[j] - a[j];
|
---|
| 844 | tmp += rd[j]*rd[j];
|
---|
| 845 | }
|
---|
| 846 | tmp = charge(i)/(tmp*sqrt(tmp));
|
---|
| 847 | for (int j=0; j<3; j++) {
|
---|
| 848 | efield[j] += rd[j] * tmp;
|
---|
| 849 | }
|
---|
| 850 | }
|
---|
| 851 | }
|
---|
| 852 | }
|
---|
| 853 |
|
---|
| 854 | int
|
---|
| 855 | Molecule::atom_at_position(double *v, double tol) const
|
---|
| 856 | {
|
---|
| 857 | SCVector3 p(v);
|
---|
| 858 | for (int i=0; i < natom(); i++) {
|
---|
| 859 | SCVector3 ai(r(i));
|
---|
| 860 | if (p.dist(ai) < tol) return i;
|
---|
| 861 | }
|
---|
| 862 | return -1;
|
---|
| 863 | }
|
---|
| 864 |
|
---|
| 865 | void
|
---|
| 866 | Molecule::symmetrize(const Ref<PointGroup> &pg, double tol)
|
---|
| 867 | {
|
---|
| 868 | pg_ = new PointGroup(pg);
|
---|
| 869 |
|
---|
| 870 | // translate to the origin of the symmetry frame
|
---|
| 871 | double r[3];
|
---|
| 872 | for (int i=0; i<3; i++) {
|
---|
| 873 | r[i] = -pg_->origin()[i];
|
---|
| 874 | pg_->origin()[i] = 0;
|
---|
| 875 | }
|
---|
| 876 | translate(r);
|
---|
| 877 |
|
---|
| 878 | symmetrize(tol);
|
---|
| 879 | }
|
---|
| 880 |
|
---|
| 881 | // We are given a molecule which may or may not have just the symmetry
|
---|
| 882 | // distinct atoms in it. We have to go through the existing set of atoms,
|
---|
| 883 | // perform each symmetry operation in the point group on each of them, and
|
---|
| 884 | // then add the new atom if it isn't in the list already
|
---|
| 885 |
|
---|
| 886 | void
|
---|
| 887 | Molecule::symmetrize(double tol)
|
---|
| 888 | {
|
---|
| 889 | // if molecule is c1, don't do anything
|
---|
| 890 | if (!strcmp(this->point_group()->symbol(),"c1")) {
|
---|
| 891 | init_symmetry_info();
|
---|
| 892 | return;
|
---|
| 893 | }
|
---|
| 894 |
|
---|
| 895 | clear_symmetry_info();
|
---|
| 896 |
|
---|
| 897 | Molecule *newmol = new Molecule(*this);
|
---|
| 898 |
|
---|
| 899 | CharacterTable ct = this->point_group()->char_table();
|
---|
| 900 |
|
---|
| 901 | SCVector3 np;
|
---|
| 902 | SymmetryOperation so;
|
---|
| 903 |
|
---|
| 904 | for (int i=0; i < natom(); i++) {
|
---|
| 905 | SCVector3 ac(r(i));
|
---|
| 906 |
|
---|
| 907 | for (int g=0; g < ct.order(); g++) {
|
---|
| 908 | so = ct.symm_operation(g);
|
---|
| 909 | for (int ii=0; ii < 3; ii++) {
|
---|
| 910 | np[ii]=0;
|
---|
| 911 | for (int jj=0; jj < 3; jj++) np[ii] += so(ii,jj) * ac[jj];
|
---|
| 912 | }
|
---|
| 913 |
|
---|
| 914 | int atom = newmol->atom_at_position(np.data(), tol);
|
---|
| 915 | if (atom < 0) {
|
---|
| 916 | newmol->add_atom(Z_[i],np[0],np[1],np[2],label(i));
|
---|
| 917 | }
|
---|
| 918 | else {
|
---|
| 919 | if (Z(i) != newmol->Z(atom)
|
---|
| 920 | || fabs(mass(i)-newmol->mass(atom))>1.0e-10) {
|
---|
| 921 | throw ToleranceExceeded("symmetrize: atom mismatch",
|
---|
| 922 | __FILE__, __LINE__,
|
---|
| 923 | 1.0e-10, fabs(mass(i)-newmol->mass(atom)),
|
---|
| 924 | class_desc());
|
---|
| 925 | }
|
---|
| 926 | }
|
---|
| 927 | }
|
---|
| 928 | }
|
---|
| 929 |
|
---|
| 930 | Ref<Units> saved_units = geometry_units_;
|
---|
| 931 | *this = *newmol;
|
---|
| 932 | geometry_units_ = saved_units;
|
---|
| 933 | delete newmol;
|
---|
| 934 |
|
---|
| 935 | init_symmetry_info();
|
---|
| 936 | }
|
---|
| 937 |
|
---|
| 938 | void
|
---|
| 939 | Molecule::translate(const double *r)
|
---|
| 940 | {
|
---|
| 941 | for (int i=0; i < natom(); i++) {
|
---|
| 942 | r_[i][0] += r[0];
|
---|
| 943 | r_[i][1] += r[1];
|
---|
| 944 | r_[i][2] += r[2];
|
---|
| 945 | }
|
---|
| 946 | }
|
---|
| 947 |
|
---|
| 948 | // move the molecule to the center of mass
|
---|
| 949 | void
|
---|
| 950 | Molecule::move_to_com()
|
---|
| 951 | {
|
---|
| 952 | SCVector3 com = -center_of_mass();
|
---|
| 953 | translate(com.data());
|
---|
| 954 | }
|
---|
| 955 |
|
---|
| 956 | // find the 3 principal coordinate axes, and rotate the molecule to be
|
---|
| 957 | // aligned along them. also rotate the symmetry frame contained in point_group
|
---|
| 958 | void
|
---|
| 959 | Molecule::transform_to_principal_axes(int trans_frame)
|
---|
| 960 | {
|
---|
| 961 | // mol_move_to_com(mol);
|
---|
| 962 |
|
---|
| 963 | double *inert[3], inert_dat[9], *evecs[3], evecs_dat[9];
|
---|
| 964 | double evals[3];
|
---|
| 965 |
|
---|
| 966 | int i,j,k;
|
---|
| 967 | for (i=0; i < 3; i++) {
|
---|
| 968 | inert[i] = &inert_dat[i*3];
|
---|
| 969 | evecs[i] = &evecs_dat[i*3];
|
---|
| 970 | }
|
---|
| 971 | memset(inert_dat,'\0',sizeof(double)*9);
|
---|
| 972 | memset(evecs_dat,'\0',sizeof(double)*9);
|
---|
| 973 |
|
---|
| 974 | for (i=0; i < natom(); i++) {
|
---|
| 975 | SCVector3 ac(r(i));
|
---|
| 976 | double m=mass(i);
|
---|
| 977 | inert[0][0] += m * (ac[1]*ac[1] + ac[2]*ac[2]);
|
---|
| 978 | inert[1][0] -= m * ac[0]*ac[1];
|
---|
| 979 | inert[1][1] += m * (ac[0]*ac[0] + ac[2]*ac[2]);
|
---|
| 980 | inert[2][0] -= m * ac[0]*ac[2];
|
---|
| 981 | inert[2][1] -= m * ac[1]*ac[2];
|
---|
| 982 | inert[2][2] += m * (ac[0]*ac[0] + ac[1]*ac[1]);
|
---|
| 983 | }
|
---|
| 984 |
|
---|
| 985 | inert[0][1] = inert[1][0];
|
---|
| 986 | inert[0][2] = inert[2][0];
|
---|
| 987 | inert[1][2] = inert[2][1];
|
---|
| 988 |
|
---|
| 989 | // cleanup inert
|
---|
| 990 | for (i=0; i < 3; i++) {
|
---|
| 991 | for (int j=0; j <= i; j++) {
|
---|
| 992 | if (fabs(inert[i][j]) < 1.0e-5) {
|
---|
| 993 | inert[i][j]=inert[j][i]=0.0;
|
---|
| 994 | }
|
---|
| 995 | }
|
---|
| 996 | }
|
---|
| 997 |
|
---|
| 998 | cmat_diag(inert, evals, evecs, 3, 1, 1e-14);
|
---|
| 999 |
|
---|
| 1000 | // cleanup evecs
|
---|
| 1001 | for (i=0; i < 3; i++) {
|
---|
| 1002 | for (int j=0; j < 3; j++) {
|
---|
| 1003 | if (fabs(evecs[i][j]) < 1.0e-5) {
|
---|
| 1004 | evecs[i][j]=0.0;
|
---|
| 1005 | }
|
---|
| 1006 | }
|
---|
| 1007 | }
|
---|
| 1008 |
|
---|
| 1009 | for (i=0; i<natom(); i++) {
|
---|
| 1010 | double a[3];
|
---|
| 1011 | a[0] = r(i,0); a[1] = r(i,1); a[2] = r(i,2);
|
---|
| 1012 | for (j=0; j<3; j++) {
|
---|
| 1013 | double e = 0.0;
|
---|
| 1014 | for (k=0; k<3; k++) {
|
---|
| 1015 | e += a[k] * evecs[k][j];
|
---|
| 1016 | }
|
---|
| 1017 | r_[i][j] = e;
|
---|
| 1018 | }
|
---|
| 1019 | }
|
---|
| 1020 |
|
---|
| 1021 | if (!trans_frame) return;
|
---|
| 1022 |
|
---|
| 1023 | SymmetryOperation tso=point_group()->symm_frame();
|
---|
| 1024 |
|
---|
| 1025 | for (i=0; i < 3; i++) {
|
---|
| 1026 | for (int j=0; j < 3; j++) {
|
---|
| 1027 | double t=0;
|
---|
| 1028 | for (int k=0; k < 3; k++) t += tso[k][j]*evecs[k][i];
|
---|
| 1029 | pg_->symm_frame()[i][j] = t;
|
---|
| 1030 | }
|
---|
| 1031 | }
|
---|
| 1032 | }
|
---|
| 1033 |
|
---|
| 1034 | void
|
---|
| 1035 | Molecule::transform_to_symmetry_frame()
|
---|
| 1036 | {
|
---|
| 1037 | int i,j,k;
|
---|
| 1038 | double t[3][3];
|
---|
| 1039 |
|
---|
| 1040 | SymmetryOperation tso=point_group()->symm_frame();
|
---|
| 1041 |
|
---|
| 1042 | for (i=0; i<3; i++) {
|
---|
| 1043 | for (j=0; j<3; j++) {
|
---|
| 1044 | t[i][j] = tso[i][j];
|
---|
| 1045 | }
|
---|
| 1046 | }
|
---|
| 1047 |
|
---|
| 1048 | for (i=0; i<natom(); i++) {
|
---|
| 1049 | double a[3];
|
---|
| 1050 | a[0] = r(i,0); a[1] = r(i,1); a[2] = r(i,2);
|
---|
| 1051 | for (j=0; j<3; j++) {
|
---|
| 1052 | double e = 0.0;
|
---|
| 1053 | for (k=0; k<3; k++) {
|
---|
| 1054 | e += a[k] * t[k][j];
|
---|
| 1055 | }
|
---|
| 1056 | r_[i][j] = e;
|
---|
| 1057 | }
|
---|
| 1058 | }
|
---|
| 1059 |
|
---|
| 1060 | for (i=0; i<3; i++) {
|
---|
| 1061 | for (j=0; j<3; j++) {
|
---|
| 1062 | double e=0;
|
---|
| 1063 | for (k=0; k<3; k++) e += tso[k][j]*t[k][i];
|
---|
| 1064 | pg_->symm_frame()[i][j] = e;
|
---|
| 1065 | }
|
---|
| 1066 | }
|
---|
| 1067 | }
|
---|
| 1068 |
|
---|
| 1069 | // given a molecule, make sure that equivalent centers have coordinates
|
---|
| 1070 | // that really map into each other
|
---|
| 1071 |
|
---|
| 1072 | void
|
---|
| 1073 | Molecule::cleanup_molecule(double tol)
|
---|
| 1074 | {
|
---|
| 1075 | // if symmetry is c1, do nothing else
|
---|
| 1076 | if (!strcmp(point_group()->symbol(),"c1")) return;
|
---|
| 1077 |
|
---|
| 1078 | int i;
|
---|
| 1079 | SCVector3 up,np,ap;
|
---|
| 1080 | SymmetryOperation so;
|
---|
| 1081 | CharacterTable ct = point_group()->char_table();
|
---|
| 1082 |
|
---|
| 1083 | // first clean up the unique atoms by replacing each coordinate with the
|
---|
| 1084 | // average of coordinates obtained by applying all symmetry operations to
|
---|
| 1085 | // the original atom, iff the new atom ends up near the original atom
|
---|
| 1086 | for (i=0; i < nunique(); i++) {
|
---|
| 1087 | // up will store the original coordinates of unique atom i
|
---|
| 1088 | up = r(unique(i));
|
---|
| 1089 | // ap will hold the average coordinate (times the number of coordinates)
|
---|
| 1090 | // initialize it to the E result
|
---|
| 1091 | ap = up;
|
---|
| 1092 | int ncoor = 1;
|
---|
| 1093 | // loop through all sym ops except E
|
---|
| 1094 | for (int g=1; g < ct.order(); g++) {
|
---|
| 1095 | so = ct.symm_operation(g);
|
---|
| 1096 | for (int ii=0; ii < 3; ii++) {
|
---|
| 1097 | np[ii]=0;
|
---|
| 1098 | for (int jj=0; jj < 3; jj++) np[ii] += so(ii,jj) * up[jj];
|
---|
| 1099 | }
|
---|
| 1100 | if (np.dist(up) < 0.1) {
|
---|
| 1101 | for (int jj=0; jj < 3; jj++) ap[jj] += np[jj];
|
---|
| 1102 | ncoor++;
|
---|
| 1103 | }
|
---|
| 1104 | }
|
---|
| 1105 | // replace the unique coordinate with the average coordinate
|
---|
| 1106 | r_[unique(i)][0] = ap[0] / ncoor;
|
---|
| 1107 | r_[unique(i)][1] = ap[1] / ncoor;
|
---|
| 1108 | r_[unique(i)][2] = ap[2] / ncoor;
|
---|
| 1109 | }
|
---|
| 1110 |
|
---|
| 1111 | // find the atoms equivalent to each unique atom and eliminate
|
---|
| 1112 | // numerical errors that may be in the equivalent atom's coordinates
|
---|
| 1113 |
|
---|
| 1114 | // loop through unique atoms
|
---|
| 1115 | for (i=0; i < nunique(); i++) {
|
---|
| 1116 | // up will store the coordinates of unique atom i
|
---|
| 1117 | up = r(unique(i));
|
---|
| 1118 |
|
---|
| 1119 | // loop through all sym ops except E
|
---|
| 1120 | for (int g=1; g < ct.order(); g++) {
|
---|
| 1121 | so = ct.symm_operation(g);
|
---|
| 1122 | for (int ii=0; ii < 3; ii++) {
|
---|
| 1123 | np[ii]=0;
|
---|
| 1124 | for (int jj=0; jj < 3; jj++) np[ii] += so(ii,jj) * up[jj];
|
---|
| 1125 | }
|
---|
| 1126 |
|
---|
| 1127 | // loop through equivalent atoms
|
---|
| 1128 | int found = 0;
|
---|
| 1129 | for (int j=0; j < natom(); j++) {
|
---|
| 1130 | // see if j is generated from i
|
---|
| 1131 | if (np.dist(SCVector3(r(j))) < tol) {
|
---|
| 1132 | r_[j][0] = np[0];
|
---|
| 1133 | r_[j][1] = np[1];
|
---|
| 1134 | r_[j][2] = np[2];
|
---|
| 1135 | found = 1;
|
---|
| 1136 | }
|
---|
| 1137 | }
|
---|
| 1138 | if (!found) {
|
---|
| 1139 | SCException ex("cleanup: couldn't find atom",
|
---|
| 1140 | __FILE__, __LINE__, class_desc());
|
---|
| 1141 | try {
|
---|
| 1142 | ex.elaborate()
|
---|
| 1143 | << "couldn't find atom at " << np << endl
|
---|
| 1144 | << "transforming uniq atom " << i << " at " << up << endl
|
---|
| 1145 | << "with symmetry op " << g << ":" << endl;
|
---|
| 1146 | so.print(ex.elaborate());
|
---|
| 1147 | }
|
---|
| 1148 | catch (...) {}
|
---|
| 1149 | throw ex;
|
---|
| 1150 | }
|
---|
| 1151 | }
|
---|
| 1152 | }
|
---|
| 1153 |
|
---|
| 1154 | }
|
---|
| 1155 |
|
---|
| 1156 | ///////////////////////////////////////////////////////////////////
|
---|
| 1157 | // Compute the principal axes and the principal moments of inertia
|
---|
| 1158 | ///////////////////////////////////////////////////////////////////
|
---|
| 1159 |
|
---|
| 1160 | void
|
---|
| 1161 | Molecule::principal_moments_of_inertia(double *evals, double **evecs) const
|
---|
| 1162 | {
|
---|
| 1163 |
|
---|
| 1164 | // The principal moments of inertia are computed in amu*angstrom^2
|
---|
| 1165 | // evals: principal moments of inertia
|
---|
| 1166 | // evecs: principal axes (optional argument)
|
---|
| 1167 |
|
---|
| 1168 | Ref<Units> units = new Units("angstroms * angstroms");
|
---|
| 1169 | double au_to_angs = units->from_atomic_units();
|
---|
| 1170 |
|
---|
| 1171 | double *inert[3]; // inertia tensor
|
---|
| 1172 |
|
---|
| 1173 | int i, j;
|
---|
| 1174 | int delete_evecs = 0;
|
---|
| 1175 |
|
---|
| 1176 | // (allocate and) initialize evecs, evals, and inert
|
---|
| 1177 | if (!evecs) {
|
---|
| 1178 | evecs = new double*[3];
|
---|
| 1179 | for (i=0; i<3; i++) evecs[i] = new double[3];
|
---|
| 1180 | delete_evecs = 1;
|
---|
| 1181 | }
|
---|
| 1182 | for (i=0; i<3; i++) {
|
---|
| 1183 | inert[i] = new double[3];
|
---|
| 1184 | memset(inert[i],'\0',sizeof(double)*3);
|
---|
| 1185 | memset(evecs[i],'\0',sizeof(double)*3);
|
---|
| 1186 | }
|
---|
| 1187 | memset(evals,'\0',sizeof(double)*3);
|
---|
| 1188 |
|
---|
| 1189 | SCVector3 com = center_of_mass();
|
---|
| 1190 |
|
---|
| 1191 | // compute inertia tensor
|
---|
| 1192 | SCVector3 ac;
|
---|
| 1193 | for (i=0; i<natom(); i++) {
|
---|
| 1194 | ac = r(i);
|
---|
| 1195 | // compute moments of inertia wrt center of mass
|
---|
| 1196 | for (j=0; j<3; j++) ac(j) -= com(j);
|
---|
| 1197 | double m=au_to_angs*mass(i);
|
---|
| 1198 | inert[0][0] += m * (ac[1]*ac[1] + ac[2]*ac[2]);
|
---|
| 1199 | inert[1][0] -= m * ac[0]*ac[1];
|
---|
| 1200 | inert[1][1] += m * (ac[0]*ac[0] + ac[2]*ac[2]);
|
---|
| 1201 | inert[2][0] -= m * ac[0]*ac[2];
|
---|
| 1202 | inert[2][1] -= m * ac[1]*ac[2];
|
---|
| 1203 | inert[2][2] += m * (ac[0]*ac[0] + ac[1]*ac[1]);
|
---|
| 1204 | }
|
---|
| 1205 | inert[0][1] = inert[1][0];
|
---|
| 1206 | inert[0][2] = inert[2][0];
|
---|
| 1207 | inert[1][2] = inert[2][1];
|
---|
| 1208 |
|
---|
| 1209 | cmat_diag(inert, evals, evecs, 3, 1, 1e-14);
|
---|
| 1210 |
|
---|
| 1211 | if (delete_evecs) {
|
---|
| 1212 | for (i=0; i<3; i++) delete[] evecs[i];
|
---|
| 1213 | delete[] evecs;
|
---|
| 1214 | }
|
---|
| 1215 | for (i=0; i<3; i++) {
|
---|
| 1216 | delete[] inert[i];
|
---|
| 1217 | }
|
---|
| 1218 | }
|
---|
| 1219 |
|
---|
| 1220 | int
|
---|
| 1221 | Molecule::n_core_electrons()
|
---|
| 1222 | {
|
---|
| 1223 | int i,n=0;
|
---|
| 1224 | for (i=0; i<natom(); i++) {
|
---|
| 1225 | if (charge(i) == 0.0) continue;
|
---|
| 1226 | int z = Z_[i];
|
---|
| 1227 | if (z > 2) n += 2;
|
---|
| 1228 | if (z > 10) n += 8;
|
---|
| 1229 | if (z > 18) n += 8;
|
---|
| 1230 | if (z > 30) n += 10;
|
---|
| 1231 | if (z > 36) n += 8;
|
---|
| 1232 | if (z > 48) n += 10;
|
---|
| 1233 | if (z > 54) n += 8;
|
---|
| 1234 | if (z > 72) {
|
---|
| 1235 | throw LimitExceeded<int>("n_core_electrons: atomic number too large",
|
---|
| 1236 | __FILE__, __LINE__, 72, z, class_desc());
|
---|
| 1237 | }
|
---|
| 1238 | }
|
---|
| 1239 | return n;
|
---|
| 1240 | }
|
---|
| 1241 |
|
---|
| 1242 | int
|
---|
| 1243 | Molecule::max_z()
|
---|
| 1244 | {
|
---|
| 1245 | int i, maxz=0;
|
---|
| 1246 | for (i=0; i<natom(); i++) {
|
---|
| 1247 | int z = Z_[i];
|
---|
| 1248 | if (z>maxz) maxz = z;
|
---|
| 1249 | }
|
---|
| 1250 | return maxz;
|
---|
| 1251 | }
|
---|
| 1252 |
|
---|
| 1253 | void
|
---|
| 1254 | Molecule::read_pdb(const char *filename)
|
---|
| 1255 | {
|
---|
| 1256 | clear();
|
---|
| 1257 | ifstream in(filename);
|
---|
| 1258 | Ref<Units> units = new Units("angstrom");
|
---|
| 1259 | while (in.good()) {
|
---|
| 1260 | const int max_line = 80;
|
---|
| 1261 | char line[max_line];
|
---|
| 1262 | in.getline(line,max_line);
|
---|
| 1263 | char *endofline = (char*) memchr(line, 0, max_line);
|
---|
| 1264 | if (endofline) memset(endofline, ' ', &line[max_line-1] - endofline);
|
---|
| 1265 | if (!in.good()) break;
|
---|
| 1266 | if (strncmp(line,"ATOM ",6) == 0
|
---|
| 1267 | ||strncmp(line,"HETATM",6) == 0) {
|
---|
| 1268 |
|
---|
| 1269 | char element[3];
|
---|
| 1270 | strncpy(element,&line[76],2); element[2] = '\0';
|
---|
| 1271 | char name[5];
|
---|
| 1272 | strncpy(name,&line[12],4); name[4] = '\0';
|
---|
| 1273 | if (element[0]==' '&&element[1]==' ') {
|
---|
| 1274 | // no element was given so get the element from the atom name
|
---|
| 1275 | if (name[0]!=' '&&name[3]!=' ') {
|
---|
| 1276 | // some of the atom label may have been
|
---|
| 1277 | // pushed over into the element fields
|
---|
| 1278 | // so check the residue
|
---|
| 1279 | char resName[4];
|
---|
| 1280 | strncpy(resName,&line[17],3); resName[3] = '\0';
|
---|
| 1281 | if (strncmp(line,"ATOM ",6)==0&&(0
|
---|
| 1282 | ||strcmp(resName,"ALA")==0||strcmp(resName,"A ")==0
|
---|
| 1283 | ||strcmp(resName,"ARG")==0||strcmp(resName,"R ")==0
|
---|
| 1284 | ||strcmp(resName,"ASN")==0||strcmp(resName,"N ")==0
|
---|
| 1285 | ||strcmp(resName,"ASP")==0||strcmp(resName,"D ")==0
|
---|
| 1286 | ||strcmp(resName,"ASX")==0||strcmp(resName,"B ")==0
|
---|
| 1287 | ||strcmp(resName,"CYS")==0||strcmp(resName,"C ")==0
|
---|
| 1288 | ||strcmp(resName,"GLN")==0||strcmp(resName,"Q ")==0
|
---|
| 1289 | ||strcmp(resName,"GLU")==0||strcmp(resName,"E ")==0
|
---|
| 1290 | ||strcmp(resName,"GLX")==0||strcmp(resName,"Z ")==0
|
---|
| 1291 | ||strcmp(resName,"GLY")==0||strcmp(resName,"G ")==0
|
---|
| 1292 | ||strcmp(resName,"HIS")==0||strcmp(resName,"H ")==0
|
---|
| 1293 | ||strcmp(resName,"ILE")==0||strcmp(resName,"I ")==0
|
---|
| 1294 | ||strcmp(resName,"LEU")==0||strcmp(resName,"L ")==0
|
---|
| 1295 | ||strcmp(resName,"LYS")==0||strcmp(resName,"K ")==0
|
---|
| 1296 | ||strcmp(resName,"MET")==0||strcmp(resName,"M ")==0
|
---|
| 1297 | ||strcmp(resName,"PHE")==0||strcmp(resName,"F ")==0
|
---|
| 1298 | ||strcmp(resName,"PRO")==0||strcmp(resName,"P ")==0
|
---|
| 1299 | ||strcmp(resName,"SER")==0||strcmp(resName,"S ")==0
|
---|
| 1300 | ||strcmp(resName,"THR")==0||strcmp(resName,"T ")==0
|
---|
| 1301 | ||strcmp(resName,"TRP")==0||strcmp(resName,"W ")==0
|
---|
| 1302 | ||strcmp(resName,"TYR")==0||strcmp(resName,"Y ")==0
|
---|
| 1303 | ||strcmp(resName,"VAL")==0||strcmp(resName,"V ")==0
|
---|
| 1304 | ||strcmp(resName,"A ")==0
|
---|
| 1305 | ||strcmp(resName,"+A ")==0
|
---|
| 1306 | ||strcmp(resName,"C ")==0
|
---|
| 1307 | ||strcmp(resName,"+C ")==0
|
---|
| 1308 | ||strcmp(resName,"G ")==0
|
---|
| 1309 | ||strcmp(resName,"+G ")==0
|
---|
| 1310 | ||strcmp(resName,"I ")==0
|
---|
| 1311 | ||strcmp(resName,"+I ")==0
|
---|
| 1312 | ||strcmp(resName,"T ")==0
|
---|
| 1313 | ||strcmp(resName,"+T ")==0
|
---|
| 1314 | ||strcmp(resName,"U ")==0
|
---|
| 1315 | ||strcmp(resName,"+U ")==0
|
---|
| 1316 | )) {
|
---|
| 1317 | // there no two letter elements for these cases
|
---|
| 1318 | element[0] = name[0];
|
---|
| 1319 | element[1] = '\0';
|
---|
| 1320 | }
|
---|
| 1321 | }
|
---|
| 1322 | else {
|
---|
| 1323 | strncpy(element,name,2); element[2] = '\0';
|
---|
| 1324 | }
|
---|
| 1325 | }
|
---|
| 1326 | if (element[0] == ' ') {
|
---|
| 1327 | element[0] = element[1];
|
---|
| 1328 | element[1] = '\0';
|
---|
| 1329 | }
|
---|
| 1330 |
|
---|
| 1331 | int Z = atominfo_->string_to_Z(element);
|
---|
| 1332 |
|
---|
| 1333 | char field[9];
|
---|
| 1334 | strncpy(field,&line[30],8); field[8] = '\0';
|
---|
| 1335 | double x = atof(field);
|
---|
| 1336 | strncpy(field,&line[38],8); field[8] = '\0';
|
---|
| 1337 | double y = atof(field);
|
---|
| 1338 | strncpy(field,&line[46],8); field[8] = '\0';
|
---|
| 1339 | double z = atof(field);
|
---|
| 1340 | add_atom(Z,
|
---|
| 1341 | x*units->to_atomic_units(),
|
---|
| 1342 | y*units->to_atomic_units(),
|
---|
| 1343 | z*units->to_atomic_units());
|
---|
| 1344 | }
|
---|
| 1345 | else {
|
---|
| 1346 | // skip to next record
|
---|
| 1347 | }
|
---|
| 1348 | }
|
---|
| 1349 | }
|
---|
| 1350 |
|
---|
| 1351 | void
|
---|
| 1352 | Molecule::print_pdb(ostream& os, char *title) const
|
---|
| 1353 | {
|
---|
| 1354 | Ref<Units> u = new Units("angstrom");
|
---|
| 1355 | double bohr = u->from_atomic_units();
|
---|
| 1356 |
|
---|
| 1357 | if (title)
|
---|
| 1358 | os << scprintf("%-10s%-60s\n","COMPND",title);
|
---|
| 1359 | else
|
---|
| 1360 | os << scprintf("%-10s%-60s\n","COMPND","Title");
|
---|
| 1361 |
|
---|
| 1362 | if (title)
|
---|
| 1363 | os << scprintf("REMARK %s\n", title);
|
---|
| 1364 |
|
---|
| 1365 | int i;
|
---|
| 1366 | for (i=0; i < natom(); i++) {
|
---|
| 1367 | char symb[4];
|
---|
| 1368 | std::string symbol(atom_symbol(i));
|
---|
| 1369 | sprintf(symb,"%s1",symbol.c_str());
|
---|
| 1370 |
|
---|
| 1371 | os << scprintf(
|
---|
| 1372 | "HETATM%5d %-3s UNK %5d %8.3f%8.3f%8.3f 0.00 0.00 0\n",
|
---|
| 1373 | i+1, symb, 0, r(i,0)*bohr, r(i,1)*bohr, r(i,2)*bohr);
|
---|
| 1374 | }
|
---|
| 1375 |
|
---|
| 1376 | for (i=0; i < natom(); i++) {
|
---|
| 1377 | double at_rad_i = atominfo_->atomic_radius(Z_[i]);
|
---|
| 1378 | SCVector3 ai(r(i));
|
---|
| 1379 |
|
---|
| 1380 | os << scprintf("CONECT%5d",i+1);
|
---|
| 1381 |
|
---|
| 1382 | for (int j=0; j < natom(); j++) {
|
---|
| 1383 |
|
---|
| 1384 | if (j==i) continue;
|
---|
| 1385 |
|
---|
| 1386 | double at_rad_j = atominfo_->atomic_radius(Z_[j]);
|
---|
| 1387 | SCVector3 aj(r(j));
|
---|
| 1388 |
|
---|
| 1389 | if (ai.dist(aj) < 1.1*(at_rad_i+at_rad_j))
|
---|
| 1390 | os << scprintf("%5d",j+1);
|
---|
| 1391 | }
|
---|
| 1392 |
|
---|
| 1393 | os << endl;
|
---|
| 1394 | }
|
---|
| 1395 |
|
---|
| 1396 | os << "END" << endl;
|
---|
| 1397 | os.flush();
|
---|
| 1398 | }
|
---|
| 1399 |
|
---|
| 1400 | double
|
---|
| 1401 | Molecule::mass(int atom) const
|
---|
| 1402 | {
|
---|
| 1403 | if (!mass_ || mass_[atom] == 0) {
|
---|
| 1404 | return atominfo_->mass(Z_[atom]);
|
---|
| 1405 | }
|
---|
| 1406 | return mass_[atom];
|
---|
| 1407 | }
|
---|
| 1408 |
|
---|
| 1409 | const char *
|
---|
| 1410 | Molecule::label(int atom) const
|
---|
| 1411 | {
|
---|
| 1412 | if (!labels_) return 0;
|
---|
| 1413 | return labels_[atom];
|
---|
| 1414 | }
|
---|
| 1415 |
|
---|
| 1416 | std::string
|
---|
| 1417 | Molecule::atom_name(int iatom) const
|
---|
| 1418 | {
|
---|
| 1419 | return atominfo_->name(Z_[iatom]);
|
---|
| 1420 | }
|
---|
| 1421 |
|
---|
| 1422 | std::string
|
---|
| 1423 | Molecule::atom_symbol(int iatom) const
|
---|
| 1424 | {
|
---|
| 1425 | return atominfo_->symbol(Z_[iatom]);
|
---|
| 1426 | }
|
---|
| 1427 |
|
---|
| 1428 | /////////////////////////////////////////////////////////////////////////////
|
---|
| 1429 |
|
---|
| 1430 | // Local Variables:
|
---|
| 1431 | // mode: c++
|
---|
| 1432 | // c-file-style: "CLJ"
|
---|
| 1433 | // End:
|
---|