[bcf653] | 1 | /*
|
---|
| 2 | * Project: MoleCuilder
|
---|
| 3 | * Description: creates and alters molecular systems
|
---|
| 4 | * Copyright (C) 2010 University of Bonn. All rights reserved.
|
---|
| 5 | * Please see the LICENSE file or "Copyright notice" in builder.cpp for details.
|
---|
| 6 | */
|
---|
| 7 |
|
---|
[e138de] | 8 | /** \file MoleculeListClass.cpp
|
---|
| 9 | *
|
---|
| 10 | * Function implementations for the class MoleculeListClass.
|
---|
| 11 | *
|
---|
| 12 | */
|
---|
| 13 |
|
---|
[bf3817] | 14 | // include config.h
|
---|
[aafd77] | 15 | #ifdef HAVE_CONFIG_H
|
---|
| 16 | #include <config.h>
|
---|
| 17 | #endif
|
---|
| 18 |
|
---|
[ad011c] | 19 | #include "CodePatterns/MemDebug.hpp"
|
---|
[112b09] | 20 |
|
---|
[49e1ae] | 21 | #include <cstring>
|
---|
| 22 |
|
---|
[aafd77] | 23 | #include <gsl/gsl_inline.h>
|
---|
| 24 | #include <gsl/gsl_heapsort.h>
|
---|
| 25 |
|
---|
[e138de] | 26 | #include "atom.hpp"
|
---|
[129204] | 27 | #include "Bond/bond.hpp"
|
---|
[e138de] | 28 | #include "boundary.hpp"
|
---|
[583081] | 29 | #include "Box.hpp"
|
---|
| 30 | #include "CodePatterns/Assert.hpp"
|
---|
| 31 | #include "CodePatterns/Log.hpp"
|
---|
| 32 | #include "CodePatterns/Verbose.hpp"
|
---|
[e138de] | 33 | #include "config.hpp"
|
---|
| 34 | #include "element.hpp"
|
---|
[129204] | 35 | #include "Graph/BondGraph.hpp"
|
---|
[952f38] | 36 | #include "Helpers/helpers.hpp"
|
---|
[583081] | 37 | #include "LinearAlgebra/RealSpaceMatrix.hpp"
|
---|
[e138de] | 38 | #include "linkedcell.hpp"
|
---|
| 39 | #include "molecule.hpp"
|
---|
| 40 | #include "periodentafel.hpp"
|
---|
[88b400] | 41 | #include "tesselation.hpp"
|
---|
[583081] | 42 | #include "World.hpp"
|
---|
| 43 | #include "WorldTime.hpp"
|
---|
[920c70] | 44 |
|
---|
[e138de] | 45 | /*********************************** Functions for class MoleculeListClass *************************/
|
---|
| 46 |
|
---|
| 47 | /** Constructor for MoleculeListClass.
|
---|
| 48 | */
|
---|
[cbc5fb] | 49 | MoleculeListClass::MoleculeListClass(World *_world) :
|
---|
[cd5047] | 50 | Observable("MoleculeListClass"),
|
---|
[81a9bc] | 51 | MaxIndex(1),
|
---|
| 52 | world(_world)
|
---|
[97b825] | 53 | {};
|
---|
[e138de] | 54 |
|
---|
| 55 | /** Destructor for MoleculeListClass.
|
---|
| 56 | */
|
---|
| 57 | MoleculeListClass::~MoleculeListClass()
|
---|
| 58 | {
|
---|
[bd6bfa] | 59 | DoLog(4) && (Log() << Verbose(4) << "Clearing ListOfMolecules." << endl);
|
---|
| 60 | for(MoleculeList::iterator MolRunner = ListOfMolecules.begin(); MolRunner != ListOfMolecules.end(); ++MolRunner)
|
---|
| 61 | (*MolRunner)->signOff(this);
|
---|
[e138de] | 62 | ListOfMolecules.clear(); // empty list
|
---|
| 63 | };
|
---|
| 64 |
|
---|
| 65 | /** Insert a new molecule into the list and set its number.
|
---|
| 66 | * \param *mol molecule to add to list.
|
---|
| 67 | */
|
---|
| 68 | void MoleculeListClass::insert(molecule *mol)
|
---|
| 69 | {
|
---|
[2ba827] | 70 | OBSERVE;
|
---|
[e138de] | 71 | mol->IndexNr = MaxIndex++;
|
---|
| 72 | ListOfMolecules.push_back(mol);
|
---|
[520c8b] | 73 | mol->signOn(this);
|
---|
[e138de] | 74 | };
|
---|
| 75 |
|
---|
[bd6bfa] | 76 | /** Erases a molecule from the list.
|
---|
| 77 | * \param *mol molecule to add to list.
|
---|
| 78 | */
|
---|
| 79 | void MoleculeListClass::erase(molecule *mol)
|
---|
| 80 | {
|
---|
| 81 | OBSERVE;
|
---|
| 82 | mol->signOff(this);
|
---|
| 83 | ListOfMolecules.remove(mol);
|
---|
| 84 | };
|
---|
| 85 |
|
---|
[a0064e] | 86 | /** Comparison function for two values.
|
---|
| 87 | * \param *a
|
---|
| 88 | * \param *b
|
---|
| 89 | * \return <0, \a *a less than \a *b, ==0 if equal, >0 \a *a greater than \a *b
|
---|
| 90 | */
|
---|
| 91 | int CompareDoubles (const void * a, const void * b)
|
---|
| 92 | {
|
---|
| 93 | if (*(double *)a > *(double *)b)
|
---|
| 94 | return -1;
|
---|
| 95 | else if (*(double *)a < *(double *)b)
|
---|
| 96 | return 1;
|
---|
| 97 | else
|
---|
| 98 | return 0;
|
---|
| 99 | };
|
---|
| 100 |
|
---|
| 101 |
|
---|
[e138de] | 102 | /** Compare whether two molecules are equal.
|
---|
| 103 | * \param *a molecule one
|
---|
| 104 | * \param *n molecule two
|
---|
| 105 | * \return lexical value (-1, 0, +1)
|
---|
| 106 | */
|
---|
| 107 | int MolCompare(const void *a, const void *b)
|
---|
| 108 | {
|
---|
| 109 | int *aList = NULL, *bList = NULL;
|
---|
| 110 | int Count, Counter, aCounter, bCounter;
|
---|
| 111 | int flag;
|
---|
| 112 |
|
---|
| 113 | // sort each atom list and put the numbers into a list, then go through
|
---|
| 114 | //Log() << Verbose(0) << "Comparing fragment no. " << *(molecule **)a << " to " << *(molecule **)b << "." << endl;
|
---|
[ea7176] | 115 | // Yes those types are awkward... but check it for yourself it checks out this way
|
---|
| 116 | molecule *const *mol1_ptr= static_cast<molecule *const *>(a);
|
---|
| 117 | molecule *mol1 = *mol1_ptr;
|
---|
| 118 | molecule *const *mol2_ptr= static_cast<molecule *const *>(b);
|
---|
| 119 | molecule *mol2 = *mol2_ptr;
|
---|
| 120 | if (mol1->getAtomCount() < mol2->getAtomCount()) {
|
---|
[e138de] | 121 | return -1;
|
---|
| 122 | } else {
|
---|
[ea7176] | 123 | if (mol1->getAtomCount() > mol2->getAtomCount())
|
---|
[e138de] | 124 | return +1;
|
---|
| 125 | else {
|
---|
[ea7176] | 126 | Count = mol1->getAtomCount();
|
---|
[e138de] | 127 | aList = new int[Count];
|
---|
| 128 | bList = new int[Count];
|
---|
| 129 |
|
---|
| 130 | // fill the lists
|
---|
| 131 | Counter = 0;
|
---|
| 132 | aCounter = 0;
|
---|
| 133 | bCounter = 0;
|
---|
[ea7176] | 134 | molecule::const_iterator aiter = mol1->begin();
|
---|
| 135 | molecule::const_iterator biter = mol2->begin();
|
---|
| 136 | for (;(aiter != mol1->end()) && (biter != mol2->end());
|
---|
[9879f6] | 137 | ++aiter, ++biter) {
|
---|
| 138 | if ((*aiter)->GetTrueFather() == NULL)
|
---|
[e138de] | 139 | aList[Counter] = Count + (aCounter++);
|
---|
| 140 | else
|
---|
[735b1c] | 141 | aList[Counter] = (*aiter)->GetTrueFather()->getNr();
|
---|
[9879f6] | 142 | if ((*biter)->GetTrueFather() == NULL)
|
---|
[e138de] | 143 | bList[Counter] = Count + (bCounter++);
|
---|
| 144 | else
|
---|
[735b1c] | 145 | bList[Counter] = (*biter)->GetTrueFather()->getNr();
|
---|
[e138de] | 146 | Counter++;
|
---|
| 147 | }
|
---|
| 148 | // check if AtomCount was for real
|
---|
| 149 | flag = 0;
|
---|
[ea7176] | 150 | if ((aiter == mol1->end()) && (biter != mol2->end())) {
|
---|
[e138de] | 151 | flag = -1;
|
---|
| 152 | } else {
|
---|
[ea7176] | 153 | if ((aiter != mol1->end()) && (biter == mol2->end()))
|
---|
[e138de] | 154 | flag = 1;
|
---|
| 155 | }
|
---|
| 156 | if (flag == 0) {
|
---|
| 157 | // sort the lists
|
---|
| 158 | gsl_heapsort(aList, Count, sizeof(int), CompareDoubles);
|
---|
| 159 | gsl_heapsort(bList, Count, sizeof(int), CompareDoubles);
|
---|
| 160 | // compare the lists
|
---|
| 161 |
|
---|
| 162 | flag = 0;
|
---|
| 163 | for (int i = 0; i < Count; i++) {
|
---|
| 164 | if (aList[i] < bList[i]) {
|
---|
| 165 | flag = -1;
|
---|
| 166 | } else {
|
---|
| 167 | if (aList[i] > bList[i])
|
---|
| 168 | flag = 1;
|
---|
| 169 | }
|
---|
| 170 | if (flag != 0)
|
---|
| 171 | break;
|
---|
| 172 | }
|
---|
| 173 | }
|
---|
| 174 | delete[] (aList);
|
---|
| 175 | delete[] (bList);
|
---|
| 176 | return flag;
|
---|
| 177 | }
|
---|
| 178 | }
|
---|
| 179 | return -1;
|
---|
| 180 | };
|
---|
| 181 |
|
---|
| 182 | /** Output of a list of all molecules.
|
---|
| 183 | * \param *out output stream
|
---|
| 184 | */
|
---|
[24a5e0] | 185 | void MoleculeListClass::Enumerate(ostream *out)
|
---|
[e138de] | 186 | {
|
---|
[ead4e6] | 187 | periodentafel *periode = World::getInstance().getPeriode();
|
---|
| 188 | std::map<atomicNumber_t,unsigned int> counts;
|
---|
[e138de] | 189 | double size=0;
|
---|
| 190 | Vector Origin;
|
---|
| 191 |
|
---|
| 192 | // header
|
---|
[835a0f] | 193 | (*out) << "Index\tName\t\tAtoms\tFormula\tCenter\tSize" << endl;
|
---|
| 194 | (*out) << "-----------------------------------------------" << endl;
|
---|
[e138de] | 195 | if (ListOfMolecules.size() == 0)
|
---|
[835a0f] | 196 | (*out) << "\tNone" << endl;
|
---|
[e138de] | 197 | else {
|
---|
| 198 | Origin.Zero();
|
---|
| 199 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++) {
|
---|
| 200 | // count atoms per element and determine size of bounding sphere
|
---|
| 201 | size=0.;
|
---|
[9879f6] | 202 | for (molecule::const_iterator iter = (*ListRunner)->begin(); iter != (*ListRunner)->end(); ++iter) {
|
---|
[d74077] | 203 | counts[(*iter)->getType()->getNumber()]++;
|
---|
| 204 | if ((*iter)->DistanceSquared(Origin) > size)
|
---|
| 205 | size = (*iter)->DistanceSquared(Origin);
|
---|
[e138de] | 206 | }
|
---|
| 207 | // output Index, Name, number of atoms, chemical formula
|
---|
[ea7176] | 208 | (*out) << ((*ListRunner)->ActiveFlag ? "*" : " ") << (*ListRunner)->IndexNr << "\t" << (*ListRunner)->name << "\t\t" << (*ListRunner)->getAtomCount() << "\t";
|
---|
[ead4e6] | 209 |
|
---|
| 210 | std::map<atomicNumber_t,unsigned int>::reverse_iterator iter;
|
---|
| 211 | for(iter=counts.rbegin(); iter!=counts.rend();++iter){
|
---|
| 212 | atomicNumber_t Z =(*iter).first;
|
---|
| 213 | (*out) << periode->FindElement(Z)->getSymbol() << (*iter).second;
|
---|
[e138de] | 214 | }
|
---|
| 215 | // Center and size
|
---|
[1883f9] | 216 | Vector *Center = (*ListRunner)->DetermineCenterOfAll();
|
---|
| 217 | (*out) << "\t" << *Center << "\t" << sqrt(size) << endl;
|
---|
| 218 | delete(Center);
|
---|
[e138de] | 219 | }
|
---|
| 220 | }
|
---|
| 221 | };
|
---|
| 222 |
|
---|
| 223 | /** Returns the molecule with the given index \a index.
|
---|
| 224 | * \param index index of the desired molecule
|
---|
[1907a7] | 225 | * \return pointer to molecule structure, NULL if not found
|
---|
[e138de] | 226 | */
|
---|
| 227 | molecule * MoleculeListClass::ReturnIndex(int index)
|
---|
| 228 | {
|
---|
| 229 | for(MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++)
|
---|
| 230 | if ((*ListRunner)->IndexNr == index)
|
---|
| 231 | return (*ListRunner);
|
---|
| 232 | return NULL;
|
---|
| 233 | };
|
---|
| 234 |
|
---|
| 235 |
|
---|
| 236 | /** Simple output of the pointers in ListOfMolecules.
|
---|
| 237 | * \param *out output stream
|
---|
| 238 | */
|
---|
| 239 | void MoleculeListClass::Output(ofstream *out)
|
---|
| 240 | {
|
---|
[a67d19] | 241 | DoLog(1) && (Log() << Verbose(1) << "MoleculeList: ");
|
---|
[e138de] | 242 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++)
|
---|
[a67d19] | 243 | DoLog(0) && (Log() << Verbose(0) << *ListRunner << "\t");
|
---|
| 244 | DoLog(0) && (Log() << Verbose(0) << endl);
|
---|
[e138de] | 245 | };
|
---|
| 246 |
|
---|
| 247 | /** Calculates necessary hydrogen correction due to unwanted interaction between saturated ones.
|
---|
| 248 | * If for a pair of two hydrogen atoms a and b, at least is a saturated one, and a and b are not
|
---|
| 249 | * bonded to the same atom, then we add for this pair a correction term constructed from a Morse
|
---|
| 250 | * potential function fit to QM calculations with respecting to the interatomic hydrogen distance.
|
---|
[35b698] | 251 | * \param &path path to file
|
---|
[e138de] | 252 | */
|
---|
[35b698] | 253 | bool MoleculeListClass::AddHydrogenCorrection(std::string &path)
|
---|
[e138de] | 254 | {
|
---|
| 255 | bond *Binder = NULL;
|
---|
| 256 | double ***FitConstant = NULL, **correction = NULL;
|
---|
| 257 | int a, b;
|
---|
| 258 | ofstream output;
|
---|
| 259 | ifstream input;
|
---|
| 260 | string line;
|
---|
| 261 | stringstream zeile;
|
---|
| 262 | double distance;
|
---|
| 263 | char ParsedLine[1023];
|
---|
| 264 | double tmp;
|
---|
| 265 | char *FragmentNumber = NULL;
|
---|
| 266 |
|
---|
[a67d19] | 267 | DoLog(1) && (Log() << Verbose(1) << "Saving hydrogen saturation correction ... ");
|
---|
[e138de] | 268 | // 0. parse in fit constant files that should have the same dimension as the final energy files
|
---|
| 269 | // 0a. find dimension of matrices with constants
|
---|
| 270 | line = path;
|
---|
| 271 | line += "1";
|
---|
| 272 | line += FITCONSTANTSUFFIX;
|
---|
| 273 | input.open(line.c_str());
|
---|
[35b698] | 274 | if (input.fail()) {
|
---|
[a67d19] | 275 | DoLog(1) && (Log() << Verbose(1) << endl << "Unable to open " << line << ", is the directory correct?" << endl);
|
---|
[e138de] | 276 | return false;
|
---|
| 277 | }
|
---|
| 278 | a = 0;
|
---|
| 279 | b = -1; // we overcount by one
|
---|
| 280 | while (!input.eof()) {
|
---|
| 281 | input.getline(ParsedLine, 1023);
|
---|
| 282 | zeile.str(ParsedLine);
|
---|
| 283 | int i = 0;
|
---|
| 284 | while (!zeile.eof()) {
|
---|
| 285 | zeile >> distance;
|
---|
| 286 | i++;
|
---|
| 287 | }
|
---|
| 288 | if (i > a)
|
---|
| 289 | a = i;
|
---|
| 290 | b++;
|
---|
| 291 | }
|
---|
[a67d19] | 292 | DoLog(0) && (Log() << Verbose(0) << "I recognized " << a << " columns and " << b << " rows, ");
|
---|
[e138de] | 293 | input.close();
|
---|
| 294 |
|
---|
| 295 | // 0b. allocate memory for constants
|
---|
[920c70] | 296 | FitConstant = new double**[3];
|
---|
[e138de] | 297 | for (int k = 0; k < 3; k++) {
|
---|
[920c70] | 298 | FitConstant[k] = new double*[a];
|
---|
[e138de] | 299 | for (int i = a; i--;) {
|
---|
[920c70] | 300 | FitConstant[k][i] = new double[b];
|
---|
| 301 | for (int j = b; j--;) {
|
---|
| 302 | FitConstant[k][i][j] = 0.;
|
---|
| 303 | }
|
---|
[e138de] | 304 | }
|
---|
| 305 | }
|
---|
| 306 | // 0c. parse in constants
|
---|
| 307 | for (int i = 0; i < 3; i++) {
|
---|
| 308 | line = path;
|
---|
| 309 | line.append("/");
|
---|
| 310 | line += FRAGMENTPREFIX;
|
---|
| 311 | sprintf(ParsedLine, "%d", i + 1);
|
---|
| 312 | line += ParsedLine;
|
---|
| 313 | line += FITCONSTANTSUFFIX;
|
---|
| 314 | input.open(line.c_str());
|
---|
| 315 | if (input == NULL) {
|
---|
[58ed4a] | 316 | DoeLog(0) && (eLog()<< Verbose(0) << endl << "Unable to open " << line << ", is the directory correct?" << endl);
|
---|
[e359a8] | 317 | performCriticalExit();
|
---|
[e138de] | 318 | return false;
|
---|
| 319 | }
|
---|
| 320 | int k = 0, l;
|
---|
| 321 | while ((!input.eof()) && (k < b)) {
|
---|
| 322 | input.getline(ParsedLine, 1023);
|
---|
| 323 | //Log() << Verbose(0) << "Current Line: " << ParsedLine << endl;
|
---|
| 324 | zeile.str(ParsedLine);
|
---|
| 325 | zeile.clear();
|
---|
| 326 | l = 0;
|
---|
| 327 | while ((!zeile.eof()) && (l < a)) {
|
---|
| 328 | zeile >> FitConstant[i][l][k];
|
---|
| 329 | //Log() << Verbose(0) << FitConstant[i][l][k] << "\t";
|
---|
| 330 | l++;
|
---|
| 331 | }
|
---|
| 332 | //Log() << Verbose(0) << endl;
|
---|
| 333 | k++;
|
---|
| 334 | }
|
---|
| 335 | input.close();
|
---|
| 336 | }
|
---|
| 337 | for (int k = 0; k < 3; k++) {
|
---|
[a67d19] | 338 | DoLog(0) && (Log() << Verbose(0) << "Constants " << k << ":" << endl);
|
---|
[e138de] | 339 | for (int j = 0; j < b; j++) {
|
---|
| 340 | for (int i = 0; i < a; i++) {
|
---|
[a67d19] | 341 | DoLog(0) && (Log() << Verbose(0) << FitConstant[k][i][j] << "\t");
|
---|
[e138de] | 342 | }
|
---|
[a67d19] | 343 | DoLog(0) && (Log() << Verbose(0) << endl);
|
---|
[e138de] | 344 | }
|
---|
[a67d19] | 345 | DoLog(0) && (Log() << Verbose(0) << endl);
|
---|
[e138de] | 346 | }
|
---|
| 347 |
|
---|
| 348 | // 0d. allocate final correction matrix
|
---|
[920c70] | 349 | correction = new double*[a];
|
---|
[e138de] | 350 | for (int i = a; i--;)
|
---|
[920c70] | 351 | correction[i] = new double[b];
|
---|
[e138de] | 352 |
|
---|
| 353 | // 1a. go through every molecule in the list
|
---|
| 354 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++) {
|
---|
| 355 | // 1b. zero final correction matrix
|
---|
| 356 | for (int k = a; k--;)
|
---|
| 357 | for (int j = b; j--;)
|
---|
| 358 | correction[k][j] = 0.;
|
---|
| 359 | // 2. take every hydrogen that is a saturated one
|
---|
[9879f6] | 360 | for (molecule::const_iterator iter = (*ListRunner)->begin(); iter != (*ListRunner)->end(); ++iter) {
|
---|
[9d83b6] | 361 | //Log() << Verbose(1) << "(*iter): " << *(*iter) << " with first bond " << *((*iter)->getListOfBonds().begin()) << "." << endl;
|
---|
[83f176] | 362 | if (((*iter)->getType()->getAtomicNumber() == 1) && (((*iter)->father == NULL)
|
---|
| 363 | || ((*iter)->father->getType()->getAtomicNumber() != 1))) { // if it's a hydrogen
|
---|
[9879f6] | 364 | for (molecule::const_iterator runner = (*ListRunner)->begin(); runner != (*ListRunner)->end(); ++runner) {
|
---|
[9d83b6] | 365 | //Log() << Verbose(2) << "Runner: " << *(*runner) << " with first bond " << *((*iter)->getListOfBonds().begin()) << "." << endl;
|
---|
[e138de] | 366 | // 3. take every other hydrogen that is the not the first and not bound to same bonding partner
|
---|
[9d83b6] | 367 | Binder = *((*runner)->getListOfBonds().begin());
|
---|
[735b1c] | 368 | if (((*runner)->getType()->getAtomicNumber() == 1) && ((*runner)->getNr() > (*iter)->getNr()) && (Binder->GetOtherAtom((*runner)) != Binder->GetOtherAtom((*iter)))) { // (hydrogens have only one bonding partner!)
|
---|
[e138de] | 369 | // 4. evaluate the morse potential for each matrix component and add up
|
---|
[d74077] | 370 | distance = (*runner)->distance(*(*iter));
|
---|
[9879f6] | 371 | //Log() << Verbose(0) << "Fragment " << (*ListRunner)->name << ": " << *(*runner) << "<= " << distance << "=>" << *(*iter) << ":" << endl;
|
---|
[e138de] | 372 | for (int k = 0; k < a; k++) {
|
---|
| 373 | for (int j = 0; j < b; j++) {
|
---|
| 374 | switch (k) {
|
---|
| 375 | case 1:
|
---|
| 376 | case 7:
|
---|
| 377 | case 11:
|
---|
| 378 | tmp = pow(FitConstant[0][k][j] * (1. - exp(-FitConstant[1][k][j] * (distance - FitConstant[2][k][j]))), 2);
|
---|
| 379 | break;
|
---|
| 380 | default:
|
---|
| 381 | tmp = FitConstant[0][k][j] * pow(distance, FitConstant[1][k][j]) + FitConstant[2][k][j];
|
---|
| 382 | };
|
---|
| 383 | correction[k][j] -= tmp; // ground state is actually lower (disturbed by additional interaction)
|
---|
| 384 | //Log() << Verbose(0) << tmp << "\t";
|
---|
| 385 | }
|
---|
| 386 | //Log() << Verbose(0) << endl;
|
---|
| 387 | }
|
---|
| 388 | //Log() << Verbose(0) << endl;
|
---|
| 389 | }
|
---|
| 390 | }
|
---|
| 391 | }
|
---|
| 392 | }
|
---|
| 393 | // 5. write final matrix to file
|
---|
| 394 | line = path;
|
---|
| 395 | line.append("/");
|
---|
| 396 | line += FRAGMENTPREFIX;
|
---|
| 397 | FragmentNumber = FixedDigitNumber(ListOfMolecules.size(), (*ListRunner)->IndexNr);
|
---|
| 398 | line += FragmentNumber;
|
---|
[920c70] | 399 | delete[] (FragmentNumber);
|
---|
[e138de] | 400 | line += HCORRECTIONSUFFIX;
|
---|
| 401 | output.open(line.c_str());
|
---|
| 402 | output << "Time\t\tTotal\t\tKinetic\t\tNonLocal\tCorrelation\tExchange\tPseudo\t\tHartree\t\t-Gauss\t\tEwald\t\tIonKin\t\tETotal" << endl;
|
---|
| 403 | for (int j = 0; j < b; j++) {
|
---|
| 404 | for (int i = 0; i < a; i++)
|
---|
| 405 | output << correction[i][j] << "\t";
|
---|
| 406 | output << endl;
|
---|
| 407 | }
|
---|
| 408 | output.close();
|
---|
| 409 | }
|
---|
[920c70] | 410 | for (int i = a; i--;)
|
---|
| 411 | delete[](correction[i]);
|
---|
| 412 | delete[](correction);
|
---|
| 413 |
|
---|
[e138de] | 414 | line = path;
|
---|
| 415 | line.append("/");
|
---|
| 416 | line += HCORRECTIONSUFFIX;
|
---|
| 417 | output.open(line.c_str());
|
---|
| 418 | output << "Time\t\tTotal\t\tKinetic\t\tNonLocal\tCorrelation\tExchange\tPseudo\t\tHartree\t\t-Gauss\t\tEwald\t\tIonKin\t\tETotal" << endl;
|
---|
| 419 | for (int j = 0; j < b; j++) {
|
---|
| 420 | for (int i = 0; i < a; i++)
|
---|
| 421 | output << 0 << "\t";
|
---|
| 422 | output << endl;
|
---|
| 423 | }
|
---|
| 424 | output.close();
|
---|
| 425 | // 6. free memory of parsed matrices
|
---|
| 426 | for (int k = 0; k < 3; k++) {
|
---|
| 427 | for (int i = a; i--;) {
|
---|
[920c70] | 428 | delete[](FitConstant[k][i]);
|
---|
[e138de] | 429 | }
|
---|
[920c70] | 430 | delete[](FitConstant[k]);
|
---|
[e138de] | 431 | }
|
---|
[920c70] | 432 | delete[](FitConstant);
|
---|
[a67d19] | 433 | DoLog(0) && (Log() << Verbose(0) << "done." << endl);
|
---|
[e138de] | 434 | return true;
|
---|
| 435 | };
|
---|
| 436 |
|
---|
| 437 | /** Store force indices, i.e. the connection between the nuclear index in the total molecule config and the respective atom in fragment config.
|
---|
[35b698] | 438 | * \param &path path to file
|
---|
[e138de] | 439 | * \param *SortIndex Index to map from the BFS labeling to the sequence how of Ion_Type in the config
|
---|
| 440 | * \return true - file written successfully, false - writing failed
|
---|
| 441 | */
|
---|
[35b698] | 442 | bool MoleculeListClass::StoreForcesFile(std::string &path, int *SortIndex)
|
---|
[e138de] | 443 | {
|
---|
| 444 | bool status = true;
|
---|
[35b698] | 445 | string filename(path);
|
---|
| 446 | filename += FORCESFILE;
|
---|
| 447 | ofstream ForcesFile(filename.c_str());
|
---|
[ead4e6] | 448 | periodentafel *periode=World::getInstance().getPeriode();
|
---|
[e138de] | 449 |
|
---|
| 450 | // open file for the force factors
|
---|
[a67d19] | 451 | DoLog(1) && (Log() << Verbose(1) << "Saving force factors ... ");
|
---|
[35b698] | 452 | if (!ForcesFile.fail()) {
|
---|
[e138de] | 453 | //Log() << Verbose(1) << "Final AtomicForcesList: ";
|
---|
| 454 | //output << prefix << "Forces" << endl;
|
---|
| 455 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++) {
|
---|
[ead4e6] | 456 | periodentafel::const_iterator elemIter;
|
---|
| 457 | for(elemIter=periode->begin();elemIter!=periode->end();++elemIter){
|
---|
[389cc8] | 458 | if ((*ListRunner)->hasElement((*elemIter).first)) { // if this element got atoms
|
---|
[a7b761b] | 459 | for(molecule::iterator atomIter = (*ListRunner)->begin(); atomIter !=(*ListRunner)->end();++atomIter){
|
---|
[d74077] | 460 | if ((*atomIter)->getType()->getNumber() == (*elemIter).first) {
|
---|
[a7b761b] | 461 | if (((*atomIter)->GetTrueFather() != NULL) && ((*atomIter)->GetTrueFather() != (*atomIter))) {// if there is a rea
|
---|
[e138de] | 462 | //Log() << Verbose(0) << "Walker is " << *Walker << " with true father " << *( Walker->GetTrueFather()) << ", it
|
---|
[735b1c] | 463 | ForcesFile << SortIndex[(*atomIter)->GetTrueFather()->getNr()] << "\t";
|
---|
[e138de] | 464 | } else
|
---|
| 465 | // otherwise a -1 to indicate an added saturation hydrogen
|
---|
| 466 | ForcesFile << "-1\t";
|
---|
| 467 | }
|
---|
| 468 | }
|
---|
| 469 | }
|
---|
| 470 | }
|
---|
| 471 | ForcesFile << endl;
|
---|
| 472 | }
|
---|
| 473 | ForcesFile.close();
|
---|
[a67d19] | 474 | DoLog(1) && (Log() << Verbose(1) << "done." << endl);
|
---|
[e138de] | 475 | } else {
|
---|
| 476 | status = false;
|
---|
[35b698] | 477 | DoLog(1) && (Log() << Verbose(1) << "failed to open file " << filename << "." << endl);
|
---|
[e138de] | 478 | }
|
---|
| 479 | ForcesFile.close();
|
---|
| 480 |
|
---|
| 481 | return status;
|
---|
| 482 | };
|
---|
| 483 |
|
---|
| 484 | /** Writes a config file for each molecule in the given \a **FragmentList.
|
---|
| 485 | * \param *out output stream for debugging
|
---|
[35b698] | 486 | * \param &prefix path and prefix to the fragment config files
|
---|
[e138de] | 487 | * \param *SortIndex Index to map from the BFS labeling to the sequence how of Ion_Type in the config
|
---|
| 488 | * \return true - success (each file was written), false - something went wrong.
|
---|
| 489 | */
|
---|
[35b698] | 490 | bool MoleculeListClass::OutputConfigForListOfFragments(std::string &prefix, int *SortIndex)
|
---|
[e138de] | 491 | {
|
---|
| 492 | ofstream outputFragment;
|
---|
[35b698] | 493 | std::string FragmentName;
|
---|
[e138de] | 494 | char PathBackup[MAXSTRINGSIZE];
|
---|
| 495 | bool result = true;
|
---|
| 496 | bool intermediateResult = true;
|
---|
| 497 | Vector BoxDimension;
|
---|
| 498 | char *FragmentNumber = NULL;
|
---|
| 499 | char *path = NULL;
|
---|
| 500 | int FragmentCounter = 0;
|
---|
| 501 | ofstream output;
|
---|
[cca9ef] | 502 | RealSpaceMatrix cell_size = World::getInstance().getDomain().getM();
|
---|
| 503 | RealSpaceMatrix cell_size_backup = cell_size;
|
---|
[3c58f8] | 504 | int count=0;
|
---|
[e138de] | 505 |
|
---|
| 506 | // store the fragments as config and as xyz
|
---|
| 507 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++) {
|
---|
| 508 | // save default path as it is changed for each fragment
|
---|
[35b698] | 509 | path = World::getInstance().getConfig()->GetDefaultPath();
|
---|
[e138de] | 510 | if (path != NULL)
|
---|
| 511 | strcpy(PathBackup, path);
|
---|
[e359a8] | 512 | else {
|
---|
[efe516] | 513 | ELOG(0, "OutputConfigForListOfFragments: NULL default path obtained from config!");
|
---|
[e359a8] | 514 | performCriticalExit();
|
---|
| 515 | }
|
---|
[e138de] | 516 |
|
---|
| 517 | // correct periodic
|
---|
[3c58f8] | 518 | if ((*ListRunner)->ScanForPeriodicCorrection()) {
|
---|
| 519 | count++;
|
---|
| 520 | }
|
---|
[e138de] | 521 |
|
---|
[efe516] | 522 | {
|
---|
| 523 | // list atoms in fragment for debugging
|
---|
| 524 | std::stringstream output;
|
---|
| 525 | output << "Contained atoms: ";
|
---|
| 526 | for (molecule::const_iterator iter = (*ListRunner)->begin(); iter != (*ListRunner)->end(); ++iter) {
|
---|
| 527 | output << (*iter)->getName() << " ";
|
---|
| 528 | }
|
---|
| 529 | LOG(2, output.str());
|
---|
| 530 | }
|
---|
[e138de] | 531 |
|
---|
[efe516] | 532 | {
|
---|
| 533 | // output xyz file
|
---|
| 534 | FragmentNumber = FixedDigitNumber(ListOfMolecules.size(), FragmentCounter++);
|
---|
| 535 | FragmentName = prefix + FragmentNumber + ".conf.xyz";
|
---|
| 536 | outputFragment.open(FragmentName.c_str(), ios::out);
|
---|
| 537 | std::stringstream output;
|
---|
| 538 | output << "Saving bond fragment No. " << FragmentNumber << "/" << FragmentCounter - 1 << " as XYZ ... ";
|
---|
| 539 | if ((intermediateResult = (*ListRunner)->OutputXYZ(&outputFragment)))
|
---|
| 540 | output << " done.";
|
---|
| 541 | else
|
---|
| 542 | output << " failed.";
|
---|
| 543 | LOG(3, output.str());
|
---|
| 544 | result = result && intermediateResult;
|
---|
| 545 | outputFragment.close();
|
---|
| 546 | outputFragment.clear();
|
---|
[e138de] | 547 | }
|
---|
| 548 |
|
---|
| 549 | // center on edge
|
---|
| 550 | (*ListRunner)->CenterEdge(&BoxDimension);
|
---|
[4c2643] | 551 | for (int k = 0; k < NDIM; k++) // if one edge is to small, set at least to 1 angstroem
|
---|
| 552 | if (BoxDimension[k] < 1.)
|
---|
| 553 | BoxDimension[k] += 1.;
|
---|
[e138de] | 554 | (*ListRunner)->SetBoxDimension(&BoxDimension); // update Box of atoms by boundary
|
---|
| 555 | for (int k = 0; k < NDIM; k++) {
|
---|
[35b698] | 556 | BoxDimension[k] = 2.5 * (World::getInstance().getConfig()->GetIsAngstroem() ? 1. : 1. / AtomicLengthToAngstroem);
|
---|
[84c494] | 557 | cell_size.at(k,k) = BoxDimension[k] * 2.;
|
---|
[e138de] | 558 | }
|
---|
[84c494] | 559 | World::getInstance().setDomain(cell_size);
|
---|
[e138de] | 560 | (*ListRunner)->Translate(&BoxDimension);
|
---|
| 561 |
|
---|
| 562 | // change path in config
|
---|
[35b698] | 563 | FragmentName = PathBackup;
|
---|
| 564 | FragmentName += "/";
|
---|
| 565 | FragmentName += FRAGMENTPREFIX;
|
---|
| 566 | FragmentName += FragmentNumber;
|
---|
| 567 | FragmentName += "/";
|
---|
| 568 | World::getInstance().getConfig()->SetDefaultPath(FragmentName.c_str());
|
---|
[e138de] | 569 |
|
---|
[efe516] | 570 | {
|
---|
| 571 | // and save as config
|
---|
| 572 | FragmentName = prefix + FragmentNumber + ".conf";
|
---|
| 573 | std::stringstream output;
|
---|
| 574 | output << "Saving bond fragment No. " << FragmentNumber << "/" << FragmentCounter - 1 << " as config ... ";
|
---|
| 575 | if ((intermediateResult = World::getInstance().getConfig()->Save(FragmentName.c_str(), (*ListRunner)->elemente, (*ListRunner))))
|
---|
| 576 | output << " done.";
|
---|
| 577 | else
|
---|
| 578 | output << " failed.";
|
---|
| 579 | LOG(3, output.str());
|
---|
| 580 | result = result && intermediateResult;
|
---|
| 581 | }
|
---|
[e138de] | 582 |
|
---|
| 583 | // restore old config
|
---|
[35b698] | 584 | World::getInstance().getConfig()->SetDefaultPath(PathBackup);
|
---|
[e138de] | 585 |
|
---|
[efe516] | 586 | {
|
---|
| 587 | // and save as mpqc input file
|
---|
| 588 | stringstream output;
|
---|
| 589 | FragmentName = prefix + FragmentNumber + ".conf";
|
---|
| 590 | output << "Saving bond fragment No. " << FragmentNumber << "/" << FragmentCounter - 1 << " as mpqc input ... ";
|
---|
| 591 | if ((intermediateResult = World::getInstance().getConfig()->SaveMPQC(FragmentName.c_str(), (*ListRunner))))
|
---|
| 592 | output << " done.";
|
---|
| 593 | else
|
---|
| 594 | output << " failed.";
|
---|
| 595 | LOG(3, output.str());
|
---|
| 596 | }
|
---|
[e138de] | 597 |
|
---|
| 598 | result = result && intermediateResult;
|
---|
| 599 | //outputFragment.close();
|
---|
| 600 | //outputFragment.clear();
|
---|
[920c70] | 601 | delete[](FragmentNumber);
|
---|
[e138de] | 602 | }
|
---|
[efe516] | 603 | LOG(0, "STATUS: done.");
|
---|
[e138de] | 604 |
|
---|
| 605 | // printing final number
|
---|
[efe516] | 606 | LOG(2, "INFO: Final number of fragments: " << FragmentCounter << ".");
|
---|
[e138de] | 607 |
|
---|
[3c58f8] | 608 | // printing final number
|
---|
[efe516] | 609 | LOG(2, "INFO: For " << count << " fragments periodic correction would have been necessary.");
|
---|
[3c58f8] | 610 |
|
---|
[b34306] | 611 | // restore cell_size
|
---|
[84c494] | 612 | World::getInstance().setDomain(cell_size_backup);
|
---|
[e138de] | 613 |
|
---|
| 614 | return result;
|
---|
| 615 | };
|
---|
| 616 |
|
---|
| 617 | /** Counts the number of molecules with the molecule::ActiveFlag set.
|
---|
[1907a7] | 618 | * \return number of molecules with ActiveFlag set to true.
|
---|
[e138de] | 619 | */
|
---|
| 620 | int MoleculeListClass::NumberOfActiveMolecules()
|
---|
| 621 | {
|
---|
| 622 | int count = 0;
|
---|
| 623 | for (MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++)
|
---|
| 624 | count += ((*ListRunner)->ActiveFlag ? 1 : 0);
|
---|
| 625 | return count;
|
---|
| 626 | };
|
---|
| 627 |
|
---|
[568be7] | 628 | /** Count all atoms in each molecule.
|
---|
| 629 | * \return number of atoms in the MoleculeListClass.
|
---|
| 630 | * TODO: the inner loop should be done by some (double molecule::CountAtom()) function
|
---|
| 631 | */
|
---|
| 632 | int MoleculeListClass::CountAllAtoms() const
|
---|
| 633 | {
|
---|
| 634 | int AtomNo = 0;
|
---|
| 635 | for (MoleculeList::const_iterator MolWalker = ListOfMolecules.begin(); MolWalker != ListOfMolecules.end(); MolWalker++) {
|
---|
[9879f6] | 636 | AtomNo += (*MolWalker)->size();
|
---|
[568be7] | 637 | }
|
---|
| 638 | return AtomNo;
|
---|
| 639 | }
|
---|
| 640 |
|
---|
[477bb2] | 641 | /***********
|
---|
| 642 | * Methods Moved here from the menus
|
---|
| 643 | */
|
---|
[568be7] | 644 |
|
---|
[477bb2] | 645 | void MoleculeListClass::createNewMolecule(periodentafel *periode) {
|
---|
[2ba827] | 646 | OBSERVE;
|
---|
[477bb2] | 647 | molecule *mol = NULL;
|
---|
[23b547] | 648 | mol = World::getInstance().createMolecule();
|
---|
[477bb2] | 649 | insert(mol);
|
---|
| 650 | };
|
---|
| 651 |
|
---|
| 652 | void MoleculeListClass::loadFromXYZ(periodentafel *periode){
|
---|
| 653 | molecule *mol = NULL;
|
---|
| 654 | Vector center;
|
---|
| 655 | char filename[MAXSTRINGSIZE];
|
---|
| 656 | Log() << Verbose(0) << "Format should be XYZ with: ShorthandOfElement\tX\tY\tZ" << endl;
|
---|
[23b547] | 657 | mol = World::getInstance().createMolecule();
|
---|
[477bb2] | 658 | do {
|
---|
| 659 | Log() << Verbose(0) << "Enter file name: ";
|
---|
| 660 | cin >> filename;
|
---|
| 661 | } while (!mol->AddXYZFile(filename));
|
---|
| 662 | mol->SetNameFromFilename(filename);
|
---|
| 663 | // center at set box dimensions
|
---|
| 664 | mol->CenterEdge(¢er);
|
---|
[cca9ef] | 665 | RealSpaceMatrix domain;
|
---|
[84c494] | 666 | for(int i =0;i<NDIM;++i)
|
---|
| 667 | domain.at(i,i) = center[i];
|
---|
| 668 | World::getInstance().setDomain(domain);
|
---|
[477bb2] | 669 | insert(mol);
|
---|
| 670 | }
|
---|
| 671 |
|
---|
| 672 | void MoleculeListClass::setMoleculeFilename() {
|
---|
| 673 | char filename[MAXSTRINGSIZE];
|
---|
| 674 | int nr;
|
---|
| 675 | molecule *mol = NULL;
|
---|
| 676 | do {
|
---|
| 677 | Log() << Verbose(0) << "Enter index of molecule: ";
|
---|
| 678 | cin >> nr;
|
---|
| 679 | mol = ReturnIndex(nr);
|
---|
| 680 | } while (mol == NULL);
|
---|
| 681 | Log() << Verbose(0) << "Enter name: ";
|
---|
| 682 | cin >> filename;
|
---|
| 683 | mol->SetNameFromFilename(filename);
|
---|
| 684 | }
|
---|
| 685 |
|
---|
| 686 | void MoleculeListClass::parseXYZIntoMolecule(){
|
---|
| 687 | char filename[MAXSTRINGSIZE];
|
---|
| 688 | int nr;
|
---|
| 689 | molecule *mol = NULL;
|
---|
| 690 | mol = NULL;
|
---|
| 691 | do {
|
---|
| 692 | Log() << Verbose(0) << "Enter index of molecule: ";
|
---|
| 693 | cin >> nr;
|
---|
| 694 | mol = ReturnIndex(nr);
|
---|
| 695 | } while (mol == NULL);
|
---|
| 696 | Log() << Verbose(0) << "Format should be XYZ with: ShorthandOfElement\tX\tY\tZ" << endl;
|
---|
| 697 | do {
|
---|
| 698 | Log() << Verbose(0) << "Enter file name: ";
|
---|
| 699 | cin >> filename;
|
---|
| 700 | } while (!mol->AddXYZFile(filename));
|
---|
| 701 | mol->SetNameFromFilename(filename);
|
---|
| 702 | };
|
---|
| 703 |
|
---|
| 704 | void MoleculeListClass::eraseMolecule(){
|
---|
| 705 | int nr;
|
---|
| 706 | molecule *mol = NULL;
|
---|
| 707 | Log() << Verbose(0) << "Enter index of molecule: ";
|
---|
| 708 | cin >> nr;
|
---|
| 709 | for(MoleculeList::iterator ListRunner = ListOfMolecules.begin(); ListRunner != ListOfMolecules.end(); ListRunner++)
|
---|
| 710 | if (nr == (*ListRunner)->IndexNr) {
|
---|
| 711 | mol = *ListRunner;
|
---|
| 712 | ListOfMolecules.erase(ListRunner);
|
---|
[23b547] | 713 | World::getInstance().destroyMolecule(mol);
|
---|
[477bb2] | 714 | break;
|
---|
| 715 | }
|
---|
| 716 | };
|
---|
| 717 |
|
---|
[77675f] | 718 |
|
---|
[e138de] | 719 | /******************************************* Class MoleculeLeafClass ************************************************/
|
---|
| 720 |
|
---|
| 721 | /** Constructor for MoleculeLeafClass root leaf.
|
---|
| 722 | * \param *Up Leaf on upper level
|
---|
| 723 | * \param *PreviousLeaf NULL - We are the first leaf on this level, otherwise points to previous in list
|
---|
| 724 | */
|
---|
| 725 | //MoleculeLeafClass::MoleculeLeafClass(MoleculeLeafClass *Up = NULL, MoleculeLeafClass *Previous = NULL)
|
---|
[97b825] | 726 | MoleculeLeafClass::MoleculeLeafClass(MoleculeLeafClass *PreviousLeaf = NULL) :
|
---|
| 727 | Leaf(NULL),
|
---|
| 728 | previous(PreviousLeaf)
|
---|
[e138de] | 729 | {
|
---|
| 730 | // if (Up != NULL)
|
---|
| 731 | // if (Up->DownLeaf == NULL) // are we the first down leaf for the upper leaf?
|
---|
| 732 | // Up->DownLeaf = this;
|
---|
| 733 | // UpLeaf = Up;
|
---|
| 734 | // DownLeaf = NULL;
|
---|
| 735 | if (previous != NULL) {
|
---|
| 736 | MoleculeLeafClass *Walker = previous->next;
|
---|
| 737 | previous->next = this;
|
---|
| 738 | next = Walker;
|
---|
| 739 | } else {
|
---|
| 740 | next = NULL;
|
---|
| 741 | }
|
---|
| 742 | };
|
---|
| 743 |
|
---|
| 744 | /** Destructor for MoleculeLeafClass.
|
---|
| 745 | */
|
---|
| 746 | MoleculeLeafClass::~MoleculeLeafClass()
|
---|
| 747 | {
|
---|
| 748 | // if (DownLeaf != NULL) {// drop leaves further down
|
---|
| 749 | // MoleculeLeafClass *Walker = DownLeaf;
|
---|
| 750 | // MoleculeLeafClass *Next;
|
---|
| 751 | // do {
|
---|
| 752 | // Next = Walker->NextLeaf;
|
---|
| 753 | // delete(Walker);
|
---|
| 754 | // Walker = Next;
|
---|
| 755 | // } while (Walker != NULL);
|
---|
| 756 | // // Last Walker sets DownLeaf automatically to NULL
|
---|
| 757 | // }
|
---|
| 758 | // remove the leaf itself
|
---|
| 759 | if (Leaf != NULL) {
|
---|
[00b59d5] | 760 | Leaf->removeAtomsinMolecule();
|
---|
[23b547] | 761 | World::getInstance().destroyMolecule(Leaf);
|
---|
[e138de] | 762 | Leaf = NULL;
|
---|
| 763 | }
|
---|
| 764 | // remove this Leaf from level list
|
---|
| 765 | if (previous != NULL)
|
---|
| 766 | previous->next = next;
|
---|
| 767 | // } else { // we are first in list (connects to UpLeaf->DownLeaf)
|
---|
| 768 | // if ((NextLeaf != NULL) && (NextLeaf->UpLeaf == NULL))
|
---|
| 769 | // NextLeaf->UpLeaf = UpLeaf; // either null as we are top level or the upleaf of the first node
|
---|
| 770 | // if (UpLeaf != NULL)
|
---|
| 771 | // UpLeaf->DownLeaf = NextLeaf; // either null as we are only leaf or NextLeaf if we are just the first
|
---|
| 772 | // }
|
---|
| 773 | // UpLeaf = NULL;
|
---|
| 774 | if (next != NULL) // are we last in list
|
---|
| 775 | next->previous = previous;
|
---|
| 776 | next = NULL;
|
---|
| 777 | previous = NULL;
|
---|
| 778 | };
|
---|
| 779 |
|
---|
| 780 | /** Adds \a molecule leaf to the tree.
|
---|
| 781 | * \param *ptr ptr to molecule to be added
|
---|
| 782 | * \param *Previous previous MoleculeLeafClass referencing level and which on the level
|
---|
| 783 | * \return true - success, false - something went wrong
|
---|
| 784 | */
|
---|
| 785 | bool MoleculeLeafClass::AddLeaf(molecule *ptr, MoleculeLeafClass *Previous)
|
---|
| 786 | {
|
---|
| 787 | return false;
|
---|
| 788 | };
|
---|
| 789 |
|
---|
| 790 | /** Fills the root stack for sites to be used as root in fragmentation depending on order or adaptivity criteria
|
---|
| 791 | * Again, as in \sa FillBondStructureFromReference steps recursively through each Leaf in this chain list of molecule's.
|
---|
| 792 | * \param *out output stream for debugging
|
---|
| 793 | * \param *&RootStack stack to be filled
|
---|
[5309ba] | 794 | * \param *AtomMask defines true/false per global Atom::Nr to mask in/out each nuclear site
|
---|
[e138de] | 795 | * \param &FragmentCounter counts through the fragments in this MoleculeLeafClass
|
---|
| 796 | * \return true - stack is non-empty, fragmentation necessary, false - stack is empty, no more sites to update
|
---|
| 797 | */
|
---|
| 798 | bool MoleculeLeafClass::FillRootStackForSubgraphs(KeyStack *&RootStack, bool *AtomMask, int &FragmentCounter)
|
---|
| 799 | {
|
---|
[9879f6] | 800 | atom *Father = NULL;
|
---|
[e138de] | 801 |
|
---|
| 802 | if (RootStack != NULL) {
|
---|
| 803 | // find first root candidates
|
---|
| 804 | if (&(RootStack[FragmentCounter]) != NULL) {
|
---|
| 805 | RootStack[FragmentCounter].clear();
|
---|
[9879f6] | 806 | for(molecule::const_iterator iter = Leaf->begin(); iter != Leaf->end(); ++iter) {
|
---|
| 807 | Father = (*iter)->GetTrueFather();
|
---|
[735b1c] | 808 | if (AtomMask[Father->getNr()]) // apply mask
|
---|
[e138de] | 809 | #ifdef ADDHYDROGEN
|
---|
[83f176] | 810 | if ((*iter)->getType()->getAtomicNumber() != 1) // skip hydrogen
|
---|
[e138de] | 811 | #endif
|
---|
[735b1c] | 812 | RootStack[FragmentCounter].push_front((*iter)->getNr());
|
---|
[e138de] | 813 | }
|
---|
| 814 | if (next != NULL)
|
---|
| 815 | next->FillRootStackForSubgraphs(RootStack, AtomMask, ++FragmentCounter);
|
---|
| 816 | } else {
|
---|
[a67d19] | 817 | DoLog(1) && (Log() << Verbose(1) << "Rootstack[" << FragmentCounter << "] is NULL." << endl);
|
---|
[e138de] | 818 | return false;
|
---|
| 819 | }
|
---|
| 820 | FragmentCounter--;
|
---|
| 821 | return true;
|
---|
| 822 | } else {
|
---|
[a67d19] | 823 | DoLog(1) && (Log() << Verbose(1) << "Rootstack is NULL." << endl);
|
---|
[e138de] | 824 | return false;
|
---|
| 825 | }
|
---|
| 826 | };
|
---|
| 827 |
|
---|
[5309ba] | 828 | /** The indices per keyset are compared to the respective father's Atom::Nr in each subgraph and thus put into \a **&FragmentList.
|
---|
[e138de] | 829 | * \param *out output stream fro debugging
|
---|
| 830 | * \param *reference reference molecule with the bond structure to be copied
|
---|
| 831 | * \param *KeySetList list with all keysets
|
---|
| 832 | * \param ***ListOfLocalAtoms Lookup table for each subgraph and index of each atom in global molecule, may be NULL on start, then it is filled
|
---|
| 833 | * \param **&FragmentList list to be allocated and returned
|
---|
| 834 | * \param &FragmentCounter counts the fragments as we move along the list
|
---|
| 835 | * \param FreeList true - ***ListOfLocalAtoms is free'd before return, false - it is not
|
---|
| 836 | * \retuen true - success, false - failure
|
---|
| 837 | */
|
---|
| 838 | bool MoleculeLeafClass::AssignKeySetsToFragment(molecule *reference, Graph *KeySetList, atom ***&ListOfLocalAtoms, Graph **&FragmentList, int &FragmentCounter, bool FreeList)
|
---|
| 839 | {
|
---|
| 840 | bool status = true;
|
---|
| 841 | int KeySetCounter = 0;
|
---|
| 842 |
|
---|
[a67d19] | 843 | DoLog(1) && (Log() << Verbose(1) << "Begin of AssignKeySetsToFragment." << endl);
|
---|
[e138de] | 844 | // fill ListOfLocalAtoms if NULL was given
|
---|
[c6123b] | 845 | if (!Leaf->FillListOfLocalAtoms(ListOfLocalAtoms[FragmentCounter], reference->getAtomCount())) {
|
---|
[a67d19] | 846 | DoLog(1) && (Log() << Verbose(1) << "Filling of ListOfLocalAtoms failed." << endl);
|
---|
[e138de] | 847 | return false;
|
---|
| 848 | }
|
---|
| 849 |
|
---|
| 850 | // allocate fragment list
|
---|
| 851 | if (FragmentList == NULL) {
|
---|
| 852 | KeySetCounter = Count();
|
---|
[920c70] | 853 | FragmentList = new Graph*[KeySetCounter];
|
---|
| 854 | for (int i=0;i<KeySetCounter;i++)
|
---|
| 855 | FragmentList[i] = NULL;
|
---|
[e138de] | 856 | KeySetCounter = 0;
|
---|
| 857 | }
|
---|
| 858 |
|
---|
| 859 | if ((KeySetList != NULL) && (KeySetList->size() != 0)) { // if there are some scanned keysets at all
|
---|
| 860 | // assign scanned keysets
|
---|
| 861 | if (FragmentList[FragmentCounter] == NULL)
|
---|
| 862 | FragmentList[FragmentCounter] = new Graph;
|
---|
| 863 | KeySet *TempSet = new KeySet;
|
---|
| 864 | for (Graph::iterator runner = KeySetList->begin(); runner != KeySetList->end(); runner++) { // key sets contain global numbers!
|
---|
[735b1c] | 865 | if (ListOfLocalAtoms[FragmentCounter][reference->FindAtom(*((*runner).first.begin()))->getNr()] != NULL) {// as we may assume that that bond structure is unchanged, we only test the first key in each set
|
---|
[e138de] | 866 | // translate keyset to local numbers
|
---|
| 867 | for (KeySet::iterator sprinter = (*runner).first.begin(); sprinter != (*runner).first.end(); sprinter++)
|
---|
[735b1c] | 868 | TempSet->insert(ListOfLocalAtoms[FragmentCounter][reference->FindAtom(*sprinter)->getNr()]->getNr());
|
---|
[e138de] | 869 | // insert into FragmentList
|
---|
| 870 | FragmentList[FragmentCounter]->insert(GraphPair(*TempSet, pair<int, double> (KeySetCounter++, (*runner).second.second)));
|
---|
| 871 | }
|
---|
| 872 | TempSet->clear();
|
---|
| 873 | }
|
---|
| 874 | delete (TempSet);
|
---|
| 875 | if (KeySetCounter == 0) {// if there are no keysets, delete the list
|
---|
[a67d19] | 876 | DoLog(1) && (Log() << Verbose(1) << "KeySetCounter is zero, deleting FragmentList." << endl);
|
---|
[e138de] | 877 | delete (FragmentList[FragmentCounter]);
|
---|
| 878 | } else
|
---|
[a67d19] | 879 | DoLog(1) && (Log() << Verbose(1) << KeySetCounter << " keysets were assigned to subgraph " << FragmentCounter << "." << endl);
|
---|
[e138de] | 880 | FragmentCounter++;
|
---|
| 881 | if (next != NULL)
|
---|
| 882 | next->AssignKeySetsToFragment(reference, KeySetList, ListOfLocalAtoms, FragmentList, FragmentCounter, FreeList);
|
---|
| 883 | FragmentCounter--;
|
---|
| 884 | } else
|
---|
[a67d19] | 885 | DoLog(1) && (Log() << Verbose(1) << "KeySetList is NULL or empty." << endl);
|
---|
[e138de] | 886 |
|
---|
| 887 | if ((FreeList) && (ListOfLocalAtoms != NULL)) {
|
---|
| 888 | // free the index lookup list
|
---|
[920c70] | 889 | delete[](ListOfLocalAtoms[FragmentCounter]);
|
---|
[e138de] | 890 | }
|
---|
[a67d19] | 891 | DoLog(1) && (Log() << Verbose(1) << "End of AssignKeySetsToFragment." << endl);
|
---|
[e138de] | 892 | return status;
|
---|
| 893 | };
|
---|
| 894 |
|
---|
| 895 | /** Translate list into global numbers (i.e. ones that are valid in "this" molecule, not in MolecularWalker->Leaf)
|
---|
| 896 | * \param *out output stream for debugging
|
---|
| 897 | * \param **FragmentList Graph with local numbers per fragment
|
---|
| 898 | * \param &FragmentCounter counts the fragments as we move along the list
|
---|
| 899 | * \param &TotalNumberOfKeySets global key set counter
|
---|
| 900 | * \param &TotalGraph Graph to be filled with global numbers
|
---|
| 901 | */
|
---|
| 902 | void MoleculeLeafClass::TranslateIndicesToGlobalIDs(Graph **FragmentList, int &FragmentCounter, int &TotalNumberOfKeySets, Graph &TotalGraph)
|
---|
| 903 | {
|
---|
[a67d19] | 904 | DoLog(1) && (Log() << Verbose(1) << "Begin of TranslateIndicesToGlobalIDs." << endl);
|
---|
[e138de] | 905 | KeySet *TempSet = new KeySet;
|
---|
| 906 | if (FragmentList[FragmentCounter] != NULL) {
|
---|
| 907 | for (Graph::iterator runner = FragmentList[FragmentCounter]->begin(); runner != FragmentList[FragmentCounter]->end(); runner++) {
|
---|
| 908 | for (KeySet::iterator sprinter = (*runner).first.begin(); sprinter != (*runner).first.end(); sprinter++)
|
---|
[735b1c] | 909 | TempSet->insert((Leaf->FindAtom(*sprinter))->GetTrueFather()->getNr());
|
---|
[e138de] | 910 | TotalGraph.insert(GraphPair(*TempSet, pair<int, double> (TotalNumberOfKeySets++, (*runner).second.second)));
|
---|
| 911 | TempSet->clear();
|
---|
| 912 | }
|
---|
| 913 | delete (TempSet);
|
---|
| 914 | } else {
|
---|
[a67d19] | 915 | DoLog(1) && (Log() << Verbose(1) << "FragmentList is NULL." << endl);
|
---|
[e138de] | 916 | }
|
---|
| 917 | if (next != NULL)
|
---|
| 918 | next->TranslateIndicesToGlobalIDs(FragmentList, ++FragmentCounter, TotalNumberOfKeySets, TotalGraph);
|
---|
| 919 | FragmentCounter--;
|
---|
[a67d19] | 920 | DoLog(1) && (Log() << Verbose(1) << "End of TranslateIndicesToGlobalIDs." << endl);
|
---|
[e138de] | 921 | };
|
---|
| 922 |
|
---|
| 923 | /** Simply counts the number of items in the list, from given MoleculeLeafClass.
|
---|
| 924 | * \return number of items
|
---|
| 925 | */
|
---|
| 926 | int MoleculeLeafClass::Count() const
|
---|
| 927 | {
|
---|
| 928 | if (next != NULL)
|
---|
| 929 | return next->Count() + 1;
|
---|
| 930 | else
|
---|
| 931 | return 1;
|
---|
| 932 | };
|
---|
| 933 |
|
---|