CMS 3D CMS Logo

CastorTowerProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Castor
4 // Class: CastorTowerProducer
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
36 
37 // Castor Object include
41 
42 // Channel quality
46 
47 //
48 // class declaration
49 //
50 
52 public:
53  explicit CastorTowerProducer(const edm::ParameterSet&);
54  ~CastorTowerProducer() override;
55 
56 private:
57  void produce(edm::Event&, const edm::EventSetup&) override;
59  double& Ehot,
60  double& depth);
61 
62  // member data
64  typedef ROOT::Math::RhoEtaPhiPoint TowerPoint;
65  typedef ROOT::Math::RhoZPhiPoint CellPoint;
67  typedef std::vector<reco::CastorTower> CastorTowerCollection;
71  double towercut_;
72  double mintime_;
73  double maxtime_;
74 };
75 
76 //
77 // constants, enums and typedefs
78 //
79 
80 //
81 // static data member definitions
82 //
83 
84 //
85 // constructor and destructor
86 //
87 
89  : towercut_(iConfig.getParameter<double>("towercut")),
90  mintime_(iConfig.getParameter<double>("mintime")),
91  maxtime_(iConfig.getParameter<double>("maxtime")) {
92  tok_input_ = consumes<CastorRecHitCollection>(iConfig.getParameter<std::string>("inputprocess"));
93  tok_channelQuality_ = esConsumes<CastorChannelQuality, CastorChannelQualityRcd>();
94  //register your products
95  produces<CastorTowerCollection>();
96  //now do what ever other initialization is needed
97 }
98 
100  // do anything here that needs to be done at desctruction time
101  // (e.g. close files, deallocate resources etc.)
102 }
103 
104 //
105 // member functions
106 //
107 
108 // ------------ method called to produce the data ------------
110  using namespace edm;
111  using namespace reco;
112  using namespace TMath;
113 
114  // Produce CastorTowers from CastorCells
115 
117  iEvent.getByToken(tok_input_, InputRecHits);
118 
119  auto OutputTowers = std::make_unique<CastorTowerCollection>();
120 
121  // get and check input size
122  int nRecHits = InputRecHits->size();
123 
124  LogDebug("CastorTowerProducer") << "2. entering CastorTowerProducer" << std::endl;
125 
126  if (nRecHits == 0)
127  LogDebug("CastorTowerProducer") << "Warning: You are trying to run the Tower algorithm with 0 input rechits.";
128 
129  // declare castor array
130  // (0,x): Energies - (1,x): emEnergies - (2,x): hadEnergies - (3,x): phi position
131 
132  double poscastortowerarray[4][16];
133  double negcastortowerarray[4][16];
134 
135  CastorRecHitRefVector poscastorusedrechits[16];
136  CastorRecHitRefVector negcastorusedrechits[16];
137 
138  // set phi values and everything else to zero
139  for (int j = 0; j < 16; j++) {
140  poscastortowerarray[3][j] = -2.94524 + j * 0.3927;
141  poscastortowerarray[0][j] = 0.;
142  poscastortowerarray[1][j] = 0.;
143  poscastortowerarray[2][j] = 0.;
144 
145  negcastortowerarray[3][j] = -2.94524 + j * 0.3927;
146  negcastortowerarray[0][j] = 0.;
147  negcastortowerarray[1][j] = 0.;
148  negcastortowerarray[2][j] = 0.;
149  }
150 
151  // retrieve the channel quality lists from database
153  std::vector<DetId> channels = p->getAllChannels();
154 
155  // loop over rechits to build castortowerarray[4][16] and castorusedrechits[16]
156  for (unsigned int i = 0; i < InputRecHits->size(); i++) {
158 
159  HcalCastorDetId id = rechit_p->id();
160  DetId genericID = (DetId)id;
161 
162  // first check if the rechit is in the BAD channel list
163  bool bad = false;
164  for (std::vector<DetId>::iterator channel = channels.begin(); channel != channels.end(); channel++) {
165  if (channel->rawId() == genericID.rawId()) {
166  // if the rechit is found in the list, set it bad
167  bad = true;
168  break;
169  }
170  }
171  // if bad, continue the loop to the next rechit
172  if (bad)
173  continue;
174 
175  double Erechit = rechit_p->energy();
176  int module = id.module();
177  int sector = id.sector();
178  double zrechit = 0;
179  if (module < 3)
180  zrechit = -14390 - 24.75 - 49.5 * (module - 1);
181  if (module > 2)
182  zrechit = -14390 - 99 - 49.5 - 99 * (module - 3);
183  double phirechit = -100;
184  if (sector < 9)
185  phirechit = 0.19635 + (sector - 1) * 0.3927;
186  if (sector > 8)
187  phirechit = -2.94524 + (sector - 9) * 0.3927;
188 
189  // add time conditions for the rechit
190  if (rechit_p->time() > mintime_ && rechit_p->time() < maxtime_) {
191  // loop over the 16 towers possibilities
192  for (int j = 0; j < 16; j++) {
193  // phi matching condition
194  if (TMath::Abs(phirechit - poscastortowerarray[3][j]) < 0.0001) {
195  // condition over rechit z value
196  if (zrechit > 0.) {
197  poscastortowerarray[0][j] += Erechit;
198  if (module < 3) {
199  poscastortowerarray[1][j] += Erechit;
200  } else {
201  poscastortowerarray[2][j] += Erechit;
202  }
203  poscastorusedrechits[j].push_back(rechit_p);
204  } else {
205  negcastortowerarray[0][j] += Erechit;
206  if (module < 3) {
207  negcastortowerarray[1][j] += Erechit;
208  } else {
209  negcastortowerarray[2][j] += Erechit;
210  }
211  negcastorusedrechits[j].push_back(rechit_p);
212  } // end condition over rechit z value
213  } // end phi matching condition
214  } // end loop over the 16 towers possibilities
215  } // end time conditions
216 
217  } // end loop over rechits to build castortowerarray[4][16] and castorusedrechits[16]
218 
219  // make towers of the arrays
220 
221  double fem, Ehot, depth;
222  double rhoTower = 88.5;
223 
224  // loop over the 16 towers possibilities
225  for (int k = 0; k < 16; k++) {
226  Ehot = 0;
227  depth = 0;
228 
229  // select the positive towers with E > sqrt(Nusedrechits)*Ecut
230  if (poscastortowerarray[0][k] > sqrt(poscastorusedrechits[k].size()) * towercut_) {
231  fem = poscastortowerarray[1][k] / poscastortowerarray[0][k];
232  CastorRecHitRefVector usedRecHits = poscastorusedrechits[k];
233  ComputeTowerVariable(usedRecHits, Ehot, depth);
234 
235  LogDebug("CastorTowerProducer") << "tower " << k + 1 << ": fem = " << fem << " ,depth = " << depth
236  << " ,Ehot = " << Ehot << std::endl;
237 
238  TowerPoint temptowerposition(rhoTower, 5.9, poscastortowerarray[3][k]);
239  Point towerposition(temptowerposition);
240 
241  CastorTower newtower(poscastortowerarray[0][k],
242  towerposition,
243  poscastortowerarray[1][k],
244  poscastortowerarray[2][k],
245  fem,
246  depth,
247  Ehot,
248  poscastorusedrechits[k]);
249  OutputTowers->push_back(newtower);
250  } // end select the positive towers with E > Ecut
251 
252  // select the negative towers with E > sqrt(Nusedrechits)*Ecut
253  if (negcastortowerarray[0][k] > sqrt(negcastorusedrechits[k].size()) * towercut_) {
254  fem = negcastortowerarray[1][k] / negcastortowerarray[0][k];
255  CastorRecHitRefVector usedRecHits = negcastorusedrechits[k];
256  ComputeTowerVariable(usedRecHits, Ehot, depth);
257 
258  LogDebug("CastorTowerProducer") << "tower " << k + 1 << " energy = " << negcastortowerarray[0][k]
259  << "EM = " << negcastortowerarray[1][k] << "HAD = " << negcastortowerarray[2][k]
260  << "phi = " << negcastortowerarray[3][k] << ": fem = " << fem
261  << " ,depth = " << depth << " ,Ehot = " << Ehot << std::endl;
262 
263  TowerPoint temptowerposition(rhoTower, -5.9, negcastortowerarray[3][k]);
264  Point towerposition(temptowerposition);
265 
266  CastorTower newtower(negcastortowerarray[0][k],
267  towerposition,
268  negcastortowerarray[1][k],
269  negcastortowerarray[2][k],
270  fem,
271  depth,
272  Ehot,
273  negcastorusedrechits[k]);
274  OutputTowers->push_back(newtower);
275  } // end select the negative towers with E > Ecut
276 
277  } // end loop over the 16 towers possibilities
278 
279  iEvent.put(std::move(OutputTowers));
280 }
281 
283  double& Ehot,
284  double& depth) {
285  using namespace reco;
286 
287  double Etot = 0;
288 
289  // loop over the cells used in the tower k
290  for (CastorRecHitRefVector::iterator it = usedRecHits.begin(); it != usedRecHits.end(); it++) {
292 
293  double Erechit = rechit_p->energy();
294  HcalCastorDetId id = rechit_p->id();
295  int module = id.module();
296  double zrechit = 0;
297  if (module < 3)
298  zrechit = -14390 - 24.75 - 49.5 * (module - 1);
299  if (module > 2)
300  zrechit = -14390 - 99 - 49.5 - 99 * (module - 3);
301 
302  if (Erechit > Ehot)
303  Ehot = Erechit;
304  depth += Erechit * zrechit;
305  Etot += Erechit;
306  }
307 
308  depth /= Etot;
309  Ehot /= Etot;
310 }
311 
312 //define this as a plug-in
size
Write out results.
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::SortedCollection< CastorRecHit > CastorRecHitCollection
ROOT::Math::RhoZPhiPoint CellPoint
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
size_type size() const
std::vector< reco::CastorTower > CastorTowerCollection
ROOT::Math::RhoEtaPhiPoint TowerPoint
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:23
edm::RefVector< CastorRecHitCollection > CastorRecHitRefVector
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Definition: DetId.h:17
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
CastorTowerProducer(const edm::ParameterSet &)
edm::EDGetTokenT< CastorRecHitCollection > tok_input_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
virtual void ComputeTowerVariable(const edm::RefVector< edm::SortedCollection< CastorRecHit > > &usedRecHits, double &Ehot, double &depth)
fixed size matrix
HLT enums.
Structure Point Contains parameters of Gaussian fits to DMRs.
edm::ESGetToken< CastorChannelQuality, CastorChannelQualityRcd > tok_channelQuality_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)