#include <MyCode/TrackerOnlyConversionProducer/src/TrackerOnlyConversionProducer.cc>
Public Member Functions | |
TrackerOnlyConversionProducer (const edm::ParameterSet &) | |
~TrackerOnlyConversionProducer () | |
Private Types | |
typedef math::XYZPointD | Point |
typedef std::vector< Point > | PointCollection |
Private Member Functions | |
virtual void | beginJob (const edm::EventSetup &) |
virtual void | beginRun (const edm::Run &, const edm::EventSetup &) |
virtual void | endJob () |
virtual void | endRun (const edm::Run &, const edm::EventSetup &) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
bool | allowSingleLeg_ |
edm::InputTag | bcBarrelCollection_ |
edm::InputTag | bcEndcapCollection_ |
std::string | ConvertedPhotonCollection_ |
double | deltaCotTheta_ |
double | deltaPhi_ |
double | dEtaTkBC_ |
double | dPhiTkBC_ |
double | energyBC_ |
double | halfWayEta_ |
double | halfWayPhi_ |
double | maxChi2Left_ |
double | maxChi2Right_ |
double | maxHitsLeft_ |
double | maxHitsRight_ |
std::vector< edm::InputTag > | src_ |
Implementation: <Notes on="" implementation>="">
Definition at line 115 of file TrackerOnlyConversionProducer.cc.
typedef math::XYZPointD TrackerOnlyConversionProducer::Point [private] |
Definition at line 128 of file TrackerOnlyConversionProducer.cc.
typedef std::vector<Point> TrackerOnlyConversionProducer::PointCollection [private] |
Definition at line 129 of file TrackerOnlyConversionProducer.cc.
TrackerOnlyConversionProducer::TrackerOnlyConversionProducer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 163 of file TrackerOnlyConversionProducer.cc.
References allowSingleLeg_, bcBarrelCollection_, bcEndcapCollection_, ConvertedPhotonCollection_, deltaCotTheta_, deltaPhi_, dEtaTkBC_, dPhiTkBC_, energyBC_, edm::ParameterSet::getParameter(), halfWayEta_, maxChi2Left_, maxChi2Right_, maxHitsLeft_, maxHitsRight_, and src_.
00164 { 00165 src_ = iConfig.getParameter<std::vector<edm::InputTag> >("src"); 00166 00167 00168 bcBarrelCollection_ = iConfig.getParameter<edm::InputTag>("bcBarrelCollection"); 00169 bcEndcapCollection_ = iConfig.getParameter<edm::InputTag>("bcEndcapCollection"); 00170 00171 halfWayEta_ = iConfig.getParameter<double>("HalfwayEta");//open angle to search track matches with BC 00172 00173 //Track-cluster matching eta and phi cuts 00174 dEtaTkBC_ = iConfig.getParameter<double>("dEtaTrackBC");//TODO research on cut endcap/barrel 00175 dPhiTkBC_ = iConfig.getParameter<double>("dPhiTrackBC"); 00176 00177 energyBC_ = iConfig.getParameter<double>("EnergyBC");//BC energy cut 00178 00179 //Track cuts on left right track: at least one leg reaches ECAL 00180 //Left track: must exist, must reach Ecal and match BC, so loose cut on Chi2 and tight on hits 00181 //Right track: not necessary to exist (if allowSingleLeg_), not necessary to reach ECAL or match BC, so tight cut on Chi2 and loose on hits 00182 maxChi2Left_ = iConfig.getParameter<double>("MaxChi2Left"); 00183 maxChi2Right_ = iConfig.getParameter<double>("MaxChi2Right"); 00184 maxHitsLeft_ = iConfig.getParameter<int>("MaxHitsLeft"); 00185 maxHitsRight_ = iConfig.getParameter<int>("MaxHitsRight"); 00186 00187 //Track Open angle cut on delta cot(theta) and delta phi 00188 deltaCotTheta_ = iConfig.getParameter<double>("DeltaCotTheta"); 00189 deltaPhi_ = iConfig.getParameter<double>("DeltaPhi"); 00190 00191 // if allow single track collection, by default False 00192 allowSingleLeg_ = iConfig.getParameter<bool>("AllowSingleLeg"); 00193 00194 //output 00195 ConvertedPhotonCollection_ = iConfig.getParameter<std::string>("convertedPhotonCollection"); 00196 00197 produces< reco::ConversionCollection >(ConvertedPhotonCollection_); 00198 00199 }
TrackerOnlyConversionProducer::~TrackerOnlyConversionProducer | ( | ) |
Definition at line 202 of file TrackerOnlyConversionProducer.cc.
00203 { 00204 00205 // do anything here that needs to be done at desctruction time 00206 // (e.g. close files, deallocate resources etc.) 00207 00208 }
void TrackerOnlyConversionProducer::beginJob | ( | const edm::EventSetup & | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 492 of file TrackerOnlyConversionProducer.cc.
void TrackerOnlyConversionProducer::beginRun | ( | const edm::Run & | run, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 498 of file TrackerOnlyConversionProducer.cc.
void TrackerOnlyConversionProducer::endRun | ( | const edm::Run & | run, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
void TrackerOnlyConversionProducer::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
for( std::vector<edm::Ref<reco::TrackCollection> >const_iterator ll = allTracks.begin(); ll != allTracks.end(); ++ ll ) {
Implements edm::EDProducer.
Definition at line 217 of file TrackerOnlyConversionProducer.cc.
References allowSingleLeg_, python::MinBiasTracking_cff::allTracks, alongMomentum, bcBarrelCollection_, bcEndcapCollection_, edm::RefVector< C, T, F >::begin(), ConvertedPhotonCollection_, kinem::delta_eta(), kinem::delta_phi(), deltaCotTheta_, deltaPhi_, dEtaTkBC_, dPhi(), dPhiTkBC_, empty, edm::RefVector< C, T, F >::end(), energyBC_, geometryDiff::epsilon, PV3DBase< T, PVType, FrameType >::eta(), edm::EventSetup::get(), edm::Event::getByLabel(), TrajectoryStateOnSurface::globalPosition(), halfWayEta_, edm::Ref< C, T, F >::isNonnull(), TrajectoryStateOnSurface::isValid(), map_phi2(), maxChi2Left_, maxChi2Right_, maxHitsLeft_, maxHitsRight_, TrajectoryStateTransform::outerStateOnSurface(), edm::ESHandle< T >::product(), PropagatorWithMaterial::propagate(), edm::RefVector< C, T, F >::push_back(), edm::PtrVector< T >::push_back(), edm::Event::put(), edm::RefVector< C, T, F >::reserve(), rot, edm::RefVector< C, T, F >::size(), size, src_, and funct::tan().
00218 { 00219 using namespace edm; 00220 00221 reco::ConversionCollection outputConvPhotonCollection; 00222 std::auto_ptr<reco::ConversionCollection> outputConvPhotonCollection_p(new reco::ConversionCollection); 00223 00224 //Read multiple track input collections 00225 std::vector<const reco::TrackCollection*> trackCollections; 00226 std::vector<edm::Handle<reco::TrackCollection> > trackCollectionHandles; 00227 for (unsigned ii = 0; ii<src_.size(); ++ii){ 00228 edm::Handle<reco::TrackCollection> temp_handle; 00229 if(iEvent.getByLabel(src_[ii],temp_handle)){//starting from 170 00230 trackCollectionHandles.push_back(temp_handle); 00231 } else { 00232 edm::LogWarning("TrackerOnlyConversionProducer") << "Collection reco::TrackCollection with label " << src_[ii] << " cannot be found, using empty collection of same type"; 00233 } 00234 } 00235 edm::Handle<edm::View<reco::CaloCluster> > bcBarrelHandle;//TODO error handling if no collection found 00236 iEvent.getByLabel( bcBarrelCollection_, bcBarrelHandle); 00237 00238 edm::Handle<edm::View<reco::CaloCluster> > bcEndcapHandle;//TODO check cluster type if BasicCluster or PFCluster 00239 iEvent.getByLabel( bcEndcapCollection_, bcEndcapHandle); 00240 00241 edm::ESHandle<TrackerGeometry> trackerGeomHandle; 00242 edm::ESHandle<MagneticField> magFieldHandle; 00243 00244 iSetup.get<TrackerDigiGeometryRecord>().get( trackerGeomHandle ); 00245 iSetup.get<IdealMagneticFieldRecord>().get( magFieldHandle ); 00246 00247 const TrackerGeometry* trackerGeom = trackerGeomHandle.product(); 00248 const MagneticField* magField = magFieldHandle.product();; 00249 00250 //pair up tracks and record unpaired but matched tracks 00251 // 00252 //1. combining all track collections 00253 reco::TrackRefVector allTracks; 00254 int total_tracks = 0; 00255 for (unsigned int ii = 0; ii<trackCollectionHandles.size(); ++ii) 00256 total_tracks += trackCollectionHandles[ii]->size(); 00257 allTracks.reserve(total_tracks);//for many tracks to save relocation time 00258 00259 for (unsigned int ii = 0; ii<trackCollectionHandles.size(); ++ii){ 00260 if (!trackCollectionHandles[ii]->empty()){ 00261 for (unsigned int jj = 0; jj<trackCollectionHandles[ii]->size(); ++jj){ 00262 edm::Ref<reco::TrackCollection> ref(trackCollectionHandles[ii], jj); 00263 if (ref.isNonnull()){ 00264 allTracks.push_back(ref);//TODO find a way to get vector directly from handle to avoid loop 00265 } 00266 } 00267 } 00268 } 00269 //2. select track by propagating 00270 //2.0 build Basic cluster geometry map to search in eta bounds for clusters 00271 std::multimap<double, reco::CaloClusterPtr> basicClusterPtrs; 00272 edm::Handle<edm::View<reco::CaloCluster> > bcHandle = bcBarrelHandle; 00273 for (unsigned jj = 0; jj < 2; ++jj ){ 00274 for (unsigned ii = 0; ii < bcHandle->size(); ++ii ) { 00275 if (bcHandle->ptrAt(ii)->energy()>energyBC_) 00276 basicClusterPtrs.insert(std::make_pair(bcHandle->ptrAt(ii)->position().eta(), bcHandle->ptrAt(ii))); 00277 } 00278 bcHandle = bcEndcapHandle; 00279 } 00280 00281 std::vector<math::XYZPoint> trackImpactPosition; 00282 trackImpactPosition.reserve(allTracks.size());//track impact position at ECAL 00283 std::vector<bool> trackValidECAL;//Does this track reach ECAL basic cluster (reach ECAL && match with BC) 00284 trackValidECAL.assign(allTracks.size(), false); 00285 00286 std::vector<reco::CaloClusterPtr> trackMatchedBC; 00287 reco::CaloClusterPtr empty_bc; 00288 trackMatchedBC.assign(allTracks.size(), empty_bc);//TODO find a better way to avoid copy constructor 00289 00290 std::multimap<double, int> trackInnerEta;//Track innermost state Eta map to TrackRef index, to be used in track pair sorting 00291 00292 //2.1 propagate all tracks into ECAL, record its eta and phi 00293 //for (std::vector<edm::Ref<reco::TrackCollection> >::const_iterator ref = allTracks.begin(); ref != allTracks.end(); ++ref){ 00294 for (reco::TrackRefVector::const_iterator ref = allTracks.begin(); ref != allTracks.end(); ++ref){ 00295 const TrackRef& tk_ref = *ref; 00296 00297 if (tk_ref->normalizedChi2() > maxChi2Left_ || tk_ref->found() < maxHitsLeft_ //track quality cut 00298 ) continue; 00299 00300 //map TrackRef to Eta 00301 trackInnerEta.insert(std::make_pair(tk_ref->innerMomentum().eta(), ref-allTracks.begin())); 00302 00303 PropagatorWithMaterial propag( alongMomentum, 0.000511, magField ); 00304 TrajectoryStateTransform transformer; 00305 ReferenceCountingPointer<Surface> ecalWall( 00306 new BoundCylinder( GlobalPoint(0.,0.,0.), TkRotation<float>(), 00307 SimpleCylinderBounds( 129, 129, -320.5, 320.5 ) ) ); 00308 const float epsilon = 0.001; 00309 Surface::RotationType rot; // unit rotation matrix 00310 const float barrelRadius = 129.f; 00311 const float barrelHalfLength = 270.9f; 00312 const float endcapRadius = 171.1f; 00313 const float endcapZ = 320.5f; 00314 ReferenceCountingPointer<BoundCylinder> theBarrel_(new BoundCylinder( Surface::PositionType(0,0,0), rot, 00315 SimpleCylinderBounds( barrelRadius-epsilon, barrelRadius+epsilon, -barrelHalfLength, barrelHalfLength))); 00316 ReferenceCountingPointer<BoundDisk> theNegativeEtaEndcap_( 00317 new BoundDisk( Surface::PositionType( 0, 0, -endcapZ), rot, 00318 SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon))); 00319 ReferenceCountingPointer<BoundDisk> thePositiveEtaEndcap_( 00320 new BoundDisk( Surface::PositionType( 0, 0, endcapZ), rot, 00321 SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon))); 00322 00323 //const TrajectoryStateOnSurface myTSOS = transformer.innerStateOnSurface(*(*ref), *trackerGeom, magField); 00324 const TrajectoryStateOnSurface myTSOS = transformer.outerStateOnSurface(*tk_ref, *trackerGeom, magField); 00325 TrajectoryStateOnSurface stateAtECAL; 00326 stateAtECAL = propag.propagate(myTSOS, *theBarrel_); 00327 if (!stateAtECAL.isValid() || ( stateAtECAL.isValid() && fabs(stateAtECAL.globalPosition().eta() ) >1.479 ) ) { 00328 //endcap propagator 00329 if (myTSOS.globalPosition().eta() > 0.) { 00330 stateAtECAL = propag.propagate(myTSOS, *thePositiveEtaEndcap_); 00331 } else { 00332 stateAtECAL = propag.propagate(myTSOS, *theNegativeEtaEndcap_); 00333 } 00334 } 00335 math::XYZPoint ew; 00336 if (stateAtECAL.isValid()){ 00337 ew = stateAtECAL.globalPosition(); 00338 } 00339 trackImpactPosition.push_back(ew);//for invalid state, it will be count as invalid, so before read it out, check trackValidECAL[] 00340 if (!stateAtECAL.isValid()){ continue;} 00341 00342 const double track_eta = ew.eta(); 00343 const double track_phi = ew.phi(); 00344 00345 //2.2 check matching with BC 00346 reco::CaloClusterPtr closest_bc; 00347 double min_eta = 999., min_phi = 999.; 00348 for (std::multimap<double, reco::CaloClusterPtr>::iterator bc = basicClusterPtrs.lower_bound(track_eta - halfWayEta_); 00349 bc != basicClusterPtrs.upper_bound(track_eta + halfWayEta_); ++bc){//use eta map to select possible BC collection then loop in 00350 const reco::CaloClusterPtr& ebc = bc->second; 00351 const double delta_eta = track_eta-(ebc->position().eta()); 00352 const double delta_phi = map_phi2(track_phi-(ebc->position().phi())); 00353 if (fabs(delta_eta)<dEtaTkBC_ && fabs(delta_phi)<dPhiTkBC_){ 00354 if (fabs(min_eta)>fabs(delta_eta) && fabs(min_phi)>fabs(delta_phi)){//take the closest to track BC 00355 min_eta = delta_eta; 00356 min_phi = delta_phi; 00357 closest_bc = bc->second; 00358 //TODO check if min_eta>delta_eta but min_phi<delta_phi 00359 } 00360 } 00361 } 00362 if (min_eta < 999.){//this track matched a BC 00363 trackMatchedBC[ref-allTracks.begin()] = closest_bc; 00364 trackValidECAL[ref-allTracks.begin()] = true; 00365 } 00366 } 00367 //3. pair up : to select one track from matched EBC, and select the other track to fit cot theta and phi open angle cut 00368 //TODO it is k-Closest pair of point problem 00369 reco::VertexCollection vertexs; 00370 std::vector<bool> isPaired; 00371 isPaired.assign(allTracks.size(), false); 00373 for( reco::TrackRefVector::const_iterator ll = allTracks.begin(); ll != allTracks.end(); ++ ll ) { 00374 //Level 1 loop, in all tracks matched with ECAL 00375 if (!trackValidECAL[ll-allTracks.begin()] //this Left leg should have valid BC 00376 || (*ll)->d0()*(*ll)->charge()<0. //d0*charge for converted photons 00377 || isPaired[ll-allTracks.begin()]) //this track should not be in another pair 00378 continue; 00379 const double left_eta = (*ll)->innerMomentum().eta(); 00380 bool found_right = false;//check if right leg found, if no but allowSingleLeg_, go build a conversion with left leg 00381 std::map<double, int> right_candidates;//store all right legs passed the cut (theta and ref pair) 00382 00383 //select right leg candidates, which passed the cuts 00384 for (std::multimap<double, int>::iterator rr = trackInnerEta.lower_bound(left_eta - halfWayEta_); 00385 rr != trackInnerEta.upper_bound(left_eta + halfWayEta_); ++rr){//select neighbor tracks by eta 00386 //Level 2 loop 00387 //TODO find the closest one rather than the first matching 00388 const edm::Ref<reco::TrackCollection> & right = allTracks[rr->second]; 00389 00390 if (right->normalizedChi2() > maxChi2Right_ || right->found() < maxHitsRight_ //track quality cut 00391 || isPaired[rr->second] //this track should not be in another pair 00392 //|| right == (*ll) //no double counting (dummy if require opposite charge) 00393 || right->d0()*right->charge()<0. //d0*charge for converted photons 00394 || ((*ll)->charge())*(right->charge()) > 0) //opposite charge 00395 continue; 00396 00397 const double theta_l = (*ll)->innerMomentum().Theta(); 00398 const double theta_r = right->innerMomentum().Theta(); 00399 const double dCotTheta = 1./tan(theta_l) - 1./tan(theta_r) ; 00400 00401 if (fabs(dCotTheta) > deltaCotTheta_) continue;//Delta Cot(Theta) cut for pair 00402 00403 const double phi_l = (*ll)->innerMomentum().phi(); 00404 const double phi_r = right->innerMomentum().phi(); 00405 const double dPhi = phi_l - phi_r; 00406 00407 if (fabs(dPhi) > deltaPhi_) continue;//Delta Phi cut for pair 00408 00409 //TODO INSERT MORE CUTS HERE! 00410 00411 right_candidates.insert(std::pair<double, int>(theta_r, rr->second)); 00412 } 00413 //take the closest to left as right 00414 double min_theta = 999.; 00415 edm::Ref<reco::TrackCollection> right; 00416 int right_index = -1; 00417 for (std::map<double, int>::iterator rr = right_candidates.begin(); rr != right_candidates.end(); ++rr){ 00418 const double theta_l = (*ll)->innerMomentum().Theta(); 00419 const double dCotTheta = 1./tan(theta_l) - 1./tan(rr->first); 00420 if (fabs(min_theta) > fabs(dCotTheta)){ 00421 min_theta = dCotTheta; 00422 right_index = rr->second; 00423 right = allTracks[right_index]; 00424 } 00425 } 00426 if (min_theta <999.){//find a good right track 00427 //if all cuts passed, go ahead to make conversion candidates 00428 std::vector<edm::Ref<reco::TrackCollection> > trackPairRef; 00429 trackPairRef.push_back((*ll));//left track 00430 trackPairRef.push_back(right);//right track 00431 00432 reco::CaloClusterPtrVector scPtrVec; 00433 reco::Vertex theConversionVertex;//Dummy vertex, validity false by default 00434 std::vector<math::XYZPoint> trkPositionAtEcal; 00435 std::vector<reco::CaloClusterPtr> matchingBC; 00436 00437 trkPositionAtEcal.push_back(trackImpactPosition[ll-allTracks.begin()]);//left track 00438 if (trackValidECAL[right_index])//second track ECAL position may be invalid 00439 trkPositionAtEcal.push_back(trackImpactPosition[right_index]); 00440 00441 matchingBC.push_back(trackMatchedBC[ll-allTracks.begin()]);//left track 00442 scPtrVec.push_back(trackMatchedBC[ll-allTracks.begin()]);//left track 00443 if (trackValidECAL[right_index]){//second track ECAL position may be invalid 00444 matchingBC.push_back(trackMatchedBC[right_index]); 00445 scPtrVec.push_back(trackMatchedBC[right_index]); 00446 } 00447 00448 //TODO: currently, scPtrVec is assigned as matching BC; no Kalman vertex fit, so theConversionVertex validity is false by default 00449 // for first track (called left), trkPositionAtEcal and matchingBC must be valid 00450 // for second track (called right), trkPositionAtEcal and matchingBC is not necessary valid 00451 // so, BE CAREFUL check number of elements before reading them out 00452 reco::Conversion newCandidate(scPtrVec, trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC); 00453 outputConvPhotonCollection.push_back(newCandidate); 00454 00455 found_right = true; 00456 isPaired[ll-allTracks.begin()] = true;//mark this track is used in pair 00457 isPaired[right_index] = true; 00458 break; // to get another left leg and start new conversion 00459 } 00460 if (!found_right && allowSingleLeg_){ 00461 std::vector<edm::Ref<reco::TrackCollection> > trackPairRef; 00462 trackPairRef.push_back((*ll));//left track 00463 00464 reco::CaloClusterPtrVector scPtrVec; 00465 reco::Vertex theConversionVertex;//Dummy vertex, validity false by default 00466 std::vector<math::XYZPoint> trkPositionAtEcal; 00467 std::vector<reco::CaloClusterPtr> matchingBC; 00468 00469 trkPositionAtEcal.push_back(trackImpactPosition[ll-allTracks.begin()]);//left track 00470 00471 matchingBC.push_back(trackMatchedBC[ll-allTracks.begin()]);//left track 00472 scPtrVec.push_back(trackMatchedBC[ll-allTracks.begin()]);//left track 00473 00474 reco::Conversion newCandidate(scPtrVec, trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC); 00475 outputConvPhotonCollection.push_back(newCandidate); 00476 00477 isPaired[ll-allTracks.begin()] = true;//mark this track is used in pair 00478 } 00479 } 00480 00481 //output sorted as track pair then single track 00482 outputConvPhotonCollection_p->assign(outputConvPhotonCollection.begin(), outputConvPhotonCollection.end()); 00483 //outputConvPhotonTrackCollection_p->assign(outputConvPhotonTrackCollection.begin(), outputConvPhotonTrackCollection.end()); 00484 iEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_); 00485 //iEvent.put( outputConvPhotonTrackCollection_p, ConvertedPhotonCollection_); 00486 //TODO add photon object to reco::Photon 00487 00488 }
Definition at line 147 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
Definition at line 133 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
Definition at line 134 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
std::string TrackerOnlyConversionProducer::ConvertedPhotonCollection_ [private] |
Definition at line 135 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::deltaCotTheta_ [private] |
Definition at line 145 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::deltaPhi_ [private] |
Definition at line 145 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::dEtaTkBC_ [private] |
Definition at line 140 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::dPhiTkBC_ [private] |
Definition at line 140 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::energyBC_ [private] |
Definition at line 139 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::halfWayEta_ [private] |
Definition at line 137 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::halfWayPhi_ [private] |
Definition at line 137 of file TrackerOnlyConversionProducer.cc.
double TrackerOnlyConversionProducer::maxChi2Left_ [private] |
Definition at line 142 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::maxChi2Right_ [private] |
Definition at line 142 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::maxHitsLeft_ [private] |
Definition at line 143 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
double TrackerOnlyConversionProducer::maxHitsRight_ [private] |
Definition at line 143 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().
std::vector<edm::InputTag> TrackerOnlyConversionProducer::src_ [private] |
Definition at line 131 of file TrackerOnlyConversionProducer.cc.
Referenced by produce(), and TrackerOnlyConversionProducer().