ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MPIAnalyzer/src/MPIntuple.cc
Revision: 1.18
Committed: Wed Oct 27 14:21:20 2010 UTC (14 years, 6 months ago) by naodell
Content type: text/plain
Branch: MAIN
Changes since 1.17: +35 -16 lines
Log Message:
Include generator level info, i.e. ptHat values

File Contents

# User Rev Content
1 naodell 1.1 // Original Author: "Nathaniel Odell"
2     // Secondary Author: Steven Won
3 andrey 1.2 // With contributions from: Andrey Pozdnyakov
4 naodell 1.1 // Created: Thurs April 22 2010
5 naodell 1.18 // $Id: MPIntuple.cc,v 1.17 2010/10/25 09:17:02 naodell Exp $
6 naodell 1.1
7     // system include files
8     #include <memory>
9     #include <string>
10    
11     // user include files
12     #include "FWCore/Framework/interface/Frameworkfwd.h"
13     #include "FWCore/Framework/interface/EDAnalyzer.h"
14     #include "FWCore/Framework/interface/ESHandle.h"
15     #include "FWCore/Framework/interface/EventSetup.h"
16     #include "FWCore/Framework/interface/Event.h"
17     #include "FWCore/Framework/interface/MakerMacros.h"
18     #include "FWCore/ParameterSet/interface/ParameterSet.h"
19     #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
20     #include "Geometry/Records/interface/CaloGeometryRecord.h"
21 naodell 1.10 #include "DataFormats/Math/interface/deltaPhi.h"
22 naodell 1.1
23 naodell 1.15 // Libraries for objects
24 naodell 1.1 #include "DataFormats/JetReco/interface/CaloJet.h"
25     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
26     #include "DataFormats/JetReco/interface/PFJet.h"
27     #include "DataFormats/JetReco/interface/PFJetCollection.h"
28     #include "DataFormats/JetReco/interface/GenJet.h"
29     #include "DataFormats/JetReco/interface/GenJetCollection.h"
30     #include "DataFormats/JetReco/interface/Jet.h"
31     #include "DataFormats/VertexReco/interface/Vertex.h"
32     #include "DataFormats/VertexReco/interface/VertexFwd.h"
33 naodell 1.7 #include "DataFormats/BTauReco/interface/JetTag.h"
34 naodell 1.15 #include "DataFormats/TrackReco/interface/Track.h"
35 naodell 1.16 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
36 naodell 1.18 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
37    
38 devildog 1.11 //#include "RecoVertex/PrimaryVertexProducer/interface/VertexHigherPtSquared.h"
39 naodell 1.7
40 naodell 1.10 //GenParticles
41 naodell 1.13 //#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
42 naodell 1.1
43 naodell 1.4 // Need to get correctors
44     #include "JetMETCorrections/Objects/interface/JetCorrector.h"
45    
46     // ntuple storage classes
47     #include "TCJet.h"
48     #include "TCPrimaryVtx.h"
49    
50 andrey 1.2 // Need for HLT trigger info:
51     #include "FWCore/Common/interface/TriggerNames.h"
52     #include "DataFormats/Common/interface/TriggerResults.h"
53     #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
54    
55 naodell 1.4 //Root stuff
56 naodell 1.1 #include "TROOT.h"
57 naodell 1.15 #include "TH1.h"
58 naodell 1.16 #include "TH2.h"
59 naodell 1.15 #include "TProfile.h"
60 naodell 1.1 #include "TFile.h"
61     #include "TTree.h"
62     #include "TString.h"
63     #include "TObject.h"
64     #include "TObjArray.h"
65     #include "TClonesArray.h"
66     #include "TRefArray.h"
67     #include "TLorentzVector.h"
68     #include "TVector3.h"
69    
70     using namespace edm;
71     using namespace std;
72     using namespace reco;
73    
74     //
75     // class declaration
76     //
77    
78     class MPIntuple : public edm::EDAnalyzer {
79     public:
80     explicit MPIntuple(const edm::ParameterSet&);
81     ~MPIntuple();
82    
83    
84     private:
85     virtual void beginJob() ;
86     virtual void analyze(const edm::Event&, const edm::EventSetup&);
87     virtual void endJob() ;
88    
89 andrey 1.2 bool triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger);
90 devildog 1.11 double sumPtSquared(const Vertex & v);
91 andrey 1.2
92 naodell 1.1 // ----------member data ---------------------------
93    
94     edm::InputTag PFJetHandle_;
95     edm::InputTag GenJetHandle_;
96     edm::InputTag PrimaryVtxHandle_;
97 naodell 1.4 edm::InputTag triggerResultsTag_;
98 naodell 1.1 // Counting variables //
99    
100 naodell 1.8 int eventNumber, runNumber, lumiSection, bunchCross;
101 naodell 1.18 double ptHat, qScale;
102    
103 naodell 1.6 TTree* sTree;
104 naodell 1.1 TFile* ntupleFile;
105    
106 naodell 1.4 TClonesArray* jet_ak5PF;
107 naodell 1.1 TClonesArray* jetP4_ak5Gen;
108 naodell 1.5 TClonesArray* primaryVtx;
109 naodell 1.1
110     bool doGenJets_;
111     bool doPFJets_;
112 naodell 1.7 bool triggerHLT_;
113 naodell 1.8 bool isRealData;
114 naodell 1.1 string rootfilename;
115    
116 andrey 1.2 //Triggers
117 naodell 1.4 string hlTriggerResults_, hltName_, triggerName_;
118     TriggerNames triggerNames;
119 andrey 1.2 HLTConfigProvider hltConfig_;
120 naodell 1.4 vector<string> hlNames;
121 naodell 1.5 unsigned int triggerStatus;
122 andrey 1.2
123 naodell 1.15 //Histograms
124 naodell 1.16 TH1F * h1_nAssociatedTracks;
125 naodell 1.15 TH1F * h1_trackZ;
126 naodell 1.16 TH1F * h1_trackDxy;
127     TH1F * h1_allTrackDeltaZ_PV;
128     TH1F * h1_allTrackDeltaZ_LJ;
129 naodell 1.18 TH1F * h1_ptHat;
130 naodell 1.16 TH2F * h2_nAssTracksVsJetPt;
131 naodell 1.15 TProfile * p1_nVtcs;
132    
133    
134 naodell 1.1 };
135    
136 naodell 1.18 MPIntuple::MPIntuple(const edm::ParameterSet& iConfig)
137 naodell 1.1 {
138 naodell 1.18 PFJetHandle_ = iConfig.getUntrackedParameter<edm::InputTag>("PFJetTag");
139     GenJetHandle_ = iConfig.getUntrackedParameter<edm::InputTag>("GenJetTag");
140     PrimaryVtxHandle_ = iConfig.getUntrackedParameter<edm::InputTag>("PrimaryVtxTag");
141     doGenJets_ = iConfig.getUntrackedParameter<bool>("doGenJets");
142     doPFJets_ = iConfig.getUntrackedParameter<bool>("doPFJets");
143     triggerHLT_ = iConfig.getUntrackedParameter<bool>("triggerHLT");
144     rootfilename = iConfig.getUntrackedParameter<string>("rootfilename");
145     hlTriggerResults_ = iConfig.getUntrackedParameter<string>("HLTriggerResults","TriggerResults");
146     hltName_ = iConfig.getUntrackedParameter<string>("hltName","RECO");
147 naodell 1.1 }
148    
149     MPIntuple::~MPIntuple()
150     {
151    
152     }
153    
154     //
155     // member functions
156     //
157    
158     // ------------ method called to for each event ------------
159     void MPIntuple::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
160     {
161    
162 naodell 1.6 eventNumber = iEvent.id().event();
163     runNumber = iEvent.id().run();
164 naodell 1.1 lumiSection = (unsigned int)iEvent.getLuminosityBlock().luminosityBlock();
165 naodell 1.8 bunchCross = (unsigned int)iEvent.bunchCrossing();
166     isRealData = iEvent.isRealData();
167    
168 naodell 1.18
169 naodell 1.16 edm::Handle<reco::BeamSpot> beamSpotHandle;
170     iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
171     reco::BeamSpot vertexBeamSpot = *beamSpotHandle;
172    
173     int pfCount = 0;
174     int genCount = 0;
175     int vtxCount = 0;
176     float ljVertexZ = -999;
177     float primaryVertexZ = -999;
178    
179     //////////////////////////
180     //Get vertex information//
181     //////////////////////////
182    
183     Handle<reco::VertexCollection> primaryVtcs;
184     iEvent.getByLabel(PrimaryVtxHandle_, primaryVtcs);
185    
186     for(VertexCollection::const_iterator vtx_iter = primaryVtcs->begin(); vtx_iter!= primaryVtcs->end(); ++vtx_iter){
187    
188     reco::Vertex myVtx = reco::Vertex(*vtx_iter);
189     if(!myVtx.isValid() || myVtx.isFake()) continue;
190     TCPrimaryVtx* vtxCon = new ((*primaryVtx)[vtxCount]) TCPrimaryVtx;
191    
192     vtxCon->SetPosition(myVtx.x(), myVtx.y(), myVtx.z());
193     vtxCon->SetNDof(myVtx.ndof());
194     vtxCon->SetChi2(myVtx.chi2());
195     vtxCon->SetNTrks(myVtx.tracksSize());
196     vtxCon->SetSumPt2Trks(sumPtSquared(myVtx));
197    
198     if(vtxCount == 0) primaryVertexZ = myVtx.z();
199    
200     ++vtxCount;
201    
202     }
203    
204     p1_nVtcs->Fill(runNumber, vtxCount);
205    
206     ///////////////////////
207     //get jet information//
208     ///////////////////////
209 naodell 1.1
210     if(doPFJets_){
211 naodell 1.4
212 naodell 1.7 edm::Handle<reco::JetTagCollection> bTagHandle1;
213     iEvent.getByLabel("trackCountingHighEffBJetTags", bTagHandle1);
214     const reco::JetTagCollection & bTags1 = *(bTagHandle1.product());
215     reco::JetTagCollection::const_iterator jet_it_1;
216    
217     edm::Handle<reco::JetTagCollection> bTagHandle2;
218     iEvent.getByLabel("trackCountingHighPurBJetTags", bTagHandle2);
219     const reco::JetTagCollection & bTags2 = *(bTagHandle2.product());
220     reco::JetTagCollection::const_iterator jet_it_2;
221    
222 naodell 1.14
223 naodell 1.4 const JetCorrector* correctorL2 = JetCorrector::getJetCorrector ("ak5PFL2Relative",iSetup);
224     const JetCorrector* correctorL3 = JetCorrector::getJetCorrector ("ak5PFL3Absolute",iSetup);
225    
226 naodell 1.1 Handle<reco::PFJetCollection> PFJets;
227     iEvent.getByLabel(PFJetHandle_, PFJets);
228 naodell 1.16
229 naodell 1.1 for(PFJetCollection::const_iterator jet_iter = PFJets->begin(); jet_iter!= PFJets->end(); ++jet_iter){
230    
231     reco::PFJet myJet = reco::PFJet(*jet_iter);
232     if(myJet.et() > 5){
233 naodell 1.7
234     for(jet_it_1 = bTags1.begin(); jet_it_1 != bTags1.end(); jet_it_1++) {
235     if( (fabs(jet_it_1->first->eta() - myJet.eta()) < .005) && (deltaPhi(jet_it_1->first->phi(),myJet.phi()) < .005)) break;
236     }
237    
238     for(jet_it_2 = bTags2.begin(); jet_it_2 != bTags2.end(); jet_it_2++) {
239     if( (fabs(jet_it_2->first->eta() - myJet.eta()) < .005) && (deltaPhi(jet_it_2->first->phi(),myJet.phi()) < .005)) break;
240     }
241 naodell 1.6
242 naodell 1.16 TCJet* jetCon = new ((*jet_ak5PF)[pfCount]) TCJet;
243 naodell 1.6
244     jetCon->SetP4(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
245 naodell 1.15 jetCon->SetVtx(-999.0, -999.0, -999.0);
246 naodell 1.4
247 naodell 1.6 jetCon->SetChHadFrac(myJet.chargedHadronEnergyFraction());
248     jetCon->SetNeuHadFrac(myJet.neutralHadronEnergyFraction());
249     jetCon->SetChEmFrac(myJet.chargedEmEnergyFraction());
250     jetCon->SetNeuEmFrac(myJet.neutralEmEnergyFraction());
251 naodell 1.4
252 naodell 1.6 jetCon->SetNumConstit(myJet.chargedMultiplicity() + myJet.neutralMultiplicity());
253     jetCon->SetNumChPart(myJet.chargedMultiplicity());
254 naodell 1.4
255 naodell 1.6 jetCon->SetNumChPart(myJet.chargedMultiplicity());
256 naodell 1.4
257 naodell 1.7 if(jet_it_2 != bTags2.end()) jetCon->SetBDiscrTrkCountHiPure(jet_it_2->second);
258     if(jet_it_1 != bTags1.end()) jetCon->SetBDiscrTrkCountHiEff(jet_it_1->second);
259    
260     float scale2 = correctorL2->correction(myJet);
261     myJet.scaleEnergy(scale2);
262     float scale3 = correctorL3->correction(myJet);
263     myJet.scaleEnergy(scale3);
264    
265     //more corrections?
266    
267 naodell 1.6 jetCon->SetJetCorr(2, scale2);
268     jetCon->SetJetCorr(3, scale3);
269 naodell 1.7
270 naodell 1.16 /////////////////////////
271     //get associated tracks//
272     /////////////////////////
273 naodell 1.15
274     const reco::TrackRefVector &tracks = myJet.getTrackRefs();
275    
276 naodell 1.16 float sumTrackZ, sumTrackIP;
277     float meanTrackZ, meanTrackIP;
278 naodell 1.17 int nAssociatedTracks = -1;
279 naodell 1.16
280 naodell 1.17 sumTrackZ = sumTrackIP = 0;
281     meanTrackZ = meanTrackIP = -999;
282 naodell 1.15
283 naodell 1.17 if(myJet.pt() > 20 && fabs(myJet.eta()) < 2.4){
284    
285     nAssociatedTracks = 0;
286 naodell 1.15
287     for(TrackRefVector::const_iterator iTrack = tracks.begin(); iTrack != tracks.end(); ++iTrack){
288    
289     const reco::Track& myTrack = **iTrack;
290    
291     sumTrackZ += myTrack.vz();
292 naodell 1.16 sumTrackIP += myTrack.dxy(vertexBeamSpot.position());
293 naodell 1.15 ++nAssociatedTracks;
294 naodell 1.16
295 naodell 1.15 }
296 naodell 1.16
297 naodell 1.15 meanTrackZ = sumTrackZ/nAssociatedTracks;
298 naodell 1.16 h1_nAssociatedTracks->Fill(nAssociatedTracks);
299 naodell 1.17 h2_nAssTracksVsJetPt->Fill(myJet.pt(), nAssociatedTracks);
300 naodell 1.16
301 naodell 1.15 for(TrackRefVector::const_iterator iTrack = tracks.begin(); iTrack != tracks.end(); ++iTrack){
302 naodell 1.16
303 naodell 1.15 const reco::Track& myTrack = **iTrack;
304 naodell 1.16
305 naodell 1.15 h1_trackZ->Fill(meanTrackZ - myTrack.vz());
306 naodell 1.16 h1_trackDxy->Fill(myTrack.dxy(vertexBeamSpot.position()));
307    
308     }
309    
310     if(pfCount == 0){
311    
312     ljVertexZ = meanTrackZ;
313    
314     }else{
315    
316     h1_allTrackDeltaZ_LJ->Fill(ljVertexZ - meanTrackZ);
317    
318 naodell 1.15 }
319 naodell 1.16
320     h1_allTrackDeltaZ_PV->Fill(primaryVertexZ - meanTrackZ);
321    
322 naodell 1.15 }
323 naodell 1.16
324     jetCon->SetVtx(0, 0, meanTrackZ);
325    
326     ++pfCount;
327 naodell 1.4 }
328     }
329 naodell 1.1 }
330    
331 naodell 1.10 if(!isRealData){
332 naodell 1.1
333 naodell 1.18 Handle<GenEventInfoProduct> GenInfoHandle;
334     iEvent.getByLabel("generator", GenInfoHandle);
335    
336 naodell 1.1 Handle<reco::GenJetCollection> GenJets;
337     iEvent.getByLabel(GenJetHandle_, GenJets);
338 naodell 1.18
339     ptHat = qScale = -1;
340    
341     if(GenInfoHandle.isValid()){
342    
343     qScale = GenInfoHandle->qScale();
344     ptHat = (GenInfoHandle->hasBinningValues() ? GenInfoHandle->binningValues()[0] : 0.0);
345    
346     h1_ptHat->Fill(ptHat);
347    
348     }
349 naodell 1.16
350 naodell 1.1 for(GenJetCollection::const_iterator jet_iter = GenJets->begin(); jet_iter!= GenJets->end(); ++jet_iter){
351    
352     reco::GenJet myJet = reco::GenJet(*jet_iter);
353    
354 naodell 1.4 if(myJet.pt() > 5){
355 naodell 1.1
356 naodell 1.16 new ((*jetP4_ak5Gen)[genCount]) TLorentzVector(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
357 naodell 1.1
358 naodell 1.16 ++genCount;
359 naodell 1.1
360     }
361 naodell 1.13 }
362 naodell 1.1 }
363 naodell 1.10
364 naodell 1.16 ///////////////////////////
365     //get trigger information//
366     ///////////////////////////
367 andrey 1.2
368 naodell 1.7 if(triggerHLT_){
369 andrey 1.2
370 naodell 1.7 edm::Handle<TriggerResults> hltR;
371     triggerResultsTag_ = InputTag(hlTriggerResults_,"",hltName_);
372     iEvent.getByLabel(triggerResultsTag_,hltR);
373    
374     const TriggerNames & triggerNames = iEvent.triggerNames(*hltR);
375     hlNames=triggerNames.triggerNames();
376    
377 naodell 1.9 string MPI_TriggerNames[] = {"HLT_PixelTracks_Multiplicity70", "HLT_MinBiasBSC_NoBPTX", "HLT_PixelTracks_Multiplicity40","HLT_L1Tech_HCAL_HF", "HLT_IsoTrackHB_8E29", "HLT_IsoTrackHE_8E29", "HLT_L1Tech_RPC_TTU_RBst1_collisions", "HLT_L1_BscMinBiasOR_BptxPlusORMinus", "HLT_L1Tech_BSC_halo_forPhysicsBackground", "HLT_L1Tech_BSC_HighMultiplicity", "HLT_MinBiasPixel_DoubleIsoTrack5", "HLT_MinBiasPixel_DoubleTrack", "HLT_MinBiasPixel_SingleTrack", "HLT_ZeroBiasPixel_SingleTrack", "HLT_MinBiasBSC", "HLT_StoppedHSCP_8E29", "HLT_Jet15U_HcalNoiseFiltered", "HLT_QuadJet15U", "HLT_DiJetAve30U_8E29", "HLT_DiJetAve15U_8E29", "HLT_FwdJet20U", "HLT_Jet50U", "HLT_Jet30U", "HLT_Jet15U", "HLT_BTagMu_Jet10U", "HLT_DoubleJet15U_ForwardBackward", "HLT_BTagIP_Jet50U", "HLT_DoubleLooseIsoTau15", "HLT_SingleLooseIsoTau20", "HLT_HT100U", "HLT_MET100", "HLT_MET45"};
378 naodell 1.7
379     bool triggerPassed = false;
380     triggerStatus = 0x0;
381 naodell 1.4
382 naodell 1.7 for (uint i=0; i<hlNames.size(); ++i) {
383    
384     triggerPassed = triggerDecision(hltR, i);
385 naodell 1.4
386 naodell 1.7 if(triggerPassed){
387 andrey 1.3
388 naodell 1.7 for (uint j = 0; j != 32; ++j){
389    
390     if (hlNames[i] == MPI_TriggerNames[j])
391     {
392 naodell 1.16 //cout<<"trigger name: "<<hlNames[i]<<"list: "<<dec<<j+1<<endl;
393 naodell 1.7 triggerStatus |= 0x01 << j;
394    
395     }
396     }
397 naodell 1.4 }
398 naodell 1.7 }
399     }
400 naodell 1.16 //cout<< "total status: "<<hex<<triggerStatus<<endl;
401 naodell 1.18
402 naodell 1.14
403    
404 naodell 1.16 if((pfCount > 3 || genCount > 3) && vtxCount > 0) sTree -> Fill();
405 naodell 1.15
406 naodell 1.4 jet_ak5PF->Clear();
407 naodell 1.5 primaryVtx->Clear();
408 naodell 1.4 jetP4_ak5Gen->Clear();
409 andrey 1.2
410 naodell 1.1 }
411    
412    
413     // ------------ method called once each job just before starting event loop ------------
414     void MPIntuple::beginJob()
415     {
416    
417     ntupleFile = new TFile(rootfilename.c_str(), "RECREATE");
418 naodell 1.5 sTree = new TTree("mpiTree", "Tree for Jets");
419 naodell 1.15
420 naodell 1.16 h1_trackZ = new TH1F("h1_trackZ", "#Delta z for associated tracks", 60, -0.3, 0.3);
421     h1_nAssociatedTracks = new TH1F("h1_nAssociatedTracks", "Number of tracks associated to jet", 20, -0.5, 19.5);
422     h1_trackDxy = new TH1F("h1_trackDxy", "dxy for all associated tracks", 50, -0.4, 0.4);
423     h1_allTrackDeltaZ_PV = new TH1F("h1_allTrackDeltaZ_PV", "#Delta z between jet and primary vertex", 60, -0.3, 0.3);
424     h1_allTrackDeltaZ_LJ = new TH1F("h1_allTrackDeltaZ_LJ", "#Delta z between leading jet and other jet vertices", 60, -0.3, 0.3);
425 naodell 1.18 h1_ptHat = new TH1F("h1_ptHat", "ptHat", 50, 10, 100);
426 naodell 1.17 h2_nAssTracksVsJetPt = new TH2F("h2_nAssTracksVsJetPt", "nTracks vs Jet p_{T}", 50, 0, 100, 20, -0.5, 19.5);
427 naodell 1.15 p1_nVtcs = new TProfile("p1_nVtcs", "<N> per run", 5200.0, 132300.0, 137500.0, 0.0, 6.0);
428    
429 naodell 1.4 jet_ak5PF = new TClonesArray("TCJet");
430 naodell 1.1 jetP4_ak5Gen = new TClonesArray("TLorentzVector");
431 naodell 1.5 primaryVtx = new TClonesArray("TCPrimaryVtx");
432 naodell 1.1
433 naodell 1.15
434 naodell 1.4 sTree->Branch("jet_ak5PF",&jet_ak5PF, 6400, 0);
435 naodell 1.1 sTree->Branch("jetP4_ak5Gen",&jetP4_ak5Gen, 6400, 0);
436 naodell 1.5 sTree->Branch("primaryVtx",&primaryVtx, 6400, 0);
437 naodell 1.1
438     sTree->Branch("eventNumber",&eventNumber, "eventNumber/I");
439     sTree->Branch("runNumber",&runNumber, "runNumber/I");
440     sTree->Branch("lumiSection",&lumiSection, "lumiSection/I");
441 naodell 1.6 sTree->Branch("triggerStatus",&triggerStatus, "triggerStatus/i");
442 naodell 1.8 sTree->Branch("isRealData",&isRealData, "isRealData/i");
443     sTree->Branch("bunchCross",&bunchCross, "bunchCross/i");
444 naodell 1.18 sTree->Branch("ptHat",&ptHat, "ptHat/d");
445 naodell 1.1
446     }
447    
448     // ------------ method called once each job just after ending the event loop ------------
449     void MPIntuple::endJob()
450     {
451 naodell 1.15
452 naodell 1.17 ntupleFile->cd();
453    
454 naodell 1.15 h1_trackZ->Write();
455 naodell 1.16 h1_nAssociatedTracks->Write();
456     h1_trackDxy->Write();
457     h1_allTrackDeltaZ_PV->Write();
458     h1_allTrackDeltaZ_LJ->Write();
459 naodell 1.18 h1_ptHat->Write();
460 naodell 1.16 h2_nAssTracksVsJetPt->Write();
461 naodell 1.15 p1_nVtcs->Write();
462    
463 naodell 1.1 ntupleFile->Write();
464     ntupleFile->Close();
465    
466     }
467    
468    
469 andrey 1.2 bool MPIntuple::triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger)
470     {
471     bool triggerPassed = false;
472     if(hltR->wasrun(iTrigger) &&
473     hltR->accept(iTrigger) &&
474     !hltR->error(iTrigger) ){
475     triggerPassed = true;
476     }
477     return triggerPassed;
478 naodell 1.16 }
479    
480     double MPIntuple::sumPtSquared(const Vertex & v)
481     {
482     double sum = 0.;
483     double pT;
484     for (Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); it++) {
485     pT = (**it).pt();
486     double epT=(**it).ptError(); pT=pT>epT ? pT-epT : 0;
487    
488     sum += pT*pT;
489 andrey 1.2 }
490 naodell 1.16 return sum;
491     }
492 andrey 1.2
493 naodell 1.1 //define this as a plug-in
494     DEFINE_FWK_MODULE(MPIntuple);