ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/src/ProducerConversions.cc
Revision: 1.23
Committed: Mon Jan 18 14:41:33 2010 UTC (15 years, 3 months ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, Mit_013pre1, Mit_012i
Changes since 1.22: +37 -6 lines
Log Message:
Fix Track Paramater phi-related bugs and add consistent propagation to origin

File Contents

# Content
1 // $Id: ProducerConversions.cc,v 1.22 2009/12/15 23:27:34 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 std::vector<TrackParameters> trkPars2;
122 if (iStables1_ == iStables2_)
123 trkPars2 = trkPars1;
124 else for (UInt_t i = 0; i<pS2->size(); ++i) {
125 const reco::Track *t = pS2->at(i).track();
126 const reco::TransientTrack ttrk = transientTrackBuilder->build(t);
127 TrackParameters cmsTrk(ttrk);
128 TrackParameters mvfTrk = cmsTrk.mvfTrack();
129 trkPars2.push_back(mvfTrk);
130 }
131
132
133 // Create the output collection
134 auto_ptr<DecayPartCol> pD(new DecayPartCol());
135
136 ClosestApproachInRPhi helixIntersector;
137
138 int nFits = 0;
139
140 //printf("S1 size = %i\n", pS1->size());
141
142 // Simple double loop
143 for (UInt_t i = 0; i<pS1->size(); ++i) {
144 const StablePart &s1 = pS1->at(i);
145
146 const reco::Track * t1 = s1.track();
147 const TrackParameters &trkPar1 = trkPars1.at(i);
148
149 UInt_t j;
150 if (iStables1_ == iStables2_)
151 j = i+1;
152 else
153 j = 0;
154
155 TrajectoryStateTransform tsTransform;
156
157 FreeTrajectoryState initialState1 = tsTransform.initialFreeState(*s1.track(),&*magneticField);
158
159 for (; j<pS2->size(); ++j) {
160 const StablePart &s2 = pS2->at(j);
161
162 //Do fast helix fit to check if there's any hope
163 const reco::Track * t2 = s2.track();
164 const TrackParameters &trkPar2 = trkPars2.at(j);
165
166 int trackCharge = t1->charge() + t2->charge();
167
168 double dR0 = 0.0;
169
170 if (!applyChargeConstraint_ || trackCharge==0) {
171 FreeTrajectoryState initialState2 = tsTransform.initialFreeState(*s2.track(),&*magneticField);
172 helixIntersector.calculate(initialState1, initialState2);
173 if (helixIntersector.status())
174 dR0 = helixIntersector.crossingPoint().perp();
175 }
176
177 int fitStatus = 0;
178 MultiVertexFitterD fit;
179 if ( (!applyChargeConstraint_ || trackCharge==0) && (!useRhoMin_ || dR0 > rhoMin_) ) {
180
181 // Vertex fit now, possibly with conversion constraint
182 nFits++;
183
184 fit.init(bfield); // Reset to the magnetic field from the event setup
185
186 fit.addTrack(*trkPar1.pars(),*trkPar1.cMat(),1,s1.mass(),MultiVertexFitterD::VERTEX_1);
187 fit.addTrack(*trkPar2.pars(),*trkPar2.cMat(),2,s2.mass(),MultiVertexFitterD::VERTEX_1);
188 if (convConstraint3D_) {
189 fit.conversion_3d(MultiVertexFitterD::VERTEX_1);
190 //printf("applying 3d conversion constraint\n");
191 }
192 else if (convConstraint_) {
193 fit.conversion_2d(MultiVertexFitterD::VERTEX_1);
194 //printf("applying 2d conversion constraint\n");
195 }
196 //initialize primary vertex parameters in the fitter
197 if (usePVertex_) {
198 float vErr[3][3];
199 for (UInt_t vi=0; vi<3; ++vi)
200 for (UInt_t vj=0; vj<3; ++vj)
201 vErr[vi][vj] = vertex->covariance(vi,vj);
202
203 fit.setPrimaryVertex(vertex->x(),vertex->y(),vertex->z());
204 fit.setPrimaryVertexError(vErr);
205 }
206
207 fitStatus = fit.fit();
208
209 }
210
211 if (fitStatus) {
212 DecayPart *d = new DecayPart(oPid_,DecayPart::Fast);
213
214 // Update temporarily some of the quantities (prob, chi2, nDoF, mass, lxy, pt, fourMomentum)
215 d->setProb(fit.prob());
216 d->setChi2(fit.chisq());
217 d->setNdof(fit.ndof());
218
219 FourVector p4Fitted(0.,0.,0.,0.);
220 p4Fitted += fit.getTrackP4(1);
221 p4Fitted += fit.getTrackP4(2);
222 d->setFourMomentum(p4Fitted);
223 d->setPosition(fit.getVertex (MultiVertexFitterD::VERTEX_1));
224 d->setError (fit.getErrorMatrix(MultiVertexFitterD::VERTEX_1));
225 float mass, massErr;
226 const int trksIds[2] = { 1, 2 };
227 mass = fit.getMass(2,trksIds,massErr);
228
229 ThreeVector p3Fitted(p4Fitted.px(), p4Fitted.py(), p4Fitted.pz());
230
231 // Get decay length in xy plane
232 float dl, dlErr;
233 dl = fit.getDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
234 p3Fitted, dlErr);
235
236 // Get Z decay length
237 float dlz, dlzErr;
238 dlz = fit.getZDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
239 p3Fitted, dlzErr);
240
241 // Get impact parameter
242 float dxy, dxyErr;
243 dxy = fit.getImpactPar(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
244 p3Fitted, dxyErr);
245
246 BasePartPtr ptr1(hStables1,i);
247 BasePartPtr ptr2(hStables2,j);
248
249 StableData c1(fit.getTrackP4(1).px(),fit.getTrackP4(1).py(), fit.getTrackP4(1).pz(), ptr1);
250 StableData c2(fit.getTrackP4(2).px(),fit.getTrackP4(2).py(), fit.getTrackP4(2).pz(), ptr2);
251
252 const ThreeVector vtxPos = fit.getVertex(MultiVertexFitterD::VERTEX_1);
253 const ThreeVector trkMom1(fit.getTrackP4(1).px(),
254 fit.getTrackP4(1).py(),
255 fit.getTrackP4(1).pz());
256 const ThreeVector trkMom2(fit.getTrackP4(2).px(),
257 fit.getTrackP4(2).py(),
258 fit.getTrackP4(2).pz());
259
260 // Build corrected HitPattern for StableData, removing hits before the fit vertex
261 if (useHitDropper_) {
262 reco::HitPattern hits1 = dropper->CorrectedHitsAOD(s1.track(), vtxPos, trkMom1, dlErr, dlzErr);
263 reco::HitPattern hits2 = dropper->CorrectedHitsAOD(s2.track(), vtxPos, trkMom2, dlErr, dlzErr);
264
265 c1.SetHits(hits1);
266 c2.SetHits(hits2);
267 c1.SetHitsFilled();
268 c2.SetHitsFilled();
269
270 reco::HitPattern sharedHits = dropper->SharedHits(s1.track(),s2.track());
271 d->setSharedHits(sharedHits);
272 }
273
274 d->addStableChild(c1);
275 d->addStableChild(c2);
276
277
278 d->setFittedMass (mass);
279 d->setFittedMassError(massErr);
280
281 d->setLxy(dl);
282 d->setLxyError(dlErr);
283 d->setLxyToPv(dl);
284 d->setLxyToPvError(dlErr);
285
286 d->setLz(dlz);
287 d->setLzError(dlzErr);
288 d->setLzToPv(dlz);
289 d->setLzToPvError(dlzErr);
290
291 d->setDxy(dxy);
292 d->setDxyError(dxyErr);
293 d->setDxyToPv(dxy);
294 d->setDxyToPvError(dxyErr);
295
296 if (usePVertex_)
297 d->setPrimaryVertex(vPtr);
298
299 // Put the result into our collection
300 pD->push_back(*d);
301
302 delete d;
303 }
304 }
305 }
306
307 //printf("nConversionFits = %i\n",nFits);
308
309 // Write the collection even if it is empty
310 if (0) {
311 cout << " ProducerConversions::produce - " << pD->size() << " entries collection created -"
312 << " (Pid: " << oPid_ << ")\n";
313 }
314 evt.put(pD);
315 }
316
317 //define this as a plug-in
318 DEFINE_FWK_MODULE(ProducerConversions);