CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTRegionalPixelSeedGeneratorProducers.cc
Go to the documentation of this file.
1 //
2 // Package: RecoEgamma/EgammaHLTProducers
3 // Class: EgammaHLTRegionalPixelSeedGeneratorProducers
4 // Modified from TkSeedGeneratorFromTrk by Jeremy Werner, Princeton University, USA
5 // $Id: EgammaHLTRegionalPixelSeedGeneratorProducers.cc,v 1.13 2012/01/23 12:56:38 sharper Exp $
6 //
7 
8 #include <iostream>
9 #include <memory>
10 #include <string>
11 
26 
31 
34 // Math
35 #include "Math/GenVector/VectorUtil.h"
36 #include "Math/GenVector/PxPyPzE4D.h"
39 
40 using namespace std;
41 using namespace reco;
42 
44 {
45 
46  produces<TrajectorySeedCollection>();
47 
48  ptmin_ = conf.getParameter<double>("ptMin");
49  vertexz_ = conf.getParameter<double>("vertexZ");
50  originradius_= conf.getParameter<double>("originRadius");
51  halflength_ = conf.getParameter<double>("originHalfLength");
52  deltaEta_ = conf.getParameter<double>("deltaEtaRegion");
53  deltaPhi_ = conf.getParameter<double>("deltaPhiRegion");
54 
55  candTag_ = consumes<reco::RecoEcalCandidateCollection>(conf.getParameter< edm::InputTag > ("candTag"));
56  candTagEle_ = consumes<reco::ElectronCollection>(conf.getParameter< edm::InputTag > ("candTagEle"));
57  BSProducer_ = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("BSProducer"));
58 
59  useZvertex_ = conf.getParameter<bool>("UseZInVertex");
60 
61  edm::ParameterSet hitsfactoryPSet = conf.getParameter<edm::ParameterSet>("OrderedHitsFactoryPSet");
62  std::string hitsfactoryName = hitsfactoryPSet.getParameter<std::string>("ComponentName");
63 
64  // get orderd hits generator from factory
65  edm::ConsumesCollector iC = consumesCollector();
66  OrderedHitsGenerator* hitsGenerator = OrderedHitsGeneratorFactory::get()->create( hitsfactoryName, hitsfactoryPSet, iC);
67 
68  // start seed generator
69  edm::ParameterSet creatorPSet;
70  creatorPSet.addParameter<std::string>("propagator","PropagatorWithMaterial");
71 
72  combinatorialSeedGenerator.reset(new SeedGeneratorFromRegionHits( hitsGenerator, 0,
73  SeedCreatorFactory::get()->create("SeedFromConsecutiveHitsCreator", creatorPSet)
74  ));
75  // setup orderedhits setup (in order to tell seed generator to use pairs/triplets, which layers)
76 }
77 
78 // Virtual destructor needed.
80 }
81 
83 
85  desc.add<double>("ptMin", 1.5);
86  desc.add<double>("vertexZ", 0);
87  desc.add<double>("originRadius", 0.02);
88  desc.add<double>("originHalfLength", 15.0);
89  desc.add<double>("deltaEtaRegion", 0.3);
90  desc.add<double>("deltaPhiRegion", 0.3);
91  desc.add<edm::InputTag>(("candTag"), edm::InputTag("hltL1SeededRecoEcalCandidate"));
92  desc.add<edm::InputTag>(("candTagEle"), edm::InputTag("pixelMatchElectrons"));
93  desc.add<edm::InputTag>(("BSProducer"), edm::InputTag("hltOnlineBeamSpot"));
94  desc.add<bool>(("UseZInVertex"), false);
95  desc.add<std::string>("TTRHBuilder", "WithTrackAngle");
96 
97  edm::ParameterSetDescription orderedHitsPSET;
98  orderedHitsPSET.add<std::string>("ComponentName", "StandardHitPairGenerator");
99  orderedHitsPSET.add<edm::InputTag>("SeedingLayers", edm::InputTag("PixelLayerPairs"));
100  orderedHitsPSET.add<unsigned int>("maxElement", 0);
101  desc.add<edm::ParameterSetDescription>("OrderedHitsFactoryPSet", orderedHitsPSET);
102 
103  descriptions.add(("hltEgammaHLTRegionalPixelSeedGeneratorProducers"), desc);
104 }
105 
107 
108 
110 {
111 }
112 
113 // Functions that gets called by framework every event
115 {
116 
117  // resulting collection
118  std::auto_ptr<TrajectorySeedCollection> output(new TrajectorySeedCollection());
119 
120  // Get the recoEcalCandidates
122  iEvent.getByToken(candTag_,recoecalcands);
123 
124  //Get the Beam Spot position
125  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
126  iEvent.getByToken(BSProducer_,recoBeamSpotHandle);
127  // gets its position
128  const BeamSpot::Point& BSPosition = recoBeamSpotHandle->position();
129 
130  //Get the HLT electrons collection if needed
132  if(useZvertex_)
133  iEvent.getByToken(candTagEle_,electronHandle);
134 
135  reco::SuperClusterRef scRef;
136  for (reco::RecoEcalCandidateCollection::const_iterator recoecalcand= recoecalcands->begin(); recoecalcand!=recoecalcands->end(); recoecalcand++) {
137  scRef = recoecalcand->superCluster();
138  float zvertex = 0;
139  if( useZvertex_ ){
140  reco::SuperClusterRef scRefEle;
141  for(reco::ElectronCollection::const_iterator iElectron = electronHandle->begin(); iElectron != electronHandle->end(); iElectron++){
142  //Compare electron SC with EcalCandidate SC
143  scRefEle = iElectron->superCluster();
144  if(&(*scRef) == &(*scRefEle)){
145  if(iElectron->track().isNonnull()) zvertex = iElectron->track()->vz();
146  else zvertex = iElectron->gsfTrack()->vz();
147  break;
148  }
149  }
150 
151  }
152  GlobalVector dirVector((recoecalcand)->px(),(recoecalcand)->py(),(recoecalcand)->pz());
153  RectangularEtaPhiTrackingRegion etaphiRegion( dirVector,
154  GlobalPoint( BSPosition.x(), BSPosition.y(), zvertex ),
155  ptmin_,
156  originradius_,
157  halflength_,
158  deltaEta_,
159  deltaPhi_);
160 
161  // fill Trajectory seed collection
162  combinatorialSeedGenerator->run(*output, etaphiRegion, iEvent, iSetup);
163 
164  }
165 
166  iEvent.put(output);
167 }
T getParameter(std::string const &) const
virtual void produce(edm::Event &e, const edm::EventSetup &c)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void endRun(edm::Run const &run, const edm::EventSetup &es) overridefinal
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
int iEvent
Definition: GenABIO.cc:230
std::vector< TrajectorySeed > TrajectorySeedCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:142
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tuple conf
Definition: dbtoconf.py:185
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) overridefinal
void add(std::string const &label, ParameterSetDescription const &psetDescription)
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:55
Definition: Run.h:41