00001
00002 #include "FWCore/Framework/interface/Frameworkfwd.h"
00003 #include "FWCore/Framework/interface/EDProducer.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 #include "FWCore/Utilities/interface/InputTag.h"
00007
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00011 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00014 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00015 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00016
00017 #include "DataFormats/TrackReco/interface/Track.h"
00018 #include "DataFormats/VertexReco/interface/Vertex.h"
00019
00020 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00021 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00022 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00023 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00025 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00026 #include "MagneticField/Engine/interface/MagneticField.h"
00027 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00028 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00029 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00030
00031
00032 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
00033 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
00034 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00035
00036
00037 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplate.h"
00038 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplate2D.h"
00039 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplateSplit.h"
00040 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplateDefs.h"
00041 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelTemplateReco.h"
00042
00043 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripTemplate.h"
00044 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripTemplateSplit.h"
00045 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripTemplateDefs.h"
00046 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripTemplateReco.h"
00047
00048 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
00049 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00050 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00051 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00052
00053
00054 #include "SimDataFormats/TrackerDigiSimLink/interface/StripDigiSimLink.h"
00055 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00056
00057 #include <boost/range.hpp>
00058 #include <boost/foreach.hpp>
00059 #include "boost/multi_array.hpp"
00060
00061 #include <iostream>
00062 using namespace std;
00063
00064 class TrackClusterSplitter : public edm::EDProducer
00065 {
00066
00067 public:
00068 TrackClusterSplitter(const edm::ParameterSet& iConfig) ;
00069 ~TrackClusterSplitter() ;
00070 void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
00071
00072 private:
00073 edm::InputTag stripClusters_;
00074 edm::InputTag pixelClusters_;
00075
00076 bool simSplitPixel_;
00077 bool simSplitStrip_;
00078 bool tmpSplitPixel_;
00079 bool tmpSplitStrip_;
00080
00081
00082
00083
00084
00085
00086
00087 bool useStraightTracks_;
00088
00089
00090
00091
00092
00093 bool useTrajectories_;
00094
00095
00096 edm::InputTag trajectories_;
00097
00098
00099 edm::InputTag vertices_;
00100
00101
00102 std::string propagatorName_;
00103 edm::ESHandle<MagneticField> magfield_;
00104 edm::ESHandle<Propagator> propagator_;
00105 edm::ESHandle<GlobalTrackingGeometry> geometry_;
00106
00107
00108 edm::Handle< edm::DetSetVector<PixelDigiSimLink> > pixeldigisimlink;
00109
00110
00111 edm::Handle< edm::DetSetVector<StripDigiSimLink> > stripdigisimlink;
00112
00113
00114
00115
00116 mutable SiPixelTemplate templ_ ;
00117 mutable SiPixelTemplate2D templ2D_;
00118
00119 mutable SiStripTemplate strip_templ_ ;
00120
00121
00122 struct TrackAndState
00123 {
00124 TrackAndState(const reco::Track *aTrack, TrajectoryStateOnSurface aState) :
00125 track(aTrack), state(aState) {}
00126 const reco::Track* track;
00127 TrajectoryStateOnSurface state;
00128 };
00129
00130
00131 template<typename Cluster>
00132 struct ClusterWithTracks
00133 {
00134 ClusterWithTracks(const Cluster &c) : cluster(&c) {}
00135 const Cluster* cluster;
00136 std::vector<TrackAndState> tracks;
00137 };
00138
00139 typedef ClusterWithTracks<SiPixelCluster> SiPixelClusterWithTracks;
00140 typedef ClusterWithTracks<SiStripCluster> SiStripClusterWithTracks;
00141
00142
00143
00144 typedef boost::sub_range<std::vector<SiPixelClusterWithTracks> > SiPixelClustersWithTracks;
00145 typedef boost::sub_range<std::vector<SiStripClusterWithTracks> > SiStripClustersWithTracks;
00146
00147
00148
00149 typedef std::pair<uint32_t, EncodedEventId> SimHitIdpr;
00150 TrackerHitAssociator* hitAssociator;
00151
00152 template<typename C>
00153 static const C* getCluster(const TrackingRecHit* hit) ;
00154
00155 template<typename C>
00156 static const C* equalClusters(const C &c1, const C &c2)
00157 {
00158 return false;
00159 }
00160
00161
00162 template<typename Cluster> class FindCluster
00163 {
00164
00165 public:
00166
00167 FindCluster(const TrackingRecHit* hit) : toFind_( getCluster<Cluster>(hit) ) { }
00168
00169 bool operator()(const ClusterWithTracks<Cluster> &test) const
00170 {
00171 assert(test.cluster);
00172 return test.cluster == toFind_ || equalClusters<Cluster>(*test.cluster, *toFind_);
00173 }
00174
00175 private:
00176
00177 const Cluster* toFind_;
00178
00179 };
00180
00181
00182 template<typename Cluster>
00183 void markClusters(std::map<uint32_t, boost::sub_range<std::vector<ClusterWithTracks<Cluster> > > >& map,
00184 const TrackingRecHit* hit,
00185 const reco::Track* track,
00186 const TrajectoryStateOnSurface& tsos) const ;
00187
00188 template<typename Cluster>
00189 std::auto_ptr<edmNew::DetSetVector<Cluster> >
00190 splitClusters(const std::map<uint32_t,
00191 boost::sub_range<std::vector<ClusterWithTracks<Cluster> > > > &input,
00192 const reco::Vertex &vtx) const ;
00193
00194 template<typename Cluster>
00195 void splitCluster(const ClusterWithTracks<Cluster> &cluster,
00196 const GlobalVector &dir,
00197 typename edmNew::DetSetVector<Cluster>::FastFiller &output,
00198 DetId detId) const ;
00199
00201 std::vector<SiPixelClusterWithTracks> allSiPixelClusters;
00202 std::map<uint32_t, SiPixelClustersWithTracks> siPixelDetsWithClusters;
00203
00204 std::vector<SiStripClusterWithTracks> allSiStripClusters;
00205 std::map<uint32_t, SiStripClustersWithTracks> siStripDetsWithClusters;
00206
00207
00208 };
00209
00210 template<> const SiPixelCluster * TrackClusterSplitter::getCluster<SiPixelCluster>(const TrackingRecHit *hit) ;
00211
00212 template<>
00213 void TrackClusterSplitter::splitCluster<SiPixelCluster> (const SiPixelClusterWithTracks &cluster,
00214 const GlobalVector &dir,
00215 edmNew::DetSetVector<SiPixelCluster>::FastFiller &output,
00216 DetId detId
00217 ) const ;
00218
00219 template<> const SiStripCluster * TrackClusterSplitter::getCluster<SiStripCluster>(const TrackingRecHit *hit) ;
00220
00221 template<>
00222 void TrackClusterSplitter::splitCluster<SiStripCluster> (const SiStripClusterWithTracks &cluster,
00223 const GlobalVector &dir,
00224 edmNew::DetSetVector<SiStripCluster>::FastFiller &output,
00225 DetId detId
00226 ) const ;
00227
00228
00229 #define foreach BOOST_FOREACH
00230
00231 TrackClusterSplitter::TrackClusterSplitter(const edm::ParameterSet& iConfig):
00232 stripClusters_(iConfig.getParameter<edm::InputTag>("stripClusters")),
00233 pixelClusters_(iConfig.getParameter<edm::InputTag>("pixelClusters")),
00234 useTrajectories_(iConfig.getParameter<bool>("useTrajectories")),
00235 trajectories_(iConfig.getParameter<edm::InputTag>(useTrajectories_ ? "trajTrackAssociations" : "tracks")),
00236 vertices_(iConfig.getParameter<edm::InputTag>("vertices"))
00237 {
00238 if ( !useTrajectories_ )
00239 propagatorName_ = iConfig.getParameter<std::string>("propagator");
00240
00241 produces< edmNew::DetSetVector<SiPixelCluster> >();
00242
00243 produces< edmNew::DetSetVector<SiStripCluster> >();
00244
00245 simSplitPixel_ = (iConfig.getParameter<bool>("simSplitPixel"));
00246 simSplitStrip_ = (iConfig.getParameter<bool>("simSplitStrip"));
00247 tmpSplitPixel_ = (iConfig.getParameter<bool>("tmpSplitPixel"));
00248 tmpSplitStrip_ = (iConfig.getParameter<bool>("tmpSplitStrip"));
00249
00250 useStraightTracks_ = (iConfig.getParameter<bool>("useStraightTracks"));
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270 templ_.pushfile( 40 );
00271 templ_.pushfile( 41 );
00272 templ2D_.pushfile( 40 );
00273 templ2D_.pushfile( 41 );
00274
00275
00276 strip_templ_.pushfile( 11 );
00277 strip_templ_.pushfile( 12 );
00278 strip_templ_.pushfile( 13 );
00279 strip_templ_.pushfile( 14 );
00280 strip_templ_.pushfile( 15 );
00281 strip_templ_.pushfile( 16 );
00282
00283 }
00284
00285
00286 template<>
00287 const SiStripCluster*
00288 TrackClusterSplitter::getCluster<SiStripCluster>(const TrackingRecHit* hit)
00289 {
00290 if ( typeid(*hit) == typeid(SiStripRecHit2D) )
00291 {
00292 return (static_cast<const SiStripRecHit2D &>(*hit)).cluster().get();
00293 }
00294 else if ( typeid(*hit) == typeid(SiStripRecHit1D) )
00295 {
00296 return (static_cast<const SiStripRecHit1D &>(*hit)).cluster().get();
00297 }
00298 else
00299 throw cms::Exception("Unsupported") << "Detid of type " << typeid(*hit).name() << " not supported.\n";
00300 }
00301
00302 template<>
00303 const SiPixelCluster*
00304 TrackClusterSplitter::getCluster<SiPixelCluster>(const TrackingRecHit* hit)
00305 {
00306 if ( typeid(*hit) == typeid(SiPixelRecHit) )
00307 {
00308 return (static_cast<const SiPixelRecHit&>(*hit)).cluster().get();
00309 }
00310 else
00311 throw cms::Exception("Unsupported") << "Detid of type " << typeid(*hit).name() << " not supported.\n";
00312 }
00313
00314 TrackClusterSplitter::~TrackClusterSplitter()
00315 {
00316 }
00317
00318 void
00319 TrackClusterSplitter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00320 {
00321 using namespace edm;
00322
00323 iSetup.get<GlobalTrackingGeometryRecord>().get(geometry_);
00324
00325 if ( !useTrajectories_ )
00326 {
00327 iSetup.get<IdealMagneticFieldRecord>().get( magfield_ );
00328 iSetup.get<TrackingComponentsRecord>().get( "AnalyticalPropagator", propagator_ );
00329 }
00330
00331 Handle<edmNew::DetSetVector<SiPixelCluster> > inputPixelClusters;
00332 Handle<edmNew::DetSetVector<SiStripCluster> > inputStripClusters;
00333
00334 iEvent.getByLabel(pixelClusters_, inputPixelClusters);
00335
00336 iEvent.getByLabel(stripClusters_, inputStripClusters);
00337
00338 if(simSplitStrip_)
00339 hitAssociator = new TrackerHitAssociator(iEvent);
00340
00341
00342 allSiPixelClusters.clear(); siPixelDetsWithClusters.clear();
00343 allSiStripClusters.clear(); siStripDetsWithClusters.clear();
00344
00345
00346 allSiPixelClusters.reserve(inputPixelClusters->dataSize());
00347 allSiStripClusters.reserve(inputStripClusters->dataSize());
00348
00349
00350
00351 foreach(const edmNew::DetSet<SiPixelCluster> &ds, *inputPixelClusters)
00352 {
00353 std::vector<SiPixelClusterWithTracks>::iterator start = allSiPixelClusters.end();
00354 allSiPixelClusters.insert(start, ds.begin(), ds.end());
00355
00356 std::vector<SiPixelClusterWithTracks>::iterator end = allSiPixelClusters.end();
00357 siPixelDetsWithClusters[ds.detId()] = SiPixelClustersWithTracks(start,end);
00358 }
00359
00360 foreach(const edmNew::DetSet<SiStripCluster> &ds, *inputStripClusters)
00361 {
00362 std::vector<SiStripClusterWithTracks>::iterator start = allSiStripClusters.end();
00363 allSiStripClusters.insert(start, ds.begin(), ds.end());
00364
00365 std::vector<SiStripClusterWithTracks>::iterator end = allSiStripClusters.end();
00366 siStripDetsWithClusters[ds.detId()] = SiStripClustersWithTracks(start,end);
00367 }
00368
00369 if ( useTrajectories_ )
00370 {
00371
00372
00373 Handle<TrajTrackAssociationCollection> trajectories;
00374 iEvent.getByLabel(trajectories_, trajectories);
00375 for ( TrajTrackAssociationCollection::const_iterator it = trajectories->begin(),
00376 ed = trajectories->end(); it != ed; ++it )
00377 {
00378 const Trajectory & traj = *it->key;
00379 const reco::Track * tk = &*it->val;
00380
00381 if ( traj.measurements().size() != tk->recHitsSize() )
00382 throw cms::Exception("Aargh") << "Sizes don't match: traj " << traj.measurements().size()
00383 << ", tk " << tk->recHitsSize() << "\n";
00384
00385 trackingRecHit_iterator it_hit = tk->recHitsBegin(), ed_hit = tk->recHitsEnd();
00386
00387 const Trajectory::DataContainer & tms = traj.measurements();
00388
00389 size_t i_hit = 0, last_hit = tms.size()-1;
00390
00391 bool first = true, reversed = false;
00392
00393 for (; it_hit != ed_hit; ++it_hit, ++i_hit)
00394 {
00395
00396
00397 if ((*it_hit)->geographicalId().rawId() == 0)
00398 {
00399
00400 continue;
00401 }
00402
00403
00404 if (first)
00405 {
00406
00407 if ((*it_hit)->geographicalId() == tms[i_hit].recHit()->hit()->geographicalId())
00408 {
00409 reversed = false;
00410 }
00411 else if ((*it_hit)->geographicalId() == tms[last_hit-i_hit].recHit()->hit()->geographicalId())
00412 {
00413 reversed = true;
00414 }
00415 else
00416 {
00417 throw cms::Exception("Aargh") << "DetIDs don't match either way :-( \n";
00418 }
00419 }
00420
00421 const TrackingRecHit *hit = it_hit->get();
00422 if ( hit == 0 || !hit->isValid() )
00423 continue;
00424
00425 int subdet = hit->geographicalId().subdetId();
00426
00427 if (subdet >= 3)
00428 {
00429 markClusters<SiStripCluster>(siStripDetsWithClusters, hit, tk, tms[reversed ? last_hit-i_hit : i_hit].updatedState());
00430 }
00431 else if (subdet >= 1)
00432 {
00433 markClusters<SiPixelCluster>(siPixelDetsWithClusters, hit, tk, tms[reversed ? last_hit-i_hit : i_hit].updatedState());
00434 }
00435 else
00436 {
00437 edm::LogWarning("HitNotFound") << "Hit of type " << typeid(*hit).name() << ", detid "
00438 << hit->geographicalId().rawId() << ", subdet " << subdet;
00439 }
00440 }
00441 }
00442 }
00443 else
00444 {
00445
00446
00447 Handle<std::vector<reco::Track> > tracks;
00448 iEvent.getByLabel(trajectories_, tracks);
00449
00450 foreach (const reco::Track &track, *tracks)
00451 {
00452 FreeTrajectoryState atVtx = trajectoryStateTransform::innerFreeState(track, &*magfield_);
00453 trackingRecHit_iterator it_hit = track.recHitsBegin(), ed_hit = track.recHitsEnd();
00454 for (; it_hit != ed_hit; ++it_hit)
00455 {
00456 const TrackingRecHit *hit = it_hit->get();
00457 if ( hit == 0 || !hit->isValid() )
00458 continue;
00459
00460 int subdet = hit->geographicalId().subdetId();
00461
00462 if ( subdet == 0 )
00463 continue;
00464
00465 const GeomDet *det = geometry_->idToDet( hit->geographicalId() );
00466
00467 if ( det == 0 )
00468 {
00469 edm::LogError("MissingDetId") << "DetIDs " << (int)(hit->geographicalId()) << " is not in geometry.\n";
00470 continue;
00471 }
00472
00473 TrajectoryStateOnSurface prop = propagator_->propagate(atVtx, det->surface());
00474 if ( subdet >= 3 )
00475 {
00476 markClusters<SiStripCluster>(siStripDetsWithClusters, hit, &track, prop);
00477 }
00478 else if (subdet >= 1)
00479 {
00480 markClusters<SiPixelCluster>(siPixelDetsWithClusters, hit, &track, prop);
00481 }
00482 else
00483 {
00484 edm::LogWarning("HitNotFound") << "Hit of type " << typeid(*hit).name() << ", detid "
00485 << hit->geographicalId().rawId() << ", subdet " << subdet;
00486 }
00487 }
00488 }
00489 }
00490
00491 Handle<std::vector<reco::Vertex> > vertices;
00492 iEvent.getByLabel(vertices_, vertices);
00493
00494
00495 if ( simSplitPixel_ )
00496 iEvent.getByLabel("simSiPixelDigis", pixeldigisimlink);
00497
00498
00499 if ( simSplitStrip_ )
00500 iEvent.getByLabel("simSiStripDigis", stripdigisimlink);
00501
00502
00503 std::auto_ptr<edmNew::DetSetVector<SiPixelCluster> > newPixelClusters( splitClusters( siPixelDetsWithClusters, vertices->front() ) );
00504 std::auto_ptr<edmNew::DetSetVector<SiStripCluster> > newStripClusters( splitClusters( siStripDetsWithClusters, vertices->front() ) );
00505
00506 if ( simSplitStrip_ )
00507 delete hitAssociator;
00508
00509 iEvent.put(newPixelClusters);
00510 iEvent.put(newStripClusters);
00511
00512 allSiPixelClusters.clear(); siPixelDetsWithClusters.clear();
00513 allSiStripClusters.clear(); siStripDetsWithClusters.clear();
00514
00515 }
00516
00517 template<typename Cluster>
00518 void TrackClusterSplitter::markClusters( std::map<uint32_t, boost::sub_range<std::vector<ClusterWithTracks<Cluster> > > >& map,
00519 const TrackingRecHit* hit,
00520 const reco::Track* track,
00521 const TrajectoryStateOnSurface& tsos) const
00522 {
00523 boost::sub_range<std::vector<ClusterWithTracks<Cluster> > >& range = map[hit->geographicalId().rawId()];
00524
00525 typedef typename std::vector<ClusterWithTracks<Cluster> >::iterator IT;
00526 IT match = std::find_if(range.begin(), range.end(), FindCluster<Cluster>(hit));
00527
00528 if ( match != range.end() )
00529 {
00530 match->tracks.push_back( TrackAndState(track,tsos) );
00531 }
00532 else
00533 {
00534 edm::LogWarning("ClusterNotFound") << "Cluster of type " << typeid(Cluster).name() << " on detid "
00535 << hit->geographicalId().rawId() << " from hit of type " << typeid(*hit).name();
00536 }
00537 }
00538
00539 template<typename Cluster>
00540 std::auto_ptr<edmNew::DetSetVector<Cluster> >
00541 TrackClusterSplitter::splitClusters(const std::map<uint32_t, boost::sub_range<std::vector<ClusterWithTracks<Cluster> > > > &input,
00542 const reco::Vertex &vtx) const
00543 {
00544 std::auto_ptr<edmNew::DetSetVector<Cluster> > output(new edmNew::DetSetVector<Cluster>());
00545 typedef std::pair<uint32_t, boost::sub_range<std::vector<ClusterWithTracks<Cluster> > > > pair;
00546
00547 foreach(const pair &p, input)
00548 {
00549 const GeomDet* det = geometry_->idToDet( DetId(p.first) );
00550
00551 if ( det == 0 )
00552 {
00553 edm::LogError("MissingDetId") << "DetIDs " << p.first << " is not in geometry.\n";
00554 continue;
00555 }
00556
00557
00558
00559 GlobalVector primary_vtx( vtx.x(), vtx.y(), vtx.z() );
00560
00561
00562 typename edmNew::DetSetVector<Cluster>::FastFiller detset(*output, p.first);
00563
00564
00565 foreach(const ClusterWithTracks<Cluster> &c, p.second)
00566 {
00567 splitCluster(c, primary_vtx, detset, DetId(p.first) );
00568 }
00569 }
00570
00571 return output;
00572 }
00573
00574 template<typename Cluster>
00575 void TrackClusterSplitter::splitCluster(const ClusterWithTracks<Cluster> &cluster,
00576 const GlobalVector &dir,
00577 typename edmNew::DetSetVector<Cluster>::FastFiller &output,
00578 DetId detId) const
00579 {
00580
00581 throw cms::Exception("LogicError", "This should not be called");
00582 }
00583
00584 template<>
00585 void TrackClusterSplitter::splitCluster<SiStripCluster> (const SiStripClusterWithTracks& c,
00586 const GlobalVector &vtx,
00587 edmNew::DetSetVector<SiStripCluster>::FastFiller &output,
00588 DetId detId
00589 ) const
00590 {
00591 if ( simSplitStrip_ )
00592 {
00593 bool cluster_was_successfully_split = false;
00594
00595 const SiStripCluster* clust = static_cast<const SiStripCluster*>(c.cluster);
00596
00597 std::vector<SimHitIdpr> associatedIdpr;
00598
00599 hitAssociator->associateSimpleRecHitCluster(clust, detId, associatedIdpr);
00600
00601 size_t splittableClusterSize = 0;
00602 splittableClusterSize = associatedIdpr.size();
00603 std::vector<uint8_t> amp = clust->amplitudes();
00604 int clusiz = amp.size();
00605 associatedIdpr.clear();
00606
00607 SiStripDetId ssdid( detId );
00608
00609
00610
00611 if ( ( splittableClusterSize > 1 && amp.size() > 2 ) &&
00612 ( (int)ssdid.moduleGeometry() == 1 ||
00613 (int)ssdid.moduleGeometry() == 2 ||
00614 (int)ssdid.moduleGeometry() == 3 ||
00615 (int)ssdid.moduleGeometry() == 4 ) )
00616 {
00617
00618 edm::DetSetVector<StripDigiSimLink>::const_iterator isearch = stripdigisimlink->find(detId);
00619
00620 int first = clust->firstStrip();
00621 int last = first + clusiz;
00622 uint16_t rawAmpl = 0, currentAmpl = 0;
00623
00624 std::vector<uint16_t> tmp1, tmp2;
00625
00626 std::vector<int> firstStrip;
00627 std::vector<bool> trackInStrip;
00628 std::vector<unsigned int> trackID;
00629 std::vector<float> trackFraction;
00630 std::vector< std::vector<uint16_t> > trackAmp;
00631 unsigned int currentChannel( 9999 );
00632 unsigned int thisTrackID = 0;
00633
00634 if ( isearch != stripdigisimlink->end() )
00635 {
00636 edm::DetSet<StripDigiSimLink> link_detset = (*isearch);
00637
00638 for ( edm::DetSet<StripDigiSimLink>::const_iterator linkiter = link_detset.data.begin();
00639 linkiter != link_detset.data.end(); linkiter++)
00640 {
00641 if ( (int)(linkiter->channel()) >= first && (int)(linkiter->channel()) < last )
00642 {
00643 int stripIdx = (int)linkiter->channel()-first;
00644 rawAmpl = (uint16_t)(amp[stripIdx]);
00645
00646
00647
00648 if ( linkiter->channel() != currentChannel )
00649 {
00650
00651 uint16_t thisAmpl;
00652
00653 for (size_t i=0; i<trackID.size(); ++i)
00654 {
00655 if ( trackInStrip[i] )
00656 {
00657 if ( ( thisAmpl=currentAmpl ) < 254 )
00658 thisAmpl = min( uint16_t(253), max(uint16_t(0), (uint16_t)(currentAmpl*trackFraction[i]+0.5)) );
00659 trackAmp[i].push_back( thisAmpl );
00660 }
00661
00662 trackFraction[i] = 0;
00663 trackInStrip[i] = false;
00664 }
00665
00666 currentChannel = linkiter->channel();
00667 currentAmpl = rawAmpl;
00668 }
00669
00670
00671 thisTrackID = linkiter->SimTrackId();
00672
00673
00674 bool newTrack = true;
00675 int thisTrackIdx = 9999;
00676
00677 for (size_t i=0; i<trackID.size(); ++i)
00678 {
00679 if ( trackID[i] == thisTrackID )
00680 {
00681 thisTrackIdx = i;
00682 newTrack = false;
00683 }
00684 }
00685
00686 if ( newTrack )
00687 {
00688 trackInStrip.push_back(false);
00689 trackID.push_back(thisTrackID);
00690 firstStrip.push_back(currentChannel);
00691 std::vector<uint16_t> ampTmp;
00692 trackAmp.push_back(ampTmp);
00693 trackFraction.push_back(0);
00694 thisTrackIdx = trackID.size()-1;
00695 }
00696
00697 trackInStrip[thisTrackIdx] = true;
00698 trackFraction[thisTrackIdx] += linkiter->fraction();
00699 currentAmpl = rawAmpl;
00700
00701 }
00702
00703 }
00704
00705
00706
00707 std::vector<SiStripCluster> newCluster;
00708
00709 uint16_t thisAmpl;
00710
00711 for (size_t i=0; i < trackID.size(); ++i)
00712 {
00713 if ( trackInStrip[i] )
00714 {
00715 if ( ( thisAmpl=rawAmpl ) < 254 )
00716 thisAmpl = min(uint16_t(253), max(uint16_t(0), (uint16_t)(rawAmpl*trackFraction[i]+0.5)));
00717
00718 if ( thisAmpl > 0 )
00719 trackAmp[i].push_back( thisAmpl );
00720
00721
00722 }
00723
00724 newCluster.push_back( SiStripCluster( detId,
00725 firstStrip[i],
00726 trackAmp[i].begin(),
00727 trackAmp[i].end() ) );
00728
00729 }
00730
00731
00732 for ( size_t i=0; i<newCluster.size(); ++i )
00733 {
00734
00735
00736 float clusterAmp = 0.0;
00737 for (size_t j=0; j<trackAmp[i].size(); ++j )
00738 clusterAmp += (float)(trackAmp[i])[j];
00739
00740 if ( clusterAmp > 0.0 && firstStrip[i] != 9999 && trackAmp[i].size() > 0 )
00741 {
00742
00743 output.push_back( newCluster[i] );
00744
00745
00746
00747
00748
00749
00750 cluster_was_successfully_split = true;
00751 }
00752 else
00753 {
00754
00755
00756
00757
00758 }
00759 }
00760
00761 }
00762 else
00763 {
00764
00765 }
00766 }
00767
00768
00769 if ( !cluster_was_successfully_split )
00770 output.push_back( *c.cluster );
00771
00772 }
00773
00774 else if ( tmpSplitStrip_ )
00775 {
00776 bool cluster_was_successfully_split = false;
00777
00778 const SiStripCluster* theStripCluster = static_cast<const SiStripCluster*>(c.cluster);
00779
00780 if ( theStripCluster )
00781 {
00782
00783 SiStripDetId ssdid( detId.rawId() );
00784
00785
00786
00787
00788 if ( (int)theStripCluster->amplitudes().size() <= 1 ||
00789 ( (int)ssdid.moduleGeometry() != 1 &&
00790 (int)ssdid.moduleGeometry() != 2 &&
00791 (int)ssdid.moduleGeometry() != 3 &&
00792 (int)ssdid.moduleGeometry() != 4 ) )
00793 {
00794
00795
00796
00797
00798
00799
00800 }
00801 else
00802 {
00803
00804
00805
00806
00807 int ID = -99999;
00808
00809 int is_stereo = (int)( ssdid.stereo() );
00810
00811 if ( ssdid.moduleGeometry() == 1 )
00812 {
00813 if ( !is_stereo )
00814 ID = 11;
00815 else
00816 ID = 12;
00817 }
00818 else if ( ssdid.moduleGeometry() == 2 )
00819 {
00820 ID = 13;
00821 }
00822 else if ( ssdid.moduleGeometry() == 3 )
00823 {
00824 ID = 16;
00825 }
00826 else if ( ssdid.moduleGeometry() == 4 )
00827 {
00828 if ( !is_stereo )
00829 ID = 14;
00830 else
00831 ID = 15;
00832 }
00833 else
00834 {
00835 throw cms::Exception("TrackClusterSplitter::splitCluster")
00836 << "\nERROR: Wrong strip teplate ID. Should only use templates for IB1, IB2, OB1 and OB2 !!!" << "\n\n";
00837 }
00838
00839
00840
00841
00842
00843
00844 float cotalpha_ = -99999.9;
00845 float cotbeta_ = -99999.9;
00846
00847
00848
00849
00850
00851 float xcenter = theStripCluster->barycenter();
00852
00853 const GeomDetUnit* theDet = geometry_->idToDetUnit( detId );
00854 const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>( theDet );
00855
00856 if ( !stripDet )
00857 {
00858 throw cms::Exception("TrackClusterSplitter : ")
00859 << "\nERROR: Wrong stripDet !!! " << "\n\n";
00860 }
00861
00862
00863 const StripTopology* theTopol = &( stripDet->specificTopology() );
00864
00865
00866
00867
00868
00869
00870 LocalPoint lp = theTopol->localPosition( xcenter );
00871
00872
00873 GlobalPoint gp0 = theDet->surface().toGlobal( lp );
00874
00875
00876 GlobalPoint gp(gp0.x()-vtx.x(), gp0.y()-vtx.y(), gp0.z()-vtx.z() );
00877
00878
00879 float gp_mod = sqrt( gp.x()*gp.x() + gp.y()*gp.y() + gp.z()*gp.z() );
00880 float gpx = gp.x()/gp_mod;
00881 float gpy = gp.y()/gp_mod;
00882 float gpz = gp.z()/gp_mod;
00883 GlobalVector gv(gpx, gpy, gpz);
00884
00885
00886 const Local3DVector lvx(1.0, 0.0, 0.0);
00887 GlobalVector gvx = theDet->surface().toGlobal( lvx );
00888 const Local3DVector lvy(0.0, 1.0, 0.0);
00889 GlobalVector gvy = theDet->surface().toGlobal( lvy );
00890 const Local3DVector lvz(0.0, 0.0, 1.0);
00891 GlobalVector gvz = theDet->surface().toGlobal( lvz );
00892
00893
00894 float gv_dot_gvx = gv.x()*gvx.x() + gv.y()*gvx.y() + gv.z()*gvx.z();
00895 float gv_dot_gvy = gv.x()*gvy.x() + gv.y()*gvy.y() + gv.z()*gvy.z();
00896 float gv_dot_gvz = gv.x()*gvz.x() + gv.y()*gvz.y() + gv.z()*gvz.z();
00897
00898
00899
00900
00901
00902 cotalpha_ = gv_dot_gvx / gv_dot_gvz;
00903 cotbeta_ = gv_dot_gvy / gv_dot_gvz;
00904
00905
00906 if ( !useStraightTracks_ )
00907 {
00908
00909
00910
00911
00912
00913
00914
00915 std::vector<TrackAndState> vec_tracks_states = c.tracks;
00916
00917 if ( (int)vec_tracks_states.size() > 0 )
00918 {
00919
00920
00921
00922 int index_max_pt = -99999;
00923 float max_pt = -99999.9;
00924
00925 for (int i=0; i<(int)vec_tracks_states.size(); ++i )
00926 {
00927 const reco::Track* one_track = vec_tracks_states[i].track;
00928
00929 if ( one_track->pt() > max_pt )
00930 {
00931 index_max_pt = i;
00932 max_pt = one_track->pt();
00933 }
00934 }
00935
00936
00937
00938 TrajectoryStateOnSurface one_tsos = vec_tracks_states[index_max_pt].state;
00939
00940 LocalTrajectoryParameters ltp = one_tsos.localParameters();
00941
00942 LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
00943
00944 float locx = localDir.x();
00945 float locy = localDir.y();
00946 float locz = localDir.z();
00947
00948
00949
00950
00951 cotalpha_ = locx/locz;
00952 cotbeta_ = locy/locz;
00953
00954 }
00955
00956 }
00957
00958
00959
00960
00961
00962
00963
00964
00965
00966
00967 float strip_cluster_charge = 0.0;
00968 std::vector<float> vec_cluster_charge;
00969 vec_cluster_charge.clear();
00970 int cluster_size = (int)( (theStripCluster->amplitudes()).size() );
00971
00972 int cluster_charge = 0;
00973 for (int i=0; i<cluster_size; ++i)
00974 {
00975 float current_strip_charge = (float)( (theStripCluster->amplitudes())[i] );
00976
00977 strip_cluster_charge += current_strip_charge;
00978 vec_cluster_charge.push_back( current_strip_charge );
00979
00980 cluster_charge +=current_strip_charge;
00981 }
00982
00983
00984
00985
00986 int strip_templQbin_ = strip_templ_.qbin( ID, cotalpha_, cotbeta_, strip_cluster_charge );
00987
00988 if ( strip_templQbin_ < 0 || strip_templQbin_ > 5 )
00989 {
00990
00991
00992 }
00993 else
00994 {
00995 if ( strip_templQbin_ != 0 )
00996 {
00997
00998
00999
01000
01001 }
01002 else
01003 {
01004
01005
01006
01007
01008 uint16_t first_strip = theStripCluster->firstStrip();
01009
01010 LocalVector lbfield = ( stripDet->surface() ).toLocal( magfield_->inTesla( stripDet->surface().position() ) );
01011 float locBy = lbfield.y();
01012
01013
01014 float stripTemplXrec1_ = -99999.9;
01015 float stripTemplXrec2_ = -99999.9;
01016 float stripTemplSigmaX_ = -99999.9;
01017 float stripTemplProbX_ = -99999.9;
01018 int stripTemplQbin_ = -99999;
01019 float stripTemplProbQ_ = -99999.9;
01020
01021
01022
01023
01024
01025
01026
01027
01028
01029
01030
01031
01032
01033
01034
01035
01036
01037 int ierr =
01038 SiStripTemplateSplit::StripTempSplit(ID,
01039 cotalpha_, cotbeta_,
01040 locBy,
01041 vec_cluster_charge,
01042 strip_templ_,
01043 stripTemplXrec1_,
01044 stripTemplXrec2_,
01045 stripTemplSigmaX_,
01046 stripTemplProbX_,
01047 stripTemplQbin_,
01048 stripTemplProbQ_ );
01049
01050
01051
01052 stripTemplXrec1_ += 2*strip_templ_.xsize();
01053 stripTemplXrec2_ += 2*strip_templ_.xsize();
01054
01055
01056
01057 if ( ierr != 0 )
01058 {
01059
01060
01061 }
01062 else
01063 {
01064
01065
01066
01067
01068
01069
01070 cluster_was_successfully_split = true;
01071
01072 std::vector<float> strip_cluster1;
01073 std::vector<float> strip_cluster2;
01074
01075 strip_cluster1.clear();
01076 strip_cluster2.clear();
01077
01078
01079 for (int i=0; i<BXSIZE; ++i)
01080 {
01081 strip_cluster1.push_back(0.0);
01082 strip_cluster2.push_back(0.0);
01083 }
01084
01085
01086
01087
01088 strip_templ_.interpolate(ID, cotalpha_, cotbeta_, locBy);
01089 strip_templ_.sxtemp(stripTemplXrec1_, strip_cluster1);
01090 strip_templ_.sxtemp(stripTemplXrec2_, strip_cluster2);
01091
01092
01093
01094 vector<SiStripDigi> vecSiStripDigi1;
01095 vecSiStripDigi1.clear();
01096 int strip_cluster1_size = (int)strip_cluster1.size();
01097 for (int i=2; i<strip_cluster1_size; ++i)
01098 {
01099 if ( strip_cluster1[i] > 0.0 )
01100 {
01101 SiStripDigi current_digi1( first_strip + i-2, strip_cluster1[i] );
01102
01103 vecSiStripDigi1.push_back( current_digi1 );
01104 }
01105 }
01106
01107
01108
01109 vector<SiStripDigi> vecSiStripDigi2;
01110 vecSiStripDigi2.clear();
01111 int strip_cluster2_size = (int)strip_cluster2.size();
01112 for (int i=2; i<strip_cluster2_size; ++i)
01113 {
01114 if ( strip_cluster2[i] > 0.0 )
01115 {
01116 SiStripDigi current_digi2( first_strip + i-2, strip_cluster2[i] );
01117
01118 vecSiStripDigi2.push_back( current_digi2 );
01119 }
01120 }
01121
01122
01123
01124
01125 std::vector<SiStripDigi>::const_iterator SiStripDigiIterBegin1 = vecSiStripDigi1.begin();
01126 std::vector<SiStripDigi>::const_iterator SiStripDigiIterEnd1 = vecSiStripDigi1.end();
01127 std::pair<std::vector<SiStripDigi>::const_iterator,
01128 std::vector<SiStripDigi>::const_iterator> SiStripDigiRange1
01129 = make_pair(SiStripDigiIterBegin1, SiStripDigiIterEnd1);
01130
01131
01132
01133
01134
01135
01136
01137 std::vector<SiStripDigi>::const_iterator SiStripDigiIterBegin2 = vecSiStripDigi2.begin();
01138 std::vector<SiStripDigi>::const_iterator SiStripDigiIterEnd2 = vecSiStripDigi2.end();
01139 std::pair<std::vector<SiStripDigi>::const_iterator,
01140 std::vector<SiStripDigi>::const_iterator> SiStripDigiRange2
01141 = make_pair(SiStripDigiIterBegin2, SiStripDigiIterEnd2);
01142
01143
01144
01145
01146
01147
01148
01149
01150
01151
01152
01153 if ( SiStripDigiIterBegin1 != SiStripDigiIterEnd1 )
01154 {
01155
01156 SiStripCluster cl1( detId.rawId(), SiStripDigiRange1 );
01157
01158 cl1.setSplitClusterError( stripTemplSigmaX_ );
01159
01160 output.push_back( cl1 );
01161
01162 if ( (int)cl1.amplitudes().size() <= 0 )
01163 {
01164 throw cms::Exception("TrackClusterSplitter : ")
01165 << "\nERROR: (1) Wrong split cluster of size = " << (int)cl1.amplitudes().size() << "\n\n";
01166 }
01167
01168 }
01169
01170 if ( SiStripDigiIterBegin2 != SiStripDigiIterEnd2 )
01171 {
01172
01173 SiStripCluster cl2( detId.rawId(), SiStripDigiRange2 );
01174 cl2.setSplitClusterError( stripTemplSigmaX_ );
01175 output.push_back( cl2 );
01176
01177 if ( (int)cl2.amplitudes().size() <= 0 )
01178 {
01179 throw cms::Exception("TrackClusterSplitter : ")
01180 << "\nERROR: (2) Wrong split cluster of size = " << (int)cl2.amplitudes().size() << "\n\n";
01181 }
01182
01183 }
01184
01185
01186
01187 }
01188
01189 }
01190
01191 }
01192
01193 }
01194
01195
01196 if ( !cluster_was_successfully_split )
01197 output.push_back( *c.cluster );
01198
01199 }
01200 else
01201 {
01202 throw cms::Exception("TrackClusterSplitter : ")
01203 << "\nERROR: This is not a SiStripCluster !!!" << "\n\n";
01204 }
01205
01206 }
01207 else
01208 {
01209
01210 output.push_back( *c.cluster );
01211 }
01212 }
01213
01214 template<>
01215 void TrackClusterSplitter::splitCluster<SiPixelCluster> (const SiPixelClusterWithTracks &c,
01216 const GlobalVector &vtx,
01217 edmNew::DetSetVector<SiPixelCluster>::FastFiller &output,
01218 DetId detId
01219 ) const
01220 {
01221
01222 if ( simSplitPixel_ )
01223 {
01224
01225
01226 int minPixelRow = (*c.cluster).minPixelRow();
01227 int maxPixelRow = (*c.cluster).maxPixelRow();
01228 int minPixelCol = (*c.cluster).minPixelCol();
01229 int maxPixelCol = (*c.cluster).maxPixelCol();
01230 int dsl = 0;
01231
01232 edm::DetSetVector<PixelDigiSimLink>::const_iterator isearch = pixeldigisimlink->find(output.id());
01233 edm::DetSet<PixelDigiSimLink> digiLink = (*isearch);
01234
01235 edm::DetSet<PixelDigiSimLink>::const_iterator linkiter = digiLink.data.begin();
01236
01237 std::vector<int> simTrackIdV;
01238 simTrackIdV.clear();
01239
01240 std::vector<SiPixelCluster> splittedCluster;
01241 splittedCluster.clear();
01242
01243 for ( ; linkiter != digiLink.data.end(); linkiter++)
01244 {
01245 dsl++;
01246 std::pair<int,int> pixel_coord = PixelDigi::channelToPixel(linkiter->channel());
01247
01248
01249 if ( pixel_coord.first <= maxPixelRow &&
01250 pixel_coord.first >= minPixelRow &&
01251 pixel_coord.second <= maxPixelCol &&
01252 pixel_coord.second >= minPixelCol )
01253 {
01254 bool inStock(false);
01255
01256 SiPixelCluster::PixelPos newPixelPos(pixel_coord.first, pixel_coord.second);
01257
01258
01259 int newPixelCharge(0);
01260
01261 const std::vector<SiPixelCluster::Pixel>& pixvector = (*c.cluster).pixels();
01262
01263 for(std::vector<SiPixelCluster::Pixel>::const_iterator itPix = pixvector.begin(); itPix != pixvector.end(); itPix++)
01264 {
01265 if (((int) itPix->x) == ((int) pixel_coord.first)&&(((int) itPix->y) == ((int) pixel_coord.second)))
01266 {
01267 newPixelCharge = (int) (linkiter->fraction()*itPix->adc);
01268 }
01269 }
01270
01271 if ( newPixelCharge < 2500 )
01272 continue;
01273
01274
01275 int clusVecPos = 0;
01276 std::vector<int>::const_iterator sTIter = simTrackIdV.begin();
01277
01278 for ( ; sTIter < simTrackIdV.end(); sTIter++)
01279 {
01280 if (((*sTIter)== (int) linkiter->SimTrackId()))
01281 {
01282 inStock=true;
01283
01284
01285 splittedCluster.at(clusVecPos).add(newPixelPos,newPixelCharge);
01286 }
01287 clusVecPos++;
01288 }
01289
01290
01291
01292 if ( !inStock )
01293 {
01294
01295 simTrackIdV.push_back(linkiter->SimTrackId());
01296 splittedCluster.push_back(SiPixelCluster(newPixelPos,newPixelCharge));
01297 }
01298 }
01299 }
01300
01301
01302
01303 if ( ( ( (int)simTrackIdV.size() ) == 1 ) || ( *c.cluster).size()==1 )
01304 {
01305
01306 output.push_back(*c.cluster );
01307
01308 }
01309 else
01310 {
01311 for (std::vector<SiPixelCluster>::const_iterator cIter = splittedCluster.begin(); cIter != splittedCluster.end(); cIter++ )
01312 {
01313 output.push_back( (*cIter) );
01314 }
01315 }
01316
01317 simTrackIdV.clear();
01318 splittedCluster.clear();
01319 }
01320 else if ( tmpSplitPixel_ )
01321 {
01322 bool cluster_was_successfully_split = false;
01323
01324 const SiPixelCluster* thePixelCluster = static_cast<const SiPixelCluster*>(c.cluster);
01325
01326 if ( thePixelCluster )
01327 {
01328
01329 if ( (int)thePixelCluster->size() <= 1 )
01330 {
01331
01332
01333
01334 }
01335 else
01336 {
01337
01338 int ID = -99999;
01339 if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
01340 {
01341
01342 ID = 40;
01343 }
01344 else if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelEndcap )
01345 {
01346
01347 ID = 41;
01348 }
01349 else
01350 {
01351
01352 }
01353
01354
01355
01356
01357 float cotalpha_ = -99999.9;
01358 float cotbeta_ = -99999.9;
01359
01360
01361
01362
01363 float xcenter = thePixelCluster->x();
01364 float ycenter = thePixelCluster->y();
01365
01366 const GeomDetUnit* theDet = geometry_->idToDetUnit( detId );
01367 const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>( theDet );
01368
01369 const PixelTopology* theTopol = (&(pixDet->specificTopology()));
01370
01371
01372 LocalPoint lp = theTopol->localPosition( MeasurementPoint(xcenter, ycenter) );
01373
01374
01375 GlobalPoint gp0 = theDet->surface().toGlobal( lp );
01376
01377
01378 GlobalPoint gp(gp0.x()-vtx.x(), gp0.y()-vtx.y(), gp0.z()-vtx.z() );
01379
01380
01381 float gp_mod = sqrt( gp.x()*gp.x() + gp.y()*gp.y() + gp.z()*gp.z() );
01382 float gpx = gp.x()/gp_mod;
01383 float gpy = gp.y()/gp_mod;
01384 float gpz = gp.z()/gp_mod;
01385 GlobalVector gv(gpx, gpy, gpz);
01386
01387
01388 const Local3DVector lvx(1.0, 0.0, 0.0);
01389 GlobalVector gvx = theDet->surface().toGlobal( lvx );
01390 const Local3DVector lvy(0.0, 1.0, 0.0);
01391 GlobalVector gvy = theDet->surface().toGlobal( lvy );
01392 const Local3DVector lvz(0.0, 0.0, 1.0);
01393 GlobalVector gvz = theDet->surface().toGlobal( lvz );
01394
01395
01396 float gv_dot_gvx = gv.x()*gvx.x() + gv.y()*gvx.y() + gv.z()*gvx.z();
01397 float gv_dot_gvy = gv.x()*gvy.x() + gv.y()*gvy.y() + gv.z()*gvy.z();
01398 float gv_dot_gvz = gv.x()*gvz.x() + gv.y()*gvz.y() + gv.z()*gvz.z();
01399
01400
01401
01402
01403 cotalpha_ = gv_dot_gvx / gv_dot_gvz;
01404 cotbeta_ = gv_dot_gvy / gv_dot_gvz;
01405
01406
01407
01408
01409
01410 if ( !useStraightTracks_ )
01411 {
01412
01413
01414
01415
01416 std::vector<TrackAndState> vec_tracks_states = c.tracks;
01417
01418
01419
01420 if ( (int)vec_tracks_states.size() > 0 )
01421 {
01422
01423
01424
01425 int index_max_pt = -99999;
01426 float max_pt = -99999.9;
01427
01428 for (int i=0; i<(int)vec_tracks_states.size(); ++i )
01429 {
01430 const reco::Track* one_track = vec_tracks_states[i].track;
01431
01432 if ( one_track->pt() > max_pt )
01433 {
01434 index_max_pt = i;
01435 max_pt = one_track->pt();
01436 }
01437 }
01438
01439
01440
01441 TrajectoryStateOnSurface one_tsos = vec_tracks_states[index_max_pt].state;
01442
01443 LocalTrajectoryParameters ltp = one_tsos.localParameters();
01444
01445 LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
01446
01447 float locx = localDir.x();
01448 float locy = localDir.y();
01449 float locz = localDir.z();
01450
01451
01452
01453
01454 cotalpha_ = locx/locz;
01455 cotbeta_ = locy/locz;
01456
01457
01458
01459 }
01460
01461 }
01462
01463
01464
01465
01466
01467
01468 int templQbin_ = templ_.qbin( ID, cotalpha_, cotbeta_, thePixelCluster->charge() );
01469
01470 if ( templQbin_ < 0 || templQbin_ > 5 )
01471 {
01472
01473
01474
01475 }
01476 else
01477 {
01478
01479
01480
01481 if ( templQbin_ != 0 )
01482 {
01483
01484
01485
01486 }
01487 else
01488 {
01489
01490
01491
01492
01493
01494
01495
01496
01497 const std::vector<SiPixelCluster::Pixel> & pixVec = thePixelCluster->pixels();
01498 std::vector<SiPixelCluster::Pixel>::const_iterator pixIter = pixVec.begin(), pixEnd = pixVec.end();
01499
01500 const int cluster_matrix_size_x = 13;
01501 const int cluster_matrix_size_y = 21;
01502
01503 boost::multi_array<float, 2> clust_array_2d(boost::extents[cluster_matrix_size_x][cluster_matrix_size_y]);
01504
01505 int row_offset = thePixelCluster->minPixelRow();
01506 int col_offset = thePixelCluster->minPixelCol();
01507
01508
01509
01510 for ( ; pixIter != pixEnd; ++pixIter )
01511 {
01512 int irow = int(pixIter->x) - row_offset;
01513 int icol = int(pixIter->y) - col_offset;
01514
01515 if ( irow<cluster_matrix_size_x && icol<cluster_matrix_size_y )
01516 {
01517 clust_array_2d[irow][icol] = (float)pixIter->adc;
01518 }
01519 }
01520
01521
01522 std::vector<bool> ydouble(cluster_matrix_size_y), xdouble(cluster_matrix_size_x);
01523
01524
01525 for (int irow = 0; irow < cluster_matrix_size_x; ++irow)
01526 {
01527 xdouble[irow] = theTopol->isItBigPixelInX( irow+row_offset );
01528 }
01529
01530
01531 for (int icol = 0; icol < cluster_matrix_size_y; ++icol)
01532 {
01533 ydouble[icol] = theTopol->isItBigPixelInY( icol+col_offset );
01534 }
01535
01536
01537 float templYrec1_ = -99999.9;
01538 float templYrec2_ = -99999.9;
01539 float templSigmaY_ = -99999.9;
01540 float templProbY_ = -99999.9;
01541 float templXrec1_ = -99999.9;
01542 float templXrec2_ = -99999.9;
01543 float templSigmaX_ = -99999.9;
01544 float templProbX_ = -99999.9;
01545 float dchisq = -99999.9;
01546 float templProbQ_ = -99999.9;
01547
01548 int ierr =
01549 SiPixelTemplateSplit::PixelTempSplit( ID,
01550 cotalpha_, cotbeta_,
01551 clust_array_2d,
01552 ydouble, xdouble,
01553 templ_,
01554 templYrec1_, templYrec2_, templSigmaY_, templProbY_,
01555 templXrec1_, templXrec2_, templSigmaX_, templProbX_,
01556 templQbin_,
01557 templProbQ_,
01558 true,
01559 dchisq,
01560 templ2D_ );
01561
01562 if ( ierr != 0 )
01563 {
01564
01565 }
01566 else
01567 {
01568
01569
01570
01571
01572
01573
01574
01575 float xsize = templ_.xsize();
01576 float ysize = templ_.ysize();
01577
01578
01579 float yrecp1 = -99999.9;
01580 float yrecp2 = -99999.9;
01581 float xrecp1 = -99999.9;
01582 float xrecp2 = -99999.9;
01583
01584 if ( ydouble[0] )
01585 {
01586 yrecp1 = templYrec1_ + ysize;
01587 yrecp2 = templYrec2_ + ysize;
01588 }
01589 else
01590 {
01591 yrecp1 = templYrec1_ + ysize/2.0;
01592 yrecp2 = templYrec2_ + ysize/2.0;
01593 }
01594
01595 if ( xdouble[0] )
01596 {
01597 xrecp1 = templXrec1_ + xsize;
01598 xrecp2 = templXrec2_ + xsize;
01599 }
01600 else
01601 {
01602 xrecp1 = templXrec1_ + xsize/2.0;
01603 xrecp2 = templXrec2_ + xsize/2.0;
01604 }
01605
01606
01607
01608 float template2d1[BXM2][BYM2];
01609 float template2d2[BXM2][BYM2];
01610
01611 for ( int j=0; j<BXM2; ++j )
01612 {
01613 for ( int i=0; i<BYM2; ++i )
01614 {
01615 template2d1[j][i] = 0.0;
01616 template2d2[j][i] = 0.0;
01617 }
01618 }
01619
01620
01621 bool template_OK
01622 = templ2D_.xytemp(ID, cotalpha_, cotbeta_,
01623 xrecp1, yrecp1,
01624 ydouble, xdouble,
01625 template2d1);
01626
01627 template_OK
01628 = template_OK &&
01629 templ2D_.xytemp(ID, cotalpha_, cotbeta_,
01630 xrecp2, yrecp2,
01631 ydouble, xdouble,
01632 template2d2);
01633
01634 if ( !template_OK )
01635 {
01636
01637
01638
01639 for ( int j=0; j<BXM2; ++j )
01640 {
01641 for ( int i=0; i<BYM2; ++i )
01642 {
01643 template2d1[j][i] = 0.0;
01644 template2d2[j][i] = 0.0;
01645 }
01646 }
01647
01648 if ( !templ_.simpletemplate2D(xrecp1, yrecp1,
01649 xdouble, ydouble,
01650 template2d1) )
01651 {
01652 cluster_was_successfully_split = false;
01653 }
01654
01655 if ( !templ_.simpletemplate2D(xrecp2, yrecp2,
01656 xdouble, ydouble,
01657 template2d2) )
01658 {
01659 cluster_was_successfully_split = false;
01660 }
01661
01662 }
01663 else
01664 {
01665 cluster_was_successfully_split = true;
01666
01667
01668
01669
01670 float q50 = templ_.s50();
01671
01672
01673 float cluster1[TXSIZE][TYSIZE];
01674 float cluster2[TXSIZE][TYSIZE];
01675
01676 for ( int j=0; j<TXSIZE; ++j )
01677 {
01678 for ( int i=0; i<TYSIZE; ++i )
01679 {
01680 cluster1[j][i] = template2d1[j+1][i+1];
01681
01682 if ( cluster1[j][i] < q50 )
01683 cluster1[j][i] = 0.0;
01684
01685 cluster2[j][i] = template2d2[j+1][i+1];
01686
01687 if ( cluster2[j][i] < q50 )
01688 cluster2[j][i] = 0.0;
01689
01690
01691
01692 }
01693 }
01694
01695
01696
01697 int i1 = -999;
01698 int j1 = -999;
01699 int i2 = -999;
01700 int j2 = -999;
01701
01702 bool done_searching = false;
01703 for ( int i=0; i<13 && !done_searching; ++i)
01704 {
01705 for (int j=0; j<21 && !done_searching; ++j)
01706 {
01707 if ( cluster1[i][j] > 0 )
01708 {
01709 i1 = i;
01710 j1 = j;
01711 done_searching = true;
01712 }
01713 }
01714 }
01715
01716 done_searching = false;
01717 for ( int i=0; i<13 && !done_searching; ++i)
01718 {
01719 for (int j=0; j<21 && !done_searching; ++j)
01720 {
01721 if ( cluster2[i][j] > 0 )
01722 {
01723 i2 = i;
01724 j2 = j;
01725 done_searching = true;
01726 }
01727 }
01728 }
01729
01730
01731
01732
01733 SiPixelCluster cl1( SiPixelCluster::PixelPos( i1 + row_offset, j1 + col_offset),
01734 cluster1[i1][j1] );
01735
01736 SiPixelCluster cl2( SiPixelCluster::PixelPos( i2 + row_offset, j2 + col_offset),
01737 cluster2[i2][j2] );
01738
01739
01740
01741
01742 for ( int i=0; i<13; ++i)
01743 {
01744 for (int j=0; j<21; ++j)
01745 {
01746
01747 if ( cluster1[i][j] > 0.0 && (i!=i1 || j!=j1 ) )
01748 {
01749 cl1.add( SiPixelCluster::PixelPos( i + row_offset, j + col_offset),
01750 cluster1[i][j] );
01751
01752
01753 }
01754
01755
01756 if ( cluster2[i][j] > 0.0 && (i!=i2 || j!=j2 ) )
01757 {
01758 cl2.add( SiPixelCluster::PixelPos( i + row_offset, j + col_offset),
01759 cluster2[i][j] );
01760
01761
01762 }
01763 }
01764 }
01765
01766
01767
01768
01769
01770 cl1.setSplitClusterErrorX( templSigmaX_ );
01771 cl1.setSplitClusterErrorY( templSigmaY_ );
01772
01773
01774 cl2.setSplitClusterErrorX( templSigmaX_ );
01775 cl2.setSplitClusterErrorY( templSigmaY_ );
01776
01777
01778
01779
01780 output.push_back( cl1 );
01781 output.push_back( cl2 );
01782
01783
01784
01785 if ( (int)cl1.size() <= 0 )
01786 {
01787 edm::LogError("TrackClusterSplitter : ")
01788 << "1) Cluster of size = " << (int)cl1.size() << " !!! " << endl;
01789 }
01790
01791 if ( (int)cl2.size() <= 0 )
01792 {
01793 edm::LogError("TrackClusterSplitter : ")
01794 << "2) Cluster of size = " << (int)cl2.size() << " !!! " << endl;
01795 }
01796
01797 if ( cl1.charge() <= 0 )
01798 {
01799 edm::LogError("TrackClusterSplitter : ")
01800 << "1) Cluster of charge = " << (int)cl1.charge() << " !!! " << endl;
01801
01802 }
01803
01804 if ( cl2.charge() <= 0 )
01805 {
01806 edm::LogError("TrackClusterSplitter : ")
01807 << "2) Cluster of charge = " << (int)cl2.charge() << " !!!!!!!!!!!!!!!!!!!!!!!!!! " << endl;
01808 }
01809
01810
01811 }
01812
01813 }
01814
01815 }
01816
01817 }
01818
01819 }
01820
01821 if ( !cluster_was_successfully_split )
01822 output.push_back(*c.cluster);
01823
01824 }
01825 else
01826 {
01827 throw cms::Exception("TrackClusterSplitter :")
01828 << "This is not a SiPixelCluster !!! " << "\n";
01829 }
01830 }
01831 else
01832 {
01833
01834
01835 output.push_back( *c.cluster );
01836 }
01837 }
01838
01839 #include "FWCore/PluginManager/interface/ModuleDef.h"
01840 #include "FWCore/Framework/interface/MakerMacros.h"
01841 DEFINE_FWK_MODULE(TrackClusterSplitter);