remade output info
[alexxy/gromacs-testing.git] / src / spacetimecorr.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 <iostream>
44 #include <chrono>
45 #include <omp.h>
46 #include <thread>
47 #include <string>
48 #include <algorithm>
49
50 #include <gromacs/trajectoryanalysis.h>
51 #include <gromacs/pbcutil/pbc.h>
52 #include <gromacs/pbcutil/rmpbc.h>
53 #include <gromacs/utility/smalloc.h>
54 #include <gromacs/math/vectypes.h>
55 #include <gromacs/math/vec.h>
56 #include <gromacs/math/do_fit.h>
57
58 #include "gromacs/analysisdata/analysisdata.h"
59 #include "gromacs/analysisdata/modules/average.h"
60 #include "gromacs/analysisdata/modules/histogram.h"
61 #include "gromacs/analysisdata/modules/plot.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/options/basicoptions.h"
64 #include "gromacs/options/filenameoption.h"
65 #include "gromacs/options/ioptionscontainer.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/selection/selection.h"
68 #include "gromacs/selection/selectionoption.h"
69 #include "gromacs/trajectory/trajectoryframe.h"
70 #include "gromacs/trajectoryanalysis/analysissettings.h"
71 #include "gromacs/utility/arrayref.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/stringutil.h"
74
75 using namespace gmx;
76
77 using gmx::RVec;
78
79 void make_pdb_trajectory(std::vector< std::vector< RVec > > trajectory, const char* file_name)
80 {
81     FILE *file;
82     file = std::fopen(file_name, "w+");
83     for (int i = 1; i < trajectory.size(); i++) {
84         std::fprintf(file, "MODEL%9d\n", i);
85         for (int j = 0; j < trajectory[i].size(); j++) {
86             std::fprintf(file, "ATOM  %5d   CA LYS     1    %8.3f%8.3f%8.3f  1.00  0.00\n", j, trajectory[i][j][0] * 10, trajectory[i][j][1] * 10, trajectory[i][j][2] * 10);
87         }
88         std::fprintf(file, "ENDMDL\n");
89     }
90 }
91 /*
92 1 -  6        Record name     "ATOM  "
93 7 - 11        Integer         Atom serial number.
94 13 - 16        Atom            Atom name.
95 17             Character       Alternate location indicator.
96 18 - 20        Residue name    Residue name.
97 22             Character       Chain identifier.
98 23 - 26        Integer         Residue sequence number.
99 27             AChar           Code for insertion of residues.
100 31 - 38        Real(8.3)       Orthogonal coordinates for X in Angstroms.
101 39 - 46        Real(8.3)       Orthogonal coordinates for Y in Angstroms.
102 47 - 54        Real(8.3)       Orthogonal coordinates for Z in Angstroms.
103 55 - 60        Real(6.2)       Occupancy.
104 61 - 66        Real(6.2)       Temperature factor (Default = 0.0).
105 73 - 76        LString(4)      Segment identifier, left-justified.
106 77 - 78        LString(2)      Element symbol, right-justified.
107 79 - 80        LString(2)      Charge on the atom.
108 */
109
110 void make_correlation_file(std::vector< std::vector< std::vector< double > > > correlations, char* file_name, int start)
111 {
112     FILE *file;
113     file = std::fopen(file_name, "w+");
114     for (int i = start; i < correlations.size(); i++) {
115         if (correlations.size() - start > 1) {
116             std::fprintf(file, "correlation between 'n' and 'n + %d' frames\n", i);
117         }
118         if (i % 50 == 0) {
119             std::cout << "correlation between 'n' and 'n + " << i << "' frames\n";
120         }
121         for (int j = 0; j < correlations[i].size(); j++) {
122             for (int f = 0; f < correlations[i][j].size(); f++) {
123                 std::fprintf(file, "%3.2f ", correlations[i][j][f]);
124             }
125             std::fprintf(file, "\n");
126         }
127     }
128 }
129
130 void make_rout_file(std::vector< int > indx, std::vector< std::vector< int > > rout, const char* file_name)
131 {
132     FILE *file;
133     file = std::fopen(file_name, "w+");
134     for (int i = 0; i < rout.size(); i++) {
135         for (int j = 0; j < rout[i].size(); j++) {
136             std::fprintf(file, "%3d ", indx[rout[i][j]]);
137         }
138         std::fprintf(file, "\n");
139     }
140 }
141
142 int dive_in(std::vector< std::vector< std::pair< int, int > > > graph, int start, int steps, int frames, int depth, int depth_go, std::vector< std::pair< int, int > > path)
143 {
144     int steps_old = 0, steps_new = 0;
145     if (depth > frames) {
146         return steps;
147     }
148     bool flag = true;
149     for (int i = 0; i < path.size(); i++) {
150         if (path[i].first == start) {
151             flag = false;
152         }
153     }
154     if (flag) {
155         path.push_back(std::make_pair(start, depth_go));
156     }
157     for (int i = 0; i < graph[start].size(); i++) {
158         flag = true;
159         for (int j = 0; j < path.size(); j++) {
160             if (path[j].first == graph[start][i].first) {
161                 flag = false;
162             }
163         }
164         if (flag) {
165             steps_new = dive_in(graph, graph[start][i].first, steps + 1, frames, depth + graph[start][i].second, graph[start][i].second, path);
166             if (steps_old < steps_new) {
167                 steps_old = steps_new;
168             }
169         }
170     }
171     return (std::max(steps_old, steps));
172 }
173
174 bool mysortfunc (std::vector<int> a, std::vector<int> b) { return (a.size() > b.size()); }
175
176 /*! \brief
177  * Class used to compute free volume in a simulations box.
178  *
179  * Inherits TrajectoryAnalysisModule and all functions from there.
180  * Does not implement any new functionality.
181  *
182  * \ingroup module_trajectoryanalysis
183  */
184
185
186 class Domains : public TrajectoryAnalysisModule
187 {
188     public:
189
190         Domains();
191         virtual ~Domains();
192
193         //! Set the options and setting
194         virtual void initOptions(IOptionsContainer          *options,
195                                  TrajectoryAnalysisSettings *settings);
196
197         //! First routine called by the analysis framework
198         // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
199         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
200                                   const TopologyInformation        &top);
201
202         virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings   &settings,
203                                          const t_trxframe                   &fr);
204
205         //! Call for each frame of the trajectory
206         // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
207         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
208                                   TrajectoryAnalysisModuleData *pdata);
209
210         //! Last routine called by the analysis framework
211         // virtual void finishAnalysis(t_pbc *pbc);
212         virtual void finishAnalysis(int nframes);
213
214         //! Routine to write output, that is additional over the built-in
215         virtual void writeOutput();
216
217     private:
218
219         std::string                                                 fnNdx_;
220         SelectionList                                               sel_;
221
222         std::vector< std::vector< int > >                           domains;
223         std::vector< std::vector< int > >                           domains_local;
224         std::vector< std::vector< RVec > >                          trajectory;
225         std::vector< std::vector< RVec > >                          frankenstein_trajectory;
226
227
228         std::vector< int >                                          index;
229         std::vector< int >                                          numbers;
230         int                                                         frames              = 0;
231         int                                                         basic_frame         = 0;
232         int                                                         tau                 = 0;
233         double                                                      crl_border          = 0;
234         real                                                        **w_rls;
235
236         int                                                         domains_ePBC;
237         // Copy and assign disallowed by base.
238 };
239
240 Domains::Domains(): TrajectoryAnalysisModule()
241 {
242 }
243
244 Domains::~Domains()
245 {
246 }
247
248 void
249 Domains::initOptions(IOptionsContainer          *options,
250                      TrajectoryAnalysisSettings *settings)
251 {
252     static const char *const desc[] = {
253         "[THISMODULE] to be done"
254     };
255     // Add the descriptive text (program help text) to the options
256     settings->setHelpText(desc);
257     // Add option for output file name
258     options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
259                             .store(&fnNdx_).defaultBasename("domains")
260                             .description("Index file from the domains"));
261     // Add option for tau constant
262     options->addOption(gmx::IntegerOption("tau")
263                             .store(&tau)
264                             .description("number of frames for time travel"));
265     // Add option for crl_border constant
266     options->addOption(DoubleOption("crl")
267                             .store(&crl_border)
268                             .description("make graph based on corrs > constant"));
269     // Add option for selection list
270     options->addOption(SelectionOption("select_domains_and_residue").storeVector(&sel_)
271                            .required().dynamicMask().multiValue()
272                            .description("Domains to form rigid skeleton"));
273     // Control input settings
274     settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
275     settings->setPBC(true);
276 }
277
278 void
279 Domains::initAnalysis(const TrajectoryAnalysisSettings &settings,
280                       const TopologyInformation        &top)
281 {
282     domains_ePBC = top.ePBC();
283 }
284
285 void
286 Domains::initAfterFirstFrame(const TrajectoryAnalysisSettings       &settings,
287                              const t_trxframe                       &fr)
288 {
289     ConstArrayRef< int > atomind  = sel_[0].atomIndices();
290     index.resize(0);
291     for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
292         index.push_back(*ai);
293     }
294
295     trajectory.resize(1);
296     trajectory.back().resize(index.size());
297
298     for (int i = 0; i < index.size(); i++) {
299         trajectory.back()[i] = fr.x[index[i]];
300     }
301
302     domains.resize(sel_.size());
303     for (int i = 1; i < sel_.size(); i++) {
304         atomind  = sel_[i].atomIndices();
305         for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
306             domains[i].push_back(*ai);
307         }
308     }
309
310     for (int i = 0; i < domains.size(); i++) {
311         for (int j = 0; j < domains[i].size(); j++) {
312             int k = 0;
313             while (index[k] != domains[i][j]) {
314                 k++;
315             }
316             domains[i][j] = k;
317         }
318     }
319
320     for (int i = 0; i < domains.size(); i++) {
321         if (domains[i].size() >= 2) {
322             domains_local.push_back(domains[i]);
323             for (int k = 0; k < domains[i].size(); k++) {
324                 for (int j = i + 1; j < domains.size(); j++) {
325                     for (int x = 0; x < domains[j].size(); x++) {
326                         if (domains[j][x] == domains[i][k]) {
327                             domains[j].erase(domains[j].begin() + x);
328                         }
329                     }
330                 }
331             }
332         }
333     }
334
335     domains.resize(0);
336     domains = domains_local;
337
338     std::vector< bool >                 flags;
339     flags.resize(index.size(), true);
340
341     for (int i = 0; i < domains.size(); i++) {
342         for (int j = 0; j < domains[i].size(); j++) {
343             flags[domains[i][j]] = false;
344         }
345     }
346
347     for (int i = 0; i < index.size(); i++) {
348         if (flags[i]) {
349             int a, b, dist = 90000001;
350             rvec temp;
351                 for (int j = 0; j < domains.size(); j++) {
352                     for (int k = 0; k < domains[j].size(); k++) {
353                         rvec_sub(trajectory.back()[i], trajectory.back()[domains[j][k]], temp);
354                         if (norm(temp) <= dist) {
355                             dist = norm(temp);
356                             a = j;
357                             b = k;
358                         }
359                     }
360                 }
361             domains_local[a].push_back(i);
362             flags[i] = false;
363         }
364     }
365
366     frankenstein_trajectory.resize(trajectory.size());
367     frankenstein_trajectory.back() = trajectory.back();
368
369     frames++;
370 }
371
372 void
373 Domains::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
374                       TrajectoryAnalysisModuleData *pdata)
375 {
376     trajectory.resize(trajectory.size() + 1);
377     trajectory.back().resize(index.size());
378
379     for (int i = 0; i < index.size(); i++) {
380         trajectory.back()[i] = fr.x[index[i]];
381     }
382
383     snew(w_rls, domains_local.size() + 1);
384     for (int i = 0; i < domains_local.size(); i++) {
385         snew(w_rls[i], index.size());
386         for (int j = 0; j < index.size(); j++) {
387             w_rls[i][j] = 0;
388         }
389         for (int j = 0; j < domains_local[i].size(); j++) {
390             w_rls[i][domains_local[i][j]] = 1;
391         }
392     }
393     snew(w_rls[domains_local.size()], index.size());
394     for (int i = 0; i < index.size(); i++) {
395         w_rls[domains_local.size()][i] = 1;
396     }
397
398     frankenstein_trajectory.resize(trajectory.size());
399     frankenstein_trajectory.back().resize(index.size());
400
401     for (int i = 0; i < domains_local.size(); i++) {
402         rvec *basic, *traj;
403         snew(basic, index.size());
404         for (int k = 0; k < index.size(); k++) {
405             copy_rvec(trajectory[basic_frame][k].as_vec(), basic[k]);
406         }
407         snew(traj, index.size());
408         for (int k = 0; k < index.size(); k++) {
409             copy_rvec(trajectory.back()[k].as_vec(), traj[k]);
410         }
411         reset_x(index.size(), NULL, index.size(), NULL, basic, w_rls[i]);
412         reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[i]);
413         do_fit(index.size(), w_rls[i], basic, traj);
414
415         for (int j = 0; j < index.size(); j++) {
416             if (w_rls[i][j] == 0) {
417                 copy_rvec(basic[j], traj[j]);
418             }
419         }
420         reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[domains_local.size()]);
421
422         for (int j = 0; j < domains_local[i].size(); j++) {
423             frankenstein_trajectory.back()[domains_local[i][j]] = traj[domains_local[i][j]];
424         }
425
426         sfree(basic);
427         sfree(traj);
428     }
429
430     frames++;
431 }
432
433 //spacetimecorr -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/test.ndx'
434
435 void
436 Domains::finishAnalysis(int nframes)
437 {
438     //make_pdb_trajectory(frankenstein_trajectory, "/home/toluk/Develop/samples/reca_rd/10k_frames.pdb");
439
440     //frankenstein_trajectory = trajectory;
441
442     std::vector< std::vector< std::vector< double > > > crltns;
443     int k = 1000, m = 0;
444     if (tau > 0) {
445         k = tau;
446     }
447     crltns.resize(k + 1);
448     for (int i = 0; i < crltns.size(); i++) {
449         crltns[i].resize(index.size());
450         for (int j = 0; j < crltns[i].size(); j++) {
451             crltns[i][j].resize(index.size(), 0);
452         }
453     }
454
455     std::cout << "\n" ;
456     #pragma omp parallel
457     {
458         #pragma omp for schedule(dynamic)
459         for (int i = m; i <= k; i += 1) {
460             std::cout << " k = " << i << " ";
461
462             rvec *medx, *medy, temp_zero;
463             snew(medx, index.size());
464             snew(medy, index.size());
465             temp_zero[0] = 0;
466             temp_zero[1] = 0;
467             temp_zero[2] = 0;
468             for (int i = 0; i < index.size(); i++) {
469                 copy_rvec(temp_zero, medx[i]);
470                 copy_rvec(temp_zero, medy[i]);
471             }
472
473             for (int j = 0; j < frankenstein_trajectory.size() - i; j++) {
474                 for (int f = 0; f < index.size(); f++) {
475                     rvec temp;
476
477                     rvec_sub(frankenstein_trajectory[basic_frame][f], frankenstein_trajectory[j][f], temp);
478                     rvec_inc(medx[f], temp);
479                     rvec_sub(frankenstein_trajectory[basic_frame][f], frankenstein_trajectory[j + i][f], temp);
480                     rvec_inc(medy[f], temp);
481                 }
482             }
483
484             for (int j = 0; j < index.size(); j++) {
485                 rvec temp;
486                 real temp2 = frankenstein_trajectory.size();
487                 temp2 -= i;
488                 temp2 = 1 / temp2;
489                 copy_rvec(medx[j], temp);
490                 svmul(temp2, temp, medx[j]);
491                 copy_rvec(medy[j], temp);
492                 svmul(temp2, temp, medy[j]);
493             }
494
495             std::vector< std::vector< double > > a, b, c;
496             a.resize(index.size());
497             b.resize(index.size());
498             c.resize(index.size());
499             for (int j = 0; j < index.size(); j++) {
500                 a[j].resize(index.size(), 0);
501                 b[j].resize(index.size(), 0);
502                 c[j].resize(index.size(), 0);
503             }
504             for (int j = 0; j < frankenstein_trajectory.size() - i; j++) {
505                 for (int f1 = 0; f1 < index.size(); f1++) {
506                     for (int f2 = f1; f2 < index.size(); f2++) {
507                         rvec temp1, temp2;
508                         rvec_sub(frankenstein_trajectory[basic_frame][f1], frankenstein_trajectory[j][f1], temp1);
509                         rvec_dec(temp1, medx[f1]);
510                         rvec_sub(frankenstein_trajectory[basic_frame][f2], frankenstein_trajectory[j + i][f2], temp2);
511                         rvec_dec(temp2, medy[f2]);
512                         a[f1][f2] += (temp1[0] * temp2[0] + temp1[1] * temp2[1] + temp1[2] * temp2[2]);
513                         b[f1][f2] += (temp1[0] * temp1[0] + temp1[1] * temp1[1] + temp1[2] * temp1[2]);
514                         c[f1][f2] += (temp2[0] * temp2[0] + temp2[1] * temp2[1] + temp2[2] * temp2[2]);
515                     }
516                 }
517             }
518             for (int j = 0; j < index.size(); j++) {
519                 for (int f = j; f < index.size(); f++) {
520                     crltns[i][j][f] = a[j][f] / (std::sqrt(b[j][f] * c[j][f]));
521                 }
522             }
523             sfree(medx);
524             sfree(medy);
525         }
526     }
527     std::cout << "\n" ;
528
529     //make_correlation_file(crltns, "corrs.txt", m);
530
531     for (int i1 = 0; i1 < 100; i1++) {
532         int i5 = 0;
533         for (int i2 = 1; i2 < crltns.size(); i2++) {
534             for (int i3 = 0; i3 < crltns[i2].size(); i3++) {
535                 for (int i4 = 0; i4 < crltns[i2][i3].size(); i4++) {
536                     if (crltns[i2][i3][i4] >= (double)i1 / 100) {
537                         i5++;
538                     }
539                 }
540             }
541         }
542         std::cout << i1 << " - " << i5 << "   ";
543         if (i % 10 == 0) {
544             std::cout << "\n" ;
545         }
546     }
547
548     std::cout << "\n\n";
549
550     std::vector< std::vector< int > > graph;
551     graph.resize(index.size());
552     for (int i = 1; i <= k; i++) {
553         for (int j = 0; j < index.size(); j++) {
554             for (int f = 0; f < index.size(); f++) {
555                 if (std::abs(crltns[i][j][f]) >= crl_border) {
556                     bool local_flag = true;
557                     for (int f2 = 0; f2 < graph[j].size(); f2++) {
558                         if (graph[j][f2] == f) {
559                             local_flag = false;
560                         }
561                     }
562                     if (local_flag) {
563                         graph[j].push_back(f);
564                     }
565                 }
566             }
567         }
568     }
569
570     /*for (int i = 0; i < index.size(); i++) {
571         if (graph[i].size() > 0) {
572             std::cout << i << " - " << graph[i].size() << "\n";
573         }
574     }*/
575
576     std::cout << "\n";
577     std::vector< std::vector < int > > rout_old, rout_new;
578     bool flag = true;
579     rout_new.resize(0);
580     rout_old.resize(index.size());
581     for (int i = 0; i < index.size(); i++) {
582         rout_old[i].resize(1, i);
583     }
584
585     std::cout << "building routs\n\n" ;
586
587     while (flag) {
588         flag = false;
589         std::vector < int > flag3;
590         for (long int i = 0; i < rout_old.size(); i++) {
591             if (graph[rout_old[i].back()].size() == 0 && rout_old[i].size() > 2) {
592                 rout_new.push_back(rout_old[i]);
593             } else {
594                 bool flag222 = true;
595                 for (long int j = 0; j < graph[rout_old[i].back()].size(); j++) {
596                     bool flag2 = true, flag22 = true;
597                     flag3 = rout_old[i];
598                     for (int f = 0; f < flag3.size(); f++) {
599                         if (flag3[f] == graph[rout_old[i].back()][j]) {
600                             flag22 = false;
601                         }
602                     }
603                     if (flag22) {
604                         flag3.push_back(graph[rout_old[i].back()][j]);
605                         for (int f = 0; f < rout_new.size(); f++) {
606                             if (rout_new[f] == flag3) {
607                                 flag2 = false;
608                             }
609                         }
610                         if (flag2) {
611                             rout_new.push_back(flag3);
612                             flag = true;
613                             flag222 = false;
614
615                         }
616                     }
617                 }
618                 if (flag222 && rout_old[i].size() > 2) {
619                     rout_new.push_back(rout_old[i]);
620                 }
621             }
622         }
623         if (rout_new.size() != 0) {
624             std::cout << rout_old.size() << " -> " << rout_new.size() << "\n";
625             rout_old = rout_new;
626             rout_new.resize(0);
627         }
628     }
629
630     std::sort(rout_old.begin(), rout_old.end(), mysortfunc);
631
632     std::cout << "\nfinished routs\n\n" ;
633
634     /*graph.resize(0);
635     graph.resize(index.size());
636     rout_new.resize(0);
637
638     for (int i = 0; i < rout_old.size(); i++) {
639         int lflag2 = false;
640         for (int j = 0; j < rout_old[i].size() - 1; j++) {
641             int lflag1 = true;
642             for (int f = 0; f < graph[rout_old[i][j]].size(); f++) {
643                 if (graph[rout_old[i][j]][f] == rout_old[i][j + 1]) {
644                     lflag1 = false;
645                 }
646             }
647             if (lflag1) {
648                 graph[rout_old[i][j]].push_back(rout_old[i][j + 1]);
649                 lflag2 = true;
650             }
651         }
652         if (lflag2) {
653             rout_new.push_back(rout_old[i]);
654         }
655     }*/
656
657
658     /*std::cout << "\n test old \n";
659     for (int i = 0; i < rout_old.size(); i++) {
660         for (int j = 0; j < rout_old[i].size(); j++) {
661             std::cout << index[rout_old[i][j]] << " ";
662         }
663         std::cout << "\n";
664     }
665
666     std::cout << "\n test new \n";*/
667
668     rout_new.resize(0);
669     for (int i = rout_old.size() - 1; i >= 0; i--) {
670         bool lflag1 = true;
671         for (int j = 0; j < i; j++) {
672             if (rout_old[i].size() < rout_old[j].size()) {
673                 int la = 0;
674                 for (int f = 0; f < rout_old[j].size(); f++) {
675                     if (rout_old[j][f] == rout_old[i][la]) {
676                         la++;
677                     }
678                 }
679                 if (la == rout_old[i].size()) {
680                     lflag1 = false;
681                 }
682             }
683         }
684         if (lflag1) {
685             rout_new.push_back(rout_old[i]);
686         }
687     }
688
689     //make_rout_file(index, rout_new, "routs.txt");
690
691     for (int i = 0; i < rout_new.size(); i++) {
692         for (int j = 0; j < rout_new[i].size(); j++) {
693             std::cout << index[rout_new[i][j]] << " ";
694         }
695         std::cout << "\n";
696     }
697
698     std::cout << "\n\n\n job done \n\n\n";
699 }
700
701 void
702 Domains::writeOutput()
703 {
704
705 }
706
707
708 /*! \brief
709  * The main function for the analysis template.
710  */
711 int
712 main(int argc, char *argv[])
713 {
714     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<Domains>(argc, argv);
715 }