CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauProducer.cc
Go to the documentation of this file.
1 /* class PFRecoTauProducer
2  * EDProducer of the PFTauCollection, starting from the PFTauTagInfoCollection,
3  * authors: Simone Gennai (simone.gennai@cern.ch), Ludovic Houchu (Ludovic.Houchu@cern.ch)
4  */
5 
9 
19 
22 
25 
26 #include "CLHEP/Random/RandGauss.h"
27 
28 #include <memory>
29 
30 using namespace reco;
31 using namespace edm;
32 using namespace std;
33 
34 class PFRecoTauProducer : public EDProducer {
35  public:
36  explicit PFRecoTauProducer(const edm::ParameterSet& iConfig);
38  virtual void produce(edm::Event&,const edm::EventSetup&);
39  private:
43  std::string Algorithm_;
47  double JetMinPt_;
49 
50 };
51 
53  PFTauTagInfoProducer_ = iConfig.getParameter<edm::InputTag>("PFTauTagInfoProducer");
54  ElectronPreIDProducer_ = iConfig.getParameter<edm::InputTag>("ElectronPreIDProducer");
55  PVProducer_ = iConfig.getParameter<edm::InputTag>("PVProducer");
56  Algorithm_ = iConfig.getParameter<std::string>("Algorithm");
57  smearedPVsigmaX_ = iConfig.getParameter<double>("smearedPVsigmaX");
58  smearedPVsigmaY_ = iConfig.getParameter<double>("smearedPVsigmaY");
59  smearedPVsigmaZ_ = iConfig.getParameter<double>("smearedPVsigmaZ");
60  JetMinPt_ = iConfig.getParameter<double>("JetPtMin");
61 
62  if(Algorithm_ =="ConeBased") {
63  PFRecoTauAlgo_=new PFRecoTauAlgorithm(iConfig);
64  }
65  else if(Algorithm_ =="HPS") {
66  PFRecoTauAlgo_=new HPSPFRecoTauAlgorithm(iConfig);
67  }
68  else { //Add inside out Algorithm here
69 
70  //If no Algorithm found throw exception
71  throw cms::Exception("") << "Unknown Algorithkm" << std::endl;
72  }
73 
74 
75  produces<PFTauCollection>();
76 }
78  delete PFRecoTauAlgo_;
79 }
80 
82  auto_ptr<PFTauCollection> resultPFTau(new PFTauCollection);
83 
84  edm::ESHandle<TransientTrackBuilder> myTransientTrackBuilder;
85  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",myTransientTrackBuilder);
86  PFRecoTauAlgo_->setTransientTrackBuilder(myTransientTrackBuilder.product());
87 
88  //edm::ESHandle<MagneticField> myMF;
89  //iSetup.get<IdealMagneticFieldRecord>().get(myMF);
90  //PFRecoTauAlgo_->setMagneticField(myMF.product());
91 
92  // Electron PreID tracks: Temporary until integrated to PFCandidate
93  /*
94  edm::Handle<PFRecTrackCollection> myPFelecTk;
95  iEvent.getByLabel(ElectronPreIDProducer_,myPFelecTk);
96  const PFRecTrackCollection theElecTkCollection=*(myPFelecTk.product());
97  */
98  // query a rec/sim PV
100  iEvent.getByLabel(PVProducer_,thePVs);
101  const VertexCollection vertCollection=*(thePVs.product());
102  Vertex thePV;
103  if(vertCollection.size()) thePV=*(vertCollection.begin());
104  else{
105  Vertex::Error SimPVError;
106  SimPVError(0,0)=smearedPVsigmaX_*smearedPVsigmaX_;
107  SimPVError(1,1)=smearedPVsigmaY_*smearedPVsigmaY_;
108  SimPVError(2,2)=smearedPVsigmaZ_*smearedPVsigmaZ_;
109  Vertex::Point SimPVPoint(CLHEP::RandGauss::shoot(0.,smearedPVsigmaX_),
110  CLHEP::RandGauss::shoot(0.,smearedPVsigmaY_),
111  CLHEP::RandGauss::shoot(0.,smearedPVsigmaZ_));
112  thePV=Vertex(SimPVPoint,SimPVError,1,1,1);
113  }
114 
115  edm::Handle<PFTauTagInfoCollection> thePFTauTagInfoCollection;
116  iEvent.getByLabel(PFTauTagInfoProducer_,thePFTauTagInfoCollection);
117  int iinfo=0;
118  for(PFTauTagInfoCollection::const_iterator i_info=thePFTauTagInfoCollection->begin();i_info!=thePFTauTagInfoCollection->end();i_info++) {
119  if((*i_info).pfjetRef()->pt()>JetMinPt_){
120  // PFTau myPFTau=PFRecoTauAlgo_->buildPFTau(Ref<PFTauTagInfoCollection>(thePFTauTagInfoCollection,iinfo),thePV,theElecTkCollection);
121  PFTau myPFTau=PFRecoTauAlgo_->buildPFTau(Ref<PFTauTagInfoCollection>(thePFTauTagInfoCollection,iinfo),thePV);
122  resultPFTau->push_back(myPFTau);
123  }
124  ++iinfo;
125  }
126  iEvent.put(resultPFTau);
127 }
128 
edm::InputTag PFTauTagInfoProducer_
T getParameter(std::string const &) const
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
PFRecoTauProducer(const edm::ParameterSet &iConfig)
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
virtual void produce(edm::Event &, const edm::EventSetup &)
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:74
edm::InputTag ElectronPreIDProducer_
PFRecoTauAlgorithmBase * PFRecoTauAlgo_
edm::InputTag PVProducer_