[3ae731] | 1 | /*
|
---|
| 2 | * Project: MoleCuilder
|
---|
| 3 | * Description: creates and alters molecular systems
|
---|
[0aa122] | 4 | * Copyright (C) 2010-2012 University of Bonn. All rights reserved.
|
---|
[5aaa43] | 5 | * Copyright (C) 2013 Frederik Heber. All rights reserved.
|
---|
[94d5ac6] | 6 | *
|
---|
| 7 | *
|
---|
| 8 | * This file is part of MoleCuilder.
|
---|
| 9 | *
|
---|
| 10 | * MoleCuilder is free software: you can redistribute it and/or modify
|
---|
| 11 | * it under the terms of the GNU General Public License as published by
|
---|
| 12 | * the Free Software Foundation, either version 2 of the License, or
|
---|
| 13 | * (at your option) any later version.
|
---|
| 14 | *
|
---|
| 15 | * MoleCuilder is distributed in the hope that it will be useful,
|
---|
| 16 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 17 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 18 | * GNU General Public License for more details.
|
---|
| 19 | *
|
---|
| 20 | * You should have received a copy of the GNU General Public License
|
---|
| 21 | * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
|
---|
[3ae731] | 22 | */
|
---|
| 23 |
|
---|
| 24 | /*
|
---|
| 25 | * PdbParser.cpp
|
---|
| 26 | *
|
---|
| 27 | * Created on: Aug 17, 2010
|
---|
| 28 | * Author: heber
|
---|
| 29 | */
|
---|
| 30 |
|
---|
| 31 | // include config.h
|
---|
| 32 | #ifdef HAVE_CONFIG_H
|
---|
| 33 | #include <config.h>
|
---|
| 34 | #endif
|
---|
| 35 |
|
---|
[ad011c] | 36 | #include "CodePatterns/MemDebug.hpp"
|
---|
[3ae731] | 37 |
|
---|
[ad011c] | 38 | #include "CodePatterns/Assert.hpp"
|
---|
| 39 | #include "CodePatterns/Log.hpp"
|
---|
| 40 | #include "CodePatterns/toString.hpp"
|
---|
| 41 | #include "CodePatterns/Verbose.hpp"
|
---|
[42127c] | 42 |
|
---|
[6f0841] | 43 | #include "Atom/atom.hpp"
|
---|
[129204] | 44 | #include "Bond/bond.hpp"
|
---|
[42127c] | 45 | #include "Descriptors/AtomIdDescriptor.hpp"
|
---|
[3bdb6d] | 46 | #include "Element/element.hpp"
|
---|
| 47 | #include "Element/periodentafel.hpp"
|
---|
[42127c] | 48 | #include "molecule.hpp"
|
---|
[4fbca9c] | 49 | #include "Parser/PdbParser.hpp"
|
---|
[073a9e4] | 50 | #include "World.hpp"
|
---|
| 51 | #include "WorldTime.hpp"
|
---|
[bb6193] | 52 |
|
---|
[105b72] | 53 | #include <algorithm>
|
---|
| 54 | #include <cmath>
|
---|
[3ae731] | 55 | #include <map>
|
---|
| 56 | #include <vector>
|
---|
| 57 |
|
---|
[bb6193] | 58 | #include <iostream>
|
---|
| 59 | #include <iomanip>
|
---|
[3ae731] | 60 |
|
---|
| 61 | using namespace std;
|
---|
| 62 |
|
---|
[765f16] | 63 | // declare specialized static variables
|
---|
| 64 | const std::string FormatParserTrait<pdb>::name = "pdb";
|
---|
| 65 | const std::string FormatParserTrait<pdb>::suffix = "pdb";
|
---|
| 66 | const ParserTypes FormatParserTrait<pdb>::type = pdb;
|
---|
| 67 |
|
---|
[3ae731] | 68 | /**
|
---|
| 69 | * Constructor.
|
---|
| 70 | */
|
---|
[765f16] | 71 | FormatParser< pdb >::FormatParser() :
|
---|
| 72 | FormatParser_common(NULL)
|
---|
| 73 | {
|
---|
[4fbca9c] | 74 | knownTokens["ATOM"] = PdbKey::Atom;
|
---|
[16462f] | 75 | knownTokens["HETATM"] = PdbKey::Atom;
|
---|
[4fbca9c] | 76 | knownTokens["TER"] = PdbKey::Filler;
|
---|
[9dba5f] | 77 | knownTokens["END"] = PdbKey::EndOfTimestep;
|
---|
[4fbca9c] | 78 | knownTokens["CONECT"] = PdbKey::Connect;
|
---|
| 79 | knownTokens["REMARK"] = PdbKey::Remark;
|
---|
[9dba5f] | 80 | knownTokens[""] = PdbKey::EndOfTimestep;
|
---|
[16462f] | 81 |
|
---|
| 82 | // argh, why can't just PdbKey::X+(size_t)i
|
---|
| 83 | PositionEnumMap[0] = PdbKey::X;
|
---|
| 84 | PositionEnumMap[1] = PdbKey::Y;
|
---|
| 85 | PositionEnumMap[2] = PdbKey::Z;
|
---|
[3ae731] | 86 | }
|
---|
| 87 |
|
---|
| 88 | /**
|
---|
| 89 | * Destructor.
|
---|
| 90 | */
|
---|
[765f16] | 91 | FormatParser< pdb >::~FormatParser()
|
---|
| 92 | {
|
---|
[873037] | 93 | PdbAtomInfoContainer::clearknownDataKeys();
|
---|
[3ae731] | 94 | additionalAtomData.clear();
|
---|
[4fbca9c] | 95 | }
|
---|
| 96 |
|
---|
| 97 |
|
---|
| 98 | /** Parses the initial word of the given \a line and returns the token type.
|
---|
| 99 | *
|
---|
| 100 | * @param line line to scan
|
---|
| 101 | * @return token type
|
---|
| 102 | */
|
---|
[955b91] | 103 | enum PdbKey::KnownTokens FormatParser< pdb >::getToken(std::string &line)
|
---|
[4fbca9c] | 104 | {
|
---|
| 105 | // look for first space
|
---|
[3dfd9c] | 106 | std::string token = line.substr(0,6);
|
---|
| 107 | const size_t space_location = token.find(' ');
|
---|
| 108 | const size_t tab_location = token.find('\t');
|
---|
[4fbca9c] | 109 | size_t location = space_location < tab_location ? space_location : tab_location;
|
---|
| 110 | if (location != string::npos) {
|
---|
[47d041] | 111 | //LOG(1, "Found space at position " << space_location);
|
---|
[3dfd9c] | 112 | token = token.substr(0,space_location);
|
---|
[4fbca9c] | 113 | }
|
---|
| 114 |
|
---|
[47d041] | 115 | //LOG(1, "Token is " << token);
|
---|
[4fbca9c] | 116 | if (knownTokens.count(token) == 0)
|
---|
| 117 | return PdbKey::NoToken;
|
---|
| 118 | else
|
---|
| 119 | return knownTokens[token];
|
---|
| 120 |
|
---|
| 121 | return PdbKey::NoToken;
|
---|
[3ae731] | 122 | }
|
---|
| 123 |
|
---|
| 124 | /**
|
---|
[4fbca9c] | 125 | * Loads atoms from a PDB-formatted file.
|
---|
[3ae731] | 126 | *
|
---|
[4fbca9c] | 127 | * \param PDB file
|
---|
[3ae731] | 128 | */
|
---|
[765f16] | 129 | void FormatParser< pdb >::load(istream* file) {
|
---|
[4fbca9c] | 130 | string line;
|
---|
| 131 | size_t linecount = 0;
|
---|
| 132 | enum PdbKey::KnownTokens token;
|
---|
| 133 |
|
---|
[c0e28c] | 134 | // reset id maps for this file (to correctly parse CONECT entries)
|
---|
| 135 | resetIdAssociations();
|
---|
[16462f] | 136 |
|
---|
[9dba5f] | 137 | bool NotEndOfFile = true;
|
---|
[4fbca9c] | 138 | molecule *newmol = World::getInstance().createMolecule();
|
---|
| 139 | newmol->ActiveFlag = true;
|
---|
[b0a2e3] | 140 | unsigned int step = 0;
|
---|
[4fbca9c] | 141 | while (NotEndOfFile) {
|
---|
[9dba5f] | 142 | bool NotEndOfTimestep = true;
|
---|
[b0a2e3] | 143 | while (NotEndOfTimestep && NotEndOfFile) {
|
---|
[9dba5f] | 144 | std::getline(*file, line, '\n');
|
---|
[b0a2e3] | 145 | if (!line.empty()) {
|
---|
| 146 | // extract first token
|
---|
| 147 | token = getToken(line);
|
---|
| 148 | switch (token) {
|
---|
| 149 | case PdbKey::Atom:
|
---|
| 150 | LOG(3,"INFO: Parsing ATOM entry for time step " << step << ".");
|
---|
| 151 | readAtomDataLine(step, line, newmol);
|
---|
| 152 | break;
|
---|
| 153 | case PdbKey::Remark:
|
---|
| 154 | LOG(3,"INFO: Parsing REM entry for time step " << step << ".");
|
---|
| 155 | break;
|
---|
| 156 | case PdbKey::Connect:
|
---|
| 157 | LOG(3,"INFO: Parsing CONECT entry for time step " << step << ".");
|
---|
| 158 | readNeighbors(step, line);
|
---|
| 159 | break;
|
---|
| 160 | case PdbKey::Filler:
|
---|
| 161 | LOG(3,"INFO: Stumbled upon Filler entry for time step " << step << ".");
|
---|
| 162 | break;
|
---|
| 163 | case PdbKey::EndOfTimestep:
|
---|
[44f53e] | 164 | LOG(1,"INFO: Parsing END entry or empty line for time step " << step << ".");
|
---|
[b0a2e3] | 165 | NotEndOfTimestep = false;
|
---|
| 166 | break;
|
---|
| 167 | default:
|
---|
| 168 | // TODO: put a throw here
|
---|
[47d041] | 169 | ELOG(2, "Unknown token: '" << line << "' for time step " << step << ".");
|
---|
[765f16] | 170 | //ASSERT(0, "FormatParser< pdb >::load() - Unknown token in line "+toString(linecount)+": "+line+".");
|
---|
[b0a2e3] | 171 | break;
|
---|
| 172 | }
|
---|
[9dba5f] | 173 | }
|
---|
| 174 | NotEndOfFile = NotEndOfFile && (file->good());
|
---|
| 175 | linecount++;
|
---|
[4fbca9c] | 176 | }
|
---|
[b0a2e3] | 177 | ++step;
|
---|
| 178 | }
|
---|
[3dfd9c] | 179 | LOG(4, "INFO: Listing all newly parsed atoms.");
|
---|
| 180 | BOOST_FOREACH(atom *_atom, *newmol)
|
---|
[48801a] | 181 | LOG(4, "INFO: Atom " << _atom->getName() << " " << *dynamic_cast<AtomInfo *>(_atom) << ".");
|
---|
[4afa46] | 182 |
|
---|
| 183 | // refresh atom::nr and atom::name
|
---|
| 184 | newmol->getAtomCount();
|
---|
[3ae731] | 185 | }
|
---|
| 186 |
|
---|
| 187 | /**
|
---|
[73916f] | 188 | * Saves the \a atoms into as a PDB file.
|
---|
[3ae731] | 189 | *
|
---|
| 190 | * \param file where to save the state
|
---|
[73916f] | 191 | * \param atoms atoms to store
|
---|
[3ae731] | 192 | */
|
---|
[fac58f] | 193 | void FormatParser< pdb >::save(
|
---|
| 194 | ostream* file,
|
---|
| 195 | const std::vector<const atom *> &AtomList)
|
---|
[73916f] | 196 | {
|
---|
[830b3e] | 197 | LOG(2, "DEBUG: Saving changes to pdb.");
|
---|
[9dba5f] | 198 |
|
---|
| 199 | // check for maximum number of time steps
|
---|
| 200 | size_t max_timesteps = 0;
|
---|
[45b45d] | 201 | for (vector<const atom *>::const_iterator atomIt = AtomList.begin();
|
---|
| 202 | atomIt != AtomList.end(); atomIt++) {
|
---|
| 203 | const atom * _atom = *atomIt;
|
---|
| 204 | LOG(4, "INFO: Atom " << _atom->getName() << " "
|
---|
| 205 | << *dynamic_cast<const AtomInfo *>(_atom) << ".");
|
---|
[9dba5f] | 206 | if (_atom->getTrajectorySize() > max_timesteps)
|
---|
| 207 | max_timesteps = _atom->getTrajectorySize();
|
---|
[bb6193] | 208 | }
|
---|
[9dba5f] | 209 | LOG(2,"INFO: Found a maximum of " << max_timesteps << " time steps to store.");
|
---|
[3ae731] | 210 |
|
---|
[9dba5f] | 211 | // re-distribute serials
|
---|
[c0e28c] | 212 | resetIdAssociations();
|
---|
[9dba5f] | 213 | // (new atoms might have been added)
|
---|
| 214 | int AtomNo = 1; // serial number starts at 1 in pdb
|
---|
[45b45d] | 215 | for (vector<const atom *>::const_iterator atomIt = AtomList.begin();
|
---|
| 216 | atomIt != AtomList.end(); atomIt++) {
|
---|
[9dba5f] | 217 | PdbAtomInfoContainer &atomInfo = getadditionalAtomData(*atomIt);
|
---|
[c0e28c] | 218 | associateLocaltoGlobalId(AtomNo, (*atomIt)->getId());
|
---|
[9dba5f] | 219 | atomInfo.set(PdbKey::serial, toString(AtomNo));
|
---|
[c0e28c] | 220 | ++AtomNo;
|
---|
[9dba5f] | 221 | }
|
---|
| 222 |
|
---|
[5c5472] | 223 | // store all time steps (always do first step)
|
---|
| 224 | for (size_t step = 0; (step == 0) || (step < max_timesteps); ++step) {
|
---|
[9dba5f] | 225 | {
|
---|
| 226 | // add initial remark
|
---|
| 227 | *file << "REMARK created by molecuilder on ";
|
---|
| 228 | time_t now = time((time_t *)NULL); // Get the system time and put it into 'now' as 'calender time'
|
---|
| 229 | // ctime ends in \n\0, we have to cut away the newline
|
---|
| 230 | std::string time(ctime(&now));
|
---|
| 231 | size_t pos = time.find('\n');
|
---|
| 232 | if (pos != 0)
|
---|
| 233 | *file << time.substr(0,pos);
|
---|
| 234 | else
|
---|
| 235 | *file << time;
|
---|
| 236 | *file << ", time step " << step;
|
---|
| 237 | *file << endl;
|
---|
[16462f] | 238 | }
|
---|
[9dba5f] | 239 |
|
---|
| 240 | {
|
---|
| 241 | std::map<size_t,size_t> MolIdMap;
|
---|
| 242 | size_t MolNo = 1; // residue number starts at 1 in pdb
|
---|
[45b45d] | 243 | for (vector<const atom *>::const_iterator atomIt = AtomList.begin();
|
---|
| 244 | atomIt != AtomList.end(); atomIt++) {
|
---|
[9dba5f] | 245 | const molecule *mol = (*atomIt)->getMolecule();
|
---|
| 246 | if ((mol != NULL) && (MolIdMap.find(mol->getId()) == MolIdMap.end())) {
|
---|
| 247 | MolIdMap[mol->getId()] = MolNo++;
|
---|
| 248 | }
|
---|
[bb6193] | 249 | }
|
---|
[9dba5f] | 250 | const size_t MaxMol = MolNo;
|
---|
| 251 |
|
---|
| 252 | // have a count per element and per molecule (0 is for all homeless atoms)
|
---|
| 253 | std::vector<int> **elementNo = new std::vector<int>*[MaxMol];
|
---|
| 254 | for (size_t i = 0; i < MaxMol; ++i)
|
---|
| 255 | elementNo[i] = new std::vector<int>(MAX_ELEMENTS,1);
|
---|
| 256 | char name[MAXSTRINGSIZE];
|
---|
| 257 | std::string ResidueName;
|
---|
| 258 |
|
---|
| 259 | // write ATOMs
|
---|
[45b45d] | 260 | for (vector<const atom *>::const_iterator atomIt = AtomList.begin();
|
---|
| 261 | atomIt != AtomList.end(); atomIt++) {
|
---|
[9dba5f] | 262 | PdbAtomInfoContainer &atomInfo = getadditionalAtomData(*atomIt);
|
---|
| 263 | // gather info about residue
|
---|
| 264 | const molecule *mol = (*atomIt)->getMolecule();
|
---|
| 265 | if (mol == NULL) {
|
---|
| 266 | MolNo = 0;
|
---|
| 267 | atomInfo.set(PdbKey::resSeq, "0");
|
---|
| 268 | } else {
|
---|
| 269 | ASSERT(MolIdMap.find(mol->getId()) != MolIdMap.end(),
|
---|
[765f16] | 270 | "FormatParser< pdb >::save() - Mol id "+toString(mol->getId())+" not present despite we set it?!");
|
---|
[9dba5f] | 271 | MolNo = MolIdMap[mol->getId()];
|
---|
| 272 | atomInfo.set(PdbKey::resSeq, toString(MolIdMap[mol->getId()]));
|
---|
| 273 | if (atomInfo.get<std::string>(PdbKey::resName) == "-")
|
---|
| 274 | atomInfo.set(PdbKey::resName, mol->getName().substr(0,3));
|
---|
| 275 | }
|
---|
| 276 | // get info about atom
|
---|
| 277 | const size_t Z = (*atomIt)->getType()->getAtomicNumber();
|
---|
| 278 | if (atomInfo.get<std::string>(PdbKey::name) == "-") { // if no name set, give it a new name
|
---|
| 279 | sprintf(name, "%2s%02d",(*atomIt)->getType()->getSymbol().c_str(), (*elementNo[MolNo])[Z]);
|
---|
| 280 | (*elementNo[MolNo])[Z] = ((*elementNo[MolNo])[Z]+1) % 100; // confine to two digits
|
---|
| 281 | atomInfo.set(PdbKey::name, name);
|
---|
| 282 | }
|
---|
| 283 | // set position
|
---|
| 284 | for (size_t i=0; i<NDIM;++i) {
|
---|
| 285 | stringstream position;
|
---|
[105b72] | 286 | position << setprecision(7) << (*atomIt)->getPositionAtStep(step).at(i);
|
---|
[9dba5f] | 287 | atomInfo.set(PositionEnumMap[i], position.str());
|
---|
| 288 | }
|
---|
| 289 | // change element and charge if changed
|
---|
[8990879] | 290 | if (atomInfo.get<std::string>(PdbKey::element) != (*atomIt)->getType()->getSymbol()) {
|
---|
| 291 | std::string symbol = (*atomIt)->getType()->getSymbol();
|
---|
| 292 | if ((symbol[1] >= 'a') && (symbol[1] <= 'z'))
|
---|
| 293 | symbol[1] = (symbol[1] - 'a') + 'A';
|
---|
| 294 | atomInfo.set(PdbKey::element, symbol);
|
---|
| 295 | }
|
---|
[797a40] | 296 | // change particlename and charge if changed
|
---|
| 297 | if (atomInfo.get<std::string>(PdbKey::name) != (*atomIt)->getParticleName()) {
|
---|
| 298 | std::string particlename = (*atomIt)->getParticleName();
|
---|
| 299 | atomInfo.set(PdbKey::name, particlename);
|
---|
| 300 | }
|
---|
[9dba5f] | 301 |
|
---|
| 302 | // finally save the line
|
---|
| 303 | saveLine(file, atomInfo);
|
---|
[16462f] | 304 | }
|
---|
[9dba5f] | 305 | for (size_t i = 0; i < MaxMol; ++i)
|
---|
| 306 | delete elementNo[i];
|
---|
[3d70e3] | 307 | delete[] elementNo;
|
---|
[3ae731] | 308 |
|
---|
[9dba5f] | 309 | // write CONECTs
|
---|
[45b45d] | 310 | for (vector<const atom *>::const_iterator atomIt = AtomList.begin();
|
---|
| 311 | atomIt != AtomList.end(); atomIt++) {
|
---|
[9dba5f] | 312 | writeNeighbors(file, 4, *atomIt);
|
---|
| 313 | }
|
---|
[bb6193] | 314 | }
|
---|
[9dba5f] | 315 | // END
|
---|
| 316 | *file << "END" << endl;
|
---|
[3ae731] | 317 | }
|
---|
| 318 |
|
---|
| 319 | }
|
---|
| 320 |
|
---|
[6bc86c] | 321 | /** Add default info, when new atom is added to World.
|
---|
| 322 | *
|
---|
| 323 | * @param id of atom
|
---|
| 324 | */
|
---|
[765f16] | 325 | void FormatParser< pdb >::AtomInserted(atomId_t id)
|
---|
[6bc86c] | 326 | {
|
---|
[765f16] | 327 | //LOG(3, "FormatParser< pdb >::AtomInserted() - notified of atom " << id << "'s insertion.");
|
---|
[6bc86c] | 328 | ASSERT(!isPresentadditionalAtomData(id),
|
---|
[765f16] | 329 | "FormatParser< pdb >::AtomInserted() - additionalAtomData already present for newly added atom "
|
---|
[6bc86c] | 330 | +toString(id)+".");
|
---|
| 331 | // don't insert here as this is our check whether we are in the first time step
|
---|
| 332 | //additionalAtomData.insert( std::make_pair(id, defaultAdditionalData) );
|
---|
| 333 | }
|
---|
| 334 |
|
---|
| 335 | /** Remove additional AtomData info, when atom has been removed from World.
|
---|
| 336 | *
|
---|
| 337 | * @param id of atom
|
---|
| 338 | */
|
---|
[765f16] | 339 | void FormatParser< pdb >::AtomRemoved(atomId_t id)
|
---|
[6bc86c] | 340 | {
|
---|
[765f16] | 341 | //LOG(3, "FormatParser< pdb >::AtomRemoved() - notified of atom " << id << "'s removal.");
|
---|
[8bf9c6] | 342 | std::map<const atomId_t, PdbAtomInfoContainer>::iterator iter = additionalAtomData.find(id);
|
---|
[6bc86c] | 343 | // as we do not insert AtomData on AtomInserted, we cannot be assured of its presence
|
---|
| 344 | // ASSERT(iter != additionalAtomData.end(),
|
---|
[765f16] | 345 | // "FormatParser< pdb >::AtomRemoved() - additionalAtomData is not present for atom "
|
---|
[6bc86c] | 346 | // +toString(id)+" to remove.");
|
---|
| 347 | if (iter != additionalAtomData.end()) {
|
---|
| 348 | additionalAtomData.erase(iter);
|
---|
| 349 | }
|
---|
| 350 | }
|
---|
| 351 |
|
---|
| 352 |
|
---|
[9dba5f] | 353 | /** Checks whether there is an entry for the given atom's \a _id.
|
---|
| 354 | *
|
---|
| 355 | * @param _id atom's id we wish to check on
|
---|
| 356 | * @return true - entry present, false - only for atom's father or no entry
|
---|
| 357 | */
|
---|
[8bf9c6] | 358 | bool FormatParser< pdb >::isPresentadditionalAtomData(const atomId_t _id) const
|
---|
[9dba5f] | 359 | {
|
---|
[8bf9c6] | 360 | std::map<const atomId_t, PdbAtomInfoContainer>::const_iterator iter = additionalAtomData.find(_id);
|
---|
| 361 | return (iter != additionalAtomData.end());
|
---|
[9dba5f] | 362 | }
|
---|
| 363 |
|
---|
| 364 |
|
---|
[93fd43e] | 365 | /** Either returns reference to present entry or creates new with default values.
|
---|
| 366 | *
|
---|
| 367 | * @param _atom atom whose entry we desire
|
---|
| 368 | * @return
|
---|
| 369 | */
|
---|
[fac58f] | 370 | PdbAtomInfoContainer& FormatParser< pdb >::getadditionalAtomData(const atom * const _atom)
|
---|
[93fd43e] | 371 | {
|
---|
| 372 | if (additionalAtomData.find(_atom->getId()) != additionalAtomData.end()) {
|
---|
[910a5d] | 373 | } else if (additionalAtomData.find(_atom->getFather()->getId()) != additionalAtomData.end()) {
|
---|
[93fd43e] | 374 | // use info from direct father
|
---|
[910a5d] | 375 | additionalAtomData[_atom->getId()] = additionalAtomData[_atom->getFather()->getId()];
|
---|
[93fd43e] | 376 | } else if (additionalAtomData.find(_atom->GetTrueFather()->getId()) != additionalAtomData.end()) {
|
---|
| 377 | // use info from topmost father
|
---|
| 378 | additionalAtomData[_atom->getId()] = additionalAtomData[_atom->GetTrueFather()->getId()];
|
---|
| 379 | } else {
|
---|
| 380 | // create new entry use default values if nothing else is known
|
---|
| 381 | additionalAtomData[_atom->getId()] = defaultAdditionalData;
|
---|
| 382 | }
|
---|
| 383 | return additionalAtomData[_atom->getId()];
|
---|
| 384 | }
|
---|
| 385 |
|
---|
[105b72] | 386 | /** Tiny helper function to print a float with a most 8 digits.
|
---|
| 387 | *
|
---|
| 388 | * A few examples best give the picture:
|
---|
| 389 | * 1234.678
|
---|
| 390 | * 1.234
|
---|
| 391 | * 0.001
|
---|
| 392 | * 0.100
|
---|
| 393 | * 1234567.
|
---|
| 394 | * 123456.7
|
---|
| 395 | * -1234.56
|
---|
| 396 | *
|
---|
| 397 | * \param value
|
---|
| 398 | * \return string representation
|
---|
| 399 | */
|
---|
| 400 | const std::string FormatParser< pdb >::printCoordinate(
|
---|
| 401 | const double value)
|
---|
| 402 | {
|
---|
| 403 | size_t meaningful_bits=7; // one for decimal dot
|
---|
| 404 | if (value < 0) //one for the minus sign
|
---|
| 405 | --meaningful_bits;
|
---|
| 406 | // count digits before dot (without minus and round towards zero!)
|
---|
[4e65af] | 407 | int full = floor(fabs(value));
|
---|
[105b72] | 408 | size_t bits_before_dot = 1;
|
---|
| 409 | {
|
---|
| 410 | int tmp = full;
|
---|
| 411 | for (;bits_before_dot < meaningful_bits;++bits_before_dot) {
|
---|
| 412 | // even if value is 0...somethingish, we still must start at one digit
|
---|
| 413 | tmp = tmp/10;
|
---|
| 414 | if (tmp == 0)
|
---|
| 415 | break;
|
---|
| 416 | }
|
---|
| 417 | }
|
---|
| 418 | // this fixes bits available after dot
|
---|
| 419 | const size_t bits_after_dot = std::min((int)meaningful_bits - (int)bits_before_dot, 3);
|
---|
| 420 | stringstream position;
|
---|
| 421 | if (bits_after_dot > 0) {
|
---|
| 422 | if (value < 0)
|
---|
| 423 | position << "-";
|
---|
[4e65af] | 424 | // truncate 999.9 to 999 and not to 1000! (hence, extra check!)
|
---|
| 425 | int remainder = round((abs(value)-full)*pow(10,bits_after_dot));
|
---|
| 426 | if (remainder >= pow(10,bits_after_dot)) {
|
---|
| 427 | remainder = 0;
|
---|
| 428 | ++full;
|
---|
| 429 | }
|
---|
[105b72] | 430 | position << full << "." << setfill('0') << setw(bits_after_dot) << remainder;
|
---|
| 431 | if (bits_after_dot == 2)
|
---|
| 432 | ELOG(2, "PdbParser is writing coordinates with just a two decimal places.");
|
---|
| 433 | if (bits_after_dot == 1)
|
---|
| 434 | ELOG(1, "PdbParser is writing coordinates with just a single decimal places.");
|
---|
| 435 | } else {
|
---|
| 436 | ELOG(0, "PdbParser is writing coordinates without any decimal places.");
|
---|
| 437 | position << full << ".";
|
---|
| 438 | }
|
---|
| 439 | return position.str();
|
---|
| 440 | }
|
---|
| 441 |
|
---|
[3ae731] | 442 | /**
|
---|
[4fbca9c] | 443 | * Writes one line of PDB-formatted data to the provided stream.
|
---|
[3ae731] | 444 | *
|
---|
| 445 | * \param stream where to write the line to
|
---|
[bb6193] | 446 | * \param *currentAtom the atom of which information should be written
|
---|
| 447 | * \param AtomNo serial number of atom
|
---|
[16462f] | 448 | * \param *name name of atom, i.e. H01
|
---|
| 449 | * \param ResidueName Name of molecule
|
---|
[bb6193] | 450 | * \param ResidueNo number of residue
|
---|
[3ae731] | 451 | */
|
---|
[765f16] | 452 | void FormatParser< pdb >::saveLine(
|
---|
[16462f] | 453 | ostream* file,
|
---|
| 454 | const PdbAtomInfoContainer &atomInfo)
|
---|
| 455 | {
|
---|
| 456 | *file << setfill(' ') << left << setw(6)
|
---|
| 457 | << atomInfo.get<std::string>(PdbKey::token);
|
---|
| 458 | *file << setfill(' ') << right << setw(5)
|
---|
[1f5e97] | 459 | << (atomInfo.get<int>(PdbKey::serial) % 100000); /* atom serial number */
|
---|
[16462f] | 460 | *file << " "; /* char 12 is empty */
|
---|
| 461 | *file << setfill(' ') << left << setw(4)
|
---|
| 462 | << atomInfo.get<std::string>(PdbKey::name); /* atom name */
|
---|
| 463 | *file << setfill(' ') << left << setw(1)
|
---|
| 464 | << atomInfo.get<std::string>(PdbKey::altLoc); /* alternate location/conformation */
|
---|
| 465 | *file << setfill(' ') << left << setw(3)
|
---|
| 466 | << atomInfo.get<std::string>(PdbKey::resName); /* residue name */
|
---|
| 467 | *file << " "; /* char 21 is empty */
|
---|
| 468 | *file << setfill(' ') << left << setw(1)
|
---|
| 469 | << atomInfo.get<std::string>(PdbKey::chainID); /* chain identifier */
|
---|
| 470 | *file << setfill(' ') << left << setw(4)
|
---|
[1f5e97] | 471 | << (atomInfo.get<int>(PdbKey::resSeq) % 10000); /* residue sequence number */
|
---|
[16462f] | 472 | *file << setfill(' ') << left << setw(1)
|
---|
| 473 | << atomInfo.get<std::string>(PdbKey::iCode); /* iCode */
|
---|
| 474 | *file << " "; /* char 28-30 are empty */
|
---|
| 475 | // have the following operate on stringstreams such that format specifiers
|
---|
| 476 | // only act on these
|
---|
| 477 | for (size_t i=0;i<NDIM;++i) {
|
---|
[105b72] | 478 | *file << setfill(' ') << right << setw(8)
|
---|
| 479 | << printCoordinate(atomInfo.get<double>(PositionEnumMap[i]));
|
---|
[16462f] | 480 | }
|
---|
| 481 | {
|
---|
| 482 | stringstream occupancy;
|
---|
| 483 | occupancy << fixed << setprecision(2) << showpoint
|
---|
| 484 | << atomInfo.get<double>(PdbKey::occupancy); /* occupancy */
|
---|
| 485 | *file << setfill(' ') << right << setw(6) << occupancy.str();
|
---|
[3ae731] | 486 | }
|
---|
[16462f] | 487 | {
|
---|
| 488 | stringstream tempFactor;
|
---|
| 489 | tempFactor << fixed << setprecision(2) << showpoint
|
---|
| 490 | << atomInfo.get<double>(PdbKey::tempFactor); /* temperature factor */
|
---|
| 491 | *file << setfill(' ') << right << setw(6) << tempFactor.str();
|
---|
| 492 | }
|
---|
| 493 | *file << " "; /* char 68-76 are empty */
|
---|
| 494 | *file << setfill(' ') << right << setw(2) << atomInfo.get<std::string>(PdbKey::element); /* element */
|
---|
| 495 | *file << setfill(' ') << right << setw(2) << atomInfo.get<int>(PdbKey::charge); /* charge */
|
---|
[3ae731] | 496 |
|
---|
| 497 | *file << endl;
|
---|
| 498 | }
|
---|
| 499 |
|
---|
| 500 | /**
|
---|
| 501 | * Writes the neighbor information of one atom to the provided stream.
|
---|
| 502 | *
|
---|
[9d83b6] | 503 | * Note that ListOfBonds of WorldTime::CurrentTime is used.
|
---|
| 504 | *
|
---|
[473237] | 505 | * Also, we fill up the CONECT line to extend over 80 chars.
|
---|
| 506 | *
|
---|
[bb6193] | 507 | * \param *file where to write neighbor information to
|
---|
| 508 | * \param MaxnumberOfNeighbors of neighbors
|
---|
| 509 | * \param *currentAtom to the atom of which to take the neighbor information
|
---|
[3ae731] | 510 | */
|
---|
[fac58f] | 511 | void FormatParser< pdb >::writeNeighbors(
|
---|
| 512 | ostream* file,
|
---|
| 513 | int MaxnumberOfNeighbors,
|
---|
| 514 | const atom * const currentAtom) {
|
---|
[4c1230] | 515 | int MaxNo = MaxnumberOfNeighbors;
|
---|
[473237] | 516 | int charsleft = 80;
|
---|
[9d83b6] | 517 | const BondList & ListOfBonds = currentAtom->getListOfBonds();
|
---|
| 518 | if (!ListOfBonds.empty()) {
|
---|
| 519 | for(BondList::const_iterator currentBond = ListOfBonds.begin(); currentBond != ListOfBonds.end(); ++currentBond) {
|
---|
[4c1230] | 520 | if (MaxNo >= MaxnumberOfNeighbors) {
|
---|
| 521 | *file << "CONECT";
|
---|
[c0e28c] | 522 | *file << setw(5) << getLocalId(currentAtom->getId());
|
---|
[473237] | 523 | charsleft = 80-6-5;
|
---|
[4c1230] | 524 | MaxNo = 0;
|
---|
[bb6193] | 525 | }
|
---|
[c0e28c] | 526 | *file << setw(5) << getLocalId((*currentBond)->GetOtherAtom(currentAtom)->getId());
|
---|
[473237] | 527 | charsleft -= 5;
|
---|
[bb6193] | 528 | MaxNo++;
|
---|
[473237] | 529 | if (MaxNo == MaxnumberOfNeighbors) {
|
---|
| 530 | for (;charsleft > 0; charsleft--)
|
---|
| 531 | *file << ' ';
|
---|
[4c1230] | 532 | *file << "\n";
|
---|
[473237] | 533 | }
|
---|
[3ae731] | 534 | }
|
---|
[473237] | 535 | if (MaxNo != MaxnumberOfNeighbors) {
|
---|
| 536 | for (;charsleft > 0; charsleft--)
|
---|
| 537 | *file << ' ';
|
---|
[4c1230] | 538 | *file << "\n";
|
---|
[473237] | 539 | }
|
---|
[3ae731] | 540 | }
|
---|
| 541 | }
|
---|
| 542 |
|
---|
[9dba5f] | 543 | /** Either returns present atom with given id or a newly created one.
|
---|
| 544 | *
|
---|
| 545 | * @param id_string
|
---|
| 546 | * @return
|
---|
| 547 | */
|
---|
[c0e28c] | 548 | atom* FormatParser< pdb >::getAtomToParse(std::string id_string)
|
---|
[9dba5f] | 549 | {
|
---|
| 550 | // get the local ID
|
---|
| 551 | ConvertTo<int> toInt;
|
---|
[c0e28c] | 552 | const unsigned int AtomID_local = toInt(id_string);
|
---|
| 553 | LOG(4, "INFO: Local id is "+toString(AtomID_local)+".");
|
---|
[9dba5f] | 554 | // get the atomic ID if present
|
---|
| 555 | atom* newAtom = NULL;
|
---|
[c0e28c] | 556 | if (getGlobalId(AtomID_local) != -1) {
|
---|
| 557 | const unsigned int AtomID_global = getGlobalId(AtomID_local);
|
---|
| 558 | LOG(4, "INFO: Global id present as " << AtomID_global << ".");
|
---|
[9dba5f] | 559 | // check if atom exists
|
---|
[c0e28c] | 560 | newAtom = World::getInstance().getAtom(AtomById(AtomID_global));
|
---|
[ffa69c] | 561 | if (DoLog(5)) {
|
---|
| 562 | LOG(5, "INFO: Listing all present atoms with id.");
|
---|
| 563 | BOOST_FOREACH(const atom *_atom, const_cast<const World &>(World::getInstance()).getAllAtoms())
|
---|
| 564 | LOG(5, "INFO: " << *_atom << " with id " << _atom->getId());
|
---|
| 565 | }
|
---|
[9dba5f] | 566 | }
|
---|
| 567 | // if not exists, create
|
---|
| 568 | if (newAtom == NULL) {
|
---|
| 569 | newAtom = World::getInstance().createAtom();
|
---|
[8bf9c6] | 570 | //const unsigned int AtomID_global = newAtom->getId();
|
---|
[9dba5f] | 571 | LOG(4, "INFO: No association to global id present, creating atom.");
|
---|
| 572 | } else {
|
---|
| 573 | LOG(4, "INFO: Existing atom found: " << *newAtom << ".");
|
---|
| 574 | }
|
---|
| 575 | return newAtom;
|
---|
| 576 | }
|
---|
| 577 |
|
---|
[5fa2ba] | 578 | /** read a line starting with key ATOM.
|
---|
| 579 | *
|
---|
| 580 | * We check for line's length and parse only up to this value.
|
---|
| 581 | *
|
---|
| 582 | * @param atomInfo container to put information in
|
---|
| 583 | * @param line line containing key ATOM
|
---|
| 584 | */
|
---|
[765f16] | 585 | void FormatParser< pdb >::readPdbAtomInfoContainer(PdbAtomInfoContainer &atomInfo, std::string &line) const
|
---|
[9dba5f] | 586 | {
|
---|
[5fa2ba] | 587 | const size_t length = line.length();
|
---|
| 588 | if (length < 80)
|
---|
[765f16] | 589 | ELOG(2, "FormatParser< pdb >::readPdbAtomInfoContainer() - pdb is mal-formed, containing less than 80 chars!");
|
---|
[5fa2ba] | 590 | if (length >= 6) {
|
---|
| 591 | LOG(4,"INFO: Parsing token from "+line.substr(0,6)+".");
|
---|
| 592 | atomInfo.set(PdbKey::token, line.substr(0,6));
|
---|
| 593 | }
|
---|
| 594 | if (length >= 11) {
|
---|
| 595 | LOG(4,"INFO: Parsing serial from "+line.substr(6,5)+".");
|
---|
| 596 | atomInfo.set(PdbKey::serial, line.substr(6,5));
|
---|
| 597 | ASSERT(atomInfo.get<int>(PdbKey::serial) != 0,
|
---|
[765f16] | 598 | "FormatParser< pdb >::readPdbAtomInfoContainer() - serial 0 is invalid (filler id for conect entries).");
|
---|
[5fa2ba] | 599 | }
|
---|
| 600 |
|
---|
| 601 | if (length >= 16) {
|
---|
| 602 | LOG(4,"INFO: Parsing name from "+line.substr(12,4)+".");
|
---|
| 603 | atomInfo.set(PdbKey::name, line.substr(12,4));
|
---|
| 604 | }
|
---|
| 605 | if (length >= 17) {
|
---|
| 606 | LOG(4,"INFO: Parsing altLoc from "+line.substr(16,1)+".");
|
---|
| 607 | atomInfo.set(PdbKey::altLoc, line.substr(16,1));
|
---|
| 608 | }
|
---|
| 609 | if (length >= 20) {
|
---|
| 610 | LOG(4,"INFO: Parsing resName from "+line.substr(17,3)+".");
|
---|
| 611 | atomInfo.set(PdbKey::resName, line.substr(17,3));
|
---|
| 612 | }
|
---|
| 613 | if (length >= 22) {
|
---|
| 614 | LOG(4,"INFO: Parsing chainID from "+line.substr(21,1)+".");
|
---|
| 615 | atomInfo.set(PdbKey::chainID, line.substr(21,1));
|
---|
| 616 | }
|
---|
| 617 | if (length >= 26) {
|
---|
| 618 | LOG(4,"INFO: Parsing resSeq from "+line.substr(22,4)+".");
|
---|
| 619 | atomInfo.set(PdbKey::resSeq, line.substr(22,4));
|
---|
| 620 | }
|
---|
| 621 | if (length >= 27) {
|
---|
| 622 | LOG(4,"INFO: Parsing iCode from "+line.substr(26,1)+".");
|
---|
| 623 | atomInfo.set(PdbKey::iCode, line.substr(26,1));
|
---|
| 624 | }
|
---|
| 625 |
|
---|
| 626 | if (length >= 60) {
|
---|
| 627 | LOG(4,"INFO: Parsing occupancy from "+line.substr(54,6)+".");
|
---|
| 628 | atomInfo.set(PdbKey::occupancy, line.substr(54,6));
|
---|
| 629 | }
|
---|
| 630 | if (length >= 66) {
|
---|
| 631 | LOG(4,"INFO: Parsing tempFactor from "+line.substr(60,6)+".");
|
---|
| 632 | atomInfo.set(PdbKey::tempFactor, line.substr(60,6));
|
---|
| 633 | }
|
---|
| 634 | if (length >= 80) {
|
---|
| 635 | LOG(4,"INFO: Parsing charge from "+line.substr(78,2)+".");
|
---|
| 636 | atomInfo.set(PdbKey::charge, line.substr(78,2));
|
---|
| 637 | }
|
---|
| 638 | if (length >= 78) {
|
---|
| 639 | LOG(4,"INFO: Parsing element from "+line.substr(76,2)+".");
|
---|
| 640 | atomInfo.set(PdbKey::element, line.substr(76,2));
|
---|
| 641 | } else {
|
---|
| 642 | LOG(4,"INFO: Trying to parse alternative element from name "+line.substr(12,4)+".");
|
---|
| 643 | atomInfo.set(PdbKey::element, line.substr(12,4));
|
---|
| 644 | }
|
---|
[9dba5f] | 645 | }
|
---|
| 646 |
|
---|
[4fbca9c] | 647 | /** Parse an ATOM line from a PDB file.
|
---|
| 648 | *
|
---|
| 649 | * Reads one data line of a pdstatus file and interprets it according to the
|
---|
| 650 | * specifications of the PDB 3.2 format: http://www.wwpdb.org/docs.html
|
---|
| 651 | *
|
---|
| 652 | * A new atom is created and filled with available information, non-
|
---|
| 653 | * standard information is placed in additionalAtomData at the atom's id.
|
---|
[3ae731] | 654 | *
|
---|
[b0a2e3] | 655 | * \param _step time step to use
|
---|
[3ae731] | 656 | * \param line to parse as an atom
|
---|
[4fbca9c] | 657 | * \param newmol molecule to add parsed atoms to
|
---|
[3ae731] | 658 | */
|
---|
[765f16] | 659 | void FormatParser< pdb >::readAtomDataLine(const unsigned int _step, std::string &line, molecule *newmol = NULL) {
|
---|
[4fbca9c] | 660 | vector<string>::iterator it;
|
---|
[9dba5f] | 661 |
|
---|
| 662 | atom* newAtom = getAtomToParse(line.substr(6,5));
|
---|
| 663 | LOG(3,"INFO: Parsing END entry or empty line.");
|
---|
| 664 | bool FirstTimestep = isPresentadditionalAtomData(newAtom->getId()) ? false : true;
|
---|
[b0a2e3] | 665 | ASSERT((FirstTimestep && (_step == 0)) || (!FirstTimestep && (_step !=0)),
|
---|
[765f16] | 666 | "FormatParser< pdb >::readAtomDataLine() - additionalAtomData present though atom is newly parsed.");
|
---|
[9dba5f] | 667 | if (FirstTimestep) {
|
---|
[c0e28c] | 668 | LOG(3,"INFO: Parsing new atom "+toString(*newAtom)+" "+toString(newAtom->getId())+".");
|
---|
[9dba5f] | 669 | } else {
|
---|
| 670 | LOG(3,"INFO: Parsing present atom "+toString(*newAtom)+".");
|
---|
| 671 | }
|
---|
[93fd43e] | 672 | PdbAtomInfoContainer &atomInfo = getadditionalAtomData(newAtom);
|
---|
[9dba5f] | 673 | LOG(4,"INFO: Information in container is "+toString(atomInfo)+".");
|
---|
| 674 |
|
---|
[4fbca9c] | 675 | string word;
|
---|
| 676 | ConvertTo<size_t> toSize_t;
|
---|
| 677 |
|
---|
| 678 | // check whether serial exists, if so, assign next available
|
---|
| 679 |
|
---|
[47d041] | 680 | // LOG(2, "Split line:"
|
---|
[4fbca9c] | 681 | // << line.substr(6,5) << "|"
|
---|
| 682 | // << line.substr(12,4) << "|"
|
---|
| 683 | // << line.substr(16,1) << "|"
|
---|
| 684 | // << line.substr(17,3) << "|"
|
---|
| 685 | // << line.substr(21,1) << "|"
|
---|
| 686 | // << line.substr(22,4) << "|"
|
---|
| 687 | // << line.substr(26,1) << "|"
|
---|
| 688 | // << line.substr(30,8) << "|"
|
---|
| 689 | // << line.substr(38,8) << "|"
|
---|
| 690 | // << line.substr(46,8) << "|"
|
---|
| 691 | // << line.substr(54,6) << "|"
|
---|
| 692 | // << line.substr(60,6) << "|"
|
---|
| 693 | // << line.substr(76,2) << "|"
|
---|
[47d041] | 694 | // << line.substr(78,2));
|
---|
[4fbca9c] | 695 |
|
---|
[9dba5f] | 696 | if (FirstTimestep) {
|
---|
| 697 | // first time step
|
---|
| 698 | // then fill info container
|
---|
| 699 | readPdbAtomInfoContainer(atomInfo, line);
|
---|
[c0e28c] | 700 | // associate local with global id
|
---|
| 701 | associateLocaltoGlobalId(toSize_t(atomInfo.get<std::string>(PdbKey::serial)), newAtom->getId());
|
---|
[9dba5f] | 702 | // set position
|
---|
| 703 | Vector tempVector;
|
---|
| 704 | LOG(4,"INFO: Parsing position from ("
|
---|
| 705 | +line.substr(30,8)+","
|
---|
| 706 | +line.substr(38,8)+","
|
---|
| 707 | +line.substr(46,8)+").");
|
---|
| 708 | PdbAtomInfoContainer::ScanKey(tempVector[0], line.substr(30,8));
|
---|
| 709 | PdbAtomInfoContainer::ScanKey(tempVector[1], line.substr(38,8));
|
---|
| 710 | PdbAtomInfoContainer::ScanKey(tempVector[2], line.substr(46,8));
|
---|
| 711 | newAtom->setPosition(tempVector);
|
---|
[797a40] | 712 | // set particle name
|
---|
| 713 | newAtom->setParticleName(atomInfo.get<std::string>(PdbKey::name));
|
---|
[9dba5f] | 714 | // set element
|
---|
[8990879] | 715 | std::string value = atomInfo.get<std::string>(PdbKey::element);
|
---|
| 716 | // make second character lower case if not
|
---|
| 717 | if ((value[1] >= 'A') && (value[1] <= 'Z'))
|
---|
| 718 | value[1] = (value[1] - 'A') + 'a';
|
---|
[9dba5f] | 719 | const element *elem = World::getInstance().getPeriode()
|
---|
[8990879] | 720 | ->FindElement(value);
|
---|
[9dba5f] | 721 | ASSERT(elem != NULL,
|
---|
[765f16] | 722 | "FormatParser< pdb >::readAtomDataLine() - element "+atomInfo.get<std::string>(PdbKey::element)+" is unknown!");
|
---|
[9dba5f] | 723 | newAtom->setType(elem);
|
---|
| 724 |
|
---|
| 725 | if (newmol != NULL)
|
---|
| 726 | newmol->AddAtom(newAtom);
|
---|
| 727 | } else {
|
---|
| 728 | // not first time step
|
---|
| 729 | // then parse into different container
|
---|
| 730 | PdbAtomInfoContainer consistencyInfo;
|
---|
| 731 | readPdbAtomInfoContainer(consistencyInfo, line);
|
---|
| 732 | // then check additional info for consistency
|
---|
| 733 | ASSERT(atomInfo.get<std::string>(PdbKey::token) == consistencyInfo.get<std::string>(PdbKey::token),
|
---|
[765f16] | 734 | "FormatParser< pdb >::readAtomDataLine() - difference in token on multiple time step for atom with id "
|
---|
[9dba5f] | 735 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 736 | ASSERT(atomInfo.get<std::string>(PdbKey::name) == consistencyInfo.get<std::string>(PdbKey::name),
|
---|
[765f16] | 737 | "FormatParser< pdb >::readAtomDataLine() - difference in name on multiple time step for atom with id "
|
---|
[9dba5f] | 738 | +atomInfo.get<std::string>(PdbKey::serial)+":"
|
---|
| 739 | +atomInfo.get<std::string>(PdbKey::name)+"!="
|
---|
| 740 | +consistencyInfo.get<std::string>(PdbKey::name)+".");
|
---|
| 741 | ASSERT(atomInfo.get<std::string>(PdbKey::altLoc) == consistencyInfo.get<std::string>(PdbKey::altLoc),
|
---|
[765f16] | 742 | "FormatParser< pdb >::readAtomDataLine() - difference in altLoc on multiple time step for atom with id "
|
---|
[9dba5f] | 743 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 744 | ASSERT(atomInfo.get<std::string>(PdbKey::resName) == consistencyInfo.get<std::string>(PdbKey::resName),
|
---|
[765f16] | 745 | "FormatParser< pdb >::readAtomDataLine() - difference in resName on multiple time step for atom with id "
|
---|
[9dba5f] | 746 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 747 | ASSERT(atomInfo.get<std::string>(PdbKey::chainID) == consistencyInfo.get<std::string>(PdbKey::chainID),
|
---|
[765f16] | 748 | "FormatParser< pdb >::readAtomDataLine() - difference in chainID on multiple time step for atom with id "
|
---|
[9dba5f] | 749 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 750 | ASSERT(atomInfo.get<std::string>(PdbKey::resSeq) == consistencyInfo.get<std::string>(PdbKey::resSeq),
|
---|
[765f16] | 751 | "FormatParser< pdb >::readAtomDataLine() - difference in resSeq on multiple time step for atom with id "
|
---|
[9dba5f] | 752 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 753 | ASSERT(atomInfo.get<std::string>(PdbKey::iCode) == consistencyInfo.get<std::string>(PdbKey::iCode),
|
---|
[765f16] | 754 | "FormatParser< pdb >::readAtomDataLine() - difference in iCode on multiple time step for atom with id "
|
---|
[9dba5f] | 755 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 756 | ASSERT(atomInfo.get<std::string>(PdbKey::occupancy) == consistencyInfo.get<std::string>(PdbKey::occupancy),
|
---|
[765f16] | 757 | "FormatParser< pdb >::readAtomDataLine() - difference in occupancy on multiple time step for atom with id "
|
---|
[9dba5f] | 758 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 759 | ASSERT(atomInfo.get<std::string>(PdbKey::tempFactor) == consistencyInfo.get<std::string>(PdbKey::tempFactor),
|
---|
[765f16] | 760 | "FormatParser< pdb >::readAtomDataLine() - difference in tempFactor on multiple time step for atom with id "
|
---|
[9dba5f] | 761 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 762 | ASSERT(atomInfo.get<std::string>(PdbKey::charge) == consistencyInfo.get<std::string>(PdbKey::charge),
|
---|
[765f16] | 763 | "FormatParser< pdb >::readAtomDataLine() - difference in charge on multiple time step for atom with id "
|
---|
[9dba5f] | 764 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 765 | ASSERT(atomInfo.get<std::string>(PdbKey::element) == consistencyInfo.get<std::string>(PdbKey::element),
|
---|
[765f16] | 766 | "FormatParser< pdb >::readAtomDataLine() - difference in element on multiple time step for atom with id "
|
---|
[9dba5f] | 767 | +atomInfo.get<std::string>(PdbKey::serial)+"!");
|
---|
| 768 | // and parse in trajectory
|
---|
| 769 | Vector tempVector;
|
---|
| 770 | LOG(4,"INFO: Parsing trajectory position from ("
|
---|
| 771 | +line.substr(30,8)+","
|
---|
| 772 | +line.substr(38,8)+","
|
---|
| 773 | +line.substr(46,8)+").");
|
---|
| 774 | PdbAtomInfoContainer::ScanKey(tempVector[0], line.substr(30,8));
|
---|
| 775 | PdbAtomInfoContainer::ScanKey(tempVector[1], line.substr(38,8));
|
---|
| 776 | PdbAtomInfoContainer::ScanKey(tempVector[2], line.substr(46,8));
|
---|
[b0a2e3] | 777 | LOG(4,"INFO: Adding trajectory point to time step "+toString(_step)+".");
|
---|
[9dba5f] | 778 | // and set position at new time step
|
---|
[b0a2e3] | 779 | newAtom->setPositionAtStep(_step, tempVector);
|
---|
[9dba5f] | 780 | }
|
---|
| 781 |
|
---|
[4fbca9c] | 782 |
|
---|
| 783 | // printAtomInfo(newAtom);
|
---|
[3ae731] | 784 | }
|
---|
| 785 |
|
---|
[4fbca9c] | 786 | /** Prints all PDB-specific information known about an atom.
|
---|
[3ae731] | 787 | *
|
---|
| 788 | */
|
---|
[765f16] | 789 | void FormatParser< pdb >::printAtomInfo(const atom * const newAtom) const
|
---|
[4fbca9c] | 790 | {
|
---|
| 791 | const PdbAtomInfoContainer &atomInfo = additionalAtomData.at(newAtom->getId()); // operator[] const does not exist
|
---|
| 792 |
|
---|
[47d041] | 793 | LOG(1, "We know about atom " << newAtom->getId() << ":");
|
---|
| 794 | LOG(1, "\ttoken is " << atomInfo.get<std::string>(PdbKey::token));
|
---|
| 795 | LOG(1, "\tserial is " << atomInfo.get<int>(PdbKey::serial));
|
---|
| 796 | LOG(1, "\tname is " << atomInfo.get<std::string>(PdbKey::name));
|
---|
| 797 | LOG(1, "\taltLoc is " << atomInfo.get<std::string>(PdbKey::altLoc));
|
---|
| 798 | LOG(1, "\tresName is " << atomInfo.get<std::string>(PdbKey::resName));
|
---|
| 799 | LOG(1, "\tchainID is " << atomInfo.get<std::string>(PdbKey::chainID));
|
---|
| 800 | LOG(1, "\tresSeq is " << atomInfo.get<int>(PdbKey::resSeq));
|
---|
| 801 | LOG(1, "\tiCode is " << atomInfo.get<std::string>(PdbKey::iCode));
|
---|
| 802 | LOG(1, "\tX is " << atomInfo.get<double>(PdbKey::X));
|
---|
| 803 | LOG(1, "\tY is " << atomInfo.get<double>(PdbKey::Y));
|
---|
| 804 | LOG(1, "\tZ is " << atomInfo.get<double>(PdbKey::Z));
|
---|
| 805 | LOG(1, "\toccupancy is " << atomInfo.get<double>(PdbKey::occupancy));
|
---|
| 806 | LOG(1, "\ttempFactor is " << atomInfo.get<double>(PdbKey::tempFactor));
|
---|
| 807 | LOG(1, "\telement is '" << *(newAtom->getType()) << "'");
|
---|
| 808 | LOG(1, "\tcharge is " << atomInfo.get<int>(PdbKey::charge));
|
---|
[3ae731] | 809 | }
|
---|
| 810 |
|
---|
| 811 | /**
|
---|
[4fbca9c] | 812 | * Reads neighbor information for one atom from the input.
|
---|
| 813 | *
|
---|
[b0a2e3] | 814 | * \param _step time step to use
|
---|
[4fbca9c] | 815 | * \param line to parse as an atom
|
---|
[3ae731] | 816 | */
|
---|
[765f16] | 817 | void FormatParser< pdb >::readNeighbors(const unsigned int _step, std::string &line)
|
---|
[4fbca9c] | 818 | {
|
---|
| 819 | const size_t length = line.length();
|
---|
| 820 | std::list<size_t> ListOfNeighbors;
|
---|
| 821 | ConvertTo<size_t> toSize_t;
|
---|
| 822 |
|
---|
| 823 | // obtain neighbours
|
---|
| 824 | // show split line for debugging
|
---|
| 825 | string output;
|
---|
| 826 | ASSERT(length >=16,
|
---|
[765f16] | 827 | "FormatParser< pdb >::readNeighbors() - CONECT entry has not enough entries: "+line+"!");
|
---|
| 828 | // output = "Split line:|";
|
---|
| 829 | // output += line.substr(6,5) + "|";
|
---|
[4fbca9c] | 830 | const size_t id = toSize_t(line.substr(6,5));
|
---|
| 831 | for (size_t index = 11; index <= 26; index+=5) {
|
---|
| 832 | if (index+5 <= length) {
|
---|
[473237] | 833 | output += line.substr(index,5) + "|";
|
---|
| 834 | // search for digits
|
---|
| 835 | int otherid = -1;
|
---|
| 836 | PdbAtomInfoContainer::ScanKey(otherid, line.substr(index,5));
|
---|
[5fa2ba] | 837 | if (otherid > 0)
|
---|
| 838 | ListOfNeighbors.push_back(otherid);
|
---|
| 839 | else
|
---|
[44f53e] | 840 | ELOG(3, "FormatParser< pdb >::readNeighbors() - discarding CONECT entry with id 0.");
|
---|
[4fbca9c] | 841 | } else {
|
---|
| 842 | break;
|
---|
| 843 | }
|
---|
| 844 | }
|
---|
[473237] | 845 | LOG(4, output);
|
---|
[4fbca9c] | 846 |
|
---|
| 847 | // add neighbours
|
---|
[f01769] | 848 | atom * const _atom = World::getInstance().getAtom(AtomById(getGlobalId(id)));
|
---|
[473237] | 849 | LOG(2, "STATUS: Atom " << _atom->getId() << " gets " << ListOfNeighbors.size() << " more neighbours.");
|
---|
[4fbca9c] | 850 | for (std::list<size_t>::const_iterator iter = ListOfNeighbors.begin();
|
---|
| 851 | iter != ListOfNeighbors.end();
|
---|
| 852 | ++iter) {
|
---|
[c0e28c] | 853 | atom * const _Otheratom = World::getInstance().getAtom(AtomById(getGlobalId(*iter)));
|
---|
[473237] | 854 | LOG(3, "INFO: Adding Bond (" << *_atom << "," << *_Otheratom << ")");
|
---|
[b0a2e3] | 855 | _atom->addBond(_step, _Otheratom);
|
---|
[4fbca9c] | 856 | }
|
---|
[3ae731] | 857 | }
|
---|
| 858 |
|
---|
[765f16] | 859 | bool FormatParser< pdb >::operator==(const FormatParser< pdb >& b) const
|
---|
[4fbca9c] | 860 | {
|
---|
| 861 | bool status = true;
|
---|
[a58c16] | 862 | World::ConstAtomComposite atoms = const_cast<const World &>(World::getInstance()).
|
---|
| 863 | getAllAtoms();
|
---|
| 864 | for (World::ConstAtomComposite::const_iterator iter = atoms.begin(); iter != atoms.end(); ++iter) {
|
---|
[4fbca9c] | 865 | if ((additionalAtomData.find((*iter)->getId()) != additionalAtomData.end())
|
---|
| 866 | && (b.additionalAtomData.find((*iter)->getId()) != b.additionalAtomData.end())) {
|
---|
| 867 | const PdbAtomInfoContainer &atomInfo = additionalAtomData.at((*iter)->getId());
|
---|
| 868 | const PdbAtomInfoContainer &OtheratomInfo = b.additionalAtomData.at((*iter)->getId());
|
---|
| 869 |
|
---|
[16462f] | 870 | status = status && (atomInfo.get<std::string>(PdbKey::serial) == OtheratomInfo.get<std::string>(PdbKey::serial));
|
---|
[47d041] | 871 | if (!status) ELOG(1, "Mismatch in serials!");
|
---|
[16462f] | 872 | status = status && (atomInfo.get<std::string>(PdbKey::name) == OtheratomInfo.get<std::string>(PdbKey::name));
|
---|
[47d041] | 873 | if (!status) ELOG(1, "Mismatch in names!");
|
---|
[16462f] | 874 | status = status && (atomInfo.get<std::string>(PdbKey::altLoc) == OtheratomInfo.get<std::string>(PdbKey::altLoc));
|
---|
[47d041] | 875 | if (!status) ELOG(1, "Mismatch in altLocs!");
|
---|
[16462f] | 876 | status = status && (atomInfo.get<std::string>(PdbKey::resName) == OtheratomInfo.get<std::string>(PdbKey::resName));
|
---|
[47d041] | 877 | if (!status) ELOG(1, "Mismatch in resNames!");
|
---|
[16462f] | 878 | status = status && (atomInfo.get<std::string>(PdbKey::chainID) == OtheratomInfo.get<std::string>(PdbKey::chainID));
|
---|
[47d041] | 879 | if (!status) ELOG(1, "Mismatch in chainIDs!");
|
---|
[16462f] | 880 | status = status && (atomInfo.get<std::string>(PdbKey::resSeq) == OtheratomInfo.get<std::string>(PdbKey::resSeq));
|
---|
[47d041] | 881 | if (!status) ELOG(1, "Mismatch in resSeqs!");
|
---|
[16462f] | 882 | status = status && (atomInfo.get<std::string>(PdbKey::iCode) == OtheratomInfo.get<std::string>(PdbKey::iCode));
|
---|
[47d041] | 883 | if (!status) ELOG(1, "Mismatch in iCodes!");
|
---|
[16462f] | 884 | status = status && (atomInfo.get<std::string>(PdbKey::occupancy) == OtheratomInfo.get<std::string>(PdbKey::occupancy));
|
---|
[47d041] | 885 | if (!status) ELOG(1, "Mismatch in occupancies!");
|
---|
[16462f] | 886 | status = status && (atomInfo.get<std::string>(PdbKey::tempFactor) == OtheratomInfo.get<std::string>(PdbKey::tempFactor));
|
---|
[47d041] | 887 | if (!status) ELOG(1, "Mismatch in tempFactors!");
|
---|
[16462f] | 888 | status = status && (atomInfo.get<std::string>(PdbKey::charge) == OtheratomInfo.get<std::string>(PdbKey::charge));
|
---|
[47d041] | 889 | if (!status) ELOG(1, "Mismatch in charges!");
|
---|
[4fbca9c] | 890 | }
|
---|
[3ae731] | 891 | }
|
---|
| 892 |
|
---|
[4fbca9c] | 893 | return status;
|
---|
[3ae731] | 894 | }
|
---|
| 895 |
|
---|