f3f1beb6974c4aca032f4168507ebc6d2b8693c3
[alexxy/gromacs-dssp.git] / src / dssptools.h
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  * Declares trajectory analysis module for secondary structure asignment.
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 #ifndef GMX_TRAJECTORYANALYSIS_MODULES_DSSPTOOLS_H
46 #define GMX_TRAJECTORYANALYSIS_MODULES_DSSPTOOLS_H
47
48 #include <algorithm>
49 #include "gromacs/math/units.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include <gromacs/trajectoryanalysis.h>
52 #include "gromacs/trajectoryanalysis/topologyinformation.h"
53 #include <set>
54 #include <fstream>
55 #include <mutex>
56
57 namespace gmx
58 {
59
60 namespace analysismodules
61 {
62
63 enum class NBSearchMethod : std::size_t
64 {
65     Classique = 0,
66     Experimental,
67     DSSP,
68     Count
69 };
70
71 const gmx::EnumerationArray<NBSearchMethod, const char*> NBSearchMethodNames = {
72     { "Classique", "Experimental", "DSSP" }
73 };
74
75 struct initParameters {
76     Selection                                            sel_;
77     real                                                 cutoff_; // = 4.0; ???
78     NBSearchMethod                                       NBS;
79     bool                                                 verbose, PPHelices, addHydrogens;
80 };
81
82 enum class backboneAtomTypes : std::size_t
83 {
84    AtomCA,
85    AtomC,
86    AtomO,
87    AtomN,
88    AtomH,
89    Count
90 };
91 //! String values corresponding to backbone atom types
92 const gmx::EnumerationArray<backboneAtomTypes, const char*> backboneAtomTypeNames = {
93    { "CA", "C", "O", "N", "H" }
94 };
95
96 struct ResInfo {
97     std::array<std::size_t, static_cast<std::size_t>(backboneAtomTypes::Count)>     _backboneIndices{ 0, 0, 0, 0, 0 }; // TODO something with zeroes
98     std::size_t                                                                     getIndex(backboneAtomTypes atomTypeName) const;
99     t_resinfo                                                                       *info{nullptr}, *donor[2]{nullptr, nullptr}, *acceptor[2]{nullptr, nullptr};
100     ResInfo                                                                         *prevResi{nullptr}, *nextResi{nullptr};
101     float                                                                           donorEnergy[2]{}, acceptorEnergy[2]{};
102     bool                                                                            is_proline{false};
103 };
104
105
106 //class ResInfo
107 //{
108 //public:
109 //   void   setIndex(backboneAtomTypes atomTypeName, std::size_t atomIndex);
110 //   std::size_t getIndex(backboneAtomTypes atomTypeName) const;
111
112 //   ResInfo                                                                    *donor[2]{nullptr, nullptr}, *acceptor[2]{nullptr, nullptr};
113 //   float                                                                        donorEnergy[2]{}, acceptorEnergy[2]{};
114 //   std::string                                                                  resiName;
115
116 //private:
117 //   std::array<std::size_t, static_cast<std::size_t>(backboneAtomTypes::Count)>  _ResInfo{ 0, 0, 0, 0, 0 };
118 //};
119
120 enum class secondaryStructureTypes : std::size_t { // TODO
121     Loop = 0, // ~
122     Break, // =
123     Bend, // S
124     Turn, // T
125     Helis_PPII, // P
126     Helix_5, // I
127     Helix_3, // G
128     Ladder, // E
129     Bridge, // B
130     Helix_4, // H
131     Count
132
133 };
134
135 enum class turnsTypes : std::size_t {
136     Turn_3 = 0,
137     Turn_4,
138     Turn_5,
139     Count
140 };
141
142 enum class HelixPositions : std::size_t {
143     None = 0,
144     Start,
145     Middle,
146     End,
147     Start_AND_End,
148     Count
149 };
150
151 enum class bridgeTypes : std::size_t {
152     None = 0,
153     AntiParallelBridge,
154     ParallelBridge,
155     Count
156 };
157
158 class secondaryStructures{ // PatterSearch Wrapper
159 public:
160     secondaryStructures();
161     void                    initiateSearch(const std::vector<ResInfo> &ResInfoMatrix, const bool PiHelicesPreferencez = true);
162     std::string             patternSearch();
163     ~secondaryStructures();
164
165     class secondaryStructuresData{ // PatternSearch Tool
166     public:
167         void                                setStatus(const secondaryStructureTypes secondaryStructureTypeName);
168         void                                setStatus(const HelixPositions helixPosition, const turnsTypes turn);
169         bool                                getStatus(const secondaryStructureTypes secondaryStructureTypeName) const, isBreakPartnerWith(const secondaryStructuresData *partner) const;
170         HelixPositions                      getStatus(const turnsTypes turn) const;
171         secondaryStructureTypes             getStatus() const;
172         void                                setBreak(secondaryStructuresData *breakPartner), setBridge(secondaryStructuresData *bridgePartner, bridgeTypes bridgeType);
173         bool                                hasBridges() const, isBridgePartnerWith(secondaryStructuresData *bridgePartner, bridgeTypes bridgeType) const;
174     private:
175         std::array<bool, static_cast<std::size_t>(secondaryStructureTypes::Count)>  SecondaryStructuresStatusArray{ true, 0, 0, 0, 0, 0, 0 };
176         secondaryStructuresData                                        *breakPartners[2]{nullptr, nullptr};
177         secondaryStructuresData                                        *bridgePartners[2]{nullptr, nullptr};
178         secondaryStructureTypes                                         SecondaryStructuresStatus {secondaryStructureTypes::Loop};
179         std::array<HelixPositions, static_cast<std::size_t>(turnsTypes::Count)>  TurnsStatusArray {HelixPositions::None, HelixPositions::None, HelixPositions::None};
180     };
181
182     std::vector<secondaryStructuresData>     SecondaryStructuresStatusMap;
183 private:
184
185     const std::vector<ResInfo>                                       *ResInfoMap;
186
187     const gmx::EnumerationArray<secondaryStructureTypes, const char> secondaryStructureTypeNames = {
188        { '~', '=', 'S', 'T', 'P', 'I', 'G', 'B', 'E', 'H'}
189     };
190     std::string     SecondaryStructuresStringLine;
191
192     bool            hasHBondBetween(std::size_t resi1, std::size_t resi2) const;
193
194     bool            NoChainBreaksBetween(std::size_t ResiStart, std::size_t ResiEnd) const, isLoop(const std::size_t resiIndex) const, PiHelixPreference;
195     bridgeTypes     calculateBridge(std::size_t i, std::size_t j) const;
196     void            analyzeBridgesAndLaddersPatterns(), analyzeTurnsAndHelicesPatterns(), analyzePPHelicesPatterns();
197
198     const float                           HBondEnergyCutOff{ -0.5 };
199 };
200
201 class alternateNeighborhoodSearch{
202 public:
203     alternateNeighborhoodSearch(){}
204     void setCutoff(const real &cutoff_init);
205     void AltPairSearch(const t_trxframe &fr, const std::vector<ResInfo>      &IndexMap);
206     bool findNextPair();
207     std::size_t getResiI() const, getResiJ() const;
208 private:
209     real                                                 cutoff = 4.0;
210     std::size_t resiIpos{0}, resiJpos{0};
211     std::vector<std::vector<bool>> PairMap;
212     std::vector<std::size_t> MiniBoxMap, partners, num_of_miniboxes;
213     std::vector<std::vector<std::size_t>> MiniBoxesReverseMap;
214     std::vector<std::vector<std::vector<std::vector<std::size_t>>>> MiniBoxesMap;
215     std::array<int, 3> MiniBoxSize{0, 0, 0}, offset{-1, -1, -1}, fixBox{0,0,0};
216     void GetMiniBoxesMap(const t_trxframe &fr, const std::vector<ResInfo>      &IndexMap);
217     void FixAtomCoordinates(real &coordinate, const real vector_length);
218     void ReCalculatePBC(int &x, const int &x_max);
219     std::vector<std::vector<bool>>::iterator ResiI;
220     std::vector<bool>::iterator ResiJ;
221     bool init = true;
222 };
223
224 class DsspTool
225 {
226 public:
227    DsspTool();
228    void initAnalysis(/*const TrajectoryAnalysisSettings &settings,*/const TopologyInformation& top, const initParameters &initParamz);
229
230    void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc);
231
232    int                                   nres{};
233
234    std::vector<std::pair<int, std::string>> getData();
235 private:
236    initParameters                        initParams;
237    AtomsDataPtr                          atoms_;
238    std::string                           filename_;
239    void                                  calculateBends(const t_trxframe& fr, const t_pbc* pbc);
240    std::vector<std::string>              ResiNames;
241    std::vector<std::size_t>              AtomResi;
242    std::vector<ResInfo>                  IndexMap;
243    std::vector<std::vector<std::size_t>> nturnsmap, Bridges, AntiBridges;
244    secondaryStructures                   PatternSearch;
245    const float                           HBondEnergyCutOff{ -0.5 },
246                                          minimalCAdistance{ 9.0 };
247    void                                  calculateHBondEnergy(ResInfo& Donor,
248                                                                    ResInfo& Acceptor,
249                                                                    const t_trxframe&          fr,
250                                                                    const t_pbc*               pbc);
251    float CalculateAtomicDistances(const int &A, const int &B, const t_trxframe &fr, const t_pbc *pbc);
252    float CalculateAtomicDistances(const rvec &A, const int &B, const t_trxframe &fr, const t_pbc *pbc);
253
254    class DsspStorage
255    {
256    public:
257        DsspStorage();
258
259        /*Storages our dirty data, duh*/
260        void storageData(int frnr, std::string data);
261
262        /*Clear data after usage for the next pdb file*/
263        void clearAll();
264
265        /*Perform several checks and fucks and shits and returns sorted sexy data*/
266        std::vector<std::pair<int, std::string>> returnData();
267    private:
268        std::vector<std::pair<int, std::string>> storaged_data;
269        static std::mutex mx;
270    };
271
272    static DsspStorage                           Storage;
273 };
274
275 }
276 //end of namespace analysismodules
277
278 }
279 //end of namespace gmx
280
281
282 #endif