-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathramview_no_index.C
132 lines (105 loc) · 3.33 KB
/
ramview_no_index.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
//
// View a region of a RAM file, don't use index, but brute force scan.
//
// Author: Jose Javier Gonzalez Ortiz, 5/7/2017
//
#include <iostream>
#include <string>
#include <TBranch.h>
#include <TTree.h>
#include <TFile.h>
#include <TStopwatch.h>
#include <TTreePerfStats.h>
#include "ramrecord.C"
void ramview_no_index(const char *file, const char *query, bool cache = false, bool perfstats = false,
const char *perfstatsfilename = "perf.root")
{
TStopwatch stopwatch;
stopwatch.Start();
// Open the file and load tree and reader
auto f = TFile::Open(file);
auto t = RAMRecord::GetTree(f);
RAMRecord *r = 0;
if (!cache) {
t->SetCacheSize(0);
}
TTreePerfStats *ps = 0;
if (perfstats) {
ps = new TTreePerfStats("ioperf", t);
}
t->SetBranchAddress("RAMRecord.", &r);
TBranch *b = t->GetBranch("RAMRecord.");
// Parse queried region string
std::string region = query;
int chrDelimiterPos = region.find(":");
TString rname = region.substr(0, chrDelimiterPos);
int rangeDelimiterPos = region.find("-");
UInt_t rangeStart = std::stoi(region.substr(chrDelimiterPos + 1, rangeDelimiterPos - chrDelimiterPos));
UInt_t rangeEnd = std::stoi(region.substr(rangeDelimiterPos + 1, region.size() - rangeDelimiterPos));
// Default values to ensure correctness
int rnameStart = -1;
int posStart = -1;
// Assume RNAME are chunked together
// We look only at the RNAME column
// We can only do this when there are columns
if (b->GetSplitLevel() > 0) {
t->SetBranchStatus("RAMRecord.*", 0);
t->SetBranchStatus("RAMRecord.v_rname", 1);
}
for (int i = 0; i < t->GetEntries(); i++) {
t->GetEntry(i);
if (rname.EqualTo(r->GetRNAME())) {
rnameStart = i;
break;
}
}
// If the RNAME was found
if (rnameStart >= 0) {
// We need to look both at the leftmost position (v_pos)
// as well as the length of sequence (v_lseq)
if (b->GetSplitLevel() > 0) {
t->SetBranchStatus("RAMRecord.v_pos", 1);
t->SetBranchStatus("RAMRecord.v_lseq", 1);
}
for (int i = rnameStart; i < t->GetEntries(); i++) {
t->GetEntry(i);
// If the RNAME region ends
if (!rname.EqualTo(r->GetRNAME())) {
break;
} else {
if (r->GetPOS() + r->GetSEQLEN() > rangeStart) {
// Register first valid position for printing
posStart = i;
break;
}
}
}
// If the position was found
if (posStart >= 0) {
// Enable all fields for printing
if (b->GetSplitLevel() > 0) {
t->SetBranchStatus("RAMRecord.*", 1);
}
for (int i = posStart; i < t->GetEntries(); i++) {
t->GetEntry(i);
// If the RNAME region ends
if (!rname.EqualTo(r->GetRNAME())) {
break;
} else {
// Within the region
if (r->GetPOS() <= rangeEnd) {
r->Print();
} else {
break;
}
}
}
}
}
stopwatch.Print();
if (perfstats) {
ps->SaveAs(perfstatsfilename);
delete ps;
printf("Reading %lld bytes in %d transactions\n", f->GetBytesRead(), f->GetReadCalls());
}
}