StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEmcOfflineCalibrationMipAnalysis.cxx
1 /*
2  * StEmcOfflineCalibrationMipAnalysis.cxx
3  * J. Kevin Adkins, University of Kentucky
4  * June 18, 2014
5  *
6  * Based on the previous mip_histogram_maker.C
7  * macro written by Matt Walker. Converted into
8  * compiled form to increase code stability when
9  * running over files and to run over new version
10  * of StEmcOfflineCalibrationEvent classes
11  */
12 
13 // User defined includes
14 #include "StEmcOfflineCalibrationMipAnalysis.h"
15 #include "StEmcOfflineCalibrationEvent.h"
16 
17 // ROOT includes
18 #include "TH1D.h"
19 #include "TH2F.h"
20 #include "TChain.h"
21 #include "TFile.h"
22 #include "TString.h"
23 
24 // STD library includes
25 #include <set>
26 #include <map>
27 
29 
30 StEmcOfflineCalibrationMipAnalysis::StEmcOfflineCalibrationMipAnalysis(const char *name, const char* outfile, TChain *calibChain):StMaker(name),nTowers(4800){
31  mCalibChain = calibChain;
32  mOutfileName = (TString)outfile;
33  mEvent = 0;
34  mVertex = 0;
35  mTrack = 0;
36  mFile = 0;
37 }
38 
39 StEmcOfflineCalibrationMipAnalysis::~StEmcOfflineCalibrationMipAnalysis(){
40  if(mCalibChain) delete mCalibChain;
41  if(mEvent) delete mEvent;
42  if(mVertex) delete mVertex;
43  if(mTrack) delete mTrack;
44  if (mFile) delete mFile;
45 }
46 
47 Int_t StEmcOfflineCalibrationMipAnalysis::Init(){
48  mFile = new TFile(mOutfileName, "RECREATE");
49  mCalibChain->SetBranchAddress("event_branch",&mEvent);
50 
51  Char_t name[100];
52  for (Int_t iTow = 0; iTow < nTowers; ++iTow){
53  sprintf(name,"tower_histo_%i",iTow+1);
54  towerHisto[iTow] = new TH1D(name,name,250,-50.5,199.5);
55  towerHisto[iTow]->Sumw2();
56  }
57 
58  mapcheck = new TH2F("mapcheck","check mapping",4800,0.5,4800.5,4800,0.5,4800.5);
59  mapcheck->GetXaxis()->SetTitle("Track Projection ID");
60  mapcheck->GetYaxis()->SetTitle("Tower hit above 5 rms");
61 
62  return StMaker::Init();
63 }
64 
66  pedSubAdc = 0.;
67  assert(mEvent);
68 
69  trackTowers.clear();
70  excludedTowers.clear();
71 
72  // Loop over tracks to get excluded towers
73  for (Int_t iTrk = 0; iTrk < mEvent->nTracks(); ++iTrk){
74  mTrack = mEvent->track(iTrk);
75  Int_t softId = mTrack->towerId(0);
76 
77  if(trackTowers.find(softId) != trackTowers.end())
78  excludedTowers.insert(softId);
79  else
80  trackTowers.insert(softId);
81  }
82 
83  // Begin loop over vertices
84  for(Int_t iVert = 0; iVert < mEvent->nVertices(); ++iVert){
85  mVertex = mEvent->vertex(iVert);
86  assert(mVertex);
87  if (mVertex->ranking() < 1e6) continue;
88  if (TMath::Abs(mVertex->z()) > 30) continue;
89 
90  // Begin loop over tracks at vertex
91  for(Int_t iTrack = 0; iTrack < mVertex->nTracks(); ++iTrack){
92  mTrack = mVertex->track(iTrack);
93  assert(mTrack);
94  pedSubAdc = mTrack->towerAdc(0) - mTrack->towerPedestal(0);
95 
96  if (excludedTowers.find(mTrack->towerId(0)) != excludedTowers.end()) continue;
97  if (mTrack->p() < 1.) continue;
98  if (mTrack->towerId(0) != mTrack->towerExitId()) continue;
99 
100  for(Int_t k = 0; k < 9; k++){
101  if(mTrack->towerAdc(k) - mTrack->towerPedestal(k) < 5*mTrack->towerPedestalRms(k)) continue;
102  mapcheck->Fill(mTrack->towerId(0),mTrack->towerId(k));
103  }
104 
105  if(mTrack->highestNeighbor() > 2.) continue;
106  if(pedSubAdc < 1.5*mTrack->towerPedestalRms(0)) continue;
107 
108  Int_t index = mTrack->towerId(0) - 1;
109  towerHisto[index]->Fill(pedSubAdc);
110  }// Tracks loop
111  }// Vertex Loop
112 
113  return kStOK;
114 }
115 
117 {
118  mFile->Write();
119  mFile->Close();
120  delete mFile;
121  return kStOk;
122 }
Definition: Stypes.h:40
Definition: Stypes.h:41