ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/src/ProducerConversions.cc
Revision: 1.25
Committed: Fri Mar 30 01:08:39 2012 UTC (13 years, 1 month ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_028, Mit_027a, Mit_027, Mit_026, HEAD
Changes since 1.24: +4 -5 lines
Log Message:
Initial 5 version.

File Contents

# Content
1 // $Id: ProducerConversions.cc,v 1.24 2010/06/08 20:17:30 bendavid Exp $
2
3 #include "MitEdm/Producers/interface/ProducerConversions.h"
4 #include "DataFormats/Common/interface/Handle.h"
5 #include "DataFormats/TrackReco/interface/Track.h"
6 #include "DataFormats/TrackReco/interface/TrackFwd.h"
7 #include "DataFormats/VertexReco/interface/Vertex.h"
8 #include "DataFormats/VertexReco/interface/VertexFwd.h"
9 #include "MagneticField/Engine/interface/MagneticField.h"
10 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
11 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
12 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
13 #include "MitEdm/Producers/interface/HitDropperRecord.h"
14 #include "MitEdm/Producers/interface/HitDropper.h"
15 #include "MitEdm/DataFormats/interface/Types.h"
16 #include "MitEdm/DataFormats/interface/Collections.h"
17 #include "MitEdm/DataFormats/interface/DecayPart.h"
18 #include "MitEdm/DataFormats/interface/StablePart.h"
19 #include "MitEdm/DataFormats/interface/StableData.h"
20 #include "MitEdm/VertexFitInterface/interface/MvfInterface.h"
21 #include "MitEdm/VertexFitInterface/interface/TrackParameters.h"
22 #include <TMath.h>
23
24 using namespace std;
25 using namespace edm;
26 using namespace reco;
27 using namespace mitedm;
28 using namespace mithep;
29
30 //--------------------------------------------------------------------------------------------------
31 ProducerConversions::ProducerConversions(const ParameterSet& cfg) :
32 BaseCandProducer (cfg),
33 iStables1_ (cfg.getUntrackedParameter<string>("iStables1","")),
34 iStables2_ (cfg.getUntrackedParameter<string>("iStables2","")),
35 iPVertexes_ (cfg.getUntrackedParameter<string>("iPVertexes","offlinePrimaryVerticesWithBS")),
36 usePVertex_ (cfg.getUntrackedParameter<bool> ("usePVertex",true)),
37 convConstraint_ (cfg.getUntrackedParameter<bool> ("convConstraint",false)),
38 convConstraint3D_(cfg.getUntrackedParameter<bool> ("convConstraint3D",true)),
39 rhoMin_ (cfg.getUntrackedParameter<double>("rhoMin",0.0)),
40 useRhoMin_ (cfg.getUntrackedParameter<bool> ("useRhoMin",true)),
41 useHitDropper_ (cfg.getUntrackedParameter<bool> ("useHitDropper",true)),
42 applyChargeConstraint_(cfg.getUntrackedParameter<bool> ("applyChargeConstraint",false))
43 {
44 // Constructor.
45
46 produces<DecayPartCol>();
47 }
48
49 //--------------------------------------------------------------------------------------------------
50 ProducerConversions::~ProducerConversions()
51 {
52 // Destructor.
53 }
54
55 //--------------------------------------------------------------------------------------------------
56 void ProducerConversions::produce(Event &evt, const EventSetup &setup)
57 {
58 // Produce our DecayPartCol.
59
60 // First input collection
61 Handle<StablePartCol> hStables1;
62 if (!GetProduct(iStables1_, hStables1, evt)) {
63 printf("Stable collection 1 not found in ProducerConversions\n");
64 return;
65 }
66 const StablePartCol *pS1 = hStables1.product();
67 // Second input collection
68 Handle<StablePartCol> hStables2;
69 if (!GetProduct(iStables2_, hStables2, evt)) {
70 printf("Stable collection 2 not found in ProducerConversions\n");
71 return;
72 }
73 const StablePartCol *pS2 = hStables2.product();
74
75 const reco::Vertex *vertex = 0;
76 mitedm::VertexPtr vPtr;
77 if (usePVertex_) {
78 // Primary vertex collection
79 Handle<reco::VertexCollection> hVertexes;
80 if (!GetProduct(iPVertexes_, hVertexes, evt))
81 return;
82 const reco::VertexCollection *pV = hVertexes.product();
83
84 //Choose the primary vertex with the largest number of tracks
85 UInt_t maxTracks=0;
86 for (UInt_t i=0; i<pV->size(); ++i) {
87 const reco::Vertex &v = pV->at(i);
88 UInt_t nTracks = v.tracksSize();
89 if (nTracks >= maxTracks) {
90 maxTracks = nTracks;
91 vertex = &v;
92 vPtr = mitedm::VertexPtr(hVertexes,i);
93 }
94 }
95 }
96
97 // Get hit dropper
98 ESHandle<HitDropper> hDropper;
99 setup.get<HitDropperRecord>().get("HitDropper",hDropper);
100 const HitDropper *dropper = hDropper.product();
101
102 //Get Magnetic Field from event setup, taking value at (0,0,0)
103 edm::ESHandle<MagneticField> magneticField;
104 setup.get<IdealMagneticFieldRecord>().get(magneticField);
105 const double bfield = magneticField->inTesla(GlobalPoint(0.,0.,0.)).z();
106
107 edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
108 setup.get<TransientTrackRecord>().get("TransientTrackBuilder",hTransientTrackBuilder);
109 const TransientTrackBuilder *transientTrackBuilder = hTransientTrackBuilder.product();
110
111 //construct intermediate collection of TrackParameters in mvf format for vertex fit
112 std::vector<TrackParameters> trkPars1;
113 for (UInt_t i = 0; i<pS1->size(); ++i) {
114 const reco::Track *t = pS1->at(i).track();
115 const reco::TransientTrack ttrk = transientTrackBuilder->build(t);
116 TrackParameters cmsTrk(ttrk);
117 TrackParameters mvfTrk = cmsTrk.mvfTrack();
118 trkPars1.push_back(mvfTrk);
119 }
120
121 bool sameCollection = iStables1_ == iStables2_;
122
123 std::vector<TrackParameters> trkPars2;
124 if (sameCollection)
125 trkPars2 = trkPars1;
126 else for (UInt_t i = 0; i<pS2->size(); ++i) {
127 const reco::Track *t = pS2->at(i).track();
128 const reco::TransientTrack ttrk = transientTrackBuilder->build(t);
129 TrackParameters cmsTrk(ttrk);
130 TrackParameters mvfTrk = cmsTrk.mvfTrack();
131 trkPars2.push_back(mvfTrk);
132 }
133
134
135 // Create the output collection
136 auto_ptr<DecayPartCol> pD(new DecayPartCol());
137
138 ClosestApproachInRPhi helixIntersector;
139
140 int nFits = 0;
141
142 //printf("S1 size = %i\n", pS1->size());
143
144 std::set<std::pair<edm::Ptr<reco::Track>,edm::Ptr<reco::Track> > > pairSet;
145 // Simple double loop
146 for (UInt_t i = 0; i<pS1->size(); ++i) {
147 const StablePart &s1 = pS1->at(i);
148
149 const reco::Track * t1 = s1.track();
150 const TrackParameters &trkPar1 = trkPars1.at(i);
151
152 UInt_t j;
153 if (sameCollection)
154 j = i+1;
155 else
156 j = 0;
157
158
159 FreeTrajectoryState initialState1 = trajectoryStateTransform::initialFreeState(*s1.track(),&*magneticField);
160
161 for (; j<pS2->size(); ++j) {
162 const StablePart &s2 = pS2->at(j);
163
164 if (!sameCollection) {
165 std::pair<edm::Ptr<reco::Track>,edm::Ptr<reco::Track> > pair1(s1.trackPtr(),s2.trackPtr());
166 std::pair<edm::Ptr<reco::Track>,edm::Ptr<reco::Track> > pair2(s2.trackPtr(),s1.trackPtr());
167 if (pairSet.find(pair1)!=pairSet.end() || pairSet.find(pair2)!=pairSet.end()) {
168 continue;
169 }
170
171 pairSet.insert(pair1);
172 pairSet.insert(pair2);
173 }
174
175 //Do fast helix fit to check if there's any hope
176 const reco::Track * t2 = s2.track();
177 const TrackParameters &trkPar2 = trkPars2.at(j);
178
179 int trackCharge = t1->charge() + t2->charge();
180
181 double dR0 = 0.0;
182
183 if (!applyChargeConstraint_ || trackCharge==0) {
184 FreeTrajectoryState initialState2 = trajectoryStateTransform::initialFreeState(*s2.track(),&*magneticField);
185 helixIntersector.calculate(initialState1, initialState2);
186 if (helixIntersector.status())
187 dR0 = helixIntersector.crossingPoint().perp();
188 }
189
190 int fitStatus = 0;
191 MultiVertexFitterD fit;
192 if ( (!applyChargeConstraint_ || trackCharge==0) && (!useRhoMin_ || dR0 > rhoMin_) ) {
193
194 // Vertex fit now, possibly with conversion constraint
195 nFits++;
196
197 fit.init(bfield); // Reset to the magnetic field from the event setup
198
199 fit.addTrack(*trkPar1.pars(),*trkPar1.cMat(),1,s1.mass(),MultiVertexFitterD::VERTEX_1);
200 fit.addTrack(*trkPar2.pars(),*trkPar2.cMat(),2,s2.mass(),MultiVertexFitterD::VERTEX_1);
201 if (convConstraint3D_) {
202 fit.conversion_3d(MultiVertexFitterD::VERTEX_1);
203 //printf("applying 3d conversion constraint\n");
204 }
205 else if (convConstraint_) {
206 fit.conversion_2d(MultiVertexFitterD::VERTEX_1);
207 //printf("applying 2d conversion constraint\n");
208 }
209 //initialize primary vertex parameters in the fitter
210 if (usePVertex_) {
211 float vErr[3][3];
212 for (UInt_t vi=0; vi<3; ++vi)
213 for (UInt_t vj=0; vj<3; ++vj)
214 vErr[vi][vj] = vertex->covariance(vi,vj);
215
216 fit.setPrimaryVertex(vertex->x(),vertex->y(),vertex->z());
217 fit.setPrimaryVertexError(vErr);
218 }
219
220 fitStatus = fit.fit();
221
222 }
223
224 if (fitStatus) {
225 DecayPart *d = new DecayPart(oPid_,DecayPart::Fast);
226
227 // Update temporarily some of the quantities (prob, chi2, nDoF, mass, lxy, pt, fourMomentum)
228 d->setProb(fit.prob());
229 d->setChi2(fit.chisq());
230 d->setNdof(fit.ndof());
231
232 FourVector p4Fitted(0.,0.,0.,0.);
233 p4Fitted += fit.getTrackP4(1);
234 p4Fitted += fit.getTrackP4(2);
235 d->setFourMomentum(p4Fitted);
236 d->setPosition(fit.getVertex (MultiVertexFitterD::VERTEX_1));
237 d->setError (fit.getErrorMatrix(MultiVertexFitterD::VERTEX_1));
238 float mass, massErr;
239 const int trksIds[2] = { 1, 2 };
240 mass = fit.getMass(2,trksIds,massErr);
241
242 ThreeVector p3Fitted(p4Fitted.px(), p4Fitted.py(), p4Fitted.pz());
243
244 // Get decay length in xy plane
245 float dl, dlErr;
246 dl = fit.getDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
247 p3Fitted, dlErr);
248
249 // Get Z decay length
250 float dlz, dlzErr;
251 dlz = fit.getZDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
252 p3Fitted, dlzErr);
253
254 // Get impact parameter
255 float dxy, dxyErr;
256 dxy = fit.getImpactPar(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
257 p3Fitted, dxyErr);
258
259 BasePartPtr ptr1(hStables1,i);
260 BasePartPtr ptr2(hStables2,j);
261
262 StableData c1(fit.getTrackP4(1).px(),fit.getTrackP4(1).py(), fit.getTrackP4(1).pz(), ptr1);
263 StableData c2(fit.getTrackP4(2).px(),fit.getTrackP4(2).py(), fit.getTrackP4(2).pz(), ptr2);
264
265 const ThreeVector vtxPos = fit.getVertex(MultiVertexFitterD::VERTEX_1);
266 const ThreeVector trkMom1(fit.getTrackP4(1).px(),
267 fit.getTrackP4(1).py(),
268 fit.getTrackP4(1).pz());
269 const ThreeVector trkMom2(fit.getTrackP4(2).px(),
270 fit.getTrackP4(2).py(),
271 fit.getTrackP4(2).pz());
272
273 // Build corrected HitPattern for StableData, removing hits before the fit vertex
274 if (useHitDropper_) {
275 std::pair<reco::HitPattern,uint> hits1 = dropper->CorrectedHitsAOD(s1.track(), vtxPos, trkMom1, dlErr, dlzErr);
276 std::pair<reco::HitPattern,uint> hits2 = dropper->CorrectedHitsAOD(s2.track(), vtxPos, trkMom2, dlErr, dlzErr);
277
278 c1.SetHits(hits1.first);
279 c2.SetHits(hits2.first);
280 c1.SetHitsFilled();
281 c2.SetHitsFilled();
282 c1.SetNWrongHits(hits1.second);
283 c2.SetNWrongHits(hits2.second);
284
285 reco::HitPattern sharedHits = dropper->SharedHits(s1.track(),s2.track());
286 d->setSharedHits(sharedHits);
287 }
288
289 d->addStableChild(c1);
290 d->addStableChild(c2);
291
292
293 d->setFittedMass (mass);
294 d->setFittedMassError(massErr);
295
296 d->setLxy(dl);
297 d->setLxyError(dlErr);
298 d->setLxyToPv(dl);
299 d->setLxyToPvError(dlErr);
300
301 d->setLz(dlz);
302 d->setLzError(dlzErr);
303 d->setLzToPv(dlz);
304 d->setLzToPvError(dlzErr);
305
306 d->setDxy(dxy);
307 d->setDxyError(dxyErr);
308 d->setDxyToPv(dxy);
309 d->setDxyToPvError(dxyErr);
310
311 if (usePVertex_)
312 d->setPrimaryVertex(vPtr);
313
314 // Put the result into our collection
315 pD->push_back(*d);
316
317 delete d;
318 }
319 }
320 }
321
322 //printf("nConversionFits = %i\n",nFits);
323
324 // Write the collection even if it is empty
325 if (0) {
326 cout << " ProducerConversions::produce - " << pD->size() << " entries collection created -"
327 << " (Pid: " << oPid_ << ")\n";
328 }
329 evt.put(pD);
330 }
331
332 //define this as a plug-in
333 DEFINE_FWK_MODULE(ProducerConversions);