CMS 3D CMS Logo

RPCRecHitBaseAlgo.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author M. Maggi -- INFN Bari
5  */
6 
7 #include "RPCRecHitBaseAlgo.h"
8 #include "RPCClusterContainer.h"
9 #include "RPCCluster.h"
10 #include "RPCClusterizer.h"
11 #include "RPCMaskReClusterizer.h"
12 
14  // theSync = RPCTTrigSyncFactory::get()->create(config.getParameter<string>("tTrigMode"),
15  //config.getParameter<ParameterSet>("tTrigModeConfig"));
16 }
17 
18 // Build all hits in the range associated to the layerId, at the 1st step.
20  const RPCDetId& rpcId,
21  const RPCDigiCollection::Range& digiRange,
22  const RollMask& mask) {
24 
25  RPCClusterizer clizer;
26  RPCClusterContainer tcls = clizer.doAction(digiRange);
27  RPCMaskReClusterizer mrclizer;
28  RPCClusterContainer cls = mrclizer.doAction(rpcId, tcls, mask);
29 
30  for (const auto& cl : cls) {
31  LocalError tmpErr;
33  float time = 0, timeErr = -1;
34 
35  // Call the compute method
36  const bool OK = this->compute(roll, cl, point, tmpErr, time, timeErr);
37  if (!OK)
38  continue;
39 
40  // Build a new pair of 1D rechit
41  const int firstClustStrip = cl.firstStrip();
42  const int clusterSize = cl.clusterSize();
43  RPCRecHit* recHit = new RPCRecHit(rpcId, cl.bx(), firstClustStrip, clusterSize, point, tmpErr);
44  recHit->setTimeAndError(time, timeErr);
45 
46  result.push_back(recHit);
47  }
48 
49  return result;
50 }
Definition: config.py:1
RPCClusterContainer doAction(const RPCDigiCollection::Range &digiRange)
RPCClusterContainer doAction(const RPCDetId &id, RPCClusterContainer &initClusters, const RollMask &mask) const
virtual bool compute(const RPCRoll &roll, const RPCCluster &cl, LocalPoint &Point, LocalError &error, float &time, float &timeErr) const =0
standard local recHit computation
std::bitset< maskSIZE > RollMask
Definition: RPCRollMask.h:7
RPCRecHitBaseAlgo(const edm::ParameterSet &config)
Constructor.
std::pair< int, edm::FunctionWithDict > OK
Definition: findMethod.cc:126
std::pair< const_iterator, const_iterator > Range
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