CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 // system include
21 #include <memory>
22 #include <vector>
23 #include <iostream>
24 #include <TMath.h>
25 #include <TRandom3.h>
26 
27 // 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&);
55 
56  private:
57  virtual void beginJob() override ;
58  virtual void produce(edm::Event&, const edm::EventSetup&) override;
59  virtual void endJob() override ;
60  virtual void ComputeTowerVariable(const edm::RefVector<edm::SortedCollection<CastorRecHit> >& usedRecHits, double& Ehot, 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;
70  double towercut_;
71  double mintime_;
72  double maxtime_;
73 };
74 
75 //
76 // constants, enums and typedefs
77 //
78 
79 const double MYR2D = 180/M_PI;
80 
81 //
82 // static data member definitions
83 //
84 
85 //
86 // constructor and destructor
87 //
88 
90  towercut_(iConfig.getParameter<double>("towercut")),
91  mintime_(iConfig.getParameter<double>("mintime")),
92  maxtime_(iConfig.getParameter<double>("maxtime"))
93 {
94  tok_input_ = consumes<CastorRecHitCollection>(iConfig.getParameter<std::string>("inputprocess"));
95  //register your products
96  produces<CastorTowerCollection>();
97  //now do what ever other initialization is needed
98 }
99 
100 
102 {
103  // do anything here that needs to be done at desctruction time
104  // (e.g. close files, deallocate resources etc.)
105 }
106 
107 
108 //
109 // member functions
110 //
111 
112 // ------------ method called to produce the data ------------
114 
115  using namespace edm;
116  using namespace reco;
117  using namespace TMath;
118 
119  // Produce CastorTowers from CastorCells
120 
122  iEvent.getByToken(tok_input_,InputRecHits);
123 
124  std::auto_ptr<CastorTowerCollection> OutputTowers (new CastorTowerCollection);
125 
126  // get and check input size
127  int nRecHits = InputRecHits->size();
128 
129  LogDebug("CastorTowerProducer")
130  <<"2. entering CastorTowerProducer"<<std::endl;
131 
132  if (nRecHits==0)
133  LogDebug("CastorTowerProducer") <<"Warning: You are trying to run the Tower algorithm with 0 input rechits.";
134 
135  // declare castor array
136  // (0,x): Energies - (1,x): emEnergies - (2,x): hadEnergies - (3,x): phi position
137 
138  double poscastortowerarray[4][16];
139  double negcastortowerarray[4][16];
140 
141  CastorRecHitRefVector poscastorusedrechits[16];
142  CastorRecHitRefVector negcastorusedrechits[16];
143 
144  // set phi values and everything else to zero
145  for (int j = 0; j < 16; j++) {
146  poscastortowerarray[3][j] = -2.94524 + j*0.3927;
147  poscastortowerarray[0][j] = 0.;
148  poscastortowerarray[1][j] = 0.;
149  poscastortowerarray[2][j] = 0.;
150 
151  negcastortowerarray[3][j] = -2.94524 + j*0.3927;
152  negcastortowerarray[0][j] = 0.;
153  negcastortowerarray[1][j] = 0.;
154  negcastortowerarray[2][j] = 0.;
155  }
156 
157  // retrieve the channel quality lists from database
159  iSetup.get<CastorChannelQualityRcd>().get(p);
160  std::vector<DetId> channels = p->getAllChannels();
161 
162  // loop over rechits to build castortowerarray[4][16] and castorusedrechits[16]
163  for (unsigned int i = 0; i < InputRecHits->size(); i++) {
164 
166 
167  HcalCastorDetId id = rechit_p->id();
168  DetId genericID=(DetId)id;
169 
170  // first check if the rechit is in the BAD channel list
171  bool bad = false;
172  for (std::vector<DetId>::iterator channel = channels.begin();channel != channels.end();channel++) {
173  if (channel->rawId() == genericID.rawId()) {
174  // if the rechit is found in the list, set it bad
175  bad = true; break;
176  }
177  }
178  // if bad, continue the loop to the next rechit
179  if (bad) continue;
180 
181  double Erechit = rechit_p->energy();
182  int module = id.module();
183  int sector = id.sector();
184  double zrechit = 0;
185  if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);
186  if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);
187  double phirechit = -100;
188  if (sector < 9) phirechit = 0.19635 + (sector-1)*0.3927;
189  if (sector > 8) phirechit = -2.94524 + (sector - 9)*0.3927;
190 
191  // add time conditions for the rechit
192  if (rechit_p->time() > mintime_ && rechit_p->time() < maxtime_) {
193 
194  // loop over the 16 towers possibilities
195  for ( int j=0;j<16;j++) {
196 
197  // phi matching condition
198  if (TMath::Abs(phirechit - poscastortowerarray[3][j]) < 0.0001) {
199 
200  // condition over rechit z value
201  if (zrechit > 0.) {
202  poscastortowerarray[0][j]+=Erechit;
203  if (module < 3) {poscastortowerarray[1][j]+=Erechit;} else {poscastortowerarray[2][j]+=Erechit;}
204  poscastorusedrechits[j].push_back(rechit_p);
205  } else {
206  negcastortowerarray[0][j]+=Erechit;
207  if (module < 3) {negcastortowerarray[1][j]+=Erechit;} else {negcastortowerarray[2][j]+=Erechit;}
208  negcastorusedrechits[j].push_back(rechit_p);
209  } // end condition over rechit z value
210  } // end phi matching condition
211  } // end loop over the 16 towers possibilities
212  } // end time conditions
213 
214  } // end loop over rechits to build castortowerarray[4][16] and castorusedrechits[16]
215 
216  // make towers of the arrays
217 
218  double fem, Ehot, depth;
219  double rhoTower = 88.5;
220 
221  // loop over the 16 towers possibilities
222  for (int k=0;k<16;k++) {
223 
224  fem = 0;
225  Ehot = 0;
226  depth = 0;
227 
228  // select the positive towers with E > sqrt(Nusedrechits)*Ecut
229  if (poscastortowerarray[0][k] > sqrt(poscastorusedrechits[k].size())*towercut_) {
230 
231  fem = poscastortowerarray[1][k]/poscastortowerarray[0][k];
232  CastorRecHitRefVector usedRecHits = poscastorusedrechits[k];
233  ComputeTowerVariable(usedRecHits,Ehot,depth);
234 
235  LogDebug("CastorTowerProducer")
236  <<"tower "<<k+1<<": fem = "<<fem<<" ,depth = "<<depth<<" ,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],towerposition,poscastortowerarray[1][k],poscastortowerarray[2][k],fem,depth,Ehot,
242  poscastorusedrechits[k]);
243  OutputTowers->push_back(newtower);
244  } // end select the positive towers with E > Ecut
245 
246  // select the negative towers with E > sqrt(Nusedrechits)*Ecut
247  if (negcastortowerarray[0][k] > sqrt(negcastorusedrechits[k].size())*towercut_) {
248 
249  fem = negcastortowerarray[1][k]/negcastortowerarray[0][k];
250  CastorRecHitRefVector usedRecHits = negcastorusedrechits[k];
251  ComputeTowerVariable(usedRecHits,Ehot,depth);
252 
253  LogDebug("CastorTowerProducer")
254  <<"tower "<<k+1 << " energy = " << negcastortowerarray[0][k] << "EM = " << negcastortowerarray[1][k] << "HAD = " << negcastortowerarray[2][k] << "phi = " << negcastortowerarray[3][k] << ": fem = "<<fem<<" ,depth = "<<depth<<" ,Ehot = "<<Ehot<<std::endl;
255 
256  TowerPoint temptowerposition(rhoTower,-5.9,negcastortowerarray[3][k]);
257  Point towerposition(temptowerposition);
258 
259  CastorTower newtower(negcastortowerarray[0][k],towerposition,negcastortowerarray[1][k],negcastortowerarray[2][k],fem,depth,Ehot,
260  negcastorusedrechits[k]);
261  OutputTowers->push_back(newtower);
262  } // end select the negative towers with E > Ecut
263 
264  } // end loop over the 16 towers possibilities
265 
266  iEvent.put(OutputTowers);
267 }
268 
269 
270 // ------------ method called once each job just before starting event loop ------------
272  LogDebug("CastorTowerProducer")
273  <<"Starting CastorTowerProducer";
274 }
275 
276 // ------------ method called once each job just after ending the event loop ------------
278  LogDebug("CastorTowerProducer")
279  <<"Ending CastorTowerProducer";
280 }
281 
283 
284  using namespace reco;
285 
286  double Etot = 0;
287 
288  // loop over the cells used in the tower k
289  for (CastorRecHitRefVector::iterator it = usedRecHits.begin(); it != usedRecHits.end(); it++) {
290  edm::Ref<CastorRecHitCollection> rechit_p = *it;
291 
292  double Erechit = rechit_p->energy();
293  HcalCastorDetId id = rechit_p->id();
294  int module = id.module();
295  double zrechit = 0;
296  if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);
297  if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);
298 
299  if(Erechit > Ehot) Ehot = Erechit;
300  depth+=Erechit*zrechit;
301  Etot+=Erechit;
302  }
303 
304  depth/=Etot;
305  Ehot/=Etot;
306 }
307 
308 //define this as a plug-in
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
edm::SortedCollection< CastorRecHit > CastorRecHitCollection
ROOT::Math::RhoZPhiPoint CellPoint
module()
Definition: vlib.cc:994
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void endJob() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< reco::CastorTower > CastorTowerCollection
ROOT::Math::RhoEtaPhiPoint TowerPoint
int bad(Items const &cont)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
edm::RefVector< CastorRecHitCollection > CastorRecHitRefVector
int j
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &) override
int k[5][pyjets_maxn]
Definition: DetId.h:18
virtual void beginJob() override
#define M_PI
Definition: BFit3D.cc:3
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
const double MYR2D
CastorTowerProducer(const edm::ParameterSet &)
edm::EDGetTokenT< CastorRecHitCollection > tok_input_
virtual void ComputeTowerVariable(const edm::RefVector< edm::SortedCollection< CastorRecHit > > &usedRecHits, double &Ehot, double &depth)
const T & get() const
Definition: EventSetup.h:55
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
Definition: vlib.h:208
tuple size
Write out results.