CMS 3D CMS Logo

CastorCellProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Castor
4 // Class: CastorCellProducer
5 //
6 
13 //
14 // Original Author: Hans Van Haevermaet, Benoit Roland
15 // Created: Wed Jul 9 14:00:40 CEST 2008
16 //
17 //
18 
19 // system include
20 #include <memory>
21 #include <vector>
22 #include <iostream>
23 #include <TMath.h>
24 #include <TRandom3.h>
25 
26 // user include
33 
34 // Castor Object include
38 
39 //
40 // class declaration
41 //
42 
44 public:
45  explicit CastorCellProducer(const edm::ParameterSet&);
46 
47 private:
48  void beginJob() override;
49  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
50  void endJob() override;
51 
52  // member data
54  typedef ROOT::Math::RhoZPhiPoint CellPoint;
55  typedef std::vector<reco::CastorCell> CastorCellCollection;
57 };
58 
59 //
60 // constants, enums and typedefs
61 //
62 namespace {
63  constexpr double MYR2D = 180 / M_PI;
64 }
65 //
66 // static data member definitions
67 //
68 
69 //
70 // constructor and destructor
71 //
72 
74  tok_input_ = consumes<CastorRecHitCollection>(
75  edm::InputTag(iConfig.getUntrackedParameter<std::string>("inputprocess", "castorreco")));
76  // register your products
77  produces<CastorCellCollection>();
78  // now do what ever other initialization is needed
79 }
80 
81 //
82 // member functions
83 //
84 
85 // ------------ method called to produce the data ------------
86 
88  using namespace edm;
89  using namespace reco;
90  using namespace TMath;
91 
92  // Produce CastorCells from CastorRecHits
93 
95  iEvent.getByToken(tok_input_, InputRecHits);
96 
97  auto OutputCells = std::make_unique<CastorCellCollection>();
98 
99  // looping over all CastorRecHits
100 
101  LogDebug("CastorCellProducer") << "1. entering CastorCellProducer ";
102 
103  for (size_t i = 0; i < InputRecHits->size(); ++i) {
104  const CastorRecHit& rh = (*InputRecHits)[i];
105  int sector = rh.id().sector();
106  int module = rh.id().module();
107  double energy = rh.energy();
108  int zside = rh.id().zside();
109 
110  // define CastorCell properties
111  double zCell = 0.;
112  double phiCell;
113  double rhoCell;
114 
115  // set z position of the cell
116  if (module < 3) {
117  // starting in EM section
118  if (module == 1)
119  zCell = 14415;
120  if (module == 2)
121  zCell = 14464;
122  } else {
123  // starting in HAD section
124  zCell = 14534 + (module - 3) * 92;
125  }
126 
127  // set phi position of the cell
128  double castorphi[16];
129  for (int j = 0; j < 16; j++) {
130  castorphi[j] = -2.94524 + j * 0.3927;
131  }
132  if (sector > 8) {
133  phiCell = castorphi[sector - 9];
134  } else {
135  phiCell = castorphi[sector + 7];
136  }
137 
138  // add condition to select in eta sides
139  if (zside <= 0)
140  zCell = -1 * zCell;
141 
142  // set rho position of the cell (inner radius 3.7cm, outer radius 14cm)
143  rhoCell = 88.5;
144 
145  // store cell position
146  CellPoint tempcellposition(rhoCell, zCell, phiCell);
147  Point cellposition(tempcellposition);
148 
149  LogDebug("CastorCellProducer") << "cell number: " << i + 1 << std::endl
150  << "rho: " << cellposition.rho() << " phi: " << cellposition.phi() * MYR2D
151  << " eta: " << cellposition.eta() << std::endl
152  << "x: " << cellposition.x() << " y: " << cellposition.y()
153  << " z: " << cellposition.z();
154 
155  if (energy > 0.) {
156  CastorCell newCell(energy, cellposition);
157  OutputCells->push_back(newCell);
158  }
159 
160  } // end loop over CastorRecHits
161 
162  LogDebug("CastorCellProducer") << "total number of cells in the event: " << InputRecHits->size();
163 
164  iEvent.put(std::move(OutputCells));
165 }
166 
167 // ------------ method called once each job just before starting event loop ------------
168 void CastorCellProducer::beginJob() { LogDebug("CastorCellProducer") << "Starting CastorCellProducer"; }
169 
170 // ------------ method called once each job just after ending the event loop ------------
171 void CastorCellProducer::endJob() { LogDebug("CastorCellProducer") << "Ending CastorCellProducer"; }
172 
173 //define this as a plug-in
size_type size() const
void endJob() override
HcalCastorDetId id() const
get the id
Definition: CastorRecHit.h:14
int zside(DetId const &)
constexpr float energy() const
Definition: CaloRecHit.h:29
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< CastorRecHitCollection > tok_input_
int iEvent
Definition: GenABIO.cc:224
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int module() const
get the module (1-2 for EM, 1-12 for HAD)
int zside() const
get the z-side of the cell (1/-1)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int sector() const
get the sector (1-16)
std::vector< reco::CastorCell > CastorCellCollection
#define M_PI
math::XYZPointD Point
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
void beginJob() override
fixed size matrix
HLT enums.
Structure Point Contains parameters of Gaussian fits to DMRs.
CastorCellProducer(const edm::ParameterSet &)
ROOT::Math::RhoZPhiPoint CellPoint
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)