Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:26

0001 // -*- C++ -*-
0002 //
0003 // Package:    WriteVHDL
0004 // Class:      WriteVHDL
0005 //
0006 /**\class WriteVHDL WriteVHDL.cc L1TriggerConfig/WriteVHDL/src/WriteVHDL.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Tomasz Maciej Frueboes
0015 //         Created:  Tue Mar 18 15:15:30 CET 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 
0031 #include "FWCore/Framework/interface/ESHandle.h"
0032 
0033 #include "CondFormats/DataRecord/interface/L1RPCConfigRcd.h"
0034 #include "CondFormats/L1TObjects/interface/L1RPCConfig.h"
0035 #include "CondFormats/RPCObjects/interface/L1RPCHwConfig.h"
0036 
0037 #include "CondFormats/DataRecord/interface/L1RPCConeBuilderRcd.h"
0038 #include "CondFormats/RPCObjects/interface/L1RPCConeBuilder.h"
0039 
0040 #include "CondFormats/L1TObjects/interface/L1RPCConeDefinition.h"
0041 
0042 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0043 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0044 
0045 #include "CondFormats/RPCObjects/interface/RPCReadOutMapping.h"
0046 
0047 #include "CondFormats/RPCObjects/interface/RPCEMap.h"
0048 #include "CondFormats/DataRecord/interface/RPCEMapRcd.h"
0049 
0050 #include <fstream>
0051 #include <bitset>
0052 
0053 //
0054 // class decleration
0055 //
0056 
0057 //need SharedResource since internally uses a static function variable
0058 class WriteVHDL : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0059 public:
0060   explicit WriteVHDL(const edm::ParameterSet&);
0061 
0062 private:
0063   virtual void analyze(const edm::Event&, const edm::EventSetup&);
0064   int getDCCNumber(int iTower, int iSec);
0065   int getDCC(int iSec);
0066   int getTBNumber(int iTower);
0067   int getDCCNumberFromTB(int tbNumber, int iSec);
0068   int m_towerBeg;
0069   int m_towerEnd;
0070   int m_sectorBeg;
0071   int m_sectorEnd;
0072   std::string m_templateName;
0073   std::string m_outdirName;
0074   edm::ESGetToken<L1RPCConfig, L1RPCConfigRcd> m_confToken;
0075   edm::ESGetToken<L1RPCConeBuilder, L1RPCConeBuilderRcd> m_coneBuilderToken;
0076   edm::ESGetToken<RPCGeometry, MuonGeometryRecord> m_rpcGeomToken;
0077   edm::ESGetToken<L1RPCConeDefinition, L1RPCConeDefinitionRcd> m_coneDefToken;
0078   edm::ESGetToken<RPCEMap, RPCEMapRcd> m_nmapToken;
0079 
0080   struct TBLoc {
0081     TBLoc(int tb, int sec) : tbNum(tb), sector(sec){};
0082     int tbNum;
0083     int sector;
0084     bool operator<(const TBLoc& c2) const {
0085       if (this->tbNum != c2.tbNum)
0086         return (this->tbNum < c2.tbNum);
0087       else
0088         return (this->sector < c2.sector);
0089     }
0090   };
0091 
0092   typedef std::map<TBLoc, std::set<int> > TBInputsMap;  /// value contains list of used TB inputs, key localizes TB
0093   TBInputsMap m_tbInputs;
0094 
0095   typedef std::map<TBLoc, std::map<int, int> > TBInputs3to4Map;
0096   TBInputs3to4Map m_tbInputs3to4;
0097 
0098   struct TDetStrip {
0099     TDetStrip(int d, int s) : detId(d), strip(s){};
0100     int detId;
0101     int strip;
0102     bool operator<(const TDetStrip& c2) const {
0103       if (this->detId != c2.detId)
0104         return (this->detId < c2.detId);
0105       else
0106         return (this->strip < c2.strip);
0107     }
0108   };
0109 
0110   struct TStripConnection {
0111     TStripConnection(int t, int l, int p) : tbInput(t), lbInTBInput(l), packedStrip(p){};
0112     TStripConnection() : tbInput(-1), lbInTBInput(-1), packedStrip(-1){};
0113     int tbInput;
0114     int lbInTBInput;
0115     int packedStrip;
0116   };
0117 
0118   typedef std::map<TDetStrip, TStripConnection> TStrip2Con;
0119   typedef std::map<TBLoc, TStrip2Con> TTB2Con;
0120   TTB2Con m_4thPlaneConnections;
0121 
0122   void writePats(const edm::EventSetup& evtSetup, int tower, int logsector);
0123 
0124   std::string writeVersion();
0125 
0126   std::string writeCNT(const edm::EventSetup& iSetup, int tower, int logsector, std::string pacT);
0127 
0128   std::string writePACandLPDef(const edm::EventSetup& iSetup, int tower, int logsector, std::string PACt);
0129 
0130   std::string writeConeDef(const edm::EventSetup& iSetup, int tower, int sector, std::string PACt);
0131 
0132   std::string writeQualTable(const edm::EventSetup& iSetup, int tower, int sector);
0133 
0134   std::string writePatterns(const edm::EventSetup& iSetup, int tower, int sector, std::string PACt);
0135 
0136   std::string writeGB(std::string PACt);
0137 
0138   void prepareEncdap4thPlaneConnections(edm::ESHandle<RPCGeometry> geom, edm::ESHandle<RPCReadOutMapping> map);
0139   // ----------member data ---------------------------
0140 };
0141 
0142 WriteVHDL::WriteVHDL(const edm::ParameterSet& iConfig)
0143 
0144 {
0145   //now do what ever initialization is needed
0146   m_towerBeg = iConfig.getParameter<int>("minTower");
0147   m_towerEnd = iConfig.getParameter<int>("maxTower");
0148 
0149   m_sectorBeg = iConfig.getParameter<int>("minSector");
0150   m_sectorEnd = iConfig.getParameter<int>("maxSector");
0151   m_templateName = iConfig.getParameter<std::string>("templateName");
0152   m_outdirName = iConfig.getParameter<std::string>("outDir");
0153 
0154   m_confToken = esConsumes();
0155   m_coneBuilderToken = esConsumes();
0156   m_rpcGeomToken = esConsumes();
0157   m_coneDefToken = esConsumes();
0158   m_nmapToken = esConsumes();
0159 }
0160 
0161 //
0162 // member functions
0163 //
0164 
0165 // ------------ method called to for each event  ------------
0166 
0167 /*
0168 XXV -- version comment
0169 XXP -- pac/logplanes def
0170 XXC -- cones
0171 XXQ -- quality table
0172 XXS -- Patterns
0173 XXG -- ghostbuster
0174 */
0175 void WriteVHDL::analyze(const edm::Event& iEvent, const edm::EventSetup& evtSetup) {
0176   for (int tw = m_towerBeg; tw <= m_towerEnd; ++tw) {
0177     for (int sec = m_sectorBeg; sec <= m_sectorEnd; ++sec) {
0178       writePats(evtSetup, tw, sec);
0179     }
0180   }
0181 }
0182 
0183 void WriteVHDL::writePats(const edm::EventSetup& evtSetup, int tower, int logsector) {
0184   std::ifstream inputfile(m_templateName.c_str());
0185   std::stringstream fname;
0186   fname << m_outdirName << "/pac_t" << tower << "_sec" << logsector << ".vhd";
0187 
0188   std::ofstream fout(fname.str().c_str());
0189 
0190   // get PAC type
0191   const L1RPCConfig* rpcconf = &evtSetup.getData(m_confToken);
0192 
0193   RPCPattern::RPCPatVec::const_iterator it = rpcconf->m_pats.begin();
0194 
0195   while (it->getTower() != std::abs(tower) && it != rpcconf->m_pats.end())
0196     ++it;
0197 
0198   if (it == rpcconf->m_pats.end())
0199     throw cms::Exception("") << " tower not found " << tower << "\n";
0200 
0201   RPCPattern::TPatternType patType = it->getPatternType();
0202 
0203   std::string pacT = "E";
0204   if (patType == RPCPattern::PAT_TYPE_T)
0205     pacT = "T";
0206 
0207   if (inputfile.fail()) {
0208     throw cms::Exception("IO") << "Cannot open file: " << m_templateName << "\n";
0209   }
0210 
0211   char ch, chNext, chCmd;
0212   while (!inputfile.eof()) {
0213     inputfile.get(ch);
0214 
0215     if (!inputfile.eof()) {
0216       if (ch == 'X') {
0217         inputfile.get(chNext);
0218         if (!inputfile.eof()) {
0219           if (chNext == 'X') {  //input command
0220             inputfile.get(chCmd);
0221             if (!inputfile.eof()) {
0222               switch (chCmd) {
0223                 case 'V':
0224                   fout << writeVersion();
0225                   break;
0226                 case 'N':
0227                   fout << writeCNT(evtSetup, tower, logsector, pacT);
0228                   break;
0229                 case 'P':
0230                   fout << writePACandLPDef(evtSetup, tower, logsector, pacT);
0231                   break;
0232                 case 'C':
0233                   fout << writeConeDef(evtSetup, tower, logsector, pacT);
0234                   break;
0235                 case 'Q':
0236                   fout << writeQualTable(evtSetup, tower, logsector);
0237                   break;
0238                 case 'S':
0239                   fout << writePatterns(evtSetup, tower, logsector, pacT);
0240                   break;
0241                 case 'G':
0242                   fout << writeGB(pacT);
0243                   break;
0244                 default:
0245                   throw cms::Exception("BadTemplate") << " Unknown command: XX" << chCmd << "\n";
0246               }
0247 
0248             } else {
0249               throw cms::Exception("BadTemplate") << " Problem when reading template \n";
0250             }
0251 
0252           } else {
0253             fout << ch << chNext;
0254           }
0255         }
0256 
0257       } else {
0258         fout << ch;
0259       }
0260     }
0261   }
0262 }
0263 
0264 std::string WriteVHDL::writeVersion() {
0265   std::stringstream ret;
0266   //Get current time
0267   time_t rawtime;
0268   struct tm* timeinfo;
0269   time(&rawtime);
0270   timeinfo = localtime(&rawtime);
0271   //asctime (timeinfo);
0272 
0273   ret << "-- WriteVHDL " << asctime(timeinfo) << std::endl;
0274   return ret.str();
0275 }
0276 
0277 std::string WriteVHDL::writeCNT(const edm::EventSetup& iSetup, int tower, int sector, std::string pacT) {
0278   std::stringstream ret;
0279   int nT = 0, nE = 0, refGrps = 0;
0280   if (pacT == "E") {
0281     nE = 12;
0282   } else if (pacT == "T") {
0283     nT = 12;
0284   } else
0285     throw cms::Exception("") << "Unknown PAC type \n";
0286 
0287   ret << "constant TT_EPACS_COUNT         :natural := " << nE << ";" << std::endl;
0288   ret << "constant TT_TPACS_COUNT         :natural := " << nT << ";" << std::endl;
0289 
0290   // calculate number of RefGruops used
0291   if (pacT == "E") {
0292     tower = std::abs(tower);
0293 
0294     edm::ESHandle<L1RPCConfig> conf = iSetup.getHandle(m_confToken);
0295 
0296     const RPCPattern::RPCPatVec* pats = &conf.product()->m_pats;
0297     int ppt = conf.product()->getPPT();
0298     int segment = 0;
0299 
0300     if (ppt == 1 || ppt == 12) {
0301       sector = 0;
0302     }
0303 
0304     const RPCPattern::RPCPatVec::const_iterator itEnd = pats->end();
0305     RPCPattern::RPCPatVec::const_iterator it;
0306 
0307     for (int iPAC = 0; iPAC < 12; ++iPAC) {
0308       if (ppt == 144 || ppt == 12)
0309         segment = iPAC;
0310 
0311       for (it = pats->begin(); it != itEnd; ++it) {
0312         // select right pac
0313         if (it->getTower() != tower || it->getLogSector() != sector || it->getLogSegment() != segment)
0314           continue;
0315 
0316         if (it->getPatternType() != RPCPattern::PAT_TYPE_E) {
0317           throw cms::Exception("WriteVHDL") << "Expected E type pattern, got different one" << std::endl;
0318         }
0319         if (refGrps < it->getRefGroup())
0320           refGrps = it->getRefGroup();
0321 
0322       }  // patsIter
0323     }    // segment iter
0324   }      // if type E
0325 
0326   ret << "constant TT_REF_GROUP_NUMBERS   :natural := " << refGrps + 1 << ";" << std::endl;
0327 
0328   return ret.str();
0329 }
0330 
0331 std::string WriteVHDL::writePACandLPDef(const edm::EventSetup& iSetup, int tower, int logsector, std::string pacT) {
0332   std::stringstream ret;
0333 
0334   tower = std::abs(tower);
0335 
0336   // get logplane size
0337   edm::ESHandle<L1RPCConeBuilder> coneBuilder = iSetup.getHandle(m_coneBuilderToken);
0338 
0339   edm::ESHandle<L1RPCConeDefinition> l1RPCConeDefinition = iSetup.getHandle(m_coneDefToken);
0340 
0341   std::string coma = "";
0342   for (int seg = 0; seg < 12; ++seg) {
0343     ret << coma << "   (" << seg << ",  " << pacT << ", (  ";
0344 
0345     std::string coma1 = "";
0346     for (int lp = 0; lp < 6; ++lp) {
0347       //int size = l1RPCConeDefinition->getLPSizes().at(tower).at(lp);
0348       int lpSize = -1;
0349       L1RPCConeDefinition::TLPSizeVec::const_iterator it = l1RPCConeDefinition->getLPSizeVec().begin();
0350       L1RPCConeDefinition::TLPSizeVec::const_iterator itEnd = l1RPCConeDefinition->getLPSizeVec().end();
0351       for (; it != itEnd; ++it) {
0352         if (it->m_tower != std::abs(tower) || it->m_LP != lp)
0353           continue;
0354         lpSize = it->m_size;
0355       }
0356 
0357       //FIXME
0358       if (lpSize == -1) {
0359         throw cms::Exception("getLogStrip") << " lpSize==-1\n";
0360       }
0361 
0362       if (lpSize == 0 || lpSize == -1)
0363         lpSize = 1;
0364       ret << coma1 << lpSize;
0365       coma1 = ", ";
0366     }
0367 
0368     ret << "))--" << std::endl;
0369     coma = ",";
0370   }
0371 
0372   ret << ");" << std::endl;
0373   return ret.str();
0374 }
0375 
0376 std::string WriteVHDL::writeQualTable(const edm::EventSetup& iSetup, int tower, int sector) {
0377   std::stringstream ret;
0378 
0379   edm::ESHandle<L1RPCConfig> conf = iSetup.getHandle(m_confToken);
0380 
0381   const RPCPattern::TQualityVec* qvec = &conf.product()->m_quals;
0382 
0383   bool first = true;
0384   RPCPattern::TQualityVec::const_iterator it = qvec->begin();
0385   RPCPattern::TQualityVec::const_iterator itEnd = qvec->end();
0386 
0387   //unsigned int ppt = conf.product()->getPPT();
0388   // if (ppt == 1) {
0389   sector = 0;
0390   // }
0391 
0392   int noOfQualitiesWritten = 0;
0393 
0394   for (; it != itEnd; ++it) {
0395     // there is only one PACCellQuality for 12 comparators!
0396     if (it->m_tower != std::abs(tower) || it->m_logsector != sector || it->m_logsegment != 0)
0397       continue;
0398 
0399     if (first) {
0400       ret << " (";
0401       first = false;
0402     } else {
0403       ret << ", " << std::endl << " (";
0404     }
0405 
0406     std::bitset<6> fp(it->m_FiredPlanes);
0407     ret << (int)it->m_QualityTabNumber
0408         << ",\""
0409         //           << (int)it->m_FiredPlanes<<"\","
0410         << fp.to_string<char, std::char_traits<char>, std::allocator<char> >() << "\"," << (int)it->m_QualityValue
0411         << ")";
0412     ++noOfQualitiesWritten;
0413   }
0414 
0415   if (noOfQualitiesWritten == 1 && std::abs(tower) == 9) {
0416   }
0417 
0418   ret << ");" << std::endl << std::endl;
0419 
0420   return ret.str();
0421 }
0422 
0423 std::string WriteVHDL::writePatterns(const edm::EventSetup& iSetup, int tower, int sector, std::string pacT) {
0424   std::stringstream ret;
0425 
0426   tower = std::abs(tower);
0427 
0428   edm::ESHandle<L1RPCConfig> conf = iSetup.getHandle(m_confToken);
0429 
0430   const RPCPattern::RPCPatVec* pats = &conf.product()->m_pats;
0431   int ppt = conf.product()->getPPT();
0432   int segment = 0;
0433 
0434   if (ppt == 1 || ppt == 12) {
0435     sector = 0;
0436   }
0437 
0438   const RPCPattern::RPCPatVec::const_iterator itEnd = pats->end();
0439   RPCPattern::RPCPatVec::const_iterator it;
0440   int to[6], globalPatNo = 0;
0441   bool firstRun = true;
0442 
0443   for (int iPAC = 0; iPAC < 12; ++iPAC) {
0444     if (ppt == 144 || ppt == 12)
0445       segment = iPAC;
0446 
0447     for (it = pats->begin(); it != itEnd; ++it) {
0448       // select right pac
0449       if (it->getTower() != tower || it->getLogSector() != sector || it->getLogSegment() != segment)
0450         continue;
0451 
0452       for (int i = 0; i < 6; ++i) {
0453         to[i] = it->getStripTo(i) - 1;
0454         if (it->getStripFrom(i) == RPCPattern::m_NOT_CONECTED)
0455           to[i] = RPCPattern::m_NOT_CONECTED;
0456       }
0457 
0458       if (!firstRun)
0459         ret << ",";
0460       firstRun = false;
0461 
0462       int sign = it->getSign();
0463 
0464       if (sign == 0)
0465         sign = 1;
0466       else if (sign == 1)
0467         sign = 0;
0468       else
0469         throw cms::Exception("BAD sign") << "Bad sign definition: " << sign << std::endl;
0470 
0471       ret << "( " << iPAC << ", " << pacT << ", " << it->getRefGroup() << ", " << it->getQualityTabNumber()
0472           << ",("                                                // planes start
0473           << "(" << it->getStripFrom(0) << "," << to[0] << ")"   // pl1
0474           << ",(" << it->getStripFrom(1) << "," << to[1] << ")"  // pl2
0475           << ",(" << it->getStripFrom(2) << "," << to[2] << ")"  // pl3
0476           << ",(" << it->getStripFrom(3) << "," << to[3] << ")"  // pl4
0477           << ",(" << it->getStripFrom(4) << "," << to[4] << ")"  // pl5
0478           << ",(" << it->getStripFrom(5) << "," << to[5] << ")"  // pl6
0479           << ")"                                                 // planes end
0480           << ", " << sign << ", " << it->getCode() << ") -- " << globalPatNo++ << std::endl;
0481 
0482     }  // patterns iteration
0483 
0484   }  // segment
0485 
0486   ret << ");" << std::endl << std::endl;
0487 
0488   return ret.str();
0489 }
0490 
0491 std::string WriteVHDL::writeGB(std::string PACt) {
0492   std::stringstream ret;
0493   bool frun = true;
0494 
0495   for (int i = 0; i < 12; ++i) {
0496     if (frun) {
0497       frun = false;
0498       ret << "(";
0499     } else {
0500       ret << std::endl << ",(";
0501     }
0502 
0503     ret << i << ", " << PACt << ", " << i << ")";
0504   }
0505 
0506   ret << ");" << std::endl;
0507 
0508   return ret.str();
0509 }
0510 
0511 void WriteVHDL::prepareEncdap4thPlaneConnections(edm::ESHandle<RPCGeometry> rpcGeom,
0512                                                  edm::ESHandle<RPCReadOutMapping> map) {
0513   static bool jobDone = true;
0514   if (jobDone)
0515     return;
0516   jobDone = true;
0517 
0518   //std::cout << "prepareEncdap4thPlaneConnections\n " ;
0519 
0520   // build map of used TB inputs
0521   for (TrackingGeometry::DetContainer::const_iterator it = rpcGeom->dets().begin(); it != rpcGeom->dets().end(); ++it) {
0522     if (dynamic_cast<const RPCRoll*>(*it) == 0)
0523       continue;
0524     RPCRoll const* roll = dynamic_cast<RPCRoll const*>(*it);
0525     int detId = roll->id().rawId();
0526 
0527     for (int strip = 1; strip <= roll->nstrips(); ++strip) {
0528       LinkBoardElectronicIndex a = {0, 0, 0, 0};
0529       std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> linkStrip =
0530           std::make_pair(a, LinkBoardPackedStrip(0, 0));
0531 
0532       std::pair<int, int> stripInDetUnit(detId, strip);
0533       std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> > aVec = map->rawDataFrame(stripInDetUnit);
0534       std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> >::const_iterator CI;
0535 
0536       for (int iSec = 0; iSec < 12; ++iSec) {
0537         int DCC = getDCC(iSec);
0538         for (int iTB = 1; iTB < 8; ++iTB) {
0539           int DCCin = getDCCNumberFromTB(iTB, iSec);
0540           int ncons = 0;
0541           for (CI = aVec.begin(); CI != aVec.end(); ++CI) {
0542             if (CI->first.dccInputChannelNum == DCCin && CI->first.dccId == DCC) {
0543               linkStrip = *CI;
0544               ++ncons;
0545             }
0546           }
0547 
0548           if (ncons > 1)
0549             std::cout << "Problem: more then one connection for given TB\n";
0550           if (ncons == 1) {
0551             TBLoc loc(iTB, iSec);
0552             int tbInput = linkStrip.first.tbLinkInputNum;
0553             //int lbInTBInput = linkStrip.first.lbNumInLink;
0554             m_tbInputs[loc].insert(tbInput);
0555             /*std::cout << " Found con for TC" << iSec << " TB" << iTB
0556                 << " - " << tbInput
0557                 << std::endl;*/
0558           }
0559         }  // tb iter
0560       }    // TC iter
0561     }      // strip iter
0562   }        // roll iter
0563 
0564   /*
0565   for(auto const& t : m_tbInputs )
0566   {
0567     std::cout 
0568         << "TB" << t.first.tbNum
0569         << " SEC" <<  t.first.sector
0570         << std::endl
0571         << "      ";
0572     for(auto const& input : t.second ){
0573       std::cout << " " << input; 
0574     }
0575     std::cout << std::endl;
0576   }
0577   */
0578 
0579   for (TrackingGeometry::DetContainer::const_iterator it = rpcGeom->dets().begin(); it != rpcGeom->dets().end(); ++it) {
0580     RPCRoll const* roll = dynamic_cast<RPCRoll const*>(*it);
0581     if (roll == 0)
0582       continue;
0583     RPCDetId d = roll->id();
0584     if (std::abs(d.region()) != 1)
0585       continue;
0586     if (d.station() != 4 && d.station() != -4)
0587       continue;
0588     int ss = 3;
0589     if (d.station() < 0)
0590       ss = -3;
0591 
0592     RPCDetId matching4stDetId(d.region(), d.ring(), ss, d.sector(), d.layer(), d.subsector(), d.roll());
0593 
0594     //std::cout << d << std::endl;
0595     int chamberMatches = 0;
0596     // check if 3d plane matches 4th plane
0597     for (TrackingGeometry::DetContainer::const_iterator it3 = rpcGeom->dets().begin(); it3 != rpcGeom->dets().end();
0598          ++it3) {
0599       RPCRoll const* roll3 = dynamic_cast<RPCRoll const*>(*it3);
0600       if (roll3 == 0)
0601         continue;
0602       RPCDetId d3 = roll3->id().rawId();
0603       if (d3 != matching4stDetId)
0604         continue;
0605       /*
0606      if (roll->id().rawId() == 637637354) {
0607        std::cout << "Match found for:  " <<  roll->id().rawId() 
0608            << " " <<d3 
0609            << "\n " << d << std::endl; 
0610    }*/
0611       ++chamberMatches;
0612       //std::cout << "    " << d3 << std::endl;
0613 
0614       if (roll->nstrips() != roll3->nstrips()) {
0615         std::cout << d << " Strips differ\n ";
0616       } else {  // check phi pos
0617         LocalPoint c = roll->centreOfStrip(roll->nstrips() / 2);
0618         LocalPoint c4 = roll3->centreOfStrip(roll3->nstrips() / 2);
0619         GlobalPoint g = roll->toGlobal(c);
0620         GlobalPoint g4 = roll3->toGlobal(c4);
0621         if (std::abs(g.phi() - g4.phi()) > 0.01)
0622           std::cout << d << " " << g.phi() << " " << g4.phi() << std::endl;
0623       }
0624 
0625       // Prepare 4th plane cons
0626       int detId3 = roll3->id().rawId();
0627 
0628       for (int strip = 1; strip <= roll3->nstrips(); ++strip) {
0629         LinkBoardElectronicIndex a = {0, 0, 0, 0};
0630         std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> linkStrip =
0631             std::make_pair(a, LinkBoardPackedStrip(0, 0));
0632 
0633         std::pair<int, int> stripInDetUnit(detId3, strip);
0634         std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> > aVec =
0635             map->rawDataFrame(stripInDetUnit);
0636         std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> >::const_iterator CI;
0637 
0638         /*
0639        if (roll->id().rawId() == 637637354 && strip == 8) {
0640          std::cout << "Enter\n";
0641         }*/
0642 
0643         for (CI = aVec.begin(); CI != aVec.end(); ++CI) {
0644           // reverse map a DCC connection to a TB
0645           int DCCin = CI->first.dccInputChannelNum;
0646           int DCC = CI->first.dccId;
0647           int tb = -1, tc = -1;
0648 
0649           int tbInput = -1;
0650           int lbInTBInput = -1;
0651           int packedStrip = -1;
0652           for (int iSec = 0; iSec < 12; ++iSec) {
0653             int matches = 0;
0654             for (int iTB = 1; iTB < 8; ++iTB) {
0655               if (DCCin == getDCCNumberFromTB(iTB, iSec)) {
0656                 if (DCC == getDCC(iSec)) {
0657                   ++matches;
0658                   tb = iTB;
0659                   tc = iSec;
0660                   tbInput = linkStrip.first.tbLinkInputNum;
0661                   lbInTBInput = linkStrip.first.lbNumInLink;
0662                   packedStrip = linkStrip.second.packedStrip();
0663                 }
0664               }
0665             }  // TB iter
0666             if (matches > 1) {
0667               std::cout << "Found more than one match!\n";
0668             }
0669 
0670             /*
0671            if (roll->id().rawId() == 637637354 && strip == 8) {
0672              std::cout << "Going through 3 cons"
0673                 << " " <<  DCC
0674                 << " " << DCCin
0675                 <<std::endl;
0676              std::cout << "DCC=" << getDCC(10)
0677                  << " DCCin=" << getDCCNumberFromTB(7, 10)
0678                  << std::endl;
0679              if (matches == 1) {
0680                std::cout << "got one\n";
0681              }
0682           }*/
0683 
0684             if (matches != 1)
0685               continue;
0686             TBLoc loc(tb, tc);
0687 
0688             // check if coresponding chamber from plane 4 hasnt have tbInput assigned, if not find new one
0689             if (m_tbInputs3to4[loc].find(tbInput) == m_tbInputs3to4[loc].end()) {
0690               bool conAdded = false;
0691               // get first free connection
0692               for (int con = 0; con < 18; ++con) {
0693                 if (m_tbInputs[loc].find(con) == m_tbInputs[loc].end()) {
0694                   m_tbInputs3to4[loc][tbInput] = con;
0695                   m_tbInputs[loc].insert(con);  // input just taken
0696                   conAdded = true;
0697                   break;
0698                 }
0699               }
0700               if (!conAdded)
0701                 std::cout << "Warning - no more connections were avaliable!!\n";
0702             }
0703 
0704             int tbInput4 = m_tbInputs3to4[loc][tbInput];
0705             m_4thPlaneConnections[loc][TDetStrip(roll->id().rawId(), strip)] =
0706                 TStripConnection(tbInput4, lbInTBInput, packedStrip);
0707             /*
0708            if (roll->id().rawId() == 637637354) {
0709              std::cout << "Connection added for:  " <<  roll->id().rawId() 
0710                  << " " << strip << std::endl; 
0711            }*/
0712 
0713           }  // TC iter
0714         }    // connection vec iter
0715       }      // strip iter
0716     }        // roll iter
0717 
0718     if (chamberMatches != 1)
0719       std::cout << d << "  -> no  of matches is  " << chamberMatches << std::endl;
0720   }
0721 }
0722 
0723 std::string WriteVHDL::writeConeDef(const edm::EventSetup& evtSetup, int tower, int sector, std::string PACt) {
0724   std::stringstream ret;
0725 
0726   edm::ESHandle<L1RPCConeBuilder> coneBuilder = evtSetup.getHandle(m_coneBuilderToken);
0727 
0728   edm::ESHandle<RPCGeometry> rpcGeom = evtSetup.getHandle(m_rpcGeomToken);
0729 
0730   edm::ESHandle<L1RPCConeDefinition> coneDef = evtSetup.getHandle(m_coneDefToken);
0731 
0732   static edm::ESHandle<RPCReadOutMapping> map;
0733   static bool isMapValid = false;
0734 
0735   if (!isMapValid) {
0736     edm::ESHandle<RPCEMap> nmap = evtSetup.getHandle(m_nmapToken);
0737     const RPCEMap* eMap = nmap.product();
0738     map = eMap->convert();  //*/
0739     isMapValid = true;
0740   }
0741 
0742   prepareEncdap4thPlaneConnections(rpcGeom, map);
0743 
0744   /*
0745    static edm::ESWatcher<RPCEMapRcd> recordWatcher;
0746    const RPCReadOutMapping* map = 0;
0747 
0748    if (recordWatcher.check(evtSetup)) {  
0749     delete map; 
0750     edm::ESHandle<RPCEMap> readoutMapping;
0751     evtSetup.get<RPCEMapRcd>().get(readoutMapping);
0752     map = readoutMapping->convert();
0753    }*/
0754 
0755   bool beg = true;
0756   for (TrackingGeometry::DetContainer::const_iterator it = rpcGeom->dets().begin(); it != rpcGeom->dets().end(); ++it) {
0757     if (dynamic_cast<const RPCRoll*>(*it) == 0)
0758       continue;
0759     RPCRoll const* roll = dynamic_cast<RPCRoll const*>(*it);
0760 
0761     int detId = roll->id().rawId();
0762     //iterate over strips
0763 
0764     for (int strip = 1; strip <= roll->nstrips(); ++strip) {
0765       std::pair<L1RPCConeBuilder::TStripConVec::const_iterator, L1RPCConeBuilder::TStripConVec::const_iterator> itPair =
0766           coneBuilder->getConVec(detId, strip);
0767 
0768       if (itPair.first != itPair.second) {
0769         throw cms::Exception("") << " FIXME found uncompressed connection. " << tower << "\n";
0770       }
0771 
0772       std::pair<L1RPCConeBuilder::TCompressedConVec::const_iterator, L1RPCConeBuilder::TCompressedConVec::const_iterator>
0773           compressedConnPair = coneBuilder->getCompConVec(detId);
0774 
0775       L1RPCConeBuilder::TCompressedConVec::const_iterator itComp = compressedConnPair.first;
0776 
0777       for (; itComp != compressedConnPair.second; ++itComp) {
0778         int logstrip = itComp->getLogStrip(strip, coneDef->getLPSizeVec());
0779         if (logstrip == -1)
0780           continue;
0781 
0782         // iterate over all PACs
0783         if (itComp->m_tower != tower)
0784           continue;
0785 
0786         int dccInputChannel = getDCCNumber(tower, sector);
0787         int PACstart = sector * 12;
0788         int PACend = PACstart + 11;
0789 
0790         for (int PAC = PACstart; PAC <= PACend; ++PAC) {
0791           if (itComp->m_PAC != PAC)
0792             continue;
0793 
0794           LinkBoardElectronicIndex a = {0, 0, 0, 0};
0795           std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> linkStrip =
0796               std::make_pair(a, LinkBoardPackedStrip(0, 0));
0797 
0798           std::pair<int, int> stripInDetUnit(detId, strip);
0799           std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> > aVec =
0800               map->rawDataFrame(stripInDetUnit);
0801           std::vector<std::pair<LinkBoardElectronicIndex, LinkBoardPackedStrip> >::const_iterator CI;
0802 
0803           //bool connectionFound = false;
0804           for (CI = aVec.begin(); CI != aVec.end(); ++CI) {
0805             if (CI->first.dccInputChannelNum == dccInputChannel)
0806               linkStrip = *CI;
0807             // connectionFound = true;
0808           }
0809           // check if it is missing plane 4 in endcap
0810 
0811           int tbLink = linkStrip.first.tbLinkInputNum;
0812           int lbNum = linkStrip.first.lbNumInLink;
0813           int packedStrip = linkStrip.second.packedStrip();
0814           bool plane4HackWorked = false;
0815           if (std::abs(roll->id().region()) == 1 && std::abs(roll->id().station()) == 4) {
0816             std::cout << "Warning, never should get here!" << std::endl;
0817             TBLoc loc(getTBNumber(tower), sector);
0818             TDetStrip ds(detId, strip);
0819             if (m_4thPlaneConnections.find(loc) == m_4thPlaneConnections.end()) {
0820               std::cout << "4thplane problem - unkown TB\n";
0821             } else if (m_4thPlaneConnections[loc].find(ds) == m_4thPlaneConnections[loc].end()) {
0822               std::cout << "4thplane problem - unkown strip " << strip << " " << detId << " " << roll->id()
0823                         << std::endl;
0824 
0825             } else {
0826               //std::cout << "4thplane fine" << std::endl;
0827               tbLink = m_4thPlaneConnections[loc][ds].tbInput;
0828               lbNum = m_4thPlaneConnections[loc][ds].lbInTBInput;
0829               packedStrip = m_4thPlaneConnections[loc][ds].packedStrip;
0830               plane4HackWorked = true;
0831             }
0832           }  // endcap 4th plane
0833 
0834           //if(!connectionFound) { test me...
0835           if (linkStrip.second.packedStrip() == -17 && !plane4HackWorked) {
0836             std::cout << "Problem" << std::endl;
0837 
0838           } else {
0839             if (!beg)
0840               ret << ",";
0841             else
0842               beg = false;
0843             ret << "(" << PAC - PACstart << ",\t " << PACt << ", \t" << (int)itComp->m_logplane - 1 << ",\t" << tbLink
0844                 << ",\t" << lbNum << ",\t" << logstrip << ",\t " << packedStrip << ", \t";
0845             ret << 1 << ") --" << roll->id() << std::endl;
0846           }
0847 
0848         }  // PAC iteration
0849 
0850       }  // cone connections interation
0851 
0852     }  // strip in roll iteration
0853 
0854   }  // roll iteration
0855 
0856   ret << "\n );";
0857 
0858   return ret.str();
0859 }
0860 
0861 // returns DCC channel for given tower, sec
0862 int WriteVHDL::getDCCNumber(int iTower, int iSec) {
0863   int tbNumber = getTBNumber(iTower);
0864   return getDCCNumberFromTB(tbNumber, iSec);  //Count DCC input channel from 1
0865 }
0866 
0867 int WriteVHDL::getDCCNumberFromTB(int tbNumber, int iSec) {
0868   int phiFactor = iSec % 4;
0869   return (tbNumber + phiFactor * 9);  //Count DCC input channel from 1
0870 }
0871 
0872 int WriteVHDL::getTBNumber(int iTower) {
0873   int tbNumber = 0;
0874   if (iTower < -12)
0875     tbNumber = 0;
0876   else if (-13 < iTower && iTower < -8)
0877     tbNumber = 1;
0878   else if (-9 < iTower && iTower < -4)
0879     tbNumber = 2;
0880   else if (-5 < iTower && iTower < -1)
0881     tbNumber = 3;
0882   else if (-2 < iTower && iTower < 2)
0883     tbNumber = 4;
0884   else if (1 < iTower && iTower < 5)
0885     tbNumber = 5;
0886   else if (4 < iTower && iTower < 9)
0887     tbNumber = 6;
0888   else if (8 < iTower && iTower < 13)
0889     tbNumber = 7;
0890   else if (12 < iTower)
0891     tbNumber = 8;
0892   return tbNumber;
0893 }
0894 
0895 int WriteVHDL::getDCC(int iSec) {
0896   if (iSec < 0 || iSec > 12)
0897     throw cms::Exception("Problem!!!\n");
0898   int DCC = 792;  // sec 0...3
0899   if (iSec > 3 && iSec < 8) {
0900     DCC = 791;
0901   } else if (iSec > 7) {
0902     DCC = 790;
0903   }
0904 
0905   return DCC;
0906 }
0907 
0908 //define this as a plug-in
0909 DEFINE_FWK_MODULE(WriteVHDL);