StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFpostQaMaker.cxx
1 /*
2  *
3  * \class StFpostQaMaker
4  *
5  */
6 
7 #include "StFpostQaMaker.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/StFmsCollection.h"
13 #include "StRoot/StEvent/StFmsHit.h"
14 #include "StRoot/StFmsDbMaker/StFmsDbMaker.h"
15 #include "StRoot/StSpinPool/StFpsRawDaqReader/StFpsRawDaqReader.h"
16 
17 #include "TH1F.h"
18 #include "TH2F.h"
19 #include "TString.h"
20 #include "TFile.h"
21 #include "TCanvas.h"
22 
23 #include <string.h>
24 #include <time.h>
25 
26 
27 //Default Constructor
28 StFpostQaMaker::StFpostQaMaker(const Char_t* name) : StMaker(name),mFmsDbMkr(0),mFmsCollectionPtr(0), mRun(0) {};
29 
30 //Defatult Destructor
31 StFpostQaMaker::~StFpostQaMaker(){};
32 
33 //Function that will initailze variables and histograms
34 Int_t StFpostQaMaker::Init()
35 {
36  //Opening the FMS database
37  mFmsDbMkr = static_cast< StFmsDbMaker*>(GetMaker("fmsDb"));
38  if( !mFmsDbMkr )
39  {
40  LOG_FATAL << "Error finding StFmsDbMaker"<< endm;
41  return kStFatal;
42  }
43 
44  int yday=mRun/1000;//The usual 'yday' variable which looks like year followed by day
45 
46  //Making the root file that data will be written to
47  sprintf(mFilename,"%d/%d.root",yday,mRun);
48  printf("StFpostQaMaker::Init - Opening %s\n",mFilename);
49  mFile=new TFile(mFilename,"RECREATE");
50 
51  char name[100];//variable to hold names of histograms
52 
53  //Histograms for the Data size
54  mDataSize[0] = new TH1F("TotalSize","TotalSize",100,-1.0,4.0);
55  mDataSize[1] = new TH1F("DataSize","DataSize",100,-1.0,3.0);
56 
57  mRccDiff[0] = new TH1F("RccDiff_Full","RccDiff_Full",100,0,120000);
58  mRccDiff[1] = new TH1F("RccDiff_Zoom","RccDiff_Zoom",100,0,10000);
59 
60  //Histograms for Crossings
61  for(int i=0; i<mNPREPOST*2+1; i++)
62  {
63  int x=i-mNPREPOST;
64  sprintf(name,"Xing=%d",x);
65  mXing[i] = new TH1F(name,name,100,0.0,3000.0);
66  }
67 
68  //Histograms for the ADCs
69  if(mPed==0)//physics
70  {
71  mAdc2[0] = new TH2F("Adc2", "Adc2", 252,0.0,252.0,64,0.0,4096.0);
72  mAdc2[1] = new TH2F("Adc2z","Adc2z",252,0.0,252.0,50,0.0,200.0);
73  for(int i=0; i<mNID; i++){
74  sprintf(name,"ADC%03d",i);
75  mAdc[i][0]=new TH1F(name,name,128,0.0,4096.0);
76  sprintf(name,"ADC%03dz",i);
77  mAdc[i][1]=new TH1F(name,name,150,0.0,300.0);
78  }
79  }
80  else//not physics (pedestal)
81  {
82  mAdc2[0] = new TH2F("Adc2", "Adc2", 252,0.0,252.0, 100,64.0,4096.0);
83  mAdc2[1] = new TH2F("Adc2z","Adc2z",252,0.0,252.0, 100,0.0,200.0);
84  for(int i=0; i<mNID; i++){
85  sprintf(name,"ADC%03d",i);
86  mAdc[i][0]=new TH1F(name,name,128,0.0,4096.0);
87  sprintf(name,"ADC%03dz",i);
88  mAdc[i][1]=new TH1F(name,name,100,0.0,200.0);
89  }
90  }
91 
92  //Histograms for number of hits and hits per slat
93  /*
94  for( int q = 0; q <= 2 )
95  {
96  for( int l = 0; l <= 2; l++ )
97  {
98  //if( q == 2 && l == 1 ){mNHit[q][l]=0;continue;}
99  sprintf(name,"NHIT_Q%1dL%1d",q+1,l+1);
100  mNHit[q][l]=new TH1F(name,name,22,0.0,22.0);
101  sprintf(name,"HIT_Q%1dL%1d",q+1,l+1);
102  mHit[q][l]=new TH1F(name,name,14,0.5,14.5);
103  }
104 
105  //Layer 4
106  sprintf(name,"NHIT_Q%1dL%1d",q+1,4);
107  mNHit[q][3]=new TH1F(name,name,22,0.0,22.0);
108  sprintf(name,"HIT_Q%1dL%1d",q+1,4);
109  mHit[q][3]=new TH1F(name,name,25,0.5,25.5);
110 
111  //Layer 5 Top
112  sprintf(name,"NHIT_Q%1dL%1dT",q+1,5);
113  mNHit[q][4]=new TH1F(name,name,22,0.0,22.0);
114  sprintf(name,"HIT_Q%1dL%1dT",q+1,5);
115  mHit[q][4]=new TH1F(name,name,21,0.5,21.5);
116 
117  //Layer 5 Bottom
118  sprintf(name,"NHIT_Q%1dL%1dB",q+1,5);
119  mNHit[q][5]=new TH1F(name,name,22,0.0,22.0);
120  sprintf(name,"HIT_Q%1dL%1dB",q+1,5);
121  mHit[q][5]=new TH1F(name,name,21,21.5,43.5);
122  }
123 
124  */
125 
126  for(int q=0; q<mNQ; q++)
127  {
128  for(int l=0; l<mNL; l++)
129  {
130  sprintf(name,"NHIT_Q%1dL%1d",q+1,l+1);
131  mNHit[q][l]=new TH1F(name,name,22,0.0,22.0);
132  sprintf(name,"HIT_Q%1dL%1d",q+1,l+1);
133  mHit[q][l]=new TH1F(name,name,43,0.5,43.5);
134  }
135  }
136 
137  //Histograms for the triggers
138  for(int t=0; t<mNTRG+1; t++)
139  {
140  sprintf(name,"NHIT_TRG%02d",t);
141  mNHitTrg[t]=new TH1F(name,name,241,0.0,241.0);
142  }
143  sprintf(name,"NHIT_TRG");
144  mNHitTrg2=new TH2F(name,name,241,0.0,241.0,64,0.0,64.0);
145 
146  return kStOK;
147 };
148 
149 //Function that will use the database and the collection pointer to fill the histograms
151 {
152  //cout << "In Make" << endl;
153  StEvent* eventPtr=0;
154  mFmsCollectionPtr=0;
155 
156  eventPtr= (StEvent*)GetInputDS("StEvent");
157  if(!eventPtr) { LOG_INFO << "No StEvent found" << endm;}
158  else{ mFmsCollectionPtr=eventPtr->fmsCollection();}
159  if(!mFmsCollectionPtr)
160  {
161  LOG_INFO << "No StFmsCollection found" << endm;
162  return kStErr;
163  }
164 
165  unsigned int nhit=mFmsCollectionPtr->numberOfHits();
166  StSPtrVecFmsHit hits = mFmsCollectionPtr->hits();
167  printf("StFpostQaMaker found %d hits\n",nhit);
168  /*
169  if( nhit == 0 )
170  {
171  printf("There were no hits skipping this iteration\n");
172  return kStOK;
173  }
174  */
175  int nfpostdata=0; //counter for number of fpost data points for xing=0
176  int nfpostdatatot=0;//counter for total number of fpost data points
177  int nh[mNQ][mNL]; memset(nh,0,sizeof(nh));//number of hits array
178  int nhtot=0;//total number of hits over all quadrands and layers
179 
180  //Loop through all the hits and fill adcs, xing, and hits
181  //cout << "Beginning for loop for hits" << endl;
182  for (unsigned int i=0; i<nhit; i++)
183  {
184  //cout << "inside for loop" << endl;
185  int det = hits[i]->detectorId();
186  //Check if the detector id is the same as for the fpost
187  //cout << "Detector ID: " << det << endl;
188  if( det == 15 )//See StRoot/StEvent/StEnumerations.h to see what is the correct ID
189  {
190  ran_fpost = true;
191  nfpostdatatot++;
192  int xing = hits[i]->tdc(); if(xing>65536/2){xing-=65536;}
193  int adc = hits[i]->adc();
194  int qt = hits[i]->qtSlot();
195  int ch = hits[i]->qtChannel();
196  int slatid = mFmsDbMkr->fpostSlatidFromQT(qt,ch);
197 
198  int q=0;int l=0;int s=0;
199  mFmsDbMkr->fpostQLSfromSlatId(slatid,&q,&l,&s);
200 
201  //cout << "Begin outside Check: SlatID="<<slatid << "|xing="<<xing << "|adc="<<adc << "|qt="<<qt << "|ch="<<ch << "|q="<<q <<",l="<<l << ",s="<<s << endl;
202 
203  if( q>0 && l>0 && s>0 && abs(xing) <= mNPREPOST )
204  {
205  mXing[xing+mNPREPOST]->Fill((float)adc);
206  if( xing == 0 || mPed != 0 )//Collision for this event
207  {
208  nfpostdata++;
209  mAdc2[0]->Fill((float)slatid,(float)adc);
210  mAdc2[1]->Fill((float)slatid,(float)adc);
211 
212  //cout << "Begin inside Check: slatID=" << slatid << "|xing="<<xing << "|adc="<<adc << "|qt="<<qt << "|ch="<<ch << "|q="<<q <<",l="<<l << ",s="<<s << endl;
213 
214  mAdc[slatid][0]->Fill((float)adc);
215  mAdc[slatid][1]->Fill((float)adc);
216  if( adc>50 )
217  {
218  nhtot++;
219  nh[q-1][l-1]++;
220  mHit[q-1][l-1]->Fill(float(s));
221  }//if( adc>50) )
222  }//if( xing == 0 )
223  }//if(q>0 && l>0 && s>0 && abs(xing)<=mNPREPOST)
224  //hits[i]->print();
225  }//if( det == 15 )
226  }//for (unsigned int i=0; i<nhit; i++)
227 
228  //cout << "Finished first for loop" << endl;
229 
230  //Filling number of hits
231  mDataSize[0]->Fill(log10(nfpostdatatot));
232  mDataSize[1]->Fill(log10(nfpostdata));
233  //Filling number of hits per quadrant and layer
234  for(int q=0; q<mNQ; q++)
235  {
236  for(int l=0; l<mNL; l++)
237  {
238  mNHit[q][l]->Fill(float(nh[q][l]));
239  }
240  }
241 
242  //cout << "Finshed second for loop" << endl;
243 
244  //total multiplicity by triggers
245  mNHitTrg[64]->Fill(float(nhtot));//cout << "mNHitTrg->Fill" << endl;
246  unsigned long long one=1;//cout << "long long one" << endl;
247  StFpsRawDaqReader* fpsraw=(StFpsRawDaqReader*)GetMaker("fpsRawDaqReader");//cout <<"pointer: " << fpsraw <<endl;
248  if( fpsraw == 0 ){return kStOK;}
249  //cout << "Checking: " << endl; cout << fpsraw->trgMask() << endl;
250  //cout << "value" << endl;
251  unsigned long long tmask = fpsraw->trgMask();//cout << "long long tmask" << endl;
252  //cout << "Beginning for loop" << endl;
253  for(int t=0; t<mNTRG; t++)
254  {
255  //cout << "Inside third for loop" << endl;
256  if( tmask & (one<<t) )
257  {
258  mNHitTrg[t]->Fill(float(nhtot));
259  mNHitTrg2->Fill(float(nhtot),float(t));
260  }
261  }
262 
263  //cout << "Finished third for loop" << endl;
264  //printf("NFMSHIT=%4d NFPSHITTOT=%4d NFPSHIT(xing=0)=%d\n",nhit,nfpsdatatot,nfpsdata);
265 
266 
267  StTriggerData* trg = fpsraw->trgdata();
268  if(trg)
269  {
270  unsigned int tcu=trg->tcuCounter();
271  unsigned int rfpo=fpsraw->rccFpost();
272  unsigned int dfpo=rfpo-tcu;
273  if(rfpo==0)
274  {
275  dfpo=0;
276  }
277  else if(rfpo < tcu)
278  {
279  const long long one=1;
280  const long long m=one<<32;
281  long long r=rfpo;
282  long long t=tcu;
283  dfpo=(unsigned int)(r+m-t);
284  }
285  mRccDiff[0]->Fill(dfpo);
286  mRccDiff[1]->Fill(dfpo);
287  //printf("FPOST RCC=%10d TCU=%10d DIFF=%10d\n",rfpo,tcu,dfpo);
288  }
289  else
290  {
291  printf("No StTriggerData found\n");
292  }
293 
294  return kStOK;
295 };
296 
298  /*
299  mDataSize[0]->Write();
300  mDataSize[1]->Write();
301  for(int i=0; i<mNPREPOST*2+1; i++) mXing[i]->Write();
302  mAdc2->Write();
303  for(int i=0; i<mNID; i++) mAdc[i]->Write();
304  */
305 
306  if( ran_fpost )//fpost was run so write the root file
307  {
308  cout << "This run has FPOST" << endl;
309  mFile->Write();
310  mFile->Close();
311  printf("StFpostQaMaker::Finish - Closing %s\n",mFilename);
312  return kStOK;
313  }
314  else
315  {
316  cout << "This run has no FPOST" << endl;
317  //This will create an empty file instead of a root file so that I don't waste time plotting
318  int yday = mRun/1000;//The usual 'yday' variable which looks like year followed by day
319  stringstream outFile_name;
320  outFile_name << "/ldaphome/dkap7827/fpost_onl_monitoring/" << yday << "/" << mRun << ".empty.txt";
321  ofstream outFile( outFile_name.str().c_str() );
322  if( !outFile.is_open() )
323  {
324  cout << "Could not open file: " << outFile_name.str() << endl;
325  exit(0);
326  }
327  outFile << endl;
328  mFile->Close();
329  outFile.close();
330  printf("StFpostQaMaker::Finish - Closing %s and %s\n",outFile_name.str().c_str(),mFilename);
331  return kStOK;
332  }
333 };
334 
335 ClassImp(StFpostQaMaker);
336 
337 /*
338  * $Id: StFpostQaMaker.cxx,v 1.1 2017/02/22 07:14:44 akio Exp $
339  * $Log: StFpostQaMaker.cxx,v $
340  * Revision 1.1 2017/02/22 07:14:44 akio
341  * Initial version from David
342  *
343  *
344  * Revision 1.1 2015/02/01 23:22:56 dkap7827
345  * new fpost qa maker
346  *
347  */
348 
virtual Int_t Make()
Definition: Stypes.h:40
virtual Int_t Finish()
Definition: Stypes.h:44