Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:54

0001 #include "FWCore/Framework/interface/ConsumesCollector.h"
0002 #include "L1Trigger/L1TMuonCPPF/interface/RecHitProcessor.h"
0003 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
0004 #include "DataFormats/Common/interface/DetSetVector.h"
0005 #include "L1Trigger/L1TMuonCPPF/src/CPPFClusterContainer.h"
0006 #include "L1Trigger/L1TMuonCPPF/src/CPPFCluster.h"
0007 #include "L1Trigger/L1TMuonCPPF/src/CPPFClusterizer.h"
0008 #include "L1Trigger/L1TMuonCPPF/src/CPPFMaskReClusterizer.h"
0009 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
0010 
0011 RecHitProcessor::RecHitProcessor() {}
0012 
0013 RecHitProcessor::~RecHitProcessor() {}
0014 
0015 void RecHitProcessor::processLook(const edm::Event &iEvent,
0016                                   const edm::EventSetup &iSetup,
0017                                   const edm::EDGetToken &recHitToken,
0018                                   const edm::EDGetToken &rpcDigiToken,
0019                                   const edm::EDGetToken &rpcDigiSimLinkToken,
0020                                   const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> &rpcGeomToken,
0021                                   std::vector<RecHitProcessor::CppfItem> &CppfVec1,
0022                                   l1t::CPPFDigiCollection &cppfDigis,
0023                                   const int MaxClusterSize) const {
0024   edm::Handle<RPCRecHitCollection> recHits;
0025   iEvent.getByToken(recHitToken, recHits);
0026 
0027   edm::Handle<RPCDigiCollection> rpcDigis;
0028   iEvent.getByToken(rpcDigiToken, rpcDigis);
0029 
0030   edm::Handle<edm::DetSetVector<RPCDigiSimLink>> theSimlinkDigis;
0031   iEvent.getByToken(rpcDigiSimLinkToken, theSimlinkDigis);
0032 
0033   const auto &rpcGeom = iSetup.getData(rpcGeomToken);
0034 
0035   for (const auto &&rpcdgIt : (*rpcDigis)) {
0036     const RPCDetId &rpcId = rpcdgIt.first;
0037     const RPCDigiCollection::Range &range = rpcdgIt.second;
0038     if (rpcId.region() == 0)
0039       continue;
0040 
0041     CPPFClusterizer clizer;
0042     CPPFClusterContainer tcls = clizer.doAction(range);
0043     CPPFMaskReClusterizer mrclizer;
0044     CPPFRollMask mask;
0045     CPPFClusterContainer cls = mrclizer.doAction(rpcId, tcls, mask);
0046 
0047     for (const auto &cl : cls) {
0048       int isValid = rpcDigis.isValid();
0049       int rawId = rpcId.rawId();
0050       int Bx = cl.bx();
0051       const int firststrip = cl.firstStrip();
0052       const int clustersize = cl.clusterSize();
0053       const int laststrip = cl.lastStrip();
0054       const RPCRoll *roll = rpcGeom.roll(rpcId);
0055       // Get Average Strip position
0056       const float fstrip = (roll->centreOfStrip(firststrip)).x();
0057       const float lstrip = (roll->centreOfStrip(laststrip)).x();
0058       const float centreOfCluster = (fstrip + lstrip) / 2;
0059       const double y = cl.hasY() ? cl.y() : 0;
0060       LocalPoint lPos(centreOfCluster, y, 0);
0061       if (roll->id().region() != 0) {
0062         const auto &topo = dynamic_cast<const TrapezoidalStripTopology &>(roll->topology());
0063         const double angle = topo.stripAngle((firststrip + laststrip) / 2.);
0064         const double x = centreOfCluster - y * std::tan(angle);
0065         lPos = LocalPoint(x, y, 0);
0066       }
0067       const BoundPlane &rollSurface = roll->surface();
0068       GlobalPoint gPos = rollSurface.toGlobal(lPos);
0069       float global_theta = emtf::rad_to_deg(gPos.theta().value());
0070       float global_phi = emtf::rad_to_deg(gPos.phi().value());
0071 
0072       // Establish the average position of the rechit
0073       int rechitstrip = firststrip;
0074 
0075       if (clustersize > 2) {
0076         int medium = 0;
0077         if (clustersize % 2 == 0)
0078           medium = 0.5 * (clustersize);
0079         else
0080           medium = 0.5 * (clustersize - 1);
0081         rechitstrip += medium;
0082       }
0083       if (clustersize > MaxClusterSize)
0084         continue;
0085       // This is just for test CPPFDigis with the RPC Geometry, It must be
0086       // "true" in the normal runs
0087       bool Geo = true;
0088       ////:::::::::::::::::::::::::::::::::::::::::::::::::
0089       // Set the EMTF Sector
0090       int EMTFsector1 = 0;
0091       int EMTFsector2 = 0;
0092 
0093       // sector 1
0094       if ((global_phi > 15.) && (global_phi <= 16.3)) {
0095         EMTFsector1 = 1;
0096         EMTFsector2 = 6;
0097       } else if ((global_phi > 16.3) && (global_phi <= 53.)) {
0098         EMTFsector1 = 1;
0099         EMTFsector2 = 0;
0100       } else if ((global_phi > 53.) && (global_phi <= 75.)) {
0101         EMTFsector1 = 1;
0102         EMTFsector2 = 2;
0103       }
0104       // sector 2
0105       else if ((global_phi > 75.) && (global_phi <= 76.3)) {
0106         EMTFsector1 = 1;
0107         EMTFsector2 = 2;
0108       } else if ((global_phi > 76.3) && (global_phi <= 113.)) {
0109         EMTFsector1 = 2;
0110         EMTFsector2 = 0;
0111       } else if ((global_phi > 113.) && (global_phi <= 135.)) {
0112         EMTFsector1 = 2;
0113         EMTFsector2 = 3;
0114       }
0115       // sector 3
0116       // less than 180
0117       else if ((global_phi > 135.) && (global_phi <= 136.3)) {
0118         EMTFsector1 = 2;
0119         EMTFsector2 = 3;
0120       } else if ((global_phi > 136.3) && (global_phi <= 173.)) {
0121         EMTFsector1 = 3;
0122         EMTFsector2 = 0;
0123       } else if ((global_phi > 173.) && (global_phi <= 180.)) {
0124         EMTFsector1 = 3;
0125         EMTFsector2 = 4;
0126       }
0127       // Greater than -180
0128       else if ((global_phi < -165.) && (global_phi >= -180.)) {
0129         EMTFsector1 = 3;
0130         EMTFsector2 = 4;
0131       }
0132       // Fourth sector
0133       else if ((global_phi > -165.) && (global_phi <= -163.7)) {
0134         EMTFsector1 = 3;
0135         EMTFsector2 = 4;
0136       } else if ((global_phi > -163.7) && (global_phi <= -127.)) {
0137         EMTFsector1 = 4;
0138         EMTFsector2 = 0;
0139       } else if ((global_phi > -127.) && (global_phi <= -105.)) {
0140         EMTFsector1 = 4;
0141         EMTFsector2 = 5;
0142       }
0143       // fifth sector
0144       else if ((global_phi > -105.) && (global_phi <= -103.7)) {
0145         EMTFsector1 = 4;
0146         EMTFsector2 = 5;
0147       } else if ((global_phi > -103.7) && (global_phi <= -67.)) {
0148         EMTFsector1 = 5;
0149         EMTFsector2 = 0;
0150       } else if ((global_phi > -67.) && (global_phi <= -45.)) {
0151         EMTFsector1 = 5;
0152         EMTFsector2 = 6;
0153       }
0154       // sixth sector
0155       else if ((global_phi > -45.) && (global_phi <= -43.7)) {
0156         EMTFsector1 = 5;
0157         EMTFsector2 = 6;
0158       } else if ((global_phi > -43.7) && (global_phi <= -7.)) {
0159         EMTFsector1 = 6;
0160         EMTFsector2 = 0;
0161       } else if ((global_phi > -7.) && (global_phi <= 15.)) {
0162         EMTFsector1 = 6;
0163         EMTFsector2 = 1;
0164       }
0165 
0166       double EMTFLink1 = 0.;
0167       double EMTFLink2 = 0.;
0168       std::vector<RecHitProcessor::CppfItem>::iterator cppf1;
0169       std::vector<RecHitProcessor::CppfItem>::iterator cppf;
0170       for (cppf1 = CppfVec1.begin(); cppf1 != CppfVec1.end(); cppf1++) {
0171         // Condition to save the CPPFDigi
0172         if (((*cppf1).rawId == rawId) && ((*cppf1).strip == rechitstrip)) {
0173           int old_strip = (*cppf1).strip;
0174           int before = 0;
0175           int after = 0;
0176 
0177           if (cppf1 != CppfVec1.begin())
0178             before = (*(cppf1 - 2)).strip;
0179           else if (cppf1 == CppfVec1.begin())
0180             before = (*cppf1).strip;
0181           if (cppf1 != CppfVec1.end())
0182             after = (*(cppf1 + 2)).strip;
0183           else if (cppf1 == CppfVec1.end())
0184             after = (*cppf1).strip;
0185           cppf = cppf1;
0186 
0187           if (clustersize == 2) {
0188             if (firststrip == 1) {
0189               if (before < after)
0190                 cppf = (cppf1 - 1);
0191               else if (before > after)
0192                 cppf = (cppf1 + 1);
0193             } else if (firststrip > 1) {
0194               if (before < after)
0195                 cppf = (cppf1 + 1);
0196               else if (before > after)
0197                 cppf = (cppf1 - 1);
0198             }
0199           }
0200           // Using the RPCGeometry
0201           if (Geo) {
0202             std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
0203                                                                             Bx,
0204                                                                             (*cppf).int_phi,
0205                                                                             (*cppf).int_theta,
0206                                                                             isValid,
0207                                                                             (*cppf).lb,
0208                                                                             (*cppf).halfchannel,
0209                                                                             EMTFsector1,
0210                                                                             EMTFLink1,
0211                                                                             old_strip,
0212                                                                             clustersize,
0213                                                                             global_phi,
0214                                                                             global_theta));
0215             std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
0216                                                                             Bx,
0217                                                                             (*cppf).int_phi,
0218                                                                             (*cppf).int_theta,
0219                                                                             isValid,
0220                                                                             (*cppf).lb,
0221                                                                             (*cppf).halfchannel,
0222                                                                             EMTFsector2,
0223                                                                             EMTFLink2,
0224                                                                             old_strip,
0225                                                                             clustersize,
0226                                                                             global_phi,
0227                                                                             global_theta));
0228 
0229             if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
0230               cppfDigis.push_back(*MainVariables1.get());
0231             } else if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
0232               cppfDigis.push_back(*MainVariables1.get());
0233               cppfDigis.push_back(*MainVariables2.get());
0234             } else if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
0235               continue;
0236             }
0237           }  // Geo is true
0238           else {
0239             global_phi = 0.;
0240             global_theta = 0.;
0241             std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
0242                                                                             Bx,
0243                                                                             (*cppf).int_phi,
0244                                                                             (*cppf).int_theta,
0245                                                                             isValid,
0246                                                                             (*cppf).lb,
0247                                                                             (*cppf).halfchannel,
0248                                                                             EMTFsector1,
0249                                                                             EMTFLink1,
0250                                                                             old_strip,
0251                                                                             clustersize,
0252                                                                             global_phi,
0253                                                                             global_theta));
0254             std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
0255                                                                             Bx,
0256                                                                             (*cppf).int_phi,
0257                                                                             (*cppf).int_theta,
0258                                                                             isValid,
0259                                                                             (*cppf).lb,
0260                                                                             (*cppf).halfchannel,
0261                                                                             EMTFsector2,
0262                                                                             EMTFLink2,
0263                                                                             old_strip,
0264                                                                             clustersize,
0265                                                                             global_phi,
0266                                                                             global_theta));
0267             if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
0268               cppfDigis.push_back(*MainVariables1.get());
0269             } else if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
0270               cppfDigis.push_back(*MainVariables1.get());
0271               cppfDigis.push_back(*MainVariables2.get());
0272             } else if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
0273               continue;
0274             }
0275           }
0276         }  // Condition to save the CPPFDigi
0277       }    // Loop over the LUTVector
0278     }      //end loop over cludters
0279   }        //end loop over digis
0280 }  //end processlook function
0281 
0282 void RecHitProcessor::process(const edm::Event &iEvent,
0283                               const edm::EventSetup &iSetup,
0284                               const edm::EDGetToken &recHitToken,
0285                               const edm::EDGetToken &rpcDigiToken,
0286                               const edm::EDGetToken &rpcDigiSimLinkToken,
0287                               const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> &rpcGeomToken,
0288                               l1t::CPPFDigiCollection &cppfDigis) const {
0289   // Get the RPC Geometry
0290   const auto &rpcGeom = iSetup.getData(rpcGeomToken);
0291 
0292   edm::Handle<RPCDigiCollection> rpcDigis;
0293   iEvent.getByToken(rpcDigiToken, rpcDigis);
0294 
0295   // Get the RecHits from the event
0296   edm::Handle<RPCRecHitCollection> recHits;
0297   iEvent.getByToken(recHitToken, recHits);
0298 
0299   for (const auto &&rpcdgIt : (*rpcDigis)) {
0300     const RPCDetId &rpcId = rpcdgIt.first;
0301     const RPCDigiCollection::Range &range = rpcdgIt.second;
0302     if (rpcId.region() == 0)
0303       continue;
0304 
0305     CPPFClusterizer clizer;
0306     CPPFClusterContainer tcls = clizer.doAction(range);
0307     CPPFMaskReClusterizer mrclizer;
0308     CPPFRollMask mask;
0309     CPPFClusterContainer cls = mrclizer.doAction(rpcId, tcls, mask);
0310 
0311     for (const auto &cl : cls) {
0312       int region = rpcId.region();
0313       int isValid = rpcDigis.isValid();
0314       //      int rawId = rpcId.rawId();
0315       int Bx = cl.bx();
0316       const int firststrip = cl.firstStrip();
0317       const int clustersize = cl.clusterSize();
0318       const int laststrip = cl.lastStrip();
0319       const RPCRoll *roll = rpcGeom.roll(rpcId);
0320       // Get Average Strip position
0321       const float fstrip = (roll->centreOfStrip(firststrip)).x();
0322       const float lstrip = (roll->centreOfStrip(laststrip)).x();
0323       const float centreOfCluster = (fstrip + lstrip) / 2;
0324       const double y = cl.hasY() ? cl.y() : 0;
0325       LocalPoint lPos(centreOfCluster, y, 0);
0326       if (roll->id().region() != 0) {
0327         const auto &topo = dynamic_cast<const TrapezoidalStripTopology &>(roll->topology());
0328         const double angle = topo.stripAngle((firststrip + laststrip) / 2.);
0329         const double x = centreOfCluster - y * std::tan(angle);
0330         lPos = LocalPoint(x, y, 0);
0331       }
0332       const BoundPlane &rollSurface = roll->surface();
0333       GlobalPoint gPos = rollSurface.toGlobal(lPos);
0334       float global_theta = emtf::rad_to_deg(gPos.theta().value());
0335       float global_phi = emtf::rad_to_deg(gPos.phi().value());
0336 
0337       // Endcap region only
0338       if (region != 0) {
0339         int int_theta =
0340             (region == -1 ? 180. * 32. / 36.5 : 0.) + (float)region * global_theta * 32. / 36.5 - 8.5 * 32 / 36.5;
0341         if (region == 1) {
0342           if (global_theta < 8.5)
0343             int_theta = 0;
0344           if (global_theta > 45.)
0345             int_theta = 31;
0346         } else if (region == -1) {
0347           if (global_theta < 135.)
0348             int_theta = 31;
0349           if (global_theta > 171.5)
0350             int_theta = 0;
0351         }
0352         // Local EMTF
0353         double local_phi = 0.;
0354         int EMTFsector1 = 0;
0355         int EMTFsector2 = 0;
0356 
0357         // sector 1
0358         if ((global_phi > 15.) && (global_phi <= 16.3)) {
0359           local_phi = global_phi - 15.;
0360           EMTFsector1 = 1;
0361           EMTFsector2 = 6;
0362         } else if ((global_phi > 16.3) && (global_phi <= 53.)) {
0363           local_phi = global_phi - 15.;
0364           EMTFsector1 = 1;
0365           EMTFsector2 = 0;
0366         } else if ((global_phi > 53.) && (global_phi <= 75.)) {
0367           local_phi = global_phi - 15.;
0368           EMTFsector1 = 1;
0369           EMTFsector2 = 2;
0370         }
0371         // sector 2
0372         else if ((global_phi > 75.) && (global_phi <= 76.3)) {
0373           local_phi = global_phi - 15.;
0374           EMTFsector1 = 1;
0375           EMTFsector2 = 2;
0376         } else if ((global_phi > 76.3) && (global_phi <= 113.)) {
0377           local_phi = global_phi - 75.;
0378           EMTFsector1 = 2;
0379           EMTFsector2 = 0;
0380         } else if ((global_phi > 113.) && (global_phi <= 135.)) {
0381           local_phi = global_phi - 75.;
0382           EMTFsector1 = 2;
0383           EMTFsector2 = 3;
0384         }
0385         // sector 3
0386         // less than 180
0387         else if ((global_phi > 135.) && (global_phi <= 136.3)) {
0388           local_phi = global_phi - 75.;
0389           EMTFsector1 = 2;
0390           EMTFsector2 = 3;
0391         } else if ((global_phi > 136.3) && (global_phi <= 173.)) {
0392           local_phi = global_phi - 135.;
0393           EMTFsector1 = 3;
0394           EMTFsector2 = 0;
0395         } else if ((global_phi > 173.) && (global_phi <= 180.)) {
0396           local_phi = global_phi - 135.;
0397           EMTFsector1 = 3;
0398           EMTFsector2 = 4;
0399         }
0400         // Greater than -180
0401         else if ((global_phi < -165.) && (global_phi >= -180.)) {
0402           local_phi = global_phi + 225.;
0403           EMTFsector1 = 3;
0404           EMTFsector2 = 4;
0405         }
0406         // Fourth sector
0407         else if ((global_phi > -165.) && (global_phi <= -163.7)) {
0408           local_phi = global_phi + 225.;
0409           EMTFsector1 = 3;
0410           EMTFsector2 = 4;
0411         } else if ((global_phi > -163.7) && (global_phi <= -127.)) {
0412           local_phi = global_phi + 165.;
0413           EMTFsector1 = 4;
0414           EMTFsector2 = 0;
0415         } else if ((global_phi > -127.) && (global_phi <= -105.)) {
0416           local_phi = global_phi + 165.;
0417           EMTFsector1 = 4;
0418           EMTFsector2 = 5;
0419         }
0420         // fifth sector
0421         else if ((global_phi > -105.) && (global_phi <= -103.7)) {
0422           local_phi = global_phi + 165.;
0423           EMTFsector1 = 4;
0424           EMTFsector2 = 5;
0425         } else if ((global_phi > -103.7) && (global_phi <= -67.)) {
0426           local_phi = global_phi + 105.;
0427           EMTFsector1 = 5;
0428           EMTFsector2 = 0;
0429         } else if ((global_phi > -67.) && (global_phi <= -45.)) {
0430           local_phi = global_phi + 105.;
0431           EMTFsector1 = 5;
0432           EMTFsector2 = 6;
0433         }
0434         // sixth sector
0435         else if ((global_phi > -45.) && (global_phi <= -43.7)) {
0436           local_phi = global_phi + 105.;
0437           EMTFsector1 = 5;
0438           EMTFsector2 = 6;
0439         } else if ((global_phi > -43.7) && (global_phi <= -7.)) {
0440           local_phi = global_phi + 45.;
0441           EMTFsector1 = 6;
0442           EMTFsector2 = 0;
0443         } else if ((global_phi > -7.) && (global_phi <= 15.)) {
0444           local_phi = global_phi + 45.;
0445           EMTFsector1 = 6;
0446           EMTFsector2 = 1;
0447         }
0448 
0449         int int_phi = int((local_phi + 22.0) * 15. + .5);
0450         double EMTFLink1 = 0.;
0451         double EMTFLink2 = 0.;
0452         double lb = 0.;
0453         double halfchannel = 0.;
0454 
0455         // Invalid hit
0456         if (isValid == 0)
0457           int_phi = 2047;
0458         // Right integers range
0459         assert(0 <= int_phi && int_phi < 1250);
0460         assert(0 <= int_theta && int_theta < 32);
0461 
0462         std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
0463                                                                         Bx,
0464                                                                         int_phi,
0465                                                                         int_theta,
0466                                                                         isValid,
0467                                                                         lb,
0468                                                                         halfchannel,
0469                                                                         EMTFsector1,
0470                                                                         EMTFLink1,
0471                                                                         firststrip,
0472                                                                         clustersize,
0473                                                                         global_phi,
0474                                                                         global_theta));
0475         std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
0476                                                                         Bx,
0477                                                                         int_phi,
0478                                                                         int_theta,
0479                                                                         isValid,
0480                                                                         lb,
0481                                                                         halfchannel,
0482                                                                         EMTFsector2,
0483                                                                         EMTFLink2,
0484                                                                         firststrip,
0485                                                                         clustersize,
0486                                                                         global_phi,
0487                                                                         global_theta));
0488         if (int_theta == 31)
0489           continue;
0490         if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
0491           cppfDigis.push_back(*MainVariables1.get());
0492         }
0493         if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
0494           cppfDigis.push_back(*MainVariables1.get());
0495           cppfDigis.push_back(*MainVariables2.get());
0496         }
0497         if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
0498           continue;
0499         }
0500       }  // No barrel hits
0501     }    //end loop over clusters
0502   }      //end loop over digis
0503 }  // End function: void RecHitProcessor::process()