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
51 #include <gromacs/trajectoryanalysis.h>
52 #include <gromacs/utility/smalloc.h>
53 #include <gromacs/math/do_fit.h>
59 void make_correlation_matrix_file(std::vector< std::vector< std::vector< float > > > correlations, const char* file_name, int start)
62 file = std::fopen(file_name, "w+");
63 for (int i = start; i < correlations.size(); i++) {
64 if (correlations.size() - start > 1) {
65 std::fprintf(file, "correlation between 'n' and 'n + %d' frames\n", i);
68 std::cout << "correlation between 'n' and 'n + " << i << "' frames\n";
70 for (int j = 0; j < correlations[i].size(); j++) {
71 for (int f = 0; f < correlations[i][j].size(); f++) {
72 std::fprintf(file, "%3.4f ", correlations[i][j][f]);
74 std::fprintf(file, "\n");
80 void make_correlation_pairs_file(std::vector< std::vector< std::vector< float > > > correlations, const char* file_name, int start)
83 file = std::fopen(file_name, "w+");
84 for (int i = 0; i < correlations.front().size(); i++) {
85 for (int j = 0; j < correlations.front().front().size(); j++) {
86 std::fprintf(file, "correlation between point '%d' and point '%d'\n", i, j);
87 for (int k = 0; k < correlations.size(); k++) {
88 std::fprintf(file, "%3.4f ", correlations[k][i][j]);
90 std::fprintf(file, "\n");
93 std::cout << "correlations in row " << i << " out of " << correlations.front().size() << " completed\n";
99 void make_rout_file(float crl_border, std::vector< int > indx, std::vector< std::vector< std::pair< int, int > > > rout, const char* file_name)
102 file = std::fopen(file_name, "w+");
103 std::fprintf(file, "correlations >= %0.2f\n\n", crl_border);
104 for (int i = 0; i < rout.size(); i++) {
105 for (int j = 0; j < rout[i].size(); j++) {
106 std::fprintf(file, "cgo_arrow (id %3d), (id %3d), radius=0.15\n", indx[rout[i][j].first]/* + 1*/, indx[rout[i][j].second]/* + 1*/);
108 std::fprintf(file, "\n\n");
113 void make_best_corrs_graphics(std::vector< std::vector< std::vector< float > > > correlations,
114 std::vector< std::vector< std::pair< int, int > > > rout_pairs,
115 std::vector< int > indx,
116 const char* file_name)
119 file = std::fopen(file_name, "w+");
120 for (int i = 0; i < rout_pairs.size(); i++) {
121 for (int j = 0; j < rout_pairs[i].size(); j++) {
122 std::fprintf(file, "%3d %3d\n", indx[rout_pairs[i][j].first]/* + 1*/, indx[rout_pairs[i][j].second]/* + 1*/);
123 for (int k = 0; k < correlations.size(); k++) {
124 std::fprintf(file, "%3.5f ", correlations[k][rout_pairs[i][j].first][rout_pairs[i][j].second]);
126 std::fprintf(file, "\n");
132 bool mysortfunc (std::vector< int > a, std::vector< int > b) {
133 return (a.size() > b.size());
136 bool isitsubset (std::vector< int > a, std::vector< int > b) {
140 std::sort(a.begin(), a.end());
141 std::sort(b.begin(), b.end());
143 for (int i = 0; i < a.size(); i++) {
155 void correlation_evaluation(std::vector< std::vector< RVec > > traj, int b_frame, std::vector< std::vector< std::vector< float > > > &crl, int tauS, int tauE) {
156 crl.resize(tauE + 1);
157 for (int i = 0; i < crl.size(); i++) {
158 crl[i].resize(traj.front().size());
159 for (int j = 0; j < crl[i].size(); j++) {
160 crl[i][j].resize(traj.front().size(), 0);
168 std::vector< float > d;
169 d.resize(traj.front().size(), 0);
171 #pragma omp parallel shared(crl, temp_zero, d)
173 #pragma omp for schedule(dynamic)
174 for (int i = tauS; i <= tauE; i += 1) {
175 std::vector< RVec > medx, medy;
178 medx.resize(traj.front().size(), temp_zero);
179 medy.resize(traj.front().size(), temp_zero);
180 for (int j = 0; j < traj.size() - i - 1; j++) {
181 for (int f = 0; f < traj.front().size(); f++) {
182 medx[f] += (traj[b_frame][f] - traj[j][f]);
183 medy[f] += (traj[b_frame][f] - traj[j + i][f]);
186 for (int j = 0; j < traj.front().size(); j++) {
187 tmp = traj.size() - 1;
190 //tmp = 1 / (traj.size() - 1 - i);
194 std::vector< std::vector< float > > a, b, c;
195 a.resize(traj.front().size(), d);
196 b.resize(traj.front().size(), d);
197 c.resize(traj.front().size(), d);
198 for (int j = 0; j < traj.size() - i - 1; j++) {
199 for (int f1 = 0; f1 < traj.front().size(); f1++) {
200 for (int f2 = 0; f2 < traj.front().size(); f2++) {
201 temp1 = traj[b_frame][f1] - traj[j][f1] - medx[f1];
202 temp2 = traj[b_frame][f2] - traj[j + i][f2] - medy[f2];
203 a[f1][f2] += (temp1[0] * temp2[0] + temp1[1] * temp2[1] + temp1[2] * temp2[2]);
204 b[f1][f2] += (temp1[0] * temp1[0] + temp1[1] * temp1[1] + temp1[2] * temp1[2]);
205 c[f1][f2] += (temp2[0] * temp2[0] + temp2[1] * temp2[1] + temp2[2] * temp2[2]);
209 for (int j = 0; j < traj.front().size(); j++) {
210 for (int f = 0; f < traj.front().size(); f++) {
211 crl[i][j][f] = a[j][f] / (std::sqrt(b[j][f] * c[j][f]));
216 std::cout << i << " corr done\n";
222 void graph_calculation(std::vector< std::vector< std::pair< float, int > > > &graph, std::vector< std::vector< int > > &s_graph, std::vector< std::vector< std::pair< int, int > > > &s_graph_rbr,
223 std::vector< std::vector< RVec > > traj, int b_frame,
224 std::vector< std::vector< std::vector< float > > > crl, float crl_b, float e_rad, int tauE) {
225 graph.resize(traj.front().size());
226 for (int i = 0; i < traj.front().size(); i++) {
227 graph[i].resize(traj.front().size(), std::make_pair(0, -1));
230 for (int i = 1; i <= tauE; i++) {
231 for (int j = 0; j < traj.front().size(); j++) {
232 for (int f = j; f < traj.front().size(); f++) {
233 temp = traj[b_frame][j] - traj[b_frame][f];
234 if (std::max(std::abs(crl[i][j][f]), std::abs(crl[i][f][j])) >= crl_b && norm(temp) <= e_rad && std::abs(graph[j][f].first) < std::max(std::abs(crl[i][j][f]), std::abs(crl[i][f][j]))) {
235 if (std::abs(crl[i][j][f]) > std::abs(crl[i][f][j])) {
236 graph[j][f].first = crl[i][j][f];
238 graph[j][f].first = crl[i][f][j];
240 graph[j][f].second = i;
245 std::cout << "crl analysed\n";
246 std::vector< bool > graph_flags;
247 graph_flags.resize(traj.front().size(), true);
248 std::vector< int > a;
250 std::vector< std::pair< int, int > > b;
252 std::vector< int > que1, que2, que3;
253 for (int i = 0; i < traj.front().size(); i++) {
254 if (graph_flags[i]) {
255 s_graph.push_back(a);
256 s_graph_rbr.push_back(b);
262 graph_flags[i] = false;
263 while(que1.size() > 0) {
265 for (int k = 0; k < que1.size(); k++) {
266 for (int j = 0; j < traj.front().size(); j++) {
267 if (graph[que1[k]][j].second > -1 && graph_flags[j]) {
269 graph_flags[j] = false;
274 for (int j = 0; j < que2.size(); j++) {
275 que3.push_back(que2[j]);
278 s_graph.back() = que3;
279 for (int j = 0; j < que3.size(); j++) {
280 for (int k = 0; k < traj.front().size(); k++) {
281 if (graph[que3[j]][k].second > -1) {
282 s_graph_rbr.back().push_back(std::make_pair(que3[j], k));
286 //std::cout << s_graph.back().size() << " ";
291 bool myfunction (const std::pair< int, float > i, const std::pair< int, float > j) {
292 return i.second < j.second;
295 void graph_back_bone_evaluation(std::vector< std::vector < std::pair< int, int > > > &rout_n, int indxSize,
296 std::vector< std::vector< std::pair< float, int > > > graph, std::vector< std::vector< int > > s_graph, std::vector< std::vector< std::pair< int, int > > > s_graph_rbr) {
297 std::vector< float > key;
298 std::vector< int > path;
299 std::vector< std::pair< int, float > > que;
300 std::vector< std::pair< int, int > > a;
302 for (int i = 0; i < s_graph.size(); i++) {
307 if (s_graph[i].size() > 2) {
308 key.resize(indxSize, 2);
309 path.resize(indxSize, -1);
310 key[s_graph[i][0]] = 0;
311 for (int j = 0; j < s_graph[i].size(); j++) {
312 que.push_back(std::make_pair(s_graph[i][j], key[s_graph[i][j]]));
314 std::sort(que.begin(), que.end(), myfunction);
315 while (!que.empty()) {
317 que.erase(que.begin());
318 for (int j = 0; j < s_graph_rbr[i].size(); j++) {
320 if (s_graph_rbr[i][j].first == v) {
321 u = s_graph_rbr[i][j].second;
322 } else if (s_graph_rbr[i][j].second == v) {
323 u = s_graph_rbr[i][j].first;
327 for (int k = 0; k < que.size(); k++) {
328 if (que[k].first == u) {
334 if (flag && key[u] > 1 - std::abs(graph[v][u].first)) {
336 key[u] = 1 - std::abs(graph[v][u].first);
337 que[pos].second = key[u];
338 sort(que.begin(), que.end(), myfunction);
344 for (int j = 0; j < indxSize; j++) {
346 rout_n.back().push_back(std::make_pair(j, path[j]));
354 * \ingroup module_trajectoryanalysis
357 class SpaceTimeCorr : public TrajectoryAnalysisModule
362 virtual ~SpaceTimeCorr();
364 //! Set the options and setting
365 virtual void initOptions(IOptionsContainer *options,
366 TrajectoryAnalysisSettings *settings);
368 //! First routine called by the analysis framework
369 // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
370 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
371 const TopologyInformation &top);
373 virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
374 const t_trxframe &fr);
376 //! Call for each frame of the trajectory
377 // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
378 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
379 TrajectoryAnalysisModuleData *pdata);
381 //! Last routine called by the analysis framework
382 // virtual void finishAnalysis(t_pbc *pbc);
383 virtual void finishAnalysis(int nframes);
385 //! Routine to write output, that is additional over the built-in
386 virtual void writeOutput();
392 std::vector< std::vector< RVec > > trajectory;
394 std::vector< int > index;
397 int tau = 0; // selectable
398 float crl_border = 0; // selectable
399 float eff_rad = 1.5; // selectable
400 std::string OutPutName; // selectable
401 // Copy and assign disallowed by base.
404 SpaceTimeCorr::SpaceTimeCorr(): TrajectoryAnalysisModule()
408 SpaceTimeCorr::~SpaceTimeCorr()
413 SpaceTimeCorr::initOptions(IOptionsContainer *options,
414 TrajectoryAnalysisSettings *settings)
416 static const char *const desc[] = {
417 "[THISMODULE] to be done"
419 // Add the descriptive text (program help text) to the options
420 settings->setHelpText(desc);
421 // Add option for output file name
422 //options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
423 // .store(&fnNdx_).defaultBasename("domains")
424 // .description("Index file from the domains"));
425 // Add option for tau constant
426 options->addOption(gmx::IntegerOption("tau")
428 .description("number of frames for time travel"));
429 // Add option for crl_border constant
430 options->addOption(FloatOption("crl")
432 .description("make graph based on corrs > constant"));
433 // Add option for eff_rad constant
434 options->addOption(FloatOption("ef_rad")
436 .description("effective radius for atoms to evaluate corrs"));
437 // Add option for selection list
438 options->addOption(SelectionOption("select_domains_and_residue").storeVector(&sel_)
439 .required().dynamicMask().multiValue()
440 .description("Domains to form rigid skeleton"));
441 options->addOption(StringOption("out_put")
443 .description("<your name here> + <local file tag>.txt"));
444 // Control input settings
445 settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
446 settings->setPBC(true);
450 SpaceTimeCorr::initAnalysis(const TrajectoryAnalysisSettings &settings,
451 const TopologyInformation &top)
453 ArrayRef< const int > atomind = sel_[0].atomIndices();
455 for (ArrayRef< const int >::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
456 index.push_back(*ai);
458 trajectory.resize(0);
462 SpaceTimeCorr::initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
463 const t_trxframe &fr)
467 // -s '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -f '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -n '/home/toluk/Develop/samples/reca_rd/test6.ndx' -sf '/home/toluk/Develop/samples/reca_rd/SelectionList5' -tau 5 -crl 0.10
468 // -s '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -n '/home/toluk/Develop/samples/reca_rd/CorrsTestDomainsNfit.ndx' -sf '/home/toluk/Develop/samples/reca_rd/SelectionListDomainsNFit' -tau 5000 -crl 0.75 -ef_rad 1
469 // -s '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -n '/home/toluk/Develop/samples/reca_rd/TestCa.ndx' -sf '/home/toluk/Develop/samples/reca_rd/SelListCa' -tau 100 -crl 0.75 -ef_rad 1
470 // -s '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -n '/home/toluk/Develop/samples/reca_rd/testCa.ndx' -sf '/home/toluk/Develop/samples/reca_rd/SelListCa' -tau 9000 -crl 0.60 -ef_rad 1
472 // -s '/home/toluk/Develop/samples/CubeTermal/CUBETermalTest1.pdb' -f '/home/toluk/Develop/samples/CubeTermal/CUBETermalTest.pdb' -n '/home/toluk/Develop/samples/CubeTermal/testCa.ndx' -sf '/home/toluk/Develop/samples/CubeTermal/SelListCa' -tau 900 -crl 0.20 -ef_rad 1
473 // cube.000.000.10k.10.3.1stfrm
474 // -s '/home/toluk/Develop/samples/JustCube/cube.000.000.10k.10.3.1stfrm.pdb' -f '/home/toluk/Develop/samples/JustCube/cube.000.000.10k.10.3.pdb' -n '/home/toluk/Develop/samples/JustCube/system.ndx' -sf '/home/toluk/Develop/samples/JustCube/SLsystem' -tau 1000 -crl 0.20 -ef_rad 9000 -out_put OLA
476 SpaceTimeCorr::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
477 TrajectoryAnalysisModuleData *pdata)
479 trajectory.resize(trajectory.size() + 1);
480 trajectory.back().resize(index.size());
481 for (int i = 0; i < index.size(); i++) {
482 trajectory.back()[i] = fr.x[index[i]];
488 SpaceTimeCorr::finishAnalysis(int nframes)
490 std::vector< std::vector< std::vector< float > > > crltns;
491 std::vector< std::vector< std::pair< float, int > > > graph;
492 std::vector< std::vector< int > > sub_graph;
493 std::vector< std::vector< std::pair< int, int > > > sub_graph_rbr, rout_new;
499 std::cout << "\nCorrelation's evaluation - start\n";
501 correlation_evaluation(trajectory, basic_frame, crltns, m, k);
503 make_correlation_matrix_file(crltns, (OutPutName + "_matrix.txt").c_str(), 0);
504 std::cout << "corelation matrix printed\n";
505 make_correlation_pairs_file(crltns, (OutPutName + "_pairs.txt").c_str(), 0);
506 std::cout << "corelation pairs printed\n";
508 std::cout << "Correlation's evaluation - end\n" << "graph evaluation: start\n";
510 graph_calculation(graph, sub_graph, sub_graph_rbr, trajectory, basic_frame, crltns, crl_border, eff_rad, k);
512 std::cout << "graph evaluation: end\n" << "routs evaluation: start\n";
514 graph_back_bone_evaluation(rout_new, index.size(), graph, sub_graph, sub_graph_rbr);
516 std::cout << "routs evaluation: end\n";
518 make_rout_file(crl_border, index, rout_new, (OutPutName + "_routs.txt").c_str());
519 std::cout << "corelation routs printed\n";
520 make_best_corrs_graphics(crltns, rout_new, index, (OutPutName + "_routs_graphics.txt").c_str());
521 std::cout << "corelation routs' pairs' graphics printed\n";
523 std::cout << "Finish Analysis - end\n\n";
527 SpaceTimeCorr::writeOutput()
533 * The main function for the analysis template.
536 main(int argc, char *argv[])
538 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<SpaceTimeCorr>(argc, argv);