votca 2024.2-dev
Loading...
Searching...
No Matches
topology.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 The VOTCA Development Team (http://www.votca.org)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 *
16 */
17
18#ifndef VOTCA_CSG_TOPOLOGY_H
19#define VOTCA_CSG_TOPOLOGY_H
20#pragma once
21
22// Standard includes
23#include <cassert>
24#include <map>
25#include <memory>
26#include <unordered_map>
27#include <vector>
28
29// Third party includes
30#include <boost/container/deque.hpp>
31
32// VOTCA includes
33#include <votca/tools/types.h>
34
35// Local VOTCA includes
36#include "bead.h"
37#include "boundarycondition.h"
38#include "exclusionlist.h"
39#include "molecule.h"
40#include "openbox.h"
41#include "orthorhombicbox.h"
42#include "residue.h"
43#include "triclinicbox.h"
44
45namespace votca {
46namespace csg {
47
48class Interaction;
49
50/* Boost deque has been chosen and contents have been replaced with objects
51 * as opposed to heap allocated types:
52 * 1. To get rid of indirection
53 * 2. Clarify ownership
54 * 3. To ensure pointers are not invalidated if the container changes size
55 * that is not a guarantee with a vector
56 * 4. To provide better contiguous memory access, not possible with std::deque
57 * or list
58 */
59typedef boost::container::deque_options<
60 boost::container::block_size<sizeof(Residue) * 4>>::type block_residue_x4_t;
61typedef boost::container::deque_options<
62 boost::container::block_size<sizeof(Molecule) * 4>>::type
64typedef boost::container::deque_options<
65 boost::container::block_size<sizeof(Bead) * 4>>::type block_bead_x4_t;
66
68 boost::container::deque<Molecule, void, block_molecule_4x_t>;
69using BeadContainer = boost::container::deque<Bead, void, block_bead_x4_t>;
71 boost::container::deque<Residue, void, block_residue_x4_t>;
72using InteractionContainer = std::vector<Interaction *>;
73
81class Topology {
82 public:
84 Topology() { bc_ = std::make_unique<OpenBox>(); }
85 ~Topology();
86
90 void Cleanup();
91
105 Bead *CreateBead(Bead::Symmetry symmetry, std::string name, std::string type,
106 Index resnr, double m, double q);
107
113 Molecule *CreateMolecule(std::string name);
114
119 void CheckMoleculeNaming(void);
120
127 Residue &CreateResidue(std::string name);
128 Residue &CreateResidue(std::string name, Index id);
129
137 void CreateMoleculesByRange(std::string name, Index first, Index nbeads,
138 Index nmolecules);
139
144 Index MoleculeCount() const { return molecules_.size(); }
145
150 Index BeadCount() const { return beads_.size(); }
151
156 Index ResidueCount() const { return residues_.size(); }
157
164
170
176 const ResidueContainer &Residues() const { return residues_; }
177
183 const MoleculeContainer &Molecules() const { return molecules_; }
184
191 return interactions_;
192 }
193
195 std::vector<Interaction *> InteractionsInGroup(const std::string &group);
196
202 bool BeadTypeExist(std::string type) const;
203
210 void RegisterBeadType(std::string type);
211
219 Index getBeadTypeId(std::string type) const;
220
227 Bead *getBead(const Index i) { return &beads_[i]; }
228 const Bead *getBead(const Index i) const { return &beads_[i]; }
229 Residue &getResidue(const Index i) { return residues_[i]; }
230 const Residue &getResidue(const Index i) const { return residues_[i]; }
231 Molecule *getMolecule(const Index i) { return &molecules_[i]; }
232 const Molecule *getMolecule(const Index i) const { return &molecules_[i]; }
233
237 void ClearMoleculeList() { molecules_.clear(); }
238
243 void CopyTopologyData(Topology *top);
244
251 void RenameMolecules(std::string range, std::string name);
252
258 void RenameBeadType(std::string name, std::string newname);
259
265 void SetBeadTypeMass(std::string name, double value);
266
272 void setBox(const Eigen::Matrix3d &box, BoundaryCondition::eBoxtype boxtype =
274 // determine box type automatically in case boxtype==typeAuto
275 if (boxtype == BoundaryCondition::typeAuto) {
276 boxtype = autoDetectBoxType(box);
277 }
278
279 switch (boxtype) {
281 bc_ = std::make_unique<TriclinicBox>();
282 break;
284 bc_ = std::make_unique<OrthorhombicBox>();
285 break;
286 default:
287 bc_ = std::make_unique<OpenBox>();
288 break;
289 }
290
291 bc_->setBox(box);
292 };
293
298 const Eigen::Matrix3d &getBox() const { return bc_->getBox(); };
299
304 assert(bc_ && "Cannot return boundary condition is null");
305 return *bc_;
306 };
311 void setTime(double t) { time_ = t; };
312
317 double getTime() const { return time_; };
318
323 void setStep(Index s) { step_ = s; };
324
329 Index getStep() const { return step_; };
330
335 void setParticleGroup(std::string particle_group) {
336 particle_group_ = particle_group;
337 };
338
343 std::string getParticleGroup() const { return particle_group_; };
344
354 Eigen::Vector3d getDist(Index bead1, Index bead2) const;
355
365 Eigen::Vector3d BCShortestConnection(const Eigen::Vector3d &r_i,
366 const Eigen::Vector3d &r_j) const;
367
374 double ShortestBoxSize() const;
375
380 double BoxVolume() const;
381
385 void RebuildExclusions();
386
392 const ExclusionList &getExclusions() const { return exclusions_; }
393
394 BoundaryCondition::eBoxtype getBoxType() const { return bc_->getBoxType(); }
395
396 template <typename iteratable>
397 void InsertExclusion(Bead *bead1, iteratable &l);
398
399 bool HasVel() { return has_vel_; }
400 void SetHasVel(const bool v) { has_vel_ = v; }
401
402 bool HasForce() { return has_force_; }
403 void SetHasForce(const bool v) { has_force_ = v; }
404
405 protected:
406 std::unique_ptr<BoundaryCondition> bc_;
407
409 const Eigen::Matrix3d &box) const;
410
412 std::unordered_map<std::string, Index> beadtypes_;
413
416
419
422
425
427
428 std::map<std::string, Index> interaction_groups_;
429
430 std::map<std::string, std::vector<Interaction *>> interactions_by_group_;
431
432 double time_ = 0.0;
434 bool has_vel_ = false;
435 bool has_force_ = false;
436
438 std::string particle_group_ = "unassigned";
439};
440
441inline Bead *Topology::CreateBead(Bead::Symmetry symmetry, std::string name,
442 std::string type, Index resnr, double m,
443 double q) {
444
445 beads_.push_back(Bead(beads_.size(), type, symmetry, name, resnr, m, q));
446 return &beads_.back();
447}
448
449inline Molecule *Topology::CreateMolecule(std::string name) {
450 molecules_.push_back(Molecule(molecules_.size(), name));
451 return &molecules_.back();
452}
453
454inline Residue &Topology::CreateResidue(std::string name, Index id) {
455 // Note that Residue constructor is intentionally private and only topology
456 // class can create it, hence emplace back will not work because the vector
457 // class does not have access to the constructor.
458 residues_.push_back(Residue(id, name));
459 return residues_.back();
460}
461
462inline Residue &Topology::CreateResidue(std::string name) {
463 // Note that Residue constructor is intentionally private and only topology
464 // class can create it, hence emplace back will not work because the vector
465 // class does not have access to the constructor.
466 residues_.push_back(Residue(residues_.size(), name));
467 return residues_.back();
468}
469
471 return &molecules_[index];
472}
473
474template <typename iteratable>
475inline void Topology::InsertExclusion(Bead *bead1, iteratable &l) {
477}
478
479} // namespace csg
480} // namespace votca
481
482#include "interaction.h"
483
484#endif // VOTCA_CSG_TOPOLOGY_H
information about a bead
Definition bead.h:50
Symmetry
get the symmetry of the bead
Definition bead.h:85
Class keeps track of how the boundaries of the system are handled.
void InsertExclusion(Bead *bead, iterable &excluded)
base class for all interactions
Definition interaction.h:40
information about molecules
Definition molecule.h:45
class for a residue
Definition residue.h:35
topology of the whole system
Definition topology.h:81
const BoundaryCondition & getBoundary() const
Return the boundary condition object.
Definition topology.h:303
double getTime() const
Definition topology.h:317
BoundaryCondition::eBoxtype getBoxType() const
Definition topology.h:394
Residue & CreateResidue(std::string name)
Create a new resiude.
Definition topology.h:462
BoundaryCondition::eBoxtype autoDetectBoxType(const Eigen::Matrix3d &box) const
Definition topology.cc:252
void setBox(const Eigen::Matrix3d &box, BoundaryCondition::eBoxtype boxtype=BoundaryCondition::typeAuto)
Definition topology.h:272
void Cleanup()
Cleans up all the stored data.
Definition topology.cc:49
std::unordered_map< std::string, Index > beadtypes_
bead types in the topology
Definition topology.h:412
void CreateMoleculesByRange(std::string name, Index first, Index nbeads, Index nmolecules)
create molecules based on blocks of atoms
Definition topology.cc:72
Index MoleculeCount() const
number of molecules in the system
Definition topology.h:144
void RenameMolecules(std::string range, std::string name)
rename all the molecules in range
Definition topology.cc:138
Index getBeadTypeId(std::string type) const
Given a bead type this method returns the id associated with the type.
Definition topology.cc:133
BeadContainer & Beads()
Definition topology.h:169
void setTime(double t)
Definition topology.h:311
std::vector< Interaction * > InteractionsInGroup(const std::string &group)
Definition topology.cc:201
const Bead * getBead(const Index i) const
Definition topology.h:228
void RenameBeadType(std::string name, std::string newname)
rename all the bead types
Definition topology.cc:150
Index ResidueCount() const
Definition topology.h:156
bool BeadTypeExist(std::string type) const
Determine if a bead type exists.
Definition topology.cc:210
MoleculeContainer molecules_
molecules in the topology
Definition topology.h:418
std::map< std::string, Index > interaction_groups_
Definition topology.h:428
double ShortestBoxSize() const
return the shortest box size
Definition topology.cc:268
const Eigen::Matrix3d & getBox() const
Definition topology.h:298
Bead * CreateBead(Bead::Symmetry symmetry, std::string name, std::string type, Index resnr, double m, double q)
Creates a new Bead.
Definition topology.h:441
void SetHasVel(const bool v)
Definition topology.h:400
std::string particle_group_
The particle group (For H5MD file format)
Definition topology.h:438
const MoleculeContainer & Molecules() const
Definition topology.h:183
const InteractionContainer & BondedInteractions() const
Definition topology.h:190
ExclusionList exclusions_
Definition topology.h:426
InteractionContainer interactions_
bonded interactions in the topology
Definition topology.h:424
void SetHasForce(const bool v)
Definition topology.h:403
void AddBondedInteraction(Interaction *ic)
Definition topology.cc:188
Index BeadCount() const
Definition topology.h:150
Index getStep() const
Definition topology.h:329
Molecule * getMolecule(const Index i)
Definition topology.h:231
ExclusionList & getExclusions()
Definition topology.h:391
void RegisterBeadType(std::string type)
Register the bead type with the topology object.
Definition topology.cc:214
void SetBeadTypeMass(std::string name, double value)
set the mass of all the beads of a certain type
Definition topology.cc:160
Molecule * CreateMolecule(std::string name)
Creates a new molecule.
Definition topology.h:449
Molecule * MoleculeByIndex(Index index)
Definition topology.h:470
const ExclusionList & getExclusions() const
Definition topology.h:392
const ResidueContainer & Residues() const
Definition topology.h:176
std::string getParticleGroup() const
Definition topology.h:343
std::unique_ptr< BoundaryCondition > bc_
Definition topology.h:406
void CheckMoleculeNaming(void)
checks weather molecules with the same name really contain the same number of beads
Definition topology.cc:169
Eigen::Vector3d getDist(Index bead1, Index bead2) const
pbc correct distance of two beads
Definition topology.cc:243
Residue & getResidue(const Index i)
Definition topology.h:229
Eigen::Vector3d BCShortestConnection(const Eigen::Vector3d &r_i, const Eigen::Vector3d &r_j) const
calculate shortest vector connecting two points
Definition topology.cc:238
const Residue & getResidue(const Index i) const
Definition topology.h:230
Bead * getBead(const Index i)
Returns a pointer to the bead with index i.
Definition topology.h:227
std::map< std::string, std::vector< Interaction * > > interactions_by_group_
Definition topology.h:430
const Molecule * getMolecule(const Index i) const
Definition topology.h:232
BeadContainer beads_
beads in the topology
Definition topology.h:415
void CopyTopologyData(Topology *top)
copy topology data of different topology
Definition topology.cc:102
void setStep(Index s)
Definition topology.h:323
void InsertExclusion(Bead *bead1, iteratable &l)
Definition topology.h:475
double BoxVolume() const
Definition topology.cc:248
MoleculeContainer & Molecules()
Definition topology.h:182
Topology()
constructor
Definition topology.h:84
ResidueContainer & Residues()
Definition topology.h:175
InteractionContainer & BondedInteractions()
Definition topology.h:189
ResidueContainer residues_
residues in the topology
Definition topology.h:421
void setParticleGroup(std::string particle_group)
Definition topology.h:335
boost::container::deque_options< boost::container::block_size< sizeof(Residue) *4 > >::type block_residue_x4_t
Definition topology.h:60
boost::container::deque< Molecule, void, block_molecule_4x_t > MoleculeContainer
Definition topology.h:67
boost::container::deque< Bead, void, block_bead_x4_t > BeadContainer
Definition topology.h:69
std::vector< Interaction * > InteractionContainer
Definition topology.h:72
boost::container::deque_options< boost::container::block_size< sizeof(Bead) *4 > >::type block_bead_x4_t
Definition topology.h:65
boost::container::deque_options< boost::container::block_size< sizeof(Molecule) *4 > >::type block_molecule_4x_t
Definition topology.h:63
boost::container::deque< Residue, void, block_residue_x4_t > ResidueContainer
Definition topology.h:70
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26