CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
NuclearInteractionEDProducer.cc
Go to the documentation of this file.
3 
5 
8 
9 
11 
12 
14 conf_(iConfig)
15 {
16  token_primaryTrack = consumes<reco::TrackCollection>(iConfig.getParameter<std::string>("primaryProducer"));
17  token_secondaryTrack = consumes<reco::TrackCollection>(iConfig.getParameter<std::string>("secondaryProducer"));
18  token_additionalSecTracks = consumes<reco::TrackCollection>(iConfig.getParameter<std::string>("additionalSecondaryProducer"));
19  token_primaryTrajectory = consumes<TrajectoryCollection>(iConfig.getParameter<std::string>("primaryProducer"));
20  token_refMapH = consumes<TrajTrackAssociationCollection>(iConfig.getParameter<std::string>("primaryProducer"));
21  token_nuclMapH = consumes<TrajectoryToSeedsMap>(iConfig.getParameter<std::string>("seedsProducer"));
22 
23  produces<reco::NuclearInteractionCollection>();
24 }
25 
27 {
28 }
29 
30 //
31 // member functions
32 //
33 
34 // ------------ method called to produce the data ------------
35 void
37 {
38  if ( magFieldWatcher_.check(iSetup) || transientTrackWatcher_.check(iSetup) ) {
41  iSetup.get<IdealMagneticFieldRecord>().get(magField);
42 
44  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
45 
46  vertexBuilder = std::auto_ptr< NuclearVertexBuilder >(new NuclearVertexBuilder( magField.product(), builder.product(), conf_) );
47  likelihoodCalculator = std::auto_ptr< NuclearLikelihood >(new NuclearLikelihood);
48  }
49 
51  edm::Handle<reco::TrackCollection> primaryTrackCollection;
52  iEvent.getByToken(token_primaryTrack, primaryTrackCollection);
53 
55  edm::Handle< TrajectoryCollection > primaryTrajectoryCollection;
56  iEvent.getByToken(token_primaryTrajectory, primaryTrajectoryCollection);
57 
60  iEvent.getByToken(token_refMapH, refMapH);
61  const TrajTrackAssociationCollection& refMap = *(refMapH.product());
62 
65  iEvent.getByToken(token_nuclMapH, nuclMapH);
66  const TrajectoryToSeedsMap& nuclMap = *(nuclMapH.product());
67 
69  edm::Handle<reco::TrackCollection> secondaryTrackCollection;
70  iEvent.getByToken(token_secondaryTrack, secondaryTrackCollection);
71 
72  // Get eventual additional secondary tracks
73  edm::Handle<reco::TrackCollection> additionalSecTracks;
74  iEvent.getByToken(token_additionalSecTracks, additionalSecTracks);
75 
77  std::auto_ptr<reco::NuclearInteractionCollection> theNuclearInteractions(new reco::NuclearInteractionCollection);
78 
79  typedef edm::Ref<TrajectoryCollection> TrajectoryRef;
80 
82  for(unsigned int i = 0; i < primaryTrajectoryCollection->size() ; i++) {
83 
84  TrajectoryRef trajRef( primaryTrajectoryCollection, i );
85 
87  TrajTrackAssociationCollection::const_iterator itPrimTrack = refMap.find( trajRef );
88  if( itPrimTrack == refMap.end() || (itPrimTrack->val).isNull() ) continue;
89  const reco::TrackRef& primary_track = itPrimTrack->val;
90 
92  TrajectoryToSeedsMap::const_iterator itSeeds = nuclMap.find( trajRef );
93  if( itSeeds == nuclMap.end() || (itSeeds->val).isNull()) continue;
94  const TrajectorySeedRefVector& seeds = itSeeds->val;
95 
97  std::vector<reco::TrackRef> secondary_tracks;
98  for( unsigned int k=0; k < secondaryTrackCollection->size(); k++) {
99  reco::TrackRef currentTrk(secondaryTrackCollection, k);
100  if( isInside( currentTrk, seeds ) ) secondary_tracks.push_back(currentTrk);
101  }
102 
104  vertexBuilder->build(primary_track, secondary_tracks);
105  likelihoodCalculator->calculate( vertexBuilder->getVertex() );
106 
108  reco::NuclearInteraction nuclInter(seeds, vertexBuilder->getVertex(), likelihoodCalculator->result() );
109 
112  if( additionalSecTracks.isValid() )
113  findAdditionalSecondaryTracks(nuclInter, additionalSecTracks);
114 
115  theNuclearInteractions->push_back( nuclInter );
116 
117  std::ostringstream str;
118  print(str, nuclInter, vertexBuilder);
119  edm::LogInfo("NuclearInteractionMaker") << str.str();
120 
121  }
122 
123 
124  LogDebug("NuclearInteractionMaker") << "End of NuclearInteractionMaker - Number of nuclear interactions found :" << theNuclearInteractions->size();
125  iEvent.put(theNuclearInteractions);
126 }
127 
128 // ------------ method called once each job just before starting event loop ------------
129 void
131 {
132 }
133 
135 
136 // ------ method used to check whether the seed of a track belong to the vector of seeds --
138  unsigned int seedKey = track->seedRef().key();
139  for (unsigned int i=0; i< seeds.size(); i++) { if( seeds[i].key() == seedKey ) return true; }
140  return false;
141 }
142 
144  const edm::Handle<reco::TrackCollection>& additionalSecTracks) {
145 
146  LogDebug("NuclearInteractionMaker") << "Check if one of the " << additionalSecTracks->size()
147  << " additional secondary track is compatible";
148  reco::TrackRefVector allSecondary;
149  for(unsigned int i=0; i< additionalSecTracks->size(); i++) {
150  reco::TrackRef sec(additionalSecTracks, i);
151  if( vertexBuilder->isCompatible( sec ) ) {
152  // check if sec is already a secondary track (with id = tkId)
153  // else add this new track
154  vertexBuilder->addSecondaryTrack( sec );
155  }
156  }
157  likelihoodCalculator->calculate( vertexBuilder->getVertex() );
158 
159  nucl = reco::NuclearInteraction(nucl.seeds(), vertexBuilder->getVertex(), likelihoodCalculator->result() );
160 }
161 
162 // -- print out
163 void print(std::ostringstream& out, const reco::NuclearInteraction& nucl, const std::auto_ptr< NuclearVertexBuilder >& builder) {
164  out<<"Nuclear Interaction with vertex position : (";
165  out<< nucl.vertex().position().x() << " , "
166  << nucl.vertex().position().y() << " , "
167  << nucl.vertex().position().z() << ")";
168  reco::TrackRef primTrack = (nucl.primaryTrack()).castTo<reco::TrackRef>();
169  out<<"\tLikelihood : " << nucl.likelihood() << std::endl;
170  out<<"\tPrimary Track : Pt = " << primTrack->pt() << " - Nhits = "
171  << primTrack->numberOfValidHits() << std::endl;
172  out << "\tNumber of seeds : " << nucl.seedsSize() << std::endl;
173  out << "\tNumber of secondary Tracks : " << nucl.secondaryTracksSize() << std::endl;
174  int it=0;
175  for( reco::NuclearInteraction::trackRef_iterator itr_=nucl.secondaryTracks_begin(); itr_ != nucl.secondaryTracks_end(); itr_++, it++) {
176  reco::TrackRef secTrack = (*itr_).castTo<reco::TrackRef>();
177  ClosestApproachInRPhi* theApproach = builder->closestApproach(primTrack, secTrack);
178  out << "\t\t Secondary track " << it << " : Pt = " << (*itr_)->pt()
179  << " - Nhits = " << (*itr_)->numberOfValidHits()
180  << " - Dist = " << theApproach->distance()
181  << " - chi2 = " << (*itr_)->normalizedChi2() << std::endl;
182  delete theApproach;
183  }
184  out << "----------------" << std::endl;
185 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &)
bool isInside(const reco::TrackRef &track, const TrajectorySeedRefVector &seeds)
trackRef_iterator secondaryTracks_begin() const
first iterator over secondary tracks
edm::ESWatcher< TransientTrackRecord > transientTrackWatcher_
NuclearInteractionEDProducer(const edm::ParameterSet &)
edm::EDGetTokenT< reco::TrackCollection > token_additionalSecTracks
const_iterator end() const
last iterator over the map (read only)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
const edm::RefToBase< reco::Track > & primaryTrack() const
return the base reference to the primary track
std::string print(const Track &, edm::Verbosity=edm::Concise)
Track print utility.
Definition: print.cc:10
const_iterator find(const key_type &k) const
find element with specified reference key
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const Point & position() const
position
Definition: Vertex.h:106
int secondaryTracksSize() const
return the number of secondary tracks
const TrajectorySeedRefVector & seeds()
return the seeds
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< TrajTrackAssociationCollection > token_refMapH
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
std::vector< NuclearInteraction > NuclearInteractionCollection
collection of NuclearInteractions
reco::Vertex::trackRef_iterator trackRef_iterator
edm::EDGetTokenT< TrajectoryCollection > token_primaryTrajectory
int seedsSize() const
return the number of seeds
bool isValid() const
Definition: HandleBase.h:76
double likelihood() const
return the likelihood ~ probability that the vertex is a real nuclear interaction ...
tuple out
Definition: dbtoconf.py:99
void findAdditionalSecondaryTracks(reco::NuclearInteraction &nucl, const edm::Handle< reco::TrackCollection > &additionalSecTracks)
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
edm::ESWatcher< IdealMagneticFieldRecord > magFieldWatcher_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
edm::EDGetTokenT< reco::TrackCollection > token_secondaryTrack
const reco::Vertex & vertex() const
return the vertex
list key
Definition: combine.py:13
std::auto_ptr< NuclearLikelihood > likelihoodCalculator
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
edm::EDGetTokenT< reco::TrackCollection > token_primaryTrack
edm::EDGetTokenT< TrajectoryToSeedsMap > token_nuclMapH
virtual float distance() const
trackRef_iterator secondaryTracks_end() const
last iterator over secondary tracks
std::auto_ptr< NuclearVertexBuilder > vertexBuilder