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
52 #include <gromacs/trajectoryanalysis.h>
53 #include <gromacs/utility/smalloc.h>
54 #include <gromacs/math/do_fit.h>
55 #include <gromacs/trajectoryanalysis/topologyinformation.h>
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/topology/index.h"
63 const std::vector< std::vector< std::vector< double > > > read_correlation_matrix_file(const char* file_name, unsigned long size, unsigned int length)
66 std::vector< std::vector< std::vector< double > > > crrlts;
67 file = std::fopen(file_name, "r+");
68 std::vector< double > c;
70 std::vector< std::vector< double > > b;
72 crrlts.resize(length + 1, b);
74 std::cout << "reading correlations from a file:\n";
75 for (unsigned int i = 0; i < length + 1; i++) {
76 //char *t1 = std::fgets(a, 100, file);
77 int t0 = 0, t1 = std::fscanf(file, "%d\n", &t0);
78 std::cout << t0 << " | ";
82 for (unsigned int j = 0; j < size; j++) {
83 for (unsigned int f = 0; f < size; f++) {
84 int t2 = std::fscanf(file, "%lf ", &crrlts[i][j][f]);
93 template< typename T >
94 unsigned int posSrch(std::vector< T > a, T b) {
95 for (unsigned i01 = 0; i01 < a.size(); i01++) {
104 void make_rout_file(double crl_border, std::vector< int > indx, std::vector< std::string > resNames, std::vector< std::vector< std::vector< std::pair< unsigned int, unsigned int > > > > rout, const char* file_name01)
107 file01 = std::fopen(file_name01, "w+");
110 std::vector< std::tuple< int, int, std::vector< int > > > table;
112 for (unsigned int i01 = 0; i01 < rout.size(); i01++) {
114 for (unsigned int i02 = i01 + 1; i02 <rout.size(); i02++) {
115 if (rout[i02] == rout[i01]) {
122 std::fprintf(file01, "tau = %d | correlations >= %0.2f\n\n", i01, crl_border);
124 std::fprintf(file01, "tau = [%d ; %d] | correlations >= %0.2f\n\n", i01, f, crl_border);
126 for (unsigned int i02 = 0; i02 < rout[i01].size(); i02++) {
127 for (unsigned int i03 = 0; i03 < rout[i01][i02].size(); i03++) {
128 std::fprintf(file01, "cgo_arrow (id %3d), (id %3d), radius=0.15\n", indx[rout[i01][i02][i03].first] + 1, indx[rout[i01][i02][i03].second] + 1);
130 std::fprintf(file01, "\n");
133 for (unsigned int i02 = 0; i02 < rout[i01].size(); i02++) {
134 for (unsigned int i03 = 0; i03 < rout[i01][i02].size(); i03++) {
136 for (unsigned int i04 = 0; i04 < table.size(); i04++) {
137 if (std::get<0>(table[i04]) == indx[rout[i01][i02][i03].first] + 1) {
138 std::get<1>(table[i04])++;
139 std::get<2>(table[i04]).push_back(indx[rout[i01][i02][i03].second] + 1);
141 } else if (std::get<0>(table[i04]) == indx[rout[i01][i02][i03].second] + 1) {
142 std::get<1>(table[i04])++;
143 std::get<2>(table[i04]).push_back(indx[rout[i01][i02][i03].first] + 1);
148 std::vector< int > a;
150 a.push_back(indx[rout[i01][i02][i03].second] + 1);
151 table.push_back(std::make_tuple(indx[rout[i01][i02][i03].first] + 1, 1, a));
153 a.push_back(indx[rout[i01][i02][i03].first] + 1);
154 table.push_back(std::make_tuple(indx[rout[i01][i02][i03].second] + 1, 1, a));
158 for (unsigned int i02 = 0; i02 < table.size(); i02++) {
159 std::fprintf(file01, "residue %s connections %d | ", (resNames[posSrch(indx, std::get<0>(table[i02]))]).c_str(), std::get<1>(table[i02]));
160 for (unsigned int i03 = 0; i03 < std::get<2>(table[i03]).size(); i03++) {
161 std::fprintf(file01, "%s ", (resNames[posSrch(indx, std::get<2>(table[i02])[i03])]).c_str());
163 std::fprintf(file01, "\n");
169 * добавить вывод разницы в плане что +, а что -
181 bool mysortfunc (std::vector< int > a, std::vector< int > b) {
182 return (a.size() > b.size());
185 bool isitsubset (std::vector< int > a, std::vector< int > b) {
189 std::sort(a.begin(), a.end());
190 std::sort(b.begin(), b.end());
192 for (unsigned int i = 0; i < a.size(); i++) {
204 void correlation_evaluation(const std::vector< RVec > ref, const std::vector< std::vector< RVec > > trajectory, const unsigned int tauE, const unsigned int window_size, const char* file_name) {
205 std::vector< std::vector< std::vector< double > > > crl;
206 std::vector< std::vector< RVec > > trj;
208 for (unsigned int istart = 0; istart < trajectory.size() - window_size - tauE; istart++) {
210 trj.resize(window_size + tauE);
211 for (unsigned int i = 0; i < window_size + tauE; i++) {
212 for (unsigned int j = 0; j < trajectory[i].size(); j++) {
213 trj[i].push_back(trajectory[i + istart][j]);
216 crl.resize(tauE + 1);
217 for (unsigned int i = 0; i < crl.size(); i++) {
218 crl[i].resize(trajectory.front().size());
219 for (unsigned int j = 0; j < trajectory.front().size(); j++) {
220 crl[i][j].resize(trajectory.front().size());
221 for (unsigned int k = 0; k < trajectory.front().size(); k++) {
226 #pragma omp parallel for ordered schedule(dynamic) shared(crl) firstprivate(trj, ref)
227 for (unsigned int i = 0; i <= tauE; i++) {
228 std::vector< std::vector< double > > a, b, c;
229 std::vector< double > d;
231 for (unsigned int j = 0; j < trajectory.front().size(); j++) {
235 for (unsigned int k = 0; k < trajectory.front().size(); k++) {
241 for (unsigned int j = 0; j < trj.size() - tauE; j++) {
242 for (unsigned int f1 = 0; f1 < trj.front().size(); f1++) {
243 for (unsigned int f2 = 0; f2 < trj.front().size(); f2++) {
245 temp1 = trj[j][f1] - ref[f1];
246 temp2 = trj[j + i][f2] - ref[f2];
247 a[f1][f2] += (static_cast<double>(temp1[0]) * static_cast<double>(temp2[0]) +
248 static_cast<double>(temp1[1]) * static_cast<double>(temp2[1]) +
249 static_cast<double>(temp1[2]) * static_cast<double>(temp2[2]));
250 b[f1][f2] += (static_cast<double>(temp1[0]) * static_cast<double>(temp1[0]) +
251 static_cast<double>(temp1[1]) * static_cast<double>(temp1[1]) +
252 static_cast<double>(temp1[2]) * static_cast<double>(temp1[2]));
253 c[f1][f2] += (static_cast<double>(temp2[0]) * static_cast<double>(temp2[0]) +
254 static_cast<double>(temp2[1]) * static_cast<double>(temp2[1]) +
255 static_cast<double>(temp2[2]) * static_cast<double>(temp2[2]));
259 for (unsigned int j = 0; j < trj.front().size(); j++) {
260 for (unsigned int f = 0; f < trj.front().size(); f++) {
261 crl[i][j][f] = a[j][f] / (std::sqrt(b[j][f] * c[j][f]));
264 /*std::cout << i << " | ";
272 // запись блока в память
273 // заход на новый круг
276 std::cout << "crl calc " << istart << "\n";
278 for (unsigned int i1 = 0; i1 < crl.size(); i1++) {
279 for (unsigned int i2 = 0; i2 < crl[i1].size(); i2++) {
280 for (unsigned int i3 = 0; i3 < crl[i1][i2].size(); i3++) {
281 crl[i1][i2][i3] = std::round(crl[i1][i2][i3] * 10000) / 10000;
288 file = std::fopen(file_name, "a");
289 for (unsigned int i = 0; i < crl.size(); i++) {
290 std::fprintf(file, "%d %d\n", istart, i);
291 for (unsigned int j = 0; j < crl[i].size(); j++) {
292 for (unsigned int f = 0; f < crl[i][j].size(); f++) {
293 std::fprintf(file, "%.6f ", crl[i][j][f]); //~16
295 std::fprintf(file, "\n");
301 std::cout << "file input " << istart << "\n";
304 // запись блока в память
305 // заход на новый круг
312 void graph_calculation( std::vector< std::vector< std::pair< double, long int > > > &graph, std::vector< std::vector< unsigned int > > &s_graph,
313 std::vector< std::vector< std::pair< unsigned int, unsigned int > > > &s_graph_rbr,
314 const std::vector< RVec > ref,
315 const std::vector< std::vector< std::vector<double > > > crl,
316 const double crl_b, const double e_rad, const unsigned int tauB, const unsigned int tauE) {
317 graph.resize(ref.size());
318 for (unsigned int i = 0; i < ref.size(); i++) {
319 graph[i].resize(ref.size(), std::make_pair(0, -1));
322 for (unsigned int i = tauB; i <= tauE; i++) {
323 for (unsigned int j = 0; j < ref.size(); j++) {
324 for (unsigned int f = j; f < ref.size(); f++) {
325 temp = ref[i] - ref[f];
326 if (std::max(std::abs(crl[i][j][f]), std::abs(crl[i][f][j])) >= crl_b &&
327 static_cast<double>(norm(temp)) <= e_rad && std::abs(graph[j][f].first) < std::max(std::abs(crl[i][j][f]), std::abs(crl[i][f][j]))) {
328 if (std::abs(crl[i][j][f]) > std::abs(crl[i][f][j])) {
329 graph[j][f].first = crl[i][j][f];
331 graph[j][f].first = crl[i][f][j];
333 graph[j][f].second = i;
338 std::cout << "crl analysed\n";
339 std::vector< bool > graph_flags;
340 graph_flags.resize(ref.size(), true);
341 std::vector< unsigned int > a;
343 std::vector< std::pair< unsigned int, unsigned int > > b;
345 std::vector< unsigned int > que1, que2, que3;
346 for (unsigned int i = 0; i < ref.size(); i++) {
347 if (graph_flags[i]) {
348 s_graph.push_back(a);
349 s_graph_rbr.push_back(b);
355 graph_flags[i] = false;
356 while(que1.size() > 0) {
358 for (unsigned int k = 0; k < que1.size(); k++) {
359 for (unsigned int j = 0; j < ref.size(); j++) {
360 if (graph[que1[k]][j].second > -1 && graph_flags[j]) {
362 graph_flags[j] = false;
367 for (unsigned int j = 0; j < que2.size(); j++) {
368 que3.push_back(que2[j]);
371 s_graph.back() = que3;
372 for (unsigned int j = 0; j < que3.size(); j++) {
373 for (unsigned int k = 0; k < ref.size(); k++) {
374 if (graph[que3[j]][k].second > -1) {
375 s_graph_rbr.back().push_back(std::make_pair(que3[j], k));
379 //std::cout << s_graph.back().size() << " ";
384 bool myfunction (const std::pair< int, double > i, const std::pair< int, double > j) {
385 return i.second < j.second;
388 bool myfuncRVecEq (const RVec a, const RVec b) {
389 return (a[0] != b[0]) || (a[1] != b[1]) || (a[2] != b[2]);
392 void graph_back_bone_evaluation(std::vector< std::vector < std::pair< unsigned int, unsigned int > > > &rout_n, const unsigned long indxSize,
393 const std::vector< std::vector< std::pair< double, long > > > graph, const std::vector< std::vector< unsigned int > > s_graph,
394 const std::vector< std::vector< std::pair< unsigned int, unsigned int > > > s_graph_rbr) {
395 std::vector< double > key;
396 std::vector< long > path;
397 std::vector< std::pair< unsigned int, double > > que;
398 std::vector< std::pair< unsigned int, unsigned int > > a;
400 for (unsigned int i = 0; i < s_graph.size(); i++) {
405 if (s_graph[i].size() > 2) {
406 key.resize(indxSize, 2);
407 path.resize(indxSize, -1);
408 key[s_graph[i][0]] = 0;
409 for (unsigned int j = 0; j < s_graph[i].size(); j++) {
410 que.push_back(std::make_pair(s_graph[i][j], key[s_graph[i][j]]));
412 std::sort(que.begin(), que.end(), myfunction);
413 while (!que.empty()) {
415 que.erase(que.begin());
416 for (unsigned int j = 0; j < s_graph_rbr[i].size(); j++) {
418 if (s_graph_rbr[i][j].first == v) {
419 u = s_graph_rbr[i][j].second;
420 } else if (s_graph_rbr[i][j].second == v) {
421 u = s_graph_rbr[i][j].first;
424 unsigned long pos = 0;
425 for (unsigned long k = 0; k < que.size(); k++) {
426 if (que[k].first == u) {
432 if (flag && key[static_cast< unsigned long >(u)] > 1 - std::abs(graph[v][static_cast< unsigned long >(u)].first)) {
433 path[static_cast< unsigned long >(u)] = v;
434 key[static_cast< unsigned long >(u)] = 1 - std::abs(graph[v][static_cast< unsigned long >(u)].first);
435 que[pos].second = key[static_cast< unsigned long >(u)];
436 sort(que.begin(), que.end(), myfunction);
442 for (unsigned int j = 0; j < indxSize; j++) {
444 rout_n.back().push_back(std::make_pair(j, path[j]));
451 gmx::RVec evaluate_com(std::vector< RVec > frame) {
456 for (unsigned int i = 0; i < frame.size(); i++) {
459 temp[0] /= frame.size();
460 temp[1] /= frame.size();
461 temp[2] /= frame.size();
466 * \ingroup module_trajectoryanalysis
469 class SpaceTimeCorr : public TrajectoryAnalysisModule
474 virtual ~SpaceTimeCorr();
476 //! Set the options and setting
477 virtual void initOptions(IOptionsContainer *options,
478 TrajectoryAnalysisSettings *settings);
480 //! First routine called by the analysis framework
481 // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
482 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
483 const TopologyInformation &top);
485 virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
486 const t_trxframe &fr);
488 //! Call for each frame of the trajectory
489 // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
490 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
491 TrajectoryAnalysisModuleData *pdata);
493 //! Last routine called by the analysis framework
494 // virtual void finishAnalysis(t_pbc *pbc);
495 virtual void finishAnalysis(int nframes);
497 //! Routine to write output, that is additional over the built-in
498 virtual void writeOutput();
504 std::vector< std::vector< RVec > > trajectory;
505 std::vector< RVec > reference;
507 std::vector< int > index;
508 std::vector< std::string > resNms;
511 int tau = 0; // selectable
513 float crl_border = 0; // selectable
514 float eff_rad = 1.5; // selectable
515 std::string OutPutName; // selectable
516 int mode = 0; // selectable
517 std::string MtrxNm; // selectable
518 // Copy and assign disallowed by base.
521 SpaceTimeCorr::SpaceTimeCorr(): TrajectoryAnalysisModule()
525 SpaceTimeCorr::~SpaceTimeCorr()
530 SpaceTimeCorr::initOptions(IOptionsContainer *options,
531 TrajectoryAnalysisSettings *settings)
533 static const char *const desc[] = {
534 "[THISMODULE] to be done"
536 // Add the descriptive text (program help text) to the options
537 settings->setHelpText(desc);
538 // Add option for output file name
539 //options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
540 // .store(&fnNdx_).defaultBasename("domains")
541 // .description("Index file from the domains"));
542 // Add option for working mode
543 options->addOption(gmx::IntegerOption("mode")
545 .description("default 0 | rdy correlation matrixes 1, need extra params"));
546 // Add option for Matrix Input file names
547 options->addOption(StringOption("Mtrx_in_put")
549 .description("mandatory if work mode == 1"));
550 // Add option for tau constant
551 options->addOption(gmx::IntegerOption("tau")
553 .description("number of frames for time travel"));
554 // Add option for window constant
555 options->addOption(gmx::IntegerOption("window")
557 .description("size of window in frames for correlation evaluation"));
558 // Add option for crl_border constant
559 options->addOption(FloatOption("crl")
561 .description("make graph based on corrs > constant"));
562 // Add option for eff_rad constant
563 options->addOption(FloatOption("ef_rad")
565 .description("effective radius for atoms to evaluate corrs"));
566 // Add option for selection list
567 options->addOption(SelectionOption("select_domains_and_residue").storeVector(&sel_)
568 .required().dynamicMask().multiValue()
569 .description("Domains to form rigid skeleton"));
570 // Add option for output file names
571 options->addOption(StringOption("out_put")
573 .description("<your name here> + <local file tag>.txt"));
574 // Control input settings
575 settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
576 settings->setFlag(TrajectoryAnalysisSettings::efUseTopX);
577 settings->setPBC(true);
581 SpaceTimeCorr::initAnalysis(const TrajectoryAnalysisSettings &settings, const TopologyInformation &top)
583 ArrayRef< const int > atomind = sel_[0].atomIndices();
585 for (ArrayRef< const int >::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
586 index.push_back(*ai);
587 resNms.push_back(*(top.atoms()->resinfo[top.atoms()->atom[*ai].resind].name) + std::to_string((top.atoms()->atom[*ai].resind + 1)));
589 trajectory.resize(0);
592 if (top.hasFullTopology()) {
593 for (unsigned int i = 0; i < index.size(); i++) {
594 reference.push_back(top.x().at(index[i]));
600 SpaceTimeCorr::initAfterFirstFrame(const TrajectoryAnalysisSettings &settings, const t_trxframe &fr)
604 // -s '*.pdb' -f '*.pdb' -n '*.ndx' -sf 'name' -tau 1000 -crl 0.20 -ef_rad 9000 -out_put OLA
605 // -s '*.tpr' -f '*.xtc' -n '*.ndx' -sf 'name' -tau 1000 -crl 0.30 -ef_rad 20 -out_put 'test_run'
608 SpaceTimeCorr::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, TrajectoryAnalysisModuleData *pdata)
610 trajectory.resize(trajectory.size() + 1);
611 trajectory.back().resize(index.size());
612 for (unsigned int i = 0; i < index.size(); i++) {
613 trajectory.back()[i] = fr.x[index[i]];
619 SpaceTimeCorr::finishAnalysis(int nframes)
621 std::vector< std::vector< std::vector< std::pair< unsigned int, unsigned int > > > > rout_new;
622 unsigned int k = 1000, m = 0;
624 k = static_cast< unsigned int>(tau);
631 std::cout << "\nCorrelation's evaluation - start\n";
633 correlation_evaluation(reference, trajectory, k, 1000, (OutPutName + "-DumpData.txt").c_str());
635 std::cout << "Correlation's evaluation - end\n";
636 } else if (mode == 1) {
637 rout_new.resize(trajectory.size() - static_cast< unsigned long>(window) - k);
639 file = std::fopen((OutPutName + "-DumpData.txt").c_str(), "r+");
640 #pragma omp parallel for ordered schedule(dynamic) shared(rout_new) firstprivate(reference, crl_border, eff_rad, m, k)
641 for(unsigned long i = 0; i < trajectory.size() - static_cast< unsigned int>(window) - k; i++) {
642 std::vector< std::vector< std::vector< double > > > crltns;
643 std::vector< std::vector< std::pair< double, long int > > > graph;
644 std::vector< std::vector< unsigned int > > sub_graph;
645 std::vector< std::vector< std::pair< unsigned int, unsigned int > > > sub_graph_rbr;
649 crltns.resize(k + 1);
650 for (unsigned int i1 = 0; i1 < crltns.size(); i1++) {
651 crltns[i1].resize(trajectory.front().size());
652 for (unsigned int i2 = 0; i2 < trajectory.front().size(); i2++) {
653 crltns[i1][i2].resize(trajectory.front().size());
654 for (unsigned int i3 = 0; i3 < trajectory.front().size(); i3++) {
655 crltns[i1][i2][i3] = 0.;
659 for (unsigned int i1 = 0; i1 < k + 1; i1++) {
660 int t0, t1, t2 = std::fscanf(file, "%d %d\n", &t0, &t1);
661 for (unsigned int i2 = 0; i2 < trajectory.front().size(); i2++) {
662 for (unsigned int i3 = 0; i3 < trajectory.front().size(); i3++) {
663 t2 = std::fscanf(file, "%lf ", &crltns[i1][i2][i3]);
667 std::cout << "red mtrxs " << i << "\n";
670 graph_calculation(graph, sub_graph, sub_graph_rbr, reference, crltns, static_cast< double >(crl_border), static_cast< double >(eff_rad), m, k);
671 graph_back_bone_evaluation(rout_new[i], reference.size()/*index.size()*/, graph, sub_graph, sub_graph_rbr);
677 make_rout_file(static_cast< double >(crl_border), index, resNms, rout_new, (OutPutName + "-routs.txt").c_str());
679 std::cout << "Finish Analysis - end\n\n";
683 SpaceTimeCorr::writeOutput()
689 * The main function for the analysis template.
692 main(int argc, char *argv[])
694 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<SpaceTimeCorr>(argc, argv);