Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "RecoLocalTracker/SiPhase2VectorHitBuilder/interface/VectorHitBuilderAlgorithmBase.h"
0002 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0003 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0004 #include "RecoLocalTracker/Records/interface/TkPhase2OTCPERecord.h"
0005 
0006 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0007 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0008 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0009 
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 
0012 VectorHitBuilderAlgorithmBase::VectorHitBuilderAlgorithmBase(
0013     const edm::ParameterSet& conf,
0014     const TrackerGeometry* tkGeomProd,
0015     const TrackerTopology* tkTopoProd,
0016     const ClusterParameterEstimator<Phase2TrackerCluster1D>* cpeProd)
0017     : tkGeom_(tkGeomProd),
0018       tkTopo_(tkTopoProd),
0019       cpe_(cpeProd),
0020       nMaxVHforeachStack_(conf.getParameter<int>("maxVectorHitsInAStack")),
0021       barrelCut_(conf.getParameter<std::vector<double> >("BarrelCut")),
0022       endcapCut_(conf.getParameter<std::vector<double> >("EndcapCut")),
0023       cpeTag_(conf.getParameter<edm::ESInputTag>("CPE")) {}
0024 
0025 double VectorHitBuilderAlgorithmBase::computeParallaxCorrection(const PixelGeomDetUnit* geomDetUnit_low,
0026                                                                 const Point3DBase<float, LocalTag>& lPosClu_low,
0027                                                                 const PixelGeomDetUnit* geomDetUnit_upp,
0028                                                                 const Point3DBase<float, LocalTag>& lPosClu_upp) const {
0029   double parallCorr = 0.0;
0030   Global3DPoint origin(0, 0, 0);
0031   Global3DPoint gPosClu_low = geomDetUnit_low->surface().toGlobal(lPosClu_low);
0032   GlobalVector gV = gPosClu_low - origin;
0033   LogTrace("VectorHitsBuilderValidation") << " global vector passing to the origin:" << gV;
0034 
0035   LocalVector lV = geomDetUnit_low->surface().toLocal(gV);
0036   LogTrace("VectorHitsBuilderValidation")
0037       << " local vector passing to the origin (in the lower detector system of reference):" << lV;
0038   LocalVector lV_norm = lV / lV.z();
0039   LogTrace("VectorHitsBuilderValidation")
0040       << " normalized local vector passing to the origin (in low the lower detector system of reference):" << lV_norm;
0041 
0042   Global3DPoint gPosClu_upp = geomDetUnit_upp->surface().toGlobal(lPosClu_upp);
0043   Local3DPoint lPosClu_uppInLow = geomDetUnit_low->surface().toLocal(gPosClu_upp);
0044   parallCorr = lV_norm.x() * lPosClu_uppInLow.z();
0045 
0046   return parallCorr;
0047 }
0048 
0049 void VectorHitBuilderAlgorithmBase::printClusters(const edmNew::DetSetVector<Phase2TrackerCluster1D>& clusters) const {
0050   int nCluster = 0;
0051   for (const auto& DSViter : clusters) {
0052     // Loop over the clusters in the detector unit
0053     for (const auto& clustIt : DSViter) {
0054       nCluster++;
0055       // get the detector unit's id
0056       const GeomDetUnit* geomDetUnit(tkGeom_->idToDetUnit(DSViter.detId()));
0057       if (!geomDetUnit)
0058         return;
0059       printCluster(geomDetUnit, &clustIt);
0060     }
0061   }
0062   LogDebug("VectorHitBuilder") << " Number of input clusters: " << nCluster << std::endl;
0063 }
0064 
0065 void VectorHitBuilderAlgorithmBase::printCluster(const GeomDet* geomDetUnit,
0066                                                  const Phase2TrackerCluster1D* clustIt) const {
0067   if (!geomDetUnit)
0068     return;
0069   const PixelGeomDetUnit* pixelGeomDetUnit = dynamic_cast<const PixelGeomDetUnit*>(geomDetUnit);
0070   const PixelTopology& topol = pixelGeomDetUnit->specificTopology();
0071   if (!pixelGeomDetUnit)
0072     return;
0073 
0074   unsigned int layer = tkTopo_->layer(geomDetUnit->geographicalId());
0075   unsigned int module = tkTopo_->module(geomDetUnit->geographicalId());
0076   LogTrace("VectorHitBuilder") << "Layer:" << layer << " and DetId: " << geomDetUnit->geographicalId().rawId()
0077                                << std::endl;
0078   TrackerGeometry::ModuleType mType = tkGeom_->getDetectorType(geomDetUnit->geographicalId());
0079   if (mType == TrackerGeometry::ModuleType::Ph2PSP)
0080     LogTrace("VectorHitBuilder") << "Pixel cluster (module:" << module << ") " << std::endl;
0081   else if (mType == TrackerGeometry::ModuleType::Ph2SS || mType == TrackerGeometry::ModuleType::Ph2PSS)
0082     LogTrace("VectorHitBuilder") << "Strip cluster (module:" << module << ") " << std::endl;
0083   else
0084     LogTrace("VectorHitBuilder") << "no module?!" << std::endl;
0085   LogTrace("VectorHitBuilder") << "with pitch:" << topol.pitch().first << " , " << topol.pitch().second << std::endl;
0086   LogTrace("VectorHitBuilder") << " and width:" << pixelGeomDetUnit->surface().bounds().width()
0087                                << " , lenght:" << pixelGeomDetUnit->surface().bounds().length() << std::endl;
0088 
0089   auto&& lparams = cpe_->localParameters(*clustIt, *pixelGeomDetUnit);
0090   Global3DPoint gparams = pixelGeomDetUnit->surface().toGlobal(lparams.first);
0091 
0092   LogTrace("VectorHitBuilder") << "\t global pos " << gparams << std::endl;
0093   LogTrace("VectorHitBuilder") << "\t local  pos " << lparams.first << "with err " << lparams.second << std::endl;
0094   LogTrace("VectorHitBuilder") << std::endl;
0095 
0096   return;
0097 }