CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EopTreeWriter Class Reference

#include <Alignment/OfflineValidation/plugins/EopTreeWriter.cc>

Inheritance diagram for EopTreeWriter:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 EopTreeWriter (const edm::ParameterSet &)
 
 ~EopTreeWriter () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
double getDistInCM (double eta1, double phi1, double eta2, double phi2)
 

Private Attributes

edm::Service< TFileServicefs_
 
TrackAssociatorParameters parameters_
 
edm::InputTag src_
 
TrackDetectorAssociator trackAssociator_
 
TTree * tree_
 
EopVariablestreeMemPtr_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 68 of file EopTreeWriter.cc.

Constructor & Destructor Documentation

EopTreeWriter::EopTreeWriter ( const edm::ParameterSet iConfig)
explicit

Definition at line 102 of file EopTreeWriter.cc.

References edm::EDConsumerBase::consumesCollector(), fs_, edm::ParameterSet::getParameter(), TrackAssociatorParameters::loadParameters(), TFileService::make(), parameters_, tree_, and treeMemPtr_.

102  :
103  src_(iConfig.getParameter<edm::InputTag>("src"))
104 {
105  //now do what ever initialization is needed
106 
107  // TrackAssociator parameters
108  edm::ParameterSet parameters = iConfig.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
110  parameters_.loadParameters( parameters, iC );
111 
112  tree_ = fs_->make<TTree>("EopTree","EopTree");
114  tree_->Branch("EopVariables", &treeMemPtr_); // address of pointer!
115 }
T getParameter(std::string const &) const
container to hold data to be written into TTree
Definition: EopVariables.h:8
void loadParameters(const edm::ParameterSet &, edm::ConsumesCollector &)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::InputTag src_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
TrackAssociatorParameters parameters_
edm::Service< TFileService > fs_
EopVariables * treeMemPtr_
EopTreeWriter::~EopTreeWriter ( )
override

Definition at line 118 of file EopTreeWriter.cc.

119 {
120 
121  // do anything here that needs to be done at destruction time
122  // (e.g. close files, deallocate resources etc.)
123 
124 }

Member Function Documentation

void EopTreeWriter::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 133 of file EopTreeWriter.cc.

References TrackDetectorAssociator::associate(), RecoEcal_EventContent_cff::ec, TrackDetMatchInfo::EcalRecHits, PV3DBase< T, PVType, FrameType >::eta(), Exception, EopVariables::fillVariables(), geometry, edm::EventSetup::get(), edm::Event::getByLabel(), getDistInCM(), TrackDetectorAssociator::getFreeTrajectoryState(), CaloGeometry::getPosition(), TrackDetMatchInfo::HcalRecHits, mps_fire::i, info(), TrackDetMatchInfo::nXnEnergy(), parameters_, PV3DBase< T, PVType, FrameType >::phi(), edm::ESHandle< T >::product(), rpcPointValidation_cfi::recHit, src_, HiIsolationCommonParameters_cff::track, trackAssociator_, l1t::tracks, tree_, treeMemPtr_, TrackDetMatchInfo::trkGlobPosAtEcal, TrackDetectorAssociator::useDefaultPropagator(), and TrackAssociatorParameters::useMuon.

134 {
135  using namespace edm;
136 
137  // get geometry
139  iSetup.get<CaloGeometryRecord>().get(geometry);
140  const CaloGeometry* geo = geometry.product();
141 // const CaloSubdetectorGeometry* towerGeometry =
142 // geo->getSubdetectorGeometry(DetId::Calo, CaloTowerDetId::SubdetId);
143 
144  // temporary collection of EB+EE recHits
145  std::auto_ptr<EcalRecHitCollection> tmpEcalRecHitCollection(new EcalRecHitCollection);
146  std::vector<edm::InputTag> ecalLabels_;
147 
149  bool ecalInAlca = iEvent.getByLabel(edm::InputTag("IsoProd","IsoTrackEcalRecHitCollection"),tmpEc);
150  bool ecalInReco = iEvent.getByLabel(edm::InputTag("ecalRecHit","EcalRecHitsEB"),tmpEc)&&
151  iEvent.getByLabel(edm::InputTag("ecalRecHit","EcalRecHitsEE"),tmpEc);
152  if(ecalInAlca)ecalLabels_.push_back(edm::InputTag("IsoProd","IsoTrackEcalRecHitCollection"));
153  else if(ecalInReco){
154  ecalLabels_.push_back(edm::InputTag("ecalRecHit","EcalRecHitsEB"));
155  ecalLabels_.push_back(edm::InputTag("ecalRecHit","EcalRecHitsEE"));
156  }
157  else throw cms::Exception("MissingProduct","can not find EcalRecHits");
158 
159  std::vector<edm::InputTag>::const_iterator i;
160  for (i=ecalLabels_.begin(); i!=ecalLabels_.end(); i++)
161  {
163  iEvent.getByLabel(*i,ec);
164  for(EcalRecHitCollection::const_iterator recHit = (*ec).begin(); recHit != (*ec).end(); ++recHit)
165  {
166  tmpEcalRecHitCollection->push_back(*recHit);
167  }
168  }
169 
171  iEvent.getByLabel(src_, tracks);
172 
175  bool pixelInAlca = iEvent.getByLabel(edm::InputTag("IsoProd","HcalIsolatedTrackCollection"),tmpPix);
176  if(pixelInAlca)iEvent.getByLabel(edm::InputTag("IsoProd","HcalIsolatedTrackCollection"),isoPixelTracks);
177 
178  Double_t trackemc1;
179  Double_t trackemc3;
180  Double_t trackemc5;
181  Double_t trackhac1;
182  Double_t trackhac3;
183  Double_t trackhac5;
184  Double_t maxPNearby;
185  Double_t dist;
186  Double_t EnergyIn;
187  Double_t EnergyOut;
188 
189  parameters_.useMuon = false;
190 
191  if(pixelInAlca)
192  if(isoPixelTracks->empty()) return;
193 
194  for(reco::TrackCollection::const_iterator track = tracks->begin();track!=tracks->end();++track){
195 
196  bool noChargedTracks = true;
197 
198  if(track->p()<9.) continue;
199 
202 
203  trackemc1 = info.nXnEnergy(TrackDetMatchInfo::EcalRecHits, 0);
204  trackemc3 = info.nXnEnergy(TrackDetMatchInfo::EcalRecHits, 1);
205  trackemc5 = info.nXnEnergy(TrackDetMatchInfo::EcalRecHits, 2);
206  trackhac1 = info.nXnEnergy(TrackDetMatchInfo::HcalRecHits, 0);
207  trackhac3 = info.nXnEnergy(TrackDetMatchInfo::HcalRecHits, 1);
208  trackhac5 = info.nXnEnergy(TrackDetMatchInfo::HcalRecHits, 2);
209 
210  if(trackhac3<5.) continue;
211 
212  double etaecal=info.trkGlobPosAtEcal.eta();
213  double phiecal=info.trkGlobPosAtEcal.phi();
214 
215  maxPNearby=-10;
216  dist=50;
217  for (reco::TrackCollection::const_iterator track1 = tracks->begin(); track1!=tracks->end(); track1++)
218  {
219  if (track == track1) continue;
220  TrackDetMatchInfo info1 = trackAssociator_.associate(iEvent, iSetup, *track1, parameters_);
221  double etaecal1=info1.trkGlobPosAtEcal.eta();
222  double phiecal1=info1.trkGlobPosAtEcal.phi();
223 
224  if (etaecal1==0&&phiecal1==0) continue;
225 
226  double ecDist=getDistInCM(etaecal,phiecal,etaecal1,phiecal1);
227 
228  if( ecDist < 40. )
229  {
230  //calculate maximum P and sum P near seed track
231  if (track1->p()>maxPNearby)
232  {
233  maxPNearby=track1->p();
234  dist = ecDist;
235  }
236 
237  //apply loose isolation criteria
238  if (track1->p()>5.)
239  {
240  noChargedTracks = false;
241  break;
242  }
243  }
244  }
245  EnergyIn=0;
246  EnergyOut=0;
247  if(noChargedTracks){
248  for (std::vector<EcalRecHit>::const_iterator ehit=tmpEcalRecHitCollection->begin(); ehit!=tmpEcalRecHitCollection->end(); ehit++)
249  {
251  // R-scheme of ECAL CLUSTERIZATION
252  const GlobalPoint& posH = geo->getPosition((*ehit).detid());
253  double phihit = posH.phi();
254  double etahit = posH.eta();
255 
256  double dHitCM=getDistInCM(etaecal,phiecal,etahit,phihit);
257 
258  if (dHitCM<9.0)
259  {
260  EnergyIn+=ehit->energy();
261  }
262  if (dHitCM>15.0&&dHitCM<35.0)
263  {
264  EnergyOut+=ehit->energy();
265  }
266 
267  }
268 
269  treeMemPtr_->fillVariables(track->charge(), track->innerOk(), track->outerRadius(),
270  track->numberOfValidHits(), track->numberOfLostHits(),
271  track->chi2(), track->normalizedChi2(),
272  track->p(), track->pt(), track->ptError(),
273  track->theta(), track->eta(), track->phi(),
274  trackemc1, trackemc3, trackemc5,
275  trackhac1, trackhac3, trackhac5,
276  maxPNearby, dist, EnergyIn, EnergyOut);
277 
278  tree_->Fill();
279  }
280  }
281 }
static const TGPicture * info(bool iBackgroundIsBlack)
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
void useDefaultPropagator()
use the default propagator
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
std::vector< EcalRecHit >::const_iterator const_iterator
edm::InputTag src_
double nXnEnergy(const DetId &, EnergyType, int gridSize=1)
get energy of the NxN shape (N = 2*gridSize + 1) around given detector element
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:74
TrackAssociatorParameters parameters_
void fillVariables(Int_t charge, Int_t innerOk, Double_t outerRadius, Int_t numberOfValidHits, Int_t numberOfLostHits, Double_t chi2, Double_t normalizedChi2, Double_t p, Double_t pt, Double_t ptError, Double_t theta, Double_t eta, Double_t phi, Double_t emc1, Double_t emc3, Double_t emc5, Double_t hac1, Double_t hac3, Double_t hac5, Double_t maxPNearby, Double_t dist, Double_t EnergyIn, Double_t EnergyOut)
fill variables into tree
Definition: EopVariables.h:25
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
TrackDetectorAssociator trackAssociator_
T eta() const
Definition: PV3DBase.h:76
ESHandle< TrackerGeometry > geometry
HLT enums.
double getDistInCM(double eta1, double phi1, double eta2, double phi2)
T get() const
Definition: EventSetup.h:63
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
T const * product() const
Definition: ESHandle.h:86
EopVariables * treeMemPtr_
void EopTreeWriter::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 286 of file EopTreeWriter.cc.

287 {
288 }
void EopTreeWriter::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 292 of file EopTreeWriter.cc.

References treeMemPtr_.

Referenced by o2olib.O2ORunMgr::executeJob().

292  {
293 
294  delete treeMemPtr_; treeMemPtr_ = nullptr;
295 
296 }
EopVariables * treeMemPtr_
double EopTreeWriter::getDistInCM ( double  eta1,
double  phi1,
double  eta2,
double  phi2 
)
private

Definition at line 299 of file EopTreeWriter.cc.

References funct::cos(), DEFINE_FWK_MODULE, hiPixelPairStep_cff::deltaPhi, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, JetChargeProducer_cfi::exp, Pi, mathSSE::sqrt(), and funct::tan().

Referenced by analyze().

300 {
301  double deltaPhi=phi1-phi2;
302  while(deltaPhi > TMath::Pi())deltaPhi-=2*TMath::Pi();
303  while(deltaPhi <= -TMath::Pi())deltaPhi+=2*TMath::Pi();
304  double dR;
305  // double Rec;
306  double theta1=2*atan(exp(-eta1));
307  double theta2=2*atan(exp(-eta2));
308  double cotantheta1;
309  if(cos(theta1)==0)cotantheta1=0;
310  else cotantheta1=1/tan(theta1);
311  double cotantheta2;
312  if(cos(theta2)==0)cotantheta2=0;
313  else cotantheta2=1/tan(theta2);
314  // if (fabs(eta1)<1.479) Rec=129; //radius of ECAL barrel
315  // else Rec=317; //distance from IP to ECAL endcap
316  //|vect| times tg of acos(scalar product)
317  // dR=fabs((Rec/sin(theta1))*tan(acos(sin(theta1)*sin(theta2)*(sin(phi1)*sin(phi2)+cos(phi1)*cos(phi2))+cos(theta1)*cos(theta2))));
318  if(fabs(eta1)<1.479)dR=129*sqrt((cotantheta1-cotantheta2)*(cotantheta1-cotantheta2)+deltaPhi*deltaPhi);
319  else dR=317*sqrt(tan(theta1)*tan(theta1)+tan(theta2)*tan(theta2)-2*tan(theta1)*tan(theta2)*cos(deltaPhi));
320  return dR;
321 }
const double Pi
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22

Member Data Documentation

edm::Service<TFileService> EopTreeWriter::fs_
private

Definition at line 84 of file EopTreeWriter.cc.

Referenced by EopTreeWriter().

TrackAssociatorParameters EopTreeWriter::parameters_
private
edm::InputTag EopTreeWriter::src_
private

Definition at line 82 of file EopTreeWriter.cc.

Referenced by analyze().

TrackDetectorAssociator EopTreeWriter::trackAssociator_
private

Definition at line 87 of file EopTreeWriter.cc.

Referenced by analyze().

TTree* EopTreeWriter::tree_
private

Definition at line 85 of file EopTreeWriter.cc.

Referenced by analyze(), and EopTreeWriter().

EopVariables* EopTreeWriter::treeMemPtr_
private

Definition at line 86 of file EopTreeWriter.cc.

Referenced by analyze(), endJob(), and EopTreeWriter().