Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-12 04:20:16

0001 //----------Author's Name: B.Fabbro, F.X.Gentit + EB table from P.Jarry  DSM/IRFU/SPP CEA-Saclay

0002 //----------Copyright:Those valid for CEA software

0003 //----------Modified:28/01/2014

0004 
0005 #include "CalibCalorimetry/EcalCorrelatedNoiseAnalysisAlgos/interface/TEcnaNumbering.h"
0006 
0007 //--------------------------------------

0008 //  TEcnaNumbering.cc

0009 //  Class creation: 30 September 2005

0010 //  Documentation: see TEcnaNumbering.h

0011 //--------------------------------------

0012 
0013 ClassImp(TEcnaNumbering);
0014 //-------------------------------------------------------------------------------------------------------

0015 //

0016 //  Building of the numbering for the Ecal channels (EB and EE)

0017 //

0018 //  Convention for the names used here and in the other "TEcna" classes:

0019 //

0020 //  Name     Number                      Reference set   Range      Comment

0021 //

0022 //  SMTow  : Tower number                in SuperModule  [1,68]     (phi,eta) progression

0023 //  SMCrys : Crystal number              in SuperModule  [1,1700]   (phi,eta) progression

0024 //  SMEcha : Electronic channel number   in SuperModule  [0,1699]   S shape data reading order

0025 //  TowEcha: Electronic channel number   in Tower        [0,24]     S shape data reading order

0026 //

0027 //  DeeSC  : super-crystal(SC) number    in Dee          [1,200]    (IY,IX) progression

0028 //  DeeCrys: Crystal number              in Dee matrix   [1,5000]   (IY,IX) progression

0029 //  DeeEcha: Electronic channel number   in Dee matrix   [0,4999]   (IY,IX) progression (starting from 0)

0030 //  SCEcha : Electronic channel number   in SC           [1,25]     Crystal numbering for construction

0031 //-------------------------------------------------------------------------------------------------------

0032 
0033 //====================================== constructors =========================================

0034 TEcnaNumbering::TEcnaNumbering() {
0035   // Constructor without argument: call to method Init()

0036 
0037   //  std::cout << "[Info Management] CLASS: TEcnaNumbering.    CREATE OBJECT: this = " << this << std::endl;

0038 
0039   Init();
0040 }
0041 
0042 TEcnaNumbering::TEcnaNumbering(TEcnaObject* pObjectManager, const TString& SubDet) {
0043   // Constructor with argument: call to methods Init() and SetEcalSubDetector(const TString&)

0044 
0045   // std::cout << "[Info Management] CLASS: TEcnaNumbering.    CREATE OBJECT: this = " << this << std::endl;

0046 
0047   Init();
0048   Long_t i_this = (Long_t)this;
0049   pObjectManager->RegisterPointer("TEcnaNumbering", i_this);
0050 
0051   //............................ fEcal  => to be changed in fParEcal

0052   fEcal = nullptr;
0053   Long_t iParEcal = pObjectManager->GetPointerValue("TEcnaParEcal");
0054   if (iParEcal == 0) {
0055     fEcal = new TEcnaParEcal(pObjectManager, SubDet.Data()); /*fCnew++*/
0056   } else {
0057     fEcal = (TEcnaParEcal*)iParEcal;
0058   }
0059 
0060   SetEcalSubDetector(SubDet.Data());
0061 }
0062 
0063 TEcnaNumbering::TEcnaNumbering(const TString& SubDet, TEcnaParEcal* pEcal) {
0064   // Constructor with argument: call to methods Init() and SetEcalSubDetector(const TString&)

0065 
0066   // std::cout << "[Info Management] CLASS: TEcnaNumbering.    CREATE OBJECT: this = " << this << std::endl;

0067 
0068   Init();
0069   SetEcalSubDetector(SubDet.Data(), pEcal);
0070 }
0071 
0072 //====================================== destructor =========================================

0073 TEcnaNumbering::~TEcnaNumbering() {
0074   //destructor

0075 
0076   //if (fEcal   != 0){delete fEcal;   fCdelete++;}

0077 
0078   //....................... Barrel

0079   if (fT2dSMCrys != nullptr) {
0080     delete[] fT2dSMCrys;
0081     fCdelete++;
0082   }
0083   if (fT1dSMCrys != nullptr) {
0084     delete[] fT1dSMCrys;
0085     fCdelete++;
0086   }
0087   if (fT1dSMTow != nullptr) {
0088     delete[] fT1dSMTow;
0089     fCdelete++;
0090   }
0091   if (fT1dTowEcha != nullptr) {
0092     delete[] fT1dTowEcha;
0093     fCdelete++;
0094   }
0095 
0096   //....................... Endcap

0097   if (fT3dDeeCrys != nullptr) {
0098     delete[] fT3dDeeCrys;
0099     fCdelete++;
0100   }
0101   if (fT2dDeeCrys != nullptr) {
0102     delete[] fT2dDeeCrys;
0103     fCdelete++;
0104   }
0105   if (fT1dDeeCrys != nullptr) {
0106     delete[] fT1dDeeCrys;
0107     fCdelete++;
0108   }
0109   if (fT2dDeeSC != nullptr) {
0110     delete[] fT2dDeeSC;
0111     fCdelete++;
0112   }
0113   if (fT1dDeeSC != nullptr) {
0114     delete[] fT1dDeeSC;
0115     fCdelete++;
0116   }
0117   if (fT2dSCEcha != nullptr) {
0118     delete[] fT2dSCEcha;
0119     fCdelete++;
0120   }
0121   if (fT1dSCEcha != nullptr) {
0122     delete[] fT1dSCEcha;
0123     fCdelete++;
0124   }
0125   if (fT2d_jch_JY != nullptr) {
0126     delete[] fT2d_jch_JY;
0127     fCdelete++;
0128   }
0129   if (fT1d_jch_JY != nullptr) {
0130     delete[] fT1d_jch_JY;
0131     fCdelete++;
0132   }
0133   if (fT2d_ich_IX != nullptr) {
0134     delete[] fT2d_ich_IX;
0135     fCdelete++;
0136   }
0137   if (fT1d_ich_IX != nullptr) {
0138     delete[] fT1d_ich_IX;
0139     fCdelete++;
0140   }
0141   if (fT2d_DS != nullptr) {
0142     delete[] fT2d_DS;
0143     fCdelete++;
0144   }
0145   if (fT1d_DS != nullptr) {
0146     delete[] fT1d_DS;
0147     fCdelete++;
0148   }
0149   if (fT2d_DSSC != nullptr) {
0150     delete[] fT2d_DSSC;
0151     fCdelete++;
0152   }
0153   if (fT1d_DSSC != nullptr) {
0154     delete[] fT1d_DSSC;
0155     fCdelete++;
0156   }
0157   if (fT2d_DeeSCCons != nullptr) {
0158     delete[] fT2d_DeeSCCons;
0159     fCdelete++;
0160   }
0161   if (fT1d_DeeSCCons != nullptr) {
0162     delete[] fT1d_DeeSCCons;
0163     fCdelete++;
0164   }
0165   if (fT2d_RecovDeeSC != nullptr) {
0166     delete[] fT2d_RecovDeeSC;
0167     fCdelete++;
0168   }
0169   if (fT1d_RecovDeeSC != nullptr) {
0170     delete[] fT1d_RecovDeeSC;
0171     fCdelete++;
0172   }
0173 
0174   // std::cout << "[Info Management] CLASS: TEcnaNumbering.    DESTROY OBJECT: this = " << this << std::endl;

0175 }
0176 //------------------------------------------------------------- Init()

0177 void TEcnaNumbering::Init() {
0178   //Set default values and build crystal numbering table

0179 
0180   //.............................. Initialisations

0181   fTTBELL = '\007';
0182   fgMaxCar = (Int_t)512;
0183 
0184   //....................... Barrel

0185   fT2dSMCrys = nullptr;
0186   fT1dSMCrys = nullptr;
0187   fT1dSMTow = nullptr;
0188   fT1dTowEcha = nullptr;
0189 
0190   fCodeChNumberingLvrbBot = "bottom";
0191   fCodeChNumberingLvrbTop = "top";
0192 
0193   //....................... Endcap

0194   fT3dDeeCrys = nullptr;
0195   fT2dDeeCrys = nullptr;
0196   fT1dDeeCrys = nullptr;
0197   fT2dDeeSC = nullptr;
0198   fT1dDeeSC = nullptr;
0199   fT2dSCEcha = nullptr;
0200   fT1dSCEcha = nullptr;
0201   fT2d_jch_JY = nullptr;
0202   fT1d_jch_JY = nullptr;
0203   fT2d_ich_IX = nullptr;
0204   fT1d_ich_IX = nullptr;
0205   fT2d_DS = nullptr;
0206   fT1d_DS = nullptr;
0207   fT2d_DSSC = nullptr;
0208   fT1d_DSSC = nullptr;
0209   fT2d_DeeSCCons = nullptr;
0210   fT1d_DeeSCCons = nullptr;
0211   fT2d_RecovDeeSC = nullptr;
0212   fT1d_RecovDeeSC = nullptr;
0213 
0214   fCodeChNumberingITP1Bot = "bottom";  // ==> Type 1 Interface plate  IPT1 (a faire)

0215   fCodeChNumberingITP2Top = "top";     // ==> Type 2 Interface plate  IPT2 (a faire)

0216 
0217   //------------------ Init pointers on the CNA objects

0218   fEcal = nullptr;
0219 }
0220 // end of Init()

0221 //------------------------------------------------------------- SetEcalSubDetector(...)

0222 void TEcnaNumbering::SetEcalSubDetector(const TString& SubDet, TEcnaParEcal* pEcal) {
0223   //Set the current subdetector flag and the current subdetector parameters

0224 
0225   fEcal = nullptr;
0226   if (pEcal == nullptr) {
0227     fEcal = new TEcnaParEcal(SubDet.Data()); /*fCnew++*/
0228   } else {
0229     fEcal = pEcal;
0230   }
0231 
0232   Int_t MaxCar = fgMaxCar;
0233   fFlagSubDet.Resize(MaxCar);
0234   fFlagSubDet = fEcal->GetEcalSubDetector();  // fFlagSubDet = "EB" or "EE"

0235 
0236   if (fFlagSubDet == "EB") {
0237     BuildBarrelCrysTable();
0238   }
0239   if (fFlagSubDet == "EE") {
0240     BuildEndcapCrysTable();
0241     BuildEndcapSCTable();
0242   }
0243 }
0244 
0245 void TEcnaNumbering::SetEcalSubDetector(const TString& SubDet) {
0246   //Set the current subdetector flag and the current subdetector parameters

0247 
0248   Int_t MaxCar = fgMaxCar;
0249   fFlagSubDet.Resize(MaxCar);
0250   fFlagSubDet = fEcal->GetEcalSubDetector();  // fFlagSubDet = "EB" or "EE"

0251 
0252   if (fFlagSubDet == "EB") {
0253     BuildBarrelCrysTable();
0254   }
0255   if (fFlagSubDet == "EE") {
0256     BuildEndcapCrysTable();
0257     BuildEndcapSCTable();
0258   }
0259 }
0260 
0261 //====================================================================================================

0262 //

0263 //

0264 //                                  B   A   R   R   E   L

0265 //

0266 //

0267 //====================================================================================================

0268 //

0269 //               SMCrys <-> (SMTow, TowEcha) correspondance table (from Patrick Jarry)

0270 //

0271 //====================================================================================================

0272 void TEcnaNumbering::BuildBarrelCrysTable() {
0273   // Build the correspondance table: SMCrys <-> (SMTow, TowEcha) for the ECAL BARREL

0274   //

0275   //  From CMS Internal Note  "CMS ECAL Barrel channel numbering"

0276   //

0277   //       Name      Number                     Reference set   Range        Comment

0278   //

0279   //       SMTow  : Tower number               in SuperModule  [1,68]      (phi,eta) progression

0280   //       SMCrys : Crystal number             in SuperModule  [1,1700]    (phi,eta) progression

0281   //       SMEcha : Electronic channel number  in SuperModule  [0,1699]    S shape data reading order

0282   //

0283   //       TowEcha: Electronic channel number  in Tower        [0,24]      S shape data reading order

0284   //

0285   //

0286   //   fill the 2D array:  fT2dSMCrys[n_SMTow][n_TowEcha]

0287   //

0288   //   and the 1d arrays:  fT1dSMTow[i_SMCrys] and fT1dTowEcha[i_SMCrys]

0289   //

0290   //-----------------------------------------------------------------------------------------------------

0291 
0292   if (fT2dSMCrys == nullptr) {
0293     Int_t MaxSMTow = fEcal->MaxTowInSM();
0294     Int_t MaxTowEcha = fEcal->MaxCrysInTow();
0295     Int_t MaxSMCrys = fEcal->MaxCrysInSM();
0296 
0297     //................... Allocation and Init CrysNumbersTable

0298 
0299     fT2dSMCrys = new Int_t*[MaxSMTow];
0300     fCnew++;
0301     fT1dSMCrys = new Int_t[MaxSMTow * MaxTowEcha];
0302     fCnew++;
0303     for (Int_t i_SMTow = 0; i_SMTow < MaxSMTow; i_SMTow++) {
0304       fT2dSMCrys[i_SMTow] = &fT1dSMCrys[0] + i_SMTow * MaxTowEcha;
0305     }
0306     for (Int_t i = 0; i < MaxSMTow; i++) {
0307       for (Int_t j = 0; j < MaxTowEcha; j++) {
0308         fT2dSMCrys[i][j] = 0;
0309       }
0310     }
0311 
0312     fT1dSMTow = new Int_t[MaxSMCrys];
0313     fCnew++;
0314     for (Int_t i = 0; i < MaxSMCrys; i++) {
0315       fT1dSMTow[i] = 0;
0316     }
0317 
0318     fT1dTowEcha = new Int_t[MaxSMCrys];
0319     fCnew++;
0320     for (Int_t i = 0; i < MaxSMCrys; i++) {
0321       fT1dTowEcha[i] = 0;
0322     }
0323 
0324     //........................ Build table

0325     Int_t m2 = (Int_t)2;
0326     Int_t m26 = (Int_t)26;
0327 
0328     //      Int_t jch_type[2][26];

0329     Int_t** jch_type = new Int_t*[m2];
0330     fCnew++;
0331     Int_t* jch_type_d1 = new Int_t[m2 * m26];
0332     fCnew++;
0333     for (Int_t i_m2 = 0; i_m2 < m2; i_m2++) {
0334       jch_type[i_m2] = &jch_type_d1[0] + i_m2 * m26;
0335     }
0336 
0337     for (Int_t k = 25; k >= 21; k--) {
0338       jch_type[0][k] = 25 - k;
0339     }  //  k = 25,24,23,22,21 -> jch_type[0][k] = 0,1,2,3,4

0340     for (Int_t k = 16; k <= 20; k++) {
0341       jch_type[0][k] = k - 16;
0342     }  //  k = 16,17,18,19,20 -> jch_type[0][k] = 0,1,2,3,4

0343     for (Int_t k = 15; k >= 11; k--) {
0344       jch_type[0][k] = 15 - k;
0345     }  //  k = 15,14,12,13,11 -> jch_type[0][k] = 0,1,2,3,4

0346     for (Int_t k = 6; k <= 10; k++) {
0347       jch_type[0][k] = k - 6;
0348     }  //  k =  6, 7, 8, 9,10 -> jch_type[0][k] = 0,1,2,3,4

0349     for (Int_t k = 5; k >= 1; k--) {
0350       jch_type[0][k] = 5 - k;
0351     }  //  k =  5, 4, 3, 2, 1 -> jch_type[0][k] = 0,1,2,3,4

0352 
0353     for (Int_t k = 1; k <= 5; k++) {
0354       jch_type[1][k] = k - 1;
0355     }  //  k =  1, 2, 3, 4, 5 -> jch_type[1][k] = 0,1,2,3,4

0356     for (Int_t k = 10; k >= 6; k--) {
0357       jch_type[1][k] = 10 - k;
0358     }  //  k = 10, 9, 8, 7, 6 -> jch_type[1][k] = 0,1,2,3,4

0359     for (Int_t k = 11; k <= 15; k++) {
0360       jch_type[1][k] = k - 11;
0361     }  //  k = 11,12,13,14,15 -> jch_type[1][k] = 0,1,2,3,4

0362     for (Int_t k = 20; k >= 16; k--) {
0363       jch_type[1][k] = 20 - k;
0364     }  //  k = 20,19,18,17,16 -> jch_type[1][k] = 0,1,2,3,4

0365     for (Int_t k = 21; k <= 25; k++) {
0366       jch_type[1][k] = k - 21;
0367     }  //  k = 21,22,23,24,25 -> jch_type[1][k] = 0,1,2,3,4

0368 
0369     //      Int_t ich_type[2][26];

0370     Int_t** ich_type = new Int_t*[m2];
0371     fCnew++;
0372     Int_t* ich_type_d1 = new Int_t[m2 * m26];
0373     fCnew++;
0374     for (Int_t i_m2 = 0; i_m2 < m2; i_m2++) {
0375       ich_type[i_m2] = &ich_type_d1[0] + i_m2 * m26;
0376     }
0377 
0378     for (Int_t k = 25; k >= 21; k--) {
0379       ich_type[0][k] = 0;
0380     }  //  k = 25,24,23,22,21 -> ich_type[0][k] = 0

0381     for (Int_t k = 16; k <= 20; k++) {
0382       ich_type[0][k] = 1;
0383     }  //  k = 16,17,18,19,20 -> ich_type[0][k] = 1

0384     for (Int_t k = 15; k >= 11; k--) {
0385       ich_type[0][k] = 2;
0386     }  //  k = 15,14,12,13,11 -> ich_type[0][k] = 2

0387     for (Int_t k = 6; k <= 10; k++) {
0388       ich_type[0][k] = 3;
0389     }  //  k =  6, 7, 8, 9,10 -> ich_type[0][k] = 3

0390     for (Int_t k = 5; k >= 1; k--) {
0391       ich_type[0][k] = 4;
0392     }  //  k =  5, 4, 3, 2, 1 -> ich_type[0][k] = 4

0393 
0394     for (Int_t k = 1; k <= 5; k++) {
0395       ich_type[1][k] = 0;
0396     }  //  k =  1, 2, 3, 4, 5 -> ich_type[1][k] = 0

0397     for (Int_t k = 10; k >= 6; k--) {
0398       ich_type[1][k] = 1;
0399     }  //  k = 10, 9, 8, 7, 6 -> ich_type[1][k] = 1

0400     for (Int_t k = 11; k <= 15; k++) {
0401       ich_type[1][k] = 2;
0402     }  //  k = 11,12,13,14,15 -> ich_type[1][k] = 2

0403     for (Int_t k = 20; k >= 16; k--) {
0404       ich_type[1][k] = 3;
0405     }  //  k = 20,19,18,17,16 -> ich_type[1][k] = 3

0406     for (Int_t k = 21; k <= 25; k++) {
0407       ich_type[1][k] = 4;
0408     }  //  k = 21,22,23,24,25 -> ich_type[1][k] = 4

0409 
0410     //     Int_t type[17]={0,0,0,1,1, 0,0,1,1, 0,0,1,1, 0,0,1,1};

0411     Int_t m17 = 17;
0412     Int_t* type = new Int_t[m17];
0413     fCnew++;
0414 
0415     //  0 -> LVRB at the bottom, 1 -> LVRB at the top

0416     type[0] = 0;  // M1

0417     type[1] = 0;
0418     type[2] = 0;
0419     type[3] = 1;
0420     type[4] = 1;
0421 
0422     type[5] = 0;  // M2

0423     type[6] = 0;
0424     type[7] = 1;
0425     type[8] = 1;
0426 
0427     type[9] = 0;  // M3

0428     type[10] = 0;
0429     type[11] = 1;
0430     type[12] = 1;
0431 
0432     type[13] = 0;  // M4

0433     type[14] = 0;
0434     type[15] = 1;
0435     type[16] = 1;
0436 
0437     for (Int_t tow = 0; tow < MaxSMTow; tow++)  //  tow  = 0 to 67   (MaxSMTow = 68)

0438     {
0439       for (Int_t ic = 1; ic <= MaxTowEcha; ic++)  //  ic   = 1 to 25   (MaxTowEcha = 25)

0440       {
0441         Int_t jtow = tow % 4;  //  jtow = 0,1,2,3

0442         Int_t itow = tow / 4;  //  itow = 0 to 16

0443 
0444         Int_t icrys = itow * 5 + ich_type[type[itow]][ic];  //  type[0->16] = 0,1 ,

0445                                                             //  ich_type[0->1][1->25] = 0,1,2,3,4

0446                                                             //  icrys = 0 to 84  (=> eta)

0447 
0448         Int_t jcrys = jtow * 5 + jch_type[type[itow]][ic];  //  type[0->16] = 0,1 ,

0449                                                             //  jch_type[0->1][1->25] = 0,1,2,3,4

0450                                                             //  jcrys = 0 to 19  (=> phi)

0451 
0452         Int_t n1SMCrys = icrys * 20 + jcrys + 1;  //  n1SMCrys = 1 to 1700

0453 
0454         fT2dSMCrys[tow][ic - 1] = n1SMCrys;  // fT2dSMCrys[]  : range = [1,1700]

0455         fT1dSMTow[n1SMCrys - 1] = tow + 1;   // fT1dSMTow[]   : range = [1,68]

0456         fT1dTowEcha[n1SMCrys - 1] = ic - 1;  // fT1dTowEcha[] : range = [0,24]

0457       }
0458     }
0459     // std::cout << "#TEcnaNumbering::TBuildBarrelCrysTable()> Crys Table Building done" << std::endl;

0460 
0461     delete[] jch_type;
0462     fCdelete++;
0463     delete[] jch_type_d1;
0464     fCdelete++;
0465     delete[] ich_type;
0466     fCdelete++;
0467     delete[] ich_type_d1;
0468     fCdelete++;
0469     delete[] type;
0470     fCdelete++;
0471   } else {
0472     // std::cout << "#TEcnaNumbering::TBuildBarrelCrysTable()> No Building of Crys Table since it is already done." << std::endl;

0473   }
0474 }
0475 
0476 //===============================================================================

0477 //

0478 //        GetSMCrysFrom1SMTowAnd0TowEcha

0479 //        GetSMCrysFromSMEcha

0480 //

0481 //===============================================================================

0482 Int_t TEcnaNumbering::Get1SMCrysFrom1SMTowAnd0TowEcha(const Int_t& n1SMTow, const Int_t& i0TowEcha) {
0483   //get crystal number in SM from tower number in SM

0484   // and from Electronic Channel number in tower

0485 
0486   Int_t n1SMCrys = 0;
0487 
0488   if (fT2dSMCrys == nullptr) {
0489     BuildBarrelCrysTable();
0490   }
0491 
0492   if (n1SMTow >= 1 && n1SMTow <= fEcal->MaxTowInSM()) {
0493     if (i0TowEcha >= 0 && i0TowEcha < fEcal->MaxCrysInTow()) {
0494       n1SMCrys = fT2dSMCrys[n1SMTow - 1][i0TowEcha];
0495     } else {
0496       n1SMCrys = -2;  // Electronic Cnannel in Tower out of range

0497       std::cout << "!TEcnaNumbering::Get1SMCrysFrom1SMTowAnd0TowEcha(...)> Electronic Channel in Tower out of range."
0498                 << " i0TowEcha = " << i0TowEcha << "(n1SMTow = " << n1SMTow << ")" << fTTBELL << std::endl;
0499     }
0500   } else {
0501     n1SMCrys = -3;  // Tower number in SM out of range

0502     std::cout << "!TEcnaNumbering::Get1SMCrysFrom1SMTowAnd0TowEcha(...)> Tower number in SM out of range."
0503               << " n1SMTow = " << n1SMTow << "(i0TowEcha = " << i0TowEcha << ")" << fTTBELL << std::endl;
0504   }
0505 
0506   return n1SMCrys;  // Range = [1,1700]

0507 }
0508 
0509 //===============================================================================

0510 //

0511 //                Get0TowEchaFrom1SMCrys, Get1SMTowFrom1SMCrys

0512 //

0513 //===============================================================================

0514 
0515 Int_t TEcnaNumbering::Get0TowEchaFrom1SMCrys(const Int_t& n1SMCrys) {
0516   // get Electronic Channel number in Tower from Crystal number in SuperModule

0517 
0518   Int_t i0TowEcha = -1;
0519 
0520   if (n1SMCrys >= 1 && n1SMCrys <= fEcal->MaxCrysInSM()) {
0521     i0TowEcha = fT1dTowEcha[n1SMCrys - 1];
0522   } else {
0523     i0TowEcha = -2;
0524     std::cout << "!TEcnaNumbering::Get0TowEchaFrom1SMCrys(...)> Crystal number in SM out of range."
0525               << " n1SMCrys = " << n1SMCrys << fTTBELL << std::endl;
0526   }
0527   return i0TowEcha;  // range = [0,24]

0528 }
0529 
0530 Int_t TEcnaNumbering::Get1SMTowFrom1SMCrys(const Int_t& n1SMCrys) {
0531   // get Tower number in SM (range [1,68]) from Crystal number in SuperModule (range [1,1700])

0532 
0533   Int_t n1SMtox = 0;
0534 
0535   if (n1SMCrys >= 1 && n1SMCrys <= fEcal->MaxCrysInSM()) {
0536     n1SMtox = fT1dSMTow[n1SMCrys - 1];
0537   } else {
0538     n1SMtox = -1;
0539     std::cout << "!TEcnaNumbering::Get1SMTowFrom1SMCrys(...)> Crystal number in SM out of range."
0540               << " n1SMCrys = " << n1SMCrys << fTTBELL << std::endl;
0541   }
0542   return n1SMtox;  // range = [1,68]

0543 }
0544 
0545 //===============================================================================

0546 //

0547 //          Get0TowEchaFrom0SMEcha

0548 //          Get1SMTowFrom0SMEcha

0549 //

0550 //===============================================================================

0551 
0552 Int_t TEcnaNumbering::Get0TowEchaFrom0SMEcha(const Int_t& i0SMEcha) {
0553   //get electronic channel number in tower from electronic channel number in SM

0554 
0555   Int_t n1SMTow = i0SMEcha / fEcal->MaxCrysInTow() + 1;
0556   Int_t i0TowEcha = i0SMEcha - fEcal->MaxCrysInTow() * (n1SMTow - 1);
0557 
0558   return i0TowEcha;  // range = [0,24]

0559 }
0560 
0561 Int_t TEcnaNumbering::Get1SMTowFrom0SMEcha(const Int_t& i0SMEcha) {
0562   //get tower number from electronic channel number in SM

0563 
0564   Int_t n1SMTow = i0SMEcha / fEcal->MaxCrysInTow() + 1;
0565 
0566   return n1SMTow;  // range = [1,68]

0567 }
0568 
0569 Int_t TEcnaNumbering::Get0SMEchaFrom1SMTowAnd0TowEcha(const Int_t& n1SMTow, const Int_t& i0TowEcha) {
0570   //get tower number from electronic channel number in SM

0571 
0572   Int_t i0SMEcha = (n1SMTow - 1) * fEcal->MaxCrysInTow() + i0TowEcha;
0573 
0574   return i0SMEcha;
0575 }
0576 //===========================================================================

0577 //

0578 //                        GetHashedNumberFromIEtaAndIPhi

0579 //

0580 //===========================================================================

0581 Int_t TEcnaNumbering::GetHashedNumberFromIEtaAndIPhi(const Int_t& IEta, const Int_t& IPhi) {
0582   Int_t Hashed = 0;
0583 
0584   if (IEta > 0) {
0585     Hashed = (85 + IEta - 1) * 360 + IPhi - 1;
0586   }
0587   if (IEta < 0) {
0588     Hashed = (85 + IEta) * 360 + IPhi - 1;
0589   }
0590 
0591   return Hashed;
0592 }
0593 
0594 Int_t TEcnaNumbering::GetIEtaFromHashed(const Int_t& Hashed, const Int_t& SMNumber) {
0595   Int_t IEta = 0;
0596 
0597   if (GetSMHalfBarrel(SMNumber) == "EB+") {
0598     IEta = Hashed / 360 - 85 + 1;
0599   }
0600   if (GetSMHalfBarrel(SMNumber) == "EB-") {
0601     IEta = 85 + Hashed / 360;
0602   }
0603 
0604   return IEta;
0605 }
0606 
0607 Int_t TEcnaNumbering::GetIPhiFromHashed(const Int_t& Hashed) {
0608   Int_t IPhi = Hashed % 360 + 1;
0609 
0610   return IPhi;
0611 }
0612 //===========================================================================

0613 //

0614 //                        GetTowerLvrbType

0615 //

0616 //===========================================================================

0617 TString TEcnaNumbering::GetStinLvrbType(const Int_t& n1SMTow) {
0618   TString lvrb_type = GetTowerLvrbType(n1SMTow);
0619   return lvrb_type;
0620 }
0621 TString TEcnaNumbering::GetTowerLvrbType(const Int_t& n1SMTow) {
0622   //gives the LVRB type of the crystal numbering of tower

0623 
0624   TString type = fCodeChNumberingLvrbTop;  // => default value

0625 
0626   if (n1SMTow >= 1 && n1SMTow <= 12) {
0627     type = fCodeChNumberingLvrbBot;
0628   }
0629   if (n1SMTow >= 21 && n1SMTow <= 28) {
0630     type = fCodeChNumberingLvrbBot;
0631   }
0632   if (n1SMTow >= 37 && n1SMTow <= 44) {
0633     type = fCodeChNumberingLvrbBot;
0634   }
0635   if (n1SMTow >= 53 && n1SMTow <= 60) {
0636     type = fCodeChNumberingLvrbBot;
0637   }
0638 
0639   return type;
0640 }
0641 
0642 //==============================================================================

0643 //

0644 //       GetEta, GetEtaMin, GetEtaMax,  GetIEtaMin, GetIEtaMax

0645 //

0646 //==============================================================================

0647 Double_t TEcnaNumbering::GetEta(const Int_t& n1EBSM, const Int_t& n1SMTow, const Int_t& i0TowEcha) {
0648   //Gives Eta for a given (n1EBSM, n1SMTow, i0TowEcha)

0649 
0650   Double_t eta = (Double_t)0.;
0651 
0652   Int_t max_crys_eta_in_tower = fEcal->MaxCrysEtaInTow();
0653   Int_t max_tow_eta_in_sm = fEcal->MaxTowEtaInSM();
0654   Int_t max_sm_in_barrel = fEcal->MaxSMInEB();
0655 
0656   if (n1EBSM >= 1 && n1EBSM <= max_sm_in_barrel) {
0657     for (Int_t i_sm_tow_eta = 0; i_sm_tow_eta < max_tow_eta_in_sm; i_sm_tow_eta++) {
0658       Int_t i_crys_eta_min = (Int_t)(1 + i_sm_tow_eta * (max_crys_eta_in_tower - 1));
0659       Int_t i_crys_eta_max = (Int_t)((i_sm_tow_eta + 1) * (max_crys_eta_in_tower - 1));
0660       Int_t i_crys_eta = (Int_t)(i_sm_tow_eta * max_crys_eta_in_tower);
0661       // = 0,..,16 -> last xtal in eta for the previous tower

0662       Double_t d_echa_eta = (Double_t)(i0TowEcha / max_crys_eta_in_tower);  // = 0,1,2,3,4

0663 
0664       if (n1SMTow >= i_crys_eta_min && n1SMTow <= i_crys_eta_max) {
0665         if (GetTowerLvrbType(n1SMTow) == fCodeChNumberingLvrbTop) {
0666           eta = (Double_t)(i_crys_eta) + d_echa_eta + 1;
0667         }
0668         if (GetTowerLvrbType(n1SMTow) == fCodeChNumberingLvrbBot) {
0669           eta = (Double_t)(i_crys_eta + max_crys_eta_in_tower) - d_echa_eta;
0670         }
0671       }
0672     }
0673     if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0674       eta = -eta;
0675     }
0676   } else {
0677     std::cout << "TEcnaNumbering::GetEta(...)> SM = " << n1EBSM << ". Out of range (range = [1," << fEcal->MaxSMInEB()
0678               << "])" << fTTBELL << std::endl;
0679   }
0680   return eta;
0681 }
0682 //-------------------------------------------------------------------------------------

0683 Double_t TEcnaNumbering::GetEtaMin(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0684   //Gives EtaMin for a given Tower

0685 
0686   Int_t max_tow_eta_in_sm = fEcal->MaxTowEtaInSM();
0687   Int_t max_crys_eta_in_tower = fEcal->MaxCrysEtaInTow();
0688 
0689   Double_t eta_min = (Double_t)0.;
0690 
0691   for (Int_t i_sm_tow_eta = 0; i_sm_tow_eta < max_tow_eta_in_sm; i_sm_tow_eta++) {
0692     Int_t i_crys_eta_min = (Int_t)(1 + i_sm_tow_eta * (max_crys_eta_in_tower - 1));
0693     Int_t i_crys_eta_max = (Int_t)((i_sm_tow_eta + 1) * (max_crys_eta_in_tower - 1));
0694     Int_t i_crys_eta = (Int_t)(i_sm_tow_eta * max_crys_eta_in_tower);
0695 
0696     if (n1SMTow >= i_crys_eta_min && n1SMTow <= i_crys_eta_max) {
0697       if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0698         eta_min = (Double_t)i_crys_eta;
0699       }
0700       if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0701         eta_min = -(Double_t)(i_crys_eta + max_crys_eta_in_tower);
0702       }
0703     }
0704   }
0705   return eta_min;
0706 }
0707 //------------------------------------------------------------------------------------

0708 Double_t TEcnaNumbering::GetEtaMax(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0709   //Gives EtaMax for a given Tower

0710 
0711   Int_t max_tow_eta_in_sm = fEcal->MaxTowEtaInSM();
0712   Int_t max_crys_eta_in_tower = fEcal->MaxCrysEtaInTow();
0713 
0714   Double_t eta_max = (max_crys_eta_in_tower - 1);
0715 
0716   for (Int_t i_sm_tow_eta = 0; i_sm_tow_eta < max_tow_eta_in_sm; i_sm_tow_eta++) {
0717     Int_t i_crys_eta_min = (Int_t)(1 + i_sm_tow_eta * (max_crys_eta_in_tower - 1));
0718     Int_t i_crys_eta_max = (Int_t)((i_sm_tow_eta + 1) * (max_crys_eta_in_tower - 1));
0719     Int_t i_crys_eta = (Int_t)(i_sm_tow_eta * max_crys_eta_in_tower);
0720 
0721     if (n1SMTow >= i_crys_eta_min && n1SMTow <= i_crys_eta_max) {
0722       if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0723         eta_max = (Double_t)(i_crys_eta + max_crys_eta_in_tower);
0724       }
0725       if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0726         eta_max = -(Double_t)i_crys_eta;
0727       }
0728     }
0729   }
0730 
0731   return eta_max;
0732 }
0733 
0734 Double_t TEcnaNumbering::GetIEtaMin(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0735   //Gives IEtaMin for a given (n1EBSM, n1SMTow)

0736 
0737   Double_t i_eta_min = (Int_t)0.;
0738 
0739   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0740     i_eta_min = (Double_t)GetEtaMin(n1EBSM, n1SMTow) + (Double_t)0.5;
0741   }
0742   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0743     i_eta_min = (Double_t)GetEtaMin(n1EBSM, n1SMTow) - (Double_t)0.5;
0744   }
0745 
0746   return i_eta_min;
0747 }
0748 
0749 Double_t TEcnaNumbering::GetIEtaMax(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0750   //Gives IEtaMax for a given (n1EBSM, n1SMTow)

0751 
0752   Double_t i_eta_max = (Int_t)0.;
0753 
0754   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0755     i_eta_max = (Double_t)GetEtaMax(n1EBSM, n1SMTow) + (Double_t)0.5;
0756   }
0757   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0758     i_eta_max = (Double_t)GetEtaMax(n1EBSM, n1SMTow) - (Double_t)0.5;
0759   }
0760 
0761   return i_eta_max;
0762 }
0763 
0764 Double_t TEcnaNumbering::GetIEtaMin(const Int_t& n1EBSM) {
0765   //Gives IEtaMin for a given (n1EBSM)

0766 
0767   Double_t i_eta_min = (Int_t)0.;
0768 
0769   Int_t n1SMTowPlus = (Int_t)1;
0770   Int_t n1SMTowMinus = (Int_t)fEcal->MaxTowInSM();
0771 
0772   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0773     i_eta_min = (Double_t)GetIEtaMin(n1EBSM, n1SMTowPlus);
0774   }
0775   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0776     i_eta_min = (Double_t)GetIEtaMin(n1EBSM, n1SMTowMinus);
0777   }
0778 
0779   return i_eta_min;
0780 }
0781 
0782 Double_t TEcnaNumbering::GetIEtaMax(const Int_t& n1EBSM) {
0783   //Gives IEtaMax for a given (n1EBSM)

0784 
0785   Double_t i_eta_max = (Int_t)0.;
0786 
0787   Int_t n1SMTowPlus = (Int_t)fEcal->MaxTowInSM();
0788   Int_t n1SMTowMinus = (Int_t)1;
0789 
0790   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0791     i_eta_max = (Double_t)GetIEtaMax(n1EBSM, n1SMTowPlus);
0792   }
0793   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0794     i_eta_max = (Double_t)GetIEtaMax(n1EBSM, n1SMTowMinus);
0795   }
0796 
0797   return i_eta_max;
0798 }
0799 
0800 //==============================================================================

0801 //

0802 //    GetSMCentralPhi, GetPhiInSM, GetPhi,

0803 //    GetPhiMin, GetPhiMax, GetIPhiMin, GetIPhiMax

0804 //

0805 //==============================================================================

0806 Double_t TEcnaNumbering::GetSMCentralPhi(const Int_t& n1EBSM) {
0807   //Gives the central phi value of the SuperModule

0808 
0809   Double_t central_phi = (Double_t)10.;  //  DEFAULT = SM1

0810 
0811   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0812     central_phi = 10. + (Double_t)20. * (n1EBSM - 1);
0813   }
0814   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0815     central_phi = 10. + (Double_t)20. * (n1EBSM - 19);
0816   }
0817 
0818   return central_phi;
0819 }
0820 //------------------------------------------------------------------------------------

0821 Double_t TEcnaNumbering::GetPhiInSM(const Int_t& n1EBSM, const Int_t& n1SMTow, const Int_t& i0TowEcha) {
0822   //Gives Phi for a given (n1EBSM, n1SMTow, i0TowEcha)

0823 
0824   Double_t phi_in_SM = (Double_t)0.;
0825 
0826   Int_t rest_temp = (Int_t)(n1SMTow % (fEcal->MaxCrysPhiInTow() - 1));  // "Phi" of the tower = 1,2,3,0

0827 
0828   if (n1EBSM >= 1 && n1EBSM <= fEcal->MaxSMInEB()) {
0829     if (rest_temp == 1) {
0830       phi_in_SM = (Double_t)15.;
0831     }
0832     if (rest_temp == 2) {
0833       phi_in_SM = (Double_t)10.;
0834     }
0835     if (rest_temp == 3) {
0836       phi_in_SM = (Double_t)5.;
0837     }
0838     if (rest_temp == 0) {
0839       phi_in_SM = (Double_t)0.;
0840     }
0841 
0842     if (GetTowerLvrbType(n1SMTow) == fCodeChNumberingLvrbTop) {
0843       if (i0TowEcha == 4 || i0TowEcha == 5 || i0TowEcha == 14 || i0TowEcha == 15 || i0TowEcha == 24) {
0844         phi_in_SM = phi_in_SM + 0;
0845       }
0846 
0847       if (i0TowEcha == 3 || i0TowEcha == 6 || i0TowEcha == 13 || i0TowEcha == 16 || i0TowEcha == 23) {
0848         phi_in_SM = phi_in_SM + 1;
0849       }
0850 
0851       if (i0TowEcha == 2 || i0TowEcha == 7 || i0TowEcha == 12 || i0TowEcha == 17 || i0TowEcha == 22) {
0852         phi_in_SM = phi_in_SM + 2;
0853       }
0854 
0855       if (i0TowEcha == 1 || i0TowEcha == 8 || i0TowEcha == 11 || i0TowEcha == 18 || i0TowEcha == 21) {
0856         phi_in_SM = phi_in_SM + 3;
0857       }
0858 
0859       if (i0TowEcha == 0 || i0TowEcha == 9 || i0TowEcha == 10 || i0TowEcha == 19 || i0TowEcha == 20) {
0860         phi_in_SM = phi_in_SM + 4;
0861       }
0862     }
0863     if (GetTowerLvrbType(n1SMTow) == fCodeChNumberingLvrbBot) {
0864       if (i0TowEcha == 20 || i0TowEcha == 19 || i0TowEcha == 10 || i0TowEcha == 9 || i0TowEcha == 0) {
0865         phi_in_SM = phi_in_SM + 0;
0866       }
0867 
0868       if (i0TowEcha == 21 || i0TowEcha == 18 || i0TowEcha == 11 || i0TowEcha == 8 || i0TowEcha == 1) {
0869         phi_in_SM = phi_in_SM + 1;
0870       }
0871 
0872       if (i0TowEcha == 22 || i0TowEcha == 17 || i0TowEcha == 12 || i0TowEcha == 7 || i0TowEcha == 2) {
0873         phi_in_SM = phi_in_SM + 2;
0874       }
0875 
0876       if (i0TowEcha == 23 || i0TowEcha == 16 || i0TowEcha == 13 || i0TowEcha == 6 || i0TowEcha == 3) {
0877         phi_in_SM = phi_in_SM + 3;
0878       }
0879 
0880       if (i0TowEcha == 24 || i0TowEcha == 15 || i0TowEcha == 14 || i0TowEcha == 5 || i0TowEcha == 4) {
0881         phi_in_SM = phi_in_SM + 4;
0882       }
0883     }
0884   } else {
0885     std::cout << "TEcnaNumbering::GetPhiInSM(...)> SM = " << n1EBSM << ". Out of range (range = [1,"
0886               << fEcal->MaxSMInEB() << "])" << fTTBELL << std::endl;
0887   }
0888   phi_in_SM = 20 - phi_in_SM;
0889   return phi_in_SM;
0890 }
0891 //---------------------------------------------------------------------------------------

0892 Double_t TEcnaNumbering::GetPhi(const Int_t& n1EBSM, const Int_t& n1SMTow, const Int_t& i0TowEcha) {
0893   //Gives Phi for a given (n1EBSM, n1SMTow, i0TowEcha)

0894 
0895   Double_t phi = (Double_t)0.;
0896 
0897   if (n1EBSM >= 1 && n1EBSM <= fEcal->MaxSMInEB()) {
0898     Double_t phiInSM = GetPhiInSM(n1EBSM, n1SMTow, i0TowEcha);
0899     Double_t phi_start = GetSMCentralPhi(n1EBSM);
0900 
0901     phi = 20 - phiInSM + phi_start - (Double_t)10.;
0902   } else {
0903     std::cout << "TEcnaNumbering::GetPhi(...)> SM = " << n1EBSM << ". Out of range (range = [1," << fEcal->MaxSMInEB()
0904               << "])" << fTTBELL << std::endl;
0905   }
0906   return phi;
0907 }
0908 
0909 //-----------------------------------------------------------------------------------------

0910 Double_t TEcnaNumbering::GetPhiMin(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0911   //Gives PhiMin for a given Tower

0912 
0913   Int_t max_crys_phi_in_tower = fEcal->MaxCrysPhiInTow();
0914 
0915   Double_t phi_min = (Double_t)0.;  // DEFAULT

0916   Double_t phi_start = GetSMCentralPhi(n1EBSM);
0917 
0918   Int_t rest_temp = (Int_t)(n1SMTow % (max_crys_phi_in_tower - 1));
0919 
0920   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0921     if (rest_temp == 1) {
0922       phi_min = phi_start + (Double_t)5.;
0923     }
0924     if (rest_temp == 2) {
0925       phi_min = phi_start + (Double_t)0.;
0926     }
0927     if (rest_temp == 3) {
0928       phi_min = phi_start - (Double_t)5.;
0929     }
0930     if (rest_temp == 0) {
0931       phi_min = phi_start - (Double_t)10.;
0932     }
0933   }
0934   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0935     if (rest_temp == 0) {
0936       phi_min = phi_start + (Double_t)5.;
0937     }
0938     if (rest_temp == 3) {
0939       phi_min = phi_start + (Double_t)0.;
0940     }
0941     if (rest_temp == 2) {
0942       phi_min = phi_start - (Double_t)5.;
0943     }
0944     if (rest_temp == 1) {
0945       phi_min = phi_start - (Double_t)10.;
0946     }
0947   }
0948   return phi_min;
0949 }
0950 //-----------------------------------------------------------------------------------------

0951 Double_t TEcnaNumbering::GetPhiMax(const Int_t& n1EBSM, const Int_t& n1SMTow) {
0952   //Gives PhiMax for a given Tower

0953 
0954   Int_t max_crys_phi_in_tower = fEcal->MaxCrysPhiInTow();
0955 
0956   Double_t phi_max = (Double_t)20.;  // DEFAULT

0957   Double_t phi_start = GetSMCentralPhi(n1EBSM);
0958 
0959   Int_t rest_temp = (Int_t)(n1SMTow % (max_crys_phi_in_tower - 1));
0960 
0961   if (GetSMHalfBarrel(n1EBSM) == "EB+") {
0962     if (rest_temp == 1) {
0963       phi_max = phi_start + (Double_t)10.;
0964     }
0965     if (rest_temp == 2) {
0966       phi_max = phi_start + (Double_t)5.;
0967     }
0968     if (rest_temp == 3) {
0969       phi_max = phi_start - (Double_t)0.;
0970     }
0971     if (rest_temp == 0) {
0972       phi_max = phi_start - (Double_t)5.;
0973     }
0974   }
0975 
0976   if (GetSMHalfBarrel(n1EBSM) == "EB-") {
0977     if (rest_temp == 0) {
0978       phi_max = phi_start + (Double_t)10.;
0979     }
0980     if (rest_temp == 3) {
0981       phi_max = phi_start + (Double_t)5.;
0982     }
0983     if (rest_temp == 2) {
0984       phi_max = phi_start - (Double_t)0.;
0985     }
0986     if (rest_temp == 1) {
0987       phi_max = phi_start - (Double_t)5.;
0988     }
0989   }
0990 
0991   return phi_max;
0992 }
0993 //-----------------------------------------------------------------------------------------

0994 Double_t TEcnaNumbering::GetPhiMin(const Int_t& n1EBSM) {
0995   //Gives PhiMin for a given SuperModule

0996 
0997   Double_t phi_min = GetSMCentralPhi(n1EBSM) - (Double_t)10.;
0998 
0999   return phi_min;
1000 }
1001 //-----------------------------------------------------------------------------------------

1002 Double_t TEcnaNumbering::GetPhiMax(const Int_t& n1EBSM) {
1003   //Gives PhiMax for a given SuperModule

1004 
1005   Double_t phi_max = GetSMCentralPhi(n1EBSM) + (Double_t)10.;
1006 
1007   return phi_max;
1008 }
1009 //-----------------------------------------------------------------------------------------

1010 Double_t TEcnaNumbering::GetJPhiMin(const Int_t& n1EBSM, const Int_t& n1SMTow) {
1011   //Gives JPhiMin for a given Tower

1012 
1013   Double_t j_phi_min = (Double_t)1.;
1014   Int_t max_crys_phi_in_tower = fEcal->MaxCrysPhiInTow();
1015   Int_t rest_temp = (Int_t)(n1SMTow % (max_crys_phi_in_tower - 1));
1016 
1017   if (rest_temp == 1) {
1018     j_phi_min = (Double_t)1. - (Double_t)0.5;
1019   }
1020   if (rest_temp == 2) {
1021     j_phi_min = (Double_t)6. - (Double_t)0.5;
1022   }
1023   if (rest_temp == 3) {
1024     j_phi_min = (Double_t)11. - (Double_t)0.5;
1025   }
1026   if (rest_temp == 0) {
1027     j_phi_min = (Double_t)16. - (Double_t)0.5;
1028   }
1029 
1030   return j_phi_min;
1031 }
1032 //-----------------------------------------------------------------------------------------

1033 Double_t TEcnaNumbering::GetJPhiMax(const Int_t& n1EBSM, const Int_t& n1SMTow) {
1034   //Gives JPhiMax for a given Tower

1035 
1036   Double_t j_phi_max = (Double_t)20.;
1037   Int_t max_crys_phi_in_tower = fEcal->MaxCrysPhiInTow();
1038   Int_t rest_temp = (Int_t)(n1SMTow % (max_crys_phi_in_tower - 1));
1039 
1040   if (rest_temp == 1) {
1041     j_phi_max = (Double_t)5. + (Double_t)0.5;
1042   }
1043   if (rest_temp == 2) {
1044     j_phi_max = (Double_t)10. + (Double_t)0.5;
1045   }
1046   if (rest_temp == 3) {
1047     j_phi_max = (Double_t)15. + (Double_t)0.5;
1048   }
1049   if (rest_temp == 0) {
1050     j_phi_max = (Double_t)20. + (Double_t)0.5;
1051   }
1052 
1053   return j_phi_max;
1054 }
1055 
1056 //-----------------------------------------------------------------------------------------

1057 Double_t TEcnaNumbering::GetJPhiMin(const Int_t& n1EBSM) {
1058   //Gives JPhiMin for a given SuperModule

1059 
1060   Double_t j_phi_min = (Double_t)1. - (Double_t)0.5;
1061 
1062   return j_phi_min;
1063 }
1064 //-----------------------------------------------------------------------------------------

1065 Double_t TEcnaNumbering::GetJPhiMax(const Int_t& n1EBSM) {
1066   //Gives JPhiMax for a given SuperModule

1067 
1068   Double_t j_phi_max = (Double_t)20. + (Double_t)0.5;
1069 
1070   return j_phi_max;
1071 }
1072 
1073 //==============================================================================

1074 //

1075 //   GetXDirectionEB, GetYDirectionEB, GetJYDirectionEB

1076 //

1077 //==============================================================================

1078 
1079 //------------------------------------------------------------------

1080 TString TEcnaNumbering::GetXDirectionEB(const Int_t& SMNumber) {
1081   TString xdirection = "x";  // DEFAULT

1082 
1083   if (GetSMHalfBarrel(SMNumber) == "EB+") {
1084     xdirection = "x";
1085   }
1086   if (GetSMHalfBarrel(SMNumber) == "EB-") {
1087     xdirection = "x";
1088   }
1089 
1090   return xdirection;
1091 }
1092 //---------------------------------------------------------

1093 TString TEcnaNumbering::GetYDirectionEB(const Int_t& SMNumber) {
1094   TString ydirection = "-x";  // DEFAULT

1095 
1096   if (GetSMHalfBarrel(SMNumber) == "EB+") {
1097     ydirection = "-x";
1098   }
1099   if (GetSMHalfBarrel(SMNumber) == "EB-") {
1100     ydirection = "-x";
1101   }
1102 
1103   return ydirection;
1104 }
1105 
1106 //---------------------------------------------------------

1107 TString TEcnaNumbering::GetJYDirectionEB(const Int_t& SMNumber) {
1108   TString jydirection = "-x";  // DEFAULT

1109 
1110   if (GetSMHalfBarrel(SMNumber) == "EB+") {
1111     jydirection = "x";
1112   }
1113   if (GetSMHalfBarrel(SMNumber) == "EB-") {
1114     jydirection = "-x";
1115   }
1116 
1117   return jydirection;
1118 }
1119 
1120 //============================================================================

1121 TString TEcnaNumbering::GetSMHalfBarrel(const Int_t& SMNumber) {
1122   //gives the half-barrel of the SM ("EB+" or "EB-")

1123 
1124   TString type = "EB-";  // => default value

1125 
1126   if (SMNumber >= 1 && SMNumber <= fEcal->MaxSMInEBPlus()) {
1127     type = "EB+";
1128   }
1129   if (SMNumber > fEcal->MaxSMInEBPlus() && SMNumber <= fEcal->MaxSMInEB()) {
1130     type = "EB-";
1131   }
1132 
1133   return type;
1134 }
1135 
1136 Int_t TEcnaNumbering::PlusMinusSMNumber(const Int_t& PlusSMNumber) {
1137   Int_t PMSMNumber = PlusSMNumber;
1138   if (PlusSMNumber > fEcal->MaxSMPhiInEB()) {
1139     PMSMNumber = -PlusSMNumber + fEcal->MaxSMPhiInEB();
1140   }
1141   return PMSMNumber;
1142 }
1143 //====================================================================================================

1144 //

1145 //

1146 //                                  E   N   D   C   A   P

1147 //

1148 //

1149 //====================================================================================================

1150 //

1151 //               DeeCrys <-> (DeeSC, SCEcha) correspondance table

1152 //               (from the barrel table given by Patrick Jarry)

1153 //

1154 //====================================================================================================

1155 void TEcnaNumbering::BuildEndcapCrysTable() {
1156   // Build the correspondance table: DeeCrys <-> (DeeSC, SCEcha) for the ECAL ENDCAP

1157   //

1158   //  From CMS Internal Note  "CMS ECAL Endcap channel numbering"

1159   //

1160   //       Name      Number                     Reference set  Range          Comment

1161   //

1162   //       DeeSC   : Super-Crystal (SC) number    in Dee      [1,200]   (IY,IX) progression

1163   //       DeeCrys : Crystal number               in Dee      [1,5000]  (IY,IX) progression

1164   //       DeeEcha : Electronic channel number    in Dee      [0,4999]  Crystal numbering for construction

1165   //

1166   //       SCEcha  : Electronic channel number    in SC       [1,25]    Crystal numbering for construction

1167   //

1168   //

1169   //   fill the 3D array:  fT3dDeeCrys[n_DeeSC][n_SCEcha][2]

1170   //

1171   //   and the 2d arrays:  fT2dDeeSC[n_DeeCrys][2] and fT2dSCEcha[n_DeeCrys][2]

1172   //

1173   //------------------------------------------------------------------------------------------------------

1174 
1175   if (fT3dDeeCrys == nullptr) {
1176     Int_t MaxDeeSC = fEcal->MaxSCEcnaInDee();  // fEcal->MaxSCEcnaInDee() = 200

1177     Int_t MaxSCEcha = fEcal->MaxCrysInSC();
1178     Int_t MaxDeeCrys = fEcal->MaxCrysEcnaInDee();
1179 
1180     Int_t MaxDirections = 2;  //  (2 directions: left and right)

1181 
1182     //................... Allocation and Init CrysNumbersTable

1183     fT3dDeeCrys = new Int_t**[MaxDeeSC];
1184     fCnew++;
1185     fT2dDeeCrys = new Int_t*[MaxDeeSC * MaxSCEcha];
1186     fCnew++;
1187     fT1dDeeCrys = new Int_t[MaxDeeSC * MaxSCEcha * MaxDirections];
1188     fCnew++;
1189 
1190     for (Int_t i_DeeSC = 0; i_DeeSC < MaxDeeSC; i_DeeSC++) {
1191       fT3dDeeCrys[i_DeeSC] = &fT2dDeeCrys[0] + i_DeeSC * MaxSCEcha;
1192       for (Int_t i_SCEcha = 0; i_SCEcha < MaxSCEcha; i_SCEcha++) {
1193         fT2dDeeCrys[i_DeeSC * MaxSCEcha + i_SCEcha] =
1194             &fT1dDeeCrys[0] + (i_DeeSC * MaxSCEcha + i_SCEcha) * MaxDirections;
1195       }
1196     }
1197     for (Int_t i = 0; i < MaxDeeSC; i++) {
1198       for (Int_t j = 0; j < MaxSCEcha; j++) {
1199         for (Int_t k = 0; k < MaxDirections; k++) {
1200           fT3dDeeCrys[i][j][k] = 0;
1201         }
1202       }
1203     }
1204 
1205     fT2dDeeSC = new Int_t*[MaxDeeCrys];
1206     fCnew++;
1207     fT1dDeeSC = new Int_t[MaxDeeCrys * MaxDirections];
1208     fCnew++;
1209     for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxDeeCrys; i_DeeCrys++) {
1210       fT2dDeeSC[i_DeeCrys] = &fT1dDeeSC[0] + i_DeeCrys * MaxDirections;
1211     }
1212     for (Int_t i = 0; i < MaxDeeCrys; i++) {
1213       for (Int_t j = 0; j < MaxDirections; j++) {
1214         fT2dDeeSC[i][j] = 0;
1215       }
1216     }
1217 
1218     fT2dSCEcha = new Int_t*[MaxDeeCrys];
1219     fCnew++;
1220     fT1dSCEcha = new Int_t[MaxDeeCrys * MaxDirections];
1221     fCnew++;
1222     for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxDeeCrys; i_DeeCrys++) {
1223       fT2dSCEcha[i_DeeCrys] = &fT1dSCEcha[0] + i_DeeCrys * MaxDirections;
1224     }
1225     for (Int_t i = 0; i < MaxDeeCrys; i++) {
1226       for (Int_t j = 0; j < MaxDirections; j++) {
1227         fT2dSCEcha[i][j] = 0;
1228       }
1229     }
1230 
1231     //................................ Build table

1232     Int_t MaxTyp = (Int_t)4;
1233     Int_t MaxCrysP1 = (Int_t)(fEcal->MaxCrysInSC() + 1);
1234 
1235     //      Int_t fT2d_jch_JY[type][Echa+1]:  JY in the SC as a function of type and Echa+1

1236     fT2d_jch_JY = new Int_t*[MaxTyp];
1237     fCnew++;
1238     fT1d_jch_JY = new Int_t[MaxTyp * MaxCrysP1];
1239     fCnew++;
1240     for (Int_t i_MaxTyp = 0; i_MaxTyp < MaxTyp; i_MaxTyp++) {
1241       fT2d_jch_JY[i_MaxTyp] = &fT1d_jch_JY[0] + i_MaxTyp * MaxCrysP1;
1242     }
1243 
1244     // type: 0=(top/right), 1=(top/left), 2=(bottom/left), 3=(bottom/right),

1245 
1246     //................. top/right

1247     for (Int_t k = 5; k >= 1; k--) {
1248       fT2d_jch_JY[0][k] = 4;
1249     }  //  k =  5, 4, 3, 2, 1 -> fT2d_jch_JY[0][k] = 4

1250     for (Int_t k = 10; k >= 6; k--) {
1251       fT2d_jch_JY[0][k] = 3;
1252     }  //  k = 10, 9, 8, 7, 6 -> fT2d_jch_JY[0][k] = 3

1253     for (Int_t k = 15; k >= 11; k--) {
1254       fT2d_jch_JY[0][k] = 2;
1255     }  //  k = 15,14,13,12,11 -> fT2d_jch_JY[0][k] = 2

1256     for (Int_t k = 20; k >= 16; k--) {
1257       fT2d_jch_JY[0][k] = 1;
1258     }  //  k = 20,19,18,17,16 -> fT2d_jch_JY[0][k] = 1

1259     for (Int_t k = 25; k >= 21; k--) {
1260       fT2d_jch_JY[0][k] = 0;
1261     }  //  k = 25,24,23,22,21 -> fT2d_jch_JY[0][k] = 0

1262     //................. top/left

1263     for (Int_t k = 5; k >= 1; k--) {
1264       fT2d_jch_JY[1][k] = 5 - k;
1265     }  //  k =  5, 4, 3, 2, 1 -> fT2d_jch_JY[1][k] = 0,1,2,3,4

1266     for (Int_t k = 10; k >= 6; k--) {
1267       fT2d_jch_JY[1][k] = 10 - k;
1268     }  //  k = 10, 9, 8, 7, 6 -> fT2d_jch_JY[1][k] = 0,1,2,3,4

1269     for (Int_t k = 15; k >= 11; k--) {
1270       fT2d_jch_JY[1][k] = 15 - k;
1271     }  //  k = 15,14,13,12,11 -> fT2d_jch_JY[1][k] = 0,1,2,3,4

1272     for (Int_t k = 20; k >= 16; k--) {
1273       fT2d_jch_JY[1][k] = 20 - k;
1274     }  //  k = 20,19,18,17,16 -> fT2d_jch_JY[1][k] = 0,1,2,3,4

1275     for (Int_t k = 25; k >= 21; k--) {
1276       fT2d_jch_JY[1][k] = 25 - k;
1277     }  //  k = 25,24,23,22,21 -> fT2d_jch_JY[1][k] = 0,1,2,3,4

1278     //................. bottom/left

1279     for (Int_t k = 1; k <= 5; k++) {
1280       fT2d_jch_JY[2][k] = 0;
1281     }  //  k =  1, 2, 3, 4, 5 -> fT2d_jch_JY[2][k] = 0

1282     for (Int_t k = 6; k <= 10; k++) {
1283       fT2d_jch_JY[2][k] = 1;
1284     }  //  k =  6, 7, 8, 9,10 -> fT2d_jch_JY[2][k] = 1

1285     for (Int_t k = 11; k <= 15; k++) {
1286       fT2d_jch_JY[2][k] = 2;
1287     }  //  k = 11,12,13,14,15 -> fT2d_jch_JY[2][k] = 2

1288     for (Int_t k = 16; k <= 20; k++) {
1289       fT2d_jch_JY[2][k] = 3;
1290     }  //  k = 16,17,18,19,20 -> fT2d_jch_JY[2][k] = 3

1291     for (Int_t k = 21; k <= 25; k++) {
1292       fT2d_jch_JY[2][k] = 4;
1293     }  //  k = 21,22,23,24,25 -> fT2d_jch_JY[2][k] = 4

1294     //................. bottom/right

1295     for (Int_t k = 1; k <= 5; k++) {
1296       fT2d_jch_JY[3][k] = k - 1;
1297     }  //  k =  1, 2, 3, 4, 5 -> fT2d_jch_JY[3][k] = 0,1,2,3,4

1298     for (Int_t k = 6; k <= 10; k++) {
1299       fT2d_jch_JY[3][k] = k - 6;
1300     }  //  k =  6, 7, 8, 9,10 -> fT2d_jch_JY[3][k] = 0,1,2,3,4

1301     for (Int_t k = 11; k <= 15; k++) {
1302       fT2d_jch_JY[3][k] = k - 11;
1303     }  //  k = 11,12,13,14,15 -> fT2d_jch_JY[3][k] = 0,1,2,3,4

1304     for (Int_t k = 16; k <= 20; k++) {
1305       fT2d_jch_JY[3][k] = k - 16;
1306     }  //  k = 16,17,18,19,20 -> fT2d_jch_JY[3][k] = 0,1,2,3,4

1307     for (Int_t k = 21; k <= 25; k++) {
1308       fT2d_jch_JY[3][k] = k - 21;
1309     }  //  k = 21,22,23,24,25 -> fT2d_jch_JY[3][k] = 0,1,2,3,4

1310 
1311     //      Int_t fT2d_ich_IX[type][Echa+1]:  IX in the SC as a function of type and Echa+1

1312     fT2d_ich_IX = new Int_t*[MaxTyp];
1313     fCnew++;
1314     fT1d_ich_IX = new Int_t[MaxTyp * MaxCrysP1];
1315     fCnew++;
1316     for (Int_t i_MaxTyp = 0; i_MaxTyp < MaxTyp; i_MaxTyp++) {
1317       fT2d_ich_IX[i_MaxTyp] = &fT1d_ich_IX[0] + i_MaxTyp * MaxCrysP1;
1318     }
1319 
1320     //................. top/right

1321     for (Int_t k = 5; k >= 1; k--) {
1322       fT2d_ich_IX[0][k] = 5 - k;
1323     }  //  k =  5, 4, 3, 2, 1 -> fT2d_ich_IX[0][k] = 0,1,2,3,4

1324     for (Int_t k = 10; k >= 6; k--) {
1325       fT2d_ich_IX[0][k] = 10 - k;
1326     }  //  k = 10, 9, 8, 7, 6 -> fT2d_ich_IX[0][k] = 0,1,2,3,4

1327     for (Int_t k = 15; k >= 11; k--) {
1328       fT2d_ich_IX[0][k] = 15 - k;
1329     }  //  k = 15,14,13,12,11 -> fT2d_ich_IX[0][k] = 0,1,2,3,4

1330     for (Int_t k = 20; k >= 16; k--) {
1331       fT2d_ich_IX[0][k] = 20 - k;
1332     }  //  k = 20,19,18,17,16 -> fT2d_ich_IX[0][k] = 0,1,2,3,4

1333     for (Int_t k = 25; k >= 21; k--) {
1334       fT2d_ich_IX[0][k] = 25 - k;
1335     }  //  k = 25,24,23,22,21 -> fT2d_ich_IX[0][k] = 0,1,2,3,4

1336     //................. top/left

1337     for (Int_t k = 5; k >= 1; k--) {
1338       fT2d_ich_IX[1][k] = 4;
1339     }  //  k =  5, 4, 3, 2, 1 -> fT2d_ich_IX[1][k] = 4

1340     for (Int_t k = 10; k >= 6; k--) {
1341       fT2d_ich_IX[1][k] = 3;
1342     }  //  k = 10, 9, 8, 7, 6 -> fT2d_ich_IX[1][k] = 3

1343     for (Int_t k = 15; k >= 11; k--) {
1344       fT2d_ich_IX[1][k] = 2;
1345     }  //  k = 15,14,13,12,11 -> fT2d_ich_IX[1][k] = 2

1346     for (Int_t k = 20; k >= 16; k--) {
1347       fT2d_ich_IX[1][k] = 1;
1348     }  //  k = 20,19,18,17,16 -> fT2d_ich_IX[1][k] = 1

1349     for (Int_t k = 25; k >= 21; k--) {
1350       fT2d_ich_IX[1][k] = 0;
1351     }  //  k = 25,24,23,22,21 -> fT2d_ich_IX[1][k] = 0

1352     //................. bottom/left

1353     for (Int_t k = 1; k <= 5; k++) {
1354       fT2d_ich_IX[2][k] = 5 - k;
1355     }  //  k =  1, 2, 3, 4, 5 -> fT2d_ich_IX[2][k] = 0,1,2,3,4

1356     for (Int_t k = 6; k <= 10; k++) {
1357       fT2d_ich_IX[2][k] = 10 - k;
1358     }  //  k =  6, 7, 8, 9,10 -> fT2d_ich_IX[2][k] = 0,1,2,3,4

1359     for (Int_t k = 11; k <= 15; k++) {
1360       fT2d_ich_IX[2][k] = 15 - k;
1361     }  //  k = 11,12,13,14,15 -> fT2d_ich_IX[2][k] = 0,1,2,3,4

1362     for (Int_t k = 16; k <= 20; k++) {
1363       fT2d_ich_IX[2][k] = 20 - k;
1364     }  //  k = 16,17,18,19,20 -> fT2d_ich_IX[2][k] = 0,1,2,3,4

1365     for (Int_t k = 21; k <= 25; k++) {
1366       fT2d_ich_IX[2][k] = 25 - k;
1367     }  //  k = 21,22,23,24,25 -> fT2d_ich_IX[2][k] = 0,1,2,3,4

1368     //................. bottom/right

1369     for (Int_t k = 1; k <= 5; k++) {
1370       fT2d_ich_IX[3][k] = 4;
1371     }  //  k =  1, 2, 3, 4, 5 -> fT2d_ich_IX[3][k] = 4

1372     for (Int_t k = 6; k <= 10; k++) {
1373       fT2d_ich_IX[3][k] = 3;
1374     }  //  k =  6, 7, 8, 9,10 -> fT2d_ich_IX[3][k] = 3

1375     for (Int_t k = 11; k <= 15; k++) {
1376       fT2d_ich_IX[3][k] = 2;
1377     }  //  k = 11,12,13,14,15 -> fT2d_ich_IX[3][k] = 2

1378     for (Int_t k = 16; k <= 20; k++) {
1379       fT2d_ich_IX[3][k] = 1;
1380     }  //  k = 16,17,18,19,20 -> fT2d_ich_IX[3][k] = 1

1381     for (Int_t k = 21; k <= 25; k++) {
1382       fT2d_ich_IX[3][k] = 0;
1383     }  //  k = 21,22,23,24,25 -> fT2d_ich_IX[3][k] = 0

1384 
1385     //............................................ type

1386     Int_t Nb_DeeSC_JY = fEcal->MaxSCIYInDee();
1387     Int_t** type = new Int_t*[Nb_DeeSC_JY];
1388     fCnew++;
1389     Int_t* type_d1 = new Int_t[Nb_DeeSC_JY * MaxDirections];
1390     fCnew++;
1391     for (Int_t i_DeeSC_JY = 0; i_DeeSC_JY < Nb_DeeSC_JY; i_DeeSC_JY++) {
1392       type[i_DeeSC_JY] = &type_d1[0] + i_DeeSC_JY * MaxDirections;
1393     }
1394 
1395     //  bottom = (0,9), top  = (10,19)

1396     //  right  = 0    , left = 1

1397     //  type = Quadrant number - 1

1398 
1399     type[10][0] = 0;  // Q1 top right

1400     type[11][0] = 0;
1401     type[12][0] = 0;
1402     type[13][0] = 0;
1403     type[14][0] = 0;
1404     type[15][0] = 0;
1405     type[16][0] = 0;
1406     type[17][0] = 0;
1407     type[18][0] = 0;
1408     type[19][0] = 0;
1409 
1410     type[10][1] = 1;  // Q2 top left

1411     type[11][1] = 1;
1412     type[12][1] = 1;
1413     type[13][1] = 1;
1414     type[14][1] = 1;
1415     type[15][1] = 1;
1416     type[16][1] = 1;
1417     type[17][1] = 1;
1418     type[18][1] = 1;
1419     type[19][1] = 1;
1420 
1421     type[0][1] = 2;  // Q3 : bottom left

1422     type[1][1] = 2;
1423     type[2][1] = 2;
1424     type[3][1] = 2;
1425     type[4][1] = 2;
1426     type[5][1] = 2;
1427     type[6][1] = 2;
1428     type[7][1] = 2;
1429     type[8][1] = 2;
1430     type[9][1] = 2;
1431 
1432     type[0][0] = 3;  // Q4 : bottom right

1433     type[1][0] = 3;
1434     type[2][0] = 3;
1435     type[3][0] = 3;
1436     type[4][0] = 3;
1437     type[5][0] = 3;
1438     type[6][0] = 3;
1439     type[7][0] = 3;
1440     type[8][0] = 3;
1441     type[9][0] = 3;
1442 
1443     Int_t Nb_SCCrys_IX = fEcal->MaxCrysIXInSC();
1444     Int_t Nb_SCCrys_JY = fEcal->MaxCrysIYInSC();
1445 
1446     for (Int_t kSC = 0; kSC < MaxDeeSC; kSC++)  //  kSC  = 0 to 199   (MaxDeeSC = 200)

1447     {
1448       for (Int_t n_Echa = 1; n_Echa <= MaxSCEcha; n_Echa++)  //  n_Echa   = 1 to 25    (MaxSCEcha = 25)

1449       {
1450         for (Int_t idir = 0; idir < 2; idir++) {
1451           Int_t ikSC = kSC / Nb_DeeSC_JY;  //  ikSC = 0 to 9

1452           Int_t jkSC = kSC % Nb_DeeSC_JY;  //  jkSC = 0,1,2,..,19

1453 
1454           Int_t icrys = ikSC * Nb_SCCrys_IX + fT2d_ich_IX[type[jkSC][idir]][n_Echa];
1455           //  type[0->9][1->2] = 0,1,2,3

1456           //  fT2d_ich_IX[0->3][1->25] = 0,1,2,3,4

1457           //  icrys = 0 to 49  (=> IX)

1458 
1459           Int_t jcrys = jkSC * Nb_SCCrys_JY + fT2d_jch_JY[type[jkSC][idir]][n_Echa];
1460           //  type[0->9][1->2] = 0,1,2,3

1461           //  fT2d_jch_JY[0->3][1->25] = 0,1,2,3,4

1462           //  jcrys = 0 to 99  (=> IY)

1463 
1464           Int_t n_DeeCrys = icrys * Nb_DeeSC_JY * Nb_SCCrys_JY + jcrys + 1;  //  n_DeeCrys = 1 to 5000

1465 
1466           fT3dDeeCrys[kSC][n_Echa - 1][idir] = n_DeeCrys;  // fT3dDeeCrys[][][] : range = [1,5000]

1467           fT2dDeeSC[n_DeeCrys - 1][idir] = kSC + 1;        // fT2dDeeSC[][]     : range = [1,200]

1468           fT2dSCEcha[n_DeeCrys - 1][idir] = n_Echa;        // fT2dSCEcha[][]    : range = [1,25]

1469         }
1470       }
1471     }
1472     // std::cout << "#TEcnaNumbering::TBuildEndcapCrysTable()> Crys Table Building done" << std::endl;

1473 
1474     delete[] type;
1475     fCdelete++;
1476     delete[] type_d1;
1477     fCdelete++;
1478   } else {
1479     // std::cout << "#TEcnaNumbering::TBuildEndcapCrysTable()> No Building of Crys Table since it is already done " << std::endl;

1480   }
1481 }
1482 
1483 void TEcnaNumbering::BuildEndcapSCTable() {
1484   // Build the correspondance table: (Dee, DeeSC) <-> (DS , DSSC) for the ECAL ENDCAP

1485   //

1486   //  From CMS Internal Note  "CMS ECAL Endcap channel numbering"

1487   //

1488   //       Name       Number                     Reference set    Range     Comment

1489   //

1490   //       Dee      : Dee number                                  [1,4]

1491   //       DeeSC    : Super-Crystal (SC) number  in Dee           [1,200]   (IY,IX) progression

1492   //       DS       : Data Sector number         in EE + or -     [1,9]     (IY,IX) progression

1493   //       DSSC     : Super-Crystal (SC) number  in Data Sector   [1,32]    Crystal numbering in data sector

1494   //       DeeSCCons: Super-Crystal (SC) number  for construction [1,297]

1495   //

1496   //   fill the 2d arrays:  fT2d_DS[4][200], fT2d_DSSC[4][200] and fT2d_DeeSCCons[4][200]

1497   //

1498   //------------------------------------------------------------------------------------------------------

1499 
1500   //.................. Allocations and Init

1501   Int_t MaxEEDee = fEcal->MaxDeeInEE();
1502   Int_t MaxDeeSC = fEcal->MaxSCEcnaInDee();
1503   Int_t MaxEESCForCons = 2 * fEcal->MaxSCForConsInDee();
1504 
1505   fT2d_DS = new Int_t*[MaxEEDee];
1506   fCnew++;  // = DS[Dee - 1, CNA_SCInDee - 1]

1507   fT1d_DS = new Int_t[MaxEEDee * MaxDeeSC];
1508   fCnew++;
1509   for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxEEDee; i_DeeCrys++) {
1510     fT2d_DS[i_DeeCrys] = &fT1d_DS[0] + i_DeeCrys * MaxDeeSC;
1511   }
1512   for (Int_t i = 0; i < MaxEEDee; i++) {
1513     for (Int_t j = 0; j < MaxDeeSC; j++) {
1514       fT2d_DS[i][j] = 0;
1515     }
1516   }
1517 
1518   fT2d_DSSC = new Int_t*[MaxEEDee];
1519   fCnew++;  // = SCInDS[Dee - 1, CNA_SCInDee - 1]

1520   fT1d_DSSC = new Int_t[MaxEEDee * MaxDeeSC];
1521   fCnew++;
1522   for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxEEDee; i_DeeCrys++) {
1523     fT2d_DSSC[i_DeeCrys] = &fT1d_DSSC[0] + i_DeeCrys * MaxDeeSC;
1524   }
1525   for (Int_t i = 0; i < MaxEEDee; i++) {
1526     for (Int_t j = 0; j < MaxDeeSC; j++) {
1527       fT2d_DSSC[i][j] = 0;
1528     }
1529   }
1530 
1531   fT2d_DeeSCCons = new Int_t*[MaxEEDee];
1532   fCnew++;  // = SCConsInDee[Dee - 1, CNA_SCInDee - 1]

1533   fT1d_DeeSCCons = new Int_t[MaxEEDee * MaxDeeSC];
1534   fCnew++;
1535   for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxEEDee; i_DeeCrys++) {
1536     fT2d_DeeSCCons[i_DeeCrys] = &fT1d_DeeSCCons[0] + i_DeeCrys * MaxDeeSC;
1537   }
1538   for (Int_t i = 0; i < MaxEEDee; i++) {
1539     for (Int_t j = 0; j < MaxDeeSC; j++) {
1540       fT2d_DeeSCCons[i][j] = 0;
1541     }
1542   }
1543 
1544   fT2d_RecovDeeSC = new Int_t*[MaxEEDee];
1545   fCnew++;  // = CNA_SCInDee[Dee - 1, SCConsInDee - 1]

1546   fT1d_RecovDeeSC = new Int_t[MaxEEDee * MaxEESCForCons];
1547   fCnew++;
1548   for (Int_t i_DeeCrys = 0; i_DeeCrys < MaxEEDee; i_DeeCrys++) {
1549     fT2d_RecovDeeSC[i_DeeCrys] = &fT1d_RecovDeeSC[0] + i_DeeCrys * MaxEESCForCons;
1550   }
1551   for (Int_t i = 0; i < MaxEEDee; i++) {
1552     for (Int_t j = 0; j < MaxEESCForCons; j++) {
1553       fT2d_RecovDeeSC[i][j] = 0;
1554     }
1555   }
1556 
1557   //.............................. Data sector (DS) numbers: fT2d_DS[][]

1558   //                               isc = ECNA numbers

1559 
1560   Int_t ids = 0;
1561 
1562   //........... (D1,S1)=(D2,S9)=(D3,S9)=(D4,S1)

1563   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1564     if (dee == 1 || dee == 4) {
1565       ids = 1;
1566     }
1567     if (dee == 2 || dee == 3) {
1568       ids = 9;
1569     }
1570     for (Int_t isc = 13; isc <= 20; isc++) [[clang::suppress]]
1571       fT2d_DS[dee - 1][isc - 1] = ids;
1572     for (Int_t isc = 33; isc <= 40; isc++)
1573       fT2d_DS[dee - 1][isc - 1] = ids;
1574     for (Int_t isc = 54; isc <= 60; isc++)
1575       fT2d_DS[dee - 1][isc - 1] = ids;
1576     for (Int_t isc = 75; isc <= 79; isc++)
1577       fT2d_DS[dee - 1][isc - 1] = ids;
1578     for (Int_t isc = 96; isc <= 99; isc++)
1579       fT2d_DS[dee - 1][isc - 1] = ids;
1580     for (Int_t isc = 118; isc <= 119; isc++)
1581       fT2d_DS[dee - 1][isc - 1] = ids;
1582   }
1583   //........... (D1,S2)=(D2,S8)=(D3,S8)=(D4,S2)

1584   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1585     if (dee == 1 || dee == 4) {
1586       ids = 2;
1587     }
1588     if (dee == 2 || dee == 3) {
1589       ids = 8;
1590     }
1591     for (Int_t isc = 32; isc <= 32; isc++)
1592       fT2d_DS[dee - 1][isc - 1] = ids;
1593     for (Int_t isc = 51; isc <= 53; isc++)
1594       fT2d_DS[dee - 1][isc - 1] = ids;
1595     for (Int_t isc = 72; isc <= 74; isc++)
1596       fT2d_DS[dee - 1][isc - 1] = ids;
1597     for (Int_t isc = 92; isc <= 95; isc++)
1598       fT2d_DS[dee - 1][isc - 1] = ids;
1599     for (Int_t isc = 112; isc <= 117; isc++)
1600       fT2d_DS[dee - 1][isc - 1] = ids;
1601     for (Int_t isc = 132; isc <= 138; isc++)
1602       fT2d_DS[dee - 1][isc - 1] = ids;
1603     for (Int_t isc = 152; isc <= 157; isc++)
1604       fT2d_DS[dee - 1][isc - 1] = ids;
1605     for (Int_t isc = 173; isc <= 176; isc++)
1606       fT2d_DS[dee - 1][isc - 1] = ids;
1607     for (Int_t isc = 193; isc <= 193; isc++)
1608       fT2d_DS[dee - 1][isc - 1] = ids;
1609   }
1610   //........... (D1,S3)=(D2,S7)=(D3,S7)=(D4,S3)

1611   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1612     if (dee == 1 || dee == 4) {
1613       ids = 3;
1614     }
1615     if (dee == 2 || dee == 3) {
1616       ids = 7;
1617     }
1618     for (Int_t isc = 50; isc <= 50; isc++)
1619       fT2d_DS[dee - 1][isc - 1] = ids;
1620     for (Int_t isc = 69; isc <= 71; isc++)
1621       fT2d_DS[dee - 1][isc - 1] = ids;
1622     for (Int_t isc = 88; isc <= 91; isc++)
1623       fT2d_DS[dee - 1][isc - 1] = ids;
1624     for (Int_t isc = 108; isc <= 111; isc++)
1625       fT2d_DS[dee - 1][isc - 1] = ids;
1626     for (Int_t isc = 127; isc <= 131; isc++)
1627       fT2d_DS[dee - 1][isc - 1] = ids;
1628     for (Int_t isc = 147; isc <= 151; isc++)
1629       fT2d_DS[dee - 1][isc - 1] = ids;
1630     for (Int_t isc = 166; isc <= 172; isc++)
1631       fT2d_DS[dee - 1][isc - 1] = ids;
1632     for (Int_t isc = 188; isc <= 192; isc++)
1633       fT2d_DS[dee - 1][isc - 1] = ids;
1634   }
1635   //........... (D1,S4)=(D2,S6)=(D3,S6)=(D4,S4)

1636   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1637     if (dee == 1 || dee == 4) {
1638       ids = 4;
1639     }
1640     if (dee == 2 || dee == 3) {
1641       ids = 6;
1642     }
1643     for (Int_t isc = 27; isc <= 29; isc++)
1644       fT2d_DS[dee - 1][isc - 1] = ids;
1645     for (Int_t isc = 44; isc <= 49; isc++)
1646       fT2d_DS[dee - 1][isc - 1] = ids;
1647     for (Int_t isc = 62; isc <= 68; isc++)
1648       fT2d_DS[dee - 1][isc - 1] = ids;
1649     for (Int_t isc = 82; isc <= 87; isc++)
1650       fT2d_DS[dee - 1][isc - 1] = ids;
1651     for (Int_t isc = 102; isc <= 107; isc++)
1652       fT2d_DS[dee - 1][isc - 1] = ids;
1653     for (Int_t isc = 123; isc <= 126; isc++)
1654       fT2d_DS[dee - 1][isc - 1] = ids;
1655     for (Int_t isc = 144; isc <= 146; isc++)
1656       fT2d_DS[dee - 1][isc - 1] = ids;
1657     for (Int_t isc = 165; isc <= 165; isc++)
1658       fT2d_DS[dee - 1][isc - 1] = ids;
1659   }
1660   //........... (D1,S5)=(D2,S5)=(D3,S5)=(D4,S5)

1661   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1662     for (Int_t isc = 1; isc <= 8; isc++)
1663       fT2d_DS[dee - 1][isc - 1] = 5;
1664     for (Int_t isc = 21; isc <= 26; isc++)
1665       fT2d_DS[dee - 1][isc - 1] = 5;
1666     for (Int_t isc = 41; isc <= 43; isc++)
1667       fT2d_DS[dee - 1][isc - 1] = 5;
1668   }
1669 
1670   //............................................ SC numbers in Data Sectors: fT2d_DSSC[][]

1671   //                                             fT2d_DSSC[dee-1][SC ECNA number - 1] = SC number in DS;

1672   for (Int_t dee = 1; dee <= 4; dee++) {
1673     for (Int_t isc = 1; isc <= MaxDeeSC; isc++) {
1674       fT2d_DSSC[dee - 1][isc - 1] = -1;
1675     }
1676     //.......................................... S1 (D1,D4), S9 (D2,D3) ; 33 SC's

1677     fT2d_DSSC[dee - 1][13 - 1] = 12;
1678     fT2d_DSSC[dee - 1][14 - 1] = 11;
1679     fT2d_DSSC[dee - 1][15 - 1] = 10;
1680     fT2d_DSSC[dee - 1][16 - 1] = 9;
1681     fT2d_DSSC[dee - 1][17 - 1] = 4;
1682     fT2d_DSSC[dee - 1][18 - 1] = 3;
1683     fT2d_DSSC[dee - 1][19 - 1] = 2;
1684     fT2d_DSSC[dee - 1][20 - 1] = 1;
1685 
1686     fT2d_DSSC[dee - 1][33 - 1] = 16;
1687     fT2d_DSSC[dee - 1][34 - 1] = 15;
1688     fT2d_DSSC[dee - 1][35 - 1] = 14;
1689     fT2d_DSSC[dee - 1][36 - 1] = 13;
1690     fT2d_DSSC[dee - 1][37 - 1] = 8;
1691     fT2d_DSSC[dee - 1][38 - 1] = 7;
1692     fT2d_DSSC[dee - 1][39 - 1] = 6;
1693     fT2d_DSSC[dee - 1][40 - 1] = 5;
1694 
1695     fT2d_DSSC[dee - 1][54 - 1] = 33;
1696     fT2d_DSSC[dee - 1][55 - 1] = 31;
1697     fT2d_DSSC[dee - 1][56 - 1] = 27;
1698     fT2d_DSSC[dee - 1][57 - 1] = 24;
1699     fT2d_DSSC[dee - 1][58 - 1] = 20;
1700     fT2d_DSSC[dee - 1][59 - 1] = 17;
1701     fT2d_DSSC[dee - 1][60 - 1] = 30;  //  (182a, 33a for construction)

1702 
1703     fT2d_DSSC[dee - 1][75 - 1] = 32;
1704     fT2d_DSSC[dee - 1][76 - 1] = 28;
1705     fT2d_DSSC[dee - 1][77 - 1] = 25;
1706     fT2d_DSSC[dee - 1][78 - 1] = 21;
1707     fT2d_DSSC[dee - 1][79 - 1] = 18;
1708 
1709     fT2d_DSSC[dee - 1][96 - 1] = 29;
1710     fT2d_DSSC[dee - 1][97 - 1] = 26;
1711     fT2d_DSSC[dee - 1][98 - 1] = 22;
1712     fT2d_DSSC[dee - 1][99 - 1] = 19;
1713 
1714     fT2d_DSSC[dee - 1][118 - 1] = 23;
1715     fT2d_DSSC[dee - 1][119 - 1] = 30;  //  (182b, 33b for construction)

1716 
1717     //.......................................... S2 (D1,D4), S8(D2,D3) ; 32 SC's

1718     fT2d_DSSC[dee - 1][32 - 1] = 25;  // also 3;  // ( (207c, 58c) also (178c, 29c) for construction)

1719 
1720     fT2d_DSSC[dee - 1][51 - 1] = 32;
1721     fT2d_DSSC[dee - 1][52 - 1] = 26;
1722     fT2d_DSSC[dee - 1][53 - 1] = 18;
1723 
1724     fT2d_DSSC[dee - 1][72 - 1] = 27;
1725     fT2d_DSSC[dee - 1][73 - 1] = 19;
1726     fT2d_DSSC[dee - 1][74 - 1] = 12;
1727 
1728     fT2d_DSSC[dee - 1][92 - 1] = 28;
1729     fT2d_DSSC[dee - 1][93 - 1] = 20;
1730     fT2d_DSSC[dee - 1][94 - 1] = 13;
1731     fT2d_DSSC[dee - 1][95 - 1] = 7;
1732 
1733     fT2d_DSSC[dee - 1][112 - 1] = 29;
1734     fT2d_DSSC[dee - 1][113 - 1] = 21;
1735     fT2d_DSSC[dee - 1][114 - 1] = 14;
1736     fT2d_DSSC[dee - 1][115 - 1] = 8;
1737     fT2d_DSSC[dee - 1][116 - 1] = 4;
1738     fT2d_DSSC[dee - 1][117 - 1] = 1;
1739 
1740     fT2d_DSSC[dee - 1][132 - 1] = 30;
1741     fT2d_DSSC[dee - 1][133 - 1] = 22;
1742     fT2d_DSSC[dee - 1][134 - 1] = 15;
1743     fT2d_DSSC[dee - 1][135 - 1] = 9;
1744     fT2d_DSSC[dee - 1][136 - 1] = 5;
1745     fT2d_DSSC[dee - 1][137 - 1] = 2;
1746     fT2d_DSSC[dee - 1][138 - 1] = 3;  // (178a, 29a for construction)

1747 
1748     fT2d_DSSC[dee - 1][152 - 1] = 31;
1749     fT2d_DSSC[dee - 1][153 - 1] = 23;
1750     fT2d_DSSC[dee - 1][154 - 1] = 16;
1751     fT2d_DSSC[dee - 1][155 - 1] = 10;
1752     fT2d_DSSC[dee - 1][156 - 1] = 6;
1753     fT2d_DSSC[dee - 1][157 - 1] = 3;  // (178b, 29b for construction)

1754 
1755     fT2d_DSSC[dee - 1][173 - 1] = 24;
1756     fT2d_DSSC[dee - 1][174 - 1] = 17;
1757     fT2d_DSSC[dee - 1][175 - 1] = 11;
1758     fT2d_DSSC[dee - 1][176 - 1] = 25;  // (207a, 58a for construction)

1759 
1760     fT2d_DSSC[dee - 1][193 - 1] = 25;  // (207b, 58b for construction)

1761 
1762     //.......................................... S3 (D1,D4), S7 (D2,D3)  ; 34 SC's

1763     fT2d_DSSC[dee - 1][50 - 1] = 10;
1764 
1765     fT2d_DSSC[dee - 1][69 - 1] = 18;
1766     fT2d_DSSC[dee - 1][70 - 1] = 11;
1767     fT2d_DSSC[dee - 1][71 - 1] = 3;
1768 
1769     fT2d_DSSC[dee - 1][88 - 1] = 25;
1770     fT2d_DSSC[dee - 1][89 - 1] = 19;
1771     fT2d_DSSC[dee - 1][90 - 1] = 12;
1772     fT2d_DSSC[dee - 1][91 - 1] = 4;
1773 
1774     fT2d_DSSC[dee - 1][108 - 1] = 26;
1775     fT2d_DSSC[dee - 1][109 - 1] = 20;
1776     fT2d_DSSC[dee - 1][110 - 1] = 13;
1777     fT2d_DSSC[dee - 1][111 - 1] = 5;
1778 
1779     fT2d_DSSC[dee - 1][127 - 1] = 31;
1780     fT2d_DSSC[dee - 1][128 - 1] = 27;
1781     fT2d_DSSC[dee - 1][129 - 1] = 21;
1782     fT2d_DSSC[dee - 1][130 - 1] = 14;
1783     fT2d_DSSC[dee - 1][131 - 1] = 6;
1784 
1785     fT2d_DSSC[dee - 1][147 - 1] = 32;
1786     fT2d_DSSC[dee - 1][148 - 1] = 28;
1787     fT2d_DSSC[dee - 1][149 - 1] = 22;
1788     fT2d_DSSC[dee - 1][150 - 1] = 15;
1789     fT2d_DSSC[dee - 1][151 - 1] = 7;
1790 
1791     fT2d_DSSC[dee - 1][166 - 1] = 33;
1792     fT2d_DSSC[dee - 1][167 - 1] = 30;
1793     fT2d_DSSC[dee - 1][168 - 1] = 29;
1794     fT2d_DSSC[dee - 1][169 - 1] = 23;
1795     fT2d_DSSC[dee - 1][170 - 1] = 16;
1796     fT2d_DSSC[dee - 1][171 - 1] = 8;
1797     fT2d_DSSC[dee - 1][172 - 1] = 1;
1798 
1799     fT2d_DSSC[dee - 1][188 - 1] = 34;  // (298a, 149a for construction)

1800     fT2d_DSSC[dee - 1][189 - 1] = 24;
1801     fT2d_DSSC[dee - 1][190 - 1] = 17;
1802     fT2d_DSSC[dee - 1][191 - 1] = 9;
1803     fT2d_DSSC[dee - 1][192 - 1] = 2;
1804 
1805     //.......................................... S4 (D1,D4), S6 (D2,D3) ; 33 SC's

1806     fT2d_DSSC[dee - 1][27 - 1] = 33;
1807     fT2d_DSSC[dee - 1][28 - 1] = 32;
1808     fT2d_DSSC[dee - 1][29 - 1] = 14;  // also 21;  //  ( (261a, 112a) also (268a, 119a) for construction)

1809 
1810     fT2d_DSSC[dee - 1][44 - 1] = 22;
1811     fT2d_DSSC[dee - 1][45 - 1] = 15;
1812     fT2d_DSSC[dee - 1][46 - 1] = 8;
1813     fT2d_DSSC[dee - 1][47 - 1] = 4;
1814     fT2d_DSSC[dee - 1][48 - 1] = 2;
1815     fT2d_DSSC[dee - 1][49 - 1] = 1;
1816 
1817     fT2d_DSSC[dee - 1][62 - 1] = 29;
1818     fT2d_DSSC[dee - 1][63 - 1] = 28;
1819     fT2d_DSSC[dee - 1][64 - 1] = 23;
1820     fT2d_DSSC[dee - 1][65 - 1] = 16;
1821     fT2d_DSSC[dee - 1][66 - 1] = 9;
1822     fT2d_DSSC[dee - 1][67 - 1] = 5;
1823     fT2d_DSSC[dee - 1][68 - 1] = 3;
1824 
1825     fT2d_DSSC[dee - 1][82 - 1] = 31;
1826     fT2d_DSSC[dee - 1][83 - 1] = 30;
1827     fT2d_DSSC[dee - 1][84 - 1] = 24;
1828     fT2d_DSSC[dee - 1][85 - 1] = 17;
1829     fT2d_DSSC[dee - 1][86 - 1] = 10;
1830     fT2d_DSSC[dee - 1][87 - 1] = 6;
1831 
1832     fT2d_DSSC[dee - 1][102 - 1] = 21;  //   (268c, 119c for construction)

1833     fT2d_DSSC[dee - 1][103 - 1] = 27;
1834     fT2d_DSSC[dee - 1][104 - 1] = 25;
1835     fT2d_DSSC[dee - 1][105 - 1] = 18;
1836     fT2d_DSSC[dee - 1][106 - 1] = 11;
1837     fT2d_DSSC[dee - 1][107 - 1] = 7;
1838 
1839     fT2d_DSSC[dee - 1][123 - 1] = 21;  //   (268b, 119b for construction)

1840     fT2d_DSSC[dee - 1][124 - 1] = 26;
1841     fT2d_DSSC[dee - 1][125 - 1] = 19;
1842     fT2d_DSSC[dee - 1][126 - 1] = 12;
1843 
1844     fT2d_DSSC[dee - 1][144 - 1] = 14;  //   (261c, 112c for construction)

1845     fT2d_DSSC[dee - 1][145 - 1] = 20;
1846     fT2d_DSSC[dee - 1][146 - 1] = 13;
1847 
1848     fT2d_DSSC[dee - 1][165 - 1] = 14;  //   (261b, 112b for construction)

1849 
1850     //.......................................... S5   (2*17 = 34 SC's)

1851     if (dee == 1 || dee == 3) {
1852       fT2d_DSSC[dee - 1][1 - 1] = 34;
1853       fT2d_DSSC[dee - 1][2 - 1] = 33;
1854       fT2d_DSSC[dee - 1][3 - 1] = 32;
1855       fT2d_DSSC[dee - 1][4 - 1] = 31;
1856       fT2d_DSSC[dee - 1][5 - 1] = 26;
1857       fT2d_DSSC[dee - 1][6 - 1] = 25;
1858       fT2d_DSSC[dee - 1][7 - 1] = 24;
1859       fT2d_DSSC[dee - 1][8 - 1] = 23;
1860 
1861       fT2d_DSSC[dee - 1][21 - 1] = 30;
1862       fT2d_DSSC[dee - 1][22 - 1] = 29;
1863       fT2d_DSSC[dee - 1][23 - 1] = 28;
1864       fT2d_DSSC[dee - 1][24 - 1] = 27;
1865       fT2d_DSSC[dee - 1][25 - 1] = 22;
1866       fT2d_DSSC[dee - 1][26 - 1] = 21;
1867 
1868       fT2d_DSSC[dee - 1][41 - 1] = 20;  //   (281a for construction)

1869       fT2d_DSSC[dee - 1][42 - 1] = 19;
1870       fT2d_DSSC[dee - 1][43 - 1] = 18;
1871     }
1872 
1873     if (dee == 2 || dee == 4) {
1874       fT2d_DSSC[dee - 1][1 - 1] = 17;
1875       fT2d_DSSC[dee - 1][2 - 1] = 16;
1876       fT2d_DSSC[dee - 1][3 - 1] = 15;
1877       fT2d_DSSC[dee - 1][4 - 1] = 14;
1878       fT2d_DSSC[dee - 1][5 - 1] = 9;
1879       fT2d_DSSC[dee - 1][6 - 1] = 8;
1880       fT2d_DSSC[dee - 1][7 - 1] = 7;
1881       fT2d_DSSC[dee - 1][8 - 1] = 6;
1882 
1883       fT2d_DSSC[dee - 1][21 - 1] = 13;
1884       fT2d_DSSC[dee - 1][22 - 1] = 12;
1885       fT2d_DSSC[dee - 1][23 - 1] = 11;
1886       fT2d_DSSC[dee - 1][24 - 1] = 10;
1887       fT2d_DSSC[dee - 1][25 - 1] = 5;
1888       fT2d_DSSC[dee - 1][26 - 1] = 4;
1889 
1890       fT2d_DSSC[dee - 1][41 - 1] = 3;  //   (132a for construction)

1891       fT2d_DSSC[dee - 1][42 - 1] = 2;
1892       fT2d_DSSC[dee - 1][43 - 1] = 1;
1893     }
1894   }
1895   //............................... Numbers for construction: fT2d_DeeSCCons[][]

1896   //                                fT2d_DeeSCCons[dee-1][SC ECNA number - 1] = SC number for construction;

1897 
1898   //............................... init to -1

1899   for (Int_t dee = 1; dee <= MaxEEDee; dee++) {
1900     for (Int_t isc = 1; isc <= MaxDeeSC; isc++) {
1901       fT2d_DeeSCCons[dee - 1][isc - 1] = -1;
1902     }
1903   }
1904 
1905   for (Int_t i_dee_type = 1; i_dee_type <= 2; i_dee_type++) {
1906     Int_t dee = -1;
1907     if (i_dee_type == 1) {
1908       dee = 1;
1909     }
1910     if (i_dee_type == 2) {
1911       dee = 3;
1912     }
1913 
1914     //.......................................... (D1,S1 or D3,S9) AND (D2,S9 or D4,S1)

1915     //  number in comment = fT2d_DSSC[dee-1][SC ECNA nb - 1] (= SC number in DS)

1916     fT2d_DeeSCCons[dee - 1][13 - 1] = 161;
1917     fT2d_DeeSCCons[dee][13 - 1] = 12;  // 12;

1918     fT2d_DeeSCCons[dee - 1][14 - 1] = 160;
1919     fT2d_DeeSCCons[dee][14 - 1] = 11;  // 11;

1920     fT2d_DeeSCCons[dee - 1][15 - 1] = 159;
1921     fT2d_DeeSCCons[dee][15 - 1] = 10;  // 10;

1922     fT2d_DeeSCCons[dee - 1][16 - 1] = 158;
1923     fT2d_DeeSCCons[dee][16 - 1] = 9;  //  9;

1924     fT2d_DeeSCCons[dee - 1][17 - 1] = 153;
1925     fT2d_DeeSCCons[dee][17 - 1] = 4;  //  4;

1926     fT2d_DeeSCCons[dee - 1][18 - 1] = 152;
1927     fT2d_DeeSCCons[dee][18 - 1] = 3;  //  3;

1928     fT2d_DeeSCCons[dee - 1][19 - 1] = 151;
1929     fT2d_DeeSCCons[dee][19 - 1] = 2;  //  2;

1930     fT2d_DeeSCCons[dee - 1][20 - 1] = 150;
1931     fT2d_DeeSCCons[dee][20 - 1] = 1;  //  1;

1932 
1933     fT2d_DeeSCCons[dee - 1][33 - 1] = 165;
1934     fT2d_DeeSCCons[dee][33 - 1] = 16;  // 16;

1935     fT2d_DeeSCCons[dee - 1][34 - 1] = 164;
1936     fT2d_DeeSCCons[dee][34 - 1] = 15;  // 15;

1937     fT2d_DeeSCCons[dee - 1][35 - 1] = 163;
1938     fT2d_DeeSCCons[dee][35 - 1] = 14;  // 14;

1939     fT2d_DeeSCCons[dee - 1][36 - 1] = 162;
1940     fT2d_DeeSCCons[dee][36 - 1] = 13;  // 13;

1941     fT2d_DeeSCCons[dee - 1][37 - 1] = 157;
1942     fT2d_DeeSCCons[dee][37 - 1] = 8;  //  8;

1943     fT2d_DeeSCCons[dee - 1][38 - 1] = 156;
1944     fT2d_DeeSCCons[dee][38 - 1] = 7;  //  7;

1945     fT2d_DeeSCCons[dee - 1][39 - 1] = 155;
1946     fT2d_DeeSCCons[dee][39 - 1] = 6;  //  6;

1947     fT2d_DeeSCCons[dee - 1][40 - 1] = 154;
1948     fT2d_DeeSCCons[dee][40 - 1] = 5;  //  5;

1949 
1950     fT2d_DeeSCCons[dee - 1][54 - 1] = 193;
1951     fT2d_DeeSCCons[dee][54 - 1] = 44;  // 33;

1952     fT2d_DeeSCCons[dee - 1][55 - 1] = 186;
1953     fT2d_DeeSCCons[dee][55 - 1] = 37;  // 31;

1954     fT2d_DeeSCCons[dee - 1][56 - 1] = 179;
1955     fT2d_DeeSCCons[dee][56 - 1] = 30;  // 27;

1956     fT2d_DeeSCCons[dee - 1][57 - 1] = 173;
1957     fT2d_DeeSCCons[dee][57 - 1] = 24;  // 24;

1958     fT2d_DeeSCCons[dee - 1][58 - 1] = 169;
1959     fT2d_DeeSCCons[dee][58 - 1] = 20;  // 20;

1960     fT2d_DeeSCCons[dee - 1][59 - 1] = 166;
1961     fT2d_DeeSCCons[dee][59 - 1] = 17;  // 17;

1962     fT2d_DeeSCCons[dee - 1][60 - 1] = 182;
1963     fT2d_DeeSCCons[dee][60 - 1] = 33;  // 30;    // 182a ;  33a;

1964 
1965     fT2d_DeeSCCons[dee - 1][75 - 1] = 187;
1966     fT2d_DeeSCCons[dee][75 - 1] = 38;  // 32;

1967     fT2d_DeeSCCons[dee - 1][76 - 1] = 180;
1968     fT2d_DeeSCCons[dee][76 - 1] = 31;  // 28;

1969     fT2d_DeeSCCons[dee - 1][77 - 1] = 174;
1970     fT2d_DeeSCCons[dee][77 - 1] = 25;  // 25;

1971     fT2d_DeeSCCons[dee - 1][78 - 1] = 170;
1972     fT2d_DeeSCCons[dee][78 - 1] = 21;  // 21;

1973     fT2d_DeeSCCons[dee - 1][79 - 1] = 167;
1974     fT2d_DeeSCCons[dee][79 - 1] = 18;  // 18;

1975 
1976     fT2d_DeeSCCons[dee - 1][96 - 1] = 181;
1977     fT2d_DeeSCCons[dee][96 - 1] = 32;  // 29;

1978     fT2d_DeeSCCons[dee - 1][97 - 1] = 175;
1979     fT2d_DeeSCCons[dee][97 - 1] = 26;  // 26;

1980     fT2d_DeeSCCons[dee - 1][98 - 1] = 171;
1981     fT2d_DeeSCCons[dee][98 - 1] = 22;  // 22;

1982     fT2d_DeeSCCons[dee - 1][99 - 1] = 168;
1983     fT2d_DeeSCCons[dee][99 - 1] = 19;  // 19;

1984 
1985     fT2d_DeeSCCons[dee - 1][118 - 1] = 172;
1986     fT2d_DeeSCCons[dee][118 - 1] = 23;  // 23;

1987     fT2d_DeeSCCons[dee - 1][119 - 1] = 182;
1988     fT2d_DeeSCCons[dee][119 - 1] = 33;  // 30;    // 182b ;  33b;

1989 
1990     //.......................................... (D1,S2 or D3,S8) AND (D2,S8 or D4,S2)

1991     fT2d_DeeSCCons[dee - 1][32 - 1] = 178;
1992     fT2d_DeeSCCons[dee][32 - 1] = 29;  // 25-3;   // 178c-207c ; 29c-58c;

1993 
1994     fT2d_DeeSCCons[dee - 1][51 - 1] = 216;
1995     fT2d_DeeSCCons[dee][51 - 1] = 67;  // 32;

1996     fT2d_DeeSCCons[dee - 1][52 - 1] = 208;
1997     fT2d_DeeSCCons[dee][52 - 1] = 59;  // 26;

1998     fT2d_DeeSCCons[dee - 1][53 - 1] = 200;
1999     fT2d_DeeSCCons[dee][53 - 1] = 51;  // 18;

2000 
2001     fT2d_DeeSCCons[dee - 1][72 - 1] = 209;
2002     fT2d_DeeSCCons[dee][72 - 1] = 60;  // 27;

2003     fT2d_DeeSCCons[dee - 1][73 - 1] = 201;
2004     fT2d_DeeSCCons[dee][73 - 1] = 52;  // 19;

2005     fT2d_DeeSCCons[dee - 1][74 - 1] = 194;
2006     fT2d_DeeSCCons[dee][74 - 1] = 45;  // 12;

2007 
2008     fT2d_DeeSCCons[dee - 1][92 - 1] = 210;
2009     fT2d_DeeSCCons[dee][92 - 1] = 61;  // 28;

2010     fT2d_DeeSCCons[dee - 1][93 - 1] = 202;
2011     fT2d_DeeSCCons[dee][93 - 1] = 53;  // 20;

2012     fT2d_DeeSCCons[dee - 1][94 - 1] = 195;
2013     fT2d_DeeSCCons[dee][94 - 1] = 46;  // 13;

2014     fT2d_DeeSCCons[dee - 1][95 - 1] = 188;
2015     fT2d_DeeSCCons[dee][95 - 1] = 39;  //  7;

2016 
2017     fT2d_DeeSCCons[dee - 1][112 - 1] = 211;
2018     fT2d_DeeSCCons[dee][112 - 1] = 62;  // 29;

2019     fT2d_DeeSCCons[dee - 1][113 - 1] = 203;
2020     fT2d_DeeSCCons[dee][113 - 1] = 54;  // 21;

2021     fT2d_DeeSCCons[dee - 1][114 - 1] = 196;
2022     fT2d_DeeSCCons[dee][114 - 1] = 47;  // 14;

2023     fT2d_DeeSCCons[dee - 1][115 - 1] = 189;
2024     fT2d_DeeSCCons[dee][115 - 1] = 40;  //  8;

2025     fT2d_DeeSCCons[dee - 1][116 - 1] = 183;
2026     fT2d_DeeSCCons[dee][116 - 1] = 34;  //  4;

2027     fT2d_DeeSCCons[dee - 1][117 - 1] = 176;
2028     fT2d_DeeSCCons[dee][117 - 1] = 27;  //  1;

2029 
2030     fT2d_DeeSCCons[dee - 1][132 - 1] = 212;
2031     fT2d_DeeSCCons[dee][132 - 1] = 63;  // 30;

2032     fT2d_DeeSCCons[dee - 1][133 - 1] = 204;
2033     fT2d_DeeSCCons[dee][133 - 1] = 55;  // 22;

2034     fT2d_DeeSCCons[dee - 1][134 - 1] = 197;
2035     fT2d_DeeSCCons[dee][134 - 1] = 48;  // 15;

2036     fT2d_DeeSCCons[dee - 1][135 - 1] = 190;
2037     fT2d_DeeSCCons[dee][135 - 1] = 41;  //  9;

2038     fT2d_DeeSCCons[dee - 1][136 - 1] = 184;
2039     fT2d_DeeSCCons[dee][136 - 1] = 35;  //  5;

2040     fT2d_DeeSCCons[dee - 1][137 - 1] = 177;
2041     fT2d_DeeSCCons[dee][137 - 1] = 28;  //  2;

2042     fT2d_DeeSCCons[dee - 1][138 - 1] = 178;
2043     fT2d_DeeSCCons[dee][138 - 1] = 29;  //  3; //  178a ;  29a;

2044 
2045     fT2d_DeeSCCons[dee - 1][152 - 1] = 213;
2046     fT2d_DeeSCCons[dee][152 - 1] = 64;  // 31;

2047     fT2d_DeeSCCons[dee - 1][153 - 1] = 205;
2048     fT2d_DeeSCCons[dee][153 - 1] = 56;  // 23;

2049     fT2d_DeeSCCons[dee - 1][154 - 1] = 198;
2050     fT2d_DeeSCCons[dee][154 - 1] = 49;  // 16;

2051     fT2d_DeeSCCons[dee - 1][155 - 1] = 191;
2052     fT2d_DeeSCCons[dee][155 - 1] = 42;  // 10;

2053     fT2d_DeeSCCons[dee - 1][156 - 1] = 185;
2054     fT2d_DeeSCCons[dee][156 - 1] = 36;  //  6;

2055     fT2d_DeeSCCons[dee - 1][157 - 1] = 178;
2056     fT2d_DeeSCCons[dee][157 - 1] = 29;  //  3; //  178b ;  29b;

2057 
2058     fT2d_DeeSCCons[dee - 1][173 - 1] = 206;
2059     fT2d_DeeSCCons[dee][173 - 1] = 57;  // 24;

2060     fT2d_DeeSCCons[dee - 1][174 - 1] = 199;
2061     fT2d_DeeSCCons[dee][174 - 1] = 50;  // 17;

2062     fT2d_DeeSCCons[dee - 1][175 - 1] = 192;
2063     fT2d_DeeSCCons[dee][175 - 1] = 43;  // 11;

2064     fT2d_DeeSCCons[dee - 1][176 - 1] = 207;
2065     fT2d_DeeSCCons[dee][176 - 1] = 58;  // 25; //  58a ;  207a;

2066 
2067     fT2d_DeeSCCons[dee - 1][193 - 1] = 207;
2068     fT2d_DeeSCCons[dee][193 - 1] = 58;  // 25; //  58b ;  207b;

2069 
2070     //.......................................... (D1,S3 or D3,S7) AND  (D2,S7 or D4,S3)

2071     fT2d_DeeSCCons[dee - 1][50 - 1] = 224;
2072     fT2d_DeeSCCons[dee][50 - 1] = 75;  // 10;

2073 
2074     fT2d_DeeSCCons[dee - 1][69 - 1] = 233;
2075     fT2d_DeeSCCons[dee][69 - 1] = 84;  // 18;

2076     fT2d_DeeSCCons[dee - 1][70 - 1] = 225;
2077     fT2d_DeeSCCons[dee][70 - 1] = 76;  // 11;

2078     fT2d_DeeSCCons[dee - 1][71 - 1] = 217;
2079     fT2d_DeeSCCons[dee][71 - 1] = 68;  //  3;

2080 
2081     fT2d_DeeSCCons[dee - 1][88 - 1] = 242;
2082     fT2d_DeeSCCons[dee][88 - 1] = 93;  // 25;

2083     fT2d_DeeSCCons[dee - 1][89 - 1] = 234;
2084     fT2d_DeeSCCons[dee][89 - 1] = 85;  // 19;

2085     fT2d_DeeSCCons[dee - 1][90 - 1] = 226;
2086     fT2d_DeeSCCons[dee][90 - 1] = 77;  // 12;

2087     fT2d_DeeSCCons[dee - 1][91 - 1] = 218;
2088     fT2d_DeeSCCons[dee][91 - 1] = 69;  //  4;

2089 
2090     fT2d_DeeSCCons[dee - 1][108 - 1] = 243;
2091     fT2d_DeeSCCons[dee][108 - 1] = 94;  // 26;

2092     fT2d_DeeSCCons[dee - 1][109 - 1] = 235;
2093     fT2d_DeeSCCons[dee][109 - 1] = 86;  // 20;

2094     fT2d_DeeSCCons[dee - 1][110 - 1] = 227;
2095     fT2d_DeeSCCons[dee][110 - 1] = 78;  // 13;

2096     fT2d_DeeSCCons[dee - 1][111 - 1] = 219;
2097     fT2d_DeeSCCons[dee][111 - 1] = 70;  //  5;

2098 
2099     fT2d_DeeSCCons[dee - 1][127 - 1] = 252;
2100     fT2d_DeeSCCons[dee][127 - 1] = 103;  // 31;

2101     fT2d_DeeSCCons[dee - 1][128 - 1] = 244;
2102     fT2d_DeeSCCons[dee][128 - 1] = 95;  // 27;

2103     fT2d_DeeSCCons[dee - 1][129 - 1] = 236;
2104     fT2d_DeeSCCons[dee][129 - 1] = 87;  // 21;

2105     fT2d_DeeSCCons[dee - 1][130 - 1] = 228;
2106     fT2d_DeeSCCons[dee][130 - 1] = 79;  // 14;

2107     fT2d_DeeSCCons[dee - 1][131 - 1] = 220;
2108     fT2d_DeeSCCons[dee][131 - 1] = 71;  //  6;

2109 
2110     fT2d_DeeSCCons[dee - 1][147 - 1] = 253;
2111     fT2d_DeeSCCons[dee][147 - 1] = 104;  // 32;

2112     fT2d_DeeSCCons[dee - 1][148 - 1] = 245;
2113     fT2d_DeeSCCons[dee][148 - 1] = 96;  // 28;

2114     fT2d_DeeSCCons[dee - 1][149 - 1] = 237;
2115     fT2d_DeeSCCons[dee][149 - 1] = 88;  // 22;

2116     fT2d_DeeSCCons[dee - 1][150 - 1] = 229;
2117     fT2d_DeeSCCons[dee][150 - 1] = 80;  // 15;

2118     fT2d_DeeSCCons[dee - 1][151 - 1] = 221;
2119     fT2d_DeeSCCons[dee][151 - 1] = 72;  //  7;

2120 
2121     fT2d_DeeSCCons[dee - 1][166 - 1] = 254;
2122     fT2d_DeeSCCons[dee][166 - 1] = 105;  // 33;

2123     fT2d_DeeSCCons[dee - 1][167 - 1] = 247;
2124     fT2d_DeeSCCons[dee][167 - 1] = 98;  // 30;

2125     fT2d_DeeSCCons[dee - 1][168 - 1] = 246;
2126     fT2d_DeeSCCons[dee][168 - 1] = 97;  // 29;

2127     fT2d_DeeSCCons[dee - 1][169 - 1] = 238;
2128     fT2d_DeeSCCons[dee][169 - 1] = 89;  // 23;

2129     fT2d_DeeSCCons[dee - 1][170 - 1] = 230;
2130     fT2d_DeeSCCons[dee][170 - 1] = 81;  // 16;

2131     fT2d_DeeSCCons[dee - 1][171 - 1] = 222;
2132     fT2d_DeeSCCons[dee][171 - 1] = 73;  //  8;

2133     fT2d_DeeSCCons[dee - 1][172 - 1] = 214;
2134     fT2d_DeeSCCons[dee][172 - 1] = 65;  //  1;

2135 
2136     fT2d_DeeSCCons[dee - 1][188 - 1] = 298;
2137     fT2d_DeeSCCons[dee][188 - 1] = 149;  // 24; //  298a ;  149a;

2138     fT2d_DeeSCCons[dee - 1][189 - 1] = 239;
2139     fT2d_DeeSCCons[dee][189 - 1] = 90;  // 24;

2140     fT2d_DeeSCCons[dee - 1][190 - 1] = 231;
2141     fT2d_DeeSCCons[dee][190 - 1] = 82;  // 17;

2142     fT2d_DeeSCCons[dee - 1][191 - 1] = 223;
2143     fT2d_DeeSCCons[dee][191 - 1] = 74;  //  9;

2144     fT2d_DeeSCCons[dee - 1][192 - 1] = 215;
2145     fT2d_DeeSCCons[dee][192 - 1] = 66;  //  2;

2146 
2147     //.......................................... (D1,S4 or D3,S6) AND  (D2,S6 or D4,S4)

2148     fT2d_DeeSCCons[dee - 1][29 - 1] = 261;
2149     fT2d_DeeSCCons[dee][29 - 1] = 112;  // 14-21;   // 261a-268a ; 112a-119a;

2150     fT2d_DeeSCCons[dee - 1][27 - 1] = 283;
2151     fT2d_DeeSCCons[dee][27 - 1] = 134;  // 33;

2152     fT2d_DeeSCCons[dee - 1][28 - 1] = 282;
2153     fT2d_DeeSCCons[dee][28 - 1] = 133;  // 32;

2154 
2155     fT2d_DeeSCCons[dee - 1][44 - 1] = 269;
2156     fT2d_DeeSCCons[dee][44 - 1] = 120;  // 22;

2157     fT2d_DeeSCCons[dee - 1][45 - 1] = 262;
2158     fT2d_DeeSCCons[dee][45 - 1] = 113;  // 15;

2159     fT2d_DeeSCCons[dee - 1][46 - 1] = 255;
2160     fT2d_DeeSCCons[dee][46 - 1] = 106;  //  8;

2161     fT2d_DeeSCCons[dee - 1][47 - 1] = 248;
2162     fT2d_DeeSCCons[dee][47 - 1] = 99;  //  4;

2163     fT2d_DeeSCCons[dee - 1][48 - 1] = 240;
2164     fT2d_DeeSCCons[dee][48 - 1] = 91;  //  2;

2165     fT2d_DeeSCCons[dee - 1][49 - 1] = 232;
2166     fT2d_DeeSCCons[dee][49 - 1] = 83;  //  1;

2167 
2168     fT2d_DeeSCCons[dee - 1][62 - 1] = 276;
2169     fT2d_DeeSCCons[dee][62 - 1] = 127;  // 29;

2170     fT2d_DeeSCCons[dee - 1][63 - 1] = 275;
2171     fT2d_DeeSCCons[dee][63 - 1] = 126;  // 28;

2172     fT2d_DeeSCCons[dee - 1][64 - 1] = 270;
2173     fT2d_DeeSCCons[dee][64 - 1] = 121;  // 23;

2174     fT2d_DeeSCCons[dee - 1][65 - 1] = 263;
2175     fT2d_DeeSCCons[dee][65 - 1] = 114;  // 16;

2176     fT2d_DeeSCCons[dee - 1][66 - 1] = 256;
2177     fT2d_DeeSCCons[dee][66 - 1] = 107;  //  9;

2178     fT2d_DeeSCCons[dee - 1][67 - 1] = 249;
2179     fT2d_DeeSCCons[dee][67 - 1] = 100;  //  5;

2180     fT2d_DeeSCCons[dee - 1][68 - 1] = 241;
2181     fT2d_DeeSCCons[dee][68 - 1] = 92;  //  3;

2182 
2183     fT2d_DeeSCCons[dee - 1][82 - 1] = 278;
2184     fT2d_DeeSCCons[dee][82 - 1] = 129;  // 31;

2185     fT2d_DeeSCCons[dee - 1][83 - 1] = 277;
2186     fT2d_DeeSCCons[dee][83 - 1] = 128;  // 30;

2187     fT2d_DeeSCCons[dee - 1][84 - 1] = 271;
2188     fT2d_DeeSCCons[dee][84 - 1] = 122;  // 24;

2189     fT2d_DeeSCCons[dee - 1][85 - 1] = 264;
2190     fT2d_DeeSCCons[dee][85 - 1] = 115;  // 17;

2191     fT2d_DeeSCCons[dee - 1][86 - 1] = 257;
2192     fT2d_DeeSCCons[dee][86 - 1] = 108;  // 10;

2193     fT2d_DeeSCCons[dee - 1][87 - 1] = 250;
2194     fT2d_DeeSCCons[dee][87 - 1] = 101;  //  6;

2195 
2196     fT2d_DeeSCCons[dee - 1][102 - 1] = 268;
2197     fT2d_DeeSCCons[dee][102 - 1] = 119;  // 21; //  268c ;  119c;

2198     fT2d_DeeSCCons[dee - 1][103 - 1] = 274;
2199     fT2d_DeeSCCons[dee][103 - 1] = 125;  // 27;

2200     fT2d_DeeSCCons[dee - 1][104 - 1] = 272;
2201     fT2d_DeeSCCons[dee][104 - 1] = 123;  // 25;

2202     fT2d_DeeSCCons[dee - 1][105 - 1] = 265;
2203     fT2d_DeeSCCons[dee][105 - 1] = 116;  // 18;

2204     fT2d_DeeSCCons[dee - 1][106 - 1] = 258;
2205     fT2d_DeeSCCons[dee][106 - 1] = 109;  // 11;

2206     fT2d_DeeSCCons[dee - 1][107 - 1] = 251;
2207     fT2d_DeeSCCons[dee][107 - 1] = 102;  //  7;

2208 
2209     fT2d_DeeSCCons[dee - 1][123 - 1] = 268;
2210     fT2d_DeeSCCons[dee][123 - 1] = 119;  // 27; //  268b ;  119b;

2211     fT2d_DeeSCCons[dee - 1][124 - 1] = 273;
2212     fT2d_DeeSCCons[dee][124 - 1] = 124;  // 26;

2213     fT2d_DeeSCCons[dee - 1][125 - 1] = 266;
2214     fT2d_DeeSCCons[dee][125 - 1] = 117;  // 19;

2215     fT2d_DeeSCCons[dee - 1][126 - 1] = 259;
2216     fT2d_DeeSCCons[dee][126 - 1] = 110;  // 12;

2217 
2218     fT2d_DeeSCCons[dee - 1][144 - 1] = 261;
2219     fT2d_DeeSCCons[dee][144 - 1] = 112;  // 27; //  261c ;  112c;

2220     fT2d_DeeSCCons[dee - 1][145 - 1] = 267;
2221     fT2d_DeeSCCons[dee][145 - 1] = 118;  // 20;

2222     fT2d_DeeSCCons[dee - 1][146 - 1] = 260;
2223     fT2d_DeeSCCons[dee][146 - 1] = 111;  // 13;

2224 
2225     fT2d_DeeSCCons[dee - 1][165 - 1] = 261;
2226     fT2d_DeeSCCons[dee][165 - 1] = 112;  // 27; //  261b ;  112b;

2227 
2228     //.......................................... D1 or D3, right half of S5

2229     fT2d_DeeSCCons[dee - 1][1 - 1] = 297;  // 34;

2230     fT2d_DeeSCCons[dee - 1][2 - 1] = 296;  // 33;

2231     fT2d_DeeSCCons[dee - 1][3 - 1] = 295;  // 32;

2232     fT2d_DeeSCCons[dee - 1][4 - 1] = 294;  // 31;

2233     fT2d_DeeSCCons[dee - 1][5 - 1] = 289;  // 26;

2234     fT2d_DeeSCCons[dee - 1][6 - 1] = 288;  // 25;

2235     fT2d_DeeSCCons[dee - 1][7 - 1] = 287;  // 24;

2236     fT2d_DeeSCCons[dee - 1][8 - 1] = 286;  // 23;

2237 
2238     fT2d_DeeSCCons[dee - 1][21 - 1] = 293;  // 30;

2239     fT2d_DeeSCCons[dee - 1][22 - 1] = 292;  // 29;

2240     fT2d_DeeSCCons[dee - 1][23 - 1] = 291;  // 28;

2241     fT2d_DeeSCCons[dee - 1][24 - 1] = 290;  // 27;

2242     fT2d_DeeSCCons[dee - 1][25 - 1] = 285;  // 22;

2243     fT2d_DeeSCCons[dee - 1][26 - 1] = 284;  // 21;

2244 
2245     fT2d_DeeSCCons[dee - 1][41 - 1] = 281;  // 20; //  281a

2246     fT2d_DeeSCCons[dee - 1][42 - 1] = 280;  // 19;

2247     fT2d_DeeSCCons[dee - 1][43 - 1] = 279;  // 18;

2248 
2249     //.......................................... D2 or D4, left half of S5

2250     fT2d_DeeSCCons[dee][1 - 1] = 148;  // 17;

2251     fT2d_DeeSCCons[dee][2 - 1] = 147;  // 16;

2252     fT2d_DeeSCCons[dee][3 - 1] = 146;  // 15;

2253     fT2d_DeeSCCons[dee][4 - 1] = 145;  // 14;

2254     fT2d_DeeSCCons[dee][5 - 1] = 140;  //  9;

2255     fT2d_DeeSCCons[dee][6 - 1] = 139;  //  8;

2256     fT2d_DeeSCCons[dee][7 - 1] = 138;  //  7;

2257     fT2d_DeeSCCons[dee][8 - 1] = 137;  //  6;

2258 
2259     fT2d_DeeSCCons[dee][21 - 1] = 144;  // 13;

2260     fT2d_DeeSCCons[dee][22 - 1] = 143;  // 12;

2261     fT2d_DeeSCCons[dee][23 - 1] = 142;  // 11;

2262     fT2d_DeeSCCons[dee][24 - 1] = 141;  // 10;

2263     fT2d_DeeSCCons[dee][25 - 1] = 136;  //  5;

2264     fT2d_DeeSCCons[dee][26 - 1] = 135;  //  4;

2265 
2266     fT2d_DeeSCCons[dee][41 - 1] = 132;  // 3; //  132a

2267     fT2d_DeeSCCons[dee][42 - 1] = 131;  // 2;

2268     fT2d_DeeSCCons[dee][43 - 1] = 130;  // 1;

2269   }
2270 
2271   //........................ ECNA numbers from numbers for constructions: fT2d_RecovDeeSC[][]

2272 
2273   for (Int_t i0EEDee = 0; i0EEDee < MaxEEDee; i0EEDee++) {
2274     for (Int_t i_ecna = 0; i_ecna < MaxDeeSC; i_ecna++) {
2275       //....... test to avoid the -1 init value in 2nd index of array fT2d_RecovDeeSC[][]

2276       //        (part of the matrix without real SC counterpart)

2277       if (fT2d_DeeSCCons[i0EEDee][i_ecna] >= 0 && fT2d_DeeSCCons[i0EEDee][i_ecna] <= MaxEESCForCons) {
2278         fT2d_RecovDeeSC[i0EEDee][fT2d_DeeSCCons[i0EEDee][i_ecna] - 1] = i_ecna + 1;
2279       }
2280     }
2281   }
2282 }
2283 //------------ (end of BuildEndcapSCTable) -------------------------

2284 
2285 //===============================================================================

2286 //

2287 //        Get1DeeCrysFrom1DeeSCEcnaAnd0SCEcha

2288 //        GetDeeCrysFromDeeEcha

2289 //

2290 //===============================================================================

2291 Int_t TEcnaNumbering::Get1DeeCrysFrom1DeeSCEcnaAnd0SCEcha(const Int_t& n1DeeSCEcna,
2292                                                           const Int_t& i0SCEcha,
2293                                                           const TString& sDeeDir) {
2294   //get crystal number in Dee from SC number in Dee

2295   // and from Electronic Channel number in super-crystal

2296 
2297   Int_t n1DeeCrys = 0;
2298   Int_t i0DeeDir = GetDeeDirIndex(sDeeDir);
2299 
2300   if (fT3dDeeCrys == nullptr) {
2301     BuildEndcapCrysTable();
2302   }
2303 
2304   if ((n1DeeSCEcna >= 1) && (n1DeeSCEcna <= fEcal->MaxSCEcnaInDee())) {
2305     if (i0SCEcha >= 0 && i0SCEcha < fEcal->MaxCrysInSC()) {
2306       n1DeeCrys = fT3dDeeCrys[n1DeeSCEcna - 1][i0SCEcha][i0DeeDir];
2307     } else {
2308       n1DeeCrys = -2;  // Electronic Cnannel in Super-Crystal out of range

2309       std::cout << "!TEcnaNumbering::Get1DeeCrysFrom1DeeSCEcnaAnd0SCEcha(...)> Electronic Channel in SuperCrystal = "
2310                 << i0SCEcha + 1 << ". Out of range (range = [1," << fEcal->MaxCrysInSC() << "])" << fTTBELL
2311                 << std::endl;
2312     }
2313   } else {
2314     n1DeeCrys = -3;  // Super-Crystal number in Dee out of range

2315     std::cout << "!TEcnaNumbering::Get1DeeCrysFrom1DeeSCEcnaAnd0SCEcha(...)> Super-Crystal number in Dee out of range."
2316               << " n1DeeSCEcna = " << n1DeeSCEcna << fTTBELL << std::endl;
2317   }
2318 
2319   return n1DeeCrys;  // Range = [1,5000]

2320 }
2321 
2322 //===============================================================================

2323 //

2324 //                Get1SCEchaFrom1DeeCrys, Get1DeeSCEcnaFrom1DeeCrys

2325 //

2326 //===============================================================================

2327 
2328 Int_t TEcnaNumbering::Get1SCEchaFrom1DeeCrys(const Int_t& n1DeeCrys, const TString& sDeeDir) {
2329   // get Electronic Channel number in Super-Crystal from Crystal ECNA number in Dee

2330 
2331   Int_t n1SCEcha = -1;
2332   Int_t iDeeDir = GetDeeDirIndex(sDeeDir);
2333 
2334   if (n1DeeCrys >= 1 && n1DeeCrys <= fEcal->MaxCrysEcnaInDee()) {
2335     n1SCEcha = fT2dSCEcha[n1DeeCrys - 1][iDeeDir];
2336   } else {
2337     n1SCEcha = -2;
2338     std::cout << "!TEcnaNumbering::Get1SCEchaFrom1DeeCrys(...)> Crystal number in Dee out of range."
2339               << " n1DeeCrys = " << n1DeeCrys << "(max = " << fEcal->MaxCrysEcnaInDee() << ")" << fTTBELL << std::endl;
2340   }
2341   return n1SCEcha;  // range = [1,25]

2342 }
2343 
2344 Int_t TEcnaNumbering::Get1DeeSCEcnaFrom1DeeCrys(const Int_t& n1DeeCrys, const TString& sDeeDir) {
2345   // get Super-Crystal number in Dee from Crystal number in Dee

2346 
2347   Int_t n1DeeSCEcna = 0;
2348   Int_t iDeeDir = GetDeeDirIndex(sDeeDir);
2349 
2350   if (n1DeeCrys >= 1 && n1DeeCrys <= fEcal->MaxCrysEcnaInDee()) {
2351     n1DeeSCEcna = fT2dDeeSC[n1DeeCrys - 1][iDeeDir];
2352   } else {
2353     n1DeeSCEcna = -1;
2354     std::cout << "!TEcnaNumbering::Get1DeeSCEcnaFrom1DeeCrys(...)> Crystal number in Dee out of range."
2355               << " n1DeeCrys = " << n1DeeCrys << "(max = " << fEcal->MaxCrysEcnaInDee() << ")" << fTTBELL << std::endl;
2356   }
2357   return n1DeeSCEcna;  // range = [1,200]

2358 }
2359 
2360 //===============================================================================

2361 //

2362 //          GetSCEchaFromDeeEcha

2363 //          Get1DeeSCEcnaFromDeeEcha

2364 //

2365 //===============================================================================

2366 
2367 Int_t TEcnaNumbering::Get1SCEchaFrom0DeeEcha(const Int_t& i0DeeEcha) {
2368   //get electronic channel number in super-crystal from electronic channel number in Dee

2369 
2370   Int_t i0DeeSC = i0DeeEcha / fEcal->MaxCrysInSC();
2371   Int_t n1SCEcha = i0DeeEcha - fEcal->MaxCrysInSC() * i0DeeSC + 1;
2372 
2373   return n1SCEcha;  //  range = [1,25]

2374 }
2375 
2376 Int_t TEcnaNumbering::Get1DeeSCEcnaFrom0DeeEcha(const Int_t& i0DeeEcha) {
2377   //get super-crystal number from electronic channel number in Dee

2378 
2379   Int_t n1DeeSC = i0DeeEcha / fEcal->MaxCrysInSC() + 1;
2380 
2381   return n1DeeSC;  //  range = [1,200]

2382 }
2383 
2384 //--------------------------------------------------------------------------------

2385 //

2386 //       Correspondance (n1DeeNumber, DeeSC)  <->  (DS, DSSC, DeeSCCons)

2387 //

2388 //       GetDSFrom1DeeSCEcna,        GetDSSCFrom1DeeSCEcna,

2389 //       GetDeeSCConsFrom1DeeSCEcna, Get1DeeSCEcnaFromDeeSCCons

2390 //

2391 //       Get the values from the relevant arrays

2392 //       with cross-check of the index values in argument

2393 //--------------------------------------------------------------------------------

2394 
2395 Int_t TEcnaNumbering::GetDSFrom1DeeSCEcna(const Int_t& n1DeeNumber, const Int_t& n1DeeSCEcna) {
2396   // Get Data Sector number from SC ECNA number in Dee

2397 
2398   Int_t data_sector = -1;
2399 
2400   if (n1DeeNumber > 0 && n1DeeNumber <= fEcal->MaxDeeInEE()) {
2401     if (n1DeeSCEcna > 0 && n1DeeSCEcna <= fEcal->MaxSCEcnaInDee()) {
2402       data_sector = fT2d_DS[n1DeeNumber - 1][n1DeeSCEcna - 1];
2403     } else {
2404       std::cout << "!TEcnaNumbering::GetDSFrom1DeeSCEcna(...)> n1DeeSCEcna = " << n1DeeSCEcna
2405                 << ". Out of range ( range = [1," << fEcal->MaxSCEcnaInDee() << "] )" << fTTBELL << std::endl;
2406     }
2407   } else {
2408     if (n1DeeNumber != 0) {
2409       std::cout << "!TEcnaNumbering::GetDSFrom1DeeSCEcna(...)> n1DeeNumber = " << n1DeeNumber
2410                 << ". Out of range ( range = [1," << fEcal->MaxDeeInEE() << "] )" << fTTBELL << std::endl;
2411     } else {
2412       std::cout << "TEcnaNumbering::GetDSFrom1DeeSCEcna(...)> Dee = " << n1DeeNumber << ". Out of range (range = [1,"
2413                 << fEcal->MaxDeeInEE() << "])" << fTTBELL << std::endl;
2414     }
2415   }
2416   return data_sector;
2417 }
2418 //..........................................................................................

2419 Int_t TEcnaNumbering::GetDSSCFrom1DeeSCEcna(const Int_t& n1DeeNumber, const Int_t& n1DeeSCEcna, const Int_t& n1SCEcha) {
2420   //.......... Get the correct SC number for the unconnected SC's (inner border)

2421   Int_t ds_sc = GetDSSCFrom1DeeSCEcna(n1DeeNumber, n1DeeSCEcna);
2422 
2423   if (n1DeeSCEcna == 29 || n1DeeSCEcna == 32) {
2424     if (n1SCEcha == 11) {
2425       if (ds_sc == 14) {
2426         ds_sc = 21;
2427       }  // 14 <=> 261/BR OR 112/BL

2428     }
2429     if (n1SCEcha == 1 || n1SCEcha == 2 || n1SCEcha == 3 || n1SCEcha == 6 || n1SCEcha == 7) {
2430       if (ds_sc == 3) {
2431         ds_sc = 25;
2432       }  // 3 <=> 178/TR OR 29/TL

2433     }
2434   }
2435   return ds_sc;
2436 }
2437 //..........................................................................................

2438 Int_t TEcnaNumbering::GetDSSCFrom1DeeSCEcna(const Int_t& n1DeeNumber, const Int_t& n1DeeSCEcna) {
2439   // Get SC number in Data Sector from SC Ecna number in Dee

2440 
2441   Int_t ds_sc = -1;
2442 
2443   if (n1DeeNumber > 0 && n1DeeNumber <= fEcal->MaxDeeInEE()) {
2444     if (n1DeeSCEcna > 0 && n1DeeSCEcna <= fEcal->MaxSCEcnaInDee()) {
2445       ds_sc = fT2d_DSSC[n1DeeNumber - 1][n1DeeSCEcna - 1];  // 25 (not 3) or 14 (not 21) if n1DeeSCEcna = 32 or 29

2446                                                             // 25 and 14 => 5 Xtals,  3 and 21 => 1 Xtal

2447     } else {
2448       std::cout << "!TEcnaNumbering::GetDSSCFrom1DeeSCEcna(...)> n1DeeSCEcna = " << n1DeeSCEcna
2449                 << ". Out of range ( range = [1," << fEcal->MaxSCEcnaInDee() << "] )" << fTTBELL << std::endl;
2450     }
2451   } else {
2452     if (n1DeeNumber != 0) {
2453       std::cout << "!TEcnaNumbering::GetDSSCFrom1DeeSCEcna(...)> n1DeeNumber = " << n1DeeNumber
2454                 << ". Out of range ( range = [1," << fEcal->MaxDeeInEE() << "] )" << fTTBELL << std::endl;
2455     } else {
2456       std::cout << "TEcnaNumbering::GetDSSCFrom1DeeSCEcna(...)> Dee = " << n1DeeNumber << ". Out of range (range = [1,"
2457                 << fEcal->MaxDeeInEE() << "])" << fTTBELL << std::endl;
2458     }
2459   }
2460   return ds_sc;
2461 }
2462 //..........................................................................................

2463 Int_t TEcnaNumbering::GetDeeSCConsFrom1DeeSCEcna(const Int_t& n1DeeNumber, const Int_t& n1DeeSCEcna) {
2464   // Get SC number for Construction in Dee from SC ECNA number in Dee

2465 
2466   Int_t dee_sc_cons = -1;
2467 
2468   if (n1DeeNumber > 0 && n1DeeNumber <= fEcal->MaxDeeInEE()) {
2469     if (n1DeeSCEcna > 0 && n1DeeSCEcna <= fEcal->MaxSCEcnaInDee()) {
2470       dee_sc_cons = fT2d_DeeSCCons[n1DeeNumber - 1][n1DeeSCEcna - 1];
2471     } else {
2472       std::cout << "!TEcnaNumbering::GetDeeSCConsFrom1DeeSCEcna(...)> *** WARNING *** n1DeeSCEcna = " << n1DeeSCEcna
2473                 << ". Out of range ( range = [1," << fEcal->MaxSCEcnaInDee() << "] ). Nb for const. forced to "
2474                 << fT2d_DeeSCCons[n1DeeNumber - 1][19] << "." << std::endl;
2475       dee_sc_cons = fT2d_DeeSCCons[n1DeeNumber - 1][19];
2476     }
2477   } else {
2478     if (n1DeeNumber != 0) {
2479       std::cout << "!TEcnaNumbering::GetDeeSCConsFrom1DeeSCEcna(...)> n1DeeNumber = " << n1DeeNumber
2480                 << ". Out of range ( range = [1," << fEcal->MaxDeeInEE() << "] )" << fTTBELL << std::endl;
2481     } else {
2482       std::cout << "TEcnaNumbering::GetDeeSCConsFrom1DeeSCEcna(...)> Dee = " << n1DeeNumber
2483                 << ". Out of range (range = [1," << fEcal->MaxDeeInEE() << "])" << fTTBELL << std::endl;
2484     }
2485   }
2486   return dee_sc_cons;
2487 }
2488 //..........................................................................................

2489 Int_t TEcnaNumbering::GetDeeSCConsFrom1DeeSCEcna(const Int_t& n1DeeNumber,
2490                                                  const Int_t& n1DeeSCEcna,
2491                                                  const Int_t& n1SCEcha) {
2492   //.......... Get the correct SC number (for cons) for the unconnected SC's (inner border)

2493   Int_t dee_sc_cons = GetDeeSCConsFrom1DeeSCEcna(n1DeeNumber, n1DeeSCEcna);
2494 
2495   if (n1DeeSCEcna == 29 || n1DeeSCEcna == 32) {
2496     if (n1SCEcha == 11) {
2497       if (dee_sc_cons == 261) {
2498         dee_sc_cons = 268;
2499       }  // 261<=>14/BR

2500       if (dee_sc_cons == 112) {
2501         dee_sc_cons = 119;
2502       }  // 112<=>14/BL

2503     }
2504     if (n1SCEcha == 1 || n1SCEcha == 2 || n1SCEcha == 3 || n1SCEcha == 6 || n1SCEcha == 7) {
2505       if (dee_sc_cons == 178) {
2506         dee_sc_cons = 207;
2507       }  // 178<=>3/TR

2508       if (dee_sc_cons == 29) {
2509         dee_sc_cons = 58;
2510       }  //  29<=>3/TL

2511     }
2512   }
2513   return dee_sc_cons;
2514 }
2515 //..........................................................................................

2516 Int_t TEcnaNumbering::Get1DeeSCEcnaFromDeeSCCons(const Int_t& n1DeeNumber, const Int_t& DeeSCCons) {
2517   // Get SC Ecna number in Dee from SC number for Construction in Dee

2518 
2519   Int_t dee_sc_ecna = -1;
2520 
2521   if (n1DeeNumber > 0 && n1DeeNumber <= fEcal->MaxDeeInEE()) {
2522     Int_t off_set_cons = 0;
2523     if (n1DeeNumber == 1 || n1DeeNumber == 3) {
2524       off_set_cons = fEcal->MaxSCForConsInDee();
2525     }
2526 
2527     if (DeeSCCons > off_set_cons && DeeSCCons <= fEcal->MaxSCForConsInDee() + off_set_cons) {
2528       dee_sc_ecna = fT2d_RecovDeeSC[n1DeeNumber - 1][DeeSCCons - 1];
2529     } else {
2530       std::cout << "!TEcnaNumbering::Get1DeeSCEcnaFromDeeSCCons(...)> DeeSCCons = " << DeeSCCons
2531                 << ". Out of range ( range = [ " << off_set_cons + 1 << "," << fEcal->MaxSCForConsInDee() + off_set_cons
2532                 << "] )" << fTTBELL << std::endl;
2533     }
2534   } else {
2535     if (n1DeeNumber != 0) {
2536       std::cout << "!TEcnaNumbering::Get1DeeSCEcnaFromDeeSCCons(...)> n1DeeNumber = " << n1DeeNumber
2537                 << ". Out of range ( range = [1," << fEcal->MaxDeeInEE() << "] )" << fTTBELL << std::endl;
2538     } else {
2539       std::cout << "TEcnaNumbering::Get1DeeSCEcnaFromDeeSCCons(...)> Dee = " << n1DeeNumber
2540                 << ". Out of range (range = [1," << fEcal->MaxDeeInEE() << "])" << fTTBELL << std::endl;
2541     }
2542   }
2543   return dee_sc_ecna;
2544 }
2545 
2546 TString TEcnaNumbering::GetSCType(const Int_t& nb_for_cons) {
2547   // gives the special not connected SC's

2548 
2549   TString SCType = "Connected";  // => default type

2550 
2551   if (nb_for_cons == 182 || nb_for_cons == 33) {
2552     SCType = "NotConnected";
2553   }  // (D1,S1) (D3,S9) || (D2,S9) (D4,S1)

2554 
2555   if (nb_for_cons == 178 || nb_for_cons == 29) {
2556     SCType = "NotConnected";
2557   }  // (D1,S2) (D3,S8) || (D2,S8) (D4,S2)

2558   if (nb_for_cons == 207 || nb_for_cons == 58) {
2559     SCType = "NotConnected";
2560   }
2561 
2562   if (nb_for_cons == 298 || nb_for_cons == 149) {
2563     SCType = "NotConnected";
2564   }  // (D1,S3) (D3,S7) || (D2,S7) (D4,S3)

2565 
2566   if (nb_for_cons == 261 || nb_for_cons == 112) {
2567     SCType = "NotConnected";
2568   }  // (D1,S4) (D3,S6) || (D2,S6) (D4,S4)

2569   if (nb_for_cons == 268 || nb_for_cons == 119) {
2570     SCType = "NotConnected";
2571   }
2572 
2573   if (nb_for_cons == 281 || nb_for_cons == 132) {
2574     SCType = "NotConnected";
2575   }  // (D1,S5) (D3,S5) || (D2,S5) (D4,S5)

2576 
2577   if (nb_for_cons == 161 || nb_for_cons == 12) {
2578     SCType = "NotComplete";
2579   }  // (D1,S1) (D3,S9) || (D2,S9) (D4,S1)

2580   if (nb_for_cons == 216 || nb_for_cons == 67) {
2581     SCType = "NotComplete";
2582   }  // (D1,S2) (D3,S8) || (D2,S8) (D4,S2)

2583   if (nb_for_cons == 224 || nb_for_cons == 75) {
2584     SCType = "NotComplete";
2585   }  // (D1,S3) (D3,S7) || (D2,S7) (D4,S3)

2586   if (nb_for_cons == 286 || nb_for_cons == 137) {
2587     SCType = "NotComplete";
2588   }  // (D1,S5) (D3,S5) || (D2,S5) (D4,S5)

2589 
2590   return SCType;
2591 }
2592 
2593 Int_t TEcnaNumbering::StexEchaForCons(const Int_t& n1DeeNumber, const Int_t& i0StexEcha) {
2594   Int_t n1StexStin = Get1StexStinFrom0StexEcha(i0StexEcha);
2595   return fT2d_DeeSCCons[n1DeeNumber - 1][n1StexStin - 1];
2596 }
2597 // return -1 if the SC does not correspond to a real SC; return the number for construction otherwise

2598 
2599 //===========================================================================

2600 //

2601 //                        GetSCQuadFrom1DeeSCEcna

2602 //

2603 //===========================================================================

2604 
2605 TString TEcnaNumbering::GetSCQuadFrom1DeeSCEcna(const Int_t& n1DeeSCEcna) {
2606   //gives the quadrant type ("top" or "bottom") from the SC number in Dee

2607 
2608   TString SCQuad = "top";  // => default value

2609 
2610   if (n1DeeSCEcna >= 1 && n1DeeSCEcna <= 10) {
2611     SCQuad = "bottom";
2612   }
2613   if (n1DeeSCEcna >= 21 && n1DeeSCEcna <= 30) {
2614     SCQuad = "bottom";
2615   }
2616   if (n1DeeSCEcna >= 41 && n1DeeSCEcna <= 50) {
2617     SCQuad = "bottom";
2618   }
2619   if (n1DeeSCEcna >= 61 && n1DeeSCEcna <= 70) {
2620     SCQuad = "bottom";
2621   }
2622   if (n1DeeSCEcna >= 81 && n1DeeSCEcna <= 90) {
2623     SCQuad = "bottom";
2624   }
2625   if (n1DeeSCEcna >= 101 && n1DeeSCEcna <= 110) {
2626     SCQuad = "bottom";
2627   }
2628   if (n1DeeSCEcna >= 121 && n1DeeSCEcna <= 130) {
2629     SCQuad = "bottom";
2630   }
2631   if (n1DeeSCEcna >= 141 && n1DeeSCEcna <= 150) {
2632     SCQuad = "bottom";
2633   }
2634   if (n1DeeSCEcna >= 161 && n1DeeSCEcna <= 170) {
2635     SCQuad = "bottom";
2636   }
2637   if (n1DeeSCEcna >= 181 && n1DeeSCEcna <= 190) {
2638     SCQuad = "bottom";
2639   }
2640 
2641   return SCQuad;
2642 }
2643 Int_t TEcnaNumbering::GetSCQuadTypeIndex(const TString& SCQuadType, const TString& sDeeDir) {
2644   //gives the index of the SC quadrant type (top/right, top/left, bottom/left, bottom/right)

2645   // = quadrant number - 1

2646 
2647   Int_t itype = 0;  // => default

2648   if (SCQuadType == "top" && sDeeDir == "right") {
2649     itype = 0;
2650   }
2651   if (SCQuadType == "top" && sDeeDir == "left") {
2652     itype = 1;
2653   }
2654   if (SCQuadType == "bottom" && sDeeDir == "left") {
2655     itype = 2;
2656   }
2657   if (SCQuadType == "bottom" && sDeeDir == "right") {
2658     itype = 3;
2659   }
2660   return itype;
2661 }
2662 //===========================================================================

2663 //

2664 //                    GetEEDeeType, GetDeeDirViewedFromIP

2665 //

2666 //===========================================================================

2667 TString TEcnaNumbering::GetEEDeeEndcap(const Int_t& n1DeeNumber) {
2668   //gives the Endcap (EE+ or EE-) of the Dee (H. Heath, CMS NOTE 2006/027)

2669 
2670   TString eetype = "EE+";  // => default

2671   if (n1DeeNumber == 1 || n1DeeNumber == 2) {
2672     eetype = "EE+";
2673   }
2674   if (n1DeeNumber == 3 || n1DeeNumber == 4) {
2675     eetype = "EE-";
2676   }
2677   return eetype;
2678 }
2679 TString TEcnaNumbering::GetEEDeeType(const Int_t& n1DeeNumber) {
2680   //gives the EE +/- and Forward/Near of the Dee (H. Heath, CMS NOTE 2006/027)

2681 
2682   TString type = "EE+F";  // => default

2683   if (n1DeeNumber == 1) {
2684     type = "EE+F";
2685   }
2686   if (n1DeeNumber == 2) {
2687     type = "EE+N";
2688   }
2689   if (n1DeeNumber == 3) {
2690     type = "EE-N";
2691   }
2692   if (n1DeeNumber == 4) {
2693     type = "EE-F";
2694   }
2695   return type;
2696 }
2697 
2698 TString TEcnaNumbering::GetDeeDirViewedFromIP(const Int_t& n1DeeNumber) {
2699   //gives the direction (left/right) of the IX axis of the Dee

2700   // looking from the interaction point

2701 
2702   TString sDeeDir = "right";  // => default

2703   if ((n1DeeNumber == 1) || (n1DeeNumber == 3)) {
2704     sDeeDir = "right";
2705   }
2706   if ((n1DeeNumber == 2) || (n1DeeNumber == 4)) {
2707     sDeeDir = "left";
2708   }
2709   return sDeeDir;
2710 }
2711 Int_t TEcnaNumbering::GetDeeDirIndex(const TString& sDeeDir) {
2712   //gives the index of the direction (left,right) of the IX axis of the Dee

2713   // looking from the interaction point (right = 0, left = 1)

2714 
2715   Int_t iDeeDir = 0;  // => default

2716   if (sDeeDir == "right") {
2717     iDeeDir = 0;
2718   }
2719   if (sDeeDir == "left") {
2720     iDeeDir = 1;
2721   }
2722   return iDeeDir;
2723 }
2724 
2725 //==============================================================================

2726 //

2727 //    GetIXCrysInSC,  GetJYCrysInSC

2728 //    GetIXSCInDee,   GetJYSCInDee

2729 //    GetIXCrysInDee, GetJYCrysInDConsFrom1DeeSCEcna(fStexNumber, StexStinEcna);

2730 //

2731 //==============================================================================

2732 Int_t TEcnaNumbering::GetIXCrysInSC(const Int_t& n1DeeNumber, const Int_t& DeeSC, const Int_t& i0SCEcha) {
2733   //Gives Crys IX in SC for a given (n1DeeNumber, DeeSC, i0SCEcha)

2734 
2735   TString SCQuadType = GetSCQuadFrom1DeeSCEcna(DeeSC);
2736   TString sDeeDir = GetDeeDirViewedFromIP(n1DeeNumber);
2737   Int_t type_index = GetSCQuadTypeIndex(SCQuadType, sDeeDir);
2738   Int_t IXCrysInSC = fT2d_ich_IX[type_index][i0SCEcha + 1] + 1;
2739   return IXCrysInSC;  // possible values: 1,2,3,4,5

2740 }
2741 
2742 Int_t TEcnaNumbering::GetIXSCInDee(const Int_t& DeeSC) {
2743   //Gives SC IX in Dee for a given (DeeSC)

2744 
2745   Int_t IXSCInDee = (DeeSC - 1) / fEcal->MaxSCIYInDee() + 1;
2746   return IXSCInDee;  //  possible values: 1,2,...,9,10

2747 }
2748 
2749 Int_t TEcnaNumbering::GetIXCrysInDee(const Int_t& n1DeeNumber, const Int_t& DeeSC, const Int_t& i0SCEcha) {
2750   //Gives Crys IX in Dee for a given (n1DeeNumber, DeeSC, i0SCEcha)

2751 
2752   Int_t IXCrysInDee = (GetIXSCInDee(DeeSC) - 1) * fEcal->MaxCrysIXInSC() + GetIXCrysInSC(n1DeeNumber, DeeSC, i0SCEcha);
2753   return IXCrysInDee;  // possible values: 1,2,...,49,50

2754 }
2755 //---------------------------------------------------------------------------------

2756 Int_t TEcnaNumbering::GetJYCrysInSC(const Int_t& n1DeeNumber, const Int_t& DeeSC, const Int_t& i0SCEcha) {
2757   //Gives Crys JY in SC  for a given (n1DeeNumber, DeeSC, i0SCEcha)

2758 
2759   TString SCQuadType = GetSCQuadFrom1DeeSCEcna(DeeSC);
2760   TString sDeeDir = GetDeeDirViewedFromIP(n1DeeNumber);
2761   Int_t type_index = GetSCQuadTypeIndex(SCQuadType, sDeeDir);
2762   Int_t JYCrysInSC = fT2d_jch_JY[type_index][i0SCEcha + 1] + 1;
2763   return JYCrysInSC;  // possible values: 1,2,3,4,5

2764 }
2765 
2766 Int_t TEcnaNumbering::GetJYSCInDee(const Int_t& DeeSC) {
2767   //Gives SC JY in Dee for a given (n1DeeNumber, DeeSC, i0SCEcha)

2768 
2769   Int_t JYSCInDee = (DeeSC - 1) % fEcal->MaxSCIYInDee() + 1;
2770   return JYSCInDee;  //  possible values: 1,2,...,19,20

2771 }
2772 
2773 Int_t TEcnaNumbering::GetJYCrysInDee(const Int_t& n1DeeNumber, const Int_t& DeeSC, const Int_t& i0SCEcha) {
2774   //Gives Crys JY in Dee for a given (n1DeeNumber, DeeSC, i0SCEcha)

2775 
2776   Int_t JYCrysInDee = (GetJYSCInDee(DeeSC) - 1) * fEcal->MaxCrysIYInSC() + GetJYCrysInSC(n1DeeNumber, DeeSC, i0SCEcha);
2777   return JYCrysInDee;  // possible values: 1,2,...,99,100

2778 }
2779 //---------------------------------------------------------------------------------

2780 Int_t TEcnaNumbering::GetMaxSCInDS(const Int_t& DeeDS) {
2781   // Gives the number of SC's in Data Sector DeeDS

2782 
2783   Int_t nb_of_sc = -1;
2784   if (DeeDS == 1 || DeeDS == 9) {
2785     nb_of_sc = 33;
2786   }
2787   if (DeeDS == 2 || DeeDS == 8) {
2788     nb_of_sc = 32;
2789   }
2790   if (DeeDS == 3 || DeeDS == 7) {
2791     nb_of_sc = 34;
2792   }
2793   if (DeeDS == 4 || DeeDS == 6) {
2794     nb_of_sc = 33;
2795   }
2796   if (DeeDS == 5) {
2797     nb_of_sc = 34;
2798   }
2799   return nb_of_sc;
2800 }
2801 
2802 //==============================================================================

2803 //

2804 //       GetIXMin, GetIXMax,  GetIIXMin, GetIIXMax

2805 //

2806 //==============================================================================

2807 Double_t TEcnaNumbering::GetIIXMin(const Int_t& DeeSC) {
2808   //Gives IIXMin for a given DeeSC , unit = crystal

2809 
2810   Double_t IX_min = (Double_t)((DeeSC - 1) / fEcal->MaxSCIYInDee()) * fEcal->MaxCrysIXInSC() + 1.;
2811   return IX_min;
2812 }
2813 
2814 Double_t TEcnaNumbering::GetIIXMax(const Int_t& DeeSC) {
2815   //Gives IIXMax for a given DeeSC , unit = crystal

2816 
2817   Double_t IX_max = ((Double_t)((DeeSC - 1) / fEcal->MaxSCIYInDee()) + 1.) * fEcal->MaxCrysIXInSC();
2818   return IX_max;
2819 }
2820 
2821 Double_t TEcnaNumbering::GetIIXMin() {
2822   //Gives IIXMin for Dee , unit = Super-crystal

2823 
2824   Double_t i_IX_min = (Int_t)1.;
2825   return i_IX_min;
2826 }
2827 
2828 Double_t TEcnaNumbering::GetIIXMax() {
2829   //Gives IIXMax for Dee , unit = Super-crystal

2830 
2831   Double_t i_IX_max = (Int_t)fEcal->MaxSCIXInDee();
2832   return i_IX_max;
2833 }
2834 
2835 //==============================================================================

2836 //

2837 //     GetIIYMin, GetIIYMax

2838 //

2839 //==============================================================================

2840 Double_t TEcnaNumbering::GetJIYMin(const Int_t& n1DeeNumber, const Int_t& DeeSC) {
2841   //Gives JIYMin for a given Super-Crystal

2842 
2843   Double_t IY_DeeSC = DeeSC % fEcal->MaxSCIYInDee();
2844   if (IY_DeeSC == 0.) {
2845     IY_DeeSC = fEcal->MaxSCIYInDee();
2846   }
2847 
2848   Double_t j_IY_min = (IY_DeeSC - 1) * fEcal->MaxCrysIYInSC() + 1.;
2849 
2850   return j_IY_min;
2851 }
2852 //-----------------------------------------------------------------------------------------

2853 Double_t TEcnaNumbering::GetJIYMax(const Int_t& n1DeeNumber, const Int_t& DeeSC) {
2854   //Gives JIYMax for a given Super-Crystal

2855 
2856   Double_t IY_DeeSC = DeeSC % fEcal->MaxSCIYInDee();
2857   if (IY_DeeSC == 0) {
2858     IY_DeeSC = fEcal->MaxSCIYInDee();
2859   }
2860 
2861   Double_t j_IY_max = IY_DeeSC * fEcal->MaxCrysIYInSC();
2862 
2863   return j_IY_max;
2864 }
2865 
2866 //-----------------------------------------------------------------------------------------

2867 Double_t TEcnaNumbering::GetJIYMin(const Int_t& n1DeeNumber) {
2868   //Gives JIYMin for a given Dee

2869 
2870   Double_t j_IY_min = (Double_t)1.;
2871 
2872   return j_IY_min;
2873 }
2874 //-----------------------------------------------------------------------------------------

2875 Double_t TEcnaNumbering::GetJIYMax(const Int_t& n1DeeNumber) {
2876   //Gives JIYMax for a given Dee

2877 
2878   Double_t j_IY_max = (Double_t)fEcal->MaxSCIYInDee();
2879 
2880   return j_IY_max;
2881 }
2882 //=====================================================================

2883 TString TEcnaNumbering::GetDeeHalfEndcap(const Int_t& n1DeeNumber) {
2884   //gives the half-endcap of the Dee ("EE+" or "EE-")

2885 
2886   TString type = "EE-";  // => default value

2887 
2888   if (n1DeeNumber == 1 || n1DeeNumber == 2) {
2889     type = "EE+";
2890   }
2891   if (n1DeeNumber == 3 || n1DeeNumber == 4) {
2892     type = "EE-";
2893   }
2894 
2895   return type;
2896 }
2897 //==============================================================================

2898 //

2899 //   GetXDirectionEE, GetYDirectionEE, GetJYDirectionEE

2900 //

2901 //==============================================================================

2902 TString TEcnaNumbering::GetXDirectionEE(const Int_t& n1DeeNumber) {
2903   TString xdirection = "x";  // DEFAULT

2904 
2905   if (GetEEDeeType(n1DeeNumber) == "EE+F") {
2906     xdirection = "-x";
2907   }  //  Dee 1

2908   if (GetEEDeeType(n1DeeNumber) == "EE+N") {
2909     xdirection = "-x";
2910   }  //  Dee 2

2911   if (GetEEDeeType(n1DeeNumber) == "EE-N") {
2912     xdirection = "x";
2913   }  //  Dee 3

2914   if (GetEEDeeType(n1DeeNumber) == "EE-F") {
2915     xdirection = "x";
2916   }  //  Dee 4

2917 
2918   return xdirection;
2919 }
2920 //---------------------------------------------------------

2921 TString TEcnaNumbering::GetYDirectionEE(const Int_t& n1DeeNumber) {
2922   TString ydirection = "-x";  // DEFAULT

2923 
2924   if (GetEEDeeType(n1DeeNumber) == "endcap+") {
2925     ydirection = "-x";
2926   }
2927   if (GetEEDeeType(n1DeeNumber) == "endcap-") {
2928     ydirection = "-x";
2929   }
2930 
2931   return ydirection;
2932 }
2933 
2934 //---------------------------------------------------------

2935 TString TEcnaNumbering::GetJYDirectionEE(const Int_t& n1DeeNumber) {
2936   TString jydirection = "x";  // ALWAYS IN THIS CASE

2937   return jydirection;
2938 }
2939 
2940 //==========================================================================================================

2941 //

2942 //

2943 //                               B A R R E L    A N D    E N D C A P

2944 //

2945 //

2946 //

2947 //==========================================================================================================

2948 //===============================================================================

2949 //

2950 //        Get1StexStinFrom0StexEcha

2951 //     // GetStinEchaFromStexEcha

2952 //        Get0StexEchaFrom1StexStinAnd0StinEcha

2953 //        Get1StexCrysFrom1StexStinAnd0StinEcha

2954 //

2955 //===============================================================================

2956 Int_t TEcnaNumbering::Get1StexStinFrom0StexEcha(const Int_t& i0StexEcha) {
2957   Int_t n1StexStin = 0;
2958 
2959   if (fFlagSubDet == "EB") {
2960     n1StexStin = Get1SMTowFrom0SMEcha(i0StexEcha);
2961   }
2962   if (fFlagSubDet == "EE") {
2963     n1StexStin = Get1DeeSCEcnaFrom0DeeEcha(i0StexEcha);
2964   }
2965 
2966   return n1StexStin;
2967 }
2968 
2969 Int_t TEcnaNumbering::Get0StexEchaFrom1StexStinAnd0StinEcha(const Int_t& n1StexStin, const Int_t& i0StinEcha) {
2970   // Electronic Channel number in Stex from Stin number in Stex

2971   // and from Electronic Channel number in Stin

2972 
2973   Int_t StexEcha = (Int_t)(-1.);
2974 
2975   if (n1StexStin > 0 && n1StexStin <= fEcal->MaxStinEcnaInStex() && i0StinEcha >= 0 &&
2976       i0StinEcha < fEcal->MaxCrysInStin()) {
2977     StexEcha = (n1StexStin - 1) * fEcal->MaxCrysInStin() + i0StinEcha;
2978   } else {
2979     std::cout << "!TEcnaNumbering::Get0StexEchaFrom1StexStinAnd0StinEcha *** ERROR ***> VALUE"
2980               << " OUT OF RANGE. Forced to -1. Argument values: n1StexStin = " << n1StexStin
2981               << ", channel = " << i0StinEcha << fTTBELL << std::endl;
2982   }
2983   return StexEcha;
2984 }
2985 Int_t TEcnaNumbering::Get1StexCrysFrom1StexStinAnd0StinEcha(const Int_t& n1StexStin,
2986                                                             const Int_t& i0StinEcha,
2987                                                             const Int_t& StexNumber) {
2988   // Crystal number in Stex from Stin number in Stex

2989   // and from Electronic Channel number in Stin (for the StexNumber_th Stex)

2990   // argument StexNumber used only in "EE" case

2991 
2992   Int_t n1StexCrys = (Int_t)0;
2993   if (fFlagSubDet == "EB") {
2994     n1StexCrys = Get1SMCrysFrom1SMTowAnd0TowEcha(n1StexStin, i0StinEcha);
2995   }
2996   if (fFlagSubDet == "EE") {
2997     TString sDeeDir = GetDeeDirViewedFromIP(StexNumber);
2998     n1StexCrys = Get1DeeCrysFrom1DeeSCEcnaAnd0SCEcha(n1StexStin, i0StinEcha, sDeeDir);
2999   }
3000 
3001   return n1StexCrys;
3002 }
3003 //===============================================================================

3004 //

3005 //        GetIHocoMin, GetIHocoMax, GetVecoMin, GetVecoMax

3006 //        GetJVecoMin, GetJVecoMax

3007 //

3008 //===============================================================================

3009 Double_t TEcnaNumbering::GetIHocoMin(const Int_t& Stex, const Int_t& StexStin) {
3010   Double_t IHocoMin = (Double_t)0.;
3011   if (fFlagSubDet == "EB") {
3012     IHocoMin = GetIEtaMin(Stex, StexStin);
3013   }
3014   if (fFlagSubDet == "EE") {
3015     IHocoMin = GetIIXMin(StexStin);
3016   }
3017   return IHocoMin;
3018 }
3019 
3020 Double_t TEcnaNumbering::GetIHocoMax(const Int_t& Stex, const Int_t& StexStin) {
3021   Double_t IHocoMax = (Double_t)0.;
3022   if (fFlagSubDet == "EB") {
3023     IHocoMax = GetIEtaMax(Stex, StexStin);
3024   }
3025   if (fFlagSubDet == "EE") {
3026     IHocoMax = GetIIXMax(StexStin);
3027   }
3028   return IHocoMax;
3029 }
3030 
3031 Double_t TEcnaNumbering::GetVecoMin(const Int_t& Stex, const Int_t& StexStin) {
3032   Double_t IVecoMin = (Double_t)0.;
3033   if (fFlagSubDet == "EB") {
3034     IVecoMin = GetPhiMin(Stex, StexStin);
3035   }
3036   if (fFlagSubDet == "EE") {
3037     IVecoMin = GetJIYMin(Stex, StexStin);
3038   }
3039   return IVecoMin;
3040 }
3041 
3042 Double_t TEcnaNumbering::GetVecoMax(const Int_t& Stex, const Int_t& StexStin) {
3043   Double_t IVecoMax = (Double_t)0.;
3044   if (fFlagSubDet == "EB") {
3045     IVecoMax = GetPhiMax(Stex, StexStin);
3046   }
3047   if (fFlagSubDet == "EE") {
3048     IVecoMax = GetJIYMax(Stex, StexStin);
3049   }
3050   return IVecoMax;
3051 }
3052 
3053 Double_t TEcnaNumbering::GetJVecoMin(const Int_t& Stex, const Int_t& StexStin) {
3054   Double_t JVecoMin = (Double_t)0.;
3055   if (fFlagSubDet == "EB") {
3056     JVecoMin = GetJPhiMin(Stex, StexStin);
3057   }
3058   if (fFlagSubDet == "EE") {
3059     JVecoMin = GetJIYMin(Stex, StexStin);
3060   }  // not used

3061   return JVecoMin;
3062 }
3063 Double_t TEcnaNumbering::GetJVecoMax(const Int_t& Stex, const Int_t& StexStin) {
3064   Double_t JVecoMax = (Double_t)0.;
3065   if (fFlagSubDet == "EB") {
3066     JVecoMax = GetJPhiMax(Stex, StexStin);
3067   }
3068   if (fFlagSubDet == "EE") {
3069     JVecoMax = GetJIYMax(Stex, StexStin);
3070   }  // not used

3071   return JVecoMax;
3072 }
3073 //===========================================================================

3074 //

3075 //                        GetStexHalfBarrel

3076 //

3077 //===========================================================================

3078 TString TEcnaNumbering::GetStexHalfStas(const Int_t& SMNumber) {
3079   TString half_stas = "EB? EE?";
3080 
3081   if (fFlagSubDet == "EB") {
3082     half_stas = GetSMHalfBarrel(SMNumber);
3083   }
3084   if (fFlagSubDet == "EE") {
3085     half_stas = GetDeeHalfEndcap(SMNumber);
3086   }
3087 
3088   return half_stas;
3089 }
3090 //===========================================================================

3091 //

3092 //                  GetSMFromFED, GetDSFromFED

3093 //

3094 //===========================================================================

3095 Int_t TEcnaNumbering::GetSMFromFED(const Int_t& FEDNumber) {
3096   Int_t EBSMNumber = 0;  // SM = Super Module

3097   if (FEDNumber >= 610 && FEDNumber <= 645) {
3098     EBSMNumber = FEDNumber - 609;
3099   }
3100   return EBSMNumber;
3101 }
3102 
3103 Int_t TEcnaNumbering::GetDSFromFED(const Int_t& FEDNumber) {
3104   Int_t EEDSNumber = 0;  // DS = Data Sector

3105 
3106   if (FEDNumber >= 600 && FEDNumber <= 609) {
3107     EEDSNumber = FEDNumber - 599;
3108   }
3109   if (FEDNumber >= 646 && FEDNumber <= 655) {
3110     EEDSNumber = FEDNumber - 645;
3111   }
3112 
3113   return EEDSNumber;
3114 }
3115 
3116 //--------------------------------------------

3117 Int_t TEcnaNumbering::MaxCrysInStinEcna(const Int_t& n1DeeNumber, const Int_t& n1DeeSCEcna, const TString& s_option) {
3118   // Number of Xtals in "Ecna SC" for not complete and not connected SC's.

3119   // Also valid for all connected and complete SC's and for towers of EB

3120 
3121   Int_t max_crys = fEcal->MaxCrysInStin();  // valid for EB and for connected and complete SC's of EE

3122 
3123   // Number of Xtals in SC for not complete and not connected SC's

3124 
3125   if (fFlagSubDet == "EE") {
3126     Int_t n_for_cons = GetDeeSCConsFrom1DeeSCEcna(n1DeeNumber, n1DeeSCEcna);
3127 
3128     //............ not complete SC's (inner border)

3129     if (n_for_cons == 12 || n_for_cons == 67 || n_for_cons == 75 || n_for_cons == 137 || n_for_cons == 161 ||
3130         n_for_cons == 216 || n_for_cons == 224 || n_for_cons == 286) {
3131       max_crys = 20;
3132     }
3133 
3134     //............ not connected SC's

3135     if ((n_for_cons == 182 || n_for_cons == 33) && (n1DeeSCEcna == 60 || n1DeeSCEcna == 119)) {
3136       max_crys = 10;
3137     }
3138 
3139     if ((n_for_cons == 178 || n_for_cons == 29) && (n1DeeSCEcna == 138 || n1DeeSCEcna == 157)) {
3140       max_crys = 10;
3141     }
3142     if ((n_for_cons == 207 || n_for_cons == 58) && (n1DeeSCEcna == 176 || n1DeeSCEcna == 193)) {
3143       max_crys = 10;
3144     }
3145 
3146     if ((n_for_cons == 298 || n_for_cons == 149) && (n1DeeSCEcna == 188)) {
3147       max_crys = 10;
3148     }
3149 
3150     if ((n_for_cons == 261 || n_for_cons == 112) && (n1DeeSCEcna == 144 || n1DeeSCEcna == 165)) {
3151       max_crys = 10;
3152     }
3153     if ((n_for_cons == 268 || n_for_cons == 119) && (n1DeeSCEcna == 102 || n1DeeSCEcna == 123)) {
3154       max_crys = 10;
3155     }
3156 
3157     if ((n_for_cons == 281 || n_for_cons == 132) && (n1DeeSCEcna == 41)) {
3158       max_crys = 10;
3159     }
3160 
3161     //............. not connected and mixed Ecna SC's

3162     if (s_option == "TEcnaRun" || s_option == "TEcnaRead") {
3163       if (s_option == "TEcnaRun") {
3164         // special translation of Xtal 11 of SCEcna 29 and 32 to respectively Xtal 11 of SCEcna 10 and 11

3165         if (n1DeeSCEcna == 29 || n1DeeSCEcna == 32) {
3166           max_crys = 5;
3167         }
3168         if (n1DeeSCEcna == 10 || n1DeeSCEcna == 11) {
3169           max_crys = 1;
3170         }
3171       }
3172       if (s_option == "TEcnaRead") {
3173         //if( n1DeeSCEcna == 29 || n1DeeSCEcna == 32 ){max_crys = 6;}

3174         if (n1DeeSCEcna == 29 || n1DeeSCEcna == 32) {
3175           max_crys = 5;
3176         }
3177         if (n1DeeSCEcna == 10 || n1DeeSCEcna == 11) {
3178           max_crys = 1;
3179         }
3180       }
3181     } else {
3182       std::cout << "!TEcnaNumbering::MaxCrysInStinEcna(...)> " << s_option << ": unknown option." << fTTBELL
3183                 << std::endl;
3184     }
3185   }
3186   return max_crys;
3187 }