Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:40

0001 #include <cmath>
0002 #include <iostream>
0003 #include <sstream>
0004 #include <cstdio>
0005 #include "L1Trigger/RPCTrigger/interface/RPCConst.h"
0006 #include "L1Trigger/RPCTrigger/interface/RPCException.h"
0007 
0008 #define m_pi 3.14159265358979
0009 
0010 int RPCConst::iptFromPt(const double pt) {
0011   if (pt == 0.)
0012     return 0;
0013   if (pt < m_pts[0]) {
0014     //edm::LogError("RPCTrigger")<<"** RPCConst ** iptFromPt called with illegal pt="<<pt;
0015     std::string msg = "[RPCConst::iptFromPt] called with illegal pt=";
0016     std::ostringstream ostr;
0017     ostr << pt;
0018     msg += ostr.str();
0019     throw RPCException(msg);
0020     return 0;
0021   }
0022   int ipt = RPCConst::IPT_MAX;
0023   while (pt < m_pts[ipt]) {
0024     ipt--;
0025   };
0026   return ipt;
0027 }
0028 
0029 double RPCConst::ptFromIpt(const int ipt) {
0030   if (ipt < 0 || ipt > RPCConst::IPT_MAX) {
0031     //edm::LogError("RPCTrigger") <<"**RPCConst::ptFromIpt** problem with ipt: "<<ipt;
0032     std::string msg = "[RPCConst::ptFromIpt] problem with ipt: ";
0033     std::ostringstream ostr;
0034     ostr << ipt;
0035     msg += ostr.str();
0036     throw RPCException(msg);
0037     return 0.;
0038   } else
0039     return m_pts[ipt];
0040 }
0041 
0042 double RPCConst::etaFromTowerNum(const int atower) {
0043   int iabsitow = (atower >= 0) ? atower : -atower;
0044   if (0 == iabsitow)
0045     return 0.;
0046   if (iabsitow > RPCConst::ITOW_MAX) {
0047     //edm::LogError("RPCTrigger") << "**RPCConst::etaFromTowerNum** iabsitow>ITOW_MAX for m_tower:"
0048     //     << atower ;
0049     std::string msg = "[RPCConst::etaFromTowerNum] iabsitow>ITOW_MAX for m_tower:";
0050     std::ostringstream ostr;
0051     ostr << atower;
0052     msg += ostr.str();
0053     throw RPCException(msg);
0054     return 0.;
0055   }
0056   double eta = (m_etas[iabsitow] + m_etas[iabsitow + 1]) / 2.;
0057   return (atower >= 0) ? eta : -eta;
0058 }
0059 
0060 int RPCConst::towerNumFromEta(const double eta) {
0061   int m_tower = 0;
0062   double abseta = (eta >= 0.) ? eta : -eta;
0063   while (m_tower <= ITOW_MAX) {
0064     if (m_etas[m_tower] <= abseta && abseta < m_etas[m_tower + 1])
0065       break;
0066     m_tower++;
0067   }
0068   if (m_tower > ITOW_MAX)
0069     m_tower = ITOW_MAX;
0070   return (eta >= 0) ? m_tower : -m_tower;
0071 }
0072 
0073 double RPCConst::phiFromSegmentNum(const int iseg) {
0074   double phi = OFFSET + 2. * m_pi * (iseg) / (double)RPCConst::NSEG;
0075   return (phi < 2. * m_pi) ? phi : phi - 2. * m_pi;
0076 }
0077 
0078 double RPCConst::phiFromLogSegSec(const int logSegment, const int logSector) {
0079   int iseg = logSegment * 12 + logSector;
0080   double phi = OFFSET + 2. * m_pi * (iseg) / (double)RPCConst::NSEG;
0081   return (phi < 2. * m_pi) ? phi : phi - 2. * m_pi;
0082 }
0083 
0084 int RPCConst::segmentNumFromPhi(const double phi) {
0085   double iphi;
0086   if (phi - OFFSET < 0) {
0087     iphi = 2 * m_pi + phi;
0088   } else {
0089     iphi = phi - OFFSET;
0090   }
0091   int iseg = (int)(iphi * RPCConst::NSEG / (2. * m_pi));
0092   return iseg;
0093 }
0094 
0095 /*
0096 int RPCConst::checkBarrel(const int atower) {
0097   int iabsitow = (atower >= 0)? atower: -atower;
0098   if(iabsitow <= RPCConst::ITOW_MAX_LOWPT) {
0099     return 1;
0100   } else if (iabsitow <= RPCConst::ITOW_MAX) {
0101     return 0;
0102   }
0103   return -1;
0104 } */
0105 
0106 double RPCConst::vxMuRate(int ptCode) {
0107   double pt_ev = RPCConst::ptFromIpt(ptCode);
0108   if (pt_ev == 0)
0109     return 0.0;
0110   const double lum = 2.0e33;  //defoult is 1.0e34;
0111   const double dabseta = 1.0;
0112   const double dpt = 1.0;
0113   const double afactor = 1.0e-34 * lum * dabseta * dpt;
0114   const double a = 2 * 1.3084E6;
0115   const double mu = -0.725;
0116   const double sigma = 0.4333;
0117   const double s2 = 2 * sigma * sigma;
0118 
0119   double ptlog10;
0120   ptlog10 = log10(pt_ev);
0121   double ex = (ptlog10 - mu) * (ptlog10 - mu) / s2;
0122   double rate = (a * exp(-ex) * afactor);
0123 
0124   //edm::LogError("RPCTrigger")<<ptCode<<" "<<rate;//<<<<<<<<<<<<<<<<<<<<<<<<
0125   return rate;
0126 }
0127 
0128 //muon rate for pt from ptCode to ptCode + 1
0129 //i.e for ptCode bin
0130 double RPCConst::vxIntegMuRate(int ptCode, double etaFrom, double etaTo) {
0131   //calkowanie metoda trapezow - nie do konca dobre
0132   double rate =
0133       0.5 * (vxMuRate(ptCode) + vxMuRate(ptCode + 1)) * (RPCConst::ptFromIpt(ptCode + 1) - RPCConst::ptFromIpt(ptCode));
0134 
0135   rate = rate * (etaTo - etaFrom);
0136 
0137   //edm::LogError("RPCTrigger")<<ptCode<<" "<<rate;//<<<<<<<<<<<<<<<<<<<<<<<<
0138   return rate;
0139 }
0140 
0141 //muon rate for pt from ptCode to ptCode + 1 in a given m_tower - only one!!! (mutliply by 2 to take oalso negative!!!)
0142 //i.e for ptCode bin
0143 double RPCConst::vxIntegMuRate(int ptCode, int m_tower) {
0144   //calkowanie metoda trapezow - nie do konca dobre
0145   double rate = vxIntegMuRate(ptCode, RPCConst::m_etas[abs(m_tower)], RPCConst::m_etas[abs(m_tower) + 1]);
0146 
0147   //edm::LogError("RPCTrigger")<<ptCode<<" "<<rate;//<<<<<<<<<<<<<<<<<<<<<<<<
0148   return rate;
0149 }
0150 
0151 /*
0152 const int RPCConst::IPT_THRESHOLD [2][RPCConst::ITOW_MAX+1]={
0153 //0   1   2   3   4     5   6   7   8   9    10  11  12  13  14  15  16  m_Tower
0154 {17, 17, 17, 17, 17,   16, 16, 15, 17, 14,   12, 11, 12, 17, 16, 15, 15}, //LOW
0155 {12, 12, 12, 12, 12,   11, 8,  11, 12, 9,    9,  8,  7,  11, 11, 11, 11} //VLOW
0156 };
0157 */
0158 
0159 const double RPCConst::m_pts[RPCConst::IPT_MAX + 1] = {
0160     0.0, 0.01,  //<<<<<<<<<<<<<<<<<dla ptCode = 1 bylo 0, ale to powoduje problemy w vxMuRate
0161     1.5, 2.0,  2.5, 3.0, 3.5, 4.0, 4.5, 5.,  6.,  7.,  8.,  10., 12.,  14.,  16.,
0162     18., 20.,  25., 30., 35., 40., 45., 50., 60., 70., 80., 90., 100., 120., 140.};
0163 
0164 // m_etas contain approximate lower egges of eta towers
0165 // 0:ITOW_MAX  +additionaly upper edge  of last m_tower
0166 const double RPCConst::m_etas[RPCConst::ITOW_MAX + 2] = {
0167     0.00, 0.07, 0.27, 0.44, 0.58, 0.72, 0.83, 0.93, 1.04, 1.14, 1.24, 1.36, 1.48, 1.61, 1.73, 1.85, 1.97, 2.10};
0168 
0169 // imported constants
0170 
0171 const std::string RPCConst::m_LOGPLANE_STR[RPCConst::m_LOGPLANES_COUNT] = {
0172     "m_LOGPLANE1", "m_LOGPLANE2", "m_LOGPLANE3", "m_LOGPLANE4", "m_LOGPLANE5", "m_LOGPLANE6"};
0173 
0174 const unsigned int RPCConst::m_LOGPLANE_SIZE[m_TOWER_COUNT][m_LOGPLANES_COUNT] = {
0175     //LOGPLANE  1,  2,  3   4   5   6
0176     {72, 56, 8, 40, 40, 24},  //TOWER 0
0177     {72, 56, 8, 40, 40, 24},  //TOWER 1
0178     {72, 56, 8, 40, 40, 24},  //TOWER 2
0179     {72, 56, 8, 40, 40, 24},  //TOWER 3
0180     {72, 56, 8, 40, 40, 24},  //TOWER 4
0181     {72, 56, 40, 8, 40, 24},  //TOWER 5
0182     {56, 72, 40, 8, 24, 0},   //TOWER 6
0183     {72, 56, 40, 8, 24, 0},   //TOWER 7
0184     {72, 24, 40, 8, 0, 0},    //TOWER 8
0185     {72, 8, 40, 0, 0, 0},     //TOWER 9
0186     {72, 8, 40, 24, 0, 0},    //TOWER 10
0187     {72, 8, 40, 24, 0, 0},    //TOWER 11
0188     {72, 8, 40, 24, 0, 0},    //TOWER 12
0189     {72, 8, 40, 24, 0, 0},    //TOWER 13
0190     {72, 8, 40, 24, 0, 0},    //TOWER 14
0191     {72, 8, 40, 24, 0, 0},    //TOWER 15
0192     {72, 8, 40, 24, 0, 0}     //TOWER 16
0193 };
0194 
0195 const int RPCConst::m_VLPT_PLANES_COUNT[m_TOWER_COUNT] = {4, 4, 4, 4, 4, 4, 4, 4, 3, 3, 3, 3, 3, 3, 3, 3, 3};
0196 
0197 const int RPCConst::m_USED_PLANES_COUNT[m_TOWER_COUNT] = {
0198     //0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16
0199     6,
0200     6,
0201     6,
0202     6,
0203     6,
0204     6,
0205     5,
0206     5,
0207     4,
0208     3,
0209     4,
0210     4,
0211     4,
0212     4,
0213     4,
0214     4,
0215     4};
0216 
0217 const int RPCConst::m_REF_PLANE[m_TOWER_COUNT] = {
0218     //     0,         1,         2,         3,         4,
0219     m_LOGPLANE3,
0220     m_LOGPLANE3,
0221     m_LOGPLANE3,
0222     m_LOGPLANE3,
0223     m_LOGPLANE3,
0224     //     5,         6,         7,         8,
0225     m_LOGPLANE4,
0226     m_LOGPLANE4,
0227     m_LOGPLANE4,
0228     m_LOGPLANE4,
0229     //     9,         10,       11,        12,        13,        14,         15,        16,
0230     m_LOGPLANE2,
0231     m_LOGPLANE2,
0232     m_LOGPLANE2,
0233     m_LOGPLANE2,
0234     m_LOGPLANE2,
0235     m_LOGPLANE2,
0236     m_LOGPLANE2,
0237     m_LOGPLANE2};
0238 
0239 /*    
0240     const int m_PT_CODE_MAX = 31; //!< Pt_code range = 0-m_PT_CODE_MAX
0241     
0242     const int m_LOGPLANE1 = 0; //!< The Logic Planes are named starting from '1', but in varoius loop indeks are from '0', that's why always use these consts 
0243     const int m_LOGPLANE2 = 1;
0244     const int m_LOGPLANE3 = 2;
0245     const int m_LOGPLANE4 = 3;
0246     const int m_LOGPLANE5 = 4;
0247     const int m_LOGPLANE6 = 5;
0248     
0249     const int m_FIRST_PLANE = m_LOGPLANE1; //!< Use ase a first index in loops.
0250     const int m_LAST_PLANE  = m_LOGPLANE6; //!< Use ase a last index in loops.
0251 */
0252 
0253 //------- imported fucntions
0254 
0255 int RPCConst::stringToInt(std::string str) {
0256   for (unsigned int i = 0; i < str.size(); i++)
0257     if (str[i] < '0' || str[i] > '9')
0258       throw RPCException("Error in stringToInt(): the string cannot be converted to a number");
0259   //edm::LogError("RPCTrigger")<< "Error in stringToInt(): the string cannot be converted to a number";
0260   return atoi(str.c_str());
0261 }
0262 
0263 //inline
0264 std::string RPCConst::intToString(int number) {
0265   std::string str;
0266   /* Some problems. AK
0267   std::ostringstream ostr;
0268   ostr<<number;
0269   str = ostr.str();
0270   edm::LogError("RPCTrigger")<<"std::string intToString(int number)";
0271   edm::LogError("RPCTrigger")<<str;
0272   */
0273   char tmp[20];
0274   sprintf(tmp, "%d", number);
0275   str.append(tmp);
0276   return str;
0277 }
0278 
0279 bool RPCConst::l1RpcConeCrdnts::operator<(const l1RpcConeCrdnts& cone) const {
0280   if (m_Tower != cone.m_Tower)
0281     return (m_Tower < cone.m_Tower);
0282   if (m_LogSector != cone.m_LogSector)
0283     return (m_LogSector < cone.m_LogSector);
0284   if (m_LogSegment != cone.m_LogSegment)
0285     return (m_LogSegment < cone.m_LogSegment);
0286 
0287   return false;
0288 }
0289 
0290 bool RPCConst::l1RpcConeCrdnts::operator==(const l1RpcConeCrdnts& cone) const {
0291   if (m_Tower != cone.m_Tower)
0292     return false;
0293   if (m_LogSector != cone.m_LogSector)
0294     return false;
0295   if (m_LogSegment != cone.m_LogSegment)
0296     return false;
0297 
0298   return true;
0299 }