ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/src/ProducerConversions.cc
Revision: 1.18
Committed: Mon Nov 2 22:56:18 2009 UTC (15 years, 6 months ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_012b, Mit_012a, Mit_012
Changes since 1.17: +9 -2 lines
Log Message:
Use magnetic field from event setup instead of hardcoded 3.8T

File Contents

# Content
1 // $Id: ProducerConversions.cc,v 1.17 2009/10/04 12:49:26 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 "MitEdm/Producers/interface/HitDropperRecord.h"
12 #include "MitEdm/Producers/interface/HitDropper.h"
13 #include "MitEdm/DataFormats/interface/Types.h"
14 #include "MitEdm/DataFormats/interface/Collections.h"
15 #include "MitEdm/DataFormats/interface/DecayPart.h"
16 #include "MitEdm/DataFormats/interface/StablePart.h"
17 #include "MitEdm/DataFormats/interface/StableData.h"
18 #include "MitEdm/VertexFitInterface/interface/MvfInterface.h"
19 #include "MitEdm/VertexFitInterface/interface/HisInterface.h"
20
21 using namespace std;
22 using namespace edm;
23 using namespace reco;
24 using namespace mitedm;
25 using namespace mithep;
26
27 //--------------------------------------------------------------------------------------------------
28 ProducerConversions::ProducerConversions(const ParameterSet& cfg) :
29 BaseCandProducer (cfg),
30 iStables1_ (cfg.getUntrackedParameter<string>("iStables1","")),
31 iStables2_ (cfg.getUntrackedParameter<string>("iStables2","")),
32 iPVertexes_ (cfg.getUntrackedParameter<string>("iPVertexes","offlinePrimaryVerticesWithBS")),
33 usePVertex_ (cfg.getUntrackedParameter<bool> ("usePVertex",true)),
34 convConstraint_ (cfg.getUntrackedParameter<bool> ("convConstraint",false)),
35 convConstraint3D_(cfg.getUntrackedParameter<bool> ("convConstraint3D",true)),
36 rhoMin_ (cfg.getUntrackedParameter<double>("rhoMin",0.0)),
37 useHitDropper_ (cfg.getUntrackedParameter<bool> ("useHitDropper",true))
38 {
39 // Constructor.
40
41 produces<DecayPartCol>();
42 }
43
44 //--------------------------------------------------------------------------------------------------
45 ProducerConversions::~ProducerConversions()
46 {
47 // Destructor.
48 }
49
50 //--------------------------------------------------------------------------------------------------
51 void ProducerConversions::produce(Event &evt, const EventSetup &setup)
52 {
53 // Produce our DecayPartCol.
54
55 // First input collection
56 Handle<StablePartCol> hStables1;
57 if (!GetProduct(iStables1_, hStables1, evt)) {
58 printf("Stable collection 1 not found in ProducerConversions\n");
59 return;
60 }
61 const StablePartCol *pS1 = hStables1.product();
62 // Second input collection
63 Handle<StablePartCol> hStables2;
64 if (!GetProduct(iStables2_, hStables2, evt)) {
65 printf("Stable collection 2 not found in ProducerConversions\n");
66 return;
67 }
68 const StablePartCol *pS2 = hStables2.product();
69
70 const reco::Vertex *vertex = 0;
71 mitedm::VertexPtr vPtr;
72 if (usePVertex_) {
73 // Primary vertex collection
74 Handle<reco::VertexCollection> hVertexes;
75 if (!GetProduct(iPVertexes_, hVertexes, evt))
76 return;
77 const reco::VertexCollection *pV = hVertexes.product();
78
79 //Choose the primary vertex with the largest number of tracks
80 UInt_t maxTracks=0;
81 for (UInt_t i=0; i<pV->size(); ++i) {
82 const reco::Vertex &v = pV->at(i);
83 UInt_t nTracks = v.tracksSize();
84 if (nTracks >= maxTracks) {
85 maxTracks = nTracks;
86 vertex = &v;
87 vPtr = mitedm::VertexPtr(hVertexes,i);
88 }
89 }
90 }
91
92 // Get hit dropper
93 ESHandle<HitDropper> hDropper;
94 setup.get<HitDropperRecord>().get("HitDropper",hDropper);
95 const HitDropper *dropper = hDropper.product();
96
97 //Get Magnetic Field from event setup, taking value at (0,0,0)
98 edm::ESHandle<MagneticField> magneticField;
99 setup.get<IdealMagneticFieldRecord>().get(magneticField);
100 const double bfield = magneticField->inTesla(GlobalPoint(0.,0.,0.)).z();
101
102 // Create the output collection
103 auto_ptr<DecayPartCol> pD(new DecayPartCol());
104
105 // Simple double loop
106 for (UInt_t i = 0; i<pS1->size(); ++i) {
107 const StablePart &s1 = pS1->at(i);
108
109 UInt_t j;
110 if (iStables1_ == iStables2_)
111 j = i+1;
112 else
113 j = 0;
114
115 for (; j<pS2->size(); ++j) {
116 const StablePart &s2 = pS2->at(j);
117
118 //Do fast helix fit to check if there's any hope
119 const reco::Track * t1 = s1.track();
120 const reco::Track * t2 = s2.track();
121
122 int trackCharge = t1->charge() + t2->charge();
123 double dR0 = 0.0;
124
125 if (trackCharge==0) {
126 HisInterface hisInt(t1,t2);
127 const mithep::HelixIntersector *his = hisInt.hISector();
128 dR0 = sqrt(his->IntersectionI(0).Location().X()*his->IntersectionI(0).Location().X()+
129 his->IntersectionI(0).Location().Y()*his->IntersectionI(0).Location().Y());
130 }
131
132
133 int fitStatus = 0;
134 MultiVertexFitterD fit;
135 if (trackCharge==0 && dR0 > rhoMin_) {
136
137 // Vertex fit now, possibly with conversion constraint
138
139 fit.init(bfield); // Reset to the magnetic field from the event setup
140 MvfInterface fitInt(&fit);
141 fitInt.addTrack(s1.track(),1,s1.mass(),MultiVertexFitterD::VERTEX_1);
142 fitInt.addTrack(s2.track(),2,s2.mass(),MultiVertexFitterD::VERTEX_1);
143 if (convConstraint3D_) {
144 fit.conversion_3d(MultiVertexFitterD::VERTEX_1);
145 //printf("applying 3d conversion constraint\n");
146 }
147 else if (convConstraint_) {
148 fit.conversion_2d(MultiVertexFitterD::VERTEX_1);
149 //printf("applying 2d conversion constraint\n");
150 }
151 //initialize primary vertex parameters in the fitter
152 if (usePVertex_) {
153 float vErr[3][3];
154 for (UInt_t vi=0; vi<3; ++vi)
155 for (UInt_t vj=0; vj<3; ++vj)
156 vErr[vi][vj] = vertex->covariance(vi,vj);
157
158 fit.setPrimaryVertex(vertex->x(),vertex->y(),vertex->z());
159 fit.setPrimaryVertexError(vErr);
160 }
161
162 fitStatus = fit.fit();
163 }
164
165 if (fitStatus) {
166 DecayPart *d = new DecayPart(oPid_,DecayPart::Fast);
167
168 // Update temporarily some of the quantities (prob, chi2, nDoF, mass, lxy, pt, fourMomentum)
169 d->setProb(fit.prob());
170 d->setChi2(fit.chisq());
171 d->setNdof(fit.ndof());
172
173 FourVector p4Fitted(0.,0.,0.,0.);
174 p4Fitted += fit.getTrackP4(1);
175 p4Fitted += fit.getTrackP4(2);
176 d->setFourMomentum(p4Fitted);
177 d->setPosition(fit.getVertex (MultiVertexFitterD::VERTEX_1));
178 d->setError (fit.getErrorMatrix(MultiVertexFitterD::VERTEX_1));
179 float mass, massErr;
180 const int trksIds[2] = { 1, 2 };
181 mass = fit.getMass(2,trksIds,massErr);
182
183 ThreeVector p3Fitted(p4Fitted.px(), p4Fitted.py(), p4Fitted.pz());
184
185 // Get decay length in xy plane
186 float dl, dlErr;
187 dl = fit.getDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
188 p3Fitted, dlErr);
189
190 // Get Z decay length
191 float dlz, dlzErr;
192 dlz = fit.getZDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
193 p3Fitted, dlzErr);
194
195 // Get impact parameter
196 float dxy, dxyErr;
197 dxy = fit.getImpactPar(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
198 p3Fitted, dxyErr);
199
200 BasePartPtr ptr1(hStables1,i);
201 BasePartPtr ptr2(hStables2,j);
202
203 StableData c1(fit.getTrackP4(1).px(),fit.getTrackP4(1).py(), fit.getTrackP4(1).pz(), ptr1);
204 StableData c2(fit.getTrackP4(2).px(),fit.getTrackP4(2).py(), fit.getTrackP4(2).pz(), ptr2);
205
206 const ThreeVector vtxPos = fit.getVertex(MultiVertexFitterD::VERTEX_1);
207 const ThreeVector trkMom1(fit.getTrackP4(1).px(),
208 fit.getTrackP4(1).py(),
209 fit.getTrackP4(1).pz());
210 const ThreeVector trkMom2(fit.getTrackP4(2).px(),
211 fit.getTrackP4(2).py(),
212 fit.getTrackP4(2).pz());
213
214 // Build corrected HitPattern for StableData, removing hits before the fit vertex
215 if (useHitDropper_) {
216 reco::HitPattern hits1 = dropper->CorrectedHitsAOD(s1.track(), vtxPos, trkMom1, dlErr, dlzErr);
217 reco::HitPattern hits2 = dropper->CorrectedHitsAOD(s2.track(), vtxPos, trkMom2, dlErr, dlzErr);
218
219 c1.SetHits(hits1);
220 c2.SetHits(hits2);
221 c1.SetHitsFilled();
222 c2.SetHitsFilled();
223 }
224
225 d->addStableChild(c1);
226 d->addStableChild(c2);
227
228
229 d->setFittedMass (mass);
230 d->setFittedMassError(massErr);
231
232 d->setLxy(dl);
233 d->setLxyError(dlErr);
234 d->setLxyToPv(dl);
235 d->setLxyToPvError(dlErr);
236
237 d->setLz(dlz);
238 d->setLzError(dlzErr);
239 d->setLzToPv(dlz);
240 d->setLzToPvError(dlzErr);
241
242 d->setDxy(dxy);
243 d->setDxyError(dxyErr);
244 d->setDxyToPv(dxy);
245 d->setDxyToPvError(dxyErr);
246
247 if (usePVertex_)
248 d->setPrimaryVertex(vPtr);
249
250 // Put the result into our collection
251 pD->push_back(*d);
252
253 delete d;
254 }
255 }
256 }
257
258 // Write the collection even if it is empty
259 if (0) {
260 cout << " ProducerConversions::produce - " << pD->size() << " entries collection created -"
261 << " (Pid: " << oPid_ << ")\n";
262 }
263 evt.put(pD);
264 }
265
266 //define this as a plug-in
267 DEFINE_FWK_MODULE(ProducerConversions);