2 * This file is part of the GROMACS molecular simulation package.
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.
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::Freevolume.
39 * \author Titov Anatoly <Wapuk-cobaka@yandex.ru>
40 * \ingroup module_trajectoryanalysis
43 #include <gromacs/trajectoryanalysis.h>
44 #include <gromacs/trajectoryanalysis/topologyinformation.h>
45 #include <gromacs/selection/nbsearch.h>
56 // структура углов для одной краски
57 struct colorLocalAngles {
62 std::vector< double > betaAngles;
65 // хрен пойми почему, но захотелось рекурсивно сделать | можно сделать и через вайл
66 long long returnBetaEnd(const std::string ¤tLine, long long pos) {
67 switch (currentLine[pos]) {
69 returnBetaEnd(currentLine, pos + 1);
72 returnBetaEnd(currentLine, pos + 1);
78 return (pos - 1); // формально логически лишняя строчка, но Qt ругается
81 // функция парсинга одной строки
82 void parseBetaListDATLine(const std::string ¤tLine, std::vector< std::vector< unsigned int > > &localInputBL) {
83 size_t equalCount = 0;
84 std::vector< unsigned int > a;
86 for (size_t i = 0; i < currentLine.size(); ++i) {
87 if (currentLine[i] == '=') { // подсчитываем число "пустых" символов
90 long long temp = returnBetaEnd(currentLine, i);
91 if (temp - static_cast< long long >(i) > 3) {
92 localInputBL.push_back(a);
93 for (size_t j = i; j <= temp; ++j) {
94 localInputBL.back().push_back(j - equalCount);
102 // функция нахождения бета-листов в структуре по файлу ДССП
103 void betaListDigestion(const std::string &inputFile, std::vector< std::vector< std::vector< unsigned int > > > &inputBL) {
105 std::ifstream file(inputFile);
107 getline(file, line); // считываем число в первой строке - кол-во осмысленных элементов в строках - нам не нужно
109 if (line.size() > 3) {
111 inputBL.resize(inputBL.size() + 1);
112 parseBetaListDATLine(line, inputBL.back());
115 } while (line.size() > 3);
117 throw "DSSP DAT FILE IS EMPTY";
122 // функция выделения индексов для бэталистов
123 inline void aminoacidsIndexation(const std::vector< size_t > &inputIndex, const gmx::TopologyInformation &top, std::vector< std::vector< size_t > > &aminoacidsIndex) {
124 aminoacidsIndex.resize(0);
125 for (size_t i = 0; i < inputIndex.size(); ++i) {
126 aminoacidsIndex.resize(std::max(aminoacidsIndex.size(), static_cast< size_t >(top.atoms()->atom[inputIndex[i]].resind + 1)));
127 aminoacidsIndex[top.atoms()->atom[inputIndex[i]].resind].push_back(inputIndex[i]);
131 // функция поиска RVec в кадре по имени->индексу
132 gmx::RVec returnRVec(const std::vector< gmx::RVec > &frame, const std::vector< std::pair< std::string, size_t > > &colorsIndex, std::string toFind) {
133 for (auto &i : colorsIndex) {
134 if (i.first == toFind) {
135 return frame[i.second];
138 throw "WRONG COLOR ATOM NAME TO EVALUATE VECTORS";
141 // функция векторного произведения двух RVec'ов
142 inline void RVecVecMultiply(gmx::RVec &n, const gmx::RVec &a, const gmx::RVec &b) {
143 n[0] = a[1] * b[2] - a[2] * b[1];
144 n[1] = a[2] * b[0] - a[0] * b[2];
145 n[2] = a[0] * b[1] - a[2] * b[0];
148 // поиск угла между двумя RVec'ами
149 inline double RVecAngle(const gmx::RVec &a, const gmx::RVec &b) {
150 return std::acos((a[0] * b[0] + a[1] * b[1] + a[2] * b[2]) / (a.norm() * b.norm()));
153 // вычисление внутренних углов в краске
154 void colorsAnglesEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< std::pair< std::string, size_t > > > &colorsIndex, std::vector< colorLocalAngles > &colorStruct) {
155 colorStruct.resize(colorsIndex.size());
156 #pragma omp parallel for ordered schedule(dynamic)
157 for (size_t i = 0; i < colorsIndex.size(); ++i) {
158 colorStruct[i].betaAngles.resize(0);
160 colorStruct[i].a1 = returnRVec(frame, colorsIndex[i], "CAF") - returnRVec(frame, colorsIndex[i], "CAJ") + returnRVec(frame, colorsIndex[i], "CAK") - returnRVec(frame, colorsIndex[i], "CAO");
161 colorStruct[i].a1 /= colorStruct[i].a1.norm();
162 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");
163 colorStruct[i].b1 /= colorStruct[i].b1.norm();
164 RVecVecMultiply(colorStruct[i].n1, colorStruct[i].a1, colorStruct[i].b1);
165 colorStruct[i].n1 /= colorStruct[i].n1.norm();
167 colorStruct[i].a2 = returnRVec(frame, colorsIndex[i], "CAB") - returnRVec(frame, colorsIndex[i], "CBQ") + returnRVec(frame, colorsIndex[i], "CBP") - returnRVec(frame, colorsIndex[i], "CBL");
168 colorStruct[i].a2 /= colorStruct[i].a2.norm();
169 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");
170 colorStruct[i].b2 /= colorStruct[i].b2.norm();
171 RVecVecMultiply(colorStruct[i].n2, colorStruct[i].a2, colorStruct[i].b2);
172 colorStruct[i].n2 *= -1; // для "сонаправленности" векторов нормали
173 colorStruct[i].n2 /= colorStruct[i].n2.norm();
174 colorStruct[i].a12 = RVecAngle(colorStruct[i].a1, colorStruct[i].a2);
175 colorStruct[i].b12 = RVecAngle(colorStruct[i].b1, colorStruct[i].b2);
176 colorStruct[i].n12 = RVecAngle(colorStruct[i].n1, colorStruct[i].n2);
181 // вычисление направляющего вектора в бэта-листе
182 inline void betaListsRVecsEvaluation(const std::vector< gmx::RVec > &frame, const std::vector< std::vector< unsigned int > > &inputBetaLists, std::vector< gmx::RVec > &temp, const std::vector< size_t > &inputCA) {
185 for (const auto &i : inputBetaLists) {
186 tempA = frame[inputCA[i[i.size() - 1]]] + frame[inputCA[i[i.size() - 2]]] - frame[inputCA[i[1]]] - frame[inputCA[i[0]]];
187 tempA /= tempA.norm();
188 temp.push_back(tempA);
192 // определение близко ли к белку находится краска
193 bool isNearPeptide(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > &inputFrame,
194 gmx::AnalysisNeighborhood &nbhood, const std::vector< size_t > &inputIndex,
195 const std::vector< std::pair< std::string, size_t > > &inputColor) {
196 gmx::AnalysisNeighborhoodSearch nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
197 gmx::AnalysisNeighborhoodPair pair;
198 for (const auto &i : inputColor) {
199 gmx::AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(inputFrame[i.second].as_vec());
200 while (pairSearch.findNextPair(&pair)) {
201 for (const auto &j : inputIndex) {
202 if (pair.refIndex() == j) {
211 // поиск ближайших к краске бэта-листов
212 inline void searchNearBetaLists(const t_trxframe &fr, const t_pbc *pbc, const std::vector< gmx::RVec > inputFrame,
213 gmx::AnalysisNeighborhood &nbhood,
214 const std::vector< std::vector< unsigned int > > inputBLists,
215 const std::vector< std::pair< std::string, size_t > > &inputColor,
216 std::vector< bool > &outputList, const std::vector< std::vector< size_t > > &inputAminoacids) {
217 outputList.resize(0);
218 outputList.resize(inputBLists.size(), false);
219 gmx::AnalysisNeighborhoodSearch nbsearch = nbhood.initSearch(pbc, gmx::AnalysisNeighborhoodPositions(fr.x, fr.natoms));
220 gmx::AnalysisNeighborhoodPair pair;
221 for (const auto &i : inputColor) {
222 while (nbsearch.startPairSearch(inputFrame[i.second].as_vec()).findNextPair(&pair)) {
223 for (size_t j = 0; j < inputBLists.size(); ++j) {
224 for (size_t k = 0; k < inputBLists[j].size(); ++k) {
225 for (size_t m = 0; m < inputAminoacids[inputBLists[j][k]].size(); ++m) {
226 if (pair.testIndex() == inputAminoacids[inputBLists[j][k]][m]) {
227 outputList[j] = true;
237 // определение углов между краской и направляющим вектором бэта-листа
238 inline void computeAnglesColorsVsBeta(const gmx::RVec &inputBetaRVec, colorLocalAngles &colorFormation) {
239 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.a1));
240 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.b1));
241 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.n1));
242 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.a2));
243 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.b2));
244 colorFormation.betaAngles.push_back(RVecAngle(inputBetaRVec, colorFormation.n2));
247 // функция записи в файл значений углов для кадра
248 void anglesFileDump(const int frameNum, const std::string &output, const std::vector< bool > &toPeptide, const std::vector< colorLocalAngles > &colorFormation) {
249 std::ofstream file(output);
251 std::vector< double > betaTemp;
253 file << "frame =" << std::setw(8) << frameNum << std::endl;
254 for (size_t i = 0; i < colorFormation.size(); ++i) {
255 file << "color #" << std::setw(3) << i;
257 file << std::setw(4) << "yes";
259 file << std::setw(4) << "no";
261 file << std::setw(7) << std::setprecision(2) << colorFormation[i].a12 << colorFormation[i].b12 << colorFormation[i].n12;
262 temp = colorFormation[i].betaAngles.size() / 6;
264 file << std::setw(4) << 0;
266 betaTemp.resize(6, 0.); // magic number, meh
267 for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
268 betaTemp[j % 6] += colorFormation[i].betaAngles[j];
270 for (size_t j = 0; j < betaTemp.size(); ++j) {
273 file << std::setw(4) << temp;
274 for (size_t j = 0; j < betaTemp.size(); ++j) {
275 file << std::setw(7) << std::setprecision(2) << betaTemp[j];
277 for (size_t j = 0; j < colorFormation[i].betaAngles.size(); ++j) {
278 file << std::setw(7) << std::setprecision(2) << colorFormation[i].betaAngles[j];
286 * \ingroup module_trajectoryanalysis
288 class colorVec : public gmx::TrajectoryAnalysisModule
295 //! Set the options and setting
296 virtual void initOptions(gmx::IOptionsContainer *options,
297 gmx::TrajectoryAnalysisSettings *settings);
299 //! First routine called by the analysis framework
300 // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
301 virtual void initAnalysis(const gmx::TrajectoryAnalysisSettings &settings,
302 const gmx::TopologyInformation &top);
304 //! Call for each frame of the trajectory
305 // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
306 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
307 gmx::TrajectoryAnalysisModuleData *pdata);
309 //! Last routine called by the analysis framework
310 // virtual void finishAnalysis(t_pbc *pbc);
311 virtual void finishAnalysis(int nframes);
313 //! Routine to write output, that is additional over the built-in
314 virtual void writeOutput();
318 gmx::SelectionList sel_;
319 std::string fnOut {"name"}; // selectable
320 std::string fnBetaListsDat; // selectable
321 float effRad {0.8}; // selectable
322 std::vector< size_t > index;
323 std::vector< size_t > indexCA;
324 std::vector< std::vector< size_t > > aminoacidsIndex;
325 std::vector< std::vector< std::pair< std::string, size_t > > > colorsNames;
326 std::vector< std::vector< std::vector< unsigned int > > > betaLists;
327 gmx::AnalysisNeighborhood nb_;
329 // Copy and assign disallowed by base.
332 colorVec::colorVec(): TrajectoryAnalysisModule()
336 colorVec::~colorVec()
349 colorVec::initOptions( gmx::IOptionsContainer *options,
350 gmx::TrajectoryAnalysisSettings *settings)
352 static const char *const desc[] = {
353 "[THISMODULE] to be done"
355 // Add the descriptive text (program help text) to the options
356 settings->setHelpText(desc);
357 // Add option for selection list
358 options->addOption(gmx::SelectionOption("sel")
360 .required().dynamicMask().multiValue()
361 .description("select pepride and colors / -sf"));
362 // Add option for input file names
363 options->addOption(gmx::StringOption("dat")
364 .store(&fnBetaListsDat)
365 .description("a file to make dynamic beta lists"));
366 // Add option for output file name
367 options->addOption(gmx::StringOption("out")
369 .description("Index file for the algorithm output."));
370 // Add option for effRad constant
371 options->addOption(gmx::FloatOption("efRad")
373 .description("max distance from colors to peptide in nm to consider to be \"near\""));
374 // Control input settings
375 settings->setFlags(gmx::TrajectoryAnalysisSettings::efNoUserPBC);
376 settings->setFlag(gmx::TrajectoryAnalysisSettings::efUseTopX);
377 //settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
378 settings->setPBC(true);
382 colorVec::initAnalysis( const gmx::TrajectoryAnalysisSettings &settings,
383 const gmx::TopologyInformation &top)
385 // считывание индекса
387 for (gmx::ArrayRef< const int >::iterator ai {sel_.front().atomIndices().begin()}; (ai < sel_.front().atomIndices().end()); ai++) {
388 index.push_back(static_cast< size_t >(*ai));
390 // считывание индекса
392 for (gmx::ArrayRef< const int >::iterator ai {sel_[1].atomIndices().begin()}; (ai < sel_[2].atomIndices().end()); ai++) {
393 index.push_back(static_cast< size_t >(*ai));
396 colorsNames.resize(sel_.size() - 2);
397 for (size_t i = 2; i < sel_.size(); ++i) {
398 for (gmx::ArrayRef< const int >::iterator ai {sel_[i].atomIndices().begin()}; (ai < sel_[i].atomIndices().end()); ai++) {
399 colorsNames[i - 2].push_back(std::make_pair(*(top.atoms()->atomname[*ai]), static_cast< size_t >(*ai)));
402 // разбор dat файла для создания бэта листов
403 betaListDigestion(fnBetaListsDat, betaLists);
404 // разбор топологии для индексации бета листов
405 aminoacidsIndexation(index, top, aminoacidsIndex);
406 // задание радиуса рассматриваемых соседей
407 nb_.setCutoff(effRad);
409 * формально можно сделать:
410 * найти соотношение между именами для углов и индексными номерами
411 * заполнить std::vector< std::vector< size_t > > nameToIndex;
412 * и в последствии считать:
413 * a1 = frame[nameToIndex[0][1] + ... - ... - ...
414 * b1 a2 b2 по аналогии только для [@<-[1, 2, 3]]
420 colorVec::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, gmx::TrajectoryAnalysisModuleData *pdata)
422 std::cout << "\tFrame analisis start ..." << std::endl;
423 std::vector< gmx::RVec > trajectoryFrame;
424 trajectoryFrame.resize(0);
425 // считывания текущего фрейма траектории
426 for (size_t i {0}; i < index.size(); ++i) {
427 trajectoryFrame.push_back(fr.x[index[i]]);
429 // подсчёт углов в красках
430 std::vector< colorLocalAngles > colorFormation;
431 std::cout << "\t\tColors' angles evaluation." << std::endl;
432 colorsAnglesEvaluation(trajectoryFrame, colorsNames, colorFormation);
433 // рассчёт положения относительно белка
435 * формально можно совместить определение близости с белком и определение ближайших бэта-листов
436 * для этого думаю нужно как-то соотнести или сделать соотношение между индексом и бэта-листами (т.е. хранить что бы за О(1) делать)
440 std::vector< bool > colorsToPeptide;
441 colorsToPeptide.resize(0);
442 colorsToPeptide.resize(colorsNames.size(), false);
443 std::cout << "\t\tWhich colors are \"near\" the pepride mass." << std::endl;
444 for (size_t i = 0; i < colorsNames.size(); ++i) {
445 colorsToPeptide[i] = isNearPeptide(fr, pbc, trajectoryFrame, nb_, index, colorsNames[i]);
447 // расчёт угла и среднего угла с ближайшими бета листами
448 std::vector< std::vector< bool > > colorsToBeta;
449 colorsToBeta.resize(0);
450 colorsToBeta.resize(colorsNames.size());
451 std::vector< std::vector< double > > colorsToBetaAngles;
452 colorsToBetaAngles.resize(0);
453 std::vector< gmx::RVec > betaListsRVecs;
454 std::cout << "\t\tBeta-lists' RVecs search." << std::endl;
455 betaListsRVecsEvaluation(trajectoryFrame, betaLists[frnr], betaListsRVecs, indexCA);
456 std::cout << "\t\tSearching nearby beta-lists." << std::endl;
457 for (size_t i = 0; i < colorsToPeptide.size(); ++i) {
458 if (colorsToPeptide[i]) {
459 searchNearBetaLists(fr, pbc, trajectoryFrame, nb_, betaLists[frnr], colorsNames[i], colorsToBeta[i], aminoacidsIndex);
462 std::cout << "\t\tComputing angles colors vs betas." << std::endl;
463 for (size_t i = 0; i < colorsToBeta.size(); ++i) {
464 for (size_t j = 0; j < colorsToBeta[i].size(); ++j) {
465 if (colorsToBeta[i][j]) {
466 computeAnglesColorsVsBeta(betaListsRVecs[j], colorFormation[i]);
470 // вывод в файл(ы) информацию для фрейма
471 std::cout << "\t\tDumping data." << std::endl;
472 anglesFileDump(frnr, fnOut, colorsToPeptide, colorFormation);
473 std::cout << "\tFrame analysed." << std::endl;
477 colorVec::finishAnalysis(int nframes)
483 colorVec::writeOutput()
487 * frame #<current frame number>
488 * 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>]
491 std::cout << "\n\t colorAngles finished successfully" << std::endl;
495 * The main function for the analysis template.
498 main(int argc, char *argv[])
500 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<colorVec>(argc, argv);