StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StBemcBeamBckgFinderMaker.cxx
1 /******************************************************************************
2  * $Id: StBemcBeamBckgFinderMaker.cxx,v 1.13 2010/04/15 19:13:29 mattheww Exp $
3  * \author Issam Qattan , IUCF, 2006
4  ******************************************************************************
5  * Description:
6  * Pattern recognition of the Barrel Beam Background on an event-by-event basis.
7  *******************************************************************************
8  */
9 
10 #include <TH2.h>
11 #include <TCanvas.h>
12 #include <TStyle.h>
13 #include <stdio.h>
14 
15 #include "StBemcBeamBckgFinderMaker.h"
16 #include "StChain.h"
17 #include "St_DataSetIter.h"
18 #include "StDAQMaker/StDAQReader.h"
19 
20 #include <StMessMgr.h>
21 #include <StEmcUtil/geometry/StEmcGeom.h>
22 #include <StEmcUtil/database/StBemcTables.h>
23 
24 #include "StEmcUtil/database/StEmcDecoder.h"
25 #include "StEmcRawMaker/defines.h"
26 
27 //MuDst
28 #include "StMuDSTMaker/COMMON/StMuDstMaker.h"
29 #include "StMuDSTMaker/COMMON/StMuDst.h"
30 #include "StMuDSTMaker/COMMON/StMuTriggerIdCollection.h"
31 #include "StMuDSTMaker/COMMON/StMuEvent.h"
32 
33 
35 
37 
38  mGeomB = 0;
39  mHList = 0;
40  mMappB = 0;
41  memset(mhisto, 0,sizeof(mhisto));
42  memset(mevtH, 0,sizeof(mevtH));
43  mAccEve=mInpEve=0;
44  mTrigId=0;
45  mAdcThreshold=0;
46  mAdcSumThreshold=0;
47  mpattern=0;
48 
49 }
50 /*=============================================================================*/
51 /*=============================================================================*/
52 
53 StBemcBeamBckgFinderMaker::~StBemcBeamBckgFinderMaker(){
54  // Nothing for now
55 }
56 
57 /*=============================================================================*/
58 /*=============================================================================*/
59 
60 Int_t StBemcBeamBckgFinderMaker::Init(){
61 
62  mGeomB = StEmcGeom::instance("bemc");
63 
64  assert(mHList);
65  mhisto[0]= new TH1F("Frequency1","Towers Frequency (status =1); Tower SoftId; Frequency",mxSoftId,0.5,mxSoftId+0.5);
66  mhisto[1]= new TH1F("Frequency2","Towers Frequency (status !=1); Tower SoftId; Frequency",mxSoftId,0.5,mxSoftId+0.5);
67  mhisto[2]= new TH1F("Frequency3","bin1= All, bin3= Accept, bin5= !Backgnd, bin7=Backgnd bin9=East bin11=Center bin13=West; Bin No",15,0.5,15.5);
68  mhisto[3]= new TH1F("StartingphiVspattern","Pattern Starting #phi Weighted By Pattern Length; #phi[deg]",360,0.5,360.5);
69  mhisto[4]= new TH1F("StartingphiVsadcsum","Pattern Starting #phi Weighted By Pattern Summed Adc; #phi[deg]",360,0.5,360.5);
70 
71  mhisto[5]= new TH1F("Allbunchcrossing7","All Bunch Crossings; Bunch Crossing No (7bit)",128,-0.50,127.5);
72  mhisto[6]= new TH1F("Trigbunchcrossing7","Trigger Passing Bunch Crossings; Bunch Crossing No (7bit)",128,-0.50,127.5);
73  mhisto[7]= new TH1F("Bckgbunchcrossing7","Background Passing Bunch Crossings; Bunch Crossing No (7bit)",128,-0.50,127.5);
74 
75  mevtH[0]= new TH2F("BtowphiVseta","Barrel Towers: ADC-PED; #eta Bin; #phi Bin",mxEta,-0.5,mxEta-0.5,mxPhi,-0.5,mxPhi-0.5);
76  mevtH[1]= new TH2F("StartingphiVseta","Pattern Starting #phi Vs. Pattern Starting #eta; #eta; #phi(deg)",10,-1,1,30,0.5,360.5);
77  mevtH[2]= new TH2F("AvgWeightedPhiVsEta","Pattern #phi Vs. Pattern Weighted Average #eta; <#eta>; #phi(deg)",20,-1,1,60,0.5,360.5);
78 
79  for(int i=0; i<mMaxH; i++){
80  if(mhisto[i])
81  mHList->Add(mhisto[i]);
82  }
83 
84  for(int i=0; i<mMaxevtH; i++){
85  if(mevtH[i])
86  mHList->Add(mevtH[i]);
87  }
88 
89  LOG_INFO<< Form("%s::Init() accept TrigID=%d\n",GetName(),mTrigId)<<endm;
90 
91  return StMaker::Init();
92 }
93 
94 /*=============================================================================*/
95 /*=============================================================================*/
96 
97 void StBemcBeamBckgFinderMaker::Clear(const Option_t* option){
98 
99  mevtH[0]->Reset(); //keep only eta-phi-adc histogram for the last event processed.
100  memset(mAdcArray, 0,sizeof(mAdcArray)); //array cleared for each event
101  memset(mPattSoftId, 0,sizeof(mPattSoftId)); //array cleared for each event
102 
103  mDecision=0;
104  metaBegin=0;
105  metaEnd=0;
106  mphiBegin=0;
107  mpatternLength=0;
108  msumAdc=0;
109  mAvgEtaAdc=0;
110  mSearchDone=0;
111 
112  StMaker::Clear(option);
113 }
114 
115 /*=============================================================================*/
116 /*=============================================================================*/
117 
118 // ***** Make: this method is called in loop for each event
119 
121 
122 
123  mInpEve++;
124  printf("\n");
125  //printf("In%s::Make by Issam Qattan: Working on Event=%d\n",GetName(),mInpEve);
126  StMuDstMaker *muMk = (StMuDstMaker*)GetMaker("MuDst");
127  assert(muMk);
128 
129  StMuEmcCollection *emc = muMk->muDst()->muEmcCollection();
130 
131  //========process 7-bit Bunch-Crossing (bx7)==================================
132  StMuDst *dst = muMk->muDst();
133  StMuEvent *muEve = dst->event();
134  StL0Trigger &trig = muEve->l0Trigger();
135  StEventInfo &info = muEve->eventInfo();
136 
137  int bx7=trig.bunchCrossingId7bit(info.runId());
138  mhisto[5]->Fill(bx7); //fill 7bit bunch crossing for every event we processed.
139 
140  //printf("eve=%d bx7=%d\n",mInpEve,bx7);
141 
142 
143  //========process trigger info==================================================
144 
145  assert(muEve);
146  StMuTriggerIdCollection ticB = muEve -> triggerIdCollection();
147  const StTriggerId *nominal = &ticB.nominal();
148 
149  mhisto[2]->Fill(1); //bin=1 gets filled once for every event we processed.
150 
151 
152  if(mTrigId<0) { // print all triggers
153  vector<unsigned int> nominalL=nominal->triggerIds();
154  printf("trigL len=%d totEve=%d\n",nominalL.size(),mInpEve);
155  uint ii;
156  for(ii=0;ii<nominalL.size();ii++)printf("trgID=%d, ",nominalL[ii]);
157  printf("\n");
158  }
159 
160  if(mTrigId>0 && !nominal->isTrigger(mTrigId)) return kStOK; // discard events
161  mAccEve++;
162 
163 
164  //...........................................................................//
165  //......................... B T O W ......................................//
166 
167  int id; // id = this is soft ID
168 
169  for (id=0; id<mxSoftId; id++) {
170 
171  int ietabin = mdb_btowetaBin[id];
172  int iphibin = mdb_btowphiBin[id];
173 
174  // int irdo = mdb_btowRdo[id];
175  int isoftid = mdb_btowSoftId[id];
176  // int istatus = mdb_btowStat[id];
177  float rawAdc = emc->getTowerADC(id+1);
178 
179  if(mdb_btowStat[id] != 1) continue; //added to skip failed towers (status != 1)
180 
181  float ped = mdb_btowPed[id];
182  float adc = (rawAdc - ped);
183 
184  //printf("counter id=%d softId=%d Rdo=%d status=%d etabin=%d phibin=%d rawadc=%d ped=%f adc=%f\n",id,isoftid,irdo,istatus,ietabin,iphibin,rawAdc,ped,adc);
185 
186  if(adc>mAdcThreshold) {
187  //mhisto[0]->Fill(id+1);
188  FillAdc(adc,ietabin,iphibin,isoftid);
189  mevtH[0]->Fill(ietabin,iphibin,adc);
190  }
191 
192  }
193 
194  mhisto[2]->Fill(3); //bin=3 gets filled once for every event we accepted (satisfied trigger condition).
195  mhisto[6]->Fill(bx7); //fill 7bit bunch crossing for every event we accepted (satisfied trigger condition).
196 
197  int myetapattern,myphipattern,myetaendpattern,mypatternlength;
198  float myadcsum,myavgetaadc;
199 
200 
201  bool IsItBackGround = CheckPatternType3(myetapattern,myphipattern,myetaendpattern,mypatternlength,myadcsum,myavgetaadc);
202 
203  mSearchDone=true; //set to true when search for background is done.
204 
205  if(IsItBackGround==1){
206  mDecision=1;
207  metaBegin=myetapattern;
208  mphiBegin=myphipattern;
209  metaEnd=myetaendpattern;
210  mpatternLength=mypatternlength;
211  msumAdc=myadcsum;
212  mAvgEtaAdc=myavgetaadc;
213 
214  //printf("Yes: Event=%d IS a Background\n",mInpEve);
215  mhisto[2]->Fill(7); //bin=7 gets filled for every background event we processed.
216  mevtH[1]->Fill((myetapattern-20.)/20.,myphipattern*3);
217  mevtH[2]->Fill((myavgetaadc-20.)/20.,myphipattern*3);
218  mhisto[3]->Fill(myphipattern*3,mypatternlength);
219  mhisto[4]->Fill(myphipattern*3,myadcsum);
220  mhisto[7]->Fill(bx7); //gets filled for every background event we processed.
221 
222 
223 
224  if(metaEnd<20){
225  mhisto[2]->Fill(9); // EastSide Background
226  sprintf(mLocation,"East");
227  }
228 
229  if((metaBegin<20)&&(metaEnd>=20)){
230  mhisto[2]->Fill(11); // CentralSide Background
231  sprintf(mLocation,"Central");
232  }
233 
234  if(metaBegin>19){
235  mhisto[2]->Fill(13); // WestSide Background
236  sprintf(mLocation,"West");
237  }
238 
239  //printf("Finally my background pattern starts at: phi=%d eta=%d pattern length=%d sum adc=%f\n", myphipattern,myetapattern, mypatternlength, myadcsum);
240  }
241 
242  else {
243  mDecision=0;
244  //printf("Finally No: Event=%d IS NOT a Background\n",mInpEve);
245  mhisto[2]->Fill(5); //bin=5 gets filled once for every event that is not background.
246  }
247 
248  //printf("**********************************************************************************\n");
249 
250  if(IsItBackGround==1) {
251  if((mMaxYesPlots--)>0) PlotOneEvent();
252  } else {
253  if((mMaxNoPlots--)>0) PlotOneEvent();
254  }
255 
256  return kStOK;
257 }
258 
259 /*=============================================================================*/
260 /*=============================================================================*/
261 
262 Int_t StBemcBeamBckgFinderMaker::InitRun(int runNo){
263 
264  LOG_INFO << GetName()<<"::InitRun() run=" <<runNo<<endm;
265  mRunNumber=runNo;
266 
267  // this is how the BTOW mapping is accesible
268  assert(mMappB==0);
269  StEmcDecoder *mMappB = new StEmcDecoder(GetDate(),GetTime());
270 
271  StBemcTables *myTable=new StBemcTables;
272  myTable->loadTables(this );
273 
274 
275  memset(mdb_btowPed, 0,sizeof(mdb_btowPed));
276  memset(mdb_btowStat, 0,sizeof(mdb_btowStat));
277 
278  memset(mdb_btowetaBin, 0,sizeof(mdb_btowetaBin));
279  memset(mdb_btowphiBin, 0,sizeof(mdb_btowphiBin));
280 
281 
282  memset(mdb_btowRdo, 0,sizeof(mdb_btowRdo));
283  memset(mdb_btowSoftId, 0,sizeof(mdb_btowSoftId));
284 
285 
286  int softID;
287 
288  for(softID=1; softID<=BTOWSIZE; softID++) {
289 
290  //...................................................//
291  //........... querry BTOW DB/geom ...................//
292 
293  int status;
294  myTable->getStatus(BTOW, softID, status);
295  int m,e,s;
296  mGeomB->getBin(softID,m,e,s);
297 
298  float etaF,phiF;
299  mGeomB->getEta(m,e,etaF);
300  mGeomB->getPhi(m,s,phiF); // -pi <= phi < pi
301  if(phiF<0) phiF+=2*C_PI; // I want phi in [0,2pi]
302 
303  float ped,sig;
304  myTable->getPedestal(BTOW, softID, 0, ped,sig);
305 
306  float gain;
307  myTable->getCalib(BTOW, softID, 1, gain);
308 
309  mdb_btowPed[softID-1]=ped;
310 
311  int Rdo;
312  assert(mMappB->GetDaqIdFromTowerId(softID,Rdo)==1);
313 
314  mdb_btowRdo[softID-1] = Rdo;
315  mdb_btowSoftId[softID-1] = softID;
316 
317  status= GetNewStatus(softID,Rdo,mRunNumber);
318  mdb_btowStat[softID-1]=status;
319 
320  // printf("back in InitRun: id=%d status=%d\n",softID,status);
321 
322  if(status != 1) mhisto[1]->Fill(softID);
323 
324  int kEta = int ((etaF+1.)*20.);
325  int kPhi = int (phiF*19.10); // New Mapping == 19.10=(120/2pi)
326  //int kPhi = 24- (int) (phiF/C_PI*60.); // Old Jan mapping
327  //if(kPhi<0) kPhi+=120; // Old Jan mapping
328 
329  mdb_btowetaBin[softID-1] = kEta;
330  mdb_btowphiBin[softID-1] = kPhi;
331 
332  mSoftId[kPhi][kEta] = softID;
333 
334  //printf("soft=%4d Rdo=%d ped=%f eta=%.2f phi(rad)=%f Etabin=%d Phibin=%d\n",softID,Rdo,ped,etaF,phiF,kEta,kPhi);
335 
336  }
337 
338  return kStOk;
339 
340 }
341 
342 /*=============================================================================*/
343 /*=============================================================================*/
344 
346  //
347  LOG_INFO << Form("%s::Finish() accepted %d eve of %d input events, TrigID=%d\n",GetName(),mAccEve, mInpEve,mTrigId)<<endm;
348  return kStOK;
349 }
350 
351 /*=============================================================================*/
352 /*=============================================================================*/
353 
354 Int_t StBemcBeamBckgFinderMaker::GetNewStatus(int id, int rdo, int run){
355 
356  const int db_btowAddMaskA[] = {4505,533,577,578,579,580,4020,4019,816,839,1612,2916,2897,1773,1774,1775,1776,2085,2105,2257,0};
357  const int db_btowAddMaskB[] = {3287,3309,3407,3674,0};
358  const int db_btowAddMaskC[] = {0};
359 
360  const int *db_btowAddMask;
361  int xMask;
362 
363  if(run>7000000) {
364  db_btowAddMask=db_btowAddMaskA;
365  xMask=20;
366  } else {
367  if(run==6169089) {
368  db_btowAddMask=db_btowAddMaskB;
369  xMask=4;
370  } else {
371  db_btowAddMask=db_btowAddMaskC;
372  xMask=0;
373  }
374  }
375 
376 
377  for(int i=0; i<xMask; i++) {
378 
379  //printf("db_btowAddMask =%d\n",db_btowAddMask[i]);
380 
381  if(db_btowAddMask[i]<=0) {
382  printf("Error Error Error: your masked array db_btowAddMask=0 for this run number\n");
383  assert(1==2);
384  }
385 
386  if(db_btowAddMask[i]==id){
387  printf("===================================================================\n");
388  printf("===================================================================\n");
389  printf("RunNo=%d :: tower with softid=%d and Rdo=%d will be masked\n",run,id,rdo);
390  return (0);
391  }
392  }
393 
394  return (1);
395 }
396 
397 /*=============================================================================*/
398 /*=============================================================================*/
399 
400 void StBemcBeamBckgFinderMaker::FillAdc(float adc,int ieta,int iphi,int isoft){
401  assert(ieta>=0);
402  assert(ieta<mxEta);
403 
404  assert(iphi>=0);
405  assert(iphi<mxPhi);
406 
407  mAdcArray[iphi][ieta] = adc;
408 
409  //printf("In FillAdcUsed(): ieta=%d iphi=%d isoftid=%d mAdcArray[iphi][ieta]=%f\n",ieta,iphi,isoft,adc);
410 
411 }
412 
413 /*=============================================================================*/
414 /*=============================================================================*/
415 
416 void StBemcBeamBckgFinderMaker::PlotOneEvent(){
417 
418  char Decision[10];
419  char HistoTitle_Out[2000]; //title of the histogram. Replacing the original
420  char FileTitle_Out[1000]; // mevtH[0] histograms will be saved for each accepted event in a ps and/or gif files.
421 
422 
423  if(mDecision==1){
424  sprintf(Decision,"Yes");
425  sprintf(HistoTitle_Out,"R=%d E=%05d T=%d Bg=%s #eta1=%d #phi=%d #eta2=%d L=%d Adc+=%f Loct=%s",mRunNumber,mInpEve,mTrigId,Decision,metaBegin,mphiBegin,metaEnd,mpatternLength,msumAdc,mLocation);
426  }
427 
428  if(mDecision==0){
429  sprintf(Decision,"No");
430  sprintf(HistoTitle_Out,"Run=%d, Eve=%05d, TrgId=%d, Bckg=%s",mRunNumber,mInpEve,mTrigId,Decision);
431  }
432 
433 
434  sprintf(FileTitle_Out,"event%05d_%s.ps",mInpEve,Decision); // name of ps file
435  //sprintf(FileTitle_Out,"event%05d_%s.gif",mInpEve,Decision); // name of gif file if needed
436  TCanvas mycanvas("aa","aa",600,800); //create a canvas
437  mevtH[0]->SetTitle(HistoTitle_Out);
438  mevtH[0]->GetXaxis()->CenterTitle();
439  mevtH[0]->GetYaxis()->CenterTitle();
440  gStyle->SetPalette(1,0);
441  gStyle->SetOptStat(0);
442  gStyle->SetTitleFontSize(0.08);
443  mevtH[0]->Draw("colz");
444  mycanvas.Print(FileTitle_Out);
445  //printf("Out of ::PlotOneEvent(), just finished printing to a file=%s\n",FileTitle_Out);
446 }
447 
448 /*=============================================================================*/
449 /*=============================================================================*/
450 
451 
452 Int_t StBemcBeamBckgFinderMaker::CheckPatternType3(int &etaBegin, int &phiBegin, int &etaEnd, int &patternLength, float &sumAdc, float &AverageWeightedEta){
453 
454  int maxLength = 0; // initial value for number of adjacent towers found (in eta direction).
455  int i,j,k,n;
456 
457  int FirstEtaFound;
458  float FirstEtaAdc,sum=0.;
459  float sumEtaAdc=0.;
460 
461  //printf("**********************************************************************************\n");
462  //printf("In ::CheckPatternType3(): These towers are used in pattern searched\n");
463 
464  for(i=0; i<mxPhi; i++){ //loop over phi
465  float *adc = mAdcArray[i];
466  n = 0; //counter over adjacent towers
467  maxLength = 0;
468  memset(mPattSoftId, 0,sizeof(mPattSoftId));
469 
470 
471  for(j=0; j<mxEta; j++){ //loop over eta
472  if(adc[j]<=0) continue;
473 
474  n = 1;
475  FirstEtaFound=j;
476  FirstEtaAdc = adc[j];
477  sum = FirstEtaAdc;
478  sumEtaAdc=(FirstEtaFound*sum);
479  mPattSoftId[n-1] = mSoftId[i][j]; //[n-1] to make sure array starts from zero
480 
481  //printf("start mPattSoftId[%d]=%d\n",n,mPattSoftId[n-1]);
482  //printf("(Start probing from: phi=%d eta=%d adc=%f n=%d\n",i,j,adc[j],n);
483 
484  for(k=j+1; k<mxEta; k++){
485 
486  if((adc[k]<mAdcThreshold)&&(adc[k+1]<mAdcThreshold)) {
487 
488  break; //exist for(k=j+1; k<mxEta; k++) loop.
489  }
490 
491  sum = sum + adc[k];
492  sumEtaAdc=sumEtaAdc+(k*adc[k]);
493  n++;
494  mPattSoftId[n-1] = mSoftId[i][k];
495 
496  //printf("next mPattSoftId[%d]=%d\n",n,mPattSoftId[n-1]);
497  //printf("(k=eta+1 i.e, looping over adjacent towers if any) In k loop: phi=%d eta=%d k=%d adc=%f n=%d\n",i,j,k,adc[k],n);
498  }
499 
500  maxLength=0;
501 
502  if(maxLength<n) maxLength=n;
503 
504  j=k;
505 
506 
507  if((maxLength>=mpattern)&&(sum>=mAdcSumThreshold)){
508  //printf("**********************************************************************************\n");
509  //printf("Found The First BACKGROUND Pattern. No Need To Search Further: Will Stop Here\n");
510  //printf("Number of adjacent towers: maxLength = %d adjacent towers >= background pattern length=%d\n",maxLength,pattern);
511  //printf("******************** SUMMARY SUMMARY SUMMARY *************************************\n");
512  //printf("phi=%d eta=%d pattern length=%d sum adc=%f\n",i,FirstEtaFound,maxLength,sum);
513  //printf("Yes %d %d %d %f\n",i,FirstEtaFound,maxLength,sum);
514  //printf("**********************************************************************************\n");
515 
516  phiBegin=i;
517  etaBegin=FirstEtaFound;
518  patternLength=maxLength;
519  etaEnd=(etaBegin + (patternLength-1));
520  sumAdc=sum;
521  AverageWeightedEta=(sumEtaAdc/sum);
522  //printf("my average weighted eta for this pattern=%f\n",AverageWeightedEta);
523  return(true);
524  break;
525  }
526  }
527  }
528  //printf("Not Background pattern\n");
529  return(false);
530 }
531 
532 /*=============================================================================*/
533 /*=============================================================================*/
534 
535 void StBemcBeamBckgFinderMaker::GetDecision(int &fDecision, int &eta1, int &phi1, int &eta2, int &patternleng, float &Adcsum){
536 
537  if((mSearchDone)&&(mDecision)){
538 
539  fDecision=mDecision;
540  eta1=metaBegin;
541  eta2=metaEnd;
542  phi1=mphiBegin;
543  patternleng=mpatternLength;
544  Adcsum=msumAdc;
545  } else {
546 
547  if((mSearchDone)&&(!mDecision)){
548  fDecision=mDecision;
549  } else {
550  fDecision=-1;
551  }
552  eta1=0;
553  eta2=0;
554  phi1=0;
555  patternleng=0;
556  Adcsum=0;
557 
558  }
559 }
560 
561 /*=============================================================================*/
562 /*=============================================================================*/
563 
564 /**********************************************************************
565  $Log: StBemcBeamBckgFinderMaker.cxx,v $
566  Revision 1.13 2010/04/15 19:13:29 mattheww
567  fixed some future gcc issues
568 
569  Revision 1.12 2009/02/04 20:33:30 ogrebeny
570  Moved the EEMC database functionality from StEEmcDbMaker to StEEmcUtil/database. See ticket http://www.star.bnl.gov/rt2/Ticket/Display.html?id=1388
571 
572  Revision 1.11 2007/08/22 14:03:50 kocolosk
573  #included some additional headers that used to be included implicitly from StBemcTables
574 
575  Revision 1.10 2006/06/27 16:54:47 qattan
576  *** empty log message ***
577 
578  Revision 1.9 2006/06/27 15:41:30 qattan
579  *** empty log message ***
580 
581  Revision 1.8 2006/06/13 21:50:36 qattan
582  *** empty log message ***
583 
584  Revision 1.7 2006/06/13 21:42:42 qattan
585  *** empty log message ***
586 
587  Revision 1.6 2006/06/13 21:26:25 qattan
588  *** empty log message ***
589 
590  Revision 1.5 2006/05/30 23:05:23 qattan
591  check5
592 
593  Revision 1.4 2006/05/30 22:40:29 qattan
594  check4
595 
596  Revision 1.3 2006/05/30 22:21:41 qattan
597  *** empty log message ***
598 
599  Revision 1.2 2006/05/30 20:12:56 qattan
600  fix 1
601 
602  Revision 1.1 2006/05/30 20:08:03 qattan
603  start1
604 
605 
606 */
StMuDst * muDst()
Definition: StMuDstMaker.h:425
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StMaker.cxx:634
void getPedestal(Int_t det, Int_t softId, Int_t cap, Float_t &ped, Float_t &rms) const
Return pedestal mean and rms.
void loadTables(StMaker *anyMaker)
load tables.
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
void getCalib(Int_t det, Int_t softId, Int_t power, Float_t &calib) const
Return calibration constant.
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:40
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
void getStatus(Int_t det, Int_t softId, Int_t &status, const char *option="") const
Return status.
int GetDaqIdFromTowerId(int softId, int &RDO) const
Get Daq Id from Software Id for towers.
Collection of trigger ids as stored in MuDst.
Definition: Stypes.h:41