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

using AlgoProductCollection = TrackProducerAlgorithm< reco::Track >::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
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 

Public Member Functions

 NuclearTrackCorrector (const edm::ParameterSet &)
 
 ~NuclearTrackCorrector ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 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
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (std::string const &iProcessName, std::string const &iModuleLabel, bool iPrint, std::vector< char const * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) 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::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 72 of file NuclearTrackCorrector.h.

Member Typedef Documentation

Definition at line 80 of file NuclearTrackCorrector.h.

Definition at line 78 of file NuclearTrackCorrector.h.

Definition at line 77 of file NuclearTrackCorrector.h.

Definition at line 76 of file NuclearTrackCorrector.h.

Definition at line 75 of file NuclearTrackCorrector.h.

Constructor & Destructor Documentation

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

Definition at line 35 of file NuclearTrackCorrector.cc.

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

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

Definition at line 58 of file NuclearTrackCorrector.cc.

59 {
60 }

Member Function Documentation

void NuclearTrackCorrector::endJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 199 of file NuclearTrackCorrector.cc.

199  {
200 }
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 310 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().

313  {
314 
315  TrajectoryStateOnSurface theInitialStateForRefitting;
316  //the starting state is the state closest to the first hit along seedDirection.
317 
318  //avoiding to use transientTrack, it should be faster;
319  TrajectoryStateOnSurface innerStateFromTrack=trajectoryStateTransform::innerStateOnSurface(*theT,*theG,theMF);
320  TrajectoryStateOnSurface outerStateFromTrack=trajectoryStateTransform::outerStateOnSurface(*theT,*theG,theMF);
321  TrajectoryStateOnSurface initialStateFromTrack =
322  ( (innerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() <
323  (outerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() ) ?
324  innerStateFromTrack: outerStateFromTrack;
325 
326  // error is rescaled, but correlation are kept.
327  initialStateFromTrack.rescaleError(100);
328  theInitialStateForRefitting = TrajectoryStateOnSurface(initialStateFromTrack.localParameters(),
329  initialStateFromTrack.localError(),
330  initialStateFromTrack.surface(),
331  theMF);
332  return theInitialStateForRefitting;
333 }
const LocalTrajectoryParameters & localParameters() const
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
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 innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
reco::TrackExtra NuclearTrackCorrector::getNewTrackExtra ( const AlgoProductCollection algoresults)
private

get a new TrackExtra from an AlgoProductCollection

Definition at line 276 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().

276  {
277  Trajectory* theTraj = algoResults[0].trajectory;
278  PropagationDirection seedDir = algoResults[0].pDir;
279 
280  TrajectoryStateOnSurface outertsos;
281  TrajectoryStateOnSurface innertsos;
282  unsigned int innerId, outerId;
283  if (theTraj->direction() == alongMomentum) {
284  outertsos = theTraj->lastMeasurement().updatedState();
285  innertsos = theTraj->firstMeasurement().updatedState();
286  outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
287  innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
288  } else {
289  outertsos = theTraj->firstMeasurement().updatedState();
290  innertsos = theTraj->lastMeasurement().updatedState();
291  outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
292  innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
293  }
294 
295  GlobalPoint v = outertsos.globalParameters().position();
296  GlobalVector p = outertsos.globalParameters().momentum();
297  math::XYZVector outmom( p.x(), p.y(), p.z() );
298  math::XYZPoint outpos( v.x(), v.y(), v.z() );
299  v = innertsos.globalParameters().position();
300  p = innertsos.globalParameters().momentum();
301  math::XYZVector inmom( p.x(), p.y(), p.z() );
302  math::XYZPoint inpos( v.x(), v.y(), v.z() );
303 
304  return reco::TrackExtra (outpos, outmom, true, inpos, inmom, true,
305  outertsos.curvilinearError(), outerId,
306  innertsos.curvilinearError(), innerId, seedDir);
307 
308 }
ConstRecHitPointer const & recHit() const
const CurvilinearTrajectoryError & curvilinearError() const
T y() const
Definition: PV3DBase.h:63
PropagationDirection
PropagationDirection const & direction() const
Definition: Trajectory.cc:125
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:228
T z() const
Definition: PV3DBase.h:64
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:241
const GlobalTrajectoryParameters & globalParameters() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
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 236 of file NuclearTrackCorrector.cc.

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

Referenced by produce().

236  {
237 
238  const Trajectory* it = &newTraj;
239 
241  it->validRecHits( hits );
242 
243 
244  float ndof=0;
245  for(unsigned int h=0 ; h<hits.size() ; h++)
246  {
247  if( hits[h]->isValid() )
248  {
249  ndof = ndof + hits[h]->dimension() * hits[h]->weight();
250  }
251  else {
252  LogDebug("NuclearSeedGenerator") << " HIT IS INVALID ???";
253  }
254  }
255 
256 
257  ndof = ndof - 5;
258  reco::TrackRef theT = m_TrajToTrackCollection->operator[]( initialTrajRef );
259  LogDebug("NuclearSeedGenerator") << " TrackCorrector - number of valid hits" << hits.size() << "\n"
260  << " - number of hits from Track " << theT->recHitsSize() << "\n"
261  << " - number of valid hits from initial track " << theT->numberOfValidHits();
262 
263 
264  if( hits.size() > 1){
265 
266  TrajectoryStateOnSurface theInitialStateForRefitting = getInitialState(&(*theT),hits,theG.product(),theMF.product()
267 );
268 
269  reco::BeamSpot bs;
270  return theAlgo->buildTrack(theFitter.product(), thePropagator.product(), algoResults, hits, theInitialStateForRefitting ,it->seed(), ndof, bs, theT->seedRef());
271  }
272 
273  return false;
274 }
#define LogDebug(id)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:330
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:117
TrackProducerAlgorithm< reco::Track > * theAlgo
edm::ESHandle< Propagator > thePropagator
edm::ESHandle< TrajectoryFitter > theFitter
T const * product() const
Definition: ESHandle.h:86
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 202 of file NuclearTrackCorrector.cc.

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

Referenced by produce().

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

References edm::RefToBase< T >::castTo(), conf_, event(), plotBeamSpotDB::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::AssociationMap< Tag >::refProd(), edm::second(), str_Input_NuclearInteraction, str_Input_Trajectory, AlCaHLTBitMon_QueryRunRegistry::string, swap_map(), theFitter, theG, theMF, thePropagator, edm::helpers::KeyVal< K, V >::val, and verbosity.

64 {
65  // Create Output Collections
66  // --------------------------------------------------------------------------------------------------
67  std::auto_ptr<TrajectoryCollection> Output_traj ( new TrajectoryCollection );
68  std::auto_ptr<TrajectoryToTrajectoryMap> Output_trajmap ( new TrajectoryToTrajectoryMap );
69 
70  std::auto_ptr<reco::TrackExtraCollection> Output_trackextra ( new reco::TrackExtraCollection );
71  std::auto_ptr<reco::TrackCollection> Output_track ( new reco::TrackCollection );
72  std::auto_ptr<TrackToTrajectoryMap> Output_trackmap ( new TrackToTrajectoryMap );
73 
74 
75 
76 
77 
78 
79  // Load Reccord
80  // --------------------------------------------------------------------------------------------------
81  std::string fitterName = conf_.getParameter<std::string>("Fitter");
82  iSetup.get<TrajectoryFitter::Record>().get(fitterName,theFitter);
83 
84  std::string propagatorName = conf_.getParameter<std::string>("Propagator");
85  iSetup.get<TrackingComponentsRecord>().get(propagatorName,thePropagator);
86 
87  iSetup.get<TrackerDigiGeometryRecord>().get(theG);
88 
90 
91  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
92 
93  // Load Inputs
94  // --------------------------------------------------------------------------------------------------
95  edm::Handle< TrajectoryCollection > temp_m_TrajectoryCollection;
96  iEvent.getByLabel( str_Input_Trajectory.c_str(), temp_m_TrajectoryCollection );
97  const TrajectoryCollection m_TrajectoryCollection = *(temp_m_TrajectoryCollection.product());
98 
99  edm::Handle< NuclearInteractionCollection > temp_m_NuclearInteractionCollection;
100  iEvent.getByLabel( str_Input_NuclearInteraction.c_str(), temp_m_NuclearInteractionCollection );
101  const NuclearInteractionCollection m_NuclearInteractionCollection = *(temp_m_NuclearInteractionCollection.product());
102 
103  edm::Handle< TrajTrackAssociationCollection > h_TrajToTrackCollection;
104  iEvent.getByLabel( str_Input_Trajectory.c_str(), h_TrajToTrackCollection );
105  m_TrajToTrackCollection = h_TrajToTrackCollection.product();
106 
107 
108  // Correct the trajectories (Remove trajectory's hits that are located after the nuclear interacion)
109  // --------------------------------------------------------------------------------------------------
110  if(verbosity>=1){
111  LogDebug("NuclearTrackCorrector")
112  <<"Number of trajectories = "<<m_TrajectoryCollection.size() <<std::endl
113  <<"Number of nuclear interactions = "<<m_NuclearInteractionCollection.size();
114  }
115 
116  std::map<reco::TrackRef,TrajectoryRef> m_TrackToTrajMap;
117  swap_map(temp_m_TrajectoryCollection, m_TrackToTrajMap);
118 
119  for(unsigned int i = 0 ; i < m_NuclearInteractionCollection.size() ; i++)
120  {
121  reco::NuclearInteraction ni = m_NuclearInteractionCollection[i];
122  if( ni.likelihood()<0.4) continue;
123 
124  reco::TrackRef primTrackRef = ni.primaryTrack().castTo<reco::TrackRef>();
125 
126  TrajectoryRef trajRef = m_TrackToTrajMap[primTrackRef];
127 
128  Trajectory newTraj;
129  if( newTrajNeeded(newTraj, trajRef, ni) ) {
130 
131  AlgoProductCollection algoResults;
132  bool isOK = getTrackFromTrajectory( newTraj , trajRef, algoResults);
133 
134  if( isOK ) {
135 
136  pair<unsigned int, unsigned int> tempory_pair;
137  tempory_pair.first = Output_track->size();
138  tempory_pair.second = i;
139  Indice_Map.push_back(tempory_pair);
140 
141  reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, i );
142  reco::TrackExtra newTrackExtra = getNewTrackExtra(algoResults);
143  (algoResults[0].track)->setExtra( teref );
144 
145  Output_track->push_back(*algoResults[0].track);
146  Output_trackextra->push_back( newTrackExtra );
147  Output_traj->push_back(newTraj);
148 
149  }
150  }
151  else {
153  Output_track->push_back(*primTrackRef);
154  Output_trackextra->push_back( *primTrackRef->extra() );
155  Output_traj->push_back(*trajRef);
156  }
157  }
158 
159  }
160  const edm::OrphanHandle<TrajectoryCollection> Handle_traj = iEvent.put(Output_traj);
161  const edm::OrphanHandle<reco::TrackCollection> Handle_tracks = iEvent.put(Output_track);
162  iEvent.put(Output_trackextra);
163 
164  // Make Maps between elements
165  // --------------------------------------------------------------------------------------------------
166  if(Handle_tracks->size() != Handle_traj->size() )
167  {
168  printf("ERROR Handle_tracks->size() != Handle_traj->size() \n");
169  return;
170  }
171 
172  std::auto_ptr<TrackToTrackMap> Output_tracktrackmap ( new TrackToTrackMap(Handle_tracks, m_TrajToTrackCollection->refProd().val) );
173 
174  for(unsigned int i = 0 ; i < Indice_Map.size() ; i++)
175  {
176  TrajectoryRef InTrajRef ( temp_m_TrajectoryCollection, Indice_Map[i].second );
177  TrajectoryRef OutTrajRef ( Handle_traj, Indice_Map[i].first );
178  reco::TrackRef TrackRef ( Handle_tracks, Indice_Map[i].first );
179 
180  Output_trajmap ->insert(OutTrajRef,InTrajRef);
181  Output_trackmap->insert(TrackRef,InTrajRef);
182 
183  try{
184  reco::TrackRef PrimaryTrackRef = m_TrajToTrackCollection->operator[]( InTrajRef );
185  Output_tracktrackmap->insert(TrackRef,PrimaryTrackRef);
186  }catch(edm::Exception event){}
187 
188  }
189  iEvent.put(Output_trajmap);
190  iEvent.put(Output_trackmap);
191  iEvent.put(Output_tracktrackmap);
192 
193 
194  if(verbosity>=3)printf("-----------------------\n");
195 }
#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:14
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:17
U second(std::pair< T, U > const &p)
const ref_type & refProd() const
return ref-prod structure
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
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
edm::AssociationMap< edm::OneToOne< reco::TrackCollection, reco::TrackCollection > > TrackToTrackMap
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:418
RefProd< PROD > getRefBeforePut()
Definition: Event.h:141
edm::ESHandle< Propagator > thePropagator
edm::ESHandle< TrajectoryFitter > theFitter
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
T const * product() const
Definition: Handle.h:81
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
REF castTo() const
Definition: RefToBase.h:271
const T & get() const
Definition: EventSetup.h:56
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
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
void swap_map(const edm::Handle< TrajectoryCollection > &trajColl, std::map< reco::TrackRef, edm::Ref< TrajectoryCollection > > &result)
TrackProducerAlgorithm< reco::Track >::AlgoProductCollection AlgoProductCollection
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 335 of file NuclearTrackCorrector.cc.

References i, m_TrajToTrackCollection, and mps_fire::result.

Referenced by produce().

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

Member Data Documentation

edm::ParameterSet NuclearTrackCorrector::conf_
private

Definition at line 125 of file NuclearTrackCorrector.h.

Referenced by produce().

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

Definition at line 118 of file NuclearTrackCorrector.h.

Referenced by produce().

int NuclearTrackCorrector::int_Input_Hit_Distance
private

Definition at line 113 of file NuclearTrackCorrector.h.

int NuclearTrackCorrector::KeepOnlyCorrectedTracks
private

Definition at line 116 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

const TrajTrackAssociationCollection* NuclearTrackCorrector::m_TrajToTrackCollection
private

Definition at line 129 of file NuclearTrackCorrector.h.

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

std::string NuclearTrackCorrector::str_Input_NuclearInteraction
private

Definition at line 112 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

std::string NuclearTrackCorrector::str_Input_Trajectory
private

Definition at line 111 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

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

Definition at line 128 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and NuclearTrackCorrector().

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

Definition at line 123 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

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

Definition at line 121 of file NuclearTrackCorrector.h.

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

TransientInitialStateEstimator* NuclearTrackCorrector::theInitialState
private

Definition at line 126 of file NuclearTrackCorrector.h.

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

Definition at line 122 of file NuclearTrackCorrector.h.

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

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

Definition at line 124 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

int NuclearTrackCorrector::verbosity
private

Definition at line 115 of file NuclearTrackCorrector.h.

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