CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFPSRecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFPSRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFPSRecHitCreator_h
3 
8 
12 
13 
19 
26 
28 
29  public:
31  PFRecHitCreatorBase(iConfig,iC)
32  {
34  }
35 
36  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
37 
38  beginEvent(iEvent,iSetup);
39 
42  iSetup.get<CaloGeometryRecord>().get(geoHandle);
43 
44  // get the ecal geometry
45  const CaloSubdetectorGeometry *psGeometry =
46  geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
47 
48  iEvent.getByToken(recHitToken_,recHitHandle);
49  for( const auto& erh : *recHitHandle ) {
50  ESDetId detid(erh.detid());
51  double energy = erh.energy();
52 
53 
55 
57 
58  switch( detid.plane() ) {
59  case 1:
60  layer = PFLayer::PS1;
61  break;
62  case 2:
63  layer = PFLayer::PS2;
64  break;
65  default:
66  throw cms::Exception("PFRecHitBadInput")
67  <<"incorrect preshower plane !! plane number "
68  <<detid.plane()<<std::endl;
69  }
70 
71 
72 
73  const CaloCellGeometry *thisCell;
74  thisCell= psGeometry->getGeometry(detid);
75 
76  // find rechit geometry
77  if(!thisCell) {
78  edm::LogError("PFPSRecHitCreator")
79  <<"warning detid "<<detid.rawId()
80  <<" not found in geometry"<<std::endl;
81  continue;
82  }
83 
84  auto const point = thisCell->getPosition();
85  position.SetCoordinates ( point.x(),
86  point.y(),
87  point.z() );
88 
89  reco::PFRecHit rh( detid.rawId(),layer,
90  energy,
91  position.x(), position.y(), position.z(),
92  0.0,0.0,0.0);
93 
94  rh.setDepth(detid.plane());
95  rh.setTime(erh.time());
96 
97  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
98  assert( corners.size() == 8 );
99 
100  rh.setNECorner( corners[0].x(), corners[0].y(), corners[0].z() );
101  rh.setSECorner( corners[1].x(), corners[1].y(), corners[1].z() );
102  rh.setSWCorner( corners[2].x(), corners[2].y(), corners[2].z() );
103  rh.setNWCorner( corners[3].x(), corners[3].y(), corners[3].z() );
104 
105 
106  bool rcleaned = false;
107  bool keep=true;
108 
109  //Apply Q tests
110  for( const auto& qtest : qualityTests_ ) {
111  if (!qtest->test(rh,erh,rcleaned)) {
112  keep = false;
113  }
114  }
115 
116  if(keep) {
117  out->push_back(rh);
118  }
119  else if (rcleaned)
120  cleaned->push_back(rh);
121  }
122  }
123 
124 
125 
126  protected:
128 
129 
130 };
131 
132 
133 #endif
void setSECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:100
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
void setDepth(int depth)
Definition: PFRecHit.h:79
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
void setNECorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:105
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
assert(m_qm.get())
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
void setSWCorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:95
const int keep
int iEvent
Definition: GenABIO.cc:230
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
void setTime(double time)
Definition: PFRecHit.h:78
void setNWCorner(double posx, double posy, double posz)
Definition: PFRecHit.cc:90
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
Layer
layer definition
Definition: PFLayer.h:31
size_type size() const
Definition: EZArrayFL.h:81
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:56
edm::EDGetTokenT< EcalRecHitCollection > recHitToken_
static int position[264][3]
Definition: ReadPGInfo.cc:509
PFPSRecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
const CornersVec & getCorners() const
Returns the corner points of this cell&#39;s volume.
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5