StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFmsEventDisplay.cxx
1 // \class StFmsEventDisplay
2 // \author Akio Ogawa
3 //
4 // $Id: StFmsEventDisplay.cxx,v 1.6 2017/05/26 18:55:53 akio Exp $
5 // $Log: StFmsEventDisplay.cxx,v $
6 // Revision 1.6 2017/05/26 18:55:53 akio
7 // added protection for a crash in MC mudst file reading
8 //
9 // Revision 1.5 2016/06/08 16:31:50 akio
10 // c++11 style initialization
11 //
12 // Revision 1.4 2016/01/20 19:56:39 akio
13 // *** empty log message ***
14 //
15 // Revision 1.2 2016/01/20 19:46:40 akio
16 // *** empty log message ***
17 //
18 // Revision 1.1 2015/10/20 19:55:51 akio
19 // Initial version of FMS event display
20 //
21 //
22 
23 #include "StFmsEventDisplay.h"
24 
25 #include "StMessMgr.h"
26 #include "Stypes.h"
27 #include "StMuDSTMaker/COMMON/StMuTypes.hh"
28 
29 #include "StThreeVectorF.hh"
30 #include "StFmsDbMaker/StFmsDbMaker.h"
31 #include "StEnumerations.h"
32 #include "StEventTypes.h"
33 #include "StEvent/StEvent.h"
34 #include "StEvent/StFmsCollection.h"
35 #include "StEvent/StFmsHit.h"
36 #include "StEvent/StFmsPoint.h"
37 #include "StEvent/StFmsPointPair.h"
38 
39 #include "TApplication.h"
40 #include "TCanvas.h"
41 #include "TFile.h"
42 #include "TH1F.h"
43 #include "TH2F.h"
44 #include "TLine.h"
45 #include "TBox.h"
46 #include "TMarker.h"
47 #include "TString.h"
48 #include "TColor.h"
49 #include "TStyle.h"
50 #include "TText.h"
51 #include "TROOT.h"
52 
53 ClassImp(StFmsEventDisplay);
54 
55 StFmsEventDisplay::StFmsEventDisplay(const Char_t* name):
56  StMaker(name),mFilename((char *)"fmsEventDisplay.pdf"){}
57 
58 StFmsEventDisplay::~StFmsEventDisplay(){}
59 
60 Int_t StFmsEventDisplay::Init(){
61  mFmsDbMaker=static_cast<StFmsDbMaker*>(GetMaker("fmsDb"));
62  if(!mFmsDbMaker){
63  LOG_ERROR << "StFmsEventDisplay::InitRun Failed to get StFmsDbMaker" << endm;
64  return kStFatal;
65  }
66 
67  mApplication = new TApplication("EvtDsp",0,0);
68  mCanvas=new TCanvas("FMSEventtDisplay","FMSEventtDisplay",10,10,900,900);
69  gStyle->SetOptStat(0);
70  return kStOK;
71 }
72 
74  mApplication->Terminate();
75  return kStOK;
76 }
77 
78 void setColor(int id, float v){ //blue->green for 0->0.5 and green->red for 0.5->1.0
79  TColor *color = gROOT->GetColor(id);
80  if(!color) color=new TColor(id,0,0,0);
81  float r=0.0,g=0.0,b=0.0;
82  if(v>1.0) v=1.0;
83  if(v<0.0) v=0.0;
84  if(v>0.5){
85  r=(v-0.5)*2.0;
86  g=1.0-(v-0.5)*2.0;
87  }else{
88  g=v*2.0;
89  b=1.0-(v*2.0);
90  }
91  color->SetRGB(r,g,b);
92  //printf("%d %f %f %f %f\n",id,v,r,g,b);
93 }
94 
95 void scale(float x=100.0, float dx=5.0, float ymin=-100.0, float ymax=100.0, int ny=40){
96  for(int i=0; i<ny; i++){
97  setColor(2000+i,float(i+0.5)/float(ny));
98  float dy=(ymax-ymin)/float(ny);
99  TBox* b=new TBox(x,ymin+dy*i,x+dx,ymin+dy*(i+1));
100  b->SetFillColor(2000+i);
101  b->Draw();
102  }
103 }
104 
106  int bunch=-1;
107  StMuDst* mudst = (StMuDst*)GetInputDS("MuDst");
108  if(mudst) {
109  if(mudst->event()->triggerData()){
110  bunch = mudst->event()->triggerData()->bunchId7Bit();
111  LOG_INFO << Form("Bunch=%3d from Mudst\n",bunch) << endm;
112  }
113  }
114 
115  StEvent* event = (StEvent*)GetInputDS("StEvent");
116  if(!event) {LOG_ERROR << "StFmsEventDisplay::Make did not find StEvent"<<endm; return kStErr;}
117  mFmsColl = event->fmsCollection();
118  if(!mFmsColl) {LOG_ERROR << "StFmsEventDisplay::Make did not find StEvent->FmsCollection"<<endm; return kStErr;}
119  if(bunch==-1){
120  if(event->triggerData()) {
121  bunch = event->triggerData()->bunchId7Bit();
122  LOG_INFO << Form("Bunch=%3d from StEvent\n",bunch) <<endm;
123  }
124  }
125 
126  StSPtrVecFmsHit& hits = mFmsColl->hits();
127  StSPtrVecFmsCluster& clusters = mFmsColl->clusters();
128  //StSPtrVecFmsPoint& points = mFmsColl->points();
129  //StSPtrVecFpsSlat& slats = mFmsColl->fpsSlats();
130  //vector<StFmsPointPair*>& pairs = mFmsColl->pointPairs();
131  int nh=mFmsColl->numberOfHits();
132  int nc=mFmsColl->numberOfClusters();
133  int np=mFmsColl->numberOfPoints();
134  //int ns=slats.size();
135  //int npair=mFmsColl->numberOfPointPairs();
136 
137  if(mNAccepted < mMaxEvents){
138  //filters
139  if(mFilter==1 && np<2) {mNEvents++; return kStOK;}
140  if(mFilter==2 && (bunch<30 || bunch>=40)) {mNEvents++; return kStOK;}
141 
142  char cc[100];
143  sprintf(cc,"FMSEventDisplayEvt=%d",mNEvents);
144  mCanvas->Clear();
145  //frame
146  TH2F* frame=new TH2F(cc,cc,1,-100.0,105.0,1,-100.0,100.0);
147  frame->Draw();
148  scale();
149  TText* t1=new TText(106, 100.0,"100GeV"); t1->SetTextSize(0.03); t1->Draw();
150  TText* t2=new TText(106, 33.3, "10GeV"); t2->SetTextSize(0.03); t2->Draw();
151  TText* t3=new TText(106,-33.3, "1GeV"); t3->SetTextSize(0.03); t3->Draw();
152  TText* t4=new TText(106,-100.0,"0.1GeV"); t4->SetTextSize(0.03); t4->Draw();
153 
154  for(int det=kFmsNorthLargeDetId; det<=kFmsSouthSmallDetId; det++){
155  for(int ch=0; ch<mFmsDbMaker->maxChannel(det); ch++){
156  if(mFmsDbMaker->getGain(det,ch)>0.0){
157  StThreeVectorF xyz = mFmsDbMaker->getStarXYZ(det,ch);
158  float wx=mFmsDbMaker->getXWidth(det);
159  float wy=mFmsDbMaker->getYWidth(det);
160  TBox* cell=new TBox(xyz.x()-wx/2.0, xyz.y()-wy/2.0, xyz.x()+wx/2.0, xyz.y()+wy/2.0);
161  cell->SetFillColor(0);
162  cell->SetLineColor(1); cell->SetLineWidth(1);
163  cell->Draw("l");
164  }
165  }
166  }
167 
168  //First draw hits
169  for(int i=0; i<nh; i++){
170  StFmsHit* hit=hits[i];
171  int det=hit->detectorId();
172  if(det>=kFmsNorthLargeDetId && det<=kFmsSouthSmallDetId){
173  //int ch=hit->channel();
174  StThreeVectorF xyz = mFmsDbMaker->getStarXYZ(hit);
175  float wx=mFmsDbMaker->getXWidth(det);
176  float wy=mFmsDbMaker->getYWidth(det);
177  TBox* cell=new TBox(xyz.x()-wx/2.0, xyz.y()-wy/2.0, xyz.x()+wx/2.0, xyz.y()+wy/2.0);
178  int color = 1000+i;
179  float ladc=log10(hit->energy());
180  if(ladc>2.0) ladc=2.0;
181  if(ladc<-1.0) ladc=-1.0;
182  setColor(color,(ladc+1.0)/3.0);
183  cell->SetFillColor(color);
184  cell->Draw();
185  }
186  }
187 
188  for(int i=0; i<nc; i++){
189  StFmsCluster* clu=clusters[i];
190  float nt=float(clu->nTowers());
191  for(int j=0; j<nt; j++){
192  StFmsHit* hit = clu->hits()[j];
193  int det=hit->detectorId();
194  StThreeVectorF xyz = mFmsDbMaker->getStarXYZ(hit);
195  float wx=mFmsDbMaker->getXWidth(det);
196  float wy=mFmsDbMaker->getYWidth(det);
197  TBox* cell=new TBox(xyz.x()-wx/2.0, xyz.y()-wy/2.0, xyz.x()+wx/2.0, xyz.y()+wy/2.0);
198  cell->SetFillStyle(0);
199  cell->SetLineColor(2+i);
200  cell->SetLineWidth(2);
201  cell->Draw("l");
202  }
203  int nphoton=clu->nPhotons();
204  for(int j=0; j<nphoton; j++){
205  StFmsPoint* point = clu->points()[j];
206  TMarker* m=new TMarker(point->XYZ().x(),point->XYZ().y(),29);
207  m->SetMarkerColor(1);
208  m->SetMarkerSize(2);
209  m->Draw();
210  }
211  }
212 
213  mCanvas->Update();
214  TString f(mFilename);
215  if(f.Contains(".pdf")){ //not working
216  if(mNEvents==0) f.Append("(");
217  else if(mNEvents==mMaxEvents-1) f.Append(")");
218  LOG_INFO << "Saving " << f.Data() << endm;
219  mCanvas->Print(f.Data(),"pdf");
220  }else if(f.Contains(".png")){
221  f.ReplaceAll("eventDisplay.png",Form("%d.eventDisplay.png",mNEvents));
222  LOG_INFO << "Saving " << f.Data() << endm;
223  mCanvas->SaveAs(f.Data());
224  }
225  mNAccepted++;
226  mNEvents++;
227  //std::cin.get();
228  }
229  return kStOK;
230 }
Float_t getXWidth(Int_t detectorId)
get the offset of the detector
StThreeVectorF getStarXYZ(Int_t detectorId, Float_t FmsX, Float_t FmsY)
get the Y width of the cell
UShort_t maxChannel(Int_t detectorId) const
number of column
Float_t getYWidth(Int_t detectorId)
get the X width of the cell
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
Definition: Stypes.h:40
Definition: Stypes.h:44