CMS 3D CMS Logo

CollectionMerger.cc
Go to the documentation of this file.
3 
4 
7 
10 
11 
15 
21 
22 
25 
28 
35 
39 
40 
41 
42 // -------- Here Tracker Merger -----------
43 template <typename T1, typename T2>
44 void CollectionMerger<T1,T2>::fill_output_obj_tracker(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections, bool print_pixel)
45 {
46 
47  std::map<uint32_t, std::vector<BaseHit> > output_map;
48  // First merge the collections with the help of the output map
49  for (auto inputCollection : inputCollections){
50  for ( typename MergeCollection::const_iterator clustSet = inputCollection->begin(); clustSet!= inputCollection->end(); ++clustSet ) {
51  DetId detIdObject( clustSet->detId() );
52  for (typename edmNew::DetSet<BaseHit>::const_iterator clustIt = clustSet->begin(); clustIt != clustSet->end(); ++clustIt ) {
53  output_map[detIdObject.rawId()].push_back(*clustIt);
54  }
55  }
56  }
57  // Now save it into the standard CMSSW format, with the standard Filler
58  for (typename std::map<uint32_t, std::vector<BaseHit> >::const_iterator outHits = output_map.begin(); outHits != output_map.end(); ++outHits ) {
59  DetId detIdObject(outHits->first);
60  typename MergeCollection::FastFiller spc(*output, detIdObject);
61  for (auto Hit : outHits->second){
62  spc.push_back(Hit);
63  }
64  }
65 
66 }
67 
68 
69 
70 
71 template <typename T1, typename T2>
72 void CollectionMerger<T1,T2>::fill_output_obj_calo(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
73 {
74  std::map<uint32_t, BaseHit > output_map;
75  // First merge the two collections again
76  for (auto inputCollection : inputCollections){
77  for ( typename MergeCollection::const_iterator recHit = inputCollection->begin(); recHit!= inputCollection->end(); ++recHit ) {
78  DetId detIdObject( recHit->detid().rawId() );
79  T2 *akt_calo_obj = &output_map[detIdObject.rawId()];
80  float new_energy = akt_calo_obj->energy() + recHit->energy();
81  T2 newRecHit(*recHit);
82  newRecHit.setEnergy(new_energy);
83  *akt_calo_obj = newRecHit;
84  }
85  // Now save it into the standard CMSSW format
86  for (typename std::map<uint32_t, BaseHit >::const_iterator outHits = output_map.begin(); outHits != output_map.end(); ++outHits ) {
87  output->push_back(outHits->second);
88  }
89  }
90  output->sort(); //Do a sort for this collection
91 }
92 
93 
94 
95 // -------- Here Muon Chamber Merger -----------
96 template <typename T1, typename T2>
97 void CollectionMerger<T1,T2>::fill_output_obj_muonchamber(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
98 {
99  std::map<uint32_t, std::vector<BaseHit> > output_map;
100  // First merge the collections with the help of the output map
101  for (auto inputCollection : inputCollections){
102  for ( typename MergeCollection::const_iterator recHit = inputCollection->begin(); recHit!= inputCollection->end(); ++recHit ) {
103  DetId detIdObject( recHit->geographicalId() );
104  output_map[detIdObject].push_back(*recHit);
105  }
106  }
107  // Now save it into the standard CMSSW format, with the standard Filler
108  for (typename std::map<uint32_t, std::vector<BaseHit> >::const_iterator outHits = output_map.begin(); outHits != output_map.end(); ++outHits ) {
109  output->put((typename T1::id_iterator::value_type) outHits->first , outHits->second.begin(), outHits->second.end()); // The DTLayerId misses the automatic type cast
110  }
111 }
112 
113 
114 
115 
116 // Here some overloaded functions, which are needed such that the right merger function is called for the indivudal Collections
117 template <typename T1, typename T2>
118 void CollectionMerger<T1,T2>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
119 {
120  assert(0); // CV: make sure general function never gets called;
121  // always use template specializations
122 }
123 
124 // Start with the Tracker collections
125 template <>
126 void CollectionMerger<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
127 {
128  fill_output_obj_tracker(output,inputCollections,true);
129 
130 }
131 
132 
133 template <>
134 void CollectionMerger<edmNew::DetSetVector<SiStripCluster>, SiStripCluster>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
135 {
136  fill_output_obj_tracker(output,inputCollections);
137 }
138 
139 
140 // Next are the Calo entries
141 template <>
142 void CollectionMerger<edm::SortedCollection<EcalRecHit>, EcalRecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
143 {
144  fill_output_obj_calo(output,inputCollections);
145 }
146 
147 
148 template <>
149 void CollectionMerger<edm::SortedCollection<HBHERecHit>, HBHERecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
150 {
151  fill_output_obj_calo(output,inputCollections);
152 }
153 
154 template <>
155 void CollectionMerger<edm::SortedCollection<HFRecHit>, HFRecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
156 {
157  fill_output_obj_calo(output,inputCollections);
158 }
159 
160 template <>
161 void CollectionMerger<edm::SortedCollection<HORecHit>, HORecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
162 {
163  fill_output_obj_calo(output,inputCollections);
164 }
165 
166 template <>
167 void CollectionMerger<edm::SortedCollection<CastorRecHit>, CastorRecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
168 {
169  fill_output_obj_calo(output,inputCollections);
170 }
171 
172 
173 template <>
174 void CollectionMerger<edm::SortedCollection<ZDCRecHit>, ZDCRecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
175 {
176  fill_output_obj_calo(output,inputCollections);
177 }
178 
179 
180 
181 // Here the Muon Chamber
182 template <>
183 void CollectionMerger<edm::RangeMap<DTLayerId, edm::OwnVector<DTRecHit1DPair> >, DTRecHit1DPair >::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
184 {
185  fill_output_obj_muonchamber(output,inputCollections);
186 }
187 
188 template <>
189 void CollectionMerger<edm::RangeMap<CSCDetId, edm::OwnVector<CSCRecHit2D> >, CSCRecHit2D >::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
190 {
191  fill_output_obj_muonchamber(output,inputCollections);
192 }
193 
194 
195 template <>
196 void CollectionMerger<edm::RangeMap<RPCDetId, edm::OwnVector<RPCRecHit> >, RPCRecHit>::fill_output_obj(std::unique_ptr<MergeCollection > & output, std::vector<edm::Handle<MergeCollection> > &inputCollections)
197 {
198  fill_output_obj_muonchamber(output,inputCollections);
199 }
200 
201 
202 
205 
212 
CollectionMerger< edm::SortedCollection< HORecHit >, HORecHit > HORecHitColMerger
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
CollectionMerger< edm::SortedCollection< CastorRecHit >, CastorRecHit > CastorRecHitColMerger
data_type const * const_iterator
Definition: DetSetNew.h:30
CollectionMerger< edmNew::DetSetVector< SiPixelCluster >, SiPixelCluster > PixelColMerger
Definition: DetId.h:18
void fill_output_obj_muonchamber(std::unique_ptr< MergeCollection > &output, std::vector< edm::Handle< MergeCollection > > &inputCollections)
void fill_output_obj_calo(std::unique_ptr< MergeCollection > &output, std::vector< edm::Handle< MergeCollection > > &inputCollections)
CollectionMerger< edm::SortedCollection< ZDCRecHit >, ZDCRecHit > ZDCRecHitColMerger
Pixel cluster – collection of neighboring pixels above threshold.
void fill_output_obj(std::unique_ptr< MergeCollection > &output, std::vector< edm::Handle< MergeCollection > > &inputCollections)
CollectionMerger< edmNew::DetSetVector< SiStripCluster >, SiStripCluster > StripColMerger
CollectionMerger< edm::RangeMap< DTLayerId, edm::OwnVector< DTRecHit1DPair > >, DTRecHit1DPair > DTRecHitColMerger
CollectionMerger< edm::RangeMap< CSCDetId, edm::OwnVector< CSCRecHit2D > >, CSCRecHit2D > CSCRecHitColMerger
CollectionMerger< edm::SortedCollection< HFRecHit >, HFRecHit > HFRecHitColMerger
CollectionMerger< edm::RangeMap< RPCDetId, edm::OwnVector< RPCRecHit > >, RPCRecHit > RPCRecHitColMerger
void fill_output_obj_tracker(std::unique_ptr< MergeCollection > &output, std::vector< edm::Handle< MergeCollection > > &inputCollections, bool print_pixel=false)
CollectionMerger< edm::SortedCollection< HBHERecHit >, HBHERecHit > HBHERecHitColMerger
CollectionMerger< edm::SortedCollection< EcalRecHit >, EcalRecHit > EcalRecHitColMerger