CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/L1Trigger/CSCTrackFinder/src/CSCSectorReceiverLUT.cc

Go to the documentation of this file.
00001 #include <L1Trigger/CSCTrackFinder/interface/CSCSectorReceiverLUT.h>
00002 #include <L1Trigger/CSCTrackFinder/interface/CSCSectorReceiverMiniLUT.h>
00003 #include <L1Trigger/CSCCommonTrigger/interface/CSCTriggerGeometry.h>
00004 #include <L1Trigger/CSCCommonTrigger/interface/CSCTriggerGeomManager.h>
00005 #include <L1Trigger/CSCCommonTrigger/interface/CSCPatternLUT.h>
00006 #include <L1Trigger/CSCCommonTrigger/interface/CSCFrontRearLUT.h>
00007 #include <DataFormats/L1CSCTrackFinder/interface/CSCBitWidths.h>
00008 #include <DataFormats/L1CSCTrackFinder/interface/CSCTFConstants.h>
00009 #include <L1Trigger/CSCCommonTrigger/interface/CSCConstants.h>
00010 
00011 #include <Geometry/CSCGeometry/interface/CSCLayerGeometry.h>
00012 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00013 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00014 
00015 #include <DataFormats/MuonDetId/interface/CSCTriggerNumbering.h>
00016 
00017 #include <FWCore/MessageLogger/interface/MessageLogger.h>
00018 
00019 #include <fstream>
00020 #include <cstring>
00021 
00022 lclphidat* CSCSectorReceiverLUT::me_lcl_phi = NULL;
00023 bool CSCSectorReceiverLUT::me_lcl_phi_loaded = false;
00024 
00025 
00026 CSCSectorReceiverLUT::CSCSectorReceiverLUT(int endcap, int sector, int subsector, int station,
00027                                            const edm::ParameterSet & pset, bool TMB07):_endcap(endcap),_sector(sector),
00028                                                                            _subsector(subsector),
00029                                                                            _station(station),isTMB07(TMB07)
00030 {
00031   LUTsFromFile = pset.getUntrackedParameter<bool>("ReadLUTs",false);
00032   useMiniLUTs = pset.getUntrackedParameter<bool>("UseMiniLUTs", true);
00033   isBinary = pset.getUntrackedParameter<bool>("Binary",false);
00034 
00035   me_global_eta = NULL;
00036   me_global_phi = NULL;
00037   mb_global_phi = NULL;
00038   if(LUTsFromFile && !useMiniLUTs)
00039     {
00040       me_lcl_phi_file = pset.getUntrackedParameter<edm::FileInPath>("LocalPhiLUT", edm::FileInPath(std::string("L1Trigger/CSCTrackFinder/LUTs/LocalPhiLUT"
00041                                                                                                                + (isBinary ? std::string(".bin") : std::string(".dat")))));
00042       me_gbl_phi_file = pset.getUntrackedParameter<edm::FileInPath>("GlobalPhiLUTME", edm::FileInPath((std::string("L1Trigger/CSCTrackFinder/LUTs/GlobalPhiME")
00043                                                                                                        + encodeFileIndex()
00044                                                                                                        + (isBinary ? std::string(".bin") : std::string(".dat")))));
00045       if(station == 1)
00046         mb_gbl_phi_file = pset.getUntrackedParameter<edm::FileInPath>("GlobalPhiLUTMB", edm::FileInPath((std::string("L1Trigger/CSCTrackFinder/LUTs/GlobalPhiMB")
00047                                                                                                          + encodeFileIndex()
00048                                                                                                          + (isBinary ? std::string(".bin") : std::string(".dat")))));
00049       me_gbl_eta_file = pset.getUntrackedParameter<edm::FileInPath>("GlobalEtaLUTME", edm::FileInPath((std::string("L1Trigger/CSCTrackFinder/LUTs/GlobalEtaME")
00050                                                                                                        + encodeFileIndex()
00051                                                                                                        + (isBinary ? std::string(".bin") : std::string(".dat")))));
00052       readLUTsFromFile();
00053     }
00054 
00055 }
00056 
00057 CSCSectorReceiverLUT::CSCSectorReceiverLUT(const CSCSectorReceiverLUT& lut):_endcap(lut._endcap),
00058                                                                             _sector(lut._sector),
00059                                                                             _subsector(lut._subsector),
00060                                                                             _station(lut._station),
00061                                                                             me_lcl_phi_file(lut.me_lcl_phi_file),
00062                                                                             me_gbl_phi_file(lut.me_gbl_phi_file),
00063                                                                             mb_gbl_phi_file(lut.mb_gbl_phi_file),
00064                                                                             me_gbl_eta_file(lut.me_gbl_eta_file),
00065                                                                             LUTsFromFile(lut.LUTsFromFile),
00066                                                                             isBinary(lut.isBinary)
00067 {
00068   if(lut.mb_global_phi)
00069     {
00070       mb_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00071       memcpy(mb_global_phi, lut.mb_global_phi, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(gblphidat));
00072     }
00073   else mb_global_phi = NULL;
00074   if(lut.me_global_phi)
00075     {
00076       me_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00077       memcpy(me_global_phi, lut.me_global_phi, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(gblphidat));
00078     }
00079   else me_global_phi = NULL;
00080   if(lut.me_global_eta)
00081     {
00082       me_global_eta = new gbletadat[1<<CSCBitWidths::kGlobalEtaAddressWidth];
00083       memcpy(me_global_eta, lut.me_global_eta, (1<<CSCBitWidths::kGlobalEtaAddressWidth)*sizeof(gbletadat));
00084     }
00085   else me_global_eta = NULL;
00086 }
00087 
00088 CSCSectorReceiverLUT& CSCSectorReceiverLUT::operator=(const CSCSectorReceiverLUT& lut)
00089 {
00090   if(this != &lut)
00091     {
00092       _endcap = lut._endcap;
00093       _sector = lut._sector;
00094       _subsector = lut._subsector;
00095       _station = lut._station;
00096       me_lcl_phi_file = lut.me_lcl_phi_file;
00097       me_gbl_phi_file = lut.me_gbl_phi_file;
00098       mb_gbl_phi_file = lut.mb_gbl_phi_file;
00099       me_gbl_eta_file = lut.me_gbl_eta_file;
00100       LUTsFromFile = lut.LUTsFromFile;
00101       isBinary = lut.isBinary;
00102 
00103       if(lut.mb_global_phi)
00104         {
00105           mb_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00106           memcpy(mb_global_phi, lut.mb_global_phi, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(gblphidat));
00107         }
00108       else mb_global_phi = NULL;
00109 
00110       if(lut.me_global_phi)
00111         {
00112           me_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00113           memcpy(me_global_phi, lut.me_global_phi, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(gblphidat));
00114         }
00115       else me_global_phi = NULL;
00116 
00117       if(lut.me_global_eta)
00118         {
00119           me_global_eta = new gbletadat[1<<CSCBitWidths::kGlobalEtaAddressWidth];
00120           memcpy(me_global_eta, lut.me_global_eta, (1<<CSCBitWidths::kGlobalEtaAddressWidth)*sizeof(gbletadat));
00121         }
00122       else me_global_eta = NULL;
00123     }
00124   return *this;
00125 }
00126 
00127 CSCSectorReceiverLUT::~CSCSectorReceiverLUT()
00128 {
00129   if(me_lcl_phi_loaded)
00130     {
00131       delete me_lcl_phi;
00132       me_lcl_phi = NULL;
00133       me_lcl_phi_loaded = false;
00134     }
00135   if(me_global_eta)
00136     {
00137       delete me_global_eta;
00138       me_global_eta = NULL;
00139     }
00140   if(me_global_phi)
00141     {
00142       delete me_global_phi;
00143       me_global_phi = NULL;
00144     }
00145   if(mb_global_phi)
00146     {
00147       delete mb_global_phi;
00148       mb_global_phi = NULL;
00149     }
00150 }
00151 
00152 lclphidat CSCSectorReceiverLUT::calcLocalPhi(const lclphiadd& theadd) const
00153 {
00154   lclphidat data;
00155 
00156   static int maxPhiL = 1<<CSCBitWidths::kLocalPhiDataBitWidth;
00157   double binPhiL = static_cast<double>(maxPhiL)/(2.*CSCConstants::MAX_NUM_STRIPS);
00158 
00159   memset(&data,0,sizeof(lclphidat));
00160 
00161   double patternOffset;
00162 
00163   if(isTMB07) patternOffset = CSCPatternLUT::get2007Position((theadd.pattern_type<<3) + theadd.clct_pattern);
00164   else patternOffset = CSCPatternLUT::getPosition(theadd.clct_pattern);
00165 
00166   // The phiL value stored is for the center of the half-/di-strip.
00167   if(theadd.strip < 2*CSCConstants::MAX_NUM_STRIPS)
00168     if(theadd.pattern_type == 1 || isTMB07) // if halfstrip (Note: no distrips in TMB 2007 patterns)
00169       data.phi_local = static_cast<unsigned>((0.5 + theadd.strip + patternOffset)*binPhiL);
00170     else // if distrip
00171       data.phi_local = static_cast<unsigned>((2 + theadd.strip + 4.*patternOffset)*binPhiL);
00172   else {
00173     throw cms::Exception("CSCSectorReceiverLUT")
00174       << "+++ Value of strip, " << theadd.strip
00175       << ", exceeds max allowed, " << 2*CSCConstants::MAX_NUM_STRIPS-1
00176       << " +++\n";
00177   }
00178 
00179   if (data.phi_local >= maxPhiL) {
00180     throw cms::Exception("CSCSectorReceiverLUT")
00181       << "+++ Value of phi_local, " << data.phi_local
00182       << ", exceeds max allowed, " << maxPhiL-1 << " +++\n";
00183   }
00184 
00185   LogDebug("CSCSectorReceiver")
00186     << "endcap = " << _endcap << " station = " << _station
00187     << " maxPhiL = " << maxPhiL << " binPhiL = " << binPhiL;
00188   LogDebug("CSCSectorReceiver")
00189     << "strip # " << theadd.strip << " hs/ds = " << theadd.pattern_type
00190     << " pattern = " << theadd.clct_pattern << " offset = " << patternOffset
00191     << " phi_local = " << data.phi_local;
00192 
00194   data.phi_bend_local = 0;
00195 
00196   return data; //return LUT result
00197 }
00198 
00199 
00200 void CSCSectorReceiverLUT::fillLocalPhiLUT()
00201 {
00202   // read data in from a file... Add this later.
00203 }
00204 
00205 lclphidat CSCSectorReceiverLUT::localPhi(const int strip, const int pattern,
00206                                          const int quality, const int lr) const
00207 {
00208   lclphiadd theadd;
00209 
00210   theadd.strip = strip;
00211   theadd.clct_pattern = pattern & 0x7;
00212   theadd.pattern_type = (pattern & 0x8) >> 3;
00213   theadd.quality = quality;
00214   theadd.lr = lr;
00215   theadd.spare = 0;
00216 
00217   return localPhi(theadd);
00218 }
00219 
00220 lclphidat CSCSectorReceiverLUT::localPhi(unsigned address) const
00221 {
00222   lclphidat result;
00223   lclphiadd theadd(address);
00224   
00225   if(useMiniLUTs && isTMB07)
00226     {
00227       result = CSCSectorReceiverMiniLUT::calcLocalPhiMini(address);
00228     }
00229   else if(LUTsFromFile) result = me_lcl_phi[address];
00230   else result = calcLocalPhi(theadd);
00231 
00232   return result;
00233 }
00234 
00235 lclphidat CSCSectorReceiverLUT::localPhi(lclphiadd address) const
00236 {
00237   lclphidat result;
00238 
00239   if(useMiniLUTs && isTMB07)
00240     {
00241       result = CSCSectorReceiverMiniLUT::calcLocalPhiMini(address.toint()); 
00242     }
00243   else if(LUTsFromFile) result = me_lcl_phi[address.toint()];
00244   else result = calcLocalPhi(address);
00245 
00246   return result;
00247 }
00248 
00249 double CSCSectorReceiverLUT::getGlobalPhiValue(const CSCLayer* thelayer, const unsigned& strip, const unsigned& wire_group) const
00250 {
00251   double result = 0.0;
00252   //CSCLayerGeometry* thegeom;
00253   //LocalPoint lp;
00254   //GlobalPoint gp;
00255 
00256   try
00257     {
00258       //thegeom = const_cast<CSCLayerGeometry*>(thelayer->geometry());
00259       //lp = thegeom->stripWireGroupIntersection(strip, wire_group);
00260       //gp = thelayer->surface().toGlobal(lp);
00261       result = thelayer->centerOfStrip(strip).phi();//gp.phi();
00262 
00263       if (result < 0.) result += 2.*M_PI;
00264     }
00265   catch(edm::Exception& e)
00266     {
00267       LogDebug("CSCSectorReceiverLUT|getGlobalPhiValue") << e.what();
00268     }
00269 
00270   return result;
00271 }
00272 
00273 gblphidat CSCSectorReceiverLUT::calcGlobalPhiME(const gblphiadd& address) const
00274 {
00275   gblphidat result(0);
00276   CSCTriggerGeomManager* thegeom = CSCTriggerGeometry::get();
00277   CSCChamber* thechamber = NULL;
00278   CSCLayer* thelayer = NULL;
00279   CSCLayerGeometry* layergeom = NULL;
00280   int cscid = address.cscid;
00281   unsigned wire_group = address.wire_group;
00282   unsigned local_phi = address.phi_local;
00283   const double sectorOffset = (CSCTFConstants::SECTOR1_CENT_RAD-CSCTFConstants::SECTOR_RAD/2.) + (_sector-1)*M_PI/3.;
00284 
00285   //Number of global phi units per radian.
00286   static int maxPhiG = 1<<CSCBitWidths::kGlobalPhiDataBitWidth;
00287   double binPhiG = static_cast<double>(maxPhiG)/CSCTFConstants::SECTOR_RAD;
00288 
00289   // We will use these to convert the local phi into radians.
00290   static unsigned int maxPhiL = 1<<CSCBitWidths::kLocalPhiDataBitWidth;
00291   const double binPhiL = static_cast<double>(maxPhiL)/(2.*CSCConstants::MAX_NUM_STRIPS);
00292 
00293   if(cscid < CSCTriggerNumbering::minTriggerCscId())
00294     {
00295       edm::LogWarning("CSCSectorReceiverLUT|getGlobalPhiValue")
00296         << " warning: cscId " << cscid << " is out of bounds ["
00297         << CSCTriggerNumbering::maxTriggerCscId() << "-"
00298         << CSCTriggerNumbering::maxTriggerCscId() << "]\n";
00299       throw cms::Exception("CSCSectorReceiverLUT")
00300         << "+++ Value of CSC ID, " << cscid
00301         << ", is out of bounds [" << CSCTriggerNumbering::minTriggerCscId() << "-"
00302         << CSCTriggerNumbering::maxTriggerCscId() << "] +++\n";
00303     }
00304 
00305   if(cscid > CSCTriggerNumbering::maxTriggerCscId())
00306     {
00307       edm::LogWarning("CSCSectorReceiverLUT|getGlobalPhiValue")
00308         << " warning: cscId " << cscid << " is out of bounds ["
00309         << CSCTriggerNumbering::maxTriggerCscId() << "-"
00310         << CSCTriggerNumbering::maxTriggerCscId() << "]\n";
00311       throw cms::Exception("CSCSectorReceiverLUT")
00312         << "+++ Value of CSC ID, " << cscid
00313         << ", is out of bounds [" << CSCTriggerNumbering::minTriggerCscId() << "-"
00314         << CSCTriggerNumbering::maxTriggerCscId() << "] +++\n";
00315     }
00316 
00317   if(wire_group >= 1<<5)
00318     {
00319       edm::LogWarning("CSCSectorReceiverLUT|getGlobalPhiValue")
00320         << "warning: wire_group" << wire_group
00321         << " is out of bounds (1-" << ((1<<5)-1) << "]\n";
00322       throw cms::Exception("CSCSectorReceiverLUT")
00323         << "+++ Value of wire_group, " << wire_group
00324         << ", is out of bounds (1-" << ((1<<5)-1) << "] +++\n";
00325     }
00326 
00327   if(local_phi >= maxPhiL)
00328     {
00329       edm::LogWarning("CSCSectorReceiverLUT|getGlobalPhiValue")
00330         << "warning: local_phi" << local_phi
00331         << " is out of bounds [0-" << maxPhiL << ")\n";
00332       throw cms::Exception("CSCSectorReceiverLUT")
00333         << "+++ Value of local_phi, " << local_phi
00334         << ", is out of bounds [0-, " << maxPhiL << ") +++\n";
00335     }
00336 
00337   try
00338     {
00339       thechamber = thegeom->chamber(_endcap,_station,_sector,_subsector,cscid);
00340       if(thechamber)
00341         {
00342           if(isTMB07)
00343             {
00344               layergeom = const_cast<CSCLayerGeometry*>(thechamber->layer(CSCConstants::KEY_CLCT_LAYER)->geometry());
00345               thelayer = const_cast<CSCLayer*>(thechamber->layer(CSCConstants::KEY_CLCT_LAYER));
00346             }
00347           else
00348             {
00349               layergeom = const_cast<CSCLayerGeometry*>(thechamber->layer(CSCConstants::KEY_CLCT_LAYER_PRE_TMB07)->geometry());
00350               thelayer = const_cast<CSCLayer*>(thechamber->layer(CSCConstants::KEY_CLCT_LAYER_PRE_TMB07));
00351             }
00352           const int nStrips = layergeom->numberOfStrips();
00353           // PhiL is the strip number converted into some units between 0 and
00354           // 1023.  When we did the conversion in fillLocalPhiTable(), we did
00355           // not know for which chamber we do it (and, therefore, how many strips
00356           // it has), and always used the maximum possible number of strips
00357           // per chamber, MAX_NUM_STRIPS=80.  Now, since we know the chamber id
00358           // and how many strips the chamber has, we can re-adjust the scale.
00359           //const double scale = static_cast<double>(CSCConstants::MAX_NUM_STRIPS)/nStrips;
00360 
00361           int strip = 0, halfstrip = 0;
00362 
00363           halfstrip = static_cast<int>(local_phi/binPhiL);
00364           strip     = halfstrip/2;
00365 
00366           // Find the phi width of the chamber and the position of its "left"
00367           // (lower phi) edge (both in radians).
00368           // Phi positions of the centers of the first and of the last strips
00369           // in the chamber.
00370           const double phi_f = getGlobalPhiValue(thelayer, 1, wire_group);
00371           const double phi_l = getGlobalPhiValue(thelayer, nStrips, wire_group);
00372           // Phi widths of the half-strips at both ends of the chamber;
00373           // surprisingly, they are not the same.
00374           const double hsWidth_f = fabs(getGlobalPhiValue(thelayer, 2, wire_group) - phi_f)/2.;
00375           const double hsWidth_l = fabs(phi_l - getGlobalPhiValue(thelayer, nStrips - 1, wire_group))/2.;
00376 
00377           // The "natural" match between the strips and phi values -- when
00378           // a larger strip number corresponds to a larger phi value, i.e. strips
00379           // are counted clockwise if we look at them from the inside of the
00380           // detector -- is reversed for some stations.  At the moment, these
00381           // are stations 3 and 4 of the 1st endcap, and stations 1 and 2 of
00382           // the 2nd endcap.  Instead of using
00383           // if ((theEndcap == 1 && theStation <= 2) ||
00384           // (theEndcap == 2 && theStation >= 3)),
00385           // we get the order from the phi values of the first and the last strip
00386           // in a chamber, just in case the counting scheme changes in the future.
00387           // Once we know how the strips are counted, we can go from the middle
00388           // of the strips to their outer edges.
00389           bool   clockwiseOrder;
00390           double leftEdge, rightEdge;
00391           if (fabs(phi_f - phi_l) < M_PI)
00392             {
00393               if (phi_f < phi_l) clockwiseOrder = true;
00394               else clockwiseOrder = false;
00395             }
00396           else
00397             { // the chamber crosses the phi = pi boundary
00398               if (phi_f < phi_l) clockwiseOrder = false;
00399               else clockwiseOrder = true;
00400             }
00401           if (clockwiseOrder)
00402             {
00403               leftEdge  = phi_f - hsWidth_f;
00404               rightEdge = phi_l + hsWidth_l;
00405             }
00406           else
00407             {
00408               leftEdge  = phi_l - hsWidth_l;
00409               rightEdge = phi_f + hsWidth_f;
00410             }
00411           if (fabs(phi_f - phi_l) >= M_PI) {rightEdge += 2.*M_PI;}
00412           //double chamberWidth = (rightEdge - leftEdge);
00413 
00414           // Chamber offset, relative to the edge of the sector.
00415           //double chamberOffset = leftEdge - sectorOffset;
00416           //if (chamberOffset < -M_PI) chamberOffset += 2*M_PI;
00417 
00418           double temp_phi = 0.0, strip_phi = 0.0, delta_phi = 0.0;
00419           double distFromHalfStripCenter = 0.0, halfstripWidth = 0.0;
00420 
00421           if (strip < nStrips)
00422             {
00423               // Approximate distance from the center of the half-strip to the center
00424               // of this phil bin, in units of half-strip width.
00425               distFromHalfStripCenter = (local_phi+0.5)/binPhiL - halfstrip - 0.5;
00426               // Half-strip width (in rad), calculated as the half-distance between
00427               // the adjacent strips.  Since in the current ORCA implementation
00428               // the half-strip width changes from strip to strip, base the choice
00429               // of the adjacent strip on the half-strip number.
00430               if ((halfstrip%2 == 0 && halfstrip != 0) || halfstrip == 2*nStrips-1) {
00431                 halfstripWidth =
00432                   fabs(getGlobalPhiValue(thelayer, strip+1, wire_group) - getGlobalPhiValue(thelayer, strip, wire_group)) / 2.;
00433               }
00434               else
00435                 {
00436                   halfstripWidth =
00437                     fabs(getGlobalPhiValue(thelayer, strip+1, wire_group) - getGlobalPhiValue(thelayer, strip+2, wire_group)) / 2.;
00438                 }
00439               // Correction for the strips crossing the 180 degree boundary.
00440               if (halfstripWidth > M_PI/2.) halfstripWidth = M_PI - halfstripWidth;
00441               // Phi at the center of the strip.
00442               strip_phi = getGlobalPhiValue(thelayer, strip+1, wire_group);
00443               // Distance between the center of the strip and the phil position.
00444               delta_phi = halfstripWidth*(((halfstrip%2)-0.5)+distFromHalfStripCenter);
00445               if (clockwiseOrder)
00446                 temp_phi = strip_phi+ delta_phi;
00447               else
00448                 temp_phi = strip_phi- delta_phi;
00449             }
00450           else
00451             {
00452               // PhiL values that do not have corresponding strips (the chamber
00453               // has less than 80 strips assumed in fillLocalPhi).  It does not
00454               // really matter what we do with these values; at the moment, just
00455               // set them to the phis of the edges of the chamber.
00456               if (clockwiseOrder) temp_phi = rightEdge;
00457               else temp_phi = leftEdge;
00458             }
00459 
00460           // Finally, subtract the sector offset and convert to the scale of
00461           // the global phi.
00462 
00463           temp_phi -= sectorOffset;
00464 
00465           if (temp_phi < 0.) temp_phi += 2.*M_PI;
00466 
00467           temp_phi *= binPhiG;
00468 
00469           if (temp_phi < 0.)
00470             {
00471               result.global_phi = 0;
00472             }
00473           else if (temp_phi >= maxPhiG)
00474             {
00475               result.global_phi = maxPhiG - 1;
00476             }
00477           else
00478             {
00479              result.global_phi = static_cast<unsigned short>(temp_phi);
00480             }
00481 
00482           LogDebug("CSCSectorReceiverLUT")
00483             << "local_phi = " << local_phi
00484             << " halfstrip = " << halfstrip << " strip = " << strip
00485             << " distFromHalfStripCenter = " << distFromHalfStripCenter
00486             << " halfstripWidth = " << halfstripWidth
00487             << " strip phi = " << strip_phi/(M_PI/180.)
00488             << " temp_phi = " << temp_phi*CSCTFConstants::SECTOR_DEG/maxPhiG
00489             << " global_phi = "    << result.global_phi
00490             << " " << result.global_phi*CSCTFConstants::SECTOR_DEG/maxPhiG;
00491 
00492         }
00493     }
00494   catch(edm::Exception& e)
00495     {
00496       edm::LogError("CSCSectorReceiverLUT|getGlobalPhiValue") << e.what();
00497     }
00498 
00499   return result;
00500 }
00501 
00502 gblphidat CSCSectorReceiverLUT::globalPhiME(int phi_local, int wire_group, int cscid) const
00503 {
00504   gblphidat result;
00505   gblphiadd theadd;
00506   theadd.phi_local = phi_local;
00507   theadd.wire_group = ((1<<5)-1)&(wire_group >> 2); // want 2-7 of wg
00508   theadd.cscid = cscid;
00509 
00510   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMEMini(_endcap, _sector, _station, _subsector, theadd.toint());
00511   else if(LUTsFromFile) result = me_global_phi[theadd.toint()];
00512   else result = calcGlobalPhiME(theadd);
00513 
00514   return result;
00515 }
00516 
00517 gblphidat CSCSectorReceiverLUT::globalPhiME(unsigned address) const
00518 {
00519   gblphidat result;
00520 
00521   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMEMini(_endcap, _sector, _station, _subsector, address);
00522   else if(LUTsFromFile) result = me_global_phi[address];
00523   else result = calcGlobalPhiME(gblphiadd(address));
00524 
00525   return result;
00526 }
00527 
00528 gblphidat CSCSectorReceiverLUT::globalPhiME(gblphiadd address) const
00529 {
00530   gblphidat result;
00531 
00532   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMEMini(_endcap, _sector, _station, _subsector, address.toint());
00533   else if(LUTsFromFile) result = me_global_phi[address.toint()];
00534   else result = calcGlobalPhiME(address);
00535 
00536   return result;
00537 }
00538 
00539 gblphidat CSCSectorReceiverLUT::calcGlobalPhiMB(const gblphidat &csclut) const
00540 {
00541   gblphidat dtlut;
00542 
00543   // The following method was ripped from D. Holmes' LUT conversion program
00544   // modifications from Darin and GP
00545   int GlobalPhiMin = (_subsector == 1) ? 0x42 : 0x800;  // (0.999023 : 31 in degrees)
00546   int GlobalPhiMax = (_subsector == 1) ? 0x7ff : 0xfbd; // (30.985 : 60.986 in degrees)
00547   double GlobalPhiShift = (1.0*GlobalPhiMin + (GlobalPhiMax - GlobalPhiMin)/2.0);
00548 
00549   double dt_out = static_cast<double>(csclut.global_phi) - GlobalPhiShift;
00550 
00551   // these numbers are 62 deg / 1 rad (CSC phi scale vs. DT phi scale)
00552   dt_out = (dt_out/1982)*2145; //CSC phi 62 degrees; DT phi 57.3 degrees
00553 
00554   if(dt_out >= 0) // msb != 1
00555     {
00556       dtlut.global_phi = 0x7ff&static_cast<unsigned>(dt_out);
00557     }
00558   else
00559     {
00560       dtlut.global_phi = static_cast<unsigned>(-dt_out);
00561       dtlut.global_phi = ~dtlut.global_phi;
00562       dtlut.global_phi |= 0x800;
00563     }
00564 
00565   return dtlut;
00566 }
00567 
00568 gblphidat CSCSectorReceiverLUT::globalPhiMB(int phi_local,int wire_group, int cscid) const
00569 {
00570   gblphiadd address;
00571   gblphidat result;
00572 
00573   address.cscid = cscid;
00574   address.wire_group = ((1<<5)-1)&(wire_group>>2);
00575   address.phi_local = phi_local;
00576 
00577   // comment for now
00578   //  if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMBMini(_endcap, _sector, _subsector, address.toint());
00579   //else 
00580   if(LUTsFromFile) result = mb_global_phi[address.toint()];
00581   else result = calcGlobalPhiMB(globalPhiME(address));
00582 
00583   return result;
00584 }
00585 
00586 gblphidat CSCSectorReceiverLUT::globalPhiMB(unsigned address) const
00587 {
00588   gblphidat result;
00589   gblphiadd theadd(address);
00590 
00591   //if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMBMini(_endcap, _sector, _subsector, address);
00592   //else 
00593   if(LUTsFromFile) result = mb_global_phi[theadd.toint()];
00594   else result = calcGlobalPhiMB(globalPhiME(address));
00595 
00596   return result;
00597 }
00598 
00599 gblphidat CSCSectorReceiverLUT::globalPhiMB(gblphiadd address) const
00600 {
00601   gblphidat result;
00602 
00603   //if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalPhiMBMini(_endcap, _sector, _subsector, address.toint());
00604   //else 
00605   if(LUTsFromFile) result = mb_global_phi[address.toint()];
00606   else result = calcGlobalPhiMB(globalPhiME(address));
00607 
00608   return result;
00609 }
00610 
00611 double CSCSectorReceiverLUT::getGlobalEtaValue(const unsigned& thecscid, const unsigned& thewire_group, const unsigned& thephi_local) const
00612 {
00613   double result = 0.0;
00614   unsigned wire_group = thewire_group;
00615   int cscid = thecscid;
00616   unsigned phi_local = thephi_local;
00617 
00618   // Flag to be set if one wants to apply phi corrections ONLY in ME1/1.
00619   // Turn it into a parameter?
00620   bool me1ir_only = false;
00621 
00622   if(cscid < CSCTriggerNumbering::minTriggerCscId() ||
00623      cscid > CSCTriggerNumbering::maxTriggerCscId()) {
00624        edm::LogWarning("CSCSectorReceiverLUT|getEtaValue")
00625          << " warning: cscId " << cscid
00626          << " is out of bounds [1-" << CSCTriggerNumbering::maxTriggerCscId()
00627          << "]\n";
00628       cscid = CSCTriggerNumbering::maxTriggerCscId();
00629     }
00630 
00631   CSCTriggerGeomManager* thegeom = CSCTriggerGeometry::get();
00632   CSCLayerGeometry* layerGeom = NULL;
00633   const unsigned numBins = 1 << 2; // 4 local phi bins
00634 
00635   if(phi_local > numBins - 1) {
00636       edm::LogWarning("CSCSectorReceiverLUT|getEtaValue")
00637         << "warning: phiL " << phi_local
00638         << " is out of bounds [0-" << numBins - 1 << "]\n";
00639       phi_local = numBins - 1;
00640   }
00641   try
00642     {
00643       const CSCChamber* thechamber = thegeom->chamber(_endcap,_station,_sector,_subsector,cscid);
00644       if(thechamber) {
00645         layerGeom = const_cast<CSCLayerGeometry*>(thechamber->layer(CSCConstants::KEY_ALCT_LAYER)->geometry());
00646         const unsigned nWireGroups = layerGeom->numberOfWireGroups();
00647 
00648         // Check wire group numbers; expect them to be counted from 0, as in
00649         // CorrelatedLCTDigi class.
00650         if (wire_group >= nWireGroups) {
00651            edm::LogWarning("CSCSectorReceiverLUT|getEtaValue")
00652              << "warning: wireGroup " << wire_group
00653             << " is out of bounds [0-" << nWireGroups << ")\n";
00654           wire_group = nWireGroups - 1;
00655         }
00656         // Convert to [1; nWireGroups] range used in geometry methods.
00657         wire_group += 1;
00658 
00659         // If me1ir_only is set, apply phi corrections only in ME1/1.
00660         if (me1ir_only &&
00661             (_station != 1 ||
00662              CSCTriggerNumbering::ringFromTriggerLabels(_station, cscid) != 1))
00663           {
00664             result = thechamber->layer(CSCConstants::KEY_ALCT_LAYER)->centerOfWireGroup(wire_group).eta();
00665           }
00666         else {
00667           const unsigned nStrips = layerGeom->numberOfStrips();
00668           const unsigned nStripsPerBin = CSCConstants::MAX_NUM_STRIPS/numBins;
00673           // Check that no strips will be left out.
00674           if (nStrips%numBins != 0 || CSCConstants::MAX_NUM_STRIPS%numBins != 0)
00675             edm::LogWarning("CSCSectorReceiverLUT")
00676               << "getGlobalEtaValue warning: number of strips "
00677               << nStrips << " (" << CSCConstants::MAX_NUM_STRIPS
00678               << ") is not divisible by numBins " << numBins
00679               << " Station " << _station << " sector " << _sector
00680               << " subsector " << _subsector << " cscid " << cscid << "\n";
00681 
00682           unsigned    maxStripPrevBin = 0, maxStripThisBin = 0;
00683           unsigned    correctionStrip;
00684           LocalPoint  lPoint;
00685           GlobalPoint gPoint;
00686           // Bins phi_local and find the the middle strip for each bin.
00687           maxStripThisBin = nStripsPerBin * (phi_local+1);
00688           if (maxStripThisBin <= nStrips) {
00689             correctionStrip = nStripsPerBin/2 * (2*phi_local+1);
00690           }
00691           else {
00692             // If the actual number of strips in the chamber is smaller than
00693             // the number of strips corresponding to the right edge of this phi
00694             // local bin, we take the middle strip between number of strips
00695             // at the left edge of the bin and the actual number of strips.
00696             maxStripPrevBin = nStripsPerBin * phi_local;
00697             correctionStrip = (nStrips+maxStripPrevBin)/2;
00698           }
00699 
00700           lPoint = layerGeom->stripWireGroupIntersection(correctionStrip, wire_group);
00701           gPoint = thechamber->layer(CSCConstants::KEY_ALCT_LAYER)->surface().toGlobal(lPoint);
00702 
00703           // end calc of eta correction.
00704           result = gPoint.eta();
00705         }
00706       }
00707     }
00708   catch (cms::Exception &e)
00709     {
00710       LogDebug("CSCSectorReceiver|OutofBoundInput") << e.what();
00711     }
00712 
00713   return std::fabs(result);
00714 }
00715 
00716 
00717 gbletadat CSCSectorReceiverLUT::calcGlobalEtaME(const gbletaadd& address) const
00718 {
00719   gbletadat result;
00720   double float_eta = getGlobalEtaValue(address.cscid, address.wire_group, address.phi_local);
00721   unsigned int_eta = 0;
00722   unsigned bend_global = 0; // not filled yet... will change when it is.
00723   const double etaPerBin = (CSCTFConstants::maxEta - CSCTFConstants::minEta)/CSCTFConstants::etaBins;
00724   const unsigned me12EtaCut = 56;
00725 
00726   if ((float_eta < CSCTFConstants::minEta) || (float_eta >= CSCTFConstants::maxEta))
00727     {
00728       edm::LogWarning("CSCSectorReceiverLUT:OutOfBounds")
00729         << "CSCSectorReceiverLUT warning: float_eta = " << float_eta
00730         << " minEta = " << CSCTFConstants::minEta << " maxEta = " << CSCTFConstants::maxEta
00731         << "   station " << _station << " sector " << _sector
00732         << " chamber "   << address.cscid << " wire group " << address.wire_group;
00733 
00734       throw cms::Exception("CSCSectorReceiverLUT")
00735         << "+++ Value of CSC ID, " << float_eta
00736         << ", is out of bounds [" << CSCTFConstants::minEta << "-"
00737         << CSCTFConstants::maxEta << ") +++\n";
00738 
00739       //if (float_eta < CSCTFConstants::minEta)
00740       //result.global_eta = 0;
00741       //else if (float_eta >= CSCTFConstants::maxEta)
00742       //result.global_eta = CSCTFConstants::etaBins - 1;
00743     }
00744   else
00745     {
00746       float_eta -= CSCTFConstants::minEta;
00747       float_eta = float_eta/etaPerBin;
00748       int_eta = static_cast<unsigned>(float_eta);
00749       /* Commented until I find out its use.
00750       // Fine-tune eta boundary between DT and CSC.
00751       if ((intEta == L1MuCSCSetup::CscEtaStart() && (L1MuCSCSetup::CscEtaStartCorr() > 0.) ) ||
00752           (intEta == L1MuCSCSetup::CscEtaStart() - 1 && (L1MuCSCSetup::CscEtaStartCorr() < 0.) ) ) {
00753         bitEta = (thisEta-minEta-L1MuCSCSetup::CscEtaStartCorr())/EtaPerBin;
00754         intEta = static_cast<int>(bitEta);
00755       }
00756       */
00757       if (_station == 1 && address.cscid >= static_cast<unsigned>(CSCTriggerNumbering::minTriggerCscId())
00758           && address.cscid <= static_cast<unsigned>(CSCTriggerNumbering::maxTriggerCscId()) )
00759         {
00760           unsigned ring = CSCTriggerNumbering::ringFromTriggerLabels(_station, address.cscid);
00761 
00762           if      (ring == 1 && int_eta <  me12EtaCut) {int_eta = me12EtaCut;}
00763           else if (ring == 2 && int_eta >= me12EtaCut) {int_eta = me12EtaCut-1;}
00764         }
00765       result.global_eta = int_eta;
00766     }
00767   result.global_bend = bend_global;
00768 
00769   return result;
00770 }
00771 
00772 gbletadat CSCSectorReceiverLUT::globalEtaME(int tphi_bend, int tphi_local, int twire_group, int tcscid) const
00773 {
00774   gbletadat result;
00775   gbletaadd theadd;
00776 
00777   theadd.phi_bend = tphi_bend;
00778   theadd.phi_local = (tphi_local>>(CSCBitWidths::kLocalPhiDataBitWidth - 2)) & 0x3; // want 2 msb of local phi
00779   theadd.wire_group = twire_group;
00780   theadd.cscid = tcscid;
00781 
00782   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalEtaMEMini(_endcap, _sector, _station, _subsector, theadd.toint());
00783   else if(LUTsFromFile) result = me_global_eta[theadd.toint()];
00784   else result = calcGlobalEtaME(theadd);
00785 
00786   return result;
00787 }
00788 
00789 gbletadat CSCSectorReceiverLUT::globalEtaME(unsigned address) const
00790 {
00791   gbletadat result;
00792   gbletaadd theadd(address);
00793 
00794   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalEtaMEMini(_endcap, _sector, _station, _subsector, address);
00795   else if(LUTsFromFile) result = me_global_eta[address];
00796   else result = calcGlobalEtaME(theadd);
00797   return result;
00798 }
00799 
00800 gbletadat CSCSectorReceiverLUT::globalEtaME(gbletaadd address) const
00801 {
00802   gbletadat result;
00803 
00804   if(useMiniLUTs && isTMB07) result = CSCSectorReceiverMiniLUT::calcGlobalEtaMEMini(_endcap, _sector, _station, _subsector, address.toint());
00805   else if(LUTsFromFile) result = me_global_eta[address.toint()];
00806   else result = calcGlobalEtaME(address);
00807   return result;
00808 }
00809 
00810 std::string CSCSectorReceiverLUT::encodeFileIndex() const {
00811   std::string fileName = "";
00812   if (_station == 1) {
00813     if (_subsector == 1) fileName += "1a";
00814     if (_subsector == 2) fileName += "1b";
00815   }
00816   else if (_station == 2) fileName += "2";
00817   else if (_station == 3) fileName += "3";
00818   else if (_station == 4) fileName += "4";
00819   fileName += "End";
00820   if (_endcap == 1) fileName += "1";
00821   else                fileName += "2";
00822   fileName += "Sec";
00823   if      (_sector == 1) fileName += "1";
00824   else if (_sector == 2) fileName += "2";
00825   else if (_sector == 3) fileName += "3";
00826   else if (_sector == 4) fileName += "4";
00827   else if (_sector == 5) fileName += "5";
00828   else if (_sector == 6) fileName += "6";
00829   fileName += "LUT";
00830   return fileName;
00831 }
00832 
00833 void CSCSectorReceiverLUT::readLUTsFromFile()
00834 {
00835   if(!me_lcl_phi_loaded)
00836     {
00837       me_lcl_phi = new lclphidat[1<<CSCBitWidths::kLocalPhiAddressWidth];
00838       memset(me_lcl_phi, 0, (1<<CSCBitWidths::kLocalPhiAddressWidth)*sizeof(short));
00839       std::string fName(me_lcl_phi_file.fullPath());
00840       std::ifstream LocalPhiLUT;
00841 
00842       edm::LogInfo("CSCSectorReceiverLUT") << "Loading SR LUT: " << fName;
00843 
00844       if(isBinary)
00845         {
00846           LocalPhiLUT.open(fName.c_str(),std::ios::binary);
00847           LocalPhiLUT.seekg(0,std::ios::end);
00848           int length = LocalPhiLUT.tellg();
00849           if(length == (1<<CSCBitWidths::kLocalPhiAddressWidth)*sizeof(short))
00850             {
00851               LocalPhiLUT.seekg(0,std::ios::beg);
00852               LocalPhiLUT.read(reinterpret_cast<char*>(me_lcl_phi),length);
00853               LocalPhiLUT.close();
00854             }
00855           else
00856             edm::LogError("CSCSectorReceiverLUT") << "File "<< fName << " is incorrect size!";
00857           LocalPhiLUT.close();
00858         }
00859       else
00860         {
00861           LocalPhiLUT.open(fName.c_str());
00862           unsigned i = 0;
00863           unsigned short temp = 0;
00864           while(!LocalPhiLUT.eof() && i < 1<<CSCBitWidths::kLocalPhiAddressWidth)
00865             {
00866               LocalPhiLUT >> temp;
00867               me_lcl_phi[i++] = (*reinterpret_cast<lclphidat*>(&temp));
00868             }
00869           LocalPhiLUT.close();
00870         }
00871     }
00872   if(!me_global_phi)
00873     {
00874       me_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00875       memset(me_global_phi, 0, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(short));
00876       std::string fName = me_gbl_phi_file.fullPath();
00877       std::ifstream GlobalPhiLUT;
00878 
00879       edm::LogInfo("CSCSectorReceiverLUT") << "Loading SR LUT: " << fName;
00880 
00881       if(isBinary)
00882         {
00883           GlobalPhiLUT.open(fName.c_str(),std::ios::binary);
00884           GlobalPhiLUT.seekg(0,std::ios::end);
00885           int length = GlobalPhiLUT.tellg();
00886           if(length == (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(short))
00887             {
00888               GlobalPhiLUT.seekg(0,std::ios::beg);
00889               GlobalPhiLUT.read(reinterpret_cast<char*>(me_global_phi),length);
00890             }
00891           else
00892             edm::LogError("CSCSectorReceiverLUT") << "File "<< fName << " is incorrect size!";
00893           GlobalPhiLUT.close();
00894         }
00895       else
00896         {
00897           GlobalPhiLUT.open( fName.c_str());
00898           unsigned short temp = 0;
00899           unsigned i = 0;
00900           while(!GlobalPhiLUT.eof() && i < 1<<CSCBitWidths::kGlobalPhiAddressWidth)
00901             {
00902               GlobalPhiLUT >> temp;
00903               me_global_phi[i++] = (*reinterpret_cast<gblphidat*>(&temp));
00904             }
00905           GlobalPhiLUT.close();
00906         }
00907     }
00908   if(!mb_global_phi && _station == 1) // MB lut only in station one.
00909     {
00910       mb_global_phi = new gblphidat[1<<CSCBitWidths::kGlobalPhiAddressWidth];
00911       memset(mb_global_phi, 0, (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(short));
00912       std::string fName = mb_gbl_phi_file.fullPath();
00913       std::ifstream GlobalPhiLUT;
00914 
00915       edm::LogInfo("CSCSectorReceiverLUT") << "Loading SR LUT: " << fName;
00916 
00917       if(isBinary)
00918         {
00919           GlobalPhiLUT.open( fName.c_str(),std::ios::binary);
00920           GlobalPhiLUT.seekg(0,std::ios::end);
00921           int length = GlobalPhiLUT.tellg();
00922           if(length == (1<<CSCBitWidths::kGlobalPhiAddressWidth)*sizeof(short))
00923             {
00924               GlobalPhiLUT.seekg(0,std::ios::beg);
00925               GlobalPhiLUT.read(reinterpret_cast<char*>(mb_global_phi),length);
00926             }
00927           else
00928             edm::LogError("CSCSectorReceiverLUT") << "File "<< fName << " is incorrect size!";
00929           GlobalPhiLUT.close();
00930         }
00931       else
00932         {
00933           GlobalPhiLUT.open(fName.c_str());
00934           unsigned short temp = 0;
00935           unsigned i = 0;
00936           while(!GlobalPhiLUT.eof() && i < 1<<CSCBitWidths::kGlobalPhiAddressWidth)
00937             {
00938               GlobalPhiLUT >> temp;
00939               mb_global_phi[i++] = (*reinterpret_cast<gblphidat*>(&temp));
00940             }
00941           GlobalPhiLUT.close();
00942         }
00943     }
00944   if(!me_global_eta)
00945     {
00946       me_global_eta = new gbletadat[1<<CSCBitWidths::kGlobalEtaAddressWidth];
00947       memset(me_global_eta, 0, (1<<CSCBitWidths::kGlobalEtaAddressWidth)*sizeof(short));
00948       std::string fName = me_gbl_eta_file.fullPath();
00949       std::ifstream GlobalEtaLUT;
00950 
00951       edm::LogInfo("CSCSectorReceiverLUT") << "Loading SR LUT: " << fName;
00952 
00953       if(isBinary)
00954         {
00955           GlobalEtaLUT.open(fName.c_str(),std::ios::binary);
00956           GlobalEtaLUT.seekg(0,std::ios::end);
00957           int length = GlobalEtaLUT.tellg();
00958           if(length == (1<<CSCBitWidths::kGlobalEtaAddressWidth)*sizeof(short))
00959             {
00960               GlobalEtaLUT.seekg(0,std::ios::beg);
00961               GlobalEtaLUT.read(reinterpret_cast<char*>(me_global_eta),length);
00962             }
00963           else
00964             edm::LogError("CSCSectorReceiverLUT") << "File "<< fName << " is incorrect size!";
00965           GlobalEtaLUT.close();
00966         }
00967       else
00968         {
00969           GlobalEtaLUT.open(fName.c_str());
00970           unsigned short temp = 0;
00971           unsigned i = 0;
00972           while(!GlobalEtaLUT.eof() && i < 1<<CSCBitWidths::kGlobalEtaAddressWidth)
00973           {
00974             GlobalEtaLUT >> temp;
00975             me_global_eta[i++] = (*reinterpret_cast<gbletadat*>(&temp));
00976           }
00977           GlobalEtaLUT.close();
00978         }
00979     }
00980 }
00981