ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/plugins/InFlightDecayRecoFinder.cc
Revision: 1.3
Committed: Thu Feb 18 22:01:48 2010 UTC (15 years, 3 months ago) by gpetrucc
Content type: text/plain
Branch: MAIN
CVS Tags: V03-00-00, HEAD
Changes since 1.2: +44 -32 lines
Log Message:
Added debug printouts, fix one bug

File Contents

# User Rev Content
1 gpetrucc 1.1 //
2 gpetrucc 1.3 // $Id: InFlightDecayRecoFinder.cc,v 1.2 2010/01/31 10:20:59 gpetrucc Exp $
3 gpetrucc 1.1 //
4    
5     /**
6     \class pat::InFlightDecayRecoFinder InFlightDecayRecoFinder.h "PhysicsTools/PatAlgos/interface/InFlightDecayRecoFinder.h"
7     \brief Use StandAlone track to define the 4-momentum of a PAT Muon (normally the global one is used)
8    
9     \author Giovanni Petrucciani
10 gpetrucc 1.3 \version $Id: InFlightDecayRecoFinder.cc,v 1.2 2010/01/31 10:20:59 gpetrucc Exp $
11 gpetrucc 1.1 */
12    
13     #include <algorithm>
14    
15     #include "FWCore/Framework/interface/EDProducer.h"
16     #include "FWCore/Framework/interface/Event.h"
17     #include "FWCore/Framework/interface/ESHandle.h"
18     #include "FWCore/ParameterSet/interface/ParameterSet.h"
19     #include "FWCore/ParameterSet/interface/InputTag.h"
20    
21     #include "DataFormats/Math/interface/deltaR.h"
22     #include "DataFormats/Common/interface/View.h"
23 gpetrucc 1.2 #include "DataFormats/MuonReco/interface/MuonFwd.h"
24 gpetrucc 1.1 #include "DataFormats/MuonReco/interface/Muon.h"
25     #include "DataFormats/TrackReco/interface/Track.h"
26 gpetrucc 1.2 #include "DataFormats/VertexReco/interface/Vertex.h"
27    
28     #include "DataFormats/PatCandidates/interface/GenericParticle.h"
29 gpetrucc 1.1
30     #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
31    
32     #include "TrackingTools/GeomPropagators/interface/Propagator.h"
33     #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
34     #include "MagneticField/Engine/interface/MagneticField.h"
35     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
36     #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
37     #include "Geometry/CommonDetUnit/interface/GeomDet.h"
38     #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
39     #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
40     #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
41     #include <TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h>
42 gpetrucc 1.2 #include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
43     #include "TrackingTools/IPTools/interface/IPTools.h"
44     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
45     #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
46     #include "TrackingTools/Records/interface/TransientTrackRecord.h"
47 gpetrucc 1.1
48     #include "Math/GenVector/Boost.h"
49    
50     #include <boost/foreach.hpp>
51     #define foreach BOOST_FOREACH
52    
53    
54     class InFlightDecayRecoFinder : public edm::EDProducer {
55     public:
56     explicit InFlightDecayRecoFinder(const edm::ParameterSet & iConfig);
57     virtual ~InFlightDecayRecoFinder() { }
58    
59     virtual void produce(edm::Event & iEvent, const edm::EventSetup & iSetup);
60    
61     private:
62     /// input muons
63     edm::InputTag muons_;
64     /// input tracker tracks
65     edm::InputTag tracks_;
66 gpetrucc 1.2 /// input primary vertex
67     edm::InputTag vertices_;
68 gpetrucc 1.1 /// preselection for muons
69     StringCutObjectSelector<reco::Muon> muonCut_;
70     /// preselection for candidate "inner" tracks (hadrons that decay)
71     StringCutObjectSelector<reco::Track> innerTrackCut_;
72     /// preselection for candidate "outer" tracks (muons from decay)
73     StringCutObjectSelector<reco::Track> outerTrackCut_;
74     /// minimal deltaR to search for tracks around the muon
75     double maxDeltaR2_;
76     /// minimal deltaPt/pt to search for tracks around the muon
77     double maxDeltaPtRel_;
78     /// maximum distance between the two tracks
79     double maxTwoTrackDist_;
80    
81 gpetrucc 1.2 /// Use only variables available in AOD
82     bool aodDataOnly_;
83    
84 gpetrucc 1.1 /// Magnetic field for extrapolation
85     edm::ESHandle<MagneticField> magfield_;
86    
87     /// Propagator for extrapolation
88     edm::ESHandle<Propagator> propagator_;
89    
90     /// Global tracking geometry
91     edm::ESHandle<GlobalTrackingGeometry> geometry_;
92    
93 gpetrucc 1.2 /// Transient track builder, for better d0 determination
94     edm::ESHandle<TransientTrackBuilder> ttBuilder_;
95    
96 gpetrucc 1.1 struct EqByTrack {
97     EqByTrack() : tk() {}
98     EqByTrack(const reco::TrackRef & ref) : tk(ref) {}
99    
100     bool operator()(const reco::Muon &mu) const { return mu.track() == tk; }
101     bool operator()(const reco::Muon &mu, const reco::TrackRef &tk) const { return mu.track() == tk; }
102     bool operator()(const reco::TrackRef &tk, const reco::Muon &mu) const { return mu.track() == tk; }
103     bool operator()(const reco::Muon &m1, const reco::Muon &m2) const { return m1.track() == m2.track(); }
104     reco::TrackRef tk;
105     };
106    
107     void collectTracks(const reco::Muon &mu, const edm::Handle<reco::TrackCollection> &tracks, std::vector<reco::TrackRef> & selected) {
108     selected.clear();
109     for (size_t i = 0, n = tracks->size(); i < n; ++i) {
110     const reco::Track &tk = (*tracks)[i];
111     if ((reco::deltaR2(mu, tk) < maxDeltaR2_ ) &&
112     (fabs(mu.pt() - tk.pt())/mu.pt() < maxDeltaPtRel_) &&
113     mu.charge() == tk.charge() ) {
114     selected.push_back(reco::TrackRef(tracks,i));
115     }
116     }
117     }
118    
119     struct DecayCandidate {
120     DecayCandidate(reco::TrackRef inner, reco::TrackRef outer) :
121     innerTrack(inner), outerTrack(outer), distance(-1) {}
122     reco::TrackRef innerTrack; // never null
123     reco::TrackRef outerTrack; // never null
124 gpetrucc 1.2 edm::Ptr<reco::Muon> innerMu; // may be null
125     edm::Ptr<reco::Muon> outerMu; // may be null
126 gpetrucc 1.1 double distance;
127     GlobalPoint decayPoint;
128     reco::Particle::LorentzVector p4mu, p4nu, p4had;
129     double theta, thetaStar;
130    
131     bool operator==(const DecayCandidate &other) {
132     return (innerTrack == other.innerTrack) && (outerTrack == other.outerTrack);
133     }
134     };
135    
136 gpetrucc 1.3 /// events to print out
137     int printOut_;
138 gpetrucc 1.1 };
139    
140     InFlightDecayRecoFinder::InFlightDecayRecoFinder(const edm::ParameterSet & iConfig) :
141     muons_(iConfig.getParameter<edm::InputTag>("muons")),
142     tracks_(iConfig.getParameter<edm::InputTag>("tracks")),
143 gpetrucc 1.2 vertices_(iConfig.getParameter<edm::InputTag>("vertices")),
144 gpetrucc 1.1 muonCut_(iConfig.getParameter<std::string>("muonCut")),
145     innerTrackCut_(iConfig.getParameter<std::string>("innerTrackCut")),
146     outerTrackCut_(iConfig.getParameter<std::string>("outerTrackCut")),
147     maxDeltaR2_(std::pow(iConfig.getParameter<double>("maxDeltaR"),2)),
148     maxDeltaPtRel_(iConfig.getParameter<double>("maxDeltaPtRel")),
149 gpetrucc 1.2 maxTwoTrackDist_(iConfig.getParameter<double>("maxTwoTrackDist")),
150 gpetrucc 1.3 aodDataOnly_(iConfig.getParameter<bool>("aodDataOnly")),
151     printOut_(iConfig.getUntrackedParameter<int32_t>("eventsToPrintOut",0))
152 gpetrucc 1.1 {
153     produces<reco::TrackCollection>("inner");
154     produces<reco::TrackCollection>("outer");
155 gpetrucc 1.2 produces<std::vector<reco::Muon> >("inner");
156     produces<std::vector<reco::Muon> >("outer");
157     produces<pat::GenericParticleCollection>();
158 gpetrucc 1.1 }
159    
160     void
161     InFlightDecayRecoFinder::produce(edm::Event & iEvent, const edm::EventSetup & iSetup) {
162     using namespace edm;
163     using namespace std;
164    
165     // Get EventSetup data
166     iSetup.get<IdealMagneticFieldRecord>().get(magfield_);
167     //iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator_);
168 gpetrucc 1.2 if (!aodDataOnly_) iSetup.get<GlobalTrackingGeometryRecord>().get(geometry_);
169     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", ttBuilder_);
170    
171 gpetrucc 1.1
172     // Get Event data
173     Handle<View<reco::Muon> > muons;
174 gpetrucc 1.2 iEvent.getByLabel(muons_, muons);
175 gpetrucc 1.1 Handle<reco::TrackCollection> tracks;
176     iEvent.getByLabel(tracks_, tracks);
177 gpetrucc 1.2 Handle<std::vector<reco::Vertex> > vertices;
178     iEvent.getByLabel(vertices_, vertices);
179 gpetrucc 1.1
180 gpetrucc 1.2 auto_ptr<reco::TrackCollection> innerSelTks(new reco::TrackCollection());
181     auto_ptr<reco::TrackCollection> outerSelTks(new reco::TrackCollection());
182     auto_ptr<reco::MuonCollection> innerSelMus(new reco::MuonCollection());
183     auto_ptr<reco::MuonCollection> outerSelMus(new reco::MuonCollection());
184     auto_ptr<pat::GenericParticleCollection> output(new pat::GenericParticleCollection());
185 gpetrucc 1.1
186 gpetrucc 1.3
187     std::vector<DecayCandidate> decayCandidates;
188 gpetrucc 1.1 typedef std::vector<reco::TrackRef> TkRefs;
189     EqByTrack eq;
190     TkRefs nearbyTracks;
191     edm::View<reco::Muon>::const_iterator muBegin = muons->begin(), muEnd = muons->end();
192     for (View<reco::Muon>::const_iterator itmu = muBegin; itmu != muEnd; ++itmu) {
193     if (!muonCut_(*itmu)) continue;
194     collectTracks(*itmu, tracks, nearbyTracks);
195 gpetrucc 1.3 if (printOut_) {
196     std::cerr << "Seeding from muon with pt = " << itmu->pt() << ", eta = " << itmu->eta() << ", phi = " << itmu->phi() <<
197     " track id " << itmu->track().key() <<
198     ( itmu->isGlobalMuon() ? " Glb" : "") << ( itmu->isTrackerMuon() ? " Trk" : "") << ( itmu->isStandAloneMuon() ? " Sta" : "") <<
199     ". found " << nearbyTracks.size() << " nearby tracks.\n";
200     }
201 gpetrucc 1.1 for (TkRefs::const_iterator it1 = nearbyTracks.begin(), ed = nearbyTracks.end(); it1 != ed; ++it1) {
202     if (!innerTrackCut_(**it1)) continue;
203 gpetrucc 1.3 if (printOut_) std::cerr << " inner track id " << it1->key() << " pt = " << (*it1)->pt() << ", eta = " << (*it1)->eta() << ", phi = " << (*it1)->phi() <<
204     ", valid hits = " << (*it1)->numberOfValidHits() <<
205     ", expected in/out = " << (*it1)->trackerExpectedHitsInner().numberOfLostHits() << " / " << (*it1)->trackerExpectedHitsOuter().numberOfLostHits() << std::endl;
206 gpetrucc 1.1 for (TkRefs::const_iterator it2 = nearbyTracks.begin(); it2 != ed; ++it2) {
207     if (it1 == it2) continue;
208     if (!outerTrackCut_(**it2)) continue;
209 gpetrucc 1.3 if (printOut_) std::cerr << " outer track id " << it2->key() << " pt = " << (*it2)->pt() << ", eta = " << (*it2)->eta() << ", phi = " << (*it2)->phi() <<
210     ", valid hits = " << (*it2)->numberOfValidHits() <<
211     ", expected in/out = " << (*it2)->trackerExpectedHitsInner().numberOfLostHits() << " / " << (*it2)->trackerExpectedHitsOuter().numberOfLostHits() << std::endl;
212 gpetrucc 1.1 if ((*it1)->charge() != (*it2)->charge()) continue;
213    
214     DecayCandidate cand(*it1,*it2);
215 gpetrucc 1.3 if (std::find(decayCandidates.begin(), decayCandidates.end(), cand) != decayCandidates.end()) {
216     if (printOut_) std::cerr << " already seeded." << std::endl;
217     continue; // can be seeded by both muons
218     }
219 gpetrucc 1.1
220 gpetrucc 1.2 edm::View<reco::Muon>::const_iterator innerMuIt = std::find_if(muBegin, muEnd, EqByTrack(*it1));
221 gpetrucc 1.1 edm::View<reco::Muon>::const_iterator outerMuIt = std::find_if(muBegin, muEnd, EqByTrack(*it2));
222 gpetrucc 1.2 if (innerMuIt != muEnd) { cand.innerMu = muons->ptrAt(innerMuIt - muBegin); }
223     if (outerMuIt != muEnd) { cand.outerMu = muons->ptrAt(outerMuIt - muBegin); }
224 gpetrucc 1.1
225     TwoTrackMinimumDistance ttmd;
226     TrajectoryStateTransform transformer;
227    
228 gpetrucc 1.2 if (!aodDataOnly_) {
229     TrajectoryStateOnSurface innerInnerTSOS = transformer.innerStateOnSurface(**it1, *geometry_, &*magfield_); // start from the two far states
230     TrajectoryStateOnSurface outerOuterTSOS = transformer.outerStateOnSurface(**it2, *geometry_, &*magfield_); // so that distance is less biased by shared hits
231     if (!ttmd.calculate(innerInnerTSOS, outerOuterTSOS)) { std::cerr << " TTMD failed :-( " << std::endl; continue; }
232     } else {
233     FreeTrajectoryState innerFTS = transformer.initialFreeState(**it1, &*magfield_); // start from the two far states
234     FreeTrajectoryState outerFTS = transformer.initialFreeState(**it2, &*magfield_); // so that distance is less biased by shared hits
235     if (!ttmd.calculate(innerFTS, outerFTS)) { std::cerr << " TTMD failed :-( " << std::endl; continue; }
236     }
237 gpetrucc 1.1
238     cand.decayPoint = ttmd.crossingPoint();
239     cand.distance = ttmd.distance();
240 gpetrucc 1.2
241    
242 gpetrucc 1.3 if (printOut_) std::cerr << " TTMD dist = " << cand.distance << ", R = " << cand.decayPoint.perp() << ", Z = " << cand.decayPoint.z() << std::endl;
243 gpetrucc 1.1 if (cand.distance > maxTwoTrackDist_) continue;
244    
245     typedef reco::Candidate::Vector Vector;
246 gpetrucc 1.2 typedef reco::Candidate::Point Point;
247 gpetrucc 1.1 typedef reco::Candidate::LorentzVector LorentzVector;
248 gpetrucc 1.2 TSCPBuilderNoMaterial tscpBuilder;
249     Vector ph, pm;
250     if (!aodDataOnly_) {
251     TrajectoryStateOnSurface innerOuterTSOS = transformer.outerStateOnSurface(**it1, *geometry_, &*magfield_); // start from the two near states
252     TrajectoryStateOnSurface outerInnerTSOS = transformer.innerStateOnSurface(**it2, *geometry_, &*magfield_); // so the material is more correct
253     TrajectoryStateClosestToPoint innerTSCP = tscpBuilder(innerOuterTSOS, ttmd.crossingPoint());
254     TrajectoryStateClosestToPoint outerTSCP = tscpBuilder(outerInnerTSOS, ttmd.crossingPoint());
255     ph = Vector(innerTSCP.momentum());
256     pm = Vector(outerTSCP.momentum());
257     } else {
258     FreeTrajectoryState innerFTS = transformer.initialFreeState(**it1, &*magfield_); // start from the two far states
259     FreeTrajectoryState outerFTS = transformer.initialFreeState(**it2, &*magfield_); // so that distance is less biased by shared hits
260     TrajectoryStateClosestToPoint innerTSCP = tscpBuilder(innerFTS, ttmd.crossingPoint());
261     TrajectoryStateClosestToPoint outerTSCP = tscpBuilder(outerFTS, ttmd.crossingPoint());
262     ph = Vector(innerTSCP.momentum());
263     pm = Vector(outerTSCP.momentum());
264     }
265 gpetrucc 1.1 Vector pn = ph-pm;
266     cand.p4mu = LorentzVector(pm.x(), pm.y(), pm.z(), hypot(pm.r(), 0.10566));
267     cand.p4nu = LorentzVector(pn.x(), pn.y(), pn.z(), pn.r());
268     cand.p4had = cand.p4mu + cand.p4nu;
269     cand.theta = asin(pm.Unit().Cross(ph.Unit()).r());
270     Vector beta = cand.p4had.BoostToCM();
271     ROOT::Math::Boost boostH(beta.x(), beta.y(), beta.z());
272     LorentzVector p1star = boostH(cand.p4mu);
273     cand.thetaStar = asin(cand.p4had.Vect().Unit().Cross(p1star.Vect().Unit()).r());
274    
275 gpetrucc 1.3 if (printOut_) {
276     if (eq(*it1,*itmu)) std::cerr << " outer track is the seed reco::Muon" << std::endl;
277     else if (count_if(muons->begin(), muons->end(), EqByTrack(*it1))) std::cerr << " inner track is another reco::Muon" << std::endl;
278     if (eq(*it2,*itmu)) std::cerr << " outer track is the seed reco::Muon" << std::endl;
279     else if (count_if(muons->begin(), muons->end(), EqByTrack(*it2))) std::cerr << " outer track is another reco::Muon" << std::endl;
280    
281    
282     std::cerr << " p(mu) = " << cand.p4mu.P() << std::endl;
283     std::cerr << " p(nu) = " << cand.p4nu.P() << std::endl;
284     std::cerr << " p(had) = " << cand.p4had.P() << std::endl;
285     std::cerr << " p4(mu) = " << cand.p4mu << std::endl;
286     std::cerr << " p4(nu) = " << cand.p4nu << std::endl;
287     std::cerr << " p4(had) = " << cand.p4had << std::endl;
288     std::cerr << " dist = " << cand.distance << std::endl;
289     std::cerr << " mass = " << (cand.p4had.mass2() > 0 ? cand.p4had.mass() : -sqrt(-cand.p4had.mass2())) << std::endl;
290     std::cerr << " theta = " << cand.theta << std::endl;
291     std::cerr << " theta* = " << cand.thetaStar << std::endl;
292     }
293 gpetrucc 1.1
294     decayCandidates.push_back(cand);
295 gpetrucc 1.2 innerSelTks->push_back(*cand.innerTrack);
296     outerSelTks->push_back(*cand.outerTrack);
297     if (cand.innerMu.isNonnull()) innerSelMus->push_back(*cand.innerMu);
298     if (cand.outerMu.isNonnull()) outerSelMus->push_back(*cand.outerMu);
299    
300     output->push_back(pat::GenericParticle());
301     pat::GenericParticle &gp = output->back();
302     /// ----- Base kinematics -----
303     gp.setCharge(cand.innerTrack->charge());
304     gp.setP4(math::PtEtaPhiMLorentzVector(cand.innerTrack->pt(),
305     cand.innerTrack->eta(),
306     cand.innerTrack->phi(),
307     (cand.p4had.mass2() > 0 ? cand.p4had.mass() : -sqrt(-cand.p4had.mass2()))));
308     gp.setVertex(cand.innerTrack->vertex());
309     /// ----- References ------
310     gp.setTrack(cand.innerTrack);
311     reco::TrackRefVector tks;
312     tks.push_back(cand.innerTrack);
313     tks.push_back(cand.outerTrack);
314     gp.setTracks(tks);
315     if (cand.outerMu.isNonnull() && cand.outerMu->isStandAloneMuon()) {
316     gp.setStandAloneMuon(cand.outerMu->standAloneMuon());
317     }
318     if (cand.outerMu.isNonnull() && cand.outerMu->isGlobalMuon()) {
319     gp.setCombinedMuon(cand.outerMu->combinedMuon());
320     }
321     gp.addUserCand("innerMu", cand.innerMu);
322     gp.addUserCand("outerMu", cand.outerMu);
323     /// ---- Decay kinematics ----
324     gp.addUserFloat("distance", cand.distance);
325     gp.addUserFloat("decayRho", cand.decayPoint.perp());
326     gp.addUserFloat("decayZ", cand.decayPoint.z());
327     gp.addUserData( "decayPoint", Point(cand.decayPoint));
328     gp.addUserFloat("theta", cand.theta);
329     gp.addUserFloat("thetaStar", cand.thetaStar);
330     /// ---- Impact parameter stuff ----
331     // first the basics
332     reco::VertexRef vtx(vertices, 0); // primary vertex
333     gp.addUserFloat("innerDxyVtx",cand.innerTrack->dxy(vtx->position()));
334     gp.addUserFloat("outerDxyVtx",cand.outerTrack->dxy(vtx->position()));
335     gp.addUserFloat("innerDszVtx",cand.innerTrack->dsz(vtx->position()));
336     gp.addUserFloat("outerDszVtx",cand.outerTrack->dsz(vtx->position()));
337     // the we get fancy
338     GlobalVector p4atVtx(cand.innerTrack->px(), cand.innerTrack->py(), cand.innerTrack->pz());
339     reco::TransientTrack innerTT = ttBuilder_->build(*cand.innerTrack);
340     reco::TransientTrack outerTT = ttBuilder_->build(*cand.outerTrack);
341     std::pair<bool,Measurement1D> innerTIP = IPTools::signedTransverseImpactParameter(innerTT, p4atVtx, *vtx);
342     std::pair<bool,Measurement1D> inner3DIP = IPTools::signedImpactParameter3D(innerTT, p4atVtx, *vtx);
343     std::pair<bool,Measurement1D> innerDL3D = IPTools::signedDecayLength3D(innerTT, p4atVtx, *vtx);
344     std::pair<bool,Measurement1D> outerTIP = IPTools::signedTransverseImpactParameter(outerTT, p4atVtx, *vtx);
345     std::pair<bool,Measurement1D> outer3DIP = IPTools::signedImpactParameter3D(outerTT, p4atVtx, *vtx);
346     std::pair<bool,Measurement1D> outerDL3D = IPTools::signedDecayLength3D(outerTT, p4atVtx, *vtx);
347     if (innerTIP.first ) { gp.addUserFloat("innerTIP", innerTIP.second.value()); gp.addUserFloat("innerTIPsig", innerTIP.second.significance()); }
348     if (inner3DIP.first) { gp.addUserFloat("inner3DIP",inner3DIP.second.value()); gp.addUserFloat("inner3DIPsig",inner3DIP.second.significance()); }
349     if (innerDL3D.first) { gp.addUserFloat("innerDL3D",innerDL3D.second.value()); gp.addUserFloat("innerDL3Dsig",innerDL3D.second.significance()); }
350     if (outerTIP.first ) { gp.addUserFloat("outerTIP", outerTIP.second.value()); gp.addUserFloat("outerTIPsig", outerTIP.second.significance()); }
351     if (outer3DIP.first) { gp.addUserFloat("outer3DIP",outer3DIP.second.value()); gp.addUserFloat("outer3DIPsig",outer3DIP.second.significance()); }
352     if (outerDL3D.first) { gp.addUserFloat("outerDL3D",outerDL3D.second.value()); gp.addUserFloat("outerDL3Dsig",outerDL3D.second.significance()); }
353 gpetrucc 1.1 }
354     }
355     }
356 gpetrucc 1.2
357     iEvent.put(innerSelTks, "inner");
358     iEvent.put(outerSelTks, "outer");
359     iEvent.put(innerSelMus, "inner");
360     iEvent.put(outerSelMus, "outer");
361     iEvent.put(output);
362 gpetrucc 1.3 if (printOut_ > 0) { printOut_--; std::cerr << "\n" << std::endl; }
363 gpetrucc 1.1 }
364    
365     #include "FWCore/Framework/interface/MakerMacros.h"
366     DEFINE_FWK_MODULE(InFlightDecayRecoFinder);