StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnalyzePedAsPhysRun.C
1 TString MikesMinimalOpts = " P2017a PicoVtxDefault -hitfilt ";
2 
3 
4 TString DanielSuggestedOpts = "ry2017a in AgML tpcDB btof Tree picoWrite UseXgeom BAna ppOpt VFPPVnoCTB beamline3D l3onl fpd trgd analysis PicoVtxDefault -hitfilt epdDb epdHit";
5 
6 
7 class StSPtrVecEpdHit;
8 
9 void DrawWheel(StSPtrVecEpdHit hits, TString outputFileName);
10 
11 void AnalyzePedAsPhysRun(char* inputDirectory, char* inputFileName)
12 {
13 
14  TString _infile = inputDirectory;
15  _infile += "/";
16  _infile += inputFileName;
17 
18  TString outfile = inputFileName;
19  outfile.ReplaceAll(".daq","_Levels.txt");
20 
21 
22  gROOT->Macro("Load.C");
23  gSystem->Load("StEpdUtil");
24 
25  gROOT->LoadMacro("bfc.C");
26 
27  // TString mychain = " epdDb epdHit "; //"epdDb epdHit";
28  // bfc(0,MikesMinimalOpts + mychain,_infile);
29  bfc(0,DanielSuggestedOpts,_infile);
30 
31  const int N = 2; // process N events
32  StSPtrVecEpdHit hits;
33  for ( int i=0; i<N; i++ ) {
34 
35  chain->Clear();
36  chain->Make(); // in principle should check return status...
37 
38  StEvent* event = (StEvent*) chain->GetInputDS("StEvent");
39  StMuDst* muDst = (StMuDst*) chain->GetInputDS("muDst");
40  if (!muDst){
41  cout << "Cant find muDst :-(\n";
42  return;}
43 
44  cout << "\n\n\n\n\n --------- here in the macro -----------\n\n\n i= " << i << "\n\n";
45  cout << "Looking at EPD hits in StEvent....\n\n";
46  if (event->epdCollection() == 0){ cout << "No StEpdCollection!\n";}
47  else {
48  StEpdCollection* epdCollection = event->epdCollection();
49  hits = epdCollection->epdHits();
50  // StPtrVecEpdHitIterator it = hits->begin(); I'd rather avoid iterators...
51  cout << "StEpdCollection size is " << hits.size() << endl;
52  }
53  cout << "\n\n\n\n\n\n ----------------------------------\n\n\n";
54  }
55  DrawWheel(hits, outfile);
56 
57  // new TBrowser();
58 
59 };
60 
61 
62 
63 // Here, we'll loop over all StEpdHits
64 // Every hit will have its outline drawn
65 void DrawWheel(StSPtrVecEpdHit hits, TString outputFileName){
66 
67 
68  int lowLimitB[4] = {-1,10,150,3001};
69  int hiLimitB[4] = {9,149,3000,5000};
70  int lowLimitC[4] = {-1,10,150,3001};
71  int hiLimitC[4] = {9,149,3000,5000};
72 
73 
74  int levelColor[4] = {4,3,2,6};
75 
76 
77  int tileLevel[12][31][2];
78  int tileAdc[12][31][2];
79 
80  for (int ew=0; ew<2; ew++){
81  for (int pp=1; pp<13; pp++){
82  for (int tt=1; tt<32; tt++){
83  tileLevel[pp-1][tt-1][ew] = 0;
84  tileAdc[pp-1][tt-1][ew] = 0;
85  }
86  }
87  }
88 
89  for (unsigned int tile=0; tile<hits.size(); tile++){
90  StEpdHit* theHit = hits[tile];
91  int pp = theHit->position();
92  int tt = theHit->tile();
93  int ew = theHit->side();
94  int ewIndex = (ew<0)?0:1;
95  int adc = theHit->adc();
96  int Level,lowCutoff,hiCutoff;
97  for (Level=0; Level<4; Level++){
98  if (tt<10){lowCutoff=lowLimitC[Level]; hiCutoff=hiLimitC[Level];}
99  else {lowCutoff=lowLimitB[Level]; hiCutoff=hiLimitB[Level];}
100  if ((adc>=lowCutoff)&&(adc<=hiCutoff)) break;
101  }
102  if (Level>3) cout << "ERROR!!\n\n";
103  tileLevel[pp-1][tt-1][ewIndex] = Level;
104  tileAdc[pp-1][tt-1][ewIndex] = adc;
105  }
106 
107  ofstream ofs;
108  ofs.open(outputFileName.Data());
109 
110  for (int ew=0; ew<2; ew++){
111  int ewSigned = (ew==0)?-1:1;
112  for (int pp=1; pp<13; pp++){
113  for (int tt=1; tt<32; tt++){
114  ofs << ewSigned << "\t" << pp << "\t" << tt << "\t" << tileLevel[pp-1][tt-1][ew] << "\t" << tileAdc[pp-1][tt-1][ew] << endl;
115  }
116  }
117  }
118 
119  ofs.close();
120 
121 
122 
123  StEpdGeom* geo = new StEpdGeom;
124  TCanvas* WheelCan = new TCanvas("EpdWheels","EpdWheels",1200,600);
125  WheelCan->Draw();
126  WheelCan->Divide(2,1);
127 
128  TPad* eastPad = WheelCan->cd(1);
129  eastPad->Range(-100,-100,100,100);
130  TPaveLabel* eastLab = new TPaveLabel(60,90,90,99,"East");
131  eastLab->Draw();
132  double labelRadius=96;
133  for (int pp=1; pp<13; pp++){
134  double phi=geo->TileCenter(pp,1,-1).Phi(); // note third argument - means "east"
135  TText* ppLabel = new TText(labelRadius*cos(phi),labelRadius*sin(phi),Form("PP%d",pp));
136  ppLabel->SetTextAlign(22);
137  ppLabel->SetTextColor(kBlue);
138  ppLabel->SetTextSize(.02);
139  ppLabel->Draw();
140  }
141 
142  TPad* westPad = WheelCan->cd(2);
143  westPad->Range(-100,-100,100,100);
144  TPaveLabel* westLab = new TPaveLabel(60,90,90,99,"West");
145  westLab->Draw();
146  for (int pp=1; pp<13; pp++){
147  double phi=geo->TileCenter(pp,1,1).Phi(); // note third argument - means "west"
148  TText* ppLabel = new TText(labelRadius*cos(phi),labelRadius*sin(phi),Form("PP%d",pp));
149  ppLabel->SetTextAlign(22);
150  ppLabel->SetTextColor(kBlue);
151  ppLabel->SetTextSize(.02);
152  ppLabel->Draw();
153  }
154 
155 
156  int nCorners;
157  double xcorn[6];
158  double ycorn[6];
159  TPad* thePad;
160 
161  for (int ew=0; ew<2; ew++){
162  for (int pp=1; pp<13; pp++){
163  for (int tt=1; tt<32; tt++){
164  int id = 100*pp+tt;
165  if (ew==0) id *= -1;
166  geo->GetCorners(id,&nCorners,xcorn,ycorn);
167  xcorn[nCorners]=xcorn[0]; ycorn[nCorners]=ycorn[0]; // need to close the polylines
168  TPolyLine* pline = new TPolyLine(nCorners+1,xcorn,ycorn);
169  pline->SetLineWidth(1);
170  pline->SetLineColor(1);
171  pline->SetFillColor(levelColor[tileLevel[pp-1][tt-1][ew]]);
172  thePad = (ew>0)?westPad:eastPad;
173  thePad->cd();
174  pline->Draw();
175  pline->Draw("f");
176  TText* ttlabel = new TText(geo->TileCenter(id).X(),geo->TileCenter(id).Y(),Form("%d",tt));
177  ttlabel->SetTextAlign(22);
178  ttlabel->SetTextSize(0.015);
179  ttlabel->Draw();
180  }
181  }
182  }
183 
184  TString imageName = outputFileName + ".png";
185  WheelCan->SaveAs(imageName.Data());
186 
187  ofs.close();
188 
189 
190 }
191 
192 
193 /*
194  Double_t x[5] = {-20,-10,0,10,-20};
195  Double_t y[5] = {10,40,50,5,10};
196  TPolyLine *pline = new TPolyLine(5,x,y);
197  pline->SetFillColor(38);
198  pline->SetLineColor(2);
199  pline->SetLineWidth(4);
200  pline->Draw("f");
201  pline->Draw();
202 
203 */
204 
205  /*
206  cout << "PP\tTT\tADC\tTDC\tTAC\tnMIP\n";
207  for (unsigned int tile=0; tile<hits.size(); tile++){
208  StEpdHit* theHit = hits[tile];
209  cout << theHit->position() << "\t"
210  << theHit->tile() << "\t"
211  << theHit->adc() << "\t"
212  << theHit->tdc() << "\t"
213  << theHit->tac() << "\t"
214  << theHit->nMIP() << "\n";
215  }
216  */
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
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
virtual Int_t Make()
Definition: StChain.cxx:110
Stores information for tiles in STAR Event Plane Detector.
Definition: StEpdHit.h:43
int tile() const
tile on the supersector [1,31]
Definition: StEpdHit.h:148