ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/src/ProducerV2SS.cc
Revision: 1.23
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.22: +4 -7 lines
Error occurred while calculating annotation data.
Log Message:
Initial 5 version.

File Contents

# Content
1 // $Id: ProducerV2SS.cc,v 1.22 2010/06/08 20:17:30 bendavid Exp $
2
3 #include "DataFormats/Common/interface/Handle.h"
4 #include "DataFormats/TrackReco/interface/Track.h"
5 #include "DataFormats/TrackReco/interface/TrackFwd.h"
6 #include "MagneticField/Engine/interface/MagneticField.h"
7 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
8 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
9 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
10 #include "MitEdm/Producers/interface/HitDropperRecord.h"
11 #include "MitEdm/Producers/interface/HitDropper.h"
12 #include "MitEdm/DataFormats/interface/Types.h"
13 #include "MitEdm/DataFormats/interface/Collections.h"
14 #include "MitEdm/DataFormats/interface/DecayPart.h"
15 #include "MitEdm/DataFormats/interface/StablePart.h"
16 #include "MitEdm/DataFormats/interface/StableData.h"
17 #include "MitEdm/VertexFitInterface/interface/MvfInterface.h"
18 #include "MitEdm/VertexFitInterface/interface/TrackParameters.h"
19 #include "MitEdm/Producers/interface/ProducerV2SS.h"
20 #include <TMath.h>
21
22 using namespace std;
23 using namespace edm;
24 using namespace mitedm;
25 using namespace mithep;
26
27 //--------------------------------------------------------------------------------------------------
28 ProducerV2SS::ProducerV2SS(const ParameterSet& cfg) :
29 ProducerD2SS(cfg),
30 rhoMin_ (cfg.getUntrackedParameter<double>("minRadius", 0.0)),
31 massMin_ (cfg.getUntrackedParameter<double>("minMass", 0.0)),
32 massMax_ (cfg.getUntrackedParameter<double>("maxMass", 3.0)),
33 dZMax_ (cfg.getUntrackedParameter<double>("maxZDistance",5.0)),
34 useHitDropper_(cfg.getUntrackedParameter<bool>("useHitDropper",true)),
35 applyChargeConstraint_(cfg.getUntrackedParameter<bool> ("applyChargeConstraint",false))
36 {
37 // Constructor.
38 }
39
40 //--------------------------------------------------------------------------------------------------
41 ProducerV2SS::~ProducerV2SS()
42 {
43 // Destructor.
44 }
45
46 //--------------------------------------------------------------------------------------------------
47 void ProducerV2SS::produce(Event &evt, const EventSetup &setup)
48 {
49 // Produce the output collection.
50
51 auto_ptr<DecayPartCol> pD(new DecayPartCol());
52
53 // First input collection
54 Handle<StablePartCol> hStables1;
55 Handle<StablePartCol> hStables2;
56 if (!GetProduct(iStables1_, hStables1, evt) ) {
57 cout << "Couldn't get in collection in Producer V2SS" << endl;
58 evt.put(pD);
59 return;
60 }
61 const StablePartCol *pS1 = hStables1.product();
62
63 // Second input collection
64 if(!GetProduct(iStables2_, hStables2, evt)) {
65 cout << "Couldn't get in collection in Producer V2SS" << endl;
66 evt.put(pD);
67 return;
68 }
69 const StablePartCol *pS2 = hStables2.product();
70
71 //get hit dropper
72 ESHandle<HitDropper> hDropper;
73 const HitDropper *dropper = 0;
74 if (useHitDropper_) {
75 setup.get<HitDropperRecord>().get("HitDropper",hDropper);
76 dropper = hDropper.product();
77 }
78
79 //Get Magnetic Field from event setup, taking value at (0,0,0)
80 edm::ESHandle<MagneticField> magneticField;
81 setup.get<IdealMagneticFieldRecord>().get(magneticField);
82 const double bfield = magneticField->inTesla(GlobalPoint(0.,0.,0.)).z();
83
84 edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
85 setup.get<TransientTrackRecord>().get("TransientTrackBuilder",hTransientTrackBuilder);
86 const TransientTrackBuilder *transientTrackBuilder = hTransientTrackBuilder.product();
87
88 //construct intermediate collection of TrackParameters in mvf format for vertex fit
89 std::vector<TrackParameters> trkPars1;
90 for (UInt_t i = 0; i<pS1->size(); ++i) {
91 const reco::Track *t = pS1->at(i).track();
92 const reco::TransientTrack ttrk = transientTrackBuilder->build(t);
93 TrackParameters cmsTrk(ttrk);
94 TrackParameters mvfTrk = cmsTrk.mvfTrack();
95 trkPars1.push_back(mvfTrk);
96 }
97
98 std::vector<TrackParameters> trkPars2;
99 if (iStables1_ == iStables2_)
100 trkPars2 = trkPars1;
101 else for (UInt_t i = 0; i<pS2->size(); ++i) {
102 const reco::Track *t = pS2->at(i).track();
103 const reco::TransientTrack ttrk = transientTrackBuilder->build(t);
104 TrackParameters cmsTrk(ttrk);
105 TrackParameters mvfTrk = cmsTrk.mvfTrack();
106 trkPars2.push_back(mvfTrk);
107 }
108
109 // -----------------------------------------------------------------------------------------------
110 // Simple double loop
111 // -----------------------------------------------------------------------------------------------
112 if (0)
113 cout << "Starting V finder loop" << endl;
114
115 ClosestApproachInRPhi helixIntersector;
116
117
118 //sX_y: X= pion or proton collection.
119 //i, j = 2 loop particles.
120 //ex.: s1_i and s2_i are same particle as pion and proton
121 for (UInt_t i=0; i<pS1->size(); ++i) {
122 const StablePart &s1 = pS1->at(i);
123 const TrackParameters &trkPar1 = trkPars1.at(i);
124
125 //const reco::Track * t1 = s1.track();
126
127 UInt_t j;
128 if (iStables1_ == iStables2_)
129 j = i+1;
130 else
131 j = 0;
132
133 FreeTrajectoryState initialState1 = trajectoryStateTransform::initialFreeState(*s1.track(),&*magneticField);
134
135 for (; j<pS2->size(); ++j) {
136 const StablePart &s2 = pS2->at(j);
137 const TrackParameters &trkPar2 = trkPars2.at(j);
138
139 if( applyChargeConstraint_ && (s1.charge() + s2.charge() != 0) ) continue;
140
141 // do fast helix fit to check if there's any hope
142 //
143 //const reco::Track * t2 = s2.track();
144
145 double dZ0 = -999;
146 double dR0 = -999;
147 double mass0 = 0.0;
148
149 FreeTrajectoryState initialState2 = trajectoryStateTransform::initialFreeState(*s2.track(),&*magneticField);
150 helixIntersector.calculate(initialState1, initialState2);
151 if (helixIntersector.status()) {
152 dZ0 = fabs(helixIntersector.points().first.z() - helixIntersector.points().second.z());
153 dR0 = helixIntersector.crossingPoint().perp();
154
155 GlobalVector v1, v2;
156 v1 = helixIntersector.trajectoryParameters().first.momentum();
157 v2 = helixIntersector.trajectoryParameters().second.momentum();
158
159 double e1 = sqrt(v1.mag2()+s1.mass()*s1.mass());
160 double x1 = v1.x();
161 double y1 = v1.y();
162 double z1 = v1.z();
163
164 double e2 = sqrt(v2.mag2()+s2.mass()*s2.mass());
165 double x2 = v2.x();
166 double y2 = v2.y();
167 double z2 = v2.z();
168
169 FourVector sum(x1+x2, y1+y2, z1+z2, e1+e2);
170
171 mass0 = sqrt(sum.M2());
172 }
173
174
175 // Basic cuts on helix intersection
176 if(mass0 > massMax_ || mass0<massMin_ || fabs(dZ0) > dZMax_ || dR0 < rhoMin_) continue;
177
178 // -------------------------------------------------------------------------------------------
179 // Do vertex fit for all pairs
180 // -------------------------------------------------------------------------------------------
181 mithep::MultiVertexFitterD fit;
182 fit.init(bfield); // Reset to the magnetic field from the event setup
183 fit.setChisqMax(100);
184
185 fit.addTrack(*trkPar1.pars(),*trkPar1.cMat(),1,s1.mass(),MultiVertexFitterD::VERTEX_1);
186 fit.addTrack(*trkPar2.pars(),*trkPar2.cMat(),2,s2.mass(),MultiVertexFitterD::VERTEX_1);
187
188
189 if (fit.fit()){
190
191 DecayPart *d = new DecayPart(oPid_,DecayPart::Fast);
192
193 // Update temporarily some of the quantities (prob, chi2, nDoF, mass, lxy, pt, fourMomentum)
194 d->setProb(fit.prob());
195 d->setChi2(fit.chisq());
196 d->setNdof(fit.ndof());
197
198 FourVector p4Fitted(0.,0.,0.,0.);
199 p4Fitted += fit.getTrackP4(1);
200 p4Fitted += fit.getTrackP4(2);
201 d->setFourMomentum(p4Fitted);
202 d->setPosition (fit.getVertex (MultiVertexFitterD::VERTEX_1));
203 d->setError (fit.getErrorMatrix(MultiVertexFitterD::VERTEX_1));
204 float mass, massErr;
205 const int trksIds[2] = { 1, 2 };
206 mass = fit.getMass(2,trksIds,massErr);
207
208 ThreeVector p3Fitted(p4Fitted.px(), p4Fitted.py(), p4Fitted.pz());
209
210 //Get decay length in xy plane
211 float dl, dlErr;
212 dl = fit.getDecayLength (MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
213 p3Fitted, dlErr);
214
215 //Get Z decay length
216 float dlz, dlzErr;
217 dlz = fit.getZDecayLength(MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
218 p3Fitted, dlzErr);
219
220 //get impact parameter
221 float dxy, dxyErr;
222 dxy = fit.getImpactPar (MultiVertexFitterD::PRIMARY_VERTEX, MultiVertexFitterD::VERTEX_1,
223 p3Fitted, dxyErr);
224
225 BasePartPtr ptr1(hStables1,i);
226 BasePartPtr ptr2(hStables2,j);
227
228 StableData c1(fit.getTrackP4(1).px(),fit.getTrackP4(1).py(), fit.getTrackP4(1).pz(), ptr1);
229 StableData c2(fit.getTrackP4(2).px(),fit.getTrackP4(2).py(), fit.getTrackP4(2).pz(), ptr2);
230
231 const ThreeVector vtxPos = fit.getVertex(MultiVertexFitterD::VERTEX_1);
232 const ThreeVector trkMom1(fit.getTrackP4(1).px(),fit.getTrackP4(1).py(),
233 fit.getTrackP4(1).pz());
234 const ThreeVector trkMom2(fit.getTrackP4(2).px(),fit.getTrackP4(2).py(),
235 fit.getTrackP4(2).pz());
236
237 //build corrected HitPattern for StableData, removing hits before the fit vertex
238 if (useHitDropper_) {
239 std::pair<reco::HitPattern,uint> hits1 = dropper->CorrectedHitsAOD(s1.track(),
240 vtxPos,
241 trkMom1,
242 dlErr,
243 dlzErr);
244 std::pair<reco::HitPattern,uint> hits2 = dropper->CorrectedHitsAOD(s2.track(),
245 vtxPos,
246 trkMom2,
247 dlErr,
248 dlzErr);
249
250 c1.SetHits(hits1.first);
251 c2.SetHits(hits2.first);
252 c1.SetHitsFilled();
253 c2.SetHitsFilled();
254 c1.SetNWrongHits(hits1.second);
255 c2.SetNWrongHits(hits2.second);
256
257 reco::HitPattern sharedHits = dropper->SharedHits(s1.track(),s2.track());
258 d->setSharedHits(sharedHits);
259 }
260
261 d->addStableChild (c1);
262 d->addStableChild (c2);
263 d->setFittedMass (mass);
264 d->setFittedMassError(massErr);
265 d->setLxy (dl);
266 d->setLxyError (dlErr);
267 d->setLxyToPv (dl);
268 d->setLxyToPvError (dlErr);
269 d->setLz (dlz);
270 d->setLzError (dlzErr);
271 d->setLzToPv (dlz);
272 d->setLzToPvError (dlzErr);
273 d->setDxy (dxy);
274 d->setDxyError (dxyErr);
275 d->setDxyToPv (dxy);
276 d->setDxyToPvError (dxyErr);
277
278 // put the result into our collection
279 pD->push_back(*d);
280 delete d;
281 } //done processing fit
282 } //end j loop
283 } //end i loop
284
285 // -----------------------------------------------------------------------------------------------
286 // Write the collection even if it is empty
287 // -----------------------------------------------------------------------------------------------
288 if (0)
289 cout << " V2SS::produce - " << pD->size() << " entries collection created -"
290 << " (Pid: " << oPid_ << ")\n";
291 evt.put(pD);
292 }
293
294 // define this as a plug-in
295 DEFINE_FWK_MODULE(ProducerV2SS);