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 
20 // system include
21 #include <memory>
22 #include <vector>
23 #include <iostream>
24 #include <TMath.h>
25 #include <TRandom3.h>
26 
27 // user include
34 
35 // Castor Object include
39 
40 //
41 // class declaration
42 //
43 
45  public:
46  explicit CastorCellProducer(const edm::ParameterSet&);
47  ~CastorCellProducer() override;
48 
49  private:
50  void beginJob() override ;
51  void produce(edm::Event&, const edm::EventSetup&) override;
52  void endJob() override ;
53 
54  // member data
56  typedef ROOT::Math::RhoZPhiPoint CellPoint;
57  typedef std::vector<reco::CastorCell> CastorCellCollection;
59 };
60 
61 //
62 // constants, enums and typedefs
63 //
64 
65 const double MYR2D = 180/M_PI;
66 
67 //
68 // static data member definitions
69 //
70 
71 //
72 // constructor and destructor
73 //
74 
76 {
77  tok_input_ = consumes<CastorRecHitCollection>(edm::InputTag(iConfig.getUntrackedParameter<std::string>("inputprocess","castorreco")));
78  // register your products
79  produces<CastorCellCollection>();
80  // now do what ever other initialization is needed
81 }
82 
83 
85 {
86  // do anything here that needs to be done at desctruction time
87  // (e.g. close files, deallocate resources etc.)
88 }
89 
90 
91 //
92 // member functions
93 //
94 
95 // ------------ method called to produce the data ------------
96 
98 
99  using namespace edm;
100  using namespace reco;
101  using namespace TMath;
102 
103  // Produce CastorCells from CastorRecHits
104 
106  iEvent.getByToken(tok_input_, InputRecHits);
107 
108  auto OutputCells = std::make_unique<CastorCellCollection>();
109 
110  // looping over all CastorRecHits
111 
112  LogDebug("CastorCellProducer")
113  <<"1. entering CastorCellProducer ";
114 
115  for (size_t i = 0; i < InputRecHits->size(); ++i) {
116  const CastorRecHit & rh = (*InputRecHits)[i];
117  int sector = rh.id().sector();
118  int module = rh.id().module();
119  double energy = rh.energy();
120  int zside = rh.id().zside();
121 
122  // define CastorCell properties
123  double zCell=0.;
124  double phiCell;
125  double rhoCell;
126 
127  // set z position of the cell
128  if (module < 3) {
129  // starting in EM section
130  if (module == 1) zCell = 14415;
131  if (module == 2) zCell = 14464;
132  } else {
133  // starting in HAD section
134  zCell = 14534 + (module - 3)*92;
135  }
136 
137  // set phi position of the cell
138  double castorphi[16];
139  for (int j = 0; j < 16; j++) {
140  castorphi[j] = -2.94524 + j*0.3927;
141  }
142  if (sector > 8) {
143  phiCell = castorphi[sector - 9];
144  } else {
145  phiCell = castorphi[sector + 7];
146  }
147 
148  // add condition to select in eta sides
149  if (zside <= 0) zCell = -1*zCell;
150 
151  // set rho position of the cell (inner radius 3.7cm, outer radius 14cm)
152  rhoCell = 88.5;
153 
154  // store cell position
155  CellPoint tempcellposition(rhoCell,zCell,phiCell);
156  Point cellposition(tempcellposition);
157 
158  LogDebug("CastorCellProducer")
159  <<"cell number: "<<i+1<<std::endl
160  <<"rho: "<<cellposition.rho()<<" phi: "<<cellposition.phi()*MYR2D<<" eta: "<<cellposition.eta()<<std::endl
161  <<"x: "<<cellposition.x()<<" y: "<<cellposition.y()<<" z: "<<cellposition.z();
162 
163  if (energy > 0.) {
164  CastorCell newCell(energy,cellposition);
165  OutputCells->push_back(newCell);
166  }
167 
168  } // end loop over CastorRecHits
169 
170  LogDebug("CastorCellProducer")
171  <<"total number of cells in the event: "<<InputRecHits->size();
172 
173  iEvent.put(std::move(OutputCells));
174 
175 
176 }
177 
178 // ------------ method called once each job just before starting event loop ------------
180  LogDebug("CastorCellProducer")
181  <<"Starting CastorCellProducer";
182 }
183 
184 // ------------ method called once each job just after ending the event loop ------------
186  LogDebug("CastorCellProducer")
187  <<"Ending CastorCellProducer";
188 }
189 
190 //define this as a plug-in
#define LogDebug(id)
constexpr float energy() const
Definition: CaloRecHit.h:31
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
int sector() const
get the sector (1-16)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void endJob() override
int module() const
get the module (1-2 for EM, 1-12 for HAD)
int zside(DetId const &)
HcalCastorDetId id() const
get the id
Definition: CastorRecHit.h:15
edm::EDGetTokenT< CastorRecHitCollection > tok_input_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int zside() const
get the z-side of the cell (1/-1)
std::vector< reco::CastorCell > CastorCellCollection
#define M_PI
~CastorCellProducer() override
math::XYZPointD Point
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
const double MYR2D
void beginJob() override
void produce(edm::Event &, const edm::EventSetup &) override
fixed size matrix
HLT enums.
size_type size() const
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
Definition: vlib.h:208
CastorCellProducer(const edm::ParameterSet &)
ROOT::Math::RhoZPhiPoint CellPoint
def move(src, dest)
Definition: eostools.py:511