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