StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEmcGeom.h
1 /***************************************************************************
2  *
3  * $Id: StEmcGeom.h,v 1.5 2007/07/16 21:24:59 kocolosk Exp $
4  *
5  * Author: Aleksei Pavlinov
6  ***************************************************************************
7  *
8  * Description:
9  *
10  ***************************************************************************
11  *
12  * $Log: StEmcGeom.h,v $
13  * Revision 1.5 2007/07/16 21:24:59 kocolosk
14  * added projection against sub == -1 case in getId(phi,eta,softId). Thanks Pibero
15  *
16  * Revision 1.4 2007/04/04 17:32:11 kocolosk
17  * Added softId-based versions of getEta, getTheta, and getPhi. Also added getId(phi,eta,&softId). Implemented const-correctness and used meaningful argument names in method declarations to improve readability
18  *
19  * Revision 1.3 2004/08/19 17:31:45 pavlinov
20  * getBin(const Float_t phi, const Float_t eta, Int_t &m,Int_t &e,Int_t &s) works bsmde too - request of Dmitry Arkhipkin
21  *
22  * Revision 1.2 2003/09/02 17:58:01 perev
23  * gcc 3.2 updates + WarnOff
24  *
25  * Revision 1.1 2003/01/23 01:30:28 suaide
26  * moving to sub directories
27  *
28  * Revision 1.13 2001/11/21 18:56:19 suaide
29  * added methods to get eta and phi bins boundaries
30  *
31  * Revision 1.12 2001/09/24 17:34:51 akio
32  * Dirty bug fix.
33  *
34  * Revision 1.11 2001/09/22 00:28:58 pavlinov
35  * No public constructor for StEmcGeom
36  *
37  * Revision 1.10 2001/08/08 00:33:15 pavlinov
38  * New Jose's ID
39  *
40  * Revision 1.9 2001/07/30 00:16:22 pavlinov
41  * Correct numbering scheme for BSMDE
42  *
43  * Revision 1.8 2001/07/13 18:35:20 pavlinov
44  * New version of methods getBin and getId for BSMDE and BSMDP
45  *
46  * Revision 1.7 2001/04/26 14:23:40 akio
47  * Quick and dirty fix for crashing non-bfc chain
48  *
49  * Revision 1.6 2001/04/25 01:03:10 pavlinov
50  * Clean up
51  *
52  * Revision 1.5 2001/03/22 21:50:35 pavlinov
53  * Clean up for mdc4
54  *
55  * Revision 1.4 2001/03/15 23:47:20 pavlinov
56  * Fixed error on SUN compiler
57  *
58  * Revision 1.3 2001/03/15 20:56:11 pavlinov
59  * Jose's scheme is default
60  *
61  * Revision 1.2 2001/03/15 00:57:28 pavlinov
62  * Created new methods getBinEnv and getIdEnv
63  *
64  * Revision 1.1 2000/06/20 17:11:25 pavlinov
65  * Move StEmcGeom.h from St_emc_Maker to StEmcUtil
66  *
67  * Revision 1.11 2000/05/23 14:35:01 pavlinov
68  * Clean up for SUN
69  *
70  * Revision 1.10 2000/05/18 17:07:30 pavlinov
71  * Fixed error for methods getXYZ(...)
72  *
73  * Revision 1.9 2000/04/27 01:39:15 pavlinov
74  * Cleanup for SUN
75  *
76  * Revision 1.8 2000/04/25 17:02:06 pavlinov
77  * Added methods for gettinng x,y,z from volume ID
78  *
79  * Revision 1.7 2000/04/21 17:43:01 pavlinov
80  * Added methods for for decoding Geant volume Id
81  *
82  * Revision 1.6 2000/04/18 20:38:09 pavlinov
83  * Added ctor from Geant geometry
84  *
85  * Revision 1.5 2000/04/11 19:48:40 pavlinov
86  * Merge versions of Pavlinov and Ogawa
87  *
88  * Revision 1.4 2000/01/29 00:04:57 akio
89  * temprary fix for endcap. need more work, but no more junk messages and crash
90  *
91  * Revision 1.3 1999/07/16 18:03:45 pavlinov
92  * Little correction for StEclMake
93  *
94  * Revision 1.2 1999/07/02 03:01:56 pavlinov
95  * Little corrections for Linux
96  *
97  * Revision 1.1 1999/07/01 16:17:57 pavlinov
98  * class StEmcGeom was created and maker was remade for new maker scheme
99  *
100  **************************************************************************/
101 #ifndef STAR_StEmcGeom
102 #define STAR_StEmcGeom
103 
105 // //
106 // StEmcGeom main class for <FONT COLOR="RED"> Geometry of BEMC for OFFline </FONT> //
107 // //
109 #include <Stiostream.h>
110 #include "math_constants.h"
111 #include <math.h>
112 #include <TArrayF.h>
113 #include <TString.h>
114 #include "tables/St_calb_calg_Table.h"
115 #include "tables/St_calb_calr_Table.h"
116 #include "StMessMgr.h"
117 
118 class StMaker;
119 class TDataSet;
120 
121 class StEmcGeom {
122 private:
123  static StEmcGeom *mGeom[8];
124 
125  StMaker* mChain;
126  TDataSet* mGeantGeom;
127  St_calb_calg* mCalg;
128  calb_calg_st* mCalg_st;
129  St_calb_calr* mCalr;
130  calb_calr_st* mCalr_st;
131 
132  void defineDefaultCommonConstants();
133  void defineCommonConstants();
134  void defineModuleGridOnPhi();
135  Float_t relativeError(Float_t, Float_t) const;
136  // static Int_t getIndex(const Float_t x, TArrayF &arr);
137  Int_t getIndex(const Float_t &x, const TArrayF &arr) const;
138 
139 
140 protected:
141  TString mMode; // Empty, "geant" or "db"
142  Int_t mDetector; // Detectors number => see emc/inc/emc_def.h
143  Int_t mNModule; // Number of modeles (120 is default)
144  Int_t mNEta; // Number of eta bins
145  Int_t mNSub; // Subdivision in phi bin
146  Int_t mNes; // mNes = mNEta * mNSub
147  Int_t mNRaw; // mNraw = mNes * mNModule
148  Float_t mRadius; // Distance to beam axis
149  Float_t mYWidth; // Half Width of module for mRadius
150  Float_t mEtaMax;
151  Float_t mEtaMin; //
152  TArrayF mPhiOffset; // Phi offset from numbering scheme for center of tower
153  TArrayF mPhiStep; // Step for transition to Star Cordinate System(-/+2pi/60)
154  TArrayF mPhiBound; // Phi offset from numbering scheme for edge of tower
155  Float_t mPhiStepHalf; // 2pi/(60*2)
156 
157  TArrayF mZlocal; // Array of z coordinates (system of single module)
158  TArrayF mYlocal; // Array of y coordinates (system of single module)
159  TArrayF mEta; // Array of eta coordinates of center modules
160  TArrayF mEtaB; // Array of eta boundaries
161  TArrayF mPhi; // Array of phi coordinates (system of single module)
162  TArrayF mPhiModule; // Angle of center of module in STAR system
163  TArrayF mYB; // 25-jul-2001
164  TArrayF mPhiB; //
165 
166  Int_t mMaxAdc; // Range of ADC for each detector
167 
168 public:
169  StEmcGeom(const Int_t );
170  StEmcGeom(const Char_t*);
171  StEmcGeom(const Int_t ,const Char_t*);
172 
173  static StEmcGeom *instance(const Int_t det);
174  static StEmcGeom *getEmcGeom(const Int_t det);
175  static StEmcGeom *instance(const Char_t* cdet);
176  static StEmcGeom *getEmcGeom(const Char_t* cdet);
177  static StEmcGeom *instance(const Int_t det, const Char_t* mode);
178  static StEmcGeom *getEmcGeom(const Int_t det, const Char_t* mode);
179 
180  static Int_t getDetNumFromName(const Char_t *cdet);
181 
182  virtual ~StEmcGeom();
183 
184  const TString* Mode() const;
185  Int_t Detector() const;
186  Int_t NModule() const;
187  Int_t NEta() const;
188  Int_t NSub() const;
189  Int_t Nes() const;
190  Int_t NRaw() const;
191  Float_t Radius() const;
192  Float_t YWidth() const;
193  Float_t EtaMax() const;
194  Float_t EtaMin() const;
195  const Float_t* PhiModule() const;
196  const Float_t* PhiOffset() const;
197  const Float_t* PhiStep() const;
198  const Float_t* PhiBound() const;
199  const Float_t* Zlocal() const;
200  const Float_t* Eta() const;
201  const Float_t* Ylocal() const;
202  const Float_t* Phi() const;
203 
204  const Float_t* EtaB() const; // Eta boundaries AASUAIDE
205  const Float_t* PhiB() const; // Phi boundaries AASUAIDE
206 
207  void setDetector(const Int_t val);
208  void setRadius(const Float_t val);
209  void setYWidth(const Float_t val);
210 
211  //check that the supplied value is in an acceptable range
212  Int_t checkModule(const Int_t m) const;
213  Int_t checkEta(const Int_t e) const;
214  Int_t checkSub(const Int_t s) const;
215  Int_t checkId(const Int_t softId) const;
216 
217  //convert (eta,phi) to (m,e,s) or softId
218  Int_t getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const;
219  Int_t getId(const Float_t phi, const Float_t eta, Int_t &softId) const;
220 
221  //convert (m,e,s) <-> softId
222  Int_t getBin(const Int_t softId, Int_t &m, Int_t &e, Int_t &s) const;
223  Int_t getId(const Int_t m, const Int_t e, const Int_t s, Int_t &softId) const;
224 
225  Int_t getVolIdBemc(const Int_t ivid, Int_t &module, Int_t &eta, Int_t &sub, Int_t &detector);
226  Int_t getVolIdBsmd(const Int_t ivid, Int_t &module, Int_t &eta, Int_t &sub, Int_t &detector);
227  Int_t getVolId(const Int_t ivid, Int_t &module, Int_t &eta, Int_t &sub, Int_t &det);
228 
229  Int_t getZlYl(const Int_t softId, Float_t &zl, Float_t &yl) const;
230  void getXYZ(const Int_t m, const Int_t e, const Int_t s, Float_t &x, Float_t &y, Float_t &z) const;
231  Int_t getXYZ(const Int_t softId, Float_t &x, Float_t &y, Float_t &z) const;
232  Int_t getXYZfromGeant(const Int_t ivid, Float_t &x, Float_t &y, Float_t &z);
233 
234  Int_t getEta(const Int_t m, const Int_t e, Float_t &eta) const;
235  Int_t getEta(const Int_t softId, Float_t &eta) const;
236 
237  Int_t getTheta(const Int_t m, const Int_t e, Float_t &theta) const;
238  Int_t getTheta(const Int_t softId, Float_t &theta) const;
239 
240  Int_t getPhi(const Int_t m, const Int_t s, Float_t &phi) const;
241  Int_t getPhi(const Int_t softId, Float_t &phi) const;
242 
243  Int_t getEtaPhi(const Int_t softId, Float_t &eta, Float_t &phi) const;
244 
245  Int_t getPhiModule(const Int_t m, Float_t &phi) const;
246 
247  Int_t getMaxAdc() const {return mMaxAdc;}
248 
249  void initGeom(const Int_t);
250  void initBEMCorBPRS();
251  void initBSMDE();
252  void initBSMDP();
253  void initEEMCorEPRS();
254  void initESMDE();
255  void initESMDP();
256  void printGeom() const;
257  void print() const {printGeom();}
258  void compare(const StEmcGeom &, Bool_t) const;
259  void compare(const StEmcGeom * const g, Bool_t key) const {compare(*g,key);};
260  void printError(Float_t) const;
261 
262  Float_t toDeg(const Float_t angR) const {return C_DEG_PER_RAD*angR;} // Service
263  Float_t toRad(const Float_t angD) const {return angD/C_DEG_PER_RAD;} // functions
264  void getGeantGeometryTable();
265 
266  ClassDef(StEmcGeom,1) // Standard Root macro;
267 };
268 
269 inline const TString* StEmcGeom::Mode() const {return &mMode;}
270 inline Int_t StEmcGeom::Detector() const {return mDetector;}
271 inline Int_t StEmcGeom::NModule() const {return mNModule;}
272 inline Int_t StEmcGeom::NEta() const {return mNEta;}
273 inline Int_t StEmcGeom::NSub() const {return mNSub;}
274 inline Int_t StEmcGeom::Nes() const {return mNes;}
275 inline Int_t StEmcGeom::NRaw() const {return mNRaw;}
276 inline Float_t StEmcGeom::Radius() const {return mRadius;}
277 inline Float_t StEmcGeom::YWidth() const {return mYWidth;}
278 inline Float_t StEmcGeom::EtaMax() const {return mEtaMax;}
279 inline Float_t StEmcGeom::EtaMin() const {return mEtaMin;}
280 inline const Float_t* StEmcGeom::PhiModule() const {return mPhiModule.GetArray();}
281 inline const Float_t* StEmcGeom::PhiOffset() const {return mPhiOffset.GetArray();}
282 inline const Float_t* StEmcGeom::PhiStep() const {return mPhiStep.GetArray();}
283 inline const Float_t* StEmcGeom::PhiBound() const {return mPhiBound.GetArray();}
284 inline const Float_t* StEmcGeom::Zlocal() const {return mZlocal.GetArray();}
285 inline const Float_t* StEmcGeom::Eta() const {return mEta.GetArray();}
286 inline const Float_t* StEmcGeom::Ylocal() const {return mYlocal.GetArray();}
287 inline const Float_t* StEmcGeom::Phi() const {return mPhi.GetArray();}
288 
289 inline const Float_t* StEmcGeom::EtaB() const {return mEtaB.GetArray();}
290 inline const Float_t* StEmcGeom::PhiB() const {return mPhiB.GetArray();}
291 
292 inline void StEmcGeom::setDetector(const Int_t val) { mDetector = val;}
293 inline void StEmcGeom::setRadius(const Float_t val) { mRadius = val;}
294 inline void StEmcGeom::setYWidth(const Float_t val) { mYWidth = val;}
295 
296 // _____________________________________________________________________
297 inline Int_t StEmcGeom::checkModule(const Int_t m) const
298 {
299  if(m>=1 && m<=mNModule) return 0;
300  else {LOG_ERROR<<" Bad module# "<<m<<"/"<<mNModule<<" in Detector "<<mDetector<<endm; return 1;}
301 }
302 // _____________________________________________________________________
303 inline Int_t StEmcGeom::checkEta(const Int_t e) const
304 {
305  if(e>=1 && e<=mNEta) return 0;
306  else {LOG_ERROR<<" Bad eta# "<<e<<endm; return 1;}
307 }
308 // _____________________________________________________________________
309 inline Int_t StEmcGeom::checkSub(const Int_t s) const
310 {
311  if(s>=1 && s<=mNSub) return 0;
312  else {LOG_ERROR<<" Bad sub# "<<s<<endm; return 1;}
313 }
314 // _____________________________________________________________________
315 inline Int_t StEmcGeom::checkId(const Int_t softId) const
316 {
317  if(softId>=1 && softId<=mNRaw) return 0;
318  else {LOG_ERROR<<" Bad raw# "<<softId<<endm; return 1;}
319 }
320 // _____________________________________________________________________
321 inline Int_t StEmcGeom::getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const
322 {
323 //
325 // 19-aug-2004 - for all detectors; request of Dmitry Arkhipkin
326 // mEtaMin=0 for BEMC, BPRS, BSMDP; mEtaMin= for BMSDE
327 //
328 
329  Float_t phiw, sw;
330  if(mEtaMin<eta && eta<=mEtaMax) { // First Barrel
331  e = getIndex(eta, mEtaB);
332  phiw = mPhiBound[0] - phi;
333  if(phiw<0.0) phiw = phiw + C_2PI; // 0<phiw<=2.*pi =>must be
334  if(phiw<0.0 || phiw>C_2PI) printf(" phi %f eta %f \n",phi,eta); // For testing
335  m = int(-phiw/mPhiStep[0]) + 1;
336 
337  sw = fmod(phiw, fabs(mPhiStep[0]));
338  sw = sw - mPhiStepHalf;
339  s = getIndex(sw, mPhiB);
340 
341  return 0;
342 
343  } else if(-mEtaMax<=eta && eta<-mEtaMin) { // Second Barrel
344  e = getIndex(fabs(eta), mEtaB);
345  phiw = mPhiBound[1] - phi;
346  if(phiw<0.0) phiw = phiw + C_2PI; // 0<phiw<=2.pi =>must be
347  if(phiw<0.0 || phiw>C_2PI) printf(" phi %f eta %f \n",phi,eta); // For testing
348  m = 120 - int(phiw/fabs(mPhiStep[1]));
349 
350  sw = fmod(phiw, fabs(mPhiStep[0]));
351  sw = -(sw - mPhiStepHalf);
352  s = getIndex(sw, mPhiB);
353 
354  return 0;
355 
356  } else return 1; // Out of Bemc
357 }
358 // _____________________________________________________________________
359 inline Int_t StEmcGeom::getId(const Float_t phi, const Float_t eta, Int_t &softId) const
360 {
361  Int_t m,e,s;
362  if(getBin(phi,eta,m,e,s) == 0 && s != -1) {
363  return getId(m,e,s,softId);
364  }
365  return 1;
366 }
367 // _____________________________________________________________________
368 inline Int_t StEmcGeom::getZlYl(const Int_t softId, Float_t &zl, Float_t &yl) const
369 {
370  Int_t m, e, s;
371  if(!getBin(softId, m, e, s)){
372  zl = mZlocal[e-1];
373  yl = mYlocal[s-1];
374  return 0;
375  }
376  else return 1;
377 }
378 // _____________________________________________________________________
379 inline void StEmcGeom::getXYZ(const Int_t m, const Int_t e, const Int_t s, Float_t &x,Float_t &y,Float_t &z) const
380 {
381  Float_t phi;
382  if(m<=60) z = mZlocal[e-1];
383  else z =-mZlocal[e-1];
384  getPhi(m,s,phi);
385  x = mRadius*cos(phi);
386  y = mRadius*sin(phi);
387 }
388 // _____________________________________________________________________
389 inline Int_t StEmcGeom::getXYZ(const Int_t softId, Float_t &x,Float_t &y,Float_t &z) const
390 {
391  Int_t m, e, s;
392  if(!getBin(softId, m, e, s)){
393  getXYZ(m,e,s, x,y,z);
394  return 0;
395  }
396  else return 1;
397 }
398 // _____________________________________________________________________
399 inline Int_t StEmcGeom::getXYZfromGeant(const Int_t ivid,Float_t &x,Float_t &y,Float_t &z)
400 {
401  Int_t m, e, s, det;
402  if(getVolId(ivid, m,e,s,det) == 0){
403  getXYZ(m,e,s, x,y,z);
404  return 0;
405  }
406  else return 1;
407 }
408 // _____________________________________________________________________
409 inline Int_t StEmcGeom::getEta(const Int_t m, const Int_t e, Float_t &eta) const
410 {
411  if(!checkModule(m) && !checkEta(e)){
412  if(m <= mNModule/2) eta = mEta[e-1]; // West part of EMC
413  else eta = -mEta[e-1]; // East part of EMC
414  return 0;
415  }
416  else return 1;
417 }
418 // _____________________________________________________________________
419 inline Int_t StEmcGeom::getEta(const Int_t softId, Float_t &eta) const
420 {
421  Int_t m,e,s;
422  if(getBin(softId,m,e,s) == 0) {
423  return getEta(m,e,eta);
424  }
425  return 1;
426 }
427 // _____________________________________________________________________
428 inline Int_t StEmcGeom::getTheta(const Int_t m, const Int_t e, Float_t &theta) const
429 {
430  Float_t etaw;
431  if(!getEta(m,e, etaw)) { theta = 2.*atan(exp(-etaw)); return 0;}
432  else return 1;
433 }
434 // _____________________________________________________________________
435 inline Int_t StEmcGeom::getTheta(const Int_t softId, Float_t &theta) const
436 {
437  Int_t m,e,s;
438  if(getBin(softId,m,e,s) == 0) {
439  return getTheta(m,e,theta);
440  }
441  return 1;
442 }
443 // _____________________________________________________________________
444 inline Int_t StEmcGeom::getPhi(const Int_t m, const Int_t s, Float_t &phi) const
445 {
446  //
447  // -pi <= phi < pi
448  //
449  Int_t iphi, im;
450  if(!checkModule(m) && !checkSub(s)){
451  Double_t phiW; // phi in system of module
452 
453  // change sign befor -mPhi
454  if(m <= mNModule/2) {phiW = -mPhi[s-1]; im = 0; iphi=m-1;} // West part of EMC
455  else {phiW = mPhi[s-1]; im = 1; iphi=m-mNModule/2-1;} // East part of EMC
456 
457  phiW += mPhiOffset[im] + mPhiStep[im]*iphi;
458 
459  while(phiW >= C_PI) phiW -= C_2PI;
460  while(phiW < -C_PI) phiW += C_2PI;
461  if(phiW > (C_PI-0.0001)) phiW = -C_PI; // -pi<=phi<phi
462 
463  phi = phiW;
464  return 0;
465  }
466  else return 1;
467 }
468 // _____________________________________________________________________
469 inline Int_t StEmcGeom::getPhi(const Int_t softId, Float_t &phi) const
470 {
471  Int_t m,e,s;
472  if(getBin(softId,m,e,s) == 0) {
473  return getPhi(m,s,phi);
474  }
475  return 1;
476 }
477 // _____________________________________________________________________
478 inline Int_t StEmcGeom::getEtaPhi(const Int_t softId, Float_t &eta, Float_t &phi) const
479 {
480  Int_t m=0, e=0, s=0;
481  if(!getBin(softId, m, e, s)) {
482  getEta(m, e, eta);
483  getPhi(m, s, phi);
484  return 0;
485  }
486  else return 1;
487 }
488 // _____________________________________________________________________
489 inline Int_t StEmcGeom::getPhiModule(const Int_t m, Float_t &phi) const
490 {
491  if(!checkModule(m)) { phi=mPhiModule[m-1]; return 0;}
492  else return 1;
493 }
494 // _____________________________________________________________________
495 inline Int_t StEmcGeom::getIndex(const Float_t &x, const TArrayF &arr) const
496 {
497  // Arr is array of boundaries
498  for(Int_t i=1; i<arr.GetSize(); i++){
499  if(x>=arr[i-1] && x<arr[i]) return i; // x is in i cell
500  }
501  return -1;
502 }
503 
504 #ifndef WSUJoseMarch
505 inline Int_t StEmcGeom::getId(const Int_t m, const Int_t e, const Int_t s, Int_t &softId) const
506 {
507  if(!checkModule(m) && !checkEta(e) && !checkSub(s)) {
508  softId = mNes*(m-1) + mNEta*(s-1) + e;
509  return 0;
510  }
511  else {
512  LOG_WARN << Form("<W> getId(2001 Aug Scheme) | Det %i bad index m %i e %i s %i ",mDetector,m,e,s) << endm;
513  return 1;
514  }
515 }
516 
517 inline Int_t
518 StEmcGeom::getBin(const Int_t softId, Int_t &m, Int_t &e, Int_t &s) const
519 {
520  static Int_t j, wid;
521  if(!checkId(softId)) {
522  wid = softId - 1; // from 0 to MAX-1
523  m = wid/mNes + 1;
524  j = wid - mNes*(m-1);
525  s = j/mNEta + 1;
526  e = j%mNEta + 1;
527  return 0;
528  }
529  else return 1;
530 }
531 #endif
532 
533 #ifdef WSUJoseMarch
534 //
535 // March 2001 - now obsolete (7 Aug 2001)
536 //
537 inline Int_t
538 StEmcGeom::getId(const Int_t m, const Int_t e, const Int_t s,Int_t &rid)
539 {
540  if(!checkModule(m) && !checkEta(e) && !checkSub(s)){
541  if(mDetector==1 || mDetector==2) { // only for bemc and bprs
542  rid = 40*(m-1) + 20*(s-1) + (21-e);
543  }
544  if(mDetector==3 || mDetector==4) { // for BSMDE and BSDDP
545  rid = mNes*(m-1) + mNEta*(s-1) + (mNEta - e) + 1;
546  }
547  return 0;
548  }
549  else {
550  printf("<W> Det %i bad index m %i e %i s %i \n", mDetector, m, e, s);
551  return 1;
552  }
553 }
554 
555 inline Int_t
556 StEmcGeom::getBin(const Int_t rid,Int_t &m,Int_t &e,Int_t &s)
557 {
558  //
559  // 15-mar-2001 for transition from environment (Jose) numeration to
560  // usual EMC'c numbering. This is standard now !!!
561  // It is work only for BEMC and BPRS (15-mar-2001)
562  //
563  // if(mDetector<1 || mDetector>3) {cout<<" Wrong number of detector "
564  // <<mDetector<<endl; return 1;}
565  Int_t idw;
566  if(mDetector==1 || mDetector==2) { // BEMC and BPRS
567  if(checkId(rid) == 1) return 1;
568  m = (rid - 1) / 40 + 1; // Module number
569  idw = (rid - 1) % 40; // change from 0 to 39
570  s = idw/20 + 1;
571  e = 20 - idw%20;
572  return 0; // zero is good
573  }
574  else if(mDetector==3 || mDetector==4) { // BSMDE and BSMDP
575  if(checkId(rid) == 1) return 1;
576  m = (rid - 1) / mNes + 1; // Module number
577  idw = (rid - 1) % mNes; // change from 0 to mNes - 1
578  s = idw/mNEta + 1;
579  e = mNEta - idw%mNEta;
580  return 0; // zero is good
581  }
582  else {
583  cout<<" Wrong number of detector "<<mDetector<<endl;
584  return 1;
585  }
586 }
587 
588 #endif
589 
590 #endif
Int_t getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const
Definition: StEmcGeom.h:321