ac28492238d6905e5f8ff1c68bd6ab268b20ca3e
[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, const int _pp_stretch){
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     pp_stretch = _pp_stretch;
92     SecondaryStructuresStatusMap.resize(ResInfoMatrix.size());
93     SecondaryStructuresStringLine.resize(ResInfoMatrix.size(), '~');
94 }
95
96 void secondaryStructures::secondaryStructuresData::setStatus(const secondaryStructureTypes secondaryStructureTypeName, bool status){
97     SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)] = status;
98     if (status){
99         SecondaryStructuresStatus = secondaryStructureTypeName;
100     }
101     else{
102         SecondaryStructuresStatus = secondaryStructureTypes::Loop;
103     }
104 }
105
106 void secondaryStructures::secondaryStructuresData::setStatus(const HelixPositions helixPosition, const turnsTypes turn){
107     TurnsStatusArray[static_cast<std::size_t>(turn)] = helixPosition;
108 }
109
110 bool secondaryStructures::secondaryStructuresData::getStatus(const secondaryStructureTypes secondaryStructureTypeName) const{
111     return SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)];
112 }
113
114 bool secondaryStructures::secondaryStructuresData::isBreakPartnerWith(const secondaryStructuresData *partner) const{
115     return breakPartners[0] == partner || breakPartners[1] == partner;
116 }
117
118 HelixPositions secondaryStructures::secondaryStructuresData::getStatus(const turnsTypes turn) const{
119     return TurnsStatusArray[static_cast<std::size_t>(turn)];
120 }
121
122 secondaryStructureTypes secondaryStructures::secondaryStructuresData::getStatus() const{
123     return  SecondaryStructuresStatus;
124 }
125
126 void secondaryStructures::secondaryStructuresData::setBreak(secondaryStructuresData *breakPartner){
127     if (breakPartners[0] == nullptr){
128         breakPartners[0] = breakPartner;
129     }
130     else{
131         breakPartners[1] = breakPartner;
132     }
133     setStatus(secondaryStructureTypes::Break);
134 }
135
136 void secondaryStructures::secondaryStructuresData::setBridge(secondaryStructuresData *bridgePartner, std::size_t bridgePartnerIndex, bridgeTypes bridgeType){
137     if(bridgeType == bridgeTypes::ParallelBridge){
138         bridgePartners[0] = bridgePartner;
139         bridgePartnersIndexes[0] = bridgePartnerIndex;
140     }
141     else if (bridgeType == bridgeTypes::AntiParallelBridge){
142         bridgePartners[1] = bridgePartner;
143         bridgePartnersIndexes[1] = bridgePartnerIndex;
144     }
145 }
146
147 bool secondaryStructures::secondaryStructuresData::hasBridges() const{
148     return bridgePartners[0] || bridgePartners[1];
149
150 }
151
152 bool secondaryStructures::secondaryStructuresData::hasBridges(bridgeTypes bridgeType) const{
153     if(bridgeType == bridgeTypes::ParallelBridge){
154         return bridgePartners[0] != nullptr;
155     }
156     else if (bridgeType == bridgeTypes::AntiParallelBridge){
157         return bridgePartners[1] != nullptr;
158     }
159
160     return false;
161
162 }
163
164 bool secondaryStructures::secondaryStructuresData::isBridgePartnerWith(secondaryStructuresData *bridgePartner, bridgeTypes bridgeType) const{
165     if(bridgeType == bridgeTypes::ParallelBridge){
166         return bridgePartners[0] == bridgePartner;
167     }
168     else if (bridgeType == bridgeTypes::AntiParallelBridge){
169         return bridgePartners[1] == bridgePartner;
170     }
171     return false;
172 }
173
174 std::size_t secondaryStructures::secondaryStructuresData::getBridgePartnerIndex(bridgeTypes bridgeType) const{
175     if (bridgeType == bridgeTypes::ParallelBridge){
176         return bridgePartnersIndexes[0];
177     }
178     return bridgePartnersIndexes[1];
179 }
180
181 secondaryStructures::secondaryStructuresData secondaryStructures::secondaryStructuresData::getBridgePartner(bridgeTypes bridgeType) const{
182     if(bridgeType == bridgeTypes::ParallelBridge){
183         return *(bridgePartners[0]);
184     }
185     return *(bridgePartners[1]);
186 }
187
188 bool secondaryStructures::hasHBondBetween(std::size_t Donor, std::size_t Acceptor) const{
189     if( (*ResInfoMap)[Donor].acceptor[0] == nullptr ||
190         (*ResInfoMap)[Donor].acceptor[1] == nullptr ||
191         (*ResInfoMap)[Acceptor].info == nullptr ){
192 //        std::cout << "Bad hbond check. Reason(s): " ;
193 //        if ( (*ResInfoMap)[Donor].acceptor[0] == nullptr ){
194 //            std::cout << "Donor has no acceptor[0]; ";
195 //        }
196 //        if ( (*ResInfoMap)[Donor].acceptor[1] == nullptr ){
197 //            std::cout << "Donor has no acceptor[1]; ";
198 //        }
199 //        if ( (*ResInfoMap)[Acceptor].info == nullptr ){
200 //            std::cout << "No info about acceptor; ";
201 //        }
202 //        std::cout << std::endl;
203         return false;
204     }
205 //    else if (!( (*ResInfoMap)[Acceptor].donor[0] == nullptr ||
206 //              (*ResInfoMap)[Acceptor].donor[1] == nullptr ||
207 //              (*ResInfoMap)[Donor].info == nullptr )) {
208 //        std::cout << "Comparing DONOR №" << (*ResInfoMap)[Donor].info->nr << " And ACCEPTOR №" << (*ResInfoMap)[Acceptor].info->nr << ": " << std::endl;
209 //        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;
210 //        std::cout << "Donor's acceptors' energy are = " << (*ResInfoMap)[Donor].acceptorEnergy[0] << ", " << (*ResInfoMap)[Donor].acceptorEnergy[1] << std::endl;
211 //        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;
212 //        std::cout << "Acceptors's donors' energy are = " << (*ResInfoMap)[Acceptor].donorEnergy[0] << ", " << (*ResInfoMap)[Acceptor].donorEnergy[1] << std::endl;
213 //        if( ( (*ResInfoMap)[Donor].acceptor[0] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[0] < HBondEnergyCutOff ) ||
214 //                ( (*ResInfoMap)[Donor].acceptor[1] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[1] < HBondEnergyCutOff ) ){
215 //            std::cout << "HBond Exist" << std::endl;
216 //        }
217 //    }
218 //    else {
219 //        std::cout << "Bad hbond check. Reason(s): " ;
220 //        if ( (*ResInfoMap)[Acceptor].donor[0] == nullptr ){
221 //            std::cout << "Acceptor has no donor[0]; ";
222 //        }
223 //        if ( (*ResInfoMap)[Acceptor].donor[1] == nullptr ){
224 //            std::cout << "Acceptor has no donor[1]; ";
225 //        }
226 //        if ( (*ResInfoMap)[Donor].info == nullptr ){
227 //            std::cout << "No info about donor; ";
228 //        }
229 //        std::cout << std::endl;
230 //    }
231
232     return ( (*ResInfoMap)[Donor].acceptor[0] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[0] < HBondEnergyCutOff ) ||
233            ( (*ResInfoMap)[Donor].acceptor[1] == (*ResInfoMap)[Acceptor].info && (*ResInfoMap)[Donor].acceptorEnergy[1] < HBondEnergyCutOff );
234
235
236 }
237
238 bool secondaryStructures::NoChainBreaksBetween(std::size_t Resi1, std::size_t Resi2) const{
239     std::size_t i{Resi1}, j{Resi2}; // From i to j → i <= j
240     if ( i > j ){
241         std::swap(i, j);
242     }
243
244     for (; i != j; ++i){
245         if ( SecondaryStructuresStatusMap[i].isBreakPartnerWith(&SecondaryStructuresStatusMap[i + 1]) && SecondaryStructuresStatusMap[i + 1].isBreakPartnerWith(&SecondaryStructuresStatusMap[i]) ){
246 //            std::cout << "Patternsearch has detected a CHAINBREAK between " << Resi1 << " and " << Resi2 << std::endl;
247             return false;
248         }
249     }
250     return true;
251 }
252
253 bridgeTypes secondaryStructures::calculateBridge(std::size_t i, std::size_t j) const{
254     if( i < 1 || j < 1 || i + 1 >= ResInfoMap->size() || j + 1 >= ResInfoMap->size() ){ // Protection from idiotz, не обязательно нужен
255         return bridgeTypes::None;
256     }
257
258     ResInfo a{(*ResInfoMap)[i]}, b{(*ResInfoMap)[j]};
259
260     if(NoChainBreaksBetween(i - 1, i + 1) && NoChainBreaksBetween(j - 1, j + 1) && a.prevResi && a.nextResi && b.prevResi && b.nextResi){
261         if((hasHBondBetween(i + 1, j) && hasHBondBetween(j, i - 1)) || (hasHBondBetween(j + 1, i) && hasHBondBetween(i, j - 1)) ){
262             return bridgeTypes::ParallelBridge;
263         }
264         else if((hasHBondBetween(i + 1, j - 1) && hasHBondBetween(j + 1, i - 1)) || (hasHBondBetween(j, i) && hasHBondBetween(i, j)) ){
265             return bridgeTypes::AntiParallelBridge;
266         }
267     }
268     return bridgeTypes::None;
269 }
270
271 void secondaryStructures::analyzeBridgesAndStrandsPatterns(){
272 //  Calculate Bridgez
273     for(std::size_t i {1}; i + 4 < SecondaryStructuresStatusMap.size(); ++i){
274         for(std::size_t j {i + 3}; j + 1 < SecondaryStructuresStatusMap.size(); ++j ){
275             switch(calculateBridge(i, j)){
276             case bridgeTypes::ParallelBridge : {
277                 SecondaryStructuresStatusMap[i].setBridge(&(SecondaryStructuresStatusMap[j]), j, bridgeTypes::ParallelBridge);
278                 SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Bridge);
279                 SecondaryStructuresStatusMap[j].setBridge(&(SecondaryStructuresStatusMap[i]), i, bridgeTypes::ParallelBridge);
280                 SecondaryStructuresStatusMap[j].setStatus(secondaryStructureTypes::Bridge);
281             }
282             case bridgeTypes::AntiParallelBridge : {
283                 SecondaryStructuresStatusMap[i].setBridge(&(SecondaryStructuresStatusMap[j]), j, bridgeTypes::AntiParallelBridge);
284                 SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Bridge);
285                 SecondaryStructuresStatusMap[j].setBridge(&(SecondaryStructuresStatusMap[i]), i, bridgeTypes::AntiParallelBridge);
286                 SecondaryStructuresStatusMap[j].setStatus(secondaryStructureTypes::Bridge);
287             }
288             default :
289                 continue;
290             }
291         }
292     }
293 //  Calculate Extended Strandz
294     for(std::size_t i {1}; i + 1 < SecondaryStructuresStatusMap.size(); ++i){
295         bool is_Estrand {false};
296         for(std::size_t j {2}; j != 0; --j){
297             for(const bridgeTypes &bridgeType : {bridgeTypes::ParallelBridge, bridgeTypes::AntiParallelBridge}){
298                 if (SecondaryStructuresStatusMap[i].hasBridges(bridgeType) && SecondaryStructuresStatusMap[i + j].hasBridges(bridgeType) ){
299                     std::size_t i_partner{SecondaryStructuresStatusMap[i].getBridgePartnerIndex(bridgeType)}, j_partner{SecondaryStructuresStatusMap[i + j].getBridgePartnerIndex(bridgeType)}, second_strand{};
300                     if ( abs(i_partner - j_partner) < 6){
301                         if (i_partner < j_partner){
302                             second_strand = i_partner;
303                         }
304                         else{
305                             second_strand = j_partner;
306                         }
307                         for(int k{abs(i_partner - j_partner)}; k >= 0; --k){
308                             if (SecondaryStructuresStatusMap[second_strand + k].getStatus(secondaryStructureTypes::Bridge)){
309                                 SecondaryStructuresStatusMap[second_strand + k].setStatus(secondaryStructureTypes::Bridge, false);
310                             }
311
312                             SecondaryStructuresStatusMap[second_strand + k].setStatus(secondaryStructureTypes::Strand);
313                         }
314                         is_Estrand = true;
315                     }
316                 }
317             }
318             if (is_Estrand){
319                 for(std::size_t k{0}; k <= j; ++k){
320                     if (SecondaryStructuresStatusMap[i + k].getStatus(secondaryStructureTypes::Bridge)){
321                         SecondaryStructuresStatusMap[i + k].setStatus(secondaryStructureTypes::Bridge, false);
322                     }
323                     SecondaryStructuresStatusMap[i + k].setStatus(secondaryStructureTypes::Strand);
324                 }
325                 break;
326             }
327         }
328     }
329
330
331
332
333
334
335
336
337
338
339
340
341 //    for (std::size_t i{ 1 }; i < HBondsMap.front().size() - 1; ++i){
342 //        for (std::size_t j{ 1 }; j < HBondsMap.front().size() - 1; ++j){
343 //            if (std::abs(static_cast<int>(i) - static_cast<int>(j)) > 2){
344 //                if ((HBondsMap[i - 1][j] && HBondsMap[j][i + 1])    ||
345 //                    (HBondsMap[j - 1][i] && HBondsMap[i][j + 1])){
346 //                    Bridge[i].push_back(j);
347 //                }
348 //                if ((HBondsMap[i][j] && HBondsMap[j][i])    ||
349 //                    (HBondsMap[i - 1][j + 1] && HBondsMap[j - 1][i + 1])){
350 //                    AntiBridge[i].push_back(j);
351 //                }
352 //            }
353 //        }
354 //    }
355 //    for (std::size_t i{ 0 }; i < HBondsMap.front().size(); ++i){
356 //        if ((!Bridge[i].empty() || !AntiBridge[i].empty())){
357 //            setStatus(i, secondaryStructureTypes::Bulge);
358 //        }
359 //    }
360 //    for (std::size_t i{ 2 }; i + 2 < HBondsMap.front().size(); ++i){
361 //        for (std::size_t j { i - 2 }; j <= (i + 2); ++j){
362 //            if (j == i){
363 //                continue;
364 //            }
365 //            else {
366 //                for (std::vector<bridgeTypes>::const_iterator bridge {Bridges.begin()}; bridge != Bridges.end(); ++bridge ){
367 //                    if (!getBridge(*bridge)[i].empty() || !getBridge(*bridge)[j].empty()){
368 //                        for (std::size_t i_resi{ 0 }; i_resi < getBridge(*bridge)[i].size(); ++i_resi){
369 //                            for (std::size_t j_resi{ 0 }; j_resi < getBridge(*bridge)[j].size(); ++j_resi){
370 //                                if (abs(static_cast<int>(getBridge(*bridge)[i][i_resi])
371 //                                        - static_cast<int>(getBridge(*bridge)[j][j_resi]))
372 //                                        && (abs(static_cast<int>(getBridge(*bridge)[i][i_resi])
373 //                                        - static_cast<int>(getBridge(*bridge)[j][j_resi]))
374 //                                        < 5)){
375 //                                    if (j < i){
376 //                                        for (std::size_t k{ 0 }; k <= i - j; ++k){
377 //                                            setStatus(i + k, secondaryStructureTypes::Ladder);
378 //                                        }
379 //                                    }
380 //                                    else{
381 //                                        for (std::size_t k{ 0 }; k <= j - i; ++k){
382 //                                            setStatus(i + k, secondaryStructureTypes::Ladder);
383 //                                        }
384 //                                    }
385 //                                }
386 //                            }
387 //                        }
388 //                    }
389 //                }
390 //            }
391 //        }
392 //    }
393 }
394
395 void secondaryStructures::analyzeTurnsAndHelicesPatterns(){
396     for(const turnsTypes &i : { turnsTypes::Turn_4, turnsTypes::Turn_3, turnsTypes::Turn_5 }){
397         std::size_t stride {static_cast<std::size_t>(i) + 3};
398         for(std::size_t j {0}; j + stride < SecondaryStructuresStatusMap.size(); ++j){
399             if ( hasHBondBetween(j + stride, j) && NoChainBreaksBetween(j, j + stride) ){
400                 SecondaryStructuresStatusMap[j + stride].setStatus(HelixPositions::End, i);
401
402                 for (std::size_t k {1}; k < stride; ++k){
403                     if( SecondaryStructuresStatusMap[j + k].getStatus(i) == HelixPositions::None ){
404                         SecondaryStructuresStatusMap[j + k].setStatus(HelixPositions::Middle, i);
405                         SecondaryStructuresStatusMap[j + k].setStatus(secondaryStructureTypes::Turn);
406                     }
407                 }
408
409                 if( SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::End ){
410                     SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start_AND_End, i);
411                 }
412                 else {
413                     SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start, i);
414                 }
415             }
416         }
417     }
418
419     for(const turnsTypes &i : { turnsTypes::Turn_4, turnsTypes::Turn_3, turnsTypes::Turn_5 }){
420         std::size_t stride {static_cast<std::size_t>(i) + 3};
421         for(std::size_t j {1}; j + stride < SecondaryStructuresStatusMap.size(); ++j){
422             if ( (SecondaryStructuresStatusMap[j - 1].getStatus(i) == HelixPositions::Start || SecondaryStructuresStatusMap[j - 1].getStatus(i) == HelixPositions::Start_AND_End ) &&
423                  (SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::Start || SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::Start_AND_End ) ){
424                 bool empty = true;
425                 secondaryStructureTypes Helix;
426                 switch(i){
427                 case turnsTypes::Turn_3:
428                     for (std::size_t k {0}; empty && k < stride; ++k){
429                         empty = SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Loop ) || SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_3);
430                     }
431                     Helix = secondaryStructureTypes::Helix_3;
432                     break;
433                 case turnsTypes::Turn_5:
434                     for (std::size_t k {0}; empty && k < stride; ++k){
435                         empty = SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Loop ) || SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_5) || (PiHelixPreference && SecondaryStructuresStatusMap[j + k].getStatus(secondaryStructureTypes::Helix_4)); //TODO
436                     }
437                     Helix = secondaryStructureTypes::Helix_5;
438                     break;
439                 default:
440                     Helix = secondaryStructureTypes::Helix_4;
441                     break;
442                 }
443                 if ( empty || Helix == secondaryStructureTypes::Helix_4 ){
444                     for(std::size_t k {0}; k < stride; ++k ){
445                         SecondaryStructuresStatusMap[j + k].setStatus(Helix);
446                     }
447                 }
448             }
449         }
450     }
451
452     /* Не знаю зач они в дссп так сделали, этож полное говно */
453
454 //    for(std::size_t i {1}; i + 1 < SecondaryStructuresStatusMap.size(); ++i){
455 //        if (SecondaryStructuresStatusMap[i].getStatus() == secondaryStructureTypes::Loop || SecondaryStructuresStatusMap[i].getStatus(secondaryStructureTypes::Bend)){
456 //            bool isTurn = false;
457 //            for(const turnsTypes &j : {turnsTypes::Turn_3, turnsTypes::Turn_4, turnsTypes::Turn_5}){
458 //                std::size_t stride {static_cast<std::size_t>(j) + 3};
459 //                for(std::size_t k {1}; k < stride and !isTurn; ++k){
460 //                    isTurn = (i >= k) && (SecondaryStructuresStatusMap[i - k].getStatus(j) == HelixPositions::Start || SecondaryStructuresStatusMap[i - k].getStatus(j) == HelixPositions::Start_AND_End) ;
461 //                }
462 //            }
463 //            if (isTurn){
464 //                SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Turn);
465 //            }
466 //        }
467 //    }
468
469 }
470
471 std::string secondaryStructures::patternSearch(){
472
473
474     analyzeBridgesAndStrandsPatterns();
475     analyzeTurnsAndHelicesPatterns();
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 float DsspTool::CalculateDihedralAngle(const int &A, const int &B, const int &C, const int &D, const t_trxframe &fr, const t_pbc *pbc){
742     float result{360}, u{}, v{};
743     gmx::RVec vBA{}, vDC{}, vBD{}, p{}, x{}, y{};
744     pbc_dx(pbc, fr.x[A], fr.x[B], vBA.as_vec());
745     pbc_dx(pbc, fr.x[C], fr.x[D], vDC.as_vec());
746     pbc_dx(pbc, fr.x[D], fr.x[B], vBD.as_vec());
747
748     for(std::size_t i{XX}, j{j + 1}, k{i + 2}; i <= ZZ; ++i, ++j, ++k){
749         if (j > 2){
750             j -= 3;
751         }
752         if (k > 2){
753             k -= 3;
754         }
755         p[i] = (vBD[j] * vBA[k]) - (vBD[k] * vBA[j]);
756         x[i] = (vBD[j] * vDC[k]) - (vBD[k] * vDC[j]);
757
758     }
759     for(std::size_t i{XX}, j{j + 1}, k{i + 2}; i <= ZZ; ++i, ++j, ++k){
760         if (j > 2){
761             j -= 3;
762         }
763         if (k > 2){
764             k -= 3;
765         }
766         y[i] = (vBD[j] * x[k]) - (vBD[k] * x[j]);
767     }
768
769     u = (x[XX] * x[XX]) + (x[YY] * x[YY]) + (x[ZZ] * x[ZZ]);
770     v = (y[XX] * y[XX]) + (y[YY] * y[YY]) + (y[ZZ] * y[ZZ]);
771
772     if (u > 0 and v > 0){
773         u = ((p[XX] * x[XX]) + (p[YY] * x[YY]) + (p[ZZ] * x[ZZ])) / std::sqrt(u);
774         v = ((p[XX] * y[XX]) + (p[YY] * y[YY]) + (p[ZZ] * y[ZZ])) / std::sqrt(v);
775         if (u != 0 or v != 0){
776             result = std::atan2(v, u) * gmx::c_rad2Deg;
777         }
778     }
779     return result;
780 }
781
782 void DsspTool::calculateDihedrals(const t_trxframe &fr, const t_pbc *pbc){
783     const float epsilon = 29;
784     const float phi_min = -75 - epsilon;
785     const float phi_max = -75 + epsilon;
786     const float psi_min = 145 - epsilon;
787     const float psi_max = 145 + epsilon;
788     std::vector<float>          phi, psi;
789     phi.resize(0);
790     psi.resize(0);
791     phi.resize(IndexMap.size(), 360);
792     psi.resize(IndexMap.size(), 360);
793
794     for (std::size_t i = 1; i + 1 < IndexMap.size(); ++i){ // TODO add index verifictaion (check if those atom indexes exist)
795         phi[i] = CalculateDihedralAngle(static_cast<int>(IndexMap[i - 1].getIndex(backboneAtomTypes::AtomC)),
796                                         static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomN)),
797                                         static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomCA)),
798                                         static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
799                                         fr,
800                                         pbc);
801         psi[i] = CalculateDihedralAngle(static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomN)),
802                                         static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomCA)),
803                                         static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
804                                         static_cast<int>(IndexMap[i + 1].getIndex(backboneAtomTypes::AtomN)),
805                                         fr,
806                                         pbc);
807     }
808
809     for (std::size_t i = 1; i + 3 < IndexMap.size(); ++i){
810         switch (initParams.pp_stretch){
811         case 2: {
812             if (phi_min > phi[i] or phi[i] > phi_max or
813                 phi_min > phi[i + 1] or phi[i + 1]> phi_max){
814                 continue;
815             }
816
817             if (psi_min > psi[i] or psi[i] > psi_max or
818                 psi_min > psi[i + 1] or psi[i + 1] > psi_max){
819                 continue;
820             }
821
822             switch (PatternSearch.SecondaryStructuresStatusMap[i].getStatus(turnsTypes::Turn_PP)){
823                 case HelixPositions::None:
824                     PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start, turnsTypes::Turn_PP);
825                     break;
826
827                 case HelixPositions::End:
828                     PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start_AND_End, turnsTypes::Turn_PP);
829                     break;
830
831                 default:
832                     break;
833             }
834
835             PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(HelixPositions::End, turnsTypes::Turn_PP);
836             /* Пропустил проверку того, что заменяемая ак - петля */
837             PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Helix_PP);
838             PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(secondaryStructureTypes::Helix_PP);
839         break;
840         }
841         case 3:{
842             if (phi_min > phi[i] or phi[i] > phi_max or
843                 phi_min > phi[i + 1] or phi[i + 1]> phi_max or
844                 phi_min > phi[i + 2] or phi[i + 2]> phi_max){
845                 continue;
846             }
847
848             if (psi_min > psi[i] or psi[i] > psi_max or
849                 psi_min > psi[i + 1] or psi[i + 1] > psi_max or
850                 psi_min > psi[i + 2] or psi[i + 2] > psi_max){
851                 continue;
852             }
853
854             switch (PatternSearch.SecondaryStructuresStatusMap[i].getStatus(turnsTypes::Turn_PP)){
855                 case HelixPositions::None:
856                     PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start, turnsTypes::Turn_PP);
857                     break;
858
859                 case HelixPositions::End:
860                     PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start_AND_End, turnsTypes::Turn_PP);
861                     break;
862
863                 default:
864                     break;
865             }
866
867             PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(HelixPositions::Middle, turnsTypes::Turn_PP);
868             PatternSearch.SecondaryStructuresStatusMap[i + 2].setStatus(HelixPositions::End, turnsTypes::Turn_PP);
869             /* Пропустил проверку того, что заменяемая ак - петля */
870             PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Helix_PP);
871             PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(secondaryStructureTypes::Helix_PP);
872             PatternSearch.SecondaryStructuresStatusMap[i + 2].setStatus(secondaryStructureTypes::Helix_PP);
873
874             break;
875         }
876         default:
877             throw std::runtime_error("Unsupported stretch length");
878         }
879     }
880
881 }
882
883 void DsspTool::calculateHBondEnergy(ResInfo& Donor,
884                        ResInfo& Acceptor,
885                        const t_trxframe&          fr,
886                        const t_pbc*               pbc)
887 {
888    /*
889     * DSSP uses eq from dssp 2.x
890     * kCouplingConstant = 27.888,  //  = 332 * 0.42 * 0.2
891     * E = k * (1/rON + 1/rCH - 1/rOH - 1/rCN) where CO comes from one AA and NH from another
892     * if R is in A
893     * Hbond if E < -0.5
894     *
895     * For the note, H-Bond Donor is N-H («Donor of H») and H-Bond Acceptor is C=O («Acceptor of H»)
896     *
897     */
898
899     if (CalculateAtomicDistances(
900                 Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc)
901         >= minimalCAdistance)
902     {
903         return void();
904     }
905
906     const float kCouplingConstant = 27.888;
907     const float minimalAtomDistance{ 0.5 },
908             minEnergy{ -9.9 };
909     float HbondEnergy{ 0 };
910     float distanceNO{ 0 }, distanceHC{ 0 }, distanceHO{ 0 }, distanceNC{ 0 };
911
912 //    std::cout << "For Donor №" << Donor.info->nr - 1 << " and Accpetor №" << Acceptor.info->nr - 1 << std::endl;
913
914     if( !(Donor.is_proline) && (Acceptor.getIndex(backboneAtomTypes::AtomC) && Acceptor.getIndex(backboneAtomTypes::AtomO)
915                                 && Donor.getIndex(backboneAtomTypes::AtomN) && ( Donor.getIndex(backboneAtomTypes::AtomH) || initParams.addHydrogens ) ) ){ // TODO
916         distanceNO = CalculateAtomicDistances(
917                Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
918         distanceNC = CalculateAtomicDistances(
919                Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
920         if (initParams.addHydrogens){
921             if (Donor.prevResi != nullptr && Donor.prevResi->getIndex(backboneAtomTypes::AtomC) && Donor.prevResi->getIndex(backboneAtomTypes::AtomO)){
922                rvec atomH{};
923                float prevCODist {CalculateAtomicDistances(Donor.prevResi->getIndex(backboneAtomTypes::AtomC), Donor.prevResi->getIndex(backboneAtomTypes::AtomO), fr, pbc)};
924                for (int i{XX}; i <= ZZ; ++i){
925                    float prevCO = fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomC)][i] - fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomO)][i];
926                    atomH[i] = fr.x[Donor.getIndex(backboneAtomTypes::AtomH)][i]; // Но на самом деле берутся координаты N
927                    atomH[i] += prevCO / prevCODist;
928                }
929                distanceHO = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
930                distanceHC = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
931             }
932             else{
933                 distanceHO = distanceNO;
934                 distanceHC = distanceNC;
935             }
936        }
937        else {
938            distanceHO = CalculateAtomicDistances(
939                    Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
940            distanceHC = CalculateAtomicDistances(
941                    Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
942        }
943        if ((distanceNO < minimalAtomDistance) || (distanceHC < minimalAtomDistance)
944         || (distanceHO < minimalAtomDistance) || (distanceNC < minimalAtomDistance))
945        {
946             HbondEnergy = minEnergy;
947        }
948        else{
949            HbondEnergy =
950                    kCouplingConstant
951                    * ((1 / distanceNO) + (1 / distanceHC) - (1 / distanceHO) - (1 / distanceNC));
952        }
953
954 //       std::cout << "CA-CA distance: " << CalculateAtomicDistances(
955 //                        Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc) << std::endl;
956 //       std::cout << "N-O distance: " << distanceNO << std::endl;
957 //       std::cout << "N-C distance: " << distanceNC << std::endl;
958 //       std::cout << "H-O distance: " << distanceHO << std::endl;
959 //       std::cout << "H-C distance: " << distanceHC << std::endl;
960
961        HbondEnergy = std::round(HbondEnergy * 1000) / 1000;
962
963 //       if ( HbondEnergy < minEnergy ){ // I don't think that this is correct
964 //            HbondEnergy = minEnergy;
965 //       }
966
967 //       std::cout << "Calculated energy = " << HbondEnergy << std::endl;
968     }
969 //    else{
970 //        std::cout << "Donor Is Proline" << std::endl;
971 //    }
972
973     if (HbondEnergy < Donor.acceptorEnergy[0]){
974            Donor.acceptor[1] = Donor.acceptor[0];
975            Donor.acceptor[0] = Acceptor.info;
976            Donor.acceptorEnergy[0] = HbondEnergy;
977     }
978     else if (HbondEnergy < Donor.acceptorEnergy[1]){
979            Donor.acceptor[1] = Acceptor.info;
980            Donor.acceptorEnergy[1] = HbondEnergy;
981     }
982
983     if (HbondEnergy < Acceptor.donorEnergy[0]){
984            Acceptor.donor[1] = Acceptor.donor[0];
985            Acceptor.donor[0] = Donor.info;
986            Acceptor.donorEnergy[0] = HbondEnergy;
987     }
988     else if (HbondEnergy < Acceptor.donorEnergy[1]){
989            Acceptor.donor[1] = Donor.info;
990            Acceptor.donorEnergy[1] = HbondEnergy;
991     }
992 }
993
994
995 /* Calculate Distance From B to A */
996 float DsspTool::CalculateAtomicDistances(const int &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
997 {
998    gmx::RVec r{ 0, 0, 0 };
999    pbc_dx(pbc, fr.x[A], fr.x[B], r.as_vec());
1000    return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
1001 }
1002
1003 /* Calculate Distance From B to A, where A is only fake coordinates */
1004 float DsspTool::CalculateAtomicDistances(const rvec &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
1005 {
1006    gmx::RVec r{ 0, 0, 0 };
1007    pbc_dx(pbc, A, fr.x[B], r.as_vec());
1008    return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
1009 }
1010
1011 void DsspTool::initAnalysis(/*const TrajectoryAnalysisSettings &settings,*/const TopologyInformation& top, const initParameters &initParamz)
1012 {
1013    initParams = initParamz;
1014    ResInfo _backboneAtoms;
1015    std::size_t                 i{ 0 };
1016    std::string proLINE;
1017    int resicompare{ top.atoms()->atom[static_cast<std::size_t>(*(initParams.sel_.atomIndices().begin()))].resind };
1018    IndexMap.resize(0);
1019    IndexMap.push_back(_backboneAtoms);
1020    IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
1021    proLINE = *(IndexMap[i].info->name);
1022    if( proLINE.compare("PRO") == 0 ){
1023        IndexMap[i].is_proline = true;
1024    }
1025
1026    for (gmx::ArrayRef<const int>::iterator ai{ initParams.sel_.atomIndices().begin() }; (ai != initParams.sel_.atomIndices().end()); ++ai){
1027        if (resicompare != top.atoms()->atom[static_cast<std::size_t>(*ai)].resind)
1028        {
1029            ++i;
1030            resicompare = top.atoms()->atom[static_cast<std::size_t>(*ai)].resind;
1031            IndexMap.emplace_back(_backboneAtoms);
1032            IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
1033            proLINE = *(IndexMap[i].info->name);
1034            if( proLINE.compare("PRO") == 0 ){
1035                IndexMap[i].is_proline = true;
1036            }
1037
1038        }
1039        std::string atomname(*(top.atoms()->atomname[static_cast<std::size_t>(*ai)]));
1040        if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA])
1041        {
1042            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomCA)] = *ai;
1043        }
1044        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC])
1045        {
1046            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomC)] = *ai;
1047        }
1048        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO])
1049        {
1050            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomO)] = *ai;
1051        }
1052        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN])
1053        {
1054            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomN)] = *ai;
1055            if (initParamz.addHydrogens == true){
1056                IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
1057            }
1058        }
1059        else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH] && initParamz.addHydrogens == false) // Юзать водород в структуре
1060        {
1061            IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
1062        }
1063
1064
1065
1066 //       if( atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO]
1067 //       || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH]){
1068 //           std::cout << "Atom " << atomname << " №" << *ai << " From Resi " << *(top.atoms()->resinfo[i].name) << " №" << resicompare << std::endl;
1069 //       }
1070    }
1071
1072    for (std::size_t j {1}; j < IndexMap.size(); ++j){
1073        IndexMap[j].prevResi = &(IndexMap[j - 1]);
1074
1075        IndexMap[j - 1].nextResi = &(IndexMap[j]);
1076
1077 //           std::cout << "Resi " << IndexMap[i].info->nr << *(IndexMap[i].info->name) << std::endl;
1078 //           std::cout << "Prev resi is " << IndexMap[i].prevResi->info->nr << *(IndexMap[i].prevResi->info->name) << std::endl;
1079 //           std::cout << "Prev resi's next resi is " << IndexMap[i - 1].nextResi->info->nr << *(IndexMap[i - 1].nextResi->info->name) << std::endl;
1080 //         std::cout << IndexMap[j].prevResi->info->nr;
1081 //         std::cout << *(IndexMap[j].prevResi->info->name) ;
1082 //         std::cout << " have CA = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomCA) ;
1083 //         std::cout << " C = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomC);
1084 //         std::cout << " O = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomO);
1085 //         std::cout << " N = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomN);
1086 //         std::cout << " H = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomH) << std::endl;
1087    }
1088
1089    nres = i + 1;
1090 }
1091
1092 void DsspTool::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc)
1093 {
1094
1095     switch(initParams.NBS){
1096     case (NBSearchMethod::Classique): {
1097
1098         // store positions of CA atoms to use them for nbSearch
1099         std::vector<gmx::RVec> positionsCA_;
1100         for (std::size_t i{ 0 }; i < IndexMap.size(); ++i)
1101         {
1102             positionsCA_.emplace_back(fr.x[IndexMap[i].getIndex(backboneAtomTypes::AtomCA)]);
1103         }
1104
1105         AnalysisNeighborhood nb_;
1106         nb_.setCutoff(initParams.cutoff_);
1107         AnalysisNeighborhoodPositions       nbPos_(positionsCA_);
1108         gmx::AnalysisNeighborhoodSearch     start      = nb_.initSearch(pbc, nbPos_);
1109         gmx::AnalysisNeighborhoodPairSearch pairSearch = start.startPairSearch(nbPos_);
1110         gmx::AnalysisNeighborhoodPair       pair;
1111         while (pairSearch.findNextPair(&pair))
1112         {
1113             if(CalculateAtomicDistances(
1114                         IndexMap[pair.refIndex()].getIndex(backboneAtomTypes::AtomCA), IndexMap[pair.testIndex()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1115                 < minimalCAdistance){
1116                 calculateHBondEnergy(IndexMap[pair.refIndex()], IndexMap[pair.testIndex()], fr, pbc);
1117                 if (IndexMap[pair.testIndex()].info != IndexMap[pair.refIndex() + 1].info){
1118                     calculateHBondEnergy(IndexMap[pair.testIndex()], IndexMap[pair.refIndex()], fr, pbc);
1119                 }
1120             }
1121         }
1122
1123         break;
1124     }
1125     case (NBSearchMethod::Experimental): { // TODO FIX
1126
1127         alternateNeighborhoodSearch as_;
1128
1129         as_.setCutoff(initParams.cutoff_);
1130
1131         as_.AltPairSearch(fr, IndexMap);
1132
1133         while (as_.findNextPair()){
1134             if(CalculateAtomicDistances(
1135                         IndexMap[as_.getResiI()].getIndex(backboneAtomTypes::AtomCA), IndexMap[as_.getResiJ()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1136                 < minimalCAdistance){
1137                 calculateHBondEnergy(IndexMap[as_.getResiI()], IndexMap[as_.getResiJ()], fr, pbc);
1138                 if (IndexMap[as_.getResiJ()].info != IndexMap[as_.getResiI() + 1].info){
1139                     calculateHBondEnergy(IndexMap[as_.getResiJ()], IndexMap[as_.getResiI()], fr, pbc);
1140                 }
1141             }
1142         }
1143
1144         break;
1145     }
1146     default: {
1147
1148         for(std::vector<ResInfo>::iterator Donor {IndexMap.begin()}; Donor != IndexMap.end() ; ++Donor){
1149             for(std::vector<ResInfo>::iterator Acceptor {Donor + 1} ; Acceptor != IndexMap.end() ; ++Acceptor){
1150                 if(CalculateAtomicDistances(
1151                             Donor->getIndex(backboneAtomTypes::AtomCA), Acceptor->getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1152                     < minimalCAdistance){
1153                     calculateHBondEnergy(*Donor, *Acceptor, fr, pbc);
1154                     if (Acceptor != Donor + 1){
1155                         calculateHBondEnergy(*Acceptor, *Donor, fr, pbc);
1156                     }
1157                 }
1158             }
1159         }
1160         break;
1161     }
1162     }
1163
1164
1165 //    for(std::size_t i {0}; i < IndexMap.size(); ++i){
1166 //        std::cout << IndexMap[i].info->nr << " " << *(IndexMap[i].info->name) << std::endl;
1167 //    }
1168
1169    PatternSearch.initiateSearch(IndexMap, initParams.PPHelices, initParams.pp_stretch);
1170    calculateBends(fr, pbc);
1171    calculateDihedrals(fr, pbc);
1172    Storage.storageData(frnr, PatternSearch.patternSearch());
1173
1174 }
1175
1176 std::vector<std::pair<int, std::string>> DsspTool::getData(){
1177     return Storage.returnData();
1178 }
1179
1180 } // namespace analysismodules
1181
1182 } // namespace gmx