CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/special/src/HLTRPCTrigNoSyncFilter.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Class:      HLTRPCTrigNoSyncFilter
00004 // 
00012 //
00013 // Original Author:  Camilo Andres Carrillo Montoya
00014 //         Created:  Thu Oct 29 11:04:22 CET 2009
00015 // $Id: HLTRPCTrigNoSyncFilter.cc,v 1.3 2010/02/23 14:43:03 carrillo Exp $
00016 //
00017 //
00018 
00019 #include "HLTrigger/special/interface/HLTRPCTrigNoSyncFilter.h"
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 //
00025 // constants, enums and typedefs
00026 //
00027 
00028 //
00029 // static data member definitions
00030 //
00031 
00032 //
00033 // constructors and destructor
00034 //
00035 
00036 typedef struct {
00037   int id;
00038   int bx;
00039   GlobalPoint gp;
00040 } RPC4DHit;
00041 
00042 bool bigmag(const RPC4DHit &Point1,const RPC4DHit &Point2){
00043   if((Point2).gp.mag() > (Point1).gp.mag()) return true;
00044   else return false;
00045 }
00046 
00047 HLTRPCTrigNoSyncFilter::HLTRPCTrigNoSyncFilter(const edm::ParameterSet& iConfig)
00048 {
00049   //now do what ever initialization is needed
00050   m_GMTInputTag =iConfig.getParameter< edm::InputTag >("GMTInputTag");
00051   rpcRecHitsLabel = iConfig.getParameter<edm::InputTag>("rpcRecHits");  
00052 }
00053 
00054 
00055 HLTRPCTrigNoSyncFilter::~HLTRPCTrigNoSyncFilter()
00056 {
00057     // do anything here that needs to be done at desctruction time
00058    // (e.g. close files, deallocate resources etc.)
00059 }
00060 
00061 
00062 //
00063 // member functions
00064 //
00065 
00066 // ------------ method called on each new Event  ------------
00067 bool HLTRPCTrigNoSyncFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00068 
00069   std::vector<RPC4DHit> GlobalRPC4DHits;
00070   std::vector<RPC4DHit> GlobalRPC4DHitsNoBx0;
00071 
00072   edm::Handle<RPCRecHitCollection> rpcRecHits;
00073   
00074   //std::cout <<"\t Getting the RPC RecHits"<<std::endl;
00075 
00076   iEvent.getByLabel(rpcRecHitsLabel,rpcRecHits);
00077   
00078   if(!rpcRecHits.isValid()){
00079     //std::cout<<"no valid RPC Collection"<<std::endl;
00080     //std::cout<<"event skipped"<<std::endl;
00081     return false;
00082   }
00083   if(rpcRecHits->begin() == rpcRecHits->end()){
00084     //std::cout<<"no RPC Hits in this event"<<std::endl;
00085     //std::cout<<"event skipped"<<std::endl;
00086     return false;
00087   }
00088   
00089   RPCRecHitCollection::const_iterator recHit;
00090 
00091   edm::ESHandle<RPCGeometry> rpcGeo;
00092   iSetup.get<MuonGeometryRecord>().get(rpcGeo);
00093 
00094   int k = 0;
00095 
00096   for (recHit = rpcRecHits->begin(); recHit != rpcRecHits->end(); recHit++){
00097     RPCDetId rollId = (RPCDetId)(*recHit).rpcId();
00098     LocalPoint recHitPos=recHit->localPosition();    
00099     const RPCRoll* rollasociated = rpcGeo->roll(rollId);
00100     const BoundPlane & RPCSurface = rollasociated->surface(); 
00101     GlobalPoint RecHitInGlobal = RPCSurface.toGlobal(recHitPos);
00102     
00103     int BX = recHit->BunchX();
00104     //std::cout<<"\t \t We have an RPC Rec Hit! bx="<<BX<<" Roll="<<rollId<<" Global Position="<<RecHitInGlobal<<std::endl;
00105     
00106     RPC4DHit ThisHit;
00107     ThisHit.bx =  BX;
00108     ThisHit.gp = RecHitInGlobal;
00109     ThisHit.id = k;
00110     GlobalRPC4DHits.push_back(ThisHit);
00111     if(BX!=0)GlobalRPC4DHitsNoBx0.push_back(ThisHit);
00112     k++;
00113   }
00114 
00115   if(GlobalRPC4DHitsNoBx0.size()==0){
00116     //std::cout<<"all RPC Hits are syncrhonized"<<std::endl;
00117     //std::cout<<"event skipped"<<std::endl;
00118     return false;
00119   }
00120 
00121   if(GlobalRPC4DHitsNoBx0.size()>100){
00122     //std::cout<<"too many rpcHits preventing HLT eternal loop"<<std::endl;
00123     //std::cout<<"event skipped"<<std::endl;
00124     return false;
00125   }
00126    
00127   edm::Handle<L1MuGMTReadoutCollection> gmtrc_handle;
00128   iEvent.getByLabel(m_GMTInputTag,gmtrc_handle);
00129   
00130   std::vector<L1MuGMTExtendedCand> gmt_candidates = (*gmtrc_handle).getRecord().getGMTCands();
00131   
00132   std::vector<L1MuGMTExtendedCand>::const_iterator candidate;
00133   //std::cout<<"The number of GMT Candidates in this event is "<<gmt_candidates.size()<<std::endl;
00134 
00135   if(gmt_candidates.size()==0){
00136     //std::cout<<"event skipped no gmt candidates"<<std::endl;
00137     return false;
00138   }
00139   
00140   for(candidate=gmt_candidates.begin(); candidate!=gmt_candidates.end(); candidate++) {
00141     int qual = candidate->quality();
00142     //std::cout<<"quality of this GMT Candidate (should be >5)= "<<qual<<std::endl;
00143     if(qual < 5) continue;
00144     float eta=candidate->etaValue();
00145     float phi=candidate->phiValue();      
00146     
00147     //std::cout<<" phi="<<phi<<" eta="<<eta<<std::endl;
00148     
00149     //Creating container in this etaphi direction;
00150     
00151     std::vector<RPC4DHit> PointsForGMT;
00152     
00153     for(std::vector<RPC4DHit>::iterator Point = GlobalRPC4DHitsNoBx0.begin(); Point!=GlobalRPC4DHitsNoBx0.end(); ++Point){ 
00154       float phiP = Point->gp.phi();
00155       float etaP = Point->gp.eta();
00156       //std::cout<<"\t \t GMT   phi="<<phi<<" eta="<<eta<<std::endl;
00157       //std::cout<<"\t \t Point phi="<<phiP<<" eta="<< etaP<<std::endl;
00158       //std::cout<<"\t \t "<<fabs(phi-phiP)<<" < 0,1? "<<fabs(eta-etaP)<<" < 0.20 ?"<<std::endl;
00159       if(fabs(phi-phiP) <=0.1 && fabs(eta-etaP)<=0.20){
00160         PointsForGMT.push_back(*Point);
00161         //std::cout<<"\t \t match!"<<std::endl;
00162       }
00163     }
00164     
00165     //std::cout<<"\t Points matching this GMT="<<PointsForGMT.size()<<std::endl;
00166 
00167     if(PointsForGMT.size()<1){
00168       //std::cout<<"\t Not enough RPCRecHits not syncrhonized for this GMT!!!"<<std::endl;
00169       continue;
00170     }
00171       
00172     std::sort(PointsForGMT.begin(), PointsForGMT.end(), bigmag);
00173 
00174     //std::cout<<"GMT candidate phi="<<phi<<std::endl;
00175     //std::cout<<"GMT candidate eta="<<eta<<std::endl;
00176 
00177     int lastbx=-7;
00178     bool outOfTime = false;
00179     bool incr = true;
00180     bool anydifferentzero = true;
00181     bool anydifferentone = true;
00182     
00183     //std::cout<<"\t \t loop on the RPCHit4D!!!"<<std::endl;
00184     for(std::vector<RPC4DHit>::iterator point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
00185       //float r=point->gp.mag();
00186       outOfTime |= (point->bx!=0); //condition 1: at least one measurement must have BX!=0
00187       incr &= (point->bx>=lastbx); //condition 2: BX must be increase when going inside-out.
00188       anydifferentzero &= (!point->bx==0); //to check one knee withoutzeros
00189       anydifferentone &= (!point->bx==1); //to check one knee withoutones
00190       lastbx = point->bx;
00191       //std::cout<<"\t \t  r="<<r<<" phi="<<point->gp.phi()<<" eta="<<point->gp.eta()<<" bx="<<point->bx<<" outOfTime"<<outOfTime<<" incr"<<incr<<" anydifferentzero"<<anydifferentzero<<std::endl;
00192     }
00193     //std::cout<<"\t \t";
00194     
00195     //for(std::vector<RPC4DHit>::iterator point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
00196       //std::cout<<point->bx;
00197     //}
00198     //std::cout<<std::endl;
00199     
00200     bool Candidate = (outOfTime&&incr);
00201     
00202     if(Candidate){ 
00203       //std::cout<<" Event Passed We found an strange trigger Candidate phi="<<phi<<" eta="<<eta<<std::endl;
00204       return true;
00205     }
00206   }
00207   
00208   //std::cout<<"event skipped rechits out of time but not related with a GMT "<<std::endl;
00209   return false;
00210 }
00211 
00212 // ------------ method called once each job just before starting event loop  ------------
00213 void 
00214 HLTRPCTrigNoSyncFilter::beginJob(){
00215 }
00216 
00217 // ------------ method called once each job just after ending the event loop  ------------
00218 void 
00219 HLTRPCTrigNoSyncFilter::endJob(){
00220 }
00221 
00222 //define this as a plug-in
00223 DEFINE_FWK_MODULE(HLTRPCTrigNoSyncFilter);