CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/SUSYBSMAnalysis/HSCP/src/BetaCalculatorRPC.cc

Go to the documentation of this file.
00001 #include "SUSYBSMAnalysis/HSCP/interface/BetaCalculatorRPC.h"
00002 
00003 using namespace susybsm;
00004 
00005 
00006 BetaCalculatorRPC::BetaCalculatorRPC(const edm::ParameterSet& iConfig){
00007 
00008   rpcRecHitsLabel = iConfig.getParameter<edm::InputTag>("rpcRecHits");
00009 
00010 }
00011 
00012 void BetaCalculatorRPC::algo(std::vector<susybsm::RPCHit4D> HSCPRPCRecHits){
00013 
00014   int lastbx=-7;
00015   bool outOfTime = false;
00016   bool increasing = true;
00017   bool anydifferentzero = true;
00018   bool anydifferentone = true;
00019   
00020   //std::cout<<"Inside BetaCalculatorRPC \t \t Preliminar loop on the RPCHit4D!!!"<<std::endl;
00021 
00022   std::sort(HSCPRPCRecHits.begin(), HSCPRPCRecHits.end()); //Organizing them
00023 
00024   for(std::vector<susybsm::RPCHit4D>::iterator point = HSCPRPCRecHits.begin(); point < HSCPRPCRecHits.end(); ++point) {
00025     outOfTime |= (point->bx!=0); //condition 1: at least one measurement must have BX!=0
00026     increasing &= (point->bx>=lastbx); //condition 2: BX must be increase when going inside-out.
00027     anydifferentzero &= (!point->bx==0); //to check one knee withoutzeros
00028     anydifferentone &= (!point->bx==1); //to check one knee withoutones
00029     lastbx = point->bx;
00030     //float r=point->gp.mag();
00031     //std::cout<<"Inside BetaCalculatorRPC \t \t  r="<<r<<" phi="<<point->gp.phi()<<" eta="<<point->gp.eta()<<" bx="<<point->bx<<" outOfTime"<<outOfTime<<" increasing"<<increasing<<" anydifferentzero"<<anydifferentzero<<std::endl;
00032   }
00033   
00034   bool Candidate = (outOfTime&&increasing);
00035 
00036   // here we should get some pattern-based estimate
00037 
00038   //Counting knees
00039 
00040   float delay=12.5;
00041   lastbx=-7; //already declared for the preliminar loop
00042   int knees=0;
00043   float maginknee = 0;
00044   float maginfirstknee = 0;
00045   for(std::vector<susybsm::RPCHit4D>::iterator point = HSCPRPCRecHits.begin(); point < HSCPRPCRecHits.end(); ++point) {
00046     if(lastbx==-7){
00047       maginfirstknee = point->gp.mag();
00048     }else if((lastbx!=point->bx)){
00049       //std::cout<<"Inside BetaCalculatorRPC \t \t \t one knee between"<<lastbx<<point->bx<<std::endl;
00050       maginknee=point->gp.mag();
00051       knees++;
00052     }
00053     lastbx=point->bx;
00054   }
00055       
00056   if(knees==0){
00057     //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t knees="<<knees<<std::endl;
00058     betavalue=maginfirstknee/(25.-delay+maginfirstknee/30.)/30.;
00059   }else if(knees==1){
00060     float betavalue1=0;
00061     float betavalue2=0;
00062     //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t knees="<<knees<<std::endl;
00063     //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t anydifferentzero="<<anydifferentzero<<" anydifferentone="<<anydifferentone<<std::endl;
00064     if(!anydifferentzero){
00065       betavalue=maginknee/(25-delay+maginknee/30.)/30.;
00066     }else if(!anydifferentone){//i.e non zeros and no ones
00067       betavalue=maginknee/(50-delay+maginknee/30.)/30.;
00068     }else{
00069       betavalue1=maginknee/(25-delay+maginknee/30.)/30.;
00070       float dr =(maginknee-maginfirstknee);
00071       betavalue2 = dr/(25.-delay+dr/30.);
00072       //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t \t not zero neither ones betavalue1="<<betavalue1<<" betavalue2="<<betavalue2<<std::endl;
00073       betavalue = (betavalue1 + betavalue2)*0.5;
00074     }
00075   }else if(knees==2){
00076     //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t knees="<<knees<<std::endl;
00077     knees=0;
00078     float betavalue1=0;
00079     float betavalue2=0;
00080     lastbx=-7;
00081     //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t looping again on the RPCRecHits4D="<<knees<<std::endl;
00082     for(std::vector<susybsm::RPCHit4D>::iterator point = HSCPRPCRecHits.begin(); point < HSCPRPCRecHits.end(); ++point) {
00083       if(lastbx==-7){
00084         maginfirstknee = point->gp.mag();
00085       }else if((lastbx!=point->bx)){
00086         //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t \t one knee between"<<lastbx<<point->bx<<std::endl;
00087         knees++;
00088         if(knees==2){
00089           float maginsecondknee=point->gp.mag();
00090           betavalue1=maginknee/(25-delay+maginknee/30.)/30.;
00091           float dr =(maginknee-maginsecondknee);
00092           betavalue2 = dr/(25.+dr/30.);
00093           //std::cout<<"Inside BetaCalculatorRPC \t \t \t \t \t betavalue1="<<betavalue1<<" betavalue2="<<betavalue2<<std::endl;
00094         }
00095       }
00096       lastbx=point->bx;
00097     }
00098     betavalue = (betavalue1 + betavalue2)*0.5;
00099   }
00100   
00101   if(Candidate){
00102     //std::cout<<"Inside BetaCalculatorRPC \t \t \t yes! We found an HSCPs let's try to estimate beta"<<std::endl;
00103   }else{
00104     //std::cout<<"Inside BetaCalculatorRPC \t \t \t seems that there is no RPC HSCP Candidate in the set of RPC4DHit"<<std::endl;
00105     betavalue = 1.;
00106   }
00107   
00108   if(HSCPRPCRecHits.size()==0){
00109     //std::cout<<"Inside BetaCalculatorRPC \t WARINNG EMPTY RPC4DRecHits CONTAINER!!!"<<std::endl;
00110     betavalue = 1.;
00111   }
00112 }
00113 
00114 
00115 
00116 void BetaCalculatorRPC::addInfoToCandidate(HSCParticle& candidate, const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00117 
00118   edm::ESHandle<RPCGeometry> rpcGeo;
00119   iSetup.get<MuonGeometryRecord>().get(rpcGeo);
00120 
00121   edm::Handle<RPCRecHitCollection> rpcHits;
00122   iEvent.getByLabel(rpcRecHitsLabel,rpcHits);
00123 
00124 
00125   // here we do basically as in RPCHSCPCANDIDATE.cc, but just for the hits on the muon of interest
00126   RPCBetaMeasurement result;
00127   std::vector<RPCHit4D> hits;
00128   // so, loop on the RPC hits of the muon
00129   trackingRecHit_iterator start,stop;
00130   reco::Track track;
00131 
00132   if(      candidate.hasMuonRef() && candidate.muonRef()->combinedMuon()  .isNonnull()){ 
00133      start = candidate.muonRef()->combinedMuon()->recHitsBegin();
00134      stop  = candidate.muonRef()->combinedMuon()->recHitsEnd();    
00135   }else if(candidate.hasMuonRef() && candidate.muonRef()->standAloneMuon().isNonnull()){ track=*(candidate.muonRef()->standAloneMuon());
00136      start = candidate.muonRef()->standAloneMuon()->recHitsBegin();
00137      stop  = candidate.muonRef()->standAloneMuon()->recHitsEnd();  
00138   }else return;
00139 /*
00140   if(candidate.hasMuonCombinedTrack()) {
00141     start = candidate.combinedTrack().recHitsBegin();
00142     stop  = candidate.combinedTrack().recHitsEnd();
00143   } else if(candidate.hasMuonStaTrack()) {
00144     start = candidate.staTrack().recHitsBegin();
00145     stop  = candidate.staTrack().recHitsEnd();
00146   } else return; 
00147 */
00148   for(trackingRecHit_iterator recHit = start; recHit != stop; ++recHit) {
00149     if ( (*recHit)->geographicalId().subdetId() != MuonSubdetId::RPC ) continue;
00150     if ( (*recHit)->geographicalId().det() != DetId::Muon  ) continue;
00151     if (!(*recHit)->isValid()) continue; //Is Valid?
00152        
00153     RPCDetId rollId = (RPCDetId)(*recHit)->geographicalId();
00154 
00155     typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
00156     rangeRecHits recHitCollection =  rpcHits->get(rollId);
00157     RPCRecHitCollection::const_iterator recHitC;
00158     int size = 0;
00159     int clusterS=0;
00160     for(recHitC = recHitCollection.first; recHitC != recHitCollection.second ; recHitC++) {
00161       clusterS=(*recHitC).clusterSize();
00162 //      RPCDetId rollId = (RPCDetId)(*recHitC).geographicalId();
00163 //      std::cout<<"\t \t \t \t"<<rollId<<" bx "<<(*recHitC).BunchX()<<std::endl;
00164       size++;
00165     }
00166     if(size>1) continue; //Is the only RecHit in this roll.?                                                                                                         
00167     if(clusterS>4) continue; //Is the Cluster Size 5 or bigger?    
00168     
00169     LocalPoint recHitPos=(*recHit)->localPosition();
00170     const RPCRoll* rollasociated = rpcGeo->roll(rollId);
00171     const BoundPlane & RPCSurface = rollasociated->surface();
00172     
00173     RPCHit4D ThisHit;
00174     ThisHit.bx = ((RPCRecHit*)(&(**recHit)))->BunchX();
00175     ThisHit.gp = RPCSurface.toGlobal(recHitPos);
00176     ThisHit.id = (RPCDetId)(*recHit)->geographicalId().rawId();
00177     hits.push_back(ThisHit);
00178     
00179   }
00180   // here we go on with the RPC procedure 
00181   std::sort(hits.begin(), hits.end());
00182   int lastbx=-7;
00183   bool increasing = true;
00184   bool outOfTime = false;
00185   for(std::vector<RPCHit4D>::iterator point = hits.begin(); point < hits.end(); ++point) {
00186     outOfTime |= (point->bx!=0); //condition 1: at least one measurement must have BX!=0
00187     increasing &= (point->bx>=lastbx); //condition 2: BX must increase when going inside-out.
00188     lastbx = point->bx;
00189   }
00190   result.isCandidate = (outOfTime&&increasing);
00191  
00192   //result.beta = 1; // here we should get some pattern-based estimate
00193   algo(hits);
00194   result.beta = beta();
00195   candidate.setRpc(result);
00196 }
00197 
00198