StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StRFEmcTrigMaker.cxx
1 // *-- Author : Renee Fatemi
2 //
3 // $Id: StRFEmcTrigMaker.cxx,v 1.8 2011/04/11 19:35:52 fisyak Exp $
4 
5 #include "StRFEmcTrigMaker.h"
6 #include "StChain.h"
7 #include "StMuDSTMaker/COMMON/StMuEvent.h"
8 #include "StMuDSTMaker/COMMON/StMuDst.h"
9 #include "StMuDSTMaker/COMMON/StMuDstMaker.h"
10 #include "StMuDSTMaker/COMMON/StMuEmcCollection.h"
11 #include "StMuDSTMaker/COMMON/StMuTrack.h"
12 #include "StMuDSTMaker/COMMON/StMuEmcPoint.h"
13 #include "StMuDSTMaker/COMMON/StMuTrack.h"
14 #include "StEmcUtil/geometry/StEmcGeom.h"
15 #include "TH1F.h"
16 #include "TH2F.h"
17 #include "StEvent.h"
18 #include "StEventMaker/StEventMaker.h"
19 #include "StEvent/StEmcDetector.h"
20 #include "StEvent/StEmcCollection.h"
21 #include "StEventTypes.h"
22 #include "StEvent/StBbcTriggerDetector.h"
23 #include "TClonesArray.h"
24 
25 
26 /*
27 http://www.star.bnl.gov/STAR/comp/pkg/dev/StRoot/StUtilities/doc/StMessMgr.html
28 gMessMgr->Message("This is an error message.","E");
29 
30 The message type is specified with a single letter:
31 "E" = "Error"
32 "I" = "Info"
33 "W" = "Warning"
34 "D" = "Debug"
35 "Q" = "QAInfo"
36 
37 gMessMgr->Message("","W") << big_num << " seems too big." << endm;
38 gMessMgr->Info() << big_num << " seems too big." << endm;
39 char* myText = "hello";
40 gMessMgr->Message(myText);
41 
42 */
43 
44 
45 ClassImp(StRFEmcTrigMaker)
46 
47 StRFEmcTrigMaker::StRFEmcTrigMaker(const char *name):StMaker(name){
48  activeBBC=false;
49  emcGeom =0;
50 }
51 
52 StRFEmcTrigMaker::~StRFEmcTrigMaker(){
53 }
54 
55 //_____________________________________________________________________________
56 //_____________________________________________________________________________
57 
58 Int_t StRFEmcTrigMaker::Init(){
59  if (!DataMode) {
60  cout << "You are loading StMuDst" << endl;
61  } else {
62  cout << "You are loading StEvent!" << endl;
63  }
64 
65  emcGeom = StEmcGeom::getEmcGeom("bemc");
66  initHisto();
67 
68  return StMaker::Init();
69 }
70 
71 //____________________________________________________
72 //____________________________________________________
73 void StRFEmcTrigMaker::Sum(int *sum,int *sumadd){
74 
75  (*sum)=(*sumadd)+(*sum);
76 }
77 //______________________________________________________________________
78 void StRFEmcTrigMaker::Max(int *max,int *maxcomp){
79 if ((*max) < (*maxcomp))
80  (*max)=(*maxcomp);
81 }
82 
83 //______________________________________________________________________
84 
85 Int_t StRFEmcTrigMaker::getBBCtrig(){
86  return bbcTrig;
87 }
88 //______________________________________________________________________
89 Int_t StRFEmcTrigMaker::getBEMCtrigHT(int thres){
90 
91  int test=0;
92  if (BHTmaxt>thres) test=1;
93  if (BHTmaxt<=thres) test=0;
94  return test;
95 }
96 //______________________________________________________________________
97 Int_t StRFEmcTrigMaker::getBEMCtrigJP(int thres){
98 
99  int test=0;
100  if (BJPmaxt>thres) test = 1;
101  if (BJPmaxt<=thres) test = 0;
102  return test;
103 }
104 //______________________________________________________________________
105 Int_t StRFEmcTrigMaker::getBEMCtrigTOT(int thres){
106 
107 
108  int test=0;
109  if (BJPsumt>thres) test = 1;
110  if (BJPsumt<=thres) test = 0;
111  return test;
112 }
113 //______________________________________________________________________
114 Int_t StRFEmcTrigMaker::getEEMCtrigHT(int thres){
115 
116  int test=0;
117  if (EHTmaxt>thres) test=1;
118  if (EHTmaxt<=thres) test=0;
119  return test;
120 
121 }
122 //______________________________________________________________________
123 Int_t StRFEmcTrigMaker::getEEMCtrigJP(int thres ){
124 
125 
126  int test=0;
127  if (EJPmaxt>thres) test = 1;
128  if (EJPmaxt<=thres) test = 0;
129  return test;
130 }
131 //______________________________________________________________________
132 Int_t StRFEmcTrigMaker::getEEMCtrigTOT(int thres ){
133 
134  int test=0;
135  if (EJPsumt>thres) test = 1;
136  if (EJPsumt<=thres) test = 0;
137  return test;
138 }
139 
140 //______________________________________________________________________
141 Int_t StRFEmcTrigMaker::getEEMC_HT_ADC(){
142 
143  return EHTmaxt;
144 
145 }
146 
147 //______________________________________________________________________
148 Int_t StRFEmcTrigMaker::getEEMC_JP_ADC(){
149 
150  return EJPmaxt;
151 
152 }
153 
154 //______________________________________________________________________
155 Int_t StRFEmcTrigMaker::getEEMC_TOT_ADC(){
156  return EJPsumt;
157 
158 }
159 
160 //______________________________________________________________________
161 Int_t StRFEmcTrigMaker::getBEMC_HT_ADC(){
162 
163  return BHTmaxt;
164 
165 }
166 
167 //______________________________________________________________________
168 Int_t StRFEmcTrigMaker::getBEMC_JP_ADC(){
169 
170  return BJPmaxt;
171 
172 }
173 
174 //______________________________________________________________________
175 Int_t StRFEmcTrigMaker::getBEMC_TOT_ADC(){
176 
177  return BJPsumt;
178 
179 }
180 
181 //_____________________________________________________
182 //_____________________________________________________
183 void StRFEmcTrigMaker:: Clear(const char *opt){
184  // printf("\n\n Clear is called, JB\n\n");
185 
186  bbcTrig=0;
187 
188  // set array values to zero before each event
189  memset(jpBsum,0,sizeof(jpBsum));
190  memset(jpBmax,0,sizeof(jpBmax));
191  memset(jpB_hit_num,0,sizeof(jpB_hit_num));
192  memset(tpBsum,0,sizeof(tpBsum));
193  memset(tpBmax,0,sizeof(tpBmax));
194 
195 
196  // set array values to zero before each event
197  memset(jpEsum,0,sizeof(jpEsum));
198  memset(jpEmax,0,sizeof(jpEmax));
199  memset(jpE_hit_num,0,sizeof(jpE_hit_num));
200  memset(tpEsum,0,sizeof(tpEsum));
201  memset(tpEmax,0,sizeof(tpEmax));
202 
203  BHTmaxt=0; //Hold HT for whole Barrel
204  BJPmaxt=0; //Holds max JP sum for whole Barrel
205  BJPsumt=0; //Holds sum of all JP in Barrel
206  EHTmaxt=0; //Hold HT for whole EEMC
207  EJPmaxt=0; //Holds max JP sum for whole EEMC
208  EJPsumt=0; //Holds sum of all JP in EEMC
209 
210  // clear pointer to data collections
211  bbcCol=0;
212  muEmcCol=0;
213  stEmcCol=0;
214 }
215 
216 //_____________________________________________________
217 //_____________________________________________________
219 
220  // access all existing collections
221 
222  switch( DataMode) {
223  case 0: { // use muDst
224  StMuDstMaker *muDstMaker = (StMuDstMaker*)GetMaker("MuDst");
225  if(!muDstMaker) {
226  gMessMgr->Warning() <<GetName()<<" no muDstMaker !!! Game Over!"<< endm;
227  return kStErr;
228  }
229  muEvent=muDstMaker->muDst()->event();
230  bbcCol=&(muEvent->bbcTriggerDetector());
231  muEmcCol=muDstMaker->muDst()->muEmcCollection();
232  unpackEmcFromMu();
233  } break;
234  case 1: {// stEvent
235  stEvent = (StEvent *) GetInputDS("StEvent");
236  if(!stEvent) {
237  gMessMgr->Warning() <<GetName()<<" no StEvent !!! Game Over!"<< endm;
238  return kStErr;
239  }
240  StTriggerDetectorCollection *TrigDet=stEvent->triggerDetectorCollection();
241  bbcCol=&(TrigDet->bbc());
242  stEmcCol =(StEmcCollection*)stEvent->emcCollection();
243  unpackEmcFromSt();
244  } break;
245  default:
246  gMessMgr->Error() <<GetName()<<" Logic error 1, kill chain"<< endm;
247  return kStFatal;
248  }
249 
250  // common
251  unpackBBC();
252 
253  if(activeBBC && !getBBCtrig()) return kStErr;
254  // accept this eve
255 
256  fillHisto();
257  return kStOK;
258 }
259 
260 
261 //_____________________________________________________
262 //_____________________________________________________
263 void StRFEmcTrigMaker::unpackEmcFromMu(){
264  if(!muEmcCol) {
265  gMessMgr->Warning() <<GetName()<<" no muEmcCollection"<< endm;
266  return;
267  }
268 
269  // According to STAR NOTE#229A all detectors should be numbered first from the +z side (West End) looking toward the interactions region
270  // If a detector needs additional numbering on the -z side then the numbers should be consecutive with the +z elements.
271  // Following this, standing on the west side looking at the interaction region, module 58 is at 12 o'clock and in JP0 with JP1 and so on
272  // proceeding in a clockwise manner. Standing on the east side looking at the interaction region module 118 is at 12 o'clock
273  // and in JP6 with JP7 and so on proceeding in a clockwise manner.
274 
278  //JP0 goes from module=53/2 to module=3/1 (TP=0+26-29,30+56-59,60+86-89,90+116-119,120+146-149)
279  //JP1 goes from module=3/2 to module=13/1
280  //JP2 goes from module=13/2 to module=23/1
281  //JP3 goes from module=23/2 to module=33/1
282  //JP4 goes from module=33/2 to module=43/1
283  //JP5 goes from module=43/2 to module=53/1 (TP=21-25,51-55,81-85,111-115,141-145)
284  //JP6 goes from module=113/2 to module=63/1
285  //JP7 goes from module=63/2 to module=73/1
286  //JP8 goes from module=73/2 to module=83/1
287  //JP9 goes from module=83/2 to module=93/1
288  //JP10 goes from module=93/2 to module=103/1
289  //JP11 goes from module=103/2 to module=113/1
290 
291  //TP(0-29) for eta bin (1-4) all modules 1-60
292  //TP(30-59) for eta bin(5-8) all modules 1-60
293  //TP(60-89) for eta bin(9-12) all modules 1-60
294  //TP(90-119) for eta bin(13-16) all modules 1-60
295  //TP(120-149) for eta bin(17-20) all modules 1-60
296  //TP(150-179) for eta bin (1-4) all modules 61-120
297  //TP(180-209) for eta bin(5-8) all modules 61-120
298  //TP(210-239) for eta bin(9-12) all modules 61-120
299  //TP(240-269) for eta bin(13-16) all modules 61-120
300  //TP(270-299) for eta bin(17-20) all modules 61-120
301 
302  det=1;
303  for (int n=1; n<=BemcTow; n++){
304  emcGeom->getBin(n,Bmod,Beta,Bsub);
305  BTowADC = muEmcCol->getTowerADC(n,det);
306  if (BTowADC>0) {
307  //printf("n=%d, mod=%d, sub=%d, Beta=%d adc=%d\n",n,Bmod,Bsub,Beta,BTowADC);
308 
309  int jpBindex=(Bmod+Bsub+5)/10;
310  if (((Bmod+Bsub+5)>=60)&&((Bmod+Bsub+5)<=66)) {
311  jpBindex=0;
312  }
313  if ((Bmod == 60)&&(Bsub==2)) jpBindex=0;
314 
315  int tpBindex=((Bmod+Bsub-3)/2) + 30*((Beta-1)/4);
316  if ((Bmod==1)&&(Bsub==1)) {
317  tpBindex=(29 + 30*((Beta-1)/4));
318  }
319  if (Bmod>60){//need to add 150 to tp# for east side
320  tpBindex=150 + ((Bmod+Bsub-3)/2) + 30*((Beta-1)/4);
321  if ((Bmod==61)&&(Bsub==1)) {
322  tpBindex=150 + (59 + 30*((Beta+-1)/4));
323  }
324  }
325  Sum(&jpBsum[jpBindex],&BTowADC);
326  Max(&jpBmax[jpBindex],&BTowADC);
327  Sum(&tpBsum[tpBindex],&BTowADC);
328  Max(&tpBmax[tpBindex],&BTowADC);
329  jpB_hit_num[jpBindex]++;
330  //printf("jpBindex=%d, jpBsum=%d, jpBmax=%d\n",jpBindex,jpBsum[jpBindex],jpBmax[jpBindex]);
331  //printf("tpBindex=%d, tpBsum=%d, tpBmax=%d\n",tpBindex,tpBsum[tpBindex],tpBmax[tpBindex]);
332  }
333  }
334  for (int q=0; q < 6; q++){
335  Sum(&BJPsumt,&jpBsum[q]);
336  Max(&BHTmaxt,&jpBmax[q]);
337  Max(&BJPmaxt,&jpBsum[q]);
338  }
339 
340 
344  //Due to changes in EEMC muDST for 2004 sec,sub and eta all start at 1 and go to 12,5,12
345  //This code will not work on MuDst before 2004 --actually it will but it will give you GARBAGE!!!
346  //Tower 01TA01 has sec=1,sub=1,eta=1
347  //Tower 01TA02 has sec=1,sub=1,eta=2
348  //Define EMC towers in terms of jet patches (start 0 instead of 1)
349  //By definition jp0=towers 11D-1C with Eid #0-119 (m=636-35)
350  // jp1=tower 1D-3C with Eid #120-239 (m=36-155)
351  // jp2=tower 3D-5C with Eid #240-359 (m=156-275)
352  // jp3=tower 3D-5C with Eid #359-479
353  // jp4=tower 3D-5C with Eid #479-599
354  // jp5=tower 3D-5C with Eid #599-719
355  //tp=0 is defined as Eids 0-2 + 12-14
356  //tp=1 is defined as Eids 24-26 + 36-38
357  //tp=30 is defined as Eid 3-6 + 15-18
358  //tp=60 is defined as Eid 7-11 + 19-23
359 
360  NumETow=muEmcCol->getNEndcapTowerADC();
361  for (int m=0; m<NumETow; m++){
362  muEmcCol->getEndcapTowerADC(m,ETowADC,Esec,Esub,Eeta);
363  int Eid = m+84;
364  if (ETowADC){
365  //printf("m=%d,ETowADC=%d,Esec=%d,Esub=%d,Eeta=%d,Eid=%d\n",m,ETowADC,Esec,Esub,Eeta,Eid);
366  if (Eid > 719) Eid=Eid-720;
367  int jpEindex=Eid/120;
368  int tpEindex=Eid/24;
369  if (Eeta>7) {
370  tpEindex+=60;
371  }
372  if ((Eeta>3)&&(Eeta<8)) {
373  tpEindex+=30;
374  }
375 
376  Sum(&jpEsum[jpEindex],&ETowADC);
377  Max(&jpEmax[jpEindex],&ETowADC);
378  Sum(&tpEsum[tpEindex],&ETowADC);
379  Max(&tpEmax[tpEindex],&ETowADC);
380  jpE_hit_num[jpEindex]++;
381 
382  //printf("Etow=%d, Esec=%d, Esub=%d, Eeta=%d ETowADC=%d\n",m,Esec,Esub,Eeta,ETowADC);
383  //printf("Eid=%d, jpEindex=%d, jpEsum=%d, jpEmax=%d\n",Eid,jpEindex,jpEsum[jpEindex],jpEmax[jpEindex]);
384  //printf("Eid=%d, tpEindex=%d, tpEsum=%d, tpEmax=%d\n",Eid,tpEindex,tpEsum[tpEindex],tpEmax[tpEindex]);
385  }
386  }
387 
388 
389  for (int q=0; q < 6; q++){
390 
391  Sum(&EJPsumt,&jpEsum[q]);
392  Max(&EHTmaxt,&jpEmax[q]);
393  Max(&EJPmaxt,&jpEsum[q]);
394  //printf("q=%d; JPmax=%d ,JPsum=%d\n",q,jpEmax[q],jpEsum[q]);
395  }
396 
397  printf("EJPsum=%d ,EHTmax=%d,EJPmax=%d\n",EJPsumt,EHTmaxt,EJPmaxt);
398  printf("BJPsum=%d ,BHTmax=%d,BJPmax=%d\n",BJPsumt,BHTmaxt,BJPmaxt);
399 }
400 
401 
402 
403 //_____________________________________________________
404 //_____________________________________________________
405 void StRFEmcTrigMaker::unpackEmcFromSt(){
406  if(!stEmcCol) {
407  gMessMgr->Warning() <<GetName()<<" no stEmcCollection"<< endm;
408  return;
409  }
410  // According to STAR NOTE#229A all detectors should be numbered first from the +z side (West End) looking toward the interactions region
411  // If a detector needs additional numbering on the -z side then the numbers should be consecutive with the +z elements.
412  // Following this, standing on the west side looking at the interaction region, module 58 is at 12 o'clock and in JP0 with JP1 and so on
413  // proceeding in a clockwise manner. Standing on the east side looking at the interaction region module 118 is at 12 o'clock
414  // and in JP6 with JP7 and so on proceeding in a clockwise manner.
415 
419  //JP0 goes from module=53/2 to module=3/1 (TP=0+26-29,30+56-59,60+86-89,90+116-119,120+146-149)
420  //JP1 goes from module=3/2 to module=13/1
421  //JP2 goes from module=13/2 to module=23/1
422  //JP3 goes from module=23/2 to module=33/1
423  //JP4 goes from module=33/2 to module=43/1
424  //JP5 goes from module=43/2 to module=53/1 (TP=21-25,51-55,81-85,111-115,141-145)
425  //JP6 goes from module=113/2 to module=63/1
426  //JP7 goes from module=63/2 to module=73/1
427  //JP8 goes from module=73/2 to module=83/1
428  //JP9 goes from module=83/2 to module=93/1
429  //JP10 goes from module=93/2 to module=103/1
430  //JP11 goes from module=103/2 to module=113/1
431 
432  //TP(0-29) for eta bin (1-4) all modules 1-60
433  //TP(30-59) for eta bin(5-8) all modules 1-60
434  //TP(60-89) for eta bin(9-12) all modules 1-60
435  //TP(90-119) for eta bin(13-16) all modules 1-60
436  //TP(120-149) for eta bin(17-20) all modules 1-60
437  //TP(150-179) for eta bin (1-4) all modules 61-120
438  //TP(180-209) for eta bin(5-8) all modules 61-120
439  //TP(210-239) for eta bin(9-12) all modules 61-120
440  //TP(240-269) for eta bin(13-16) all modules 61-120
441  //TP(270-299) for eta bin(17-20) all modules 61-120
442 
443 
444  //StDetectorId BemcId=StDetectorId(kBarrelEmcTowerId);
445  StDetectorId BemcId=StDetectorId(kBarrelEmcTowerId);
446  StEmcDetector *EmcDet = stEmcCol->detector(BemcId); //BEMC tower detector number
447  if(EmcDet) {
448  for (UInt_t mod=1;mod<=EmcDet->numberOfModules();mod++){
449  StEmcModule* module=EmcDet->module(mod);
450  StSPtrVecEmcRawHit& hit=module->hits();
451  for(UInt_t ih=0;ih<hit.size();ih++){
452  StEmcRawHit *x=hit[ih];
453  Bmod=x->module();
454  Bsub=x->sub();
455  Beta=x->eta();
456  BTowADC=x->adc();
457  //printf("ih=%d, mod=%d eta=%d sub=%d adc=%d\n",ih,x->module(),x->eta(),x->sub(),x->adc());
458  if (BTowADC>0) {
459  // printf("mod=%d, sub=%d, Beta=%d adc=%d\n",Bmod,Bsub,Beta,BTowADC);
460  int jpBindex=(Bmod+Bsub+5)/10;
461  if (((Bmod+Bsub+5)>=60)&&((Bmod+Bsub+5)<=66)) {
462  jpBindex=0;
463  }
464  if ((Bmod == 60)&&(Bsub==2)) jpBindex=0;
465 
466  int tpBindex=((Bmod+Bsub-3)/2) + 30*((Beta-1)/4);
467  if ((Bmod==1)&&(Bsub==1)) {
468  tpBindex=(29 + 30*((Beta-1)/4));
469  }
470  if (Bmod>60){//need to add 150 to tp# for east side
471  tpBindex=150 + ((Bmod+Bsub-3)/2) + 30*((Beta-1)/4);
472  if ((Bmod==61)&&(Bsub==1)) {
473  tpBindex=150 + (59 + 30*((Beta+-1)/4));
474  }
475  }
476  Sum(&jpBsum[jpBindex],&BTowADC);
477  Max(&jpBmax[jpBindex],&BTowADC);
478  Sum(&tpBsum[tpBindex],&BTowADC);
479  Max(&tpBmax[tpBindex],&BTowADC);
480  jpB_hit_num[jpBindex]++;
481  //printf("jpBindex=%d, jpBsum=%d, jpBmax=%d\n",jpBindex,jpBsum[jpBindex],jpBmax[jpBindex]);
482  //printf("tpBindex=%d, tpBsum=%d, tpBmax=%d\n",tpBindex,tpBsum[tpBindex],tpBmax[tpBindex]);
483  }
484  }
485  }
486 
487  for (int q=0; q < 6; q++){
488  Sum(&BJPsumt,&jpBsum[q]);
489  Max(&BHTmaxt,&jpBmax[q]);
490  Max(&BJPmaxt,&jpBsum[q]);
491  }
492  }
493 
496  //Tower 01TA01 has sec=1,sub=1,eta=1
497  //Tower 01TA02 has sec=1,sub=1,eta=2
498  //Define EMC towers in terms of jet patches (start 0 instead of 1)
499  //By definition jp0=towers 11D-1C with Eid #0-119 (m=636-35)
500  // jp1=tower 1D-3C with Eid #120-239 (m=36-155)
501  // jp2=tower 3D-5C with Eid #240-359 (m=156-275)
502  // jp3=tower 3D-5C with Eid #359-479
503  // jp4=tower 3D-5C with Eid #479-599
504  // jp5=tower 3D-5C with Eid #599-719
505  //tp=0 is defined as Eids 0-2 + 12-14
506  //tp=1 is defined as Eids 24-26 + 36-38
507  //tp=30 is defined as Eid 3-6 + 15-18
508  //tp=60 is defined as Eid 7-11 + 19-23
509 
510  StDetectorId EemcId=StDetectorId(kEndcapEmcTowerId);
511  EmcDet = stEmcCol->detector(EemcId); //EEMC tower detector number
512  if(EmcDet) {
513  for (UInt_t mod=1;mod<=EmcDet->numberOfModules();mod++){
514  StEmcModule* module=EmcDet->module(mod);
515  StSPtrVecEmcRawHit& hit=module->hits();
516  for(UInt_t ih=0;ih<hit.size();ih++){
517  StEmcRawHit *x=hit[ih];
518  Esec=x->module();
519  Esub=x->sub();
520  Eeta=x->eta();
521  ETowADC=x->adc();
522  //printf("ih=%d, mod=%d eta=%d sub=%d adc=%d\n",ih,x->module(),x->eta(),x->sub(),x->adc());
523  int m =60*(Esec-1) + 12*(Esub-1) + (Eeta-1);
524  int Eid = m+84;
525  if (ETowADC){
526  //printf("m=%d,ETowADC=%d,Esec=%d,Esub=%d,Eeta=%d,Eid=%d\n",m,ETowADC,Esec,Esub,Eeta,Eid);
527  if (Eid > 719) Eid=Eid-720;
528  int jpEindex=Eid/120;
529  int tpEindex=Eid/24;
530  if (Eeta>7) {
531  tpEindex+=60;
532  }
533  if ((Eeta>3)&&(Eeta<8)) {
534  tpEindex+=30;
535  }
536 
537  Sum(&jpEsum[jpEindex],&ETowADC);
538  Max(&jpEmax[jpEindex],&ETowADC);
539  Sum(&tpEsum[tpEindex],&ETowADC);
540  Max(&tpEmax[tpEindex],&ETowADC);
541  jpE_hit_num[jpEindex]++;
542 
543  //printf("Etow=%d, Esec=%d, Esub=%d, Eeta=%d ETowADC=%d\n",m,Esec,Esub,Eeta,ETowADC);
544  //printf("Eid=%d, jpEindex=%d, jpEsum=%d, jpEmax=%d\n",Eid,jpEindex,jpEsum[jpEindex],jpEmax[jpEindex]);
545  //printf("Eid=%d, tpEindex=%d, tpEsum=%d, tpEmax=%d\n",Eid,tpEindex,tpEsum[tpEindex],tpEmax[tpEindex]);
546  }
547  }
548  }
549  for (int q=0; q < 6; q++){
550 
551  Sum(&EJPsumt,&jpEsum[q]);
552  Max(&EHTmaxt,&jpEmax[q]);
553  Max(&EJPmaxt,&jpEsum[q]);
554  //printf("q=%d; JPmax=%d ,JPsum=%d\n",q,jpEmax[q],jpEsum[q]);
555  }
556  printf("EJPsum=%d ,EHTmax=%d,EJPmax=%d\n",EJPsumt,EHTmaxt,EJPmaxt);
557  printf("BJPsum=%d ,BHTmax=%d,BJPmax=%d\n",BJPsumt,BHTmaxt,BJPmaxt);
558  } // end of EEMC-det
559 }
560 
561 
562 
563 //_____________________________________________________
564 //_____________________________________________________
565 void StRFEmcTrigMaker::unpackBBC(){
566  if(!bbcCol) {
567  gMessMgr->Warning() <<GetName()<<" no BBC Collection"<< endm;
568  return;
569  }
570  int Npmt=bbcCol->numberOfPMTs();
571  // bbc->dump();
572 
573  int Wbbc=0;
574  int Ebbc=0;
575  for (int pmt=0;pmt<Npmt;pmt++){
576  BBCadc[pmt]=bbcCol->adc(pmt);
577  int bbcadc=bbcCol->adc(pmt);
578  if (bbcadc>5) {
579  if (pmt<16) {
580  //cout << "BBC EAST = true" << endl;
581  Ebbc=1;
582  }
583  if (23<pmt && pmt<40) {
584  //cout << "BBC WEST = true" << endl;
585  Wbbc=1;
586  }
587  }
588  }
589 
590  if ((Ebbc==1)&&(Wbbc==1)){
591  bbcTrig=1;
592  }
593 
594  gMessMgr->Info() <<GetName()<<" BBCtrig="<<getBBCtrig()<<endm;
595 
596 }
597 
598 //____________________________________________________
599 //____________________________________________________
600 void StRFEmcTrigMaker::initHisto() {
601  int mxAdc=300;
602  int nAdc=mxAdc/4;
603  ha[0]=new TH1F("bHT","BEMC HT ADC",nAdc,-.5, mxAdc-.5);
604  ha[1]=new TH1F("bJP","BEMC JP ADC",nAdc,-.5, mxAdc-.5);
605  ha[2]=new TH1F("bTot","BEMC Tot ADC",nAdc,-.5, mxAdc-.5);
606  ha[3]=new TH1F("eHT","EEMC HT ADC",nAdc,-.5, mxAdc-.5);
607  ha[4]=new TH1F("eJP","EEMC JP ADC",nAdc,-.5, mxAdc-.5);
608  ha[5]=new TH1F("eTot","EEMC Tot ADC",nAdc,-.5, mxAdc-.5);
609 }
610 
611 //____________________________________________________
612 //____________________________________________________
613 void StRFEmcTrigMaker::fillHisto() {
614 
615  ha[0]->Fill(getBEMC_HT_ADC());
616  ha[1]->Fill(getBEMC_JP_ADC());
617  ha[2]->Fill(getBEMC_TOT_ADC());
618  ha[3]->Fill(getEEMC_HT_ADC());
619  ha[4]->Fill(getEEMC_JP_ADC());
620  ha[5]->Fill(getEEMC_TOT_ADC());
621 }
622 
623 
624 // $Log: StRFEmcTrigMaker.cxx,v $
625 // Revision 1.8 2011/04/11 19:35:52 fisyak
626 // Replace uint by UInt_t, use TMath
627 //
628 // Revision 1.7 2004/10/21 13:31:34 balewski
629 // to match new name of emcCollection in muDst
630 //
631 // Revision 1.6 2004/08/18 19:52:49 balewski
632 // works for BBC
633 //
634 // Revision 1.5 2004/08/18 14:56:29 balewski
635 // trying to get BBC working
636 //
StMuDst * muDst()
Definition: StMuDstMaker.h:425
static StMuEmcCollection * muEmcCollection()
returns pointer to current StMuEmcCollection
Definition: StMuDst.h:389
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:40
virtual Int_t Make()
Int_t getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const
Definition: StEmcGeom.h:321
Definition: Stypes.h:44