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