Первая говноверсия говнопроливновых говноспиралей.
[alexxy/gromacs-dssp.git] / src / dssp.cpp
1 /*
2 * This file is part of the GROMACS molecular simulation package.
3 *
4 * Copyright (c) 2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
8 *
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
13 *
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17 * Lesser General Public License for more details.
18 *
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23 *
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
31 *
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
34 */
35 /*! \internal \file
36 * \brief
37 * Implements gmx::analysismodules::Trajectory.
38 *
39 * \author Sergey Gorelov <gorelov_sv@pnpi.nrcki.ru>
40 * \author Anatoly Titov <titov_ai@pnpi.nrcki.ru>
41 * \author Alexey Shvetsov <alexxyum@gmail.com>
42 * \ingroup module_trajectoryanalysis
43 */
44
45 #include "dssp.h"
46 #include "dssptools.h"
47
48 #include <algorithm>
49
50 #include <iostream> //remove
51
52 #include <gromacs/trajectoryanalysis.h>
53
54 namespace gmx
55 {
56
57 namespace analysismodules
58 {
59
60 class Dssp : public TrajectoryAnalysisModule
61 {
62 public:
63    Dssp();
64
65    void initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings) override;
66
67    void optionsFinished(TrajectoryAnalysisSettings* settings) override;
68
69
70    void initAnalysis(const TrajectoryAnalysisSettings& settings, const TopologyInformation& top) override;
71
72    void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* /* pdata */) override;
73
74    void finishAnalysis(int nframes) override;
75    void writeOutput() override;
76
77 private:
78    initParameters                                       initParams;
79    std::string                                          filename_;
80    int                                                  nres;
81    DsspTool                                             DT;
82    std::vector<std::pair<int, std::string>>             data;
83 };
84
85 Dssp::Dssp() {}
86
87
88 void Dssp::initOptions(IOptionsContainer* options, TrajectoryAnalysisSettings* settings)
89 {
90    static const char* const desc[] = {
91        "[THISMODULE] todo",
92    };
93    options->addOption(
94             StringOption("o").store(&filename_).required().defaultValue("SSP.dat").description("Filename for DSSP output"));
95    options->addOption(RealOption("cutoff").store(&initParams.cutoff_).required().defaultValue(4.0).description(
96            "cutoff for neighbour search"));
97    options->addOption(
98              SelectionOption("sel").store(&initParams.sel_).required().description(
99                    "Group for DSSP"));
100    options->addOption(
101                EnumOption<NBSearchMethod>("algo").store(&initParams.NBS).required().defaultValue(NBSearchMethod::DSSP).enumValue(NBSearchMethodNames).description("Algorithm for search residues' neighbours"));
102    options->addOption(
103                BooleanOption("addh").store(&initParams.addHydrogens).defaultValue(false).description("Add hydrogens or not"));
104    options->addOption(
105                BooleanOption("p").store(&initParams.PPHelices).defaultValue(true).description("Prefer Pi Helices"));
106    options->addOption(
107                IntegerOption("s").store(&initParams.pp_stretch).defaultValue(3).description("Stretch for your pp. 2 or 3"));
108    options->addOption(
109                BooleanOption("v").store(&initParams.verbose).defaultValue(false).description(">:("));
110    settings->setHelpText(desc);
111 }
112
113
114 void Dssp::optionsFinished(TrajectoryAnalysisSettings* /* settings */) {
115 }
116
117
118 void Dssp::initAnalysis(const TrajectoryAnalysisSettings &settings, const TopologyInformation& top)
119 {
120     DT.initAnalysis(top, initParams);
121     nres = DT.nres;
122 }
123
124 void Dssp::analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* /* pdata*/)
125 {
126     DT.analyzeFrame(frnr, fr, pbc);
127 }
128
129 void Dssp::finishAnalysis(int /*nframes*/) {
130     data.resize(0);
131     data = DT.getData();
132 }
133
134
135 void Dssp::writeOutput() {
136
137     std::ofstream dsspout;
138     dsspout.open(filename_, std::ofstream::out | std::ofstream::trunc);
139     dsspout << nres << "\n";
140     for (std::size_t i{ 0 }; i < data.size(); ++i)
141     {
142         for(std::size_t j{0}; j < data[i].second.size(); ++j){
143             dsspout << data[i].second[j];
144         }
145         dsspout << "\n";
146     }
147     dsspout.close();
148     data.resize(0);
149 }
150
151 const char DsspInfo::name[]             = "dssp";
152 const char DsspInfo::shortDescription[] = "Calculate protein secondary structure via DSSP algo";
153
154 TrajectoryAnalysisModulePointer DsspInfo::create()
155 {
156    return TrajectoryAnalysisModulePointer(new Dssp);
157 }
158
159 } // namespace analysismodules
160
161 } // namespace gmx
162
163 int main(int argc, char *argv[]) {
164     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<gmx::analysismodules::Dssp>(argc, argv);
165 }