CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ShallowSimhitClustersProducer.cc
Go to the documentation of this file.
2 
6 
18 #include "boost/foreach.hpp"
19 
21  : inputTags( iConfig.getParameter<std::vector<edm::InputTag> >("InputTags") ),
22  theClustersLabel( iConfig.getParameter<edm::InputTag>("Clusters")),
23  Prefix( iConfig.getParameter<std::string>("Prefix") )
24 {
25  produces <std::vector<unsigned> > ( Prefix + "hits" );
26  produces <std::vector<float> > ( Prefix + "strip" );
27  produces <std::vector<float> > ( Prefix + "localtheta" );
28  produces <std::vector<float> > ( Prefix + "localphi" );
29  produces <std::vector<float> > ( Prefix + "localx" );
30  produces <std::vector<float> > ( Prefix + "localy" );
31  produces <std::vector<float> > ( Prefix + "localz" );
32  produces <std::vector<float> > ( Prefix + "momentum" );
33  produces <std::vector<float> > ( Prefix + "energyloss" );
34  produces <std::vector<float> > ( Prefix + "time" );
35  produces <std::vector<int> > ( Prefix + "particle" );
36  produces <std::vector<unsigned short> > ( Prefix + "process" );
37 }
38 
42 
43  int size = clustermap.size();
44  std::auto_ptr<std::vector<unsigned> > hits ( new std::vector<unsigned> (size, 0) );
45  std::auto_ptr<std::vector<float> > strip ( new std::vector<float> (size, -100) );
46  std::auto_ptr<std::vector<float> > localtheta ( new std::vector<float> (size, -100) );
47  std::auto_ptr<std::vector<float> > localphi ( new std::vector<float> (size, -100) );
48  std::auto_ptr<std::vector<float> > localx ( new std::vector<float> (size, -100) );
49  std::auto_ptr<std::vector<float> > localy ( new std::vector<float> (size, -100) );
50  std::auto_ptr<std::vector<float> > localz ( new std::vector<float> (size, -100) );
51  std::auto_ptr<std::vector<float> > momentum ( new std::vector<float> (size, 0) );
52  std::auto_ptr<std::vector<float> > energyloss ( new std::vector<float> (size, -1) );
53  std::auto_ptr<std::vector<float> > time ( new std::vector<float> (size, -1) );
54  std::auto_ptr<std::vector<int> > particle ( new std::vector<int> (size,-500) );
55  std::auto_ptr<std::vector<unsigned short> > process ( new std::vector<unsigned short> (size,0) );
56 
57  edm::ESHandle<TrackerGeometry> theTrackerGeometry; iSetup.get<TrackerDigiGeometryRecord>().get( theTrackerGeometry );
60  edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters; iEvent.getByLabel("siStripClusters", "", clusters);
61 
62  BOOST_FOREACH( const edm::InputTag inputTag, inputTags ) { edm::Handle<std::vector<PSimHit> > simhits; iEvent.getByLabel(inputTag, simhits);
63  BOOST_FOREACH( const PSimHit hit, *simhits ) {
64 
65  const uint32_t id = hit.detUnitId();
66  const StripGeomDetUnit* theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTrackerGeometry->idToDet( id ) );
67  const LocalVector drift = shallow::drift(theStripDet, *magfield, *SiStripLorentzAngle);
68 
69  const float driftedstrip_ = theStripDet->specificTopology().strip( hit.localPosition()+0.5*drift );
70  const float hitstrip_ = theStripDet->specificTopology().strip( hit.localPosition() );
71 
72  shallow::CLUSTERMAP::const_iterator cluster = match_cluster( id, driftedstrip_, clustermap, *clusters);
73  if(cluster != clustermap.end()) {
74  unsigned i = cluster->second;
75  hits->at(i)+=1;
76  if(hits->at(i) == 1) {
77  strip->at(i) = hitstrip_;
78  localtheta->at(i) = hit.thetaAtEntry();
79  localphi->at(i) = hit.phiAtEntry();
80  localx->at(i) = hit.localPosition().x();
81  localy->at(i) = hit.localPosition().y();
82  localz->at(i) = hit.localPosition().z();
83  momentum->at(i) = hit.pabs();
84  energyloss->at(i) = hit.energyLoss();
85  time->at(i) = hit.timeOfFlight();
86  particle->at(i) = hit.particleType();
87  process->at(i) = hit.processType();
88  }
89  }
90  }
91  }
92 
93  iEvent.put( hits, Prefix + "hits" );
94  iEvent.put( strip, Prefix + "strip" );
95  iEvent.put( localtheta, Prefix + "localtheta" );
96  iEvent.put( localphi, Prefix + "localphi" );
97  iEvent.put( localx, Prefix + "localx" );
98  iEvent.put( localy, Prefix + "localy" );
99  iEvent.put( localz, Prefix + "localz" );
100  iEvent.put( momentum, Prefix + "momentum" );
101  iEvent.put( energyloss, Prefix + "energyloss" );
102  iEvent.put( time, Prefix + "time" );
103  iEvent.put( particle, Prefix + "particle" );
104  iEvent.put( process, Prefix + "process" );
105 }
106 
107 shallow::CLUSTERMAP::const_iterator ShallowSimhitClustersProducer::
108 match_cluster( const unsigned& id, const float& strip_, const shallow::CLUSTERMAP& clustermap, const edmNew::DetSetVector<SiStripCluster>& clusters) const {
109  shallow::CLUSTERMAP::const_iterator cluster = clustermap.end();
110  edmNew::DetSetVector<SiStripCluster>::const_iterator clustersDetSet = clusters.find(id);
111  if( clustersDetSet != clusters.end() ) {
112  edmNew::DetSet<SiStripCluster>::const_iterator left, right=clustersDetSet->begin();
113  while( right != clustersDetSet->end() && strip_ > right->barycenter() )
114  right++;
115  left = right-1;
116  if(right!=clustersDetSet->end() && right!=clustersDetSet->begin()) {
117  unsigned firstStrip = (right->barycenter()-strip_) < (strip_-left->barycenter()) ? right->firstStrip() : left->firstStrip();
118  cluster = clustermap.find( std::make_pair( id, firstStrip));
119  }
120  else if(right != clustersDetSet->begin())
121  cluster = clustermap.find( std::make_pair( id, left->firstStrip()));
122  else
123  cluster = clustermap.find( std::make_pair( id, right->firstStrip()));
124  }
125  return cluster;
126 }
127 
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
tuple SiStripLorentzAngle
Definition: redigi_cff.py:15
shallow::CLUSTERMAP::const_iterator match_cluster(const unsigned &, const float &, const shallow::CLUSTERMAP &, const edmNew::DetSetVector< SiStripCluster > &) const
Geom::Theta< float > thetaAtEntry() const
fast and more accurate access to momentumAtEntry().theta()
Definition: PSimHit.h:57
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:39
tuple magfield
Definition: HLT_ES_cff.py:2311
T y() const
Definition: PV3DBase.h:63
const_iterator find(id_type i, bool update=true) const
data_type const * const_iterator
Definition: DetSetNew.h:30
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
virtual float strip(const LocalPoint &) const =0
std::vector< edm::InputTag > inputTags
int iEvent
Definition: GenABIO.cc:243
float timeOfFlight() const
Definition: PSimHit.h:69
Local3DPoint localPosition() const
Definition: PSimHit.h:44
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
ShallowSimhitClustersProducer(const edm::ParameterSet &)
T z() const
Definition: PV3DBase.h:64
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
CLUSTERMAP make_cluster_map(const edm::Event &, edm::InputTag &)
Definition: ShallowTools.cc:15
const T & get() const
Definition: EventSetup.h:55
unsigned short processType() const
Definition: PSimHit.h:118
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
const_iterator end(bool update=true) const
int particleType() const
Definition: PSimHit.h:85
void produce(edm::Event &, const edm::EventSetup &)
std::map< std::pair< uint32_t, uint16_t >, unsigned int > CLUSTERMAP
Definition: ShallowTools.h:16
Geom::Phi< float > phiAtEntry() const
fast and more accurate access to momentumAtEntry().phi()
Definition: PSimHit.h:60
const_iterator begin(bool update=true) const
tuple process
Definition: LaserDQM_cfg.py:3
T x() const
Definition: PV3DBase.h:62
tuple size
Write out results.
unsigned int detUnitId() const
Definition: PSimHit.h:93