94fad4042960fffe1b6214c67c7457e9ada6fc02
[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(6, 0.); // magic number, meh
302             for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
303                 betaTemp[j % 6] += colorFormation[i].betaAngles[j];
304             }
305             for (size_t j = 0; j < betaTemp.size(); ++j) {
306                 betaTemp[j] /= temp;
307             }
308             file << std::setw(4) << temp;
309             for (size_t j = 0; j < betaTemp.size(); ++j) {
310                 file << std::setw(8) << std::setprecision(3) << betaTemp[j];
311             }
312             for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
313                 file << std::setw(8) << std::setprecision(3) << colorFormation[i].betaAngles[j];
314             }
315         }
316         file << std::endl;
317     }
318     file.close();
319 }
320
321 /*! \brief
322  * \ingroup module_trajectoryanalysis
323  */
324 class colorVec : public gmx::TrajectoryAnalysisModule
325 {
326     public:
327
328         colorVec();
329         virtual ~colorVec();
330
331         //! Set the options and setting
332         virtual void initOptions(gmx::IOptionsContainer          *options,
333                                  gmx::TrajectoryAnalysisSettings *settings);
334
335         //! First routine called by the analysis framework
336         // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
337         virtual void initAnalysis(const gmx::TrajectoryAnalysisSettings &settings,
338                                   const gmx::TopologyInformation        &top);
339
340         //! Call for each frame of the trajectory
341         // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
342         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
343                                   gmx::TrajectoryAnalysisModuleData *pdata);
344
345         //! Last routine called by the analysis framework
346         // virtual void finishAnalysis(t_pbc *pbc);
347         virtual void finishAnalysis(int nframes);
348
349         //! Routine to write output, that is additional over the built-in
350         virtual void writeOutput();
351
352     private:
353
354         gmx::SelectionList                                              sel_;
355         std::string                                                     fnOut           {"name"};   // selectable
356         std::string                                                     fnBetaListsDat;             // selectable
357         float                                                           effRad          {0.8};      // selectable
358         std::vector< size_t >                                           index;
359         std::vector< size_t >                                           indexCA;
360         std::vector< std::vector< size_t > >                            aminoacidsIndex;
361         std::vector< std::vector< std::pair< std::string, size_t > > >  colorsNames;
362         std::vector< std::vector< std::vector< unsigned int > > >       betaLists;
363         gmx::AnalysisNeighborhood                                       nb_;
364
365         // Copy and assign disallowed by base.
366 };
367
368 colorVec::colorVec(): TrajectoryAnalysisModule()
369 {
370 }
371
372 colorVec::~colorVec()
373 {
374 }
375
376 /*
377  * ndx:
378  * [peptide]
379  * [CA]
380  * [color_{i}]
381  *
382  */
383
384 void
385 colorVec::initOptions(  gmx::IOptionsContainer          *options,
386                         gmx::TrajectoryAnalysisSettings *settings)
387 {
388     static const char *const desc[] = {
389         "[THISMODULE] to be done"
390     };
391     // Add the descriptive text (program help text) to the options
392     settings->setHelpText(desc);
393     // Add option for selection list
394     options->addOption(gmx::SelectionOption("sel")
395                             .storeVector(&sel_)
396                             .required().dynamicMask().multiValue()
397                             .description("select pepride and colors / -sf"));
398     // Add option for input file names
399     options->addOption(gmx::StringOption("dat")
400                             .store(&fnBetaListsDat)
401                             .description("a file to make dynamic beta lists"));
402     // Add option for output file name
403     options->addOption(gmx::StringOption("out")
404                             .store(&fnOut)
405                             .description("Index file for the algorithm output."));
406     // Add option for effRad constant
407     options->addOption(gmx::FloatOption("efRad")
408                             .store(&effRad)
409                             .description("max distance from colors to peptide in nm to consider to be \"near\""));
410     // Control input settings
411     settings->setFlags(gmx::TrajectoryAnalysisSettings::efNoUserPBC);
412     settings->setFlag(gmx::TrajectoryAnalysisSettings::efUseTopX);
413     //settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
414     settings->setPBC(true);
415 }
416
417 void
418 colorVec::initAnalysis( const gmx::TrajectoryAnalysisSettings   &settings,
419                         const gmx::TopologyInformation          &top)
420 {
421     // считывание индекса
422     index.resize(0);
423     for (gmx::ArrayRef< const int >::iterator ai {sel_.front().atomIndices().begin()}; (ai < sel_.front().atomIndices().end()); ai++) {
424         index.push_back(static_cast< size_t >(*ai));
425     }
426     // считывание индекса
427     indexCA.resize(0);
428     for (gmx::ArrayRef< const int >::iterator ai {sel_[1].atomIndices().begin()}; (ai < sel_[1].atomIndices().end()); ai++) {
429         indexCA.push_back(static_cast< size_t >(*ai));
430     }
431     // считывание красок
432     colorsNames.resize(sel_.size() - 2);
433     for (size_t i = 2; i < sel_.size(); ++i) {
434         for (gmx::ArrayRef< const int >::iterator ai {sel_[i].atomIndices().begin()}; (ai < sel_[i].atomIndices().end()); ai++) {
435             colorsNames[i - 2].push_back(std::make_pair(*(top.atoms()->atomname[*ai]), static_cast< size_t >(*ai)));
436         }
437     }
438     // разбор dat файла для создания бэта листов
439     betaListDigestion(fnBetaListsDat, betaLists);
440     // разбор топологии для индексации бета листов
441     aminoacidsIndexation(index, top, aminoacidsIndex);
442     // задание радиуса рассматриваемых соседей
443     nb_.setCutoff(effRad);
444     /*
445      * формально можно сделать:
446      * найти соотношение между именами для углов и индексными номерами
447      * заполнить std::vector< std::vector< size_t > > nameToIndex;
448      * и в последствии считать:
449      * a1 = frame[nameToIndex[0][1] + ... - ... - ...
450      * b1 a2 b2 по аналогии только для [@<-[1, 2, 3]]
451      *
452      */
453 }
454
455 void
456 colorVec::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, gmx::TrajectoryAnalysisModuleData *pdata)
457 {
458     std::vector< gmx::RVec > trajectoryFrame;
459     trajectoryFrame.resize(0);
460     // считывания текущего фрейма траектории
461     for (size_t i {0}; i < fr.natoms; ++i) {
462         trajectoryFrame.push_back(fr.x[i]);
463     }
464     // подсчёт углов в красках
465     std::vector< colorLocalAngles > colorFormation;
466     colorsAnglesEvaluation(trajectoryFrame, colorsNames, colorFormation);
467     // рассчёт положения относительно белка
468     /*
469      * формально можно совместить определение близости с белком и определение ближайших бэта-листов
470      * для этого думаю нужно как-то соотнести или сделать соотношение между индексом и бэта-листами (т.е. хранить что бы за О(1) делать)
471      *
472      */
473
474     std::vector< bool >     colorsToPeptide;
475     colorsToPeptide.resize(0);
476     colorsToPeptide.resize(colorsNames.size(), false);
477     for (size_t i = 0; i < colorsNames.size(); ++i) {
478         colorsToPeptide[i] = isNearPeptide(trajectoryFrame, index, colorsNames[i], effRad);
479     }
480     // расчёт угла и среднего угла с ближайшими бета листами
481     std::vector< std::vector< bool > > colorsToBeta;
482     colorsToBeta.resize(0);
483     colorsToBeta.resize(colorsNames.size());
484     std::vector< std::vector< double > > colorsToBetaAngles;
485     colorsToBetaAngles.resize(0);
486     std::vector< gmx::RVec > betaListsRVecs;
487     betaListsRVecsEvaluation(trajectoryFrame, betaLists[frnr], betaListsRVecs, indexCA);
488     for (size_t i = 0; i < colorsToPeptide.size(); ++i) {
489         if (colorsToPeptide[i]) {
490             searchNearBetaLists(trajectoryFrame, betaLists[frnr], colorsNames[i], colorsToBeta[i], aminoacidsIndex, effRad);
491         }
492     }
493     for (size_t i = 0; i < colorsToBeta.size(); ++i) {
494         for (size_t j = 0; j < colorsToBeta[i].size(); ++j) {
495             if (colorsToBeta[i][j]) {
496                 computeAnglesColorsVsBeta(betaListsRVecs[j], colorFormation[i]);
497             }
498         }
499     }
500     // вывод в файл(ы) информацию для фрейма
501     anglesFileDump(frnr, fnOut, colorsToPeptide, colorFormation);
502     std::cout << std::endl << "frame number " << frnr;
503 }
504
505 void
506 colorVec::finishAnalysis(int nframes)
507 {
508
509 }
510
511 void
512 colorVec::writeOutput()
513 {
514     /*
515      *
516      * frame #<current frame number>
517      * 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>]
518      *
519      */
520 }
521
522 /*! \brief
523  * The main function for the analysis template.
524  */
525 int
526 main(int argc, char *argv[])
527 {
528     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<colorVec>(argc, argv);
529 }