2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * Implements gmx::analysismodules::Trajectory.
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
46 There's something wrong with energy calculations of redidues with E ≈ -0
50 #include "dssptools.h"
53 #include "gromacs/math/units.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include <gromacs/trajectoryanalysis.h>
57 #include "gromacs/trajectoryanalysis/topologyinformation.h"
66 namespace analysismodules
69 //void ResInfo::setIndex(backboneAtomTypes atomTypeName, std::size_t atomIndex)
71 // _ResInfo.at(static_cast<std::size_t>(atomTypeName)) = atomIndex;
74 //std::size_t ResInfo::getIndex(backboneAtomTypes atomTypeName) const
76 // return _ResInfo[static_cast<std::size_t>(atomTypeName)];
79 std::size_t ResInfo::getIndex(backboneAtomTypes atomTypeName) const{
80 return _backboneIndices[static_cast<std::size_t>(atomTypeName)];
83 secondaryStructures::secondaryStructures(){
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(), '~');
96 void secondaryStructures::secondaryStructuresData::setStatus(const secondaryStructureTypes secondaryStructureTypeName, bool status){
97 SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)] = status;
99 SecondaryStructuresStatus = secondaryStructureTypeName;
102 SecondaryStructuresStatus = secondaryStructureTypes::Loop;
106 void secondaryStructures::secondaryStructuresData::setStatus(const HelixPositions helixPosition, const turnsTypes turn){
107 TurnsStatusArray[static_cast<std::size_t>(turn)] = helixPosition;
110 bool secondaryStructures::secondaryStructuresData::getStatus(const secondaryStructureTypes secondaryStructureTypeName) const{
111 return SecondaryStructuresStatusArray[static_cast<std::size_t>(secondaryStructureTypeName)];
114 bool secondaryStructures::secondaryStructuresData::isBreakPartnerWith(const secondaryStructuresData *partner) const{
115 return breakPartners[0] == partner || breakPartners[1] == partner;
118 HelixPositions secondaryStructures::secondaryStructuresData::getStatus(const turnsTypes turn) const{
119 return TurnsStatusArray[static_cast<std::size_t>(turn)];
122 secondaryStructureTypes secondaryStructures::secondaryStructuresData::getStatus() const{
123 return SecondaryStructuresStatus;
126 void secondaryStructures::secondaryStructuresData::setBreak(secondaryStructuresData *breakPartner){
127 if (breakPartners[0] == nullptr){
128 breakPartners[0] = breakPartner;
131 breakPartners[1] = breakPartner;
133 setStatus(secondaryStructureTypes::Break);
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;
141 else if (bridgeType == bridgeTypes::AntiParallelBridge){
142 bridgePartners[1] = bridgePartner;
143 bridgePartnersIndexes[1] = bridgePartnerIndex;
147 bool secondaryStructures::secondaryStructuresData::hasBridges() const{
148 return bridgePartners[0] || bridgePartners[1];
152 bool secondaryStructures::secondaryStructuresData::hasBridges(bridgeTypes bridgeType) const{
153 if(bridgeType == bridgeTypes::ParallelBridge){
154 return bridgePartners[0] != nullptr;
156 else if (bridgeType == bridgeTypes::AntiParallelBridge){
157 return bridgePartners[1] != nullptr;
164 bool secondaryStructures::secondaryStructuresData::isBridgePartnerWith(secondaryStructuresData *bridgePartner, bridgeTypes bridgeType) const{
165 if(bridgeType == bridgeTypes::ParallelBridge){
166 return bridgePartners[0] == bridgePartner;
168 else if (bridgeType == bridgeTypes::AntiParallelBridge){
169 return bridgePartners[1] == bridgePartner;
174 std::size_t secondaryStructures::secondaryStructuresData::getBridgePartnerIndex(bridgeTypes bridgeType) const{
175 if (bridgeType == bridgeTypes::ParallelBridge){
176 return bridgePartnersIndexes[0];
178 return bridgePartnersIndexes[1];
181 secondaryStructures::secondaryStructuresData secondaryStructures::secondaryStructuresData::getBridgePartner(bridgeTypes bridgeType) const{
182 if(bridgeType == bridgeTypes::ParallelBridge){
183 return *(bridgePartners[0]);
185 return *(bridgePartners[1]);
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]; ";
196 // if ( (*ResInfoMap)[Donor].acceptor[1] == nullptr ){
197 // std::cout << "Donor has no acceptor[1]; ";
199 // if ( (*ResInfoMap)[Acceptor].info == nullptr ){
200 // std::cout << "No info about acceptor; ";
202 // std::cout << std::endl;
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;
219 // std::cout << "Bad hbond check. Reason(s): " ;
220 // if ( (*ResInfoMap)[Acceptor].donor[0] == nullptr ){
221 // std::cout << "Acceptor has no donor[0]; ";
223 // if ( (*ResInfoMap)[Acceptor].donor[1] == nullptr ){
224 // std::cout << "Acceptor has no donor[1]; ";
226 // if ( (*ResInfoMap)[Donor].info == nullptr ){
227 // std::cout << "No info about donor; ";
229 // std::cout << std::endl;
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 );
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
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;
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;
258 ResInfo a{(*ResInfoMap)[i]}, b{(*ResInfoMap)[j]};
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;
264 else if((hasHBondBetween(i + 1, j - 1) && hasHBondBetween(j + 1, i - 1)) || (hasHBondBetween(j, i) && hasHBondBetween(i, j)) ){
265 return bridgeTypes::AntiParallelBridge;
268 return bridgeTypes::None;
271 void secondaryStructures::analyzeBridgesAndStrandsPatterns(){
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);
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);
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;
305 second_strand = j_partner;
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);
312 SecondaryStructuresStatusMap[second_strand + k].setStatus(secondaryStructureTypes::Strand);
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);
323 SecondaryStructuresStatusMap[i + k].setStatus(secondaryStructureTypes::Strand);
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);
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);
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);
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){
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]))
376 // for (std::size_t k{ 0 }; k <= i - j; ++k){
377 // setStatus(i + k, secondaryStructureTypes::Ladder);
381 // for (std::size_t k{ 0 }; k <= j - i; ++k){
382 // setStatus(i + k, secondaryStructureTypes::Ladder);
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);
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);
409 if( SecondaryStructuresStatusMap[j].getStatus(i) == HelixPositions::End ){
410 SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start_AND_End, i);
413 SecondaryStructuresStatusMap[j].setStatus(HelixPositions::Start, i);
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 ) ){
425 secondaryStructureTypes Helix;
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);
431 Helix = secondaryStructureTypes::Helix_3;
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
437 Helix = secondaryStructureTypes::Helix_5;
440 Helix = secondaryStructureTypes::Helix_4;
443 if ( empty || Helix == secondaryStructureTypes::Helix_4 ){
444 for(std::size_t k {0}; k < stride; ++k ){
445 SecondaryStructuresStatusMap[j + k].setStatus(Helix);
452 /* Не знаю зач они в дссп так сделали, этож полное говно */
454 for(std::size_t i {1}; i + 1 < SecondaryStructuresStatusMap.size(); ++i){
455 if (static_cast<int>(SecondaryStructuresStatusMap[i].getStatus()) <= static_cast<int>(secondaryStructureTypes::Turn)){
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) ;
464 SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Turn);
471 std::string secondaryStructures::patternSearch(){
474 analyzeBridgesAndStrandsPatterns();
475 analyzeTurnsAndHelicesPatterns();
477 // for(std::size_t i {0}; i < ResInfoMap->size(); ++i){
478 // std::cout << (*ResInfoMap)[i].info->nr << " " << *((*ResInfoMap)[i].info->name) << std::endl;
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" ;
488 // std::cout << " has no donor[0] and" ;
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] ;
494 // std::cout << " has no acceptor[0]" ;
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" ;
501 // std::cout << " has no donor[1] and" ;
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] ;
507 // std::cout << " has no acceptor[1]" ;
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] ;
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]);
533 return SecondaryStructuresStringLine;
536 secondaryStructures::~secondaryStructures(){
537 SecondaryStructuresStatusMap.resize(0);
538 SecondaryStructuresStringLine.resize(0);
541 DsspTool::DsspStorage::DsspStorage(){
542 storaged_data.resize(0);
545 void DsspTool::DsspStorage::clearAll(){
546 storaged_data.resize(0);
549 std::mutex DsspTool::DsspStorage::mx;
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);
557 std::vector<std::pair<int, std::string>> DsspTool::DsspStorage::returnData(){
558 std::sort(storaged_data.begin(), storaged_data.end());
559 return storaged_data;
562 void alternateNeighborhoodSearch::setCutoff(const real &cutoff_init){
563 cutoff = cutoff_init;
566 void alternateNeighborhoodSearch::FixAtomCoordinates(real &coordinate, const real vector_length){
567 while (coordinate < 0) {
568 coordinate += vector_length;
570 while (coordinate >= vector_length) {
571 coordinate -= vector_length;
575 void alternateNeighborhoodSearch::ReCalculatePBC(int &x, const int &x_max) {
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;
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>(
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]);
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);
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();
619 PairMap.resize(IndexMap.size(), std::vector<bool>(IndexMap.size(), false));
620 ResiI = PairMap.begin();
621 ResiJ = ResiI->begin();
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]);
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;
643 bool alternateNeighborhoodSearch::findNextPair(){
649 for(; ResiI != PairMap.end(); ++ResiI, ResiJ = ResiI->begin() ){
650 for(; ResiJ != ResiI->end(); ++ResiJ){
652 resiIpos = ResiI - PairMap.begin();
653 resiJpos = ResiJ - ResiI->begin();
654 if ( ResiJ != ResiI->end() ){
657 else if (ResiI != PairMap.end()) {
659 ResiJ = ResiI->begin();
672 std::size_t alternateNeighborhoodSearch::getResiI() const {
676 std::size_t alternateNeighborhoodSearch::getResiJ() const {
681 DsspTool::DsspStorage DsspTool::Storage;
683 DsspTool::DsspTool(){
686 void DsspTool::calculateBends(const t_trxframe &fr, const t_pbc *pbc)
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)
693 if (CalculateAtomicDistances(static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
694 static_cast<int>(IndexMap[i + 1].getIndex(backboneAtomTypes::AtomN)),
699 PatternSearch.SecondaryStructuresStatusMap[i].setBreak(&PatternSearch.SecondaryStructuresStatusMap[i + 1]);
700 PatternSearch.SecondaryStructuresStatusMap[i + 1].setBreak(&PatternSearch.SecondaryStructuresStatusMap[i]);
702 // std::cout << "Break between " << i + 1 << " and " << i + 2 << std::endl;
705 for (std::size_t i{ 2 }; i + 2 < IndexMap.size() ; ++i)
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)
715 for (int j{ 0 }; j < 3; ++j)
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];
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),
727 * gmx::c_angstrom / gmx::c_nano
728 * CalculateAtomicDistances(IndexMap[i].getIndex(backboneAtomTypes::AtomCA),
729 IndexMap[i + 2].getIndex(backboneAtomTypes::AtomCA),
732 * gmx::c_angstrom / gmx::c_nano;
733 degree = std::acos(vdist / vprod) * gmx::c_rad2Deg;
734 if (degree > benddegree)
736 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Bend);
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 v12{}, v43{}, z{}, p{}, x{}, y{};
744 pbc_dx(pbc, fr.x[A], fr.x[B], v12.as_vec());
745 pbc_dx(pbc, fr.x[D], fr.x[C], v43.as_vec());
746 pbc_dx(pbc, fr.x[B], fr.x[C], z.as_vec());
748 for(std::size_t j {XX}; j <= ZZ; ++j){
749 v12[j] *= gmx::c_nm2A;
750 v43[j] *= gmx::c_nm2A;
753 for(std::size_t i{XX}, j{j + 1}, k{i + 2}; i <= ZZ; ++i, ++j, ++k){
760 p[i] = (z[j] * v12[k]) - (z[k] * v12[j]);
761 x[i] = (z[j] * v43[k]) - (z[k] * v43[j]);
764 for(std::size_t i{XX}, j{j + 1}, k{i + 2}; i <= ZZ; ++i, ++j, ++k){
771 y[i] = (z[j] * x[k]) - (z[k] * x[j]);
774 // std::cout << "v12 = " << v12[0] << ", " << v12[1] << ", " << v12[2] << std::endl;
775 // std::cout << "v43 = " << v43[0] << ", " << v43[1] << ", " << v43[2] << std::endl;
776 // std::cout << "z = " << z[0] << ", " << z[1] << ", " << z[2] << std::endl;
777 // std::cout << "p = " << p[0] << ", " << p[1] << ", " << p[2] << std::endl;
778 // std::cout << "x = " << x[0] << ", " << x[1] << ", " << x[2] << std::endl;
779 // std::cout << "y = " << y[0] << ", " << y[1] << ", " << y[2] << std::endl;
781 u = (x[XX] * x[XX]) + (x[YY] * x[YY]) + (x[ZZ] * x[ZZ]);
782 v = (y[XX] * y[XX]) + (y[YY] * y[YY]) + (y[ZZ] * y[ZZ]);
784 // std::cout << "u = " << u << std::endl;
785 // std::cout << "v = " << v << std::endl;
787 if (u > 0 and v > 0){
788 u = ((p[XX] * x[XX]) + (p[YY] * x[YY]) + (p[ZZ] * x[ZZ])) / std::sqrt(u);
789 v = ((p[XX] * y[XX]) + (p[YY] * y[YY]) + (p[ZZ] * y[ZZ])) / std::sqrt(v);
790 // std::cout << "new u = " << u << std::endl;
791 // std::cout << "new v = " << v << std::endl;
792 if (u != 0 or v != 0){
793 result = std::atan2(v, u) * gmx::c_rad2Deg;
794 // std::cout << "result = " << result << std::endl;
800 void DsspTool::calculateDihedrals(const t_trxframe &fr, const t_pbc *pbc){
801 const float epsilon = 29;
802 const float phi_min = -75 - epsilon; // -104
803 const float phi_max = -75 + epsilon; // -46
804 const float psi_min = 145 - epsilon; // 116
805 const float psi_max = 145 + epsilon; // 176
806 std::vector<float> phi(IndexMap.size(), 360), psi(IndexMap.size(), 360);
809 // phi.resize(IndexMap.size(), 360);
810 // psi.resize(IndexMap.size(), 360);
812 for (std::size_t i = 1; i + 1 < IndexMap.size(); ++i){ // TODO add index verifictaion (check if those atom indexes exist)
813 // std::cout << "For resi " << i << " :" << std::endl;
814 phi[i] = CalculateDihedralAngle(static_cast<int>(IndexMap[i - 1].getIndex(backboneAtomTypes::AtomC)),
815 static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomN)),
816 static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomCA)),
817 static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
820 psi[i] = CalculateDihedralAngle(static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomN)),
821 static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomCA)),
822 static_cast<int>(IndexMap[i].getIndex(backboneAtomTypes::AtomC)),
823 static_cast<int>(IndexMap[i + 1].getIndex(backboneAtomTypes::AtomN)),
826 // std::cout << "For " << i << " phi = " << phi[i] << ", psi = " << psi[i] << std::endl;
827 // std::cout << "phi[" << i << "] = " << phi[i] << std::endl;
828 // std::cout << "psi[" << i << "] = " << psi[i] << std::endl;
831 for (std::size_t i = 1; i + 3 < IndexMap.size(); ++i){
832 switch (initParams.pp_stretch){
834 if (phi_min > phi[i] or phi[i] > phi_max or
835 phi_min > phi[i + 1] or phi[i + 1]> phi_max){
839 if (psi_min > psi[i] or psi[i] > psi_max or
840 psi_min > psi[i + 1] or psi[i + 1] > psi_max){
844 switch (PatternSearch.SecondaryStructuresStatusMap[i].getStatus(turnsTypes::Turn_PP)){
845 case HelixPositions::None:
846 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start, turnsTypes::Turn_PP);
849 case HelixPositions::End:
850 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start_AND_End, turnsTypes::Turn_PP);
857 PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(HelixPositions::End, turnsTypes::Turn_PP);
858 /* Пропустил проверку того, что заменяемая ак - петля */
859 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Helix_PP);
860 PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(secondaryStructureTypes::Helix_PP);
864 if (phi_min > phi[i] or phi[i] > phi_max or
865 phi_min > phi[i + 1] or phi[i + 1]> phi_max or
866 phi_min > phi[i + 2] or phi[i + 2]> phi_max){
870 if (psi_min > psi[i] or psi[i] > psi_max or
871 psi_min > psi[i + 1] or psi[i + 1] > psi_max or
872 psi_min > psi[i + 2] or psi[i + 2] > psi_max){
876 switch (PatternSearch.SecondaryStructuresStatusMap[i].getStatus(turnsTypes::Turn_PP)){
877 case HelixPositions::None:
878 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start, turnsTypes::Turn_PP);
881 case HelixPositions::End:
882 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(HelixPositions::Start_AND_End, turnsTypes::Turn_PP);
889 PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(HelixPositions::Middle, turnsTypes::Turn_PP);
890 PatternSearch.SecondaryStructuresStatusMap[i + 2].setStatus(HelixPositions::End, turnsTypes::Turn_PP);
891 /* Пропустил проверку того, что заменяемая ак - петля */
892 PatternSearch.SecondaryStructuresStatusMap[i].setStatus(secondaryStructureTypes::Helix_PP);
893 PatternSearch.SecondaryStructuresStatusMap[i + 1].setStatus(secondaryStructureTypes::Helix_PP);
894 PatternSearch.SecondaryStructuresStatusMap[i + 2].setStatus(secondaryStructureTypes::Helix_PP);
899 throw std::runtime_error("Unsupported stretch length");
905 void DsspTool::calculateHBondEnergy(ResInfo& Donor,
907 const t_trxframe& fr,
911 * DSSP uses eq from dssp 2.x
912 * kCouplingConstant = 27.888, // = 332 * 0.42 * 0.2
913 * E = k * (1/rON + 1/rCH - 1/rOH - 1/rCN) where CO comes from one AA and NH from another
917 * For the note, H-Bond Donor is N-H («Donor of H») and H-Bond Acceptor is C=O («Acceptor of H»)
921 if (CalculateAtomicDistances(
922 Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc)
923 >= minimalCAdistance)
928 const float kCouplingConstant = 27.888;
929 const float minimalAtomDistance{ 0.5 },
931 float HbondEnergy{ 0 };
932 float distanceNO{ 0 }, distanceHC{ 0 }, distanceHO{ 0 }, distanceNC{ 0 };
934 // std::cout << "For Donor №" << Donor.info->nr - 1 << " and Accpetor №" << Acceptor.info->nr - 1 << std::endl;
936 if( !(Donor.is_proline) && (Acceptor.getIndex(backboneAtomTypes::AtomC) && Acceptor.getIndex(backboneAtomTypes::AtomO)
937 && Donor.getIndex(backboneAtomTypes::AtomN) && ( Donor.getIndex(backboneAtomTypes::AtomH) || initParams.addHydrogens ) ) ){ // TODO
938 distanceNO = CalculateAtomicDistances(
939 Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
940 distanceNC = CalculateAtomicDistances(
941 Donor.getIndex(backboneAtomTypes::AtomN), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
942 if (initParams.addHydrogens){
943 if (Donor.prevResi != nullptr && Donor.prevResi->getIndex(backboneAtomTypes::AtomC) && Donor.prevResi->getIndex(backboneAtomTypes::AtomO)){
945 float prevCODist {CalculateAtomicDistances(Donor.prevResi->getIndex(backboneAtomTypes::AtomC), Donor.prevResi->getIndex(backboneAtomTypes::AtomO), fr, pbc)};
946 for (int i{XX}; i <= ZZ; ++i){
947 float prevCO = fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomC)][i] - fr.x[Donor.prevResi->getIndex(backboneAtomTypes::AtomO)][i];
948 atomH[i] = fr.x[Donor.getIndex(backboneAtomTypes::AtomH)][i]; // Но на самом деле берутся координаты N
949 atomH[i] += prevCO / prevCODist;
951 distanceHO = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
952 distanceHC = CalculateAtomicDistances(atomH, Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
955 distanceHO = distanceNO;
956 distanceHC = distanceNC;
960 distanceHO = CalculateAtomicDistances(
961 Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomO), fr, pbc);
962 distanceHC = CalculateAtomicDistances(
963 Donor.getIndex(backboneAtomTypes::AtomH), Acceptor.getIndex(backboneAtomTypes::AtomC), fr, pbc);
965 if ((distanceNO < minimalAtomDistance) || (distanceHC < minimalAtomDistance)
966 || (distanceHO < minimalAtomDistance) || (distanceNC < minimalAtomDistance))
968 HbondEnergy = minEnergy;
973 * ((1 / distanceNO) + (1 / distanceHC) - (1 / distanceHO) - (1 / distanceNC));
976 // std::cout << "CA-CA distance: " << CalculateAtomicDistances(
977 // Donor.getIndex(backboneAtomTypes::AtomCA), Acceptor.getIndex(backboneAtomTypes::AtomCA), fr, pbc) << std::endl;
978 // std::cout << "N-O distance: " << distanceNO << std::endl;
979 // std::cout << "N-C distance: " << distanceNC << std::endl;
980 // std::cout << "H-O distance: " << distanceHO << std::endl;
981 // std::cout << "H-C distance: " << distanceHC << std::endl;
983 HbondEnergy = std::round(HbondEnergy * 1000) / 1000;
985 // if ( HbondEnergy < minEnergy ){ // I don't think that this is correct
986 // HbondEnergy = minEnergy;
989 // std::cout << "Calculated energy = " << HbondEnergy << std::endl;
992 // std::cout << "Donor Is Proline" << std::endl;
995 if (HbondEnergy < Donor.acceptorEnergy[0]){
996 Donor.acceptor[1] = Donor.acceptor[0];
997 Donor.acceptor[0] = Acceptor.info;
998 Donor.acceptorEnergy[0] = HbondEnergy;
1000 else if (HbondEnergy < Donor.acceptorEnergy[1]){
1001 Donor.acceptor[1] = Acceptor.info;
1002 Donor.acceptorEnergy[1] = HbondEnergy;
1005 if (HbondEnergy < Acceptor.donorEnergy[0]){
1006 Acceptor.donor[1] = Acceptor.donor[0];
1007 Acceptor.donor[0] = Donor.info;
1008 Acceptor.donorEnergy[0] = HbondEnergy;
1010 else if (HbondEnergy < Acceptor.donorEnergy[1]){
1011 Acceptor.donor[1] = Donor.info;
1012 Acceptor.donorEnergy[1] = HbondEnergy;
1017 /* Calculate Distance From B to A */
1018 float DsspTool::CalculateAtomicDistances(const int &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
1020 gmx::RVec r{ 0, 0, 0 };
1021 pbc_dx(pbc, fr.x[A], fr.x[B], r.as_vec());
1022 return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
1025 /* Calculate Distance From B to A, where A is only fake coordinates */
1026 float DsspTool::CalculateAtomicDistances(const rvec &A, const int &B, const t_trxframe &fr, const t_pbc *pbc)
1028 gmx::RVec r{ 0, 0, 0 };
1029 pbc_dx(pbc, A, fr.x[B], r.as_vec());
1030 return r.norm() * gmx::c_nm2A; // НЕ ТРОГАТЬ
1033 void DsspTool::initAnalysis(/*const TrajectoryAnalysisSettings &settings,*/const TopologyInformation& top, const initParameters &initParamz)
1035 initParams = initParamz;
1036 ResInfo _backboneAtoms;
1038 std::string proLINE;
1039 int resicompare{ top.atoms()->atom[static_cast<std::size_t>(*(initParams.sel_.atomIndices().begin()))].resind };
1041 IndexMap.push_back(_backboneAtoms);
1042 IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
1043 proLINE = *(IndexMap[i].info->name);
1044 if( proLINE.compare("PRO") == 0 ){
1045 IndexMap[i].is_proline = true;
1048 for (gmx::ArrayRef<const int>::iterator ai{ initParams.sel_.atomIndices().begin() }; (ai != initParams.sel_.atomIndices().end()); ++ai){
1049 if (resicompare != top.atoms()->atom[static_cast<std::size_t>(*ai)].resind)
1052 resicompare = top.atoms()->atom[static_cast<std::size_t>(*ai)].resind;
1053 IndexMap.emplace_back(_backboneAtoms);
1054 IndexMap[i].info = &(top.atoms()->resinfo[resicompare]);
1055 proLINE = *(IndexMap[i].info->name);
1056 if( proLINE.compare("PRO") == 0 ){
1057 IndexMap[i].is_proline = true;
1061 std::string atomname(*(top.atoms()->atomname[static_cast<std::size_t>(*ai)]));
1062 if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA])
1064 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomCA)] = *ai;
1066 else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC])
1068 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomC)] = *ai;
1070 else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO])
1072 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomO)] = *ai;
1074 else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN])
1076 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomN)] = *ai;
1077 if (initParamz.addHydrogens == true){
1078 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
1081 else if (atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH] && initParamz.addHydrogens == false) // Юзать водород в структуре
1083 IndexMap[i]._backboneIndices[static_cast<std::size_t>(backboneAtomTypes::AtomH)] = *ai;
1088 // if( atomname == backboneAtomTypeNames[backboneAtomTypes::AtomCA] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomC] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomO]
1089 // || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomN] || atomname == backboneAtomTypeNames[backboneAtomTypes::AtomH]){
1090 // std::cout << "Atom " << atomname << " №" << *ai << " From Resi " << *(top.atoms()->resinfo[i].name) << " №" << resicompare << std::endl;
1094 for (std::size_t j {1}; j < IndexMap.size(); ++j){
1095 IndexMap[j].prevResi = &(IndexMap[j - 1]);
1097 IndexMap[j - 1].nextResi = &(IndexMap[j]);
1099 // std::cout << "Resi " << IndexMap[i].info->nr << *(IndexMap[i].info->name) << std::endl;
1100 // std::cout << "Prev resi is " << IndexMap[i].prevResi->info->nr << *(IndexMap[i].prevResi->info->name) << std::endl;
1101 // std::cout << "Prev resi's next resi is " << IndexMap[i - 1].nextResi->info->nr << *(IndexMap[i - 1].nextResi->info->name) << std::endl;
1102 // std::cout << IndexMap[j].prevResi->info->nr;
1103 // std::cout << *(IndexMap[j].prevResi->info->name) ;
1104 // std::cout << " have CA = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomCA) ;
1105 // std::cout << " C = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomC);
1106 // std::cout << " O = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomO);
1107 // std::cout << " N = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomN);
1108 // std::cout << " H = " << IndexMap[j].prevResi->getIndex(backboneAtomTypes::AtomH) << std::endl;
1114 void DsspTool::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc)
1117 switch(initParams.NBS){
1118 case (NBSearchMethod::Classique): {
1120 // store positions of CA atoms to use them for nbSearch
1121 std::vector<gmx::RVec> positionsCA_;
1122 for (std::size_t i{ 0 }; i < IndexMap.size(); ++i)
1124 positionsCA_.emplace_back(fr.x[IndexMap[i].getIndex(backboneAtomTypes::AtomCA)]);
1127 AnalysisNeighborhood nb_;
1128 nb_.setCutoff(initParams.cutoff_);
1129 AnalysisNeighborhoodPositions nbPos_(positionsCA_);
1130 gmx::AnalysisNeighborhoodSearch start = nb_.initSearch(pbc, nbPos_);
1131 gmx::AnalysisNeighborhoodPairSearch pairSearch = start.startPairSearch(nbPos_);
1132 gmx::AnalysisNeighborhoodPair pair;
1133 while (pairSearch.findNextPair(&pair))
1135 if(CalculateAtomicDistances(
1136 IndexMap[pair.refIndex()].getIndex(backboneAtomTypes::AtomCA), IndexMap[pair.testIndex()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1137 < minimalCAdistance){
1138 calculateHBondEnergy(IndexMap[pair.refIndex()], IndexMap[pair.testIndex()], fr, pbc);
1139 if (IndexMap[pair.testIndex()].info != IndexMap[pair.refIndex() + 1].info){
1140 calculateHBondEnergy(IndexMap[pair.testIndex()], IndexMap[pair.refIndex()], fr, pbc);
1147 case (NBSearchMethod::Experimental): { // TODO FIX
1149 alternateNeighborhoodSearch as_;
1151 as_.setCutoff(initParams.cutoff_);
1153 as_.AltPairSearch(fr, IndexMap);
1155 while (as_.findNextPair()){
1156 if(CalculateAtomicDistances(
1157 IndexMap[as_.getResiI()].getIndex(backboneAtomTypes::AtomCA), IndexMap[as_.getResiJ()].getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1158 < minimalCAdistance){
1159 calculateHBondEnergy(IndexMap[as_.getResiI()], IndexMap[as_.getResiJ()], fr, pbc);
1160 if (IndexMap[as_.getResiJ()].info != IndexMap[as_.getResiI() + 1].info){
1161 calculateHBondEnergy(IndexMap[as_.getResiJ()], IndexMap[as_.getResiI()], fr, pbc);
1170 for(std::vector<ResInfo>::iterator Donor {IndexMap.begin()}; Donor != IndexMap.end() ; ++Donor){
1171 for(std::vector<ResInfo>::iterator Acceptor {Donor + 1} ; Acceptor != IndexMap.end() ; ++Acceptor){
1172 if(CalculateAtomicDistances(
1173 Donor->getIndex(backboneAtomTypes::AtomCA), Acceptor->getIndex(backboneAtomTypes::AtomCA), fr, pbc)
1174 < minimalCAdistance){
1175 calculateHBondEnergy(*Donor, *Acceptor, fr, pbc);
1176 if (Acceptor != Donor + 1){
1177 calculateHBondEnergy(*Acceptor, *Donor, fr, pbc);
1187 // for(std::size_t i {0}; i < IndexMap.size(); ++i){
1188 // std::cout << IndexMap[i].info->nr << " " << *(IndexMap[i].info->name) << std::endl;
1191 PatternSearch.initiateSearch(IndexMap, initParams.PPHelices, initParams.pp_stretch);
1192 calculateBends(fr, pbc);
1193 calculateDihedrals(fr, pbc);
1194 Storage.storageData(frnr, PatternSearch.patternSearch());
1198 std::vector<std::pair<int, std::string>> DsspTool::getData(){
1199 return Storage.returnData();
1202 } // namespace analysismodules