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