votca 2024-dev
Loading...
Searching...
No Matches
topology.cc
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// Standard includes
19#include <cassert>
20#include <memory>
21#include <regex>
22#include <stdexcept>
23#include <unordered_set>
24
25// Third party includes
26#include <boost/lexical_cast.hpp>
27
28// VOTCA includes
30
31// Local VOTCA includes
34#include "votca/csg/molecule.h"
35#include "votca/csg/openbox.h"
36#include "votca/csg/topology.h"
37
38namespace votca {
39namespace csg {
40
41using namespace std;
42
43bool is_digits(const std::string &str) {
44 return str.find_first_not_of("0123456789") == std::string::npos;
45}
46
48
50 // cleanup beads
51 beads_.clear();
52
53 // cleanup molecules
54 molecules_.clear();
55
56 // cleanup residues
57 residues_.clear();
58
59 // cleanup interactions
60 {
61 for (InteractionContainer::iterator i = interactions_.begin();
62 i < interactions_.end(); ++i) {
63 delete (*i);
64 }
65 interactions_.clear();
66 }
67 // cleanup bc_ object
68 bc_ = std::make_unique<OpenBox>();
69}
70
72void Topology::CreateMoleculesByRange(string name, Index first, Index nbeads,
73 Index nmolecules) {
74 Molecule *mol = CreateMolecule(name);
75 Index beadcount = 0;
76 Index res_offset = 0;
77
78 for (auto &bead_ : beads_) {
79 // xml numbering starts with 1
80 if (--first > 0) {
81 continue;
82 }
83 // This is not 100% correct, but let's assume for now that the resnr do
84 // increase
85 if (beadcount == 0) {
86 res_offset = bead_.getResnr();
87 }
88 stringstream bname;
89 bname << bead_.getResnr() - res_offset + 1 << ":"
90 << getResidue(bead_.getResnr()).getName() << ":" << bead_.getName();
91 mol->AddBead(&bead_, bname.str());
92 if (++beadcount == nbeads) {
93 if (--nmolecules <= 0) {
94 break;
95 }
96 mol = CreateMolecule(name);
97 beadcount = 0;
98 }
99 }
100}
101
103
104 bc_->setBox(top->getBox());
105 time_ = top->time_;
106 step_ = top->step_;
107
108 // cleanup old data
109 Cleanup();
110
111 // copy all residues
112 for (const auto &residue_ : top->residues_) {
113 CreateResidue(residue_.getName());
114 }
115
116 // create all beads
117 for (const auto &bi : top->beads_) {
118 string type = bi.getType();
119 CreateBead(bi.getSymmetry(), bi.getName(), type, bi.getResnr(),
120 bi.getMass(), bi.getQ());
121 }
122
123 // copy all molecules
124 for (const auto &molecule_ : top->molecules_) {
125 Molecule *mi = CreateMolecule(molecule_.getName());
126 for (Index i = 0; i < molecule_.BeadCount(); i++) {
127 Index beadid = molecule_.getBead(i)->getId();
128 mi->AddBead(&beads_[beadid], molecule_.getBeadName(i));
129 }
130 }
131}
132
133Index Topology::getBeadTypeId(string type) const {
134 assert(beadtypes_.count(type));
135 return beadtypes_.at(type);
136}
137
138void Topology::RenameMolecules(string range, string name) {
140 rp.Parse(range);
141 for (Index i : rp) {
142 if (i > Index(molecules_.size())) {
143 throw runtime_error(
144 string("RenameMolecules: num molecules smaller than"));
145 }
146 getMolecule(i - 1)->setName(name);
147 }
148}
149
150void Topology::RenameBeadType(string name, string newname) {
151
152 for (auto &bead : beads_) {
153 string type = bead.getType();
154 if (tools::wildcmp(name, type)) {
155 bead.setType(newname);
156 }
157 }
158}
159
160void Topology::SetBeadTypeMass(string name, double value) {
161 for (auto &bead : beads_) {
162 string type = bead.getType();
163 if (tools::wildcmp(name, type)) {
164 bead.setMass(value);
165 }
166 }
167}
168
170 map<string, Index> nbeads;
171
172 for (const auto &mol : molecules_) {
173 map<string, Index>::iterator entry = nbeads.find(mol.getName());
174 if (entry != nbeads.end()) {
175 if (entry->second != mol.BeadCount()) {
176 throw runtime_error(
177 "There are molecules which have the same name but different number "
178 "of bead "
179 "please check the section manual topology handling in the votca "
180 "manual");
181 }
182 continue;
183 }
184 nbeads[mol.getName()] = mol.BeadCount();
185 }
186}
187
189 map<string, Index>::iterator iter = interaction_groups_.find(ic->getGroup());
190 if (iter != interaction_groups_.end()) {
191 ic->setGroupId(iter->second);
192 } else {
193 Index i = interaction_groups_.size();
194 interaction_groups_[ic->getGroup()] = i;
195 ic->setGroupId(i);
196 }
197 interactions_.push_back(ic);
198 interactions_by_group_[ic->getGroup()].push_back(ic);
199}
200
201std::vector<Interaction *> Topology::InteractionsInGroup(const string &group) {
202 map<string, vector<Interaction *>>::iterator iter =
203 interactions_by_group_.find(group);
204 if (iter == interactions_by_group_.end()) {
205 return vector<Interaction *>();
206 }
207 return iter->second;
208}
209
210bool Topology::BeadTypeExist(string type) const {
211 return beadtypes_.count(type);
212}
213
214void Topology::RegisterBeadType(string type) {
215 unordered_set<Index> ids;
216 for (pair<const string, Index> type_and_id : beadtypes_) {
217 ids.insert(type_and_id.second);
218 }
219
220 Index id = 0;
221 // If the type is also a number use it as the id as well provided it is not
222 // already taken
223 if (is_digits(type)) {
224 id = boost::lexical_cast<Index>(type);
225 assert(!ids.count(id) &&
226 "The type passed in is a number and has already"
227 " been registered. It is likely that you are passing in numbers as "
228 "bead types as well as strings, choose one or the other do not mix "
229 "between using numbers and strings ");
230 }
231
232 while (ids.count(id)) {
233 ++id;
234 }
235 beadtypes_[type] = id;
236}
237
239 const Eigen::Vector3d &r_i, const Eigen::Vector3d &r_j) const {
240 return bc_->BCShortestConnection(r_i, r_j);
241}
242
243Eigen::Vector3d Topology::getDist(Index bead1, Index bead2) const {
244 return BCShortestConnection(getBead(bead1)->getPos(),
245 getBead(bead2)->getPos());
246}
247
248double Topology::BoxVolume() const { return bc_->BoxVolume(); }
249
251
253 const Eigen::Matrix3d &box) const {
254 // set the box type to OpenBox in case "box" is the zero matrix,
255 // to OrthorhombicBox in case "box" is a diagonal matrix,
256 // or to TriclinicBox otherwise
257 if (box.isApproxToConstant(0)) {
259 } else if ((box - Eigen::Matrix3d(box.diagonal().asDiagonal()))
260 .isApproxToConstant(0)) {
262 } else {
264 }
266}
267
269 return bc_->getShortestBoxDimension();
270}
271
272} // namespace csg
273} // namespace votca
void CreateExclusions(Topology *top)
base class for all interactions
Definition interaction.h:40
const std::string & getGroup() const
Definition interaction.h:53
void setGroupId(Index id)
Definition interaction.h:65
information about molecules
Definition molecule.h:45
void AddBead(Bead *bead, const std::string &name)
Add a bead to the molecule.
Definition molecule.cc:29
void setName(const std::string &name)
set the name of the molecule
Definition molecule.h:54
const std::string & getName() const
get the name of the residue
Definition residue.h:52
topology of the whole system
Definition topology.h:81
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 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
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
std::vector< Interaction * > InteractionsInGroup(const std::string &group)
Definition topology.cc:201
void RenameBeadType(std::string name, std::string newname)
rename all the bead types
Definition topology.cc:150
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
ExclusionList exclusions_
Definition topology.h:426
InteractionContainer interactions_
bonded interactions in the topology
Definition topology.h:424
void AddBondedInteraction(Interaction *ic)
Definition topology.cc:188
Molecule * getMolecule(const Index i)
Definition topology.h:231
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
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
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
BeadContainer beads_
beads in the topology
Definition topology.h:415
void CopyTopologyData(Topology *top)
copy topology data of different topology
Definition topology.cc:102
double BoxVolume() const
Definition topology.cc:248
ResidueContainer residues_
residues in the topology
Definition topology.h:421
void Parse(std::string str)
STL namespace.
bool is_digits(const std::string &str)
Definition topology.cc:43
int wildcmp(const char *wild, const char *string)
Definition tokenizer.cc:28
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26