73db4033b27f5ca3603bef75692636e5f98ea812
[alexxy/gromacs-dssp.git] / src / dssptools.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 /*
46     There's something wrong with energy calculations of redidues with E ≈ -0
47 */
48
49
50 #include "dssptools.h"
51
52 #include <algorithm>
53 #include "gromacs/math/units.h"
54
55 #include "gromacs/pbcutil/pbc.h"
56 #include <gromacs/trajectoryanalysis.h>
57 #include "gromacs/trajectoryanalysis/topologyinformation.h"
58 #include <set>
59 #include <fstream>
60 #include <mutex>
61 #include <iostream>
62
63 namespace gmx
64 {
65
66 namespace analysismodules
67 {
68
69 //void ResInfo::setIndex(backboneAtomTypes atomTypeName, std::size_t atomIndex)
70 //{
71 //   _ResInfo.at(static_cast<std::size_t>(atomTypeName)) = atomIndex;
72 //}
73
74 //std::size_t ResInfo::getIndex(backboneAtomTypes atomTypeName) const
75 //{
76 //   return _ResInfo[static_cast<std::size_t>(atomTypeName)];
77 //}
78
79 std::size_t ResInfo::getIndex(backboneAtomTypes atomTypeName) const{
80     return _backboneIndices[static_cast<std::size_t>(atomTypeName)];
81 }
82
83 secondaryStructures::secondaryStructures(){
84 }
85 void secondaryStructures::initiateSearch(const std::vector<ResInfo> &ResInfoMatrix, const bool PiHelicesPreferencez){
86     SecondaryStructuresStatusMap.resize(0);
87     SecondaryStructuresStringLine.resize(0);
88     std::vector<std::size_t> temp; temp.resize(0),
89     PiHelixPreference = PiHelicesPreferencez;
90     ResInfoMap = &ResInfoMatrix;
91     SecondaryStructuresStatusMap.resize(ResInfoMatrix.size());
92     SecondaryStructuresStringLine.resize(ResInfoMatrix.size(), '~');
93 }
94
95 void secondaryStructures::secondaryStructuresData::setStatus(const secondaryStructureTypes secondaryStructureTypeName){
96     SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)] = true;
97     SecondaryStructuresStatus = secondaryStructureTypeName;
98 }
99
100 void secondaryStructures::secondaryStructuresData::setStatus(const HelixPositions helixPosition, const turnsTypes turn){
101     TurnsStatusArray[static_cast<std::size_t>(turn)] = helixPosition;
102 }
103
104 bool secondaryStructures::secondaryStructuresData::getStatus(const secondaryStructureTypes secondaryStructureTypeName) const{
105     return SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)];
106 }
107
108 bool secondaryStructures::secondaryStructuresData::isBreakPartnerWith(const secondaryStructuresData *partner) const{
109     return breakPartners[0] == partner || breakPartners[1] == partner;
110 }
111
112 HelixPositions secondaryStructures::secondaryStructuresData::getStatus(const turnsTypes turn) const{
113     return TurnsStatusArray[static_cast<std::size_t>(turn)];
114 }
115
116 secondaryStructureTypes secondaryStructures::secondaryStructuresData::getStatus() const{
117     return  SecondaryStructuresStatus;
118 }
119
120 void secondaryStructures::secondaryStructuresData::setBreak(secondaryStructuresData *breakPartner){
121     if (breakPartners[0] == nullptr){
122         breakPartners[0] = breakPartner;
123     }
124     else{
125         breakPartners[1] = breakPartner;
126     }
127     setStatus(secondaryStructureTypes::Break);
128 }
129
130 bool secondaryStructures::hasHBondBetween(std::size_t Donor, std::size_t Acceptor) const{
131     if( (*ResInfoMap)[Donor].acceptor[0] == nullptr ||
132         (*ResInfoMap)[Donor].acceptor[1] == nullptr ||
133         (*ResInfoMap)[Acceptor].info == nullptr ){
134         std::cout << "Bad hbond check. Reason(s): " ;
135         if ( (*ResInfoMap)[Donor].acceptor[0] == nullptr ){
136             std::cout << "Donor has no acceptor[0]; ";
137         }
138         if ( (*ResInfoMap)[Donor].acceptor[1] == nullptr ){
139             std::cout << "Donor has no acceptor[1]; ";
140         }
141         if ( (*ResInfoMap)[Acceptor].info == nullptr ){
142             std::cout << "No info about acceptor; ";
143         }
144         std::cout << std::endl;
145         return false;
146     }
147 //    else if (!( (*ResInfoMap)[Acceptor].donor[0] == nullptr ||
148 //              (*ResInfoMap)[Acceptor].donor[1] == nullptr ||
149 //              (*ResInfoMap)[Donor].info == nullptr )) {
150 //        std::cout << "Comparing DONOR №" << (*ResInfoMap)[Donor].info->nr << " And ACCEPTOR №" << (*ResInfoMap)[Acceptor].info->nr << ": " << std::endl;
151 //        std::cout << "Donor's acceptors' nr are = " << (*ResInfoMap)[Donor].acceptor[0]->nr << " (chain " << (*ResInfoMap)[Donor].acceptor[0]->chainid << ") , " << (*ResInfoMap)[Donor].acceptor[1]->nr << " (chain " << (*ResInfoMap)[Donor].acceptor[1]->chainid << ")" << std::endl;
152 //        std::cout << "Donor's acceptors' energy are = " << (*ResInfoMap)[Donor].acceptorEnergy[0] << ", " << (*ResInfoMap)[Donor].acceptorEnergy[1] << std::endl;
153 //        std::cout << "Acceptors's donors' nr are = " << (*ResInfoMap)[Acceptor].donor[0]->nr << " (chain " << (*ResInfoMap)[Acceptor].donor[0]->chainid << ") , " << (*ResInfoMap)[Acceptor].donor[1]->nr << " (chain " << (*ResInfoMap)[Acceptor].donor[1]->chainid << ")" << std::endl;
154 //        std::cout << "Acceptors's donors' energy are = " << (*ResInfoMap)[Acceptor].donorEnergy[0] << ", " << (*ResInfoMap)[Acceptor].donorEnergy[1] << std::endl;
155 //        if( ( (*ResInfoMap)[Donor].acceptor[0] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[0] < HBondEnergyCutOff ) ||
156 //                ( (*ResInfoMap)[Donor].acceptor[1] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[1] < HBondEnergyCutOff ) ){
157 //            std::cout << "HBond Exist" << std::endl;
158 //        }
159 //    }
160 //    else {
161 //        std::cout << "Bad hbond check. Reason(s): " ;
162 //        if ( (*ResInfoMap)[Acceptor].donor[0] == nullptr ){
163 //            std::cout << "Acceptor has no donor[0]; ";
164 //        }
165 //        if ( (*ResInfoMap)[Acceptor].donor[1] == nullptr ){
166 //            std::cout << "Acceptor has no donor[1]; ";
167 //        }
168 //        if ( (*ResInfoMap)[Donor].info == nullptr ){
169 //            std::cout << "No info about donor; ";
170 //        }
171 //        std::cout << std::endl;
172 //    }
173
174     return ( (*ResInfoMap)[Donor].acceptor[0] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[0] < HBondEnergyCutOff ) ||
175            ( (*ResInfoMap)[Donor].acceptor[1] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[1] < HBondEnergyCutOff );
176
177
178 }
179
180 bool secondaryStructures::NoChainBreaksBetween(std::size_t Resi1, std::size_t Resi2) const{
181     std::size_t i{Resi1}, j{Resi2}; // From i to j → i <= j
182     if ( i > j ){
183         std::swap(i, j);
184     }
185
186     for (; i != j; ++i){
187         if ( SecondaryStructuresStatusMap[i].isBreakPartnerWith(&SecondaryStructuresStatusMap[i + 1]) && SecondaryStructuresStatusMap[i + 1].isBreakPartnerWith(&SecondaryStructuresStatusMap[i]) ){
188             std::cout << "Patternsearch has detected a CHAINBREAK between " << Resi1 << " and " << Resi2 << std::endl;
189             return false;
190         }
191     }
192     return true;
193 }
194
195 bridgeTypes secondaryStructures::calculateBridge(std::size_t i, std::size_t j) const{
196     if( i < 1 || j < 1 || i + 1 >= ResInfoMap->size() || j + 1 >= ResInfoMap->size() ){
197         return bridgeTypes::None;
198     }
199     if(NoChainBreaksBetween(i - 1, i + 1) && NoChainBreaksBetween(j - 1, j + 1)){
200         if((hasHBondBetween(i + 1, j) && hasHBondBetween(j, i - 1)) || (hasHBondBetween(j + 1, i) && hasHBondBetween(i, j - 1)) ){ //possibly swap
201             return bridgeTypes::ParallelBridge;
202         }
203         else if((hasHBondBetween(i + 1, j - 1) && hasHBondBetween(j + 1, i - 1)) || (hasHBondBetween(j, i) && hasHBondBetween(i, j)) ){ //possibly swap
204             return bridgeTypes::AntiParallelBridge;
205         }
206     }
207     return bridgeTypes::None;
208 }
209
210 void secondaryStructures::analyzeBridgesAndLaddersPatterns(){
211     for(std::size_t i {1}; i + 4 < SecondaryStructuresStatusMap.size(); ++i){
212         for(std::size_t j {i + 3}; j + 1 < SecondaryStructuresStatusMap.size(); ++j ){
213             bridgeTypes type {calculateBridge(i, j)};
214             if (type == bridgeTypes::None){
215                 continue;
216             }
217             bool found {false};
218         }
219     }
220
221
222
223
224
225
226
227
228
229
230
231
232 //    for (std::size_t i{ 1 }; i < HBondsMap.front().size() - 1; ++i){
233 //        for (std::size_t j{ 1 }; j < HBondsMap.front().size() - 1; ++j){
234 //            if (std::abs(static_cast<int>(i) - static_cast<int>(j)) > 2){
235 //                if ((HBondsMap[i - 1][j] && HBondsMap[j][i + 1])    ||
236 //                    (HBondsMap[j - 1][i] && HBondsMap[i][j + 1])){
237 //                    Bridge[i].push_back(j);
238 //                }
239 //                if ((HBondsMap[i][j] && HBondsMap[j][i])    ||
240 //                    (HBondsMap[i - 1][j + 1] && HBondsMap[j - 1][i + 1])){
241 //                    AntiBridge[i].push_back(j);
242 //                }
243 //            }
244 //        }
245 //    }
246 //    for (std::size_t i{ 0 }; i < HBondsMap.front().size(); ++i){
247 //        if ((!Bridge[i].empty() || !AntiBridge[i].empty())){
248 //            setStatus(i, secondaryStructureTypes::Bulge);
249 //        }
250 //    }
251 //    for (std::size_t i{ 2 }; i + 2 < HBondsMap.front().size(); ++i){
252 //        for (std::size_t j { i - 2 }; j <= (i + 2); ++j){
253 //            if (j == i){
254 //                continue;
255 //            }
256 //            else {
257 //                for (std::vector<bridgeTypes>::const_iterator bridge {Bridges.begin()}; bridge != Bridges.end(); ++bridge ){
258 //                    if (!getBridge(*bridge)[i].empty() || !getBridge(*bridge)[j].empty()){
259 //                        for (std::size_t i_resi{ 0 }; i_resi < getBridge(*bridge)[i].size(); ++i_resi){
260 //                            for (std::size_t j_resi{ 0 }; j_resi < getBridge(*bridge)[j].size(); ++j_resi){
261 //                                if (abs(static_cast<int>(getBridge(*bridge)[i][i_resi])
262 //                                        - static_cast<int>(getBridge(*bridge)[j][j_resi]))
263 //                                        && (abs(static_cast<int>(getBridge(*bridge)[i][i_resi])
264 //                                        - static_cast<int>(getBridge(*bridge)[j][j_resi]))
265 //                                        < 5)){
266 //                                    if (j < i){
267 //                                        for (std::size_t k{ 0 }; k <= i - j; ++k){
268 //                                            setStatus(i + k, secondaryStructureTypes::Ladder);
269 //                                        }
270 //                                    }
271 //                                    else{
272 //                                        for (std::size_t k{ 0 }; k <= j - i; ++k){
273 //                                            setStatus(i + k, secondaryStructureTypes::Ladder);
274 //                                        }
275 //                                    }
276 //                                }
277 //                            }
278 //                        }
279 //                    }
280 //                }
281 //            }
282 //        }
283 //    }
284 }
285
286 void secondaryStructures::analyzeTurnsAndHelicesPatterns(){
287     for(const turnsTypes &i : { turnsTypes::Turn_4, turnsTypes::Turn_3, turnsTypes::Turn_5 }){
288         std::size_t stride {static_cast<std::size_t>(i) + 3};
289         std::cout << "Testing Helix_" << stride << std::endl;
290         for(std::size_t j {0}; j + stride < SecondaryStructuresStatusMap.size(); ++j){
291             std::cout << "Testing " << j << " and " << j + stride << std::endl;
292             if ( hasHBondBetween(j + stride, j) && NoChainBreaksBetween(j, j + stride) ){
293                 std::cout << j << " and " << j + stride << " has hbond!" << std::endl;
294                 SecondaryStructuresStatusMap[j + stride].setStatus(HelixPositions::End, i);
295
296                 for (std::size_t k {1}; k < stride; ++k){
297                     if( SecondaryStructuresStatusMap[j + k].getStatus(i) == HelixPositions::None ){
298                         SecondaryStructuresStatusMap[j + k].setStatus(HelixPositions::Middle, i);
299                         SecondaryStructuresStatusMap[j + k].setStatus(secondaryStructureTypes::Turn);
300                     }
301
302                 }
303
304                 if( SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::End ){
305                     SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start_AND_End, i);
306                 }
307                 else {
308                     SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start, i);
309                 }
310             }
311         }
312     }
313
314     for(const turnsTypes &i : { turnsTypes::Turn_4, turnsTypes::Turn_3, turnsTypes::Turn_5 }){
315         std::size_t stride {static_cast<std::size_t>(i) + 3};
316         for(std::size_t j {1}; j + stride < SecondaryStructuresStatusMap.size(); ++j){
317             if ( (SecondaryStructuresStatusMap[j - 1].getStatus(i) == HelixPositions::Start || SecondaryStructuresStatusMap[j - 1].getStatus(i) == HelixPositions::Start_AND_End ) &&
318                  (SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::Start || SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::Start_AND_End ) ){
319                 bool empty = true;
320                 secondaryStructureTypes Helix;
321                 switch(i){
322                 case turnsTypes::Turn_3:
323                     for (std::size_t k {0}; empty && k < stride; ++k){
324                         empty = SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Loop ) || SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_3);
325                     }
326                     Helix = secondaryStructureTypes::Helix_3;
327                     break;
328                 case turnsTypes::Turn_5:
329                     for (std::size_t k {0}; empty && k < stride; ++k){
330                         empty = SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Loop ) || SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_5) || (PiHelixPreference && SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_4)); //TODO
331                     }
332                     Helix = secondaryStructureTypes::Helix_5;
333                     break;
334                 default:
335                     Helix = secondaryStructureTypes::Helix_4;
336                     break;
337                 }
338                 if ( empty || Helix == secondaryStructureTypes::Helix_4 ){
339                     for(std::size_t k {0}; k < stride; ++k ){
340                         SecondaryStructuresStatusMap[j + k].setStatus(Helix);
341                     }
342                 }
343             }
344         }
345     }
346
347 //    for(std::size_t i {1}; i + 1 < SecondaryStructuresStatusMap.size(); ++i){
348 //        if (SecondaryStructuresStatusMap[i].getStatus() == secondaryStructureTypes::Loop || SecondaryStructuresStatusMap[i].getStatus(secondaryStructureTypes::Bend)){
349 //            bool isTurn = false;
350 //            for(const turnsTypes &j : {turnsTypes::Turn_3, turnsTypes::Turn_4, turnsTypes::Turn_5}){
351 //                std::size_t stride {static_cast<std::size_t>(i) + 3};
352 //                for(std::size_t k {1}; k < stride; ++k){
353 //                    isTurn = (i >= k) && (SecondaryStructuresStatusMap[i - k].getStatus(j) == HelixPositions::Start || SecondaryStructuresStatusMap[i - k].getStatus(j) == HelixPositions::Start_AND_End) ;
354 //                }
355 //            }
356
357 //            if (isTurn){
358 //                SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Turn);
359 //            }
360 //        }
361 //    }
362
363 }
364
365 void secondaryStructures::analyzePPHelicesPatterns(){}
366
367 std::string secondaryStructures::patternSearch(){
368
369
370 //    analyzeBridgesAndLaddersPatterns();
371     analyzeTurnsAndHelicesPatterns();
372 //    analyzePPHelicesPatterns();
373
374 //    for(std::size_t i {0}; i < ResInfoMap->size(); ++i){
375 //        std::cout << (*ResInfoMap)[i].info->nr << " " << *((*ResInfoMap)[i].info->name) << std::endl;
376 //    }
377
378 //    std::cout.precision(5);
379 //    for(std::size_t i{0}; i < ResInfoMap->size(); ++i, std::cout << std::endl << std::endl){
380 //        std::cout << (*ResInfoMap)[i].info->nr << " " << *((*ResInfoMap)[i].info->name) ;
381 //        if ( (*ResInfoMap)[i].donor[0] != nullptr ){
382 //            std::cout << " has donor[0] = " << (*ResInfoMap)[i].donor[0]->nr << " " << *((*ResInfoMap)[i].donor[0]->name) << " with E = " << (*ResInfoMap)[i].donorEnergy[0] << " and" ;
383 //        }
384 //        else {
385 //            std::cout << " has no donor[0] and" ;
386 //        }
387 //        if ( (*ResInfoMap)[i].acceptor[0] != nullptr ){
388 //            std::cout << " has acceptor[0] = " << (*ResInfoMap)[i].acceptor[0]->nr << " " << *((*ResInfoMap)[i].acceptor[0]->name) << " with E = " << (*ResInfoMap)[i].acceptorEnergy[0] ;
389 //        }
390 //        else {
391 //            std::cout << " has no acceptor[0]" ;
392 //        }
393 //        std::cout << std::endl << "Also, " << (*ResInfoMap)[i].info->nr << " " << *((*ResInfoMap)[i].info->name);
394 //        if ( (*ResInfoMap)[i].donor[1] != nullptr ){
395 //            std::cout << " has donor[1] = " << (*ResInfoMap)[i].donor[1]->nr << " " << *((*ResInfoMap)[i].donor[1]->name) << " with E = " << (*ResInfoMap)[i].donorEnergy[1] << " and" ;
396 //        }
397 //        else {
398 //            std::cout << " has no donor[1] and" ;
399 //        }
400 //        if ( (*ResInfoMap)[i].acceptor[1] != nullptr ){
401 //            std::cout << " has acceptor[1] = " << (*ResInfoMap)[i].acceptor[1]->nr << " " << *((*ResInfoMap)[i].acceptor[1]->name) << " with E = " << (*ResInfoMap)[i].acceptorEnergy[1] ;
402 //        }
403 //        else {
404 //            std::cout << " has no acceptor[1]" ;
405 //        }
406 //    }
407
408     /*Write Data*/
409
410     for(std::size_t i {static_cast<std::size_t>(secondaryStructureTypes::Bend)}; i != static_cast<std::size_t>(secondaryStructureTypes::Count); ++i){
411         for(std::size_t j {0}; j < SecondaryStructuresStatusMap.size(); ++j){
412             if (SecondaryStructuresStatusMap[j].getStatus(static_cast<secondaryStructureTypes>(i))){
413                 SecondaryStructuresStringLine[j] = secondaryStructureTypeNames[i] ;
414             }
415         }
416     }
417
418     /*Add breaks*/
419
420     if(SecondaryStructuresStatusMap.size() > 1){
421         for(std::size_t i {0}, linefactor{1}; i + 1 < SecondaryStructuresStatusMap.size(); ++i){
422             if( SecondaryStructuresStatusMap[i].getStatus(secondaryStructureTypes::Break) && SecondaryStructuresStatusMap[i + 1].getStatus(secondaryStructureTypes::Break) ){
423                 if(SecondaryStructuresStatusMap[i].isBreakPartnerWith(&SecondaryStructuresStatusMap[i + 1]) && SecondaryStructuresStatusMap[i + 1].isBreakPartnerWith(&SecondaryStructuresStatusMap[i]) ){
424                     SecondaryStructuresStringLine.insert(SecondaryStructuresStringLine.begin() + i + linefactor, secondaryStructureTypeNames[secondaryStructureTypes::Break]);
425                     ++linefactor;
426                 }
427             }
428         }
429     }
430     return SecondaryStructuresStringLine;
431 }
432
433 secondaryStructures::~secondaryStructures(){
434     SecondaryStructuresStatusMap.resize(0);
435     SecondaryStructuresStringLine.resize(0);
436 }
437
438 DsspTool::DsspStorage::DsspStorage(){
439     storaged_data.resize(0);
440 }
441
442 void DsspTool::DsspStorage::clearAll(){
443     storaged_data.resize(0);
444 }
445
446 std::mutex DsspTool::DsspStorage::mx;
447
448 void DsspTool::DsspStorage::storageData(int frnr, std::string data){
449     std::lock_guard<std::mutex> guardian(mx);
450     std::pair<int, std::string> datapair(frnr, data);
451     storaged_data.push_back(datapair);
452 }
453
454 std::vector<std::pair<int, std::string>> DsspTool::DsspStorage::returnData(){
455     std::sort(storaged_data.begin(), storaged_data.end());
456     return storaged_data;
457 }
458
459 void alternateNeighborhoodSearch::setCutoff(const real &cutoff_init){
460     cutoff = cutoff_init;
461 }
462
463 void alternateNeighborhoodSearch::FixAtomCoordinates(real &coordinate, const real vector_length){
464     while (coordinate < 0) {
465         coordinate += vector_length;
466     }
467     while (coordinate >= vector_length) {
468         coordinate -= vector_length;
469     }
470 }
471
472 void alternateNeighborhoodSearch::ReCalculatePBC(int &x, const int &x_max) {
473     if (x < 0) {
474         x += x_max;
475     }
476     if (x >= x_max) {
477         x -= x_max;
478     }
479 }
480
481 void alternateNeighborhoodSearch::GetMiniBoxesMap(const t_trxframe &fr, const std::vector<ResInfo> &IndexMap){
482     rvec coordinates, box_vector_length;
483     num_of_miniboxes.resize(0);
484     num_of_miniboxes.resize(3);
485     for (std::size_t i{XX}; i <= ZZ; ++i) {
486         box_vector_length[i] = std::sqrt(
487                     std::pow(fr.box[i][XX], 2) + std::pow(fr.box[i][YY], 2) + std::pow(fr.box[i][ZZ], 2));
488         num_of_miniboxes[i] = std::floor((box_vector_length[i] / cutoff)) + 1;
489     }
490     MiniBoxesMap.resize(0);
491     MiniBoxesReverseMap.resize(0);
492     MiniBoxesMap.resize(num_of_miniboxes[XX], std::vector<std::vector<std::vector<std::size_t> > >(
493                              num_of_miniboxes[YY], std::vector<std::vector<std::size_t> >(
494                              num_of_miniboxes[ZZ], std::vector<std::size_t>(
495                              0))));
496     MiniBoxesReverseMap.resize(IndexMap.size(), std::vector<std::size_t>(3));
497     for (std::vector<ResInfo>::const_iterator i {IndexMap.begin()}; i != IndexMap.end(); ++i) {
498         for (std::size_t j{XX}; j <= ZZ; ++j) {
499             coordinates[j] = fr.x[i->getIndex(backboneAtomTypes::AtomCA)][j];
500             FixAtomCoordinates(coordinates[j], box_vector_length[j]);
501         }
502         MiniBoxesMap[std::floor(coordinates[XX] / cutoff)][std::floor(coordinates[YY] / cutoff)][std::floor(
503                           coordinates[ZZ] / cutoff)].push_back(i - IndexMap.begin());
504         for (std::size_t j{XX}; j <= ZZ; ++j){
505             MiniBoxesReverseMap[i - IndexMap.begin()][j] = std::floor(coordinates[j] / cutoff);
506         }
507     }
508 }
509
510 void alternateNeighborhoodSearch::AltPairSearch(const t_trxframe &fr, const std::vector<ResInfo> &IndexMap){
511     GetMiniBoxesMap(fr, IndexMap);
512     MiniBoxSize[XX] = MiniBoxesMap.size();
513     MiniBoxSize[YY] = MiniBoxesMap.front().size();
514     MiniBoxSize[ZZ] = MiniBoxesMap.front().front().size();
515     PairMap.resize(0);
516     PairMap.resize(IndexMap.size(), std::vector<bool>(IndexMap.size(), false));
517     ResiI = PairMap.begin();
518     ResiJ = ResiI->begin();
519
520     for (std::vector<ResInfo>::const_iterator i = IndexMap.begin(); i != IndexMap.end(); ++i){
521         for (offset[XX] = -1; offset[XX] <= 1; ++offset[XX]) {
522             for (offset[YY] = -1; offset[YY] <= 1; ++offset[YY]) {
523                 for (offset[ZZ] = -1; offset[ZZ] <= 1; ++offset[ZZ]) {
524                     for (std::size_t k{XX}; k <= ZZ; ++k) {
525                         fixBox[k] = MiniBoxesReverseMap[i - IndexMap.begin()][k] + offset[k];
526                         ReCalculatePBC(fixBox[k], MiniBoxSize[k]);
527                     }
528                     for (std::size_t j{0}; j < MiniBoxesMap[fixBox[XX]][fixBox[YY]][fixBox[ZZ]].size(); ++j) {
529                         if ( (i - IndexMap.begin()) != MiniBoxesMap[fixBox[XX]][fixBox[YY]][fixBox[ZZ]][j]){
530                             PairMap[i - IndexMap.begin()][MiniBoxesMap[fixBox[XX]][fixBox[YY]][fixBox[ZZ]][j]] = true;
531                             PairMap[MiniBoxesMap[fixBox[XX]][fixBox[YY]][fixBox[ZZ]][j]][i - IndexMap.begin()] = true;
532                         }
533                     }
534                 }
535             }
536         }
537     }
538 }
539
540 bool alternateNeighborhoodSearch::findNextPair(){
541
542     if(!PairMap.size()){
543         return false;
544     }
545
546     for(; ResiI != PairMap.end(); ++ResiI, ResiJ = ResiI->begin() ){
547         for(; ResiJ != ResiI->end(); ++ResiJ){
548             if(*ResiJ){
549                 resiIpos = ResiI - PairMap.begin();
550                 resiJpos = ResiJ - ResiI->begin();
551                 if ( ResiJ != ResiI->end() ){
552                     ++ResiJ;
553                 }
554                 else if (ResiI != PairMap.end()) {
555                     ++ResiI;
556                     ResiJ = ResiI->begin();
557                 }
558                 else {
559                     return false; // ???
560                 }
561                 return true;
562             }
563         }
564     }
565
566     return false;
567 }
568
569 std::size_t alternateNeighborhoodSearch::getResiI() const {
570     return resiIpos;
571 }
572
573 std::size_t alternateNeighborhoodSearch::getResiJ() const {
574     return resiJpos;
575 }
576
577
578 DsspTool::DsspStorage DsspTool::Storage;
579
580 DsspTool::DsspTool(){
581 }
582
583 void DsspTool::calculateBends(const t_trxframe &fr, const t_pbc *pbc)
584 {
585    const float benddegree{ 70.0 }, maxdist{ 2.5 };
586    float       degree{ 0 }, vdist{ 0 }, vprod{ 0 };
587    gmx::RVec   a{ 0, 0, 0 }, b{ 0, 0, 0 };
588    for (std::size_t i{ 0 }; i + 1 < IndexMap.size(); ++i)
589    {
590        if (CalculateAtomicDistances(static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
591                                     static_cast<int>(IndexMap[i + 1].getIndex(backboneAtomTypes::AtomN)),
592                                     fr,
593                                     pbc)
594            > maxdist)
595        {
596            PatternSearch.SecondaryStructuresStatusMap[i].setBreak(&PatternSearch.SecondaryStructuresStatusMap[i + 1]);
597            PatternSearch.SecondaryStructuresStatusMap[i + 1].setBreak(&PatternSearch.SecondaryStructuresStatusMap[i]);
598
599 //           std::cout << "Break between " << i + 1 << " and " << i + 2 << std::endl;
600        }
601    }
602    for (std::size_t i{ 2 }; i + 2 < IndexMap.size() ; ++i)
603    {
604        if (PatternSearch.SecondaryStructuresStatusMap[i - 2].getStatus(secondaryStructureTypes::Break) ||
605            PatternSearch.SecondaryStructuresStatusMap[i - 1].getStatus(secondaryStructureTypes::Break) ||
606            PatternSearch.SecondaryStructuresStatusMap[i].getStatus(secondaryStructureTypes::Break) ||
607            PatternSearch.SecondaryStructuresStatusMap[i + 1].getStatus(secondaryStructureTypes::Break)
608           )
609        {
610            continue;
611        }
612        for (int j{ 0 }; j < 3; ++j)
613        {
614            a[j] = fr.x[IndexMap[i].getIndex(backboneAtomTypes::AtomCA)][j]
615                   - fr.x[IndexMap[i - 2].getIndex(backboneAtomTypes::AtomCA)][j];
616            b[j] = fr.x[IndexMap[i + 2].getIndex(backboneAtomTypes::AtomCA)][j]
617                   - fr.x[IndexMap[i].getIndex(backboneAtomTypes::AtomCA)][j];
618        }
619        vdist = (a[0] * b[0]) + (a[1] * b[1]) + (a[2] * b[2]);
620        vprod = CalculateAtomicDistances(IndexMap[i - 2].getIndex(backboneAtomTypes::AtomCA),
621                                         IndexMap[i].getIndex(backboneAtomTypes::AtomCA),
622                                         fr,
623                                         pbc)
624                * gmx::c_angstrom / gmx::c_nano
625                * CalculateAtomicDistances(IndexMap[i].getIndex(backboneAtomTypes::AtomCA),
626                                           IndexMap[i + 2].getIndex(backboneAtomTypes::AtomCA),
627                                           fr,
628                                           pbc)
629                * gmx::c_angstrom / gmx::c_nano;
630        degree = std::acos(vdist / vprod) * gmx::c_rad2Deg;
631        if (degree > benddegree)
632        {
633            PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Bend);
634        }
635    }
636 }
637
638 void DsspTool::calculateHBondEnergy(ResInfo& Donor,
639                        ResInfo& Acceptor,
640                        const t_trxframe&          fr,
641                        const t_pbc*               pbc)
642 {
643    /*
644     * DSSP uses eq from dssp 2.x
645     * kCouplingConstant = 27.888,  //  = 332 * 0.42 * 0.2
646     * E = k * (1/rON + 1/rCH - 1/rOH - 1/rCN) where CO comes from one AA and NH from another
647     * if R is in A
648     * Hbond if E < -0.5
649     *
650     * For the note, H-Bond Donor is N-H («Donor of H») and H-Bond Acceptor is C=O («Acceptor of H»)
651     *
652     */
653
654     if (CalculateAtomicDistances(
655                 Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc)
656         >= minimalCAdistance)
657     {
658         return void();
659     }
660
661     const float kCouplingConstant = 27.888;
662     const float minimalAtomDistance{ 0.5 },
663             minEnergy{ -9.9 };
664     float HbondEnergy{ 0 };
665     float distanceNO{ 0 }, distanceHC{ 0 }, distanceHO{ 0 }, distanceNC{ 0 };
666
667 //    std::cout << "For Donor №" << Donor.info->nr - 1 << " and Accpetor №" << Acceptor.info->nr - 1 << std::endl;
668
669     if( !(Donor.is_proline) && (Acceptor.getIndex(backboneAtomTypes::AtomC) && Acceptor.getIndex(backboneAtomTypes::AtomO)
670                                 && Donor.getIndex(backboneAtomTypes::AtomN) && ( Donor.getIndex(backboneAtomTypes::AtomH) || initParams.addHydrogens ) ) ){ // TODO
671         distanceNO = CalculateAtomicDistances(
672                Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
673         distanceNC = CalculateAtomicDistances(
674                Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
675         if (initParams.addHydrogens){
676             if (Donor.prevResi != nullptr && Donor.prevResi->getIndex(backboneAtomTypes::AtomC) && Donor.prevResi->getIndex(backboneAtomTypes::AtomO)){
677                rvec atomH{};
678                float prevCODist {CalculateAtomicDistances(Donor.prevResi->getIndex(backboneAtomTypes::AtomC), Donor.prevResi->getIndex(backboneAtomTypes::AtomO), fr, pbc)};
679                for (int i{XX}; i <= ZZ; ++i){
680                    float prevCO = fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomC)][i] - fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomO)][i];
681                    atomH[i] = fr.x[Donor.getIndex(backboneAtomTypes::AtomH)][i]; // Но на самом деле берутся координаты N
682                    atomH[i] += prevCO / prevCODist;
683                }
684                distanceHO = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
685                distanceHC = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
686             }
687             else{
688                 distanceHO = distanceNO;
689                 distanceHC = distanceNC;
690             }
691        }
692        else {
693            distanceHO = CalculateAtomicDistances(
694                    Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
695            distanceHC = CalculateAtomicDistances(
696                    Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
697        }
698        if ((distanceNO < minimalAtomDistance) || (distanceHC < minimalAtomDistance)
699         || (distanceHO < minimalAtomDistance) || (distanceNC < minimalAtomDistance))
700        {
701             HbondEnergy = minEnergy;
702        }
703        else{
704            HbondEnergy =
705                    kCouplingConstant
706                    * ((1 / distanceNO) + (1 / distanceHC) - (1 / distanceHO) - (1 / distanceNC));
707        }
708
709 //       std::cout << "CA-CA distance: " << CalculateAtomicDistances(
710 //                        Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc) << std::endl;
711 //       std::cout << "N-O distance: " << distanceNO << std::endl;
712 //       std::cout << "N-C distance: " << distanceNC << std::endl;
713 //       std::cout << "H-O distance: " << distanceHO << std::endl;
714 //       std::cout << "H-C distance: " << distanceHC << std::endl;
715
716        HbondEnergy = std::round(HbondEnergy * 1000) / 1000;
717
718 //       if ( HbondEnergy < minEnergy ){ // I don't think that this is correct
719 //            HbondEnergy = minEnergy;
720 //       }
721
722 //       std::cout << "Calculated energy = " << HbondEnergy << std::endl;
723     }
724 //    else{
725 //        std::cout << "Donor Is Proline" << std::endl;
726 //    }
727
728     if (HbondEnergy < Donor.acceptorEnergy[0]){
729            Donor.acceptor[1] = Donor.acceptor[0];
730            Donor.acceptor[0] = Acceptor.info;
731            Donor.acceptorEnergy[0] = HbondEnergy;
732     }
733     else if (HbondEnergy < Donor.acceptorEnergy[1]){
734            Donor.acceptor[1] = Acceptor.info;
735            Donor.acceptorEnergy[1] = HbondEnergy;
736     }
737
738     if (HbondEnergy < Acceptor.donorEnergy[0]){
739            Acceptor.donor[1] = Acceptor.donor[0];
740            Acceptor.donor[0] = Donor.info;
741            Acceptor.donorEnergy[0] = HbondEnergy;
742     }
743     else if (HbondEnergy < Acceptor.donorEnergy[1]){
744            Acceptor.donor[1] = Donor.info;
745            Acceptor.donorEnergy[1] = HbondEnergy;
746     }
747 }
748
749
750 /* Calculate Distance From B to A */
751 float DsspTool::CalculateAtomicDistances(const int &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
752 {
753    gmx::RVec r{ 0, 0, 0 };
754    pbc_dx(pbc, fr.x[A], fr.x[B], r.as_vec());
755    return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
756 }
757
758 /* Calculate Distance From B to A, where A is only fake coordinates */
759 float DsspTool::CalculateAtomicDistances(const rvec &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
760 {
761    gmx::RVec r{ 0, 0, 0 };
762    pbc_dx(pbc, A, fr.x[B], r.as_vec());
763    return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
764 }
765
766 void DsspTool::initAnalysis(/*const TrajectoryAnalysisSettings &settings,*/const TopologyInformation& top, const initParameters &initParamz)
767 {
768    initParams = initParamz;
769    ResInfo _backboneAtoms;
770    std::size_t                 i{ 0 };
771    std::string proLINE;
772    int resicompare{ top.atoms()->atom[static_cast<std::size_t>(*(initParams.sel_.atomIndices().begin()))].resind };
773    IndexMap.resize(0);
774    IndexMap.push_back(_backboneAtoms);
775    IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
776    proLINE = *(IndexMap[i].info->name);
777    if( proLINE.compare("PRO") == 0 ){
778        IndexMap[i].is_proline = true;
779    }
780
781    for (gmx::ArrayRef<const int>::iterator ai{ initParams.sel_.atomIndices().begin() }; (ai != initParams.sel_.atomIndices().end()); ++ai){
782        if (resicompare != top.atoms()->atom[static_cast<std::size_t>(*ai)].resind)
783        {
784            ++i;
785            resicompare = top.atoms()->atom[static_cast<std::size_t>(*ai)].resind;
786            IndexMap.emplace_back(_backboneAtoms);
787            IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
788            proLINE = *(IndexMap[i].info->name);
789            if( proLINE.compare("PRO") == 0 ){
790                IndexMap[i].is_proline = true;
791            }
792
793        }
794        std::string atomname(*(top.atoms()->atomname[static_cast<std::size_t>(*ai)]));
795        if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA])
796        {
797            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomCA)] = *ai;
798        }
799        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC])
800        {
801            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomC)] = *ai;
802        }
803        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO])
804        {
805            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomO)] = *ai;
806        }
807        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN])
808        {
809            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomN)] = *ai;
810            if (initParamz.addHydrogens == true){
811                IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
812            }
813        }
814        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH] && initParamz.addHydrogens == false) // Юзать водород в структуре
815        {
816            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
817        }
818
819
820
821 //       if( atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO]
822 //       || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH]){
823 //           std::cout << "Atom " << atomname << " №" << *ai << " From Resi " << *(top.atoms()->resinfo[i].name) << " №" << resicompare << std::endl;
824 //       }
825    }
826
827    for (std::size_t j {1}; j < IndexMap.size(); ++j){
828        IndexMap[j].prevResi = &(IndexMap[j - 1]);
829
830        IndexMap[j - 1].nextResi = &(IndexMap[j]);
831
832 //           std::cout << "Resi " << IndexMap[i].info->nr << *(IndexMap[i].info->name) << std::endl;
833 //           std::cout << "Prev resi is " << IndexMap[i].prevResi->info->nr << *(IndexMap[i].prevResi->info->name) << std::endl;
834 //           std::cout << "Prev resi's next resi is " << IndexMap[i - 1].nextResi->info->nr << *(IndexMap[i - 1].nextResi->info->name) << std::endl;
835 //         std::cout << IndexMap[j].prevResi->info->nr;
836 //         std::cout << *(IndexMap[j].prevResi->info->name) ;
837 //         std::cout << " have CA = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomCA) ;
838 //         std::cout << " C = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomC);
839 //         std::cout << " O = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomO);
840 //         std::cout << " N = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomN);
841 //         std::cout << " H = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomH) << std::endl;
842    }
843
844    nres = i + 1;
845 }
846
847 void DsspTool::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc)
848 {
849
850     switch(initParams.NBS){
851     case (NBSearchMethod::Classique): {
852
853         // store positions of CA atoms to use them for nbSearch
854         std::vector<gmx::RVec> positionsCA_;
855         for (std::size_t i{ 0 }; i < IndexMap.size(); ++i)
856         {
857             positionsCA_.emplace_back(fr.x[IndexMap[i].getIndex(backboneAtomTypes::AtomCA)]);
858         }
859
860         AnalysisNeighborhood nb_;
861         nb_.setCutoff(initParams.cutoff_);
862         AnalysisNeighborhoodPositions       nbPos_(positionsCA_);
863         gmx::AnalysisNeighborhoodSearch     start      = nb_.initSearch(pbc, nbPos_);
864         gmx::AnalysisNeighborhoodPairSearch pairSearch = start.startPairSearch(nbPos_);
865         gmx::AnalysisNeighborhoodPair       pair;
866         while (pairSearch.findNextPair(&pair))
867         {
868             if(CalculateAtomicDistances(
869                         IndexMap[pair.refIndex()].getIndex(backboneAtomTypes::AtomCA), IndexMap[pair.testIndex()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
870                 < minimalCAdistance){
871                 calculateHBondEnergy(IndexMap[pair.refIndex()], IndexMap[pair.testIndex()], fr, pbc);
872                 if (IndexMap[pair.testIndex()].info != IndexMap[pair.refIndex() + 1].info){
873                     calculateHBondEnergy(IndexMap[pair.testIndex()], IndexMap[pair.refIndex()], fr, pbc);
874                 }
875             }
876         }
877
878         break;
879     }
880     case (NBSearchMethod::Experimental): { // TODO FIX
881
882         alternateNeighborhoodSearch as_;
883
884         as_.setCutoff(initParams.cutoff_);
885
886         as_.AltPairSearch(fr, IndexMap);
887
888         while (as_.findNextPair()){
889             if(CalculateAtomicDistances(
890                         IndexMap[as_.getResiI()].getIndex(backboneAtomTypes::AtomCA), IndexMap[as_.getResiJ()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
891                 < minimalCAdistance){
892                 calculateHBondEnergy(IndexMap[as_.getResiI()], IndexMap[as_.getResiJ()], fr, pbc);
893                 if (IndexMap[as_.getResiJ()].info != IndexMap[as_.getResiI() + 1].info){
894                     calculateHBondEnergy(IndexMap[as_.getResiJ()], IndexMap[as_.getResiI()], fr, pbc);
895                 }
896             }
897         }
898
899         break;
900     }
901     default: {
902
903         for(std::vector<ResInfo>::iterator Donor {IndexMap.begin()}; Donor != IndexMap.end() ; ++Donor){
904             for(std::vector<ResInfo>::iterator Acceptor {Donor + 1} ; Acceptor != IndexMap.end() ; ++Acceptor){
905                 if(CalculateAtomicDistances(
906                             Donor->getIndex(backboneAtomTypes::AtomCA), Acceptor->getIndex(backboneAtomTypes::AtomCA), fr, pbc)
907                     < minimalCAdistance){
908                     calculateHBondEnergy(*Donor, *Acceptor, fr, pbc);
909                     if (Acceptor != Donor + 1){
910                         calculateHBondEnergy(*Acceptor, *Donor, fr, pbc);
911                     }
912                 }
913             }
914         }
915         break;
916     }
917     }
918
919
920 //    for(std::size_t i {0}; i < IndexMap.size(); ++i){
921 //        std::cout << IndexMap[i].info->nr << " " << *(IndexMap[i].info->name) << std::endl;
922 //    }
923
924    PatternSearch.initiateSearch(IndexMap, initParams.PPHelices);
925    calculateBends(fr, pbc);
926    Storage.storageData(frnr, PatternSearch.patternSearch());
927
928 }
929
930 std::vector<std::pair<int, std::string>> DsspTool::getData(){
931     return Storage.returnData();
932 }
933
934 } // namespace analysismodules
935
936 } // namespace gmx