4d807ad681520b49354bd93d42673c0bc605f34e
[alexxy/gromacs-colorvec.git] / src / colorvec.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015, 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::Freevolume.
38  *
39  * \author Titov Anatoly <Wapuk-cobaka@yandex.ru>
40  * \ingroup module_trajectoryanalysis
41  */
42
43 #include <gromacs/trajectoryanalysis.h>
44 #include <gromacs/trajectoryanalysis/topologyinformation.h>
45 #include <gromacs/selection/nbsearch.h>
46
47 #include <iostream>
48 #include <fstream>
49 //#include <chrono>
50 #include <cfloat>
51 #include <omp.h>
52 #include <vector>
53 #include <cstdio>
54 #include <iomanip>
55
56 // структура углов для одной краски
57 struct colorLocalAngles {
58     gmx::RVec   a1{0., 0., 0.}, a2{0., 0., 0.};
59     gmx::RVec   b1{0., 0., 0.}, b2{0., 0., 0.};
60     gmx::RVec   n1{0., 0., 0.}, n2{0., 0., 0.};
61     double      a12 {0}, b12 {0}, n12 {0};
62     std::vector< double > betaAngles;
63 };
64
65 // функция парсинга одной строки
66 void parseBetaListDATLine(const std::string &currentLine, std::vector< std::vector< unsigned int > > &localInputBL) {
67     size_t                      equalCount = 0;
68     std::vector< unsigned int > a;
69     a.resize(0);
70     for (size_t i = 0; i < currentLine.size(); ++i) {
71         if (currentLine[i] == '=') { // подсчитываем число "пустых" символов
72             ++equalCount;
73         } else {
74             size_t temp {i};
75             for (size_t j = i + 1; j < currentLine.size(); ++j) {
76                 if (currentLine[j] == 'B' || currentLine[j] == 'E') {
77                     ++temp;
78                 } else {
79                     break;
80                 }
81             }
82             if (temp - i > 3) {
83                 localInputBL.push_back(a);
84                 for (size_t j = i; j <= temp; ++j) {
85                     localInputBL.back().push_back(j - equalCount);
86                 }
87                 i = temp;
88             }
89         }
90     }
91 }
92
93 // функция нахождения бета-листов в структуре по файлу ДССП
94 void betaListDigestion(const std::string &inputFile, std::vector< std::vector< std::vector< unsigned int > > > &inputBL) {
95     inputBL.resize(0);
96     std::ifstream   file(inputFile);
97     std::string     line;
98     getline(file, line); // считываем число в первой строке - кол-во осмысленных элементов в строках - нам не нужно
99     getline(file, line);
100     if (line.size() > 3) {
101         do {
102             inputBL.resize(inputBL.size() + 1);
103             parseBetaListDATLine(line, inputBL.back());
104             line = "";
105             getline(file, line);
106         } while (line.size() > 3);
107     } else {
108         throw "DSSP DAT FILE IS EMPTY";
109     }
110     file.close();
111 }
112
113 // функция выделения индексов для бэталистов
114 inline void aminoacidsIndexation(const std::vector< size_t > &inputIndex, const gmx::TopologyInformation &top, std::vector< std::vector< size_t > > &aminoacidsIndex) {
115     aminoacidsIndex.resize(0);
116     for (size_t i = 0; i < inputIndex.size(); ++i) {
117         aminoacidsIndex.resize(std::max(aminoacidsIndex.size(), static_cast< size_t >(top.atoms()->atom[inputIndex[i]].resind + 1)));
118         aminoacidsIndex[top.atoms()->atom[inputIndex[i]].resind].push_back(inputIndex[i]);
119     }
120 }
121
122 // функция поиска RVec в кадре по имени->индексу
123 gmx::RVec returnRVec(const std::vector< gmx::RVec > &frame, const std::vector< std::pair< std::string, size_t > > &colorsIndex, const std::string &toFind) {
124     for (auto &i : colorsIndex) {
125         if (i.first == toFind) {
126             return frame[i.second];
127         }
128     }
129     throw "WRONG COLOR ATOM NAME TO EVALUATE VECTORS";
130 }
131
132 // функция векторного произведения двух RVec'ов
133 inline void RVecVecMultiply(gmx::RVec &n, const gmx::RVec &a, const gmx::RVec &b) {
134     n[0] = a[1] * b[2] - a[2] * b[1];
135     n[1] = a[2] * b[0] - a[0] * b[2];
136     n[2] = a[0] * b[1] - a[2] * b[0];
137 }
138
139 // поиск угла между двумя RVec'ами
140 inline double RVecAngle(const gmx::RVec &a, const gmx::RVec &b) {
141     return std::acos((a[0] * b[0] + a[1] * b[1] + a[2] * b[2]) / (a.norm() * b.norm())) * 180.0 / 3.14159265;
142 }
143
144 // вычисление внутренних углов в краске
145 void colorsAnglesEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< std::pair< std::string, size_t > > > &colorsIndex,
146                             std::vector< colorLocalAngles > &colorStruct) {
147     colorStruct.resize(colorsIndex.size());
148     #pragma omp parallel for ordered schedule(dynamic)
149     for (size_t i = 0; i < colorsIndex.size(); ++i) {
150         colorStruct[i].betaAngles.resize(0);
151         // "левый" блок
152         colorStruct[i].a1 = returnRVec(frame, colorsIndex[i], "CAF") - returnRVec(frame, colorsIndex[i], "CAJ") + returnRVec(frame, colorsIndex[i], "CAK") - returnRVec(frame, colorsIndex[i], "CAO");
153         colorStruct[i].a1 /= colorStruct[i].a1.norm();
154         colorStruct[i].b1 = returnRVec(frame, colorsIndex[i], "CAO") - returnRVec(frame, colorsIndex[i], "CAJ") + returnRVec(frame, colorsIndex[i], "CAM") - returnRVec(frame, colorsIndex[i], "CAH") + returnRVec(frame, colorsIndex[i], "CAK") - returnRVec(frame, colorsIndex[i], "CAF");
155         colorStruct[i].b1 /= colorStruct[i].b1.norm();
156         RVecVecMultiply(colorStruct[i].n1, colorStruct[i].a1, colorStruct[i].b1);
157         colorStruct[i].n1 /= colorStruct[i].n1.norm();
158         // "правый" блок
159         colorStruct[i].a2 = returnRVec(frame, colorsIndex[i], "CAB") - returnRVec(frame, colorsIndex[i], "CBQ") + returnRVec(frame, colorsIndex[i], "CBP") - returnRVec(frame, colorsIndex[i], "CBL");
160         colorStruct[i].a2 /= colorStruct[i].a2.norm();
161         colorStruct[i].b2 = returnRVec(frame, colorsIndex[i], "CBL") - returnRVec(frame, colorsIndex[i], "CBQ") + returnRVec(frame, colorsIndex[i], "CBN") - returnRVec(frame, colorsIndex[i], "CBS") + returnRVec(frame, colorsIndex[i], "CBP") - returnRVec(frame, colorsIndex[i], "CAB");
162         colorStruct[i].b2 /= colorStruct[i].b2.norm();
163         RVecVecMultiply(colorStruct[i].n2, colorStruct[i].a2, colorStruct[i].b2);
164         colorStruct[i].n2 *= -1; // для "сонаправленности" векторов нормали
165         colorStruct[i].n2 /= colorStruct[i].n2.norm();
166         colorStruct[i].a12 = RVecAngle(colorStruct[i].a1, colorStruct[i].a2);
167         colorStruct[i].b12 = RVecAngle(colorStruct[i].b1, colorStruct[i].b2);
168         colorStruct[i].n12 = RVecAngle(colorStruct[i].n1, colorStruct[i].n2);
169     }
170     #pragma omp barrier
171 }
172
173 // вычисление направляющего вектора в бэта-листе
174 inline void betaListsRVecsEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< unsigned int > > &inputBetaLists,
175                                      std::vector< gmx::RVec > &temp, const std::vector< size_t > &inputCA) {
176     temp.resize(0);
177     gmx::RVec tempA;
178     for (const auto &i : inputBetaLists) {
179         tempA = frame[inputCA[i[i.size() - 1]]] + frame[inputCA[i[i.size() - 2]]] - frame[inputCA[i[1]]] - frame[inputCA[i[0]]];
180         tempA /= tempA.norm();
181         temp.push_back(tempA);
182     }
183 }
184
185 //// определение близко ли к белку находится краска
186 //bool isNearPeptide(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > &inputFrame,
187 //                          /*gmx::AnalysisNeighborhood &nbhood, */const std::vector< size_t > &inputIndex,
188 //                          const std::vector< std::pair< std::string, size_t > > &inputColor) {
189 //    /*gmx::AnalysisNeighborhoodSearch     nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
190 //    gmx::AnalysisNeighborhoodPair       pair;*/
191
192 //    for (size_t i = 0; i < inputColor.size(); ++i) {
193 //        std::cout << inputColor[i].first << " " << i << " / "  << inputColor.size() << std::endl;
194 //        gmx::AnalysisNeighborhood           nbhood;
195 //        nbhood.setCutoff(0.8);
196 //        gmx::AnalysisNeighborhoodSearch     nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
197 //        gmx::AnalysisNeighborhoodPair       pair;
198 //        gmx::AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(inputFrame[inputColor[i].second].as_vec());
199 //        int count1 {0};
200 //        while (pairSearch.findNextPair(&pair)) {
201 //            std::cout << ++count1 << " ";
202 //            for (size_t j = 0; j < inputIndex.size(); ++j) {
203 //                if (pair.refIndex() == inputIndex[j]) {
204 //                    return true;
205 //                }
206 //            }
207 //        }
208 //        std::cout << inputColor[i].first << " atom done" << std::endl;
209 //    }
210 //    return false;
211 //}
212
213 // определение близко ли к белку находится краска
214 bool isNearPeptide( const std::vector< gmx::RVec > &inputFrame, const std::vector< size_t > &inputIndex,
215                     const std::vector< std::pair< std::string, size_t > > &inputColor, const double cutOff) {
216     for (size_t i {0}; i < inputColor.size(); ++i) {
217         for (size_t j {0}; j < inputIndex.size(); ++j) {
218             if ((inputFrame[inputIndex[j]] - inputFrame[inputColor[i].second]).norm() <= cutOff) {
219                 return true;
220             }
221         }
222     }
223     return false;
224 }
225
226 // поиск ближайших к краске бэта-листов
227 //inline void searchNearBetaLists(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > &inputFrame,
228 //                                gmx::AnalysisNeighborhood &nbhood,
229 //                                const std::vector< std::vector< unsigned int > > &inputBLists,
230 //                                const std::vector< std::pair< std::string, size_t > > &inputColor,
231 //                                std::vector< bool > &outputList, const std::vector< std::vector< size_t > > &inputAminoacids) {
232 //    outputList.resize(0);
233 //    outputList.resize(inputBLists.size(), false);
234 //    gmx::AnalysisNeighborhoodSearch      nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
235 //    gmx::AnalysisNeighborhoodPair        pair;
236 //    for (const auto &i : inputColor) {
237 //        while (nbsearch.startPairSearch(inputFrame[i.second].as_vec()).findNextPair(&pair)) {
238 //            for (size_t j = 0; j < inputBLists.size(); ++j) {
239 //                for (size_t k = 0; k < inputBLists[j].size(); ++k) {
240 //                    for (size_t m = 0; m < inputAminoacids[inputBLists[j][k]].size(); ++m) {
241 //                        if (pair.testIndex() == inputAminoacids[inputBLists[j][k]][m]) {
242 //                            outputList[j] = true;
243 //                        }
244 //                    }
245 //                }
246 //            }
247 //        }
248 //    }
249 //}
250
251 inline void searchNearBetaLists(const std::vector< gmx::RVec > &inputFrame, const std::vector< std::vector< unsigned int > > &inputBLists,
252                                 const std::vector< std::pair< std::string, size_t > > &inputColor, std::vector< bool > &outputList,
253                                 const std::vector< std::vector< size_t > > &inputAminoacids, const double cutOff) {
254     outputList.resize(0);
255     outputList.resize(inputBLists.size(), false);
256     for (size_t i {0}; i < inputColor.size(); ++i) {
257         for (size_t j = 0; j < inputBLists.size(); ++j) {
258             for (size_t k = 0; k < inputBLists[j].size(); ++k) {
259                 for (size_t m = 0; m < inputAminoacids[inputBLists[j][k]].size(); ++m) {
260                     if ((inputFrame[inputAminoacids[inputBLists[j][k]][m]] - inputFrame[inputColor[i].second]).norm() <= cutOff) {
261                         outputList[j] = true;
262                         break;
263                     }
264                 }
265             }
266         }
267     }
268 }
269
270 // определение углов между краской и направляющим вектором бэта-листа
271 inline void computeAnglesColorsVsBeta(const gmx::RVec &inputBetaRVec, colorLocalAngles &colorFormation) {
272     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.a1));
273     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.b1));
274     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.n1));
275     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.a2));
276     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.b2));
277     colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.n2));
278 }
279
280 // функция записи в файл значений углов для кадра
281 void anglesFileDump(const int frameNum, const std::string &output, const std::vector< bool > &toPeptide, const std::vector< colorLocalAngles > &colorFormation) {
282     std::ofstream   file(output, std::ofstream::app);
283     int temp = 0;
284     std::vector< double > betaTemp;
285     betaTemp.resize(0);
286     file << "frame =" << std::setw(8) << frameNum << std::endl;
287     for (size_t i = 0; i < colorFormation.size(); ++i) {
288         file << "color #" << std::setw(3) << i;
289         if (toPeptide[i]) {
290             file << std::setw(4) << "yes";
291         } else {
292             file << std::setw(4) << "no";
293         }
294         file << std::setw(8) << std::setprecision(3) << colorFormation[i].a12;
295         file << std::setw(8) << std::setprecision(3) << colorFormation[i].b12;
296         file << std::setw(8) << std::setprecision(3) << colorFormation[i].n12;
297         temp = colorFormation[i].betaAngles.size() / 6;
298         if (temp == 0) {
299             file << std::setw(4) << 0;
300         } else {
301             betaTemp.resize(0);
302             betaTemp.resize(6, 0.); // magic number, meh
303             for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
304                 betaTemp[j % 6] += colorFormation[i].betaAngles[j];
305             }
306             for (size_t j = 0; j < betaTemp.size(); ++j) {
307                 betaTemp[j] /= temp;
308             }
309             file << std::setw(4) << temp;
310             for (size_t j = 0; j < betaTemp.size(); ++j) {
311                 file << std::setw(8) << std::setprecision(3) << betaTemp[j];
312             }
313             for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
314                 file << std::setw(8) << std::setprecision(3) << colorFormation[i].betaAngles[j];
315             }
316         }
317         file << std::endl;
318     }
319     file.close();
320 }
321
322 /*! \brief
323  * \ingroup module_trajectoryanalysis
324  */
325 class colorVec : public gmx::TrajectoryAnalysisModule
326 {
327     public:
328
329         colorVec();
330         virtual ~colorVec();
331
332         //! Set the options and setting
333         virtual void initOptions(gmx::IOptionsContainer          *options,
334                                  gmx::TrajectoryAnalysisSettings *settings);
335
336         //! First routine called by the analysis framework
337         // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
338         virtual void initAnalysis(const gmx::TrajectoryAnalysisSettings &settings,
339                                   const gmx::TopologyInformation        &top);
340
341         //! Call for each frame of the trajectory
342         // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
343         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
344                                   gmx::TrajectoryAnalysisModuleData *pdata);
345
346         //! Last routine called by the analysis framework
347         // virtual void finishAnalysis(t_pbc *pbc);
348         virtual void finishAnalysis(int nframes);
349
350         //! Routine to write output, that is additional over the built-in
351         virtual void writeOutput();
352
353     private:
354
355         gmx::SelectionList                                              sel_;
356         std::string                                                     fnOut           {"name"};   // selectable
357         std::string                                                     fnBetaListsDat;             // selectable
358         float                                                           effRad          {0.8};      // selectable
359         std::vector< size_t >                                           index;
360         std::vector< size_t >                                           indexCA;
361         std::vector< std::vector< size_t > >                            aminoacidsIndex;
362         std::vector< std::vector< std::pair< std::string, size_t > > >  colorsNames;
363         std::vector< std::vector< std::vector< unsigned int > > >       betaLists;
364         gmx::AnalysisNeighborhood                                       nb_;
365
366         // Copy and assign disallowed by base.
367 };
368
369 colorVec::colorVec(): TrajectoryAnalysisModule()
370 {
371 }
372
373 colorVec::~colorVec()
374 {
375 }
376
377 /*
378  * ndx:
379  * [peptide]
380  * [CA]
381  * [color_{i}]
382  *
383  */
384
385 void
386 colorVec::initOptions(  gmx::IOptionsContainer          *options,
387                         gmx::TrajectoryAnalysisSettings *settings)
388 {
389     static const char *const desc[] = {
390         "[THISMODULE] to be done"
391     };
392     // Add the descriptive text (program help text) to the options
393     settings->setHelpText(desc);
394     // Add option for selection list
395     options->addOption(gmx::SelectionOption("sel")
396                             .storeVector(&sel_)
397                             .required().dynamicMask().multiValue()
398                             .description("select pepride and colors / -sf"));
399     // Add option for input file names
400     options->addOption(gmx::StringOption("dat")
401                             .store(&fnBetaListsDat)
402                             .description("a file to make dynamic beta lists"));
403     // Add option for output file name
404     options->addOption(gmx::StringOption("out")
405                             .store(&fnOut)
406                             .description("Index file for the algorithm output."));
407     // Add option for effRad constant
408     options->addOption(gmx::FloatOption("efRad")
409                             .store(&effRad)
410                             .description("max distance from colors to peptide in nm to consider to be \"near\""));
411     // Control input settings
412     settings->setFlags(gmx::TrajectoryAnalysisSettings::efNoUserPBC);
413     settings->setFlag(gmx::TrajectoryAnalysisSettings::efUseTopX);
414     //settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
415     settings->setPBC(true);
416 }
417
418 void
419 colorVec::initAnalysis( const gmx::TrajectoryAnalysisSettings   &settings,
420                         const gmx::TopologyInformation          &top)
421 {
422     // считывание индекса
423     index.resize(0);
424     for (gmx::ArrayRef< const int >::iterator ai {sel_.front().atomIndices().begin()}; (ai < sel_.front().atomIndices().end()); ai++) {
425         index.push_back(static_cast< size_t >(*ai));
426     }
427     // считывание индекса
428     indexCA.resize(0);
429     for (gmx::ArrayRef< const int >::iterator ai {sel_[1].atomIndices().begin()}; (ai < sel_[1].atomIndices().end()); ai++) {
430         indexCA.push_back(static_cast< size_t >(*ai));
431     }
432     // считывание красок
433     colorsNames.resize(sel_.size() - 2);
434     for (size_t i = 2; i < sel_.size(); ++i) {
435         for (gmx::ArrayRef< const int >::iterator ai {sel_[i].atomIndices().begin()}; (ai < sel_[i].atomIndices().end()); ai++) {
436             colorsNames[i - 2].push_back(std::make_pair(*(top.atoms()->atomname[*ai]), static_cast< size_t >(*ai)));
437         }
438     }
439     // разбор dat файла для создания бэта листов
440     betaListDigestion(fnBetaListsDat, betaLists);
441     // разбор топологии для индексации бета листов
442     aminoacidsIndexation(index, top, aminoacidsIndex);
443     // задание радиуса рассматриваемых соседей
444     nb_.setCutoff(effRad);
445     /*
446      * формально можно сделать:
447      * найти соотношение между именами для углов и индексными номерами
448      * заполнить std::vector< std::vector< size_t > > nameToIndex;
449      * и в последствии считать:
450      * a1 = frame[nameToIndex[0][1] + ... - ... - ...
451      * b1 a2 b2 по аналогии только для [@<-[1, 2, 3]]
452      *
453      */
454 }
455
456 void
457 colorVec::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, gmx::TrajectoryAnalysisModuleData *pdata)
458 {
459     std::vector< gmx::RVec > trajectoryFrame;
460     trajectoryFrame.resize(0);
461     // считывания текущего фрейма траектории
462     for (size_t i {0}; i < fr.natoms; ++i) {
463         trajectoryFrame.push_back(fr.x[i]);
464     }
465     // подсчёт углов в красках
466     std::vector< colorLocalAngles > colorFormation;
467     colorsAnglesEvaluation(trajectoryFrame, colorsNames, colorFormation);
468     // рассчёт положения относительно белка
469     /*
470      * формально можно совместить определение близости с белком и определение ближайших бэта-листов
471      * для этого думаю нужно как-то соотнести или сделать соотношение между индексом и бэта-листами (т.е. хранить что бы за О(1) делать)
472      *
473      */
474
475     std::vector< bool >     colorsToPeptide;
476     colorsToPeptide.resize(0);
477     colorsToPeptide.resize(colorsNames.size(), false);
478     for (size_t i = 0; i < colorsNames.size(); ++i) {
479         colorsToPeptide[i] = isNearPeptide(trajectoryFrame, index, colorsNames[i], effRad);
480     }
481     // расчёт угла и среднего угла с ближайшими бета листами
482     std::vector< std::vector< bool > > colorsToBeta;
483     colorsToBeta.resize(0);
484     colorsToBeta.resize(colorsNames.size());
485     std::vector< std::vector< double > > colorsToBetaAngles;
486     colorsToBetaAngles.resize(0);
487     std::vector< gmx::RVec > betaListsRVecs;
488     betaListsRVecsEvaluation(trajectoryFrame, betaLists[frnr], betaListsRVecs, indexCA);
489     for (size_t i = 0; i < colorsToPeptide.size(); ++i) {
490         if (colorsToPeptide[i]) {
491             searchNearBetaLists(trajectoryFrame, betaLists[frnr], colorsNames[i], colorsToBeta[i], aminoacidsIndex, effRad);
492         }
493     }
494     for (size_t i = 0; i < colorsToBeta.size(); ++i) {
495         for (size_t j = 0; j < colorsToBeta[i].size(); ++j) {
496             if (colorsToBeta[i][j]) {
497                 computeAnglesColorsVsBeta(betaListsRVecs[j], colorFormation[i]);
498             }
499         }
500     }
501     // вывод в файл(ы) информацию для фрейма
502     anglesFileDump(frnr, fnOut, colorsToPeptide, colorFormation);
503     std::cout << std::endl << "frame number " << frnr;
504 }
505
506 void
507 colorVec::finishAnalysis(int nframes)
508 {
509
510 }
511
512 void
513 colorVec::writeOutput()
514 {
515     /*
516      *
517      * frame #<current frame number>
518      * color #<color number> <yes/no for "close" to peptide> a12 b12 n12 <number of beta lists> <6 avg beta angles> [<6 angles for each beta list>]
519      *
520      */
521 }
522
523 /*! \brief
524  * The main function for the analysis template.
525  */
526 int
527 main(int argc, char *argv[])
528 {
529     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<colorVec>(argc, argv);
530 }