votca 2024-dev
Loading...
Searching...
No Matches
csg_gmxtopol.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 <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;
30using boost::format;
31
33 public:
34 string ProgramName() override { return "csg_gmxtopol"; }
35 void HelpText(ostream &out) override {
36 out << "Create skeleton for gromacs topology based on atomistic topology\n"
37 "and a mapping file. File still needs to be modified by the user.";
38 }
39
40 bool DoMapping(void) override { return true; }
41
42 void Initialize(void) override;
43 bool EvaluateOptions(void) override {
45 CheckRequired("out", "no output topology specified");
46 return true;
47 }
48 bool EvaluateTopology(Topology *top, Topology *top_ref) override;
49
50 protected:
51 void WriteAtoms(ostream &out, Molecule &cg);
52 void WriteInteractions(ostream &out, const Topology &top, Molecule &cg);
53 void WriteMolecule(ostream &out, const Topology &top, Molecule &cg);
54};
55
59 "out", boost::program_options::value<string>(),
60 "output topology (will create .top and in future also .itp)");
61}
62
64 if (top->MoleculeCount() > 1) {
65 cout << "WARNING: cannot create topology for topology with"
66 "multiple molecules, using only first molecule\n";
67 }
68 ofstream fl;
69 fl.open((OptionsMap()["out"].as<string>() + ".top"));
70 WriteMolecule(fl, *top, *(top->MoleculeByIndex(0)));
71 fl.close();
72 return true;
73}
74
75void GmxTopolApp::WriteAtoms(ostream &out, Molecule &cg) {
76 out << "[atoms]\n";
77 out << "; nr type resnr residue atom cgnr charge mass\n";
78 for (votca::Index i = 0; i < cg.BeadCount(); ++i) {
79 Bead *b = cg.getBead(i);
80 out << format("%d %s 1 RES %s %d %f %f\n") % (i + 1) % b->getType() %
81 b->getName() % (i + 1) % b->getQ() % b->getMass();
82 }
83 out << endl;
84}
85
86void GmxTopolApp::WriteInteractions(ostream &out, const Topology &top,
87 Molecule &cg) {
88 votca::Index nb = -1;
89
90 for (const Interaction *ic : top.BondedInteractions()) {
91 if (ic->getMolecule() != cg.getId()) {
92 continue;
93 }
94 if (nb != ic->BeadCount()) {
95 nb = ic->BeadCount();
96 switch (nb) {
97 case 2:
98 out << "\n[ bonds ]\n";
99 break;
100 case 3:
101 out << "\n[ angles ]\n";
102 break;
103 case 4:
104 out << "\n[ dihedrals ]\n";
105 break;
106 default:
107 string err = "cannot handle number of beads in interaction:";
108 err += to_string(ic->getMolecule() + 1) + ":" + ic->getGroup();
109 err += ":" + to_string(ic->getIndex() + 1);
110 throw runtime_error(err);
111 }
112 }
113 for (votca::Index i = 0; i < nb; ++i) {
114 out << ic->getBeadId(i) + 1 << " ";
115 }
116 out << " 1 ; ";
117 out << to_string(ic->getMolecule() + 1);
118 out << ":" + ic->getGroup();
119 out << ":" + to_string(ic->getIndex() + 1) << endl;
120 }
121}
122
123void GmxTopolApp::WriteMolecule(ostream &out, const Topology &top,
124 Molecule &cg) {
125 out << "[ moleculetype ]\n";
126 out << cg.getName() << " 3\n\n";
127
128 WriteAtoms(out, cg);
129 WriteInteractions(out, top, cg);
130}
131
132int main(int argc, char **argv) {
133 GmxTopolApp app;
134 return app.Exec(argc, argv);
135}
void WriteInteractions(ostream &out, const Topology &top, Molecule &cg)
bool EvaluateOptions(void) override
Process command line options.
string ProgramName() override
program name
void WriteMolecule(ostream &out, const Topology &top, Molecule &cg)
void HelpText(ostream &out) override
help text of application without version information
bool DoMapping(void) override
overload and return true to enable mapping command line options
void Initialize(void) override
Initialize application data.
void WriteAtoms(ostream &out, Molecule &cg)
bool EvaluateTopology(Topology *top, Topology *top_ref) override
called after topology was loaded
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
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
Index getId() const
get the molecule ID
Definition molecule.h:48
topology of the whole system
Definition topology.h:81
Index MoleculeCount() const
number of molecules in the system
Definition topology.h:144
Molecule * MoleculeByIndex(Index index)
Definition topology.h:470
InteractionContainer & BondedInteractions()
Definition topology.h:189
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.
Eigen::Index Index
Definition types.h:26