GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
bufferclean.cxx
Go to the documentation of this file.
1// Author: Ryan Dunlop
2
3// g++ offsetfind.cxx `root-config --cflags --libs` -I${GRSISYS}/include -L${GRSISYS}/libraries -lMidasFormat
4// -lXMLParser -ooffsetfind
5#include "TMidasFile.h"
6#include "TMidasEvent.h"
7#include "TStopwatch.h"
8#include "TSystem.h"
9#include <map>
10#include <iostream>
11#include <fstream>
12#include <chrono>
13#include <utility>
14
15UInt_t chanId_threshold = 100;
16
17Bool_t CheckEvent(const std::shared_ptr<TMidasEvent>& evt)
18{
19 // This function does not work if a Midas event contains multiple fragments
20 static std::map<Int_t, Bool_t> triggermap; // Map of Digitizer vs have we had a triggerId < threshold yet?
21 // First parse the Midas Event.
22 evt->SetBankList();
23
24 // Need to put something in that says "if not a Griffin fragment (ie epics) return true"
25 void* ptr = nullptr;
26 int banksize = evt->LocateBank(nullptr, "GRF2", &ptr);
27 int bank = 2;
28
29 if(banksize == 0) {
30 banksize = evt->LocateBank(nullptr, "GRF1", &ptr);
31 bank = 1;
32 }
33 uint32_t type = 0xffffffff;
34 uint32_t value = 0xffffffff;
35
36 UInt_t chanadd = 0;
37 UInt_t trigId = 0;
38 // UInt_t dettype = 0;
39 // Int_t dignum = -1;
40
41 for(int x = 0; x < banksize; x++) {
42 value = *(reinterpret_cast<int*>(ptr) + x);
43 type = value & 0xf0000000;
44
45 switch(type) {
46 case 0x80000000:
47 switch(bank) {
48 case 1: // header format from before May 2015 experiments
49 // Sets:
50 // The number of filters
51 // The Data Type
52 // Number of Pileups
53 // Channel Address
54 // Detector Type
55 if((value & 0xf0000000) != 0x80000000) {
56 return false;
57 }
58 chanadd = (value & 0x0003fff0) >> 4;
59 // dettype = (value &0x0000000f);
60
61 // if(frag-DetectorType==2)
62 // frag->ChannelAddress += 0x8000;
63 break;
64 case 2:
65 // Sets:
66 // The number of filters
67 // The Data Type
68 // Number of Pileups
69 // Channel Address
70 // Detector Type
71 if((value & 0xf0000000) != 0x80000000) {
72 return false;
73 }
74 chanadd = (value & 0x000ffff0) >> 4;
75 // dettype = (value &0x0000000f);
76
77 // if(frag-DetectorType==2)
78 // frag->ChannelAddress += 0x8000;
79 break;
80 default: printf("This bank not yet defined.\n"); break;
81 };
82 // dignum = chanadd&0x0000ff00;
83 break;
84 case 0x90000000: trigId = value & 0x0fffffff;
85 };
86 }
87 if(triggermap.find(chanadd) == triggermap.end()) {
88 triggermap[chanadd] = false; // initialize the new digitizer number to false.
89 }
90 // Check to make sure we aren't getting any triggerId's = 0. I think these are corrupt events RD.
91 if(trigId == 0) {
92 return false;
93 }
94
95 // Check against map of trigger Id's to see if we have hit the elusive < threshold mark yet.
96 if(triggermap.find(chanadd)->second) {
97 return true;
98 }
99 if(trigId < chanId_threshold) {
100 triggermap.find(chanadd)->second = true;
101 return true;
102 }
103 return false;
104}
105
106void Write(const std::shared_ptr<TMidasEvent>& evt, TMidasFile* outfile)
107{
108 outfile->FillBuffer(evt);
109}
110
111int main(int argc, char** argv)
112{
113 if(argc < 2) {
114 printf("Usage: ./bufferclear <input.mid> <output.mid>\n");
115 return 1;
116 }
117
118 if(argv[1] == argv[2]) {
119 printf("input.mid and output.mid must have different names\n");
120 return 1;
121 }
122
123 auto* file = new TMidasFile;
124 // TMidasFile *outfile = new TMidasFile;
125 file->Open(argv[1]);
126 int runnumber = file->GetRunNumber();
127 int subrunnumber = file->GetSubRunNumber();
128 if(argc < 3) {
129 file->OutOpen(Form("cleanrun%05d_%03d.mid", runnumber, subrunnumber));
130 } else {
131 file->OutOpen(argv[2]);
132 }
133
134 std::ifstream in(file->GetFilename(), std::ifstream::in | std::ifstream::binary);
135 in.seekg(0, std::ifstream::end);
136 int64_t filesize = in.tellg();
137 in.close();
138
139 int bytes = 0;
140 int64_t bytesread = 0;
141
142 TStopwatch w;
143 w.Start();
144
145 UInt_t num_bad_evt = 0;
146 UInt_t num_evt = 0;
147 std::shared_ptr<TMidasEvent> event = std::make_shared<TMidasEvent>(); // need to new for each event
148
149 while(true) {
150 bytes = file->Read(event);
151 if(bytes == 0) {
152 printf(DMAGENTA "\tfile: %s ended on %s" RESET_COLOR "\n", file->GetFilename(), file->GetLastError());
153 if(file->GetLastErrno() == -1) { // try to read some more...
154 continue;
155 }
156 break;
157 }
158 bytesread += bytes;
159
160 switch(event->GetEventId()) {
161 case 0x8000:
162 printf("start of run\n");
163 Write(event, file);
164 printf(DGREEN);
165 event->Print();
166 printf(RESET_COLOR);
167 break;
168 case 0x8001:
169 printf("end of run\n");
170 printf(DRED);
171 event->Print();
172 printf(RESET_COLOR);
173 Write(event, file);
174 break;
175 case 0x0001: // This is a GRIFFIN digitizer event
176 if(CheckEvent(event)) {
177 Write(event, file);
178 } else {
179 num_bad_evt++;
180 }
181 num_evt++;
182 break;
183 default: // Probably epics
184 Write(event, file);
185 break;
186 };
187
188 if(num_evt % 5000 == 0) {
189 gSystem->ProcessEvents();
190 std::streamsize precision = std::cout.precision();
191 std::cout.precision(2);
192 std::cout << HIDE_CURSOR << " bad events " << num_bad_evt << "/" << num_evt << " have processed " << static_cast<double>(bytesread) / 1000000. << "MB/" << static_cast<double>(filesize) / 1000000. << " MB => " << std::setprecision(1) << static_cast<double>(bytesread) / 1000000 / w.RealTime() << " MB/s " << SHOW_CURSOR << "\r";
193 std::cout.precision(precision);
194 w.Continue();
195 }
196 }
197 file->WriteBuffer();
198
199 printf("\n");
200
201 file->Close();
202 file->OutClose();
203 delete file;
204 // delete outfile;
205}
206// void WriteEventToFile(TMidasFile*,std::shared_ptr<TMidasEvent> Option_t);
#define SHOW_CURSOR
Definition Globals.h:33
#define DRED
Definition Globals.h:18
#define DGREEN
Definition Globals.h:17
#define HIDE_CURSOR
Definition Globals.h:32
#define RESET_COLOR
Definition Globals.h:5
#define DMAGENTA
Definition Globals.h:20
void Write(const std::shared_ptr< TMidasEvent > &evt, TMidasFile *outfile)
UInt_t chanId_threshold
int main(int argc, char **argv)
Bool_t CheckEvent(const std::shared_ptr< TMidasEvent > &evt)
Reader for MIDAS .mid files.
Definition TMidasFile.h:32
bool Open(const char *filename) override
Open input file.
void FillBuffer(const std::shared_ptr< TMidasEvent > &midasEvent, Option_t *opt="")