CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTracker/ConversionSeedGenerators/src/HitQuadrupletGeneratorFromLayerPairForPhotonConversion.cc

Go to the documentation of this file.
00001 #include "RecoTracker/ConversionSeedGenerators/interface/HitQuadrupletGeneratorFromLayerPairForPhotonConversion.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00005 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00006 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00007 
00008 #include "RecoTracker/TkTrackingRegions/interface/HitRZCompatibility.h"
00009 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
00010 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionBase.h"
00011 #include "RecoTracker/TkHitPairs/interface/OrderedHitPairs.h"
00012 #include "RecoTracker/TkHitPairs/src/InnerDeltaPhi.h"
00013 
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "RecoTracker/TkSeedingLayers/interface/SeedingLayer.h"
00017 
00018 using namespace GeomDetEnumerators;
00019 using namespace ctfseeding;
00020 using namespace std;
00021 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00022 
00023 //#define mydebug_QSeed
00024 
00025 typedef PixelRecoRange<float> Range;
00026 template<class T> inline T sqr( T t) {return t*t;}
00027 
00028 
00029 
00030 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00031 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00032 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00033 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00034 #include "FWCore/Framework/interface/ESHandle.h"
00035 #include "DataFormats/Math/interface/deltaPhi.h"
00036 
00037 HitQuadrupletGeneratorFromLayerPairForPhotonConversion::HitQuadrupletGeneratorFromLayerPairForPhotonConversion(
00038                                                              const Layer& inner, 
00039                                                              const Layer& outer, 
00040                                                              LayerCacheType* layerCache,
00041                                                              unsigned int nSize,
00042                                                              unsigned int max)
00043   : HitPairGenerator(nSize),
00044     theLayerCache(*layerCache), theOuterLayer(outer), theInnerLayer(inner)
00045 {
00046   theMaxElement=max;
00047   ss = new std::stringstream;
00048 }
00049 
00050 void HitQuadrupletGeneratorFromLayerPairForPhotonConversion::hitPairs(const TrackingRegion & region, OrderedHitPairs & result,
00051                                                                       const edm::Event& event, const edm::EventSetup& es)
00052 {
00053 
00054   ss->str("");
00055 
00056   typedef OrderedHitPair::InnerRecHit InnerHit;
00057   typedef OrderedHitPair::OuterRecHit OuterHit;
00058   typedef RecHitsSortedInPhi::Hit Hit;
00059 
00060   size_t totCountP2=0, totCountP1=0, totCountM2=0, totCountM1=0, selCount=0;
00061 #ifdef mydebug_QSeed
00062   (*ss) << "In " << theInnerLayer.name() << " Out " << theOuterLayer.name() << std::endl;
00063 #endif
00064 
00065   /*get hit sorted in phi for each layer: NB: doesn't apply any region cut*/
00066   const RecHitsSortedInPhi & innerHitsMap = theLayerCache(&theInnerLayer, region, event, es);
00067   if (innerHitsMap.empty()) return;
00068  
00069   const RecHitsSortedInPhi& outerHitsMap = theLayerCache(&theOuterLayer, region, event, es);
00070   if (outerHitsMap.empty()) return;
00071   /*----------------*/
00072 
00073   /*This object will check the compatibility of the his in phi among the two layers. */
00074   //InnerDeltaPhi deltaPhi(*theInnerLayer.detLayer(), region, es);
00075 
00076   vector<RecHitsSortedInPhi::Hit> innerHits;
00077   //  float outerPhimin, outerPhimax;
00078   float innerPhimin, innerPhimax;
00079   float maxDeltaPhi=1.; //sguazz
00080   //float maxDeltaPhi=1.;
00081 
00082   RecHitsSortedInPhi::Range outerHits = outerHitsMap.all();
00083 
00084   RecHitsSortedInPhi::HitIter nextoh;
00085   for (RecHitsSortedInPhi::HitIter oh = outerHits.first; oh!= outerHits.second; ++oh) { 
00086     RecHitsSortedInPhi::Hit ohit = (*oh).hit();
00087     GlobalPoint oPos = ohit->globalPosition();  
00088 
00089     totCountP2++;
00090     const HitRZCompatibility *checkRZ = region.checkRZ(theOuterLayer.detLayer(), ohit, es);    
00091     for(nextoh=oh+1;nextoh!=outerHits.second; ++nextoh){
00092       
00093       RecHitsSortedInPhi::Hit nohit = (*nextoh).hit();
00094       GlobalPoint noPos = nohit->globalPosition();  
00095       
00096     
00097 #ifdef mydebug_QSeed
00098       (*ss) << "\toPos " << oPos << " r " << oPos.perp() << " phi " << oPos.phi() <<  " cotTheta " << oPos.z()/oPos.perp() << std::endl;
00099       (*ss) << "\tnoPos " << noPos << " r " << noPos.perp() << " phi " << noPos.phi() <<  " cotTheta " << noPos.z()/noPos.perp() << std::endl;
00100       (*ss) << "\tDeltaPhi " << reco::deltaPhi(noPos.phi(),oPos.phi()) << std::endl;
00101 #endif
00102     
00103     if(fabs(reco::deltaPhi(noPos.phi(),oPos.phi()))>maxDeltaPhi)  break;
00104 
00105     totCountM2++;
00106     
00107     /*Check the compatibility of the ohit with the eta of the seeding track*/
00108     if(failCheckRZCompatibility(nohit,*theOuterLayer.detLayer(),checkRZ,region)) continue;
00109 
00110     /*  
00111     //Do I need this? it uses a compatibility that probably I wouldn't 
00112     //Removing for the time being
00113 
00114     PixelRecoRange<float> phiRange = deltaPhi( oPos.perp(), oPos.phi(), oPos.z(), nSigmaPhi*(ohit->errorGlobalRPhi()));    
00115     if (phiRange.empty()) continue;
00116     */
00117 
00118   
00119     /*Get only the inner hits compatible with the conversion region*/
00120     innerPhimin=ohit->globalPosition().phi();
00121     innerPhimax=nohit->globalPosition().phi();
00122     // checkPhiRange(innerPhimin,innerPhimax);
00123 
00124     innerHits.clear();    
00125     innerHitsMap.hits(innerPhimin, innerPhimax, innerHits);
00126 
00127 #ifdef mydebug_QSeed
00128     (*ss) << "\tiphimin, iphimax " << innerPhimin << " " << innerPhimax << std::endl;
00129 #endif    
00130 
00131     const HitRZCompatibility *checkRZb = region.checkRZ(theInnerLayer.detLayer(),  ohit, es);
00132     const HitRZCompatibility *checkRZc = region.checkRZ(theInnerLayer.detLayer(), nohit, es);
00133 
00134     /*Loop on inner hits*/
00135     vector<RecHitsSortedInPhi::Hit>::const_iterator ieh = innerHits.end();
00136     for ( vector<RecHitsSortedInPhi::Hit>::const_iterator ih=innerHits.begin(); ih < ieh; ++ih) {  
00137       RecHitsSortedInPhi::Hit ihit = *ih;
00138       
00139 #ifdef mydebug_QSeed
00140       GlobalPoint innPos = (*ih)->globalPosition();
00141       (*ss) << "\toPos " << oPos << " r " << oPos.perp() << " phi " << oPos.phi() <<  " cotTheta " << oPos.z()/oPos.perp() << std::endl;
00142       (*ss) << "\tnoPos " << noPos << " r " << noPos.perp() << " phi " << noPos.phi() <<  " cotTheta " << noPos.z()/noPos.perp() << std::endl;
00143       (*ss) << "\tinnPos " << innPos <<  " r " << innPos.perp() << " phi " << innPos.phi() << " cotTheta " << innPos.z()/innPos.perp() <<  std::endl;
00144 #endif
00145 
00146       totCountP1++;
00147       
00148       /*Check the compatibility of the ihit with the two outer hits*/
00149       if(failCheckRZCompatibility(ihit,*theInnerLayer.detLayer(),checkRZb,region) 
00150          || 
00151          failCheckRZCompatibility(ihit,*theInnerLayer.detLayer(),checkRZc,region) ) continue;
00152       
00153       
00154       for ( vector<RecHitsSortedInPhi::Hit>::const_iterator nextih=ih+1; nextih != ieh; ++nextih) {  
00155         RecHitsSortedInPhi::Hit nihit = *nextih;
00156 
00157 
00158 #ifdef mydebug_QSeed
00159         GlobalPoint ninnPos = (*nextih)->globalPosition();
00160         (*ss) << "\toPos " << oPos << " r " << oPos.perp() << " phi " << oPos.phi() <<  " cotTheta " << oPos.z()/oPos.perp() << std::endl;
00161         (*ss) << "\tnoPos " << noPos << " r " << noPos.perp() << " phi " << noPos.phi() <<  " cotTheta " << noPos.z()/noPos.perp() << std::endl;
00162         (*ss) << "\tinnPos " << innPos <<  " r " << innPos.perp() << " phi " << innPos.phi() << " cotTheta " << innPos.z()/innPos.perp() <<  std::endl;
00163         (*ss) << "\tninnPos " << ninnPos <<  " r " << ninnPos.perp() << " phi " << ninnPos.phi() << " cotTheta " << ninnPos.z()/ninnPos.perp() <<  std::endl;
00164 #endif
00165 
00166         totCountM1++;
00167 
00168         /*Check the compatibility of the nihit with the two outer hits*/
00169         if(failCheckRZCompatibility(nihit,*theInnerLayer.detLayer(),checkRZb,region) 
00170            || 
00171            failCheckRZCompatibility(nihit,*theInnerLayer.detLayer(),checkRZc,region) ) continue;
00172         
00173         /*Sguazz modifica qui*/
00174         if(failCheckSlopeTest(ohit,nohit,ihit,nihit,region)) continue;
00175 
00176         if (theMaxElement!=0 && result.size() >= theMaxElement){
00177           result.clear();
00178           edm::LogError("TooManyQuads")<<"number of Quad combinations exceed maximum, no quads produced";
00179           delete checkRZ;
00180           delete checkRZb;
00181           delete checkRZc;        
00182 #ifdef mydebug_QSeed
00183           (*ss) << "In " << theInnerLayer.name() << " Out " << theOuterLayer.name()
00184                 << "\tP2 " << totCountP2 
00185                 << "\tM2 " << totCountM2 
00186                 << "\tP1 " << totCountP1 
00187                 << "\tM1 " << totCountM1 
00188                 << "\tsel " << selCount
00189                 << std::endl;
00190           std::cout << (*ss).str();
00191 #endif
00192           return;
00193         }
00194         selCount++;
00195         result.push_back( OrderedHitPair( ihit, ohit) );
00196         result.push_back( OrderedHitPair( nihit, nohit) );
00197         //#ifdef mydebug_QSeed
00198         //(*ss) << "sizeOfresul " << result.size() << std::endl;
00199         //#endif
00200       }
00201     }
00202     delete checkRZb;
00203     delete checkRZc;
00204     }
00205     delete checkRZ;
00206   }
00207 #ifdef mydebug_QSeed
00208   (*ss) << "In " << theInnerLayer.name() << " Out " << theOuterLayer.name() 
00209         << "\tP2 " << totCountP2 
00210         << "\tM2 " << totCountM2 
00211         << "\tP1 " << totCountP1 
00212         << "\tM1 " << totCountM1 
00213         << "\tsel " << selCount
00214         << std::endl;
00215   std::cout << (*ss).str();
00216 #endif
00217 }
00218 
00219 
00220 
00221 
00222 bool HitQuadrupletGeneratorFromLayerPairForPhotonConversion::
00223 failCheckRZCompatibility(const RecHitsSortedInPhi::Hit& hit,const DetLayer& layer,const HitRZCompatibility *checkRZ,const TrackingRegion & region){
00224 
00225   if(!checkRZ) {
00226 #ifdef mydebug_QSeed
00227     (*ss) << "*******\nNo valid checkRZ\n*******" << std::endl;
00228 #endif
00229     return true;
00230   }
00231 
00232   static const float nSigmaRZ = std::sqrt(12.f);
00233   float r_reduced = std::sqrt( sqr(hit->globalPosition().x()-region.origin().x())+sqr(hit->globalPosition().y()-region.origin().y()));
00234   Range allowed;
00235   Range hitRZ;
00236   if (layer.location() == barrel) {
00237     allowed = checkRZ->range(r_reduced);
00238     float zErr = nSigmaRZ * hit->errorGlobalZ();
00239     hitRZ = Range(hit->globalPosition().z()-zErr, hit->globalPosition().z()+zErr);
00240   } else {
00241     allowed = checkRZ->range(hit->globalPosition().z());
00242     float rErr = nSigmaRZ * hit->errorGlobalR();
00243     hitRZ = Range(r_reduced-rErr, r_reduced+rErr);
00244   }
00245   Range crossRange = allowed.intersection(hitRZ);
00246 
00247 #ifdef mydebug_QSeed      
00248   (*ss) 
00249     << "\n\t\t allowed Range " << allowed.min() << " \t, " << allowed.max() 
00250     << "\n\t\t hitRz   Range " << hitRZ.min()   << " \t, " << hitRZ.max() 
00251     << "\n\t\t Cross   Range " << crossRange.min()   << " \t, " << crossRange.max() 
00252     << std::endl;
00253   if( !crossRange.empty())
00254     (*ss) << "\n\t\t !!!!ACCEPTED!!! \n\n";
00255 #endif
00256   
00257   return crossRange.empty();
00258 }
00259 
00260 
00261 bool HitQuadrupletGeneratorFromLayerPairForPhotonConversion::
00262 failCheckSlopeTest(const RecHitsSortedInPhi::Hit & ohit, const RecHitsSortedInPhi::Hit & nohit, const RecHitsSortedInPhi::Hit & ihit, const RecHitsSortedInPhi::Hit & nihit, const TrackingRegion & region){
00263 
00264   double r[5], z[5], ez[5];
00265   //  double pr[2], pz[2], e2pz[2], mr[2], mz[2], e2mz[2];
00266 
00267   //
00268   //Hits
00269   r[0] = ohit->globalPosition().perp();
00270   z[0] = ohit->globalPosition().z();
00271   ez[0] = getEffectiveErrorOnZ(ohit, region);
00272   //
00273   r[1] = nohit->globalPosition().perp();
00274   z[1] = nohit->globalPosition().z();
00275   ez[1] = getEffectiveErrorOnZ(nohit, region);
00276   //
00277   r[2] = nihit->globalPosition().perp();
00278   z[2] = nihit->globalPosition().z();
00279   ez[2] = getEffectiveErrorOnZ(nihit, region);
00280   //
00281   r[3] = ihit->globalPosition().perp();
00282   z[3] = ihit->globalPosition().z();
00283   ez[3] = getEffectiveErrorOnZ(ihit, region);
00284   //
00285   //R (r) ordering of the 4 hit arrays
00286   bubbleSortVsR(4, r, z, ez);
00287   //
00288   //Vertex
00289   r[4] = region.origin().perp();
00290   z[4] = region.origin().z();
00291   double vError = region.originZBound();
00292   if ( vError > 15. ) vError = 1.;
00293   ez[4] = 3.*vError;
00294 
00295   //
00296   //Sequence of checks
00297   //
00298   //Inner segment == vertex
00299   double rInn = r[4];
00300   double zInnMin = z[4]-ez[4];
00301   double zInnMax = z[4]+ez[4];
00302   //
00303   // Int == 2, Out == 3
00304   double rOut = r[3];
00305   double zOutMin = z[3]-ez[3];
00306   double zOutMax = z[3]+ez[3];
00307   double rInt = r[2];
00308   double zIntMin = z[2]-ez[2];
00309   double zIntMax = z[2]+ez[2];
00310   if ( failCheckSegmentZCompatibility(rInn, zInnMin, zInnMax, 
00311                                       rInt, zIntMin, zIntMax,
00312                                       rOut, zOutMin, zOutMax) ) return true;
00313   //
00314   // Int == 1, Out == 2 (with updated limits)
00315   rOut = rInt;
00316   zOutMin = zIntMin;
00317   zOutMax = zIntMax;
00318   rInt    = r[1];
00319   zIntMin = z[1]-ez[1];
00320   zIntMax = z[1]+ez[1];
00321   if ( failCheckSegmentZCompatibility(rInn, zInnMin, zInnMax, 
00322                                       rInt, zIntMin, zIntMax,
00323                                       rOut, zOutMin, zOutMax) ) return true;
00324   //
00325   // Int == 0, Out == 1 (with updated limits)
00326   rOut = rInt;
00327   zOutMin = zIntMin;
00328   zOutMax = zIntMax;
00329   rInt    = r[0];
00330   zIntMin = z[0]-ez[0];
00331   zIntMax = z[0]+ez[0];
00332   if ( failCheckSegmentZCompatibility(rInn, zInnMin, zInnMax, 
00333                                       rInt, zIntMin, zIntMax,
00334                                       rOut, zOutMin, zOutMax) ) return true;
00335 
00336   //
00337   // Test is ok!!!
00338   return false;
00339 
00340 }
00341 
00342 
00343 double HitQuadrupletGeneratorFromLayerPairForPhotonConversion::verySimpleFit(int size, double* ax, double* ay, double* e2y, double& p0, double& e2p0, double& p1){
00344   
00345   //#include "RecoTracker/ConversionSeedGenerators/interface/verySimpleFit.icc"
00346   return 0;
00347 }
00348 
00349 
00350 double HitQuadrupletGeneratorFromLayerPairForPhotonConversion::getSqrEffectiveErrorOnZ(const RecHitsSortedInPhi::Hit & hit, const TrackingRegion & region){
00351 
00352   //
00353   //Fit-wise the effective error on Z is approximately the sum in quadrature of the error on Z 
00354   //and the error on R correctly projected by using hit-vertex direction
00355 
00356   double sqrProjFactor = sqr((hit->globalPosition().z()-region.origin().z())/(hit->globalPosition().perp()-region.origin().perp()));
00357   return (hit->globalPositionError().czz()+sqrProjFactor*hit->globalPositionError().rerr(hit->globalPosition()));
00358 
00359 }
00360 
00361 double HitQuadrupletGeneratorFromLayerPairForPhotonConversion::getEffectiveErrorOnZ(const RecHitsSortedInPhi::Hit & hit, const TrackingRegion & region){
00362 
00363   //
00364   //Fit-wise the effective error on Z is approximately the sum in quadrature of the error on Z 
00365   //and the error on R correctly projected by using hit-vertex direction
00366   double sqrProjFactor = sqr((hit->globalPosition().z()-region.origin().z())/(hit->globalPosition().perp()-region.origin().perp()));
00367   double effErr = sqrt(hit->globalPositionError().czz()+sqrProjFactor*hit->globalPositionError().rerr(hit->globalPosition())); 
00368   if ( effErr>2. ) {
00369     effErr*=1.8; //Single Side error is strip length * sqrt( 12.) = 3.46
00370     //empirically found that the error on ss hits value it is already twice as large (two sigmas)!
00371     //Multiply by sqrt(12.)/2.=1.73 to have effErr equal to the strip lenght (1.8 to allow for some margin)
00372     //effErr*=2.5; //Used in some tests
00373   } else {
00374     effErr*=2.; //Tight //Double side... allowing for 2 sigma variation 
00375     //effErr*=5.; //Loose //Double side... allowing for 2 sigma variation 
00376   }
00377   return effErr;
00378 
00379 }
00380 
00381 void HitQuadrupletGeneratorFromLayerPairForPhotonConversion::bubbleSortVsR(int n, double* ar, double* az, double* aez){
00382   bool swapped = true;
00383   int j = 0;
00384   double tmpr, tmpz, tmpez;
00385   while (swapped) {
00386     swapped = false;
00387     j++;
00388     for (int i = 0; i < n - j; i++) {
00389       if (  ar[i] > ar[i+1] ) {
00390         tmpr = ar[i];
00391         ar[i] = ar[i + 1];
00392         ar[i + 1] = tmpr;
00393         tmpz = az[i];
00394         az[i] = az[i + 1];
00395         az[i + 1] = tmpz;
00396         tmpez = aez[i];
00397         aez[i] = aez[i + 1];
00398         aez[i + 1] = tmpez;
00399         swapped = true;
00400       }
00401     }
00402   }
00403 }
00404 
00405 bool HitQuadrupletGeneratorFromLayerPairForPhotonConversion::
00406 failCheckSegmentZCompatibility(double &rInn, double &zInnMin, double &zInnMax,
00407                                double &rInt, double &zIntMin, double &zIntMax,
00408                                double &rOut, double &zOutMin, double &zOutMax){
00409   //
00410   // Check the compatibility in z of an INTermediate segment between an INNer segment and an OUTer segment;
00411   // when true is returned zIntMin and zIntMax are replaced with allowed range values
00412   
00413   //Left side
00414   double zLeft = getZAtR(rInn, zInnMin, rInt, rOut, zOutMin);
00415   if ( zIntMax < zLeft ) return true; 
00416   //Right side
00417   double zRight = getZAtR(rInn, zInnMax, rInt, rOut, zOutMax);
00418   if ( zIntMin > zRight ) return true; 
00419   if ( zIntMin < zLeft && zIntMax < zRight ) {
00420     zIntMax = zLeft;
00421     return false;
00422   }
00423   if ( zIntMin > zLeft && zIntMax > zRight ) {
00424     zIntMax = zRight;
00425     return false;
00426   }
00427   
00428   //Segment is fully contained
00429   return false;
00430 
00431 }
00432 
00433 double HitQuadrupletGeneratorFromLayerPairForPhotonConversion::
00434 getZAtR(double &rInn, double &zInn,
00435         double &r, 
00436         double &rOut, double &zOut){
00437 
00438   //    z - zInn      r - rInn                                  r - rInn 
00439   //  ----------- = ----------- ==> z = zInn + (zOut - zInn) * -----------
00440   //  zOut - zInn   rOut - rInn                                rOut - rInn
00441 
00442   
00443   return zInn + (zOut - zInn)*(r - rInn)/(rOut - rInn);
00444 
00445 }