StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
emchist.C
1 #include <stdio.h>
2 #include <unistd.h>
3 #include <getopt.h>
4 #include <sys/types.h>
5 #include <stdlib.h>
6 #include <string.h>
7 
8 #include <rtsLog.h> // for my LOG() call
9 #include <rtsSystems.h>
10 
11 // this needs to be always included
12 #include <DAQ_READER/daqReader.h>
13 #include <DAQ_READER/daq_dta.h>
14 
15 #include <trgDataDefs.h>
16 
17 // only the detectors we will use need to be included
18 // for their structure definitions...
19 #include <DAQ_BSMD/daq_bsmd.h>
20 #include <DAQ_BTOW/daq_btow.h>
21 #include <DAQ_EMC/daq_emc.h>
22 #include <DAQ_ESMD/daq_esmd.h>
23 #include <DAQ_ETOW/daq_etow.h>
24 #include <DAQ_TRG/daq_trg.h>
25 
26 static int bsmd_doer(daqReader *rdr, const char *do_print) ;
27 static int esmd_doer(daqReader *rdr, const char *do_print) ;
28 static int btow_doer(daqReader *rdr, const char *do_print) ;
29 static int etow_doer(daqReader *rdr, const char *do_print) ;
30 static int trg_doer(daqReader *rdr, const char *do_print) ;
31 
32 #include "TH1.h"
33 #include "TH2.h"
34 #include "TH3.h"
35 #include "TFile.h"
36 
37 void histinit();
38 void histwrite();
39 void histdelete();
40 
41 using namespace std;
42 
43 TH2F* etow[ETOW_MAXFEE];
44 TH2F* btow[BTOW_MAXFEE];
45 TH2F* bsmd[BSMD_FIBERS];
46 TH2F* bsmd_zs[BSMD_FIBERS];
47 TH2F* esmd[ESMD_MAXFEE];
48 
49 TFile* outfile;
50 
51 static int good ;
52 
53 int main(int argc, char *argv[])
54 {
55  outfile=0;
56  char dir[FILENAME_MAX];
57  getcwd(dir,FILENAME_MAX);
58  printf("%s\n",dir);
59 
60  extern char *optarg ;
61  extern int optind ;
62  int c ;
63  const char *print_det = "" ;
64  char _mountpoint[256];
65  char *mountpoint = NULL;
66 
67  rtsLogOutput(RTS_LOG_STDERR) ;
68  rtsLogLevel((char* )WARN) ;
69 
70  while((c = getopt(argc, argv, "D:d:m:h")) != EOF) {
71  switch(c) {
72  case 'd' :
73  rtsLogLevel(optarg) ;
74  break ;
75  case 'D' :
76  print_det = optarg ;
77  break ;
78  case 'm' :
79  mountpoint = _mountpoint;
80  strcpy(mountpoint, optarg);
81  break;
82 
83  default :
84  break ;
85  }
86  }
87 
88  class daqReader *evp ; // tha main guy
89  evp = new daqReader(argv[optind]) ; // create it with the filename argument..
90  if(mountpoint) {
91  evp->setEvpDisk(mountpoint);
92  }
93 
94  good=0;
95  int bad=0;
96  int lastrun = -1;
97  int byteTot = 0;
98 
99  // "DAQ trigger ID" for triggers of interest (usually MB)
100  int trigID[2]={2048,128}; //BBCMB=2048 and JP0=128 //200 GeV Run 12 emc-check
101  int eventTrig[2]={0,0}; //counter for events satisfying triggerxs
102 
103  for(;;) {
104  char *ret = evp->get(0,EVP_TYPE_ANY);
105 
106  if(ret) {
107  if(evp->status) {
108  LOG(ERR,"evp status is non-null [0x08X, %d dec]",evp->status,evp->status) ;
109  continue ;
110  }
111  good++;
112  if(good%1000 == 0) LOG(INFO,"%d good events processed",good);
113  }
114  else { // something other than an event, check what.
115  switch(evp->status) {
116  case EVP_STAT_OK: // just a burp!
117  continue;
118  case EVP_STAT_EOR:
119  LOG(DBG, "End of Run/File");
120  if(evp->IsEvp()) { // but event pool, keep trying...
121  LOG(DBG, "Wait a second...");
122  sleep(1);
123  continue;
124  }
125  break; // file, we're done...
126  case EVP_STAT_EVT:
127  bad++;
128  LOG(WARN, "Problem getting event - skipping [good %d, bad %d]",good,bad);
129  sleep(1);
130  continue;
131  case EVP_STAT_CRIT:
132  LOG(CRIT,"evp->status CRITICAL (?)") ;
133  return -1;
134  }
135  }
136 
137 
138  if(evp->status == EVP_STAT_EOR) {
139  LOG(INFO,"End of File reached...%i good events processed",good) ;
140  LOG(INFO,"bit %i = %i events processed",trigID[0],eventTrig[0]) ;
141  LOG(INFO,"bit %i = %i events processed",trigID[1],eventTrig[1]) ;
142  break ; // of the for() loop...
143  }
144 
145  if(evp->run == 0)continue;
146  if(lastrun != evp->run){
147  printf("new run started\n");
148  lastrun = evp->run;
149  if(outfile){
150  histwrite();
151  outfile->Close();
152  delete outfile;
153  outfile = 0;
154  printf("old file closed, deleted\n");
155  histdelete();
156  printf("hists cleared deleted\n");
157  }
158  if(!outfile){
159  printf("need to create new file\n");
160  string infile(argv[optind]);
161  char name[100];
162  string daq(".daq");
163  if(infile.find(daq)!=string::npos){
164  string phys("st_physics");
165  string w("st_W");
166  size_t found = infile.find(phys);
167  if(found==string::npos)found=infile.find(w);
168  string base = infile.substr(found);
169  size_t fd = base.find(daq);
170  base.replace(fd,4,"");
171  //sprintf(name,"%s/%i.ushist.root",dir,evp->run);
172  sprintf(name,"%s/%s.ushist.root",dir,base.c_str());
173  }else{
174  sprintf(name,"%s/%i.ushist.root",dir,evp->run);
175  }
176  printf("%s\n",name);
177  outfile = new TFile(name,"RECREATE");
178  histinit();
179  printf("new hists and file created\n");
180  }
181 
182  }
183 
184  daq_dta *dd ; // generic data pointer; reused all the time
185 
186  bool goodTrig[2]={0,0};
187 
188  //try trigger filter (from daqFileChopper)
189  unsigned int bits = evp->daqbits;
190  //try daq trigger bits
191  //cout<<"event ="<<good<<", bits ="<<bits;
192  if(bits & trigID[0]) goodTrig[0]=1;
193  if(bits & trigID[1]) goodTrig[1]=1;
194 
195  if(goodTrig[1]) eventTrig[1]++;
196  if(!goodTrig[0]) continue;
197  eventTrig[0]++;
198 
199  /***************** EMCs ************************/
200 
201  if(btow_doer(evp, print_det)) ;//LOG(INFO,"BTOW found") ;
202 
203  //if(bsmd_doer(evp,print_det)) ;//LOG(INFO,"BSMD found (any bank)") ;
204 
205  if(etow_doer(evp, print_det)) ;//LOG(INFO,"ETOW found") ;
206 
207  if(esmd_doer(evp, print_det)) ;//LOG(INFO,"ESMD found") ;
208 
209 
210  }
211  if(outfile){
212  histwrite();
213  outfile->Close();
214  delete outfile;
215  outfile = 0;
216  printf("old file closed, deleted\n");
217  histdelete();
218  printf("hists cleared deleted\n");
219  }
220  return 0 ;
221 }
222 
223 void histinit(){
224 
225  for(int i = 0; i < ETOW_MAXFEE; i++){
226  char name[100];
227  char title[100];
228  sprintf(name,"ETOW_%i",i+1);
229  sprintf(title,"ETOW FEE %i",i+1);
230  etow[i] = new TH2F(name,title,ETOW_DATSIZE,-0.5,ETOW_DATSIZE-0.5,4096,-0.5,4095.5);
231  }
232  for(int i = 0; i < BTOW_MAXFEE; i++){
233  char name[100];
234  char title[100];
235  sprintf(name,"BTOW_%i",i+1);
236  sprintf(title,"BTOW FEE %i",i+1);
237  btow[i] = new TH2F(name,title,BTOW_DATSIZE,-0.5,BTOW_DATSIZE-0.5,4096,-0.5,4095.5);
238  }
239  for(int i = 0; i < ESMD_MAXFEE; i++){
240  char name[100];
241  char title[100];
242  sprintf(name,"ESMD_%i",i+1);
243  sprintf(title,"ESMD FEE %i",i+1);
244  esmd[i] = new TH2F(name,title,ESMD_DATSIZE,-0.5,ESMD_DATSIZE-0.5,4096,-0.5,4095.5);
245  }
246 
247 #if 0 //not using BSMD at the moment
248  for(int i = 0; i < BSMD_FIBERS; i++){
249  char name[100];
250  char title[100];
251  sprintf(name,"BSMD_%i",i+1);
252  sprintf(title,"BSMD FIBER %i",i+1);
253  bsmd[i] = new TH2F(name,title,BSMD_DATSIZE,-0.5,BSMD_DATSIZE-0.5,1024,-0.5,1023.5);
254  }
255  for(int i = 0; i < BSMD_FIBERS; i++){
256  char name[100];
257  char title[100];
258  sprintf(name,"BSMD_zs_%i",i+1);
259  sprintf(title,"BSMD FIBER %i",i+1);
260  bsmd_zs[i] = new TH2F(name,title,BSMD_DATSIZE,-0.5,BSMD_DATSIZE-0.5,1024,-0.5,1023.5);
261  }
262 #endif
263 
264 }
265 
266 void histwrite()
267 {
268  outfile->Write();
269 }
270 void histdelete()
271 {
272  for(int i = 0; i < ETOW_MAXFEE; i++){
273  etow[i] = 0;
274  }
275  for(int i = 0; i < BTOW_MAXFEE; i++){
276  btow[i] = 0;
277  }
278  for(int i = 0; i < ESMD_MAXFEE; i++){
279  esmd[i] = 0;
280  }
281 #if 0 //not using BSMD at the moment
282  for(int i = 0; i < BSMD_FIBERS; i++){
283  bsmd[i] = 0;
284  }
285 #endif
286 
287 }
288 static int trg_doer(daqReader *rdr, const char *do_print)
289 {
290  int found = 0 ;
291  daq_dta *dd ;
292 
293  if(strcasestr(do_print,"trg")) ;
294  else do_print = 0 ;
295 
296  // If i want decoded data ala old evpReader I call "legacy" ;
297 
298  dd = rdr->det("trg")->get("legacy") ;
299  if(dd) {
300  if(dd->iterate()) {
301  trg_t *trg_p = (trg_t *) dd->Void ;
302 
303  if(do_print) { // print something...
304  printf("Trigger: daqbits 0x%08X, trg_word 0x%04X\n",trg_p->daqbits,trg_p->trg_word) ;
305  }
306 
307  }
308  }
309 
310 
311  // if you need the void black box of untouched trigger data:
312  dd = rdr->det("trg")->get("raw") ;
313  if(dd) {
314  if(dd->iterate()) {
315  found = 1 ;
316 
317 
318  u_char *trg_raw = dd->Byte;
319 
320  if(do_print) { // I have no clue but let me print first few words...
321 
322 
323  // simple start of trig desc; the way it should be...
324  struct simple_desc {
325  short len ;
326  char evt_desc ;
327  char ver ;
328  } *desc ;
329 
330  desc = (simple_desc *) trg_raw ;
331 
332 
333  printf("Trigger: raw bank has %d bytes: ver 0x%02X, desc %d, len %d\n",dd->ncontent,desc->ver, desc->evt_desc, desc->len) ;
334 
335  // dump first 10 ints...
336  u_int *i32 = (u_int *) trg_raw ;
337  for(int i=0;i<10;i++) {
338  printf("Trigger: word %d: 0x%08X\n",i,i32[i]) ;
339  }
340  }
341  }
342  }
343 
344  return found ;
345 }
346 
347 #if 0 //not using BSMD at the moment
348 static int bsmd_doer(daqReader *rdr, const char *do_print)
349 {
350  int found = 0 ;
351  daq_dta *dd ;
352 
353  if(strcasestr(do_print,"bsmd")) ; // leave as is...
354  else do_print = 0 ;
355 
356 
357  // do I see the non-zero-suppressed bank? let's do this by fiber...
358  for(int f=1;f<=12;f++) {
359  dd = rdr->det("bsmd")->get("adc_non_zs",0,f) ; // sector is ignored (=0)
360  if(dd) {
361  while(dd->iterate()) {
362  found = 1 ;
363 
364  bsmd_t *d = (bsmd_t *) dd->Void ;
365 
366  if(do_print) printf("BSMD non-ZS: fiber %2d, capacitor %d:\n",dd->rdo,d->cap) ;
367 
368 
369  for(int i=0;i<BSMD_DATSIZE;i++) {
370  short cap = (short)d->cap;
371  if(cap!=124 && cap!=125)bsmd[f-1]->Fill(i,d->adc[i]);
372  if(do_print) printf(" %4d = %4d\n",i,d->adc[i]) ;
373  }
374  }
375  }
376  }
377 
378  // do I see the zero suppressed bank?
379  for(int f=1;f<=12;f++) {
380  dd = rdr->det("bsmd")->get("adc",0,f) ;
381  int ngood = 0;
382  if(dd) {
383  while(dd->iterate()) {
384  found = 1 ;
385 
386  bsmd_t *d = (bsmd_t *) dd->Void ;
387 
388  if(do_print) printf("BSMD ZS: fiber %2d, capacitor %d:\n",dd->rdo,d->cap) ;
389 
390  for(int i=0;i<BSMD_DATSIZE;i++) {
391  short cap = (short)d->cap;
392  if(cap!=124 && cap!=125 && d->adc[i]!=0)bsmd_zs[f-1]->Fill(i,d->adc[i]);
393  if(d->adc[i]!= 0)ngood++;
394  // since this is zero-suppressed, I'll skip zeros...
395  if(do_print) if(d->adc[i]) printf(" %4d = %4d\n",i,d->adc[i]) ;
396  }
397  }
398  }
399  if(ngood == 0){
400  printf("bsmd rdo %i had no hits in %i event %i\n",f,rdr->run,rdr->event_number);
401  }
402  }
403 
404 
405  return found ;
406 }
407 #endif
408 
409 static int esmd_doer(daqReader *rdr, const char *do_print)
410 {
411  int found = 0 ;
412  daq_dta *dd ;
413 
414  if(strcasestr(do_print,"esmd")) ; // leave as is...
415  else do_print = 0 ;
416 
417 
418  dd = rdr->det("esmd")->get("adc") ;
419  if(dd) {
420  while(dd->iterate()) {
421  found = 1 ;
422 
423  esmd_t *d = (esmd_t *) dd->Void ;
424 
425  for(int i=0;i<ESMD_MAXFEE;i++) {
426  for(int j=0;j<ESMD_PRESIZE;j++) {
427  if(do_print) printf("ESMD: fee %2d: preamble %d: 0x%04X [%d dec]\n",i,j,d->preamble[i][j], d->preamble[i][j]) ;
428  }
429  for(int j=0;j<ESMD_DATSIZE;j++) {
430  esmd[i]->Fill(j,d->adc[i][j]);
431  if(do_print) printf("ESMD: fee %2d: data %d: 0x%04X [%d dec]\n",i,j,d->adc[i][j], d->adc[i][j]) ;
432  }
433 
434  }
435  }
436  }
437 
438  return found ;
439 }
440 
441 
442 static int etow_doer(daqReader *rdr, const char *do_print)
443 {
444  int found = 0 ;
445  daq_dta *dd ;
446 
447  if(strcasestr(do_print,"etow")) ; // leave as is...
448  else do_print = 0 ;
449 
450 
451  dd = rdr->det("etow")->get("adc") ;
452  if(dd) {
453  while(dd->iterate()) {
454  found = 1 ;
455 
456  etow_t *d = (etow_t *) dd->Void ;
457 
458  for(int i=0;i<ETOW_MAXFEE;i++) {
459  for(int j=0;j<ETOW_PRESIZE;j++) {
460  if(do_print) printf("ETOW: fee %2d: preamble %d: 0x%04X [%d dec]\n",i,j,d->preamble[i][j], d->preamble[i][j]) ;
461  }
462  for(int j=0;j<ETOW_DATSIZE;j++) {
463  etow[i]->Fill(j,d->adc[i][j]);
464  if(do_print) printf("ETOW: fee %2d: data %d: 0x%04X [%d dec]\n",i,j,d->adc[i][j], d->adc[i][j]) ;
465  }
466 
467  }
468  }
469  }
470 
471  return found ;
472 }
473 
474 static int btow_doer(daqReader *rdr, const char *do_print)
475 {
476  int found = 0 ;
477  daq_dta *dd ;
478 
479  if(strcasestr(do_print,"btow")) ; // leave as is...
480  else do_print = 0 ;
481 
482 
483  dd = rdr->det("btow")->get("adc") ;
484  if(dd) {
485  while(dd->iterate()) {
486  found = 1 ;
487 
488  btow_t *d = (btow_t *) dd->Void ;
489 
490  for(int i=0;i<BTOW_MAXFEE;i++) {
491  for(int j=0;j<BTOW_PRESIZE;j++) {
492  if(do_print) printf("BTOW: fee %2d: preamble %d: 0x%04X [%d dec]\n",i,j,d->preamble[i][j], d->preamble[i][j]) ;
493  }
494  for(int j=0;j<BTOW_DATSIZE;j++) {
495  btow[i]->Fill(j,d->adc[i][j]);
496  if(do_print) printf("BTOW: fee %2d: data %d: 0x%04X [%d dec]\n",i,j,d->adc[i][j], d->adc[i][j]) ;
497  }
498 
499  }
500  }
501  }
502 
503  return found ;
504 }
505 
Definition: daq_etow.h:9
Definition: daq_btow.h:9
Definition: daq_trg.h:9