CMS 3D CMS Logo

HLTRPCTrigNoSyncFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Class: HLTRPCTrigNoSyncFilter
4 //
12 //
13 // Original Author: Camilo Andres Carrillo Montoya
14 // Created: Thu Oct 29 11:04:22 CET 2009
15 //
16 //
17 
18 #include "HLTRPCTrigNoSyncFilter.h"
19 
20 // system include files
21 #include <memory>
22 
23 //
24 // constants, enums and typedefs
25 //
26 
27 //
28 // static data member definitions
29 //
30 
31 //
32 // constructors and destructor
33 //
34 
35 typedef struct {
36  int id;
37  int bx;
39 } RPC4DHit;
40 
41 bool bigmag(const RPC4DHit &Point1,const RPC4DHit &Point2){
42  if((Point2).gp.mag() > (Point1).gp.mag()) return true;
43  else return false;
44 }
45 
47 {
48  //now do what ever initialization is needed
49  m_GMTInputTag =iConfig.getParameter< edm::InputTag >("GMTInputTag");
50  rpcRecHitsLabel = iConfig.getParameter<edm::InputTag>("rpcRecHits");
51  m_GMTInputToken = consumes<L1MuGMTReadoutCollection>(m_GMTInputTag);
52  rpcRecHitsToken = consumes<RPCRecHitCollection>(rpcRecHitsLabel);
53 
54 }
55 
56 
58 {
59  // do anything here that needs to be done at desctruction time
60  // (e.g. close files, deallocate resources etc.)
61 }
62 
63 void
67  desc.add<edm::InputTag>("GMTInputTag",edm::InputTag("hltGtDigis"));
68  desc.add<edm::InputTag>("rpcRecHits",edm::InputTag("hltRpcRecHits"));
69  descriptions.add("hltRPCTrigNoSyncFilter",desc);
70 }
71 
72 //
73 // member functions
74 //
75 
76 // ------------ method called on each new Event ------------
78 
79  std::vector<RPC4DHit> GlobalRPC4DHits;
80  std::vector<RPC4DHit> GlobalRPC4DHitsNoBx0;
81 
83 
84  //std::cout <<"\t Getting the RPC RecHits"<<std::endl;
85 
86  iEvent.getByToken(rpcRecHitsToken,rpcRecHits);
87 
88  if(!rpcRecHits.isValid()){
89  //std::cout<<"no valid RPC Collection"<<std::endl;
90  //std::cout<<"event skipped"<<std::endl;
91  return false;
92  }
93  if(rpcRecHits->begin() == rpcRecHits->end()){
94  //std::cout<<"no RPC Hits in this event"<<std::endl;
95  //std::cout<<"event skipped"<<std::endl;
96  return false;
97  }
98 
100 
102  iSetup.get<MuonGeometryRecord>().get(rpcGeo);
103 
104  int k = 0;
105 
106  for (recHit = rpcRecHits->begin(); recHit != rpcRecHits->end(); recHit++){
107  RPCDetId rollId = (RPCDetId)(*recHit).rpcId();
108  LocalPoint recHitPos=recHit->localPosition();
109  const RPCRoll* rollasociated = rpcGeo->roll(rollId);
110  const BoundPlane & RPCSurface = rollasociated->surface();
111  GlobalPoint RecHitInGlobal = RPCSurface.toGlobal(recHitPos);
112 
113  int BX = recHit->BunchX();
114  //std::cout<<"\t \t We have an RPC Rec Hit! bx="<<BX<<" Roll="<<rollId<<" Global Position="<<RecHitInGlobal<<std::endl;
115 
116  RPC4DHit ThisHit;
117  ThisHit.bx = BX;
118  ThisHit.gp = RecHitInGlobal;
119  ThisHit.id = k;
120  GlobalRPC4DHits.push_back(ThisHit);
121  if(BX!=0)GlobalRPC4DHitsNoBx0.push_back(ThisHit);
122  k++;
123  }
124 
125  if(GlobalRPC4DHitsNoBx0.empty()){
126  //std::cout<<"all RPC Hits are syncrhonized"<<std::endl;
127  //std::cout<<"event skipped"<<std::endl;
128  return false;
129  }
130 
131  if(GlobalRPC4DHitsNoBx0.size()>100){
132  //std::cout<<"too many rpcHits preventing HLT eternal loop"<<std::endl;
133  //std::cout<<"event skipped"<<std::endl;
134  return false;
135  }
136 
138  iEvent.getByToken(m_GMTInputToken,gmtrc_handle);
139 
140  std::vector<L1MuGMTExtendedCand> gmt_candidates = (*gmtrc_handle).getRecord().getGMTCands();
141 
142  std::vector<L1MuGMTExtendedCand>::const_iterator candidate;
143  //std::cout<<"The number of GMT Candidates in this event is "<<gmt_candidates.size()<<std::endl;
144 
145  if(gmt_candidates.empty()){
146  //std::cout<<"event skipped no gmt candidates"<<std::endl;
147  return false;
148  }
149 
150  for(candidate=gmt_candidates.begin(); candidate!=gmt_candidates.end(); candidate++) {
151  int qual = candidate->quality();
152  //std::cout<<"quality of this GMT Candidate (should be >5)= "<<qual<<std::endl;
153  if(qual < 5) continue;
154  float eta=candidate->etaValue();
155  float phi=candidate->phiValue();
156 
157  //std::cout<<" phi="<<phi<<" eta="<<eta<<std::endl;
158 
159  //Creating container in this etaphi direction;
160 
161  std::vector<RPC4DHit> PointsForGMT;
162 
163  for(auto & Point : GlobalRPC4DHitsNoBx0){
164  float phiP = Point.gp.phi();
165  float etaP = Point.gp.eta();
166  //std::cout<<"\t \t GMT phi="<<phi<<" eta="<<eta<<std::endl;
167  //std::cout<<"\t \t Point phi="<<phiP<<" eta="<< etaP<<std::endl;
168  //std::cout<<"\t \t "<<fabs(phi-phiP)<<" < 0,1? "<<fabs(eta-etaP)<<" < 0.20 ?"<<std::endl;
169  if(fabs(phi-phiP) <=0.1 && fabs(eta-etaP)<=0.20){
170  PointsForGMT.push_back(Point);
171  //std::cout<<"\t \t match!"<<std::endl;
172  }
173  }
174 
175  //std::cout<<"\t Points matching this GMT="<<PointsForGMT.size()<<std::endl;
176 
177  if(PointsForGMT.empty()){
178  //std::cout<<"\t Not enough RPCRecHits not syncrhonized for this GMT!!!"<<std::endl;
179  continue;
180  }
181 
182  std::sort(PointsForGMT.begin(), PointsForGMT.end(), bigmag);
183 
184  //std::cout<<"GMT candidate phi="<<phi<<std::endl;
185  //std::cout<<"GMT candidate eta="<<eta<<std::endl;
186 
187  int lastbx=-7;
188  bool outOfTime = false;
189  bool incr = true;
190 
191  //std::cout<<"\t \t loop on the RPCHit4D!!!"<<std::endl;
192  for(auto point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
193  //float r=point->gp.mag();
194  outOfTime |= (point->bx!=0); //condition 1: at least one measurement must have BX!=0
195  incr &= (point->bx>=lastbx); //condition 2: BX must be increase when going inside-out.
196  lastbx = point->bx;
197  }
198  //std::cout<<"\t \t";
199 
200  //for(std::vector<RPC4DHit>::iterator point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
201  //std::cout<<point->bx;
202  //}
203  //std::cout<<std::endl;
204 
205  bool Candidate = (outOfTime&&incr);
206 
207  if(Candidate){
208  //std::cout<<" Event Passed We found an strange trigger Candidate phi="<<phi<<" eta="<<eta<<std::endl;
209  return true;
210  }
211  }
212 
213  //std::cout<<"event skipped rechits out of time but not related with a GMT "<<std::endl;
214  return false;
215 }
216 
217 // ------------ method called once each job just before starting event loop ------------
218 void
220 }
221 
222 // ------------ method called once each job just after ending the event loop ------------
223 void
225 }
226 
227 // declare this class as a framework plugin
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HLTRPCTrigNoSyncFilter(const edm::ParameterSet &)
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsToken
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
edm::EDGetTokenT< L1MuGMTReadoutCollection > m_GMTInputToken
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool bigmag(const RPC4DHit &Point1, const RPC4DHit &Point2)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
int k[5][pyjets_maxn]
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
T get() const
Definition: EventSetup.h:71
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:75