/usr/include/votca/csg/topology.h is in libvotca-csg2-dev 1.2.4-2.1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 | /*
* Copyright 2009-2011 The VOTCA Development Team (http://www.votca.org)
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*/
#ifndef _topology_H
#define _topology_H
#include <vector>
#include <map>
#include <list>
#include <assert.h>
#include <votca/tools/types.h>
#include <votca/tools/vec.h>
#include <votca/tools/matrix.h>
#include "exclusionlist.h"
#include "bead.h"
#include "molecule.h"
#include "residue.h"
#include "beadtype.h"
#include "boundarycondition.h"
#include "triclinicbox.h"
#include "orthorhombicbox.h"
#include "openbox.h"
namespace votca { namespace csg {
using namespace votca::tools;
class Interaction;
class ExclusionList;
typedef vector<Molecule *> MoleculeContainer;
typedef vector<Bead *> BeadContainer;
typedef vector<BeadType *> BeadTypeContainer;
typedef vector<Residue *> ResidueContainer;
typedef vector<Interaction *> InteractionContainer;
using namespace std;
/**
\brief topology of the whole system
The Topology class stores the topology of the system like the beads, bonds, molecules and residues.
\todo internal management for ids and indices
**/
class Topology
{
public:
/// constructor
Topology() { _bc = new OpenBox(); }
virtual ~Topology();
/**
* \brief cleans up all the stored data
*/
virtual void Cleanup();
/**
* \brief creates a new Bead
*
* \param symmetry symmetry of the bead, 1: spherical 3: ellipsoidal
* \param name name of the bead
* \param type bead type
* \param resnr residue number
* \param m mass
* \param q charge
* \return pointer to created bead
*
* The function creates a new bead and adds it to the list of beads.
*/
virtual Bead *CreateBead(byte_t symmetry, string name, BeadType *type, int resnr, double m, double q);
/**
* \brief get bead type or create it
* \param name typename
* \return pointer to bead type
*
* Returns an existing bead type or creates one if it doesn't exist yet
*/
virtual BeadType *GetOrCreateBeadType(string name);
/**
* \brief creates a new molecule
* \param name name of the molecule
* \return pointer to created molecule
*/
virtual Molecule *CreateMolecule(string name);
/**
* \brief checks weather molecules with the same name really contain the same number of beads
*/
void CheckMoleculeNaming(void);
/**
* \brief create a new resiude
* @param name residue name
* @return created residue
*/
virtual Residue *CreateResidue(string name);
/**
* \brief create molecules based on the residue
*
* This function scans the topology and creates molecules based on the resiude id.
* All beads with the same resid are put int one molecule.
*/
void CreateMoleculesByResidue();
/**
* \brief put the whole topology in one molecule
* \param name name of the new molecule
*
* This function creates one big molecule for all beads in the topology.
*/
void CreateOneBigMolecule(string name);
/**
* \brief create molecules based on blocks of atoms
* \param name molecule name
* \param first first bead
* \param nbeads number of beads per molecule
* \param nmolecules number of molecules
*/
void CreateMoleculesByRange(string name, int first, int nbeads, int nmolecules);
/**
* \brief number of molecules in the system
* @return number of molecule in topology
*/
int MoleculeCount() { return _molecules.size(); }
/**
* number of beads in the system
* @return number of beads in the system
*/
int BeadCount() { return _beads.size(); }
/**
* number of residues in the system
* \return number of residues
*/
int ResidueCount() { return _residues.size(); }
/**
* get molecule by index
* @param index molecule number
* @return pointer to molecule
*/
Molecule *MoleculeByIndex(int index);
/**
* access containter with all beads
* @return bead container
*/
BeadContainer &Beads() { return _beads; }
/**
* access containter with all molecules
* @return molecule container
*/
MoleculeContainer &Molecules() { return _molecules; }
/**
* access containter with all bonded interactions
* @return bonded interaction container
*/
InteractionContainer &BondedInteractions() { return _interactions; }
void AddBondedInteraction(Interaction *ic);
std::list<Interaction *> InteractionsInGroup(const string &group);
BeadType *getBeadType(const int i) const { return _beadtypes[i]; }
Bead *getBead(const int i) const { return _beads[i]; }
Residue *getResidue(const int i) const { return _residues[i]; }
Molecule *getMolecule(const int i) const { return _molecules[i]; }
/**
* delete all molecule information
*/
void ClearMoleculeList(){
_molecules.clear();
}
/**
* \brief adds all the beads+molecules+residues from other topology
* \param top topology to add
*/
void Add(Topology *top);
/**
* \brief copy topology data of different topology
* \param top topology to copy from
*/
void CopyTopologyData(Topology *top);
/**
* \brief rename all the molecules in range
* \param range range string of type 1:2:10 = 1, 3, 5, 7, ...
* \param name new name of molecule
* range is a string which is parsed by RangeParser,
*/
void RenameMolecules(string range, string name);
/**
* set the simulation box
* \param box triclinic box matrix
*/
void setBox(const matrix &box, BoundaryCondition::eBoxtype boxtype=BoundaryCondition::typeAuto) {
// determine box type automatically in case boxtype==typeAuto
if(boxtype==BoundaryCondition::typeAuto) {
boxtype = autoDetectBoxType(box);
}
if(_bc) {
delete (_bc);
}
switch(boxtype) {
case BoundaryCondition::typeTriclinic:
_bc = new TriclinicBox();
break;
case BoundaryCondition::typeOrthorhombic:
_bc = new OrthorhombicBox();
break;
default:
_bc = new OpenBox();
break;
}
_bc->setBox(box);
};
/**
* get the simulation box
* \return triclinic box matrix
*/
const matrix &getBox() {
return _bc->getBox();
};
/**
* set the time of current frame
* \param t simulation time in ns
*/
void setTime(double t) { _time = t; };
/**
* get the time of current frame
* \return simulation time in ns
*/
double getTime() { return _time; };
/**
* set the step number of current frame
* \param s step number
*/
void setStep(int s) { _step = s; };
/**
* get the step number of current frame
* \return step number
*/
int getStep() { return _step; };
/**
* \brief pbc correct distance of two beads
* \param bead1 index of first bead
* \param bead2 index of second bead
* \return distance vector
*
* calculates the smallest distance between two beads with correct treatment
* of pbc
*/
vec getDist(int bead1, int bead2) const;
/**
* \brief calculate shortest vector connecting two points
* \param r1 first point
* \param r2 second point
* \return distance vector
*
* calculates the smallest distance between two points with correct treatment
* of pbc
*/
vec BCShortestConnection(const vec &r1, const vec &r2) const;
/**
* \brief return the shortest box size
* \return shortest size
*
* Calculates the shortest length to connect two sides of the box
*/
double ShortestBoxSize();
/**
* calculates the box volume
* \return box volume
*/
double BoxVolume();
/**
* rebuild exclusion list
*/
void RebuildExclusions();
/**
* access exclusion list
* \return exclusion list
*/
ExclusionList &getExclusions() { return _exclusions; }
BoundaryCondition::eBoxtype getBoxType() {
return _bc->getBoxType();
}
void InsertExclusion(int i, list<int> l);
protected:
BoundaryCondition *_bc;
BoundaryCondition::eBoxtype autoDetectBoxType(const matrix &box);
/// bead types in the topology
BeadTypeContainer _beadtypes;
/// beads in the topology
BeadContainer _beads;
/// molecules in the topology
MoleculeContainer _molecules;
/// residues in the topology
ResidueContainer _residues;
/// bonded interactions in the topology
InteractionContainer _interactions;
ExclusionList _exclusions;
map<string, int> _interaction_groups;
map<string, int> _beadtype_map;
map<string, list<Interaction *> > _interactions_by_group;
double _time;
int _step;
};
inline Bead *Topology::CreateBead(byte_t symmetry, string name, BeadType *type, int resnr, double m, double q)
{
Bead *b = new Bead(this, _beads.size(), type, symmetry, name, resnr, m, q);
_beads.push_back(b);
return b;
}
inline Molecule *Topology::CreateMolecule(string name)
{
Molecule *mol = new Molecule(this, _molecules.size(), name);
_molecules.push_back(mol);
return mol;
}
inline Residue *Topology::CreateResidue(string name)
{
Residue *res = new Residue(this, _molecules.size(), name);
_residues.push_back(res);
return res;
}
inline Molecule *Topology::MoleculeByIndex(int index)
{
return _molecules[index];
}
inline void Topology::InsertExclusion(int i, list<int> l) {
_exclusions.InsertExclusion(i,l);
}
}}
#include "interaction.h"
#endif /* _topology_H */
|