CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCRecHitBaseAlgo.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2008/12/04 16:13:16 $
5  * $Revision: 1.8 $
6  * \author M. Maggi -- INFN Bari
7  */
8 
9 
10 
16 
21 
22 
24  // theSync = RPCTTrigSyncFactory::get()->create(config.getParameter<string>("tTrigMode"),
25  //config.getParameter<ParameterSet>("tTrigModeConfig"));
26 }
27 
29 
30 
31 // Build all hits in the range associated to the layerId, at the 1st step.
33  const RPCDetId& rpcId,
34  const RPCDigiCollection::Range& digiRange,
35  const RollMask& mask) {
37 
38 
39  RPCClusterizer clizer;
40  RPCClusterContainer tcls = clizer.doAction(digiRange);
41  RPCMaskReClusterizer mrclizer;
42  RPCClusterContainer cls = mrclizer.doAction(rpcId,tcls,mask);
43 
44 
45  for (RPCClusterContainer::const_iterator cl = cls.begin();
46  cl != cls.end(); cl++){
47 
48  LocalError tmpErr;
50  // Call the compute method
51  bool OK = this->compute(roll, *cl, point, tmpErr);
52  if (!OK) continue;
53 
54  // Build a new pair of 1D rechit
55  int firstClustStrip= cl->firstStrip();
56  int clusterSize=cl->clusterSize();
57  RPCRecHit* recHit = new RPCRecHit(rpcId,cl->bx(),firstClustStrip,clusterSize,point,tmpErr);
58 
59 
60  result.push_back(recHit);
61  }
62  return result;
63 }
pair< int, edm::FunctionWithDict > OK
Definition: findMethod.cc:70
RPCClusterContainer doAction(const RPCDigiCollection::Range &digiRange)
virtual ~RPCRecHitBaseAlgo()
Destructor.
void push_back(D *&d)
Definition: OwnVector.h:273
tuple result
Definition: query.py:137
std::bitset< maskSIZE > RollMask
Definition: RPCRollMask.h:8
RPCRecHitBaseAlgo(const edm::ParameterSet &config)
Constructor.
virtual bool compute(const RPCRoll &roll, const RPCCluster &cl, LocalPoint &Point, LocalError &error) const =0
standard local recHit computation
std::pair< const_iterator, const_iterator > Range
RPCClusterContainer doAction(const RPCDetId &, RPCClusterContainer &, const RollMask &)
virtual edm::OwnVector< RPCRecHit > reconstruct(const RPCRoll &roll, const RPCDetId &rpcId, const RPCDigiCollection::Range &digiRange, const RollMask &mask)
Build all hits in the range associated to the rpcId, at the 1st step.
std::set< RPCCluster > RPCClusterContainer
*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