votca 2025.1-dev
Loading...
Searching...
No Matches
csg_dlptopol.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2023 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 <fstream>
20#include <iostream>
21
22// Third party includes
23#include <boost/format.hpp>
24
25// Local VOTCA includes
27
28using namespace votca::csg;
29using namespace std;
30
37
39 public:
40 string ProgramName() override { return "csg_dlptopol"; }
41 void HelpText(ostream &out) override {
42 out << "Create a dlpoly topology template based on an existing (atomistic) "
43 "topology and \n"
44 << "a mapping xml-file. The created template file needs to be "
45 "inspected and amended by the user!\n\n"
46 << "Examples:\n"
47 << "* csg_dlptopol --top .dlpf --out .dlpf --cg cg-map.xml\n convert "
48 "FIELD to FIELD_CGV using cg-map.xml\n"
49 << "* csg_dlptopol --top FA-dlpoly.dlpf --out CG-dlpoly.dlpf --cg "
50 "cg-map.xml\n"
51 << "* csg_dlptopol --top FA-gromacs.tpr --out FA-dlpoly.dlpf "
52 "--no-map\n";
53 }
54
55 bool DoMapping(void) override { return true; }
56
57 void Initialize(void) override;
58 bool EvaluateOptions(void) override {
60 CheckRequired("out", "no output topology specified");
61 return true;
62 }
63 bool EvaluateTopology(Topology *top, Topology *top_ref) override;
64
65 protected:
66 void WriteMoleculeAtoms(ostream &out, const Molecule &cg);
67 void WriteMoleculeInteractions(ostream &out, const Molecule &cg);
68 void WriteVDWInteractions(ostream &out, const Molecule &cg);
69 void WriteMolecularType(ostream &out, const Molecule &cg,
70 votca::Index nummol);
71};
72
75 // AddProgramOptions()
76 //("top", boost::program_options::value<string>(),
77 //" input topology in any known format:\n <name>.dlpf for dlpoly, <name>.tpr
78 // for gromacs\n (convention: '.dlpf'='use FIELD')");
79 AddProgramOptions()("out", boost::program_options::value<string>(),
80 " output topology in dlpoly format");
81}
82
84 // check the file names from options
85
86 string fname = OptionsMap()["top"].as<string>();
87
88 if (fname == ".dlpf") {
89 fname = "FIELD";
90 }
91
92 cout << "input file-name: " << fname << endl;
93
94 fname = OptionsMap()["out"].as<string>();
95
96#ifndef NDEBUG
97 cout << "output file-name given: " << fname << endl;
98#endif
99
100 if (fname == ".dlpf") {
101 fname = "FIELD_CGV";
102 }
103
104#ifndef NDEBUG
105 cout << "output file-name actual: " << fname << endl;
106#else
107 cout << "output file-name: " << fname << endl;
108#endif
109
110 // do CG mapping
111
112 MoleculeContainer &mols = top->Molecules();
113 std::vector<const Molecule *> MolecularTypes;
114
115 votca::Index prv_mol_number = 1;
116 string prv_mol_name;
117 vector<votca::Index> nummols;
118
119 vector<string> vdw_pairs;
120
121 // find all unique molecular types
122
123 for (const auto &mol : mols) {
124 // molecules are ignored during the mapping stage
125 // i.e. the ignored ones do not enter the CG topology (*top) - ?
126 // if( IsIgnored(mol->getName()) ) continue;
127
128 if (mol.getName() == prv_mol_name) {
129 prv_mol_number++;
130 continue;
131 }
132
133 nummols.push_back(prv_mol_number);
134 prv_mol_number = 1;
135 prv_mol_name = mol.getName();
136
137 cout << "'" << mol.getName() << "' added to CG molecular types" << endl;
138
139 MolecularTypes.push_back(&mol);
140
141 // collect unique bead pairs over all molecular types found
142
143 for (votca::Index ib1 = 0; ib1 < mol.BeadCount(); ib1++) {
144 string bead_name1 = mol.getBead(ib1)->getType();
145 bead_name1 = bead_name1.substr(
146 0,
147 bead_name1.find_first_of("#")); // skip #index of atom from its name
148
149 for (const auto &MolecularType : MolecularTypes) {
150
151 for (votca::Index ib2 = 0; ib2 < MolecularType->BeadCount(); ib2++) {
152
153 string bead_name2 = MolecularType->getBead(ib2)->getType();
154 bead_name2 = bead_name2.substr(
155 0, bead_name2.find_first_of("#")); // skip #index of atom from
156 // its name
157
158 stringstream ss_bp1, ss_bp2;
159
160 ss_bp1 << boost::format("%8s%8s") % bead_name1 % bead_name2;
161 ss_bp2 << boost::format("%8s%8s") % bead_name2 % bead_name1;
162
163 bool is_new_pair = true;
164
165 for (const auto &vdw_pair : vdw_pairs) {
166 if (ss_bp1.str() == vdw_pair || ss_bp2.str() == vdw_pair) {
167 is_new_pair = false;
168 break;
169 }
170 }
171 if (is_new_pair) {
172 vdw_pairs.push_back(ss_bp1.str());
173#ifndef NDEBUG
174 cout << "'" << ss_bp1.str() << "' added to CG vdw pairs" << endl;
175#endif
176 }
177 }
178 }
179 }
180 }
181 nummols.push_back(prv_mol_number);
182
183 if (MolecularTypes.size() > 1) {
184 cout << "WARNING: creation of topology for multiple molecular types "
185 "is experimental at this stage\n";
186 }
187
188 ofstream fl;
189 fl.open(fname);
190
191 fl << "From VOTCA with love"
192 << " # check/amend this file if needed!\n";
193 fl << "units kJ\n";
194 fl << "molecular types " << MolecularTypes.size() << endl;
195
196 for (votca::Index i = 0; i < votca::Index(MolecularTypes.size()); i++) {
197 WriteMolecularType(fl, *(MolecularTypes[i]), nummols[i + 1]);
198 }
199
200 // vdw seciton (pairwise vdw/CG interactions)
201
202 if (vdw_pairs.size() > 0) {
203
204 fl << "vdw " << vdw_pairs.size() << endl;
205
206 for (const auto &vdw_pair : vdw_pairs) {
207 fl << vdw_pair << " tab 1.00000 0.00000\n";
208 }
209 }
210
211 fl << "close" << endl;
212
213 cout << "Created template for dlpoly topology - please, check & amend if "
214 "needed!"
215 << endl;
216
217 fl.close();
218 return true;
219}
220
221void DLPTopolApp::WriteMoleculeAtoms(ostream &out, const Molecule &cg) {
222 out << "atoms " << cg.BeadCount() << endl;
223 out << "# name mass charge nrept ifrozen (optional: ngroup, index, "
224 "name/type, type/residue, index/res-ID) \n";
225 for (votca::Index i = 0; i < cg.BeadCount(); ++i) {
226 const Bead *b = cg.getBead(i);
227
228 string bname = b->getName();
229 string btype = b->getType();
230
231 bname = bname.substr(
232 0, bname.find_first_of("#")); // skip #index of atom from its name
233 btype = btype.substr(
234 0, btype.find_first_of("#")); // skip #index of atom from its type
235
236 out << boost::format("%8s %10f %10f 1 0 1 %10d %8s %8s %10d \n") %
237 btype % b->getMass() % b->getQ() % (i + 1) % btype % bname %
238 (i + 1);
239 //% b->getType()->getName() % b->getMass() % b->getQ() % (i+1) %
240 // b->getType()->getName() % b->getName() % (i+1);
241 }
242}
243
245
246 stringstream sout;
247
248 votca::Index n_entries = 0;
249 votca::Index nb = -1;
250
251 for (const Interaction *ic : cg.Interactions()) {
252 if (nb != ic->BeadCount()) {
253
254 if (sout.str() != "") {
255 out << n_entries << endl << sout.str();
256 }
257
258 sout.str("");
259 n_entries = 0;
260
261 nb = ic->BeadCount();
262 switch (nb) {
263 case 2:
264 out << "bonds ";
265 break;
266 case 3:
267 out << "angles ";
268 break;
269 case 4:
270 out << "dihedrals ";
271 break;
272 default:
273 string err = "cannot handle number of beads in interaction:";
274 err += to_string(ic->getMolecule() + 1) + ":" + ic->getGroup();
275 err += ":" + to_string(ic->getIndex() + 1);
276 throw runtime_error(err);
277 }
278 }
279 n_entries++;
280 // to do: is it possible to use bond/angle/dihedral function types for 1:1
281 // mapping? (CG map overwrites ic->Group anyway)
282 // sout << ic->getInteractionFunc(); // something like that (only for 1:1
283 // mapping!)
284 sout << " tab ";
285 for (votca::Index i = 0; i < nb; ++i) {
286 sout << ic->getBeadId(i) + 1 << " ";
287 }
288 sout << " 1.00000 0.00000"
289 << " # ";
290 sout << to_string(ic->getMolecule() + 1);
291 sout << ":" + ic->getGroup();
292 sout << ":" + to_string(ic->getIndex() + 1) << endl;
293 }
294 if (sout.str() != "") {
295 out << n_entries << endl << sout.str();
296 }
297}
298
299void DLPTopolApp::WriteMolecularType(ostream &out, const Molecule &cg,
300 votca::Index nummol) {
301 out << cg.getName() << endl;
302 out << "nummols " << nummol << endl;
303
304 WriteMoleculeAtoms(out, cg);
306
307 out << "finish" << endl;
308}
309
310int main(int argc, char **argv) {
311 DLPTopolApp app;
312 return app.Exec(argc, argv);
313}
class for writing dlpoly topology files
void HelpText(ostream &out) override
help text of application without version information
void WriteMolecularType(ostream &out, const Molecule &cg, votca::Index nummol)
bool DoMapping(void) override
overload and return true to enable mapping command line options
string ProgramName() override
program name
bool EvaluateOptions(void) override
Process command line options.
void WriteMoleculeAtoms(ostream &out, const Molecule &cg)
bool EvaluateTopology(Topology *top, Topology *top_ref) override
called after topology was loaded
void WriteMoleculeInteractions(ostream &out, const Molecule &cg)
void WriteVDWInteractions(ostream &out, const Molecule &cg)
void Initialize(void) override
Initialize application data.
virtual const std::string getType() const noexcept
Definition basebead.h:84
std::string getName() const
Gets the name of the bead.
Definition basebead.h:58
virtual const double & getMass() const noexcept
Definition basebead.h:104
information about a bead
Definition bead.h:50
virtual const double & getQ() const
Definition bead.h:67
bool EvaluateOptions() override
Process command line options.
void Initialize() override
Initialize application data.
base class for all interactions
Definition interaction.h:40
information about molecules
Definition molecule.h:45
Bead * getBead(Index bead)
get the id of a bead in the molecule
Definition molecule.h:59
std::vector< Interaction * > Interactions()
Definition molecule.h:79
const std::string & getName() const
get the name of the molecule
Definition molecule.h:51
Index BeadCount() const
get the number of beads in the molecule
Definition molecule.h:65
topology of the whole system
Definition topology.h:81
MoleculeContainer & Molecules()
Definition topology.h:182
int Exec(int argc, char **argv)
executes the program
boost::program_options::variables_map & OptionsMap()
get available program options & descriptions
boost::program_options::options_description_easy_init AddProgramOptions(const std::string &group="")
add option for command line
void CheckRequired(const std::string &option_name, const std::string &error_msg="")
Check weather required option is set.
int main(int argc, char **argv)
STL namespace.
boost::container::deque< Molecule, void, block_molecule_4x_t > MoleculeContainer
Definition topology.h:67
Eigen::Index Index
Definition types.h:26