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