StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFcsEpdQaMaker.cxx
1 /*
2  *
3  * \class StFcsEpdQaMaker
4  *
5  */
6 
7 #include "StFcsEpdQaMaker.h"
8 
9 #include "StRoot/StEvent/StEvent.h"
10 #include "StRoot/St_base/StMessMgr.h"
11 #include "StRoot/StEvent/StTriggerData.h"
12 #include "StRoot/StEvent/StFcsCollection.h"
13 #include "StRoot/StEvent/StEpdCollection.h"
14 #include "StRoot/StEvent/StFcsHit.h"
15 #include "StRoot/StEvent/StEpdHit.h"
16 #include "StRoot/StFcsDbMaker/StFcsDb.h"
17 #include "StRoot/StSpinPool/StFcsRawDaqReader/StFcsRawDaqReader.h"
18 
19 #include "TH1F.h"
20 #include "TH2F.h"
21 #include "TString.h"
22 #include "TFile.h"
23 #include "TCanvas.h"
24 
25 #include <string.h>
26 #include <time.h>
27 
28 StFcsEpdQaMaker::StFcsEpdQaMaker(const Char_t* name) : StMaker(name) {
29  sprintf(mFilename,"0.epdqa.root");
30 }
31 
32 StFcsEpdQaMaker::~StFcsEpdQaMaker(){};
33 
34 Int_t StFcsEpdQaMaker::Init(){
35  mFcsDb = static_cast<StFcsDb*>(GetDataSet("fcsDb"));
36  if(!mFcsDb){
37  LOG_FATAL << "Error finding StFcsDb"<< endm;
38  return kStFatal;
39  }
40 
41  if(mFilename[0]==0 && mRun>0){
42  int yday=mRun/1000;
43  sprintf(mFilename,"%d/%d.epdqa.root",yday,mRun);
44  printf("StFcsEpdQaMaker::Init - Opening %s\n",mFilename);
45  }
46  mFile=new TFile(mFilename,"RECREATE");
47 
48  char t[100], tt[100];
49  for(int i=0; i<16; i++){
50  sprintf(t,"QtDepAdcCh%d",i);
51  sprintf(tt,"Dep01Ch%d-PP10TT%d; QTADC; DEP Fit Integral",i,i*2);
52  mQtDepA[i] = new TH2F(t,tt,256,0,1024,256,0,1024*4);
53  sprintf(t,"QtDepTacCh%d",i);
54  sprintf(tt,"Dep01Ch%d-PP10TT%d; QTTAC; DEP Peak Timebin",i,i*2);
55  mQtDepT[i] = new TH2F(t,tt,100,0,3000,100,45,56);
56  sprintf(t,"QtDepRatCh%d",i);
57  sprintf(tt,"Dep01Ch%d-PP10TT%d; DEP Peak Timebin; QTADC/DEPIntg;",i,i*2);
58  mQtDepR[i] = new TH2F(t,tt,100,44,57,100,0.0,0.8);
59  }
60 
61  return kStOK;
62 };
63 
65  mFcsCollection=0;
66  StTriggerData* trg=0;
67 
68  //Getting StFcsRawDaqReader and TriggerData
69  StFcsRawDaqReader* fcsraw=(StFcsRawDaqReader*)GetMaker("daqReader");
70  StEvent* event= (StEvent*)GetInputDS("StEvent");
71  if(fcsraw){
72  //Getting trigger data (if daq file)
73  trg = fcsraw->trgdata();
74  if(!trg){
75  LOG_INFO << "Canot find Trigger Data from StFcsRawDaqReader" << endm;
76  }
77  }else if(event){
78  trg=event->triggerData();
79  if(!trg){
80  LOG_INFO << "Canot find Trigger Data from StEvent" << endm;
81  }
82  }
83 
84  //tof multiplicity from trigger data
85  int tofmult = 0;
86  //check if FCS was readout for this event
87  if(trg){
88  tofmult = trg->tofMultiplicity();
89  //unsigned short detmask=trg->getTrgDetMask();
90  //printf("TrgDetMask = %4x\n",detmask);
91  //if(! ((detmask >> 30) & 0x1)){ //FCS_ID=30 but detmask is 16bit:O
92  //printf("No FCS readout for this event detmask=%x\n",detmask);
93  //return kStOK;
94  //}
95  //unsigned short lastdsm4 = trg->lastDSM(4);
96  //unsigned short fcs2019 = (lastdsm4 >> 10) & 0x1;
97  //printf("fcs2019=%1d\n",fcs2019);
98  unsigned short lastdsm5 = trg->lastDSM(5);
99  printf("lastdsm5=%04x tofmult=%d\n",lastdsm5,tofmult);
100  }
101 
102  if(!event) {
103  LOG_INFO << "No StEvent found" << endm;
104  }else{
105  mFcsCollection=event->fcsCollection();
106  mEpdCollection=event->epdCollection();
107  }
108  if(!mFcsCollection){
109  LOG_INFO << "No StFcsCollection found" << endm;
110  return kStErr;
111  }
112  if(!mEpdCollection){
113  LOG_INFO << "No StEpdCollection found" << endm;
114  return kStErr;
115  }
116 
117  for(int det=kFcsPresNorthDetId; det<kFcsNDet; det++){
118  int nhit=mFcsCollection->numberOfHits(det);
119  printf("StFcsEpdQaMaker found %d hits for det=%d\n",nhit,det);
120  if(nhit<=0) continue;
121  StSPtrVecFcsHit& hits = mFcsCollection->hits(det);
122  for (int i=0; i<nhit; i++){
123  int id = hits[i]->id();
124  int ehp = hits[i]->ehp();
125  int ns = hits[i]->ns();
126  int dep = hits[i]->dep();
127  int ch = hits[i]->channel();
128  int ntb = hits[i]->nTimeBin();
129  int pp,tt;
130  mFcsDb->getEPDfromId(det,id,pp,tt);
131  if(pp==11) pp=10; //HACK for run21 map
132  if(pp!=10 || ch>=16 || ch<0) continue;
133 
134  //from fits
135  float fititeg=0;
136  float fitpeak=0;
137  fititeg = hits[i]->adcSum();
138  fitpeak = hits[i]->fitPeak();
139  // printf("Dep=%02d Ch=%02d PP=%02d TT=%02d DEP=%6d PEAK=%f",
140  // dep,ch,pp,tt,fititeg,fitpeak);
141 
142  //get EPD data from epd collection
143  int nepd = mEpdCollection->epdHits().size();
144  int adc=0, tac=0;
145  for(int j=0; j<nepd; j++){
146  StEpdHit* epd = mEpdCollection->epdHits()[j];
147  int ew = epd->side();
148  int epp = epd->position();
149  int ett = epd->tile();
150  if(ew==1 && epp==pp && ett==tt){
151  adc = epd->adc();
152  tac = epd->tac();
153  break;
154  }
155  }
156  printf(" Dep=%02d Ch=%02d PP=%02d TT=%02d QT=%4d DEP=%6.1f TAC=%4d PEAK=%4.2f\n",
157  dep,ch,pp,tt,adc,fititeg,tac,fitpeak);
158  mQtDepA[ch]->Fill(adc,fititeg);
159  mQtDepT[ch]->Fill(tac,fitpeak);
160  if(fititeg>100){
161  mQtDepR[ch]->Fill(fitpeak,float(adc)/fititeg);
162  }
163  }
164  }
165 
166  return kStOK;
167 };
168 
170  mFile->Write();
171  mFile->Close();
172  printf("StFcsEpdQaMaker::Finish - Closing %s\n",mFilename);
173  return kStOK;
174 };
175 
176 ClassImp(StFcsEpdQaMaker)
177 
178 /*
179  * $Id: StFcsEpdQaMaker.cxx,v 1.1 2021/05/30 21:33:05 akio Exp $
180  * $Log: StFcsEpdQaMaker.cxx,v $
181  * Revision 1.1 2021/05/30 21:33:05 akio
182  * QA for EPD West from DEP and QT comparison
183  *
184  */
virtual Int_t Make()
short side() const
+1 if tile is on West side; -1 if on East side
Definition: StEpdHit.h:145
int position() const
position of supersector on a wheel [1,12]
Definition: StEpdHit.h:147
int adc() const
ADC value [0,4095].
Definition: StEpdHit.h:149
void getEPDfromId(int det, int id, int &pp, int &tt)
Get FCS&#39;s EPD map foom EPD mapping.
Definition: StFcsDb.cxx:1484
Stores information for tiles in STAR Event Plane Detector.
Definition: StEpdHit.h:43
virtual Int_t Finish()
Definition: Stypes.h:40
int tile() const
tile on the supersector [1,31]
Definition: StEpdHit.h:148
Definition: Stypes.h:44
int tac() const
TAC value [0,4095].
Definition: StEpdHit.h:150