- changed the way data stored while work in progress
authorAnatoly Titov <toluk@omrb.pnpi.spb.ru>
Mon, 20 Feb 2017 23:02:59 +0000 (02:02 +0300)
committerAnatoly Titov <toluk@omrb.pnpi.spb.ru>
Mon, 20 Feb 2017 23:02:59 +0000 (02:02 +0300)
- added extra comments while work in progress

src/domains.cpp

index 6cbb6dfa4100a9462cbbdc9f9fb7418fc9e8a8a6..0ffedb2459d7264505fe7081b6c6189fd553f207 100644 (file)
@@ -360,8 +360,6 @@ Domains::finishAnalysis(int nframes)
     int bone = index.size() - domain_min_size + 1;
     graph.resize(bone);
     real **w_rls;
-    rvec **traj, **etalon;
-    snew(traj, bone);
 
     snew(w_rls, bone);
     for (int i = 0; i < bone; i++) {
@@ -375,6 +373,8 @@ Domains::finishAnalysis(int nframes)
         }
     }
 
+    /*
+    rvec **etalon;
     snew(etalon, bone);
     for (int i = 0; i < bone; i++) {
         snew(etalon[i], index.size());
@@ -382,15 +382,26 @@ Domains::finishAnalysis(int nframes)
             copy_rvec(trajectory[etalon_frame][j].as_vec(), etalon[i][j]);
         }
     }
+    */
 
     for (int i = 0; i < bone; i++) {
-        reset_x(index.size(), NULL, index.size(), NULL, etalon[i], w_rls[i]);
-        make_graph(index.size(), etalon[i], graph[i]);
+        rvec *etalon;
+        snew(etalon, index.size());
+        for (int j = 0; j < index.size(); j++) {
+            copy_rvec(trajectory[etalon_frame][j].as_vec(), etalon[j]);
+        }
+        reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[i]);
+        make_graph(index.size(), etalon, graph[i]);
+        sfree(etalon);
     }
 
     std::chrono::time_point<std::chrono::system_clock> start1, end1;
     start1 = std::chrono::system_clock::now();
 
+    /*
+    rvec **traj;
+    snew(traj, bone);
+
     for (int i = 1; i < frames; i++) {
         for (int j = 0; j < bone; j++) {
             sfree(traj[j]);
@@ -425,10 +436,50 @@ Domains::finishAnalysis(int nframes)
             std::cout << "foretold time: " << hours << " hour(s) " << minutes << " minute(s) " << seconds << "second(s)\n";
         }
     }
+    */
+
+    for (int i = 1; i < frames; i++) {
+        #pragma omp parallel
+        {
+            #pragma omp for schedule(dynamic)
+            for (int j = 0; j < bone; j++) {
+                rvec *etalon, *traj;
+                snew(etalon, index.size());
+                for (int k = 0; k < index.size(); k++) {
+                    copy_rvec(trajectory[etalon_frame][k].as_vec(), etalon[k]);
+                }
+                snew(traj, index.size());
+                for (int k = 0; k < index.size(); k++) {
+                    copy_rvec(trajectory[i][k].as_vec(), traj[k]);
+                }
+                reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[j]);
+                reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[j]);
+                do_fit(index.size(), w_rls[j], etalon, traj);
+                update_graph(graph[j], traj, epsi);
+                sfree(etalon);
+                sfree(traj);
+            }
+        }
+        std::cout << i << " / " << frames << " processed\n";
+        if (i % 25 == 0) {
+            end1 = std::chrono::system_clock::now();
+            int seconds = 0, minutes = 0, hours = 0;
+            seconds = std::chrono::duration_cast<std::chrono::seconds>(end1-start1).count();
+            seconds = (float)seconds * (float)( (float)(frames - i) / (float)i );
+            hours = seconds / 3600;
+            seconds -= hours * 3600;
+            minutes = seconds / 60;
+            seconds %= 60;
+            std::cout << "foretold time: " << hours << " hour(s) " << minutes << " minute(s) " << seconds << "second(s)\n";
+        }
+    }
+    std::cout << "finale cheking\n";
     check_domains(delta, frames, graph);
+    std::cout << "finding domains' sizes\n";
     find_domain_sizes(graph, domsizes);
     std::vector< int > a;
     a.resize(0);
+    std::cout << "finding domains\n";
     while (check_domsizes(domsizes, domain_min_size)) {
         domains.push_back(a);
         get_maxsized_domain(domains.back(), graph, domsizes);
@@ -438,17 +489,18 @@ Domains::finishAnalysis(int nframes)
     }
     for (int i = 0; i < bone; i++) {
         sfree(w_rls[i]);
-        sfree(etalon[i]);
-        sfree(traj[i]);
+        //sfree(etalon[i]);
+        //sfree(traj[i]);
     }
     sfree(w_rls);
-    sfree(etalon);
-    sfree(traj);
+    //sfree(etalon);
+    //sfree(traj);
 }
 
 void
 Domains::writeOutput()
 {
+    std::cout << "making output file\n";
     print_domains(domains, index, fnNdx_); // see function for details | numbers from index
     std::cout << "\n END \n";
 }