CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
NuclearTrackCorrector Class Reference

#include <RecoTracker/NuclearSeedGenerator/plugin/NuclearTrackCorrector.cc>

Inheritance diagram for NuclearTrackCorrector:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Types

typedef TrackProducerAlgorithm
< reco::Track >
::AlgoProductCollection 
AlgoProductCollection
 
typedef
TransientTrackingRecHit::ConstRecHitContainer 
ConstRecHitContainer
 
typedef edm::Ref
< TrackCandidateCollection
TrackCandidateRef
 
typedef edm::Ref
< TrajectoryCollection
TrajectoryRef
 
typedef edm::RefVector
< TrajectorySeedCollection
TrajectorySeedRefVector
 
- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 

Public Member Functions

 NuclearTrackCorrector (const edm::ParameterSet &)
 
 ~NuclearTrackCorrector ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void endJob ()
 
TrajectoryStateOnSurface getInitialState (const reco::Track *theT, TransientTrackingRecHit::RecHitContainer &hits, const TrackingGeometry *theG, const MagneticField *theMF)
 Calculate the inital state to be used to buil the track. More...
 
reco::TrackExtra getNewTrackExtra (const AlgoProductCollection &algoresults)
 get a new TrackExtra from an AlgoProductCollection More...
 
bool getTrackFromTrajectory (const Trajectory &newTraj, const TrajectoryRef &initialTrajRef, AlgoProductCollection &algoResults)
 Get the refitted track from the Trajectory. More...
 
bool newTrajNeeded (Trajectory &newtrajectory, const TrajectoryRef &trajRef, const reco::NuclearInteraction &ni)
 check if the trajectory has to be refitted and get the new trajectory More...
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 
void swap_map (const edm::Handle< TrajectoryCollection > &trajColl, std::map< reco::TrackRef, edm::Ref< TrajectoryCollection > > &result)
 

Private Attributes

edm::ParameterSet conf_
 
std::vector< std::pair
< unsigned int, unsigned int > > 
Indice_Map
 
int int_Input_Hit_Distance
 
int KeepOnlyCorrectedTracks
 
const
TrajTrackAssociationCollection
m_TrajToTrackCollection
 
std::string str_Input_NuclearInteraction
 
std::string str_Input_Trajectory
 
TrackProducerAlgorithm
< reco::Track > * 
theAlgo
 
edm::ESHandle< TrajectoryFittertheFitter
 
edm::ESHandle< TrackerGeometrytheG
 
TransientInitialStateEstimatortheInitialState
 
edm::ESHandle< MagneticFieldtheMF
 
edm::ESHandle< PropagatorthePropagator
 
int verbosity
 

Additional Inherited Members

- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 75 of file NuclearTrackCorrector.h.

Member Typedef Documentation

Definition at line 82 of file NuclearTrackCorrector.h.

Definition at line 81 of file NuclearTrackCorrector.h.

Definition at line 80 of file NuclearTrackCorrector.h.

Definition at line 79 of file NuclearTrackCorrector.h.

Definition at line 78 of file NuclearTrackCorrector.h.

Constructor & Destructor Documentation

NuclearTrackCorrector::NuclearTrackCorrector ( const edm::ParameterSet iConfig)
explicit

Definition at line 36 of file NuclearTrackCorrector.cc.

References edm::ParameterSet::getParameter(), KeepOnlyCorrectedTracks, str_Input_NuclearInteraction, str_Input_Trajectory, AlCaHLTBitMon_QueryRunRegistry::string, theAlgo, and verbosity.

36  :
37 conf_(iConfig),
39 {
40  str_Input_Trajectory = iConfig.getParameter<std::string> ("InputTrajectory");
41  str_Input_NuclearInteraction = iConfig.getParameter<std::string> ("InputNuclearInteraction");
42  verbosity = iConfig.getParameter<int> ("Verbosity");
43  KeepOnlyCorrectedTracks = iConfig.getParameter<bool> ("KeepOnlyCorrectedTracks");
44 
45 
47 
48  produces< TrajectoryCollection >();
49  produces< TrajectoryToTrajectoryMap >();
50 
51  produces< reco::TrackExtraCollection >();
52  produces< reco::TrackCollection >();
53  produces< TrackToTrajectoryMap >();
54 
55  produces< TrackToTrackMap >();
56 }
T getParameter(std::string const &) const
TransientInitialStateEstimator * theInitialState
TrackProducerAlgorithm< reco::Track > * theAlgo
std::string str_Input_NuclearInteraction
NuclearTrackCorrector::~NuclearTrackCorrector ( )

Definition at line 59 of file NuclearTrackCorrector.cc.

60 {
61 }

Member Function Documentation

void NuclearTrackCorrector::endJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 201 of file NuclearTrackCorrector.cc.

201  {
202 }
TrajectoryStateOnSurface NuclearTrackCorrector::getInitialState ( const reco::Track theT,
TransientTrackingRecHit::RecHitContainer hits,
const TrackingGeometry theG,
const MagneticField theMF 
)
private

Calculate the inital state to be used to buil the track.

Definition at line 312 of file NuclearTrackCorrector.cc.

References TrajectoryStateOnSurface::globalPosition(), trajectoryStateTransform::innerStateOnSurface(), TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), mag2(), trajectoryStateTransform::outerStateOnSurface(), TrajectoryStateOnSurface::rescaleError(), TrajectoryStateOnSurface::surface(), and theMF.

Referenced by getTrackFromTrajectory().

315  {
316 
317  TrajectoryStateOnSurface theInitialStateForRefitting;
318  //the starting state is the state closest to the first hit along seedDirection.
319 
320  //avoiding to use transientTrack, it should be faster;
321  TrajectoryStateOnSurface innerStateFromTrack=trajectoryStateTransform::innerStateOnSurface(*theT,*theG,theMF);
322  TrajectoryStateOnSurface outerStateFromTrack=trajectoryStateTransform::outerStateOnSurface(*theT,*theG,theMF);
323  TrajectoryStateOnSurface initialStateFromTrack =
324  ( (innerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() <
325  (outerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() ) ?
326  innerStateFromTrack: outerStateFromTrack;
327 
328  // error is rescaled, but correlation are kept.
329  initialStateFromTrack.rescaleError(100);
330  theInitialStateForRefitting = TrajectoryStateOnSurface(initialStateFromTrack.localParameters(),
331  initialStateFromTrack.localError(),
332  initialStateFromTrack.surface(),
333  theMF);
334  return theInitialStateForRefitting;
335 }
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
const LocalTrajectoryParameters & localParameters() const
edm::ESHandle< MagneticField > theMF
GlobalPoint globalPosition() const
const SurfaceType & surface() const
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
const LocalTrajectoryError & localError() const
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
reco::TrackExtra NuclearTrackCorrector::getNewTrackExtra ( const AlgoProductCollection algoresults)
private

get a new TrackExtra from an AlgoProductCollection

Definition at line 278 of file NuclearTrackCorrector.cc.

References alongMomentum, TrajectoryStateOnSurface::curvilinearError(), Trajectory::direction(), Trajectory::firstMeasurement(), TrajectoryStateOnSurface::globalParameters(), Trajectory::lastMeasurement(), GlobalTrajectoryParameters::momentum(), AlCaHLTBitMon_ParallelJobs::p, GlobalTrajectoryParameters::position(), TrajectoryMeasurement::recHit(), TrajectoryMeasurement::updatedState(), findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by produce().

278  {
279  Trajectory* theTraj = algoResults[0].first;
280  PropagationDirection seedDir = algoResults[0].second.second;
281 
282  TrajectoryStateOnSurface outertsos;
283  TrajectoryStateOnSurface innertsos;
284  unsigned int innerId, outerId;
285  if (theTraj->direction() == alongMomentum) {
286  outertsos = theTraj->lastMeasurement().updatedState();
287  innertsos = theTraj->firstMeasurement().updatedState();
288  outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
289  innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
290  } else {
291  outertsos = theTraj->firstMeasurement().updatedState();
292  innertsos = theTraj->lastMeasurement().updatedState();
293  outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
294  innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
295  }
296 
297  GlobalPoint v = outertsos.globalParameters().position();
298  GlobalVector p = outertsos.globalParameters().momentum();
299  math::XYZVector outmom( p.x(), p.y(), p.z() );
300  math::XYZPoint outpos( v.x(), v.y(), v.z() );
301  v = innertsos.globalParameters().position();
302  p = innertsos.globalParameters().momentum();
303  math::XYZVector inmom( p.x(), p.y(), p.z() );
304  math::XYZPoint inpos( v.x(), v.y(), v.z() );
305 
306  return reco::TrackExtra (outpos, outmom, true, inpos, inmom, true,
307  outertsos.curvilinearError(), outerId,
308  innertsos.curvilinearError(), innerId, seedDir);
309 
310 }
ConstRecHitPointer const & recHit() const
const CurvilinearTrajectoryError & curvilinearError() const
T y() const
Definition: PV3DBase.h:63
PropagationDirection
PropagationDirection const & direction() const
Definition: Trajectory.cc:196
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:193
T z() const
Definition: PV3DBase.h:64
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:206
const GlobalTrajectoryParameters & globalParameters() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
TrajectoryStateOnSurface const & updatedState() const
T x() const
Definition: PV3DBase.h:62
bool NuclearTrackCorrector::getTrackFromTrajectory ( const Trajectory newTraj,
const TrajectoryRef initialTrajRef,
AlgoProductCollection algoResults 
)
private

Get the refitted track from the Trajectory.

Definition at line 238 of file NuclearTrackCorrector.cc.

References TrackProducerAlgorithm< T >::buildTrack(), getInitialState(), h, LogDebug, m_TrajToTrackCollection, edm::ESHandle< class >::product(), Trajectory::seed(), theAlgo, theFitter, theG, theMF, thePropagator, and Trajectory::validRecHits().

Referenced by produce().

238  {
239 
240  const Trajectory* it = &newTraj;
241 
243  it->validRecHits( hits );
244 
245 
246  float ndof=0;
247  for(unsigned int h=0 ; h<hits.size() ; h++)
248  {
249  if( hits[h]->isValid() )
250  {
251  ndof = ndof + hits[h]->dimension() * hits[h]->weight();
252  }
253  else {
254  LogDebug("NuclearSeedGenerator") << " HIT IS INVALID ???";
255  }
256  }
257 
258 
259  ndof = ndof - 5;
260  reco::TrackRef theT = m_TrajToTrackCollection->operator[]( initialTrajRef );
261  LogDebug("NuclearSeedGenerator") << " TrackCorrector - number of valid hits" << hits.size() << "\n"
262  << " - number of hits from Track " << theT->recHitsSize() << "\n"
263  << " - number of valid hits from initial track " << theT->numberOfValidHits();
264 
265 
266  if( hits.size() > 1){
267 
268  TrajectoryStateOnSurface theInitialStateForRefitting = getInitialState(&(*theT),hits,theG.product(),theMF.product()
269 );
270 
271  reco::BeamSpot bs;
272  return theAlgo->buildTrack(theFitter.product(), thePropagator.product(), algoResults, hits, theInitialStateForRefitting ,it->seed(), ndof, bs, theT->seedRef());
273  }
274 
275  return false;
276 }
#define LogDebug(id)
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:277
edm::ESHandle< MagneticField > theMF
std::vector< ConstRecHitPointer > RecHitContainer
bool buildTrack(const TrajectoryFitter *, const Propagator *, AlgoProductCollection &, TransientTrackingRecHit::RecHitContainer &, TrajectoryStateOnSurface &, const TrajectorySeed &, float, const reco::BeamSpot &, SeedRef seedRef=SeedRef(), int qualityMask=0, signed char nLoops=0)
Construct Tracks to be put in the event.
edm::ESHandle< TrackerGeometry > theG
void validRecHits(ConstRecHitContainer &cont) const
Definition: Trajectory.cc:188
TrackProducerAlgorithm< reco::Track > * theAlgo
edm::ESHandle< Propagator > thePropagator
edm::ESHandle< TrajectoryFitter > theFitter
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
T const * product() const
Definition: ESHandle.h:62
const TrajTrackAssociationCollection * m_TrajToTrackCollection
TrajectoryStateOnSurface getInitialState(const reco::Track *theT, TransientTrackingRecHit::RecHitContainer &hits, const TrackingGeometry *theG, const MagneticField *theMF)
Calculate the inital state to be used to buil the track.
bool NuclearTrackCorrector::newTrajNeeded ( Trajectory newtrajectory,
const TrajectoryRef trajRef,
const reco::NuclearInteraction ni 
)
private

check if the trajectory has to be refitted and get the new trajectory

Definition at line 204 of file NuclearTrackCorrector.cc.

References alongMomentum, LogDebug, m, PV3DBase< T, PVType, FrameType >::mag(), reco::Vertex::position(), Trajectory::push(), mathSSE::sqrt(), theG, verbosity, and reco::NuclearInteraction::vertex().

Referenced by produce().

204  {
205 
206  bool needNewTraj=false;
207  reco::Vertex::Point vtx_pos = ni.vertex().position();
208  double vtx_pos_mag = sqrt (vtx_pos.X()*vtx_pos.X()+vtx_pos.Y()*vtx_pos.Y()+vtx_pos.Z()*vtx_pos.Z());
209  if(verbosity>=2) printf("Nuclear Interaction pos = %f\n",vtx_pos_mag );
210 
211 
212  newtrajectory = Trajectory(trajRef->seed(), alongMomentum);
213 
214  // Look all the Hits of the trajectory and keep only Hits before seeds
215  Trajectory::DataContainer Measurements = trajRef->measurements();
216  if(verbosity>=2)LogDebug("NuclearTrackCorrector")<<"Size of Measurements = "<<Measurements.size();
217 
218  for(unsigned int m=Measurements.size()-1 ;m!=(unsigned int)-1 ; m--){
219 
220  if(!Measurements[m].recHit()->isValid() )continue;
221  GlobalPoint hit_pos = theG->idToDet(Measurements[m].recHit()->geographicalId())->surface().toGlobal(Measurements[m].recHit()->localPosition());
222 
223  if(verbosity>=2)printf("Hit pos = %f",hit_pos.mag() );
224 
225  if(hit_pos.mag()>vtx_pos_mag){
226  if(verbosity>=2)printf(" X ");
227  needNewTraj=true;
228  }else{
229  newtrajectory.push(Measurements[m]);
230  }
231  if(verbosity>=2)printf("\n");
232  }
233 
234  return needNewTraj;
235 }
#define LogDebug(id)
const Point & position() const
position
Definition: Vertex.h:93
DataContainer const & measurements() const
Definition: Trajectory.h:215
std::vector< TrajectoryMeasurement > DataContainer
Definition: Trajectory.h:42
T mag() const
Definition: PV3DBase.h:67
edm::ESHandle< TrackerGeometry > theG
T sqrt(T t)
Definition: SSEVec.h:48
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
const reco::Vertex & vertex() const
return the vertex
void push(const TrajectoryMeasurement &tm)
Definition: Trajectory.cc:35
void NuclearTrackCorrector::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 64 of file NuclearTrackCorrector.cc.

References edm::RefToBase< T >::castTo(), conf_, event(), first, edm::EventSetup::get(), edm::Event::getByLabel(), getNewTrackExtra(), edm::ParameterSet::getParameter(), edm::Event::getRefBeforePut(), getTrackFromTrajectory(), i, Indice_Map, edm::AssociationMap< Tag >::insert(), KeepOnlyCorrectedTracks, reco::NuclearInteraction::likelihood(), LogDebug, m_TrajToTrackCollection, newTrajNeeded(), reco::NuclearInteraction::primaryTrack(), edm::Handle< T >::product(), edm::Event::put(), edm::second(), str_Input_NuclearInteraction, str_Input_Trajectory, AlCaHLTBitMon_QueryRunRegistry::string, swap_map(), theFitter, theG, theMF, thePropagator, and verbosity.

65 {
66  // Create Output Collections
67  // --------------------------------------------------------------------------------------------------
68  std::auto_ptr<TrajectoryCollection> Output_traj ( new TrajectoryCollection );
69  std::auto_ptr<TrajectoryToTrajectoryMap> Output_trajmap ( new TrajectoryToTrajectoryMap );
70 
71  std::auto_ptr<reco::TrackExtraCollection> Output_trackextra ( new reco::TrackExtraCollection );
72  std::auto_ptr<reco::TrackCollection> Output_track ( new reco::TrackCollection );
73  std::auto_ptr<TrackToTrajectoryMap> Output_trackmap ( new TrackToTrajectoryMap );
74 
75  std::auto_ptr<TrackToTrackMap> Output_tracktrackmap ( new TrackToTrackMap );
76 
77 
78 
79 
80 
81  // Load Reccord
82  // --------------------------------------------------------------------------------------------------
83  std::string fitterName = conf_.getParameter<std::string>("Fitter");
84  iSetup.get<TrajectoryFitter::Record>().get(fitterName,theFitter);
85 
86  std::string propagatorName = conf_.getParameter<std::string>("Propagator");
87  iSetup.get<TrackingComponentsRecord>().get(propagatorName,thePropagator);
88 
89  iSetup.get<TrackerDigiGeometryRecord>().get(theG);
90 
92 
93  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
94 
95  // Load Inputs
96  // --------------------------------------------------------------------------------------------------
97  edm::Handle< TrajectoryCollection > temp_m_TrajectoryCollection;
98  iEvent.getByLabel( str_Input_Trajectory.c_str(), temp_m_TrajectoryCollection );
99  const TrajectoryCollection m_TrajectoryCollection = *(temp_m_TrajectoryCollection.product());
100 
101  edm::Handle< NuclearInteractionCollection > temp_m_NuclearInteractionCollection;
102  iEvent.getByLabel( str_Input_NuclearInteraction.c_str(), temp_m_NuclearInteractionCollection );
103  const NuclearInteractionCollection m_NuclearInteractionCollection = *(temp_m_NuclearInteractionCollection.product());
104 
105  edm::Handle< TrajTrackAssociationCollection > h_TrajToTrackCollection;
106  iEvent.getByLabel( str_Input_Trajectory.c_str(), h_TrajToTrackCollection );
107  m_TrajToTrackCollection = h_TrajToTrackCollection.product();
108 
109 
110  // Correct the trajectories (Remove trajectory's hits that are located after the nuclear interacion)
111  // --------------------------------------------------------------------------------------------------
112  if(verbosity>=1){
113  LogDebug("NuclearTrackCorrector")
114  <<"Number of trajectories = "<<m_TrajectoryCollection.size() <<std::endl
115  <<"Number of nuclear interactions = "<<m_NuclearInteractionCollection.size();
116  }
117 
118  std::map<reco::TrackRef,TrajectoryRef> m_TrackToTrajMap;
119  swap_map(temp_m_TrajectoryCollection, m_TrackToTrajMap);
120 
121  for(unsigned int i = 0 ; i < m_NuclearInteractionCollection.size() ; i++)
122  {
123  reco::NuclearInteraction ni = m_NuclearInteractionCollection[i];
124  if( ni.likelihood()<0.4) continue;
125 
126  reco::TrackRef primTrackRef = ni.primaryTrack().castTo<reco::TrackRef>();
127 
128  TrajectoryRef trajRef = m_TrackToTrajMap[primTrackRef];
129 
130  Trajectory newTraj;
131  if( newTrajNeeded(newTraj, trajRef, ni) ) {
132 
133  AlgoProductCollection algoResults;
134  bool isOK = getTrackFromTrajectory( newTraj , trajRef, algoResults);
135 
136  if( isOK ) {
137 
138  pair<unsigned int, unsigned int> tempory_pair;
139  tempory_pair.first = Output_track->size();
140  tempory_pair.second = i;
141  Indice_Map.push_back(tempory_pair);
142 
143  reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, i );
144  reco::TrackExtra newTrackExtra = getNewTrackExtra(algoResults);
145  (algoResults[0].second.first)->setExtra( teref );
146 
147  Output_track->push_back(*algoResults[0].second.first);
148  Output_trackextra->push_back( newTrackExtra );
149  Output_traj->push_back(newTraj);
150 
151  }
152  }
153  else {
155  Output_track->push_back(*primTrackRef);
156  Output_trackextra->push_back( *primTrackRef->extra() );
157  Output_traj->push_back(*trajRef);
158  }
159  }
160 
161  }
162  const edm::OrphanHandle<TrajectoryCollection> Handle_traj = iEvent.put(Output_traj);
163  const edm::OrphanHandle<reco::TrackCollection> Handle_tracks = iEvent.put(Output_track);
164  iEvent.put(Output_trackextra);
165 
166  // Make Maps between elements
167  // --------------------------------------------------------------------------------------------------
168  if(Handle_tracks->size() != Handle_traj->size() )
169  {
170  printf("ERROR Handle_tracks->size() != Handle_traj->size() \n");
171  return;
172  }
173 
174 
175 
176  for(unsigned int i = 0 ; i < Indice_Map.size() ; i++)
177  {
178  TrajectoryRef InTrajRef ( temp_m_TrajectoryCollection, Indice_Map[i].second );
179  TrajectoryRef OutTrajRef ( Handle_traj, Indice_Map[i].first );
180  reco::TrackRef TrackRef ( Handle_tracks, Indice_Map[i].first );
181 
182  Output_trajmap ->insert(OutTrajRef,InTrajRef);
183  Output_trackmap->insert(TrackRef,InTrajRef);
184 
185  try{
186  reco::TrackRef PrimaryTrackRef = m_TrajToTrackCollection->operator[]( InTrajRef );
187  Output_tracktrackmap->insert(TrackRef,PrimaryTrackRef);
188  }catch(edm::Exception event){}
189 
190  }
191  iEvent.put(Output_trajmap);
192  iEvent.put(Output_trackmap);
193  iEvent.put(Output_tracktrackmap);
194 
195 
196  if(verbosity>=3)printf("-----------------------\n");
197 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
edm::ESHandle< MagneticField > theMF
const edm::RefToBase< reco::Track > & primaryTrack() const
return the base reference to the primary track
std::vector< std::pair< unsigned int, unsigned int > > Indice_Map
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool getTrackFromTrajectory(const Trajectory &newTraj, const TrajectoryRef &initialTrajRef, AlgoProductCollection &algoResults)
Get the refitted track from the Trajectory.
edm::Ref< TrackExtraCollection > TrackExtraRef
persistent reference to a TrackExtra
Definition: TrackExtraFwd.h:13
U second(std::pair< T, U > const &p)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
edm::ESHandle< TrackerGeometry > theG
std::vector< NuclearInteraction > NuclearInteractionCollection
collection of NuclearInteractions
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool first
Definition: L1TdeRCT.cc:94
double likelihood() const
return the likelihood ~ probability that the vertex is a real nuclear interaction ...
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
RefProd< PROD > getRefBeforePut()
Definition: Event.h:106
edm::ESHandle< Propagator > thePropagator
edm::ESHandle< TrajectoryFitter > theFitter
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
void insert(const key_type &k, const data_type &v)
insert an association
reco::TrackExtra getNewTrackExtra(const AlgoProductCollection &algoresults)
get a new TrackExtra from an AlgoProductCollection
const T & get() const
Definition: EventSetup.h:55
REF castTo() const
cast to a concrete type
Definition: RefToBase.h:241
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
T const * product() const
Definition: Handle.h:74
edm::Ref< TrajectoryCollection > TrajectoryRef
bool newTrajNeeded(Trajectory &newtrajectory, const TrajectoryRef &trajRef, const reco::NuclearInteraction &ni)
check if the trajectory has to be refitted and get the new trajectory
std::vector< Trajectory > TrajectoryCollection
const TrajTrackAssociationCollection * m_TrajToTrackCollection
TrackProducerAlgorithm< reco::Track >::AlgoProductCollection AlgoProductCollection
void swap_map(const edm::Handle< TrajectoryCollection > &trajColl, std::map< reco::TrackRef, edm::Ref< TrajectoryCollection > > &result)
std::string str_Input_NuclearInteraction
void NuclearTrackCorrector::swap_map ( const edm::Handle< TrajectoryCollection > &  trajColl,
std::map< reco::TrackRef, edm::Ref< TrajectoryCollection > > &  result 
)
private

Definition at line 337 of file NuclearTrackCorrector.cc.

References i, m_TrajToTrackCollection, and query::result.

Referenced by produce().

337  {
338  for(unsigned int i = 0 ; i < trajColl->size() ; i++)
339  {
340  TrajectoryRef InTrajRef ( trajColl, i);
341  reco::TrackRef PrimaryTrackRef = m_TrajToTrackCollection->operator[]( InTrajRef );
342  result[ PrimaryTrackRef ] = InTrajRef;
343  }
344 }
int i
Definition: DBlmapReader.cc:9
tuple result
Definition: query.py:137
edm::Ref< TrajectoryCollection > TrajectoryRef
const TrajTrackAssociationCollection * m_TrajToTrackCollection

Member Data Documentation

edm::ParameterSet NuclearTrackCorrector::conf_
private

Definition at line 127 of file NuclearTrackCorrector.h.

Referenced by produce().

std::vector< std::pair<unsigned int, unsigned int> > NuclearTrackCorrector::Indice_Map
private

Definition at line 120 of file NuclearTrackCorrector.h.

Referenced by produce().

int NuclearTrackCorrector::int_Input_Hit_Distance
private

Definition at line 115 of file NuclearTrackCorrector.h.

int NuclearTrackCorrector::KeepOnlyCorrectedTracks
private

Definition at line 118 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

const TrajTrackAssociationCollection* NuclearTrackCorrector::m_TrajToTrackCollection
private

Definition at line 131 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), produce(), and swap_map().

std::string NuclearTrackCorrector::str_Input_NuclearInteraction
private

Definition at line 114 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

std::string NuclearTrackCorrector::str_Input_Trajectory
private

Definition at line 113 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

TrackProducerAlgorithm<reco::Track>* NuclearTrackCorrector::theAlgo
private

Definition at line 130 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and NuclearTrackCorrector().

edm::ESHandle<TrajectoryFitter> NuclearTrackCorrector::theFitter
private

Definition at line 125 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

edm::ESHandle<TrackerGeometry> NuclearTrackCorrector::theG
private

Definition at line 123 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), newTrajNeeded(), and produce().

TransientInitialStateEstimator* NuclearTrackCorrector::theInitialState
private

Definition at line 128 of file NuclearTrackCorrector.h.

edm::ESHandle<MagneticField> NuclearTrackCorrector::theMF
private

Definition at line 124 of file NuclearTrackCorrector.h.

Referenced by getInitialState(), getTrackFromTrajectory(), and produce().

edm::ESHandle<Propagator> NuclearTrackCorrector::thePropagator
private

Definition at line 126 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

int NuclearTrackCorrector::verbosity
private

Definition at line 117 of file NuclearTrackCorrector.h.

Referenced by newTrajNeeded(), NuclearTrackCorrector(), and produce().