CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHeavyIon.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <istream>
4 #include <fstream>
5 #include <iomanip>
6 #include <string>
7 #include <cmath>
8 #include <functional>
9 #include <stdlib.h>
10 #include <string.h>
11 
14 
15 // L1 related
17 
19 using namespace std;
20 
22 
23  //set parameter defaults
24  _Monte = false;
25  _Debug=false;
26  _OR_BXes=false;
27  UnpackBxInEvent=1;
28 
29 }
30 
32 
33 
34  bool changed(true);
35  if (hltConfig_.init(run,c,processName_,changed)) {
36  // if init returns TRUE, initialisation has succeeded!
37  if (changed) {
38  // The HLT config has actually changed wrt the previous Run, hence rebook your
39  // histograms or do anything else dependent on the revised HLT config
40  cout << "Initalizing HLTConfigProvider" << endl;
41  }
42  } else {
43  // if init returns FALSE, initialisation has NOT succeeded, which indicates a problem
44  // with the file and/or code and needs to be investigated!
45  cout << " HLT config extraction failure with process name " << processName_ << endl;
46  // In this case, all access methods will return empty values!
47  }
48 
49 }
50 
51 /* Setup the analysis to put the branch-variables into the tree. */
52 void HLTHeavyIon::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
53 
54 
55  processName_ = pSet.getParameter<std::string>("HLTProcessName") ;
56 
57  edm::ParameterSet myHltParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
58  vector<std::string> parameterNames = myHltParams.getParameterNames() ;
59 
60  for ( vector<std::string>::iterator iParam = parameterNames.begin();
61  iParam != parameterNames.end(); iParam++ ){
62  if ( (*iParam) == "Debug" ) _Debug = myHltParams.getParameter<bool>( *iParam );
63  if ( (*iParam) == "Monte" ) _Monte = myHltParams.getParameter<bool>( *iParam );
64 
65  }
66 
67  HltEvtCnt = 0;
68  const int kMaxEvtPlanes = 1000;
69 
70  fNpart = -1;
71  fNcoll = -1;
72  fNhard = -1;
73  fPhi0 = -1;
74  fb = -1;
75  fNcharged = -1;
76  fNchargedMR = -1;
77  fMeanPt = -1;
78  fMeanPtMR = -1;
79 
80  fEtMR = -1;
81  fNchargedPtCut = -1;
82  fNchargedPtCutMR = -1;
83 
84  nEvtPlanes = 0;
85  hiBin = -1;
86  hiEvtPlane = new float[kMaxEvtPlanes];
87 
88  HltTree->Branch("Npart",&fNpart,"Npart/F");
89  HltTree->Branch("Ncoll",&fNcoll,"Ncoll/F");
90  HltTree->Branch("Nhard",&fNhard,"Nhard/F");
91  HltTree->Branch("phi0",&fPhi0,"NPhi0/F");
92  HltTree->Branch("b",&fb,"b/F");
93  HltTree->Branch("Ncharged",&fNcharged,"Ncharged/I");
94  HltTree->Branch("NchargedMR",&fNchargedMR,"NchargedMR/I");
95  HltTree->Branch("MeanPt",&fMeanPt,"MeanPt/F");
96  HltTree->Branch("MeanPtMR",&fMeanPtMR,"MeanPtMR/F");
97  HltTree->Branch("EtMR",&fEtMR,"EtMR/F");
98  HltTree->Branch("NchargedPtCut",&fNchargedPtCut,"NchargedPtCut/I");
99  HltTree->Branch("NchargedPtCutMR",&fNchargedPtCutMR,"NchargedPtCutMR/I");
100  HltTree->Branch("hiBin",&hiBin,"hiBin/I");
101  HltTree->Branch("hiHF",&hiHF,"hiHF/F");
102  HltTree->Branch("hiHFplus",&hiHFplus,"hiHFplus/F");
103  HltTree->Branch("hiHFminus",&hiHFminus,"hiHFminus/F");
104  HltTree->Branch("hiZDC",&hiZDC,"hiZDC/F");
105  HltTree->Branch("hiZDCplus",&hiZDCplus,"hiZDCplus/F");
106  HltTree->Branch("hiZDCminus",&hiZDCminus,"hiZDCminus/F");
107 
108  HltTree->Branch("hiHFhit",&hiHFhit,"hiHFhit/F");
109  HltTree->Branch("hiHFhitPlus",&hiHFhitPlus,"hiHFhitPlus/F");
110  HltTree->Branch("hiHFhitMinus",&hiHFhitMinus,"hiHFhitMinus/F");
111 
112  HltTree->Branch("hiET",&hiET,"hiET/F");
113  HltTree->Branch("hiEE",&hiEE,"hiEE/F");
114  HltTree->Branch("hiEB",&hiEB,"hiEB/F");
115  HltTree->Branch("hiEEplus",&hiEEplus,"hiEEplus/F");
116  HltTree->Branch("hiEEminus",&hiEEminus,"hiEEminus/F");
117  HltTree->Branch("hiNpix",&hiNpix,"hiNpix/I");
118  HltTree->Branch("hiNpixelTracks",&hiNpixelTracks,"hiNpixelTracks/I");
119  HltTree->Branch("hiNtracks",&hiNtracks,"hiNtracks/I");
120  HltTree->Branch("hiNevtPlane",&nEvtPlanes,"hiNevtPlane/I");
121  HltTree->Branch("hiEvtPlanes",hiEvtPlane,"hiEvtPlanes/F");
122  HltTree->Branch("hiNtracksPtCut",&hiNtracksPtCut,"hiNtracksPtCut/I");
123  HltTree->Branch("hiNtracksEtaCut",&hiNtracksEtaCut,"hiNtracksEtaCut/I");
124  HltTree->Branch("hiNtracksEtaPtCut",&hiNtracksEtaPtCut,"hiNtracksEtaPtCut/I");
125 
126 }
127 
128 /* **Analyze the event** */
130  const edm::Handle<reco::Centrality> & centrality,
131  const edm::Handle<reco::EvtPlaneCollection> & evtPlanes,
132  const edm::Handle<edm::GenHIEvent> & mc,
133  edm::EventSetup const& eventSetup,
134  edm::Event const& iEvent,
135  TTree* HltTree) {
136 
137  std::cout << " Beginning HLTHeavyIon " << std::endl;
138 
139  if(_Monte){
140  fb = mc->b();
141  fNpart = mc->Npart();
142  fNcoll = mc->Ncoll();
143  fNhard = mc->Nhard();
144  fPhi0 = mc->evtPlane();
145  fNcharged = mc->Ncharged();
146  fNchargedMR = mc->NchargedMR();
147  fMeanPt = mc->MeanPt();
148  fMeanPtMR = mc->MeanPtMR();
149  fEtMR = mc->EtMR();
150  fNchargedPtCut = mc->NchargedPtCut();
151  fNchargedPtCutMR = mc->NchargedPtCutMR();
152  }
153 
154  edm::Handle<int> binHandle;
155  iEvent.getByLabel("centralityBin",binHandle);
156  hiBin = *binHandle;
157 
158  hiNpix = centrality->multiplicityPixel();
159  hiNpixelTracks = centrality->NpixelTracks();
160  hiNtracks = centrality->Ntracks();
161  hiNtracksPtCut = centrality->NtracksPtCut();
162  hiNtracksEtaCut = centrality->NtracksEtaCut();
163  hiNtracksEtaPtCut = centrality->NtracksEtaPtCut();
164 
165  hiHF = centrality->EtHFtowerSum();
166  hiHFplus = centrality->EtHFtowerSumPlus();
167  hiHFminus = centrality->EtHFtowerSumMinus();
168  hiHFhit = centrality->EtHFhitSum();
169  hiHFhitPlus = centrality->EtHFhitSumPlus();
170  hiHFhitMinus = centrality->EtHFhitSumMinus();
171 
172  hiZDC = centrality->zdcSum();
173  hiZDCplus = centrality->zdcSumPlus();
174  hiZDCminus = centrality->zdcSumMinus();
175 
176  hiEEplus = centrality->EtEESumPlus();
177  hiEEminus = centrality->EtEESumMinus();
178  hiEE = centrality->EtEESum();
179  hiEB = centrality->EtEBSum();
180  hiET = centrality->EtMidRapiditySum();
181 
182 
183  if(evtPlanes.isValid()){
184  nEvtPlanes = evtPlanes->size();
185  for(unsigned int i = 0; i < evtPlanes->size(); ++i){
186  hiEvtPlane[i] = (*evtPlanes)[i].angle();
187  }
188  }
189 
190 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void analyze(const edm::Handle< edm::TriggerResults > &hltresults, const edm::Handle< reco::Centrality > &centrality, const edm::Handle< reco::EvtPlaneCollection > &evtPlanes, const edm::Handle< edm::GenHIEvent > &hiMC, edm::EventSetup const &eventSetup, edm::Event const &iEvent, TTree *tree)
Definition: HLTHeavyIon.cc:129
void beginRun(const edm::Run &, const edm::EventSetup &)
Definition: HLTHeavyIon.cc:31
int iEvent
Definition: GenABIO.cc:243
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
void setup(const edm::ParameterSet &pSet, TTree *tree)
Definition: HLTHeavyIon.cc:52
tuple cout
Definition: gather_cfg.py:121
Definition: Run.h:36