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
 
- 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 ()
 
ModuleDescription const & moduleDescription () const
 
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
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 (const std::string &iProcessName, std::vector< const char * > &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::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 74 of file NuclearTrackCorrector.h.

Member Typedef Documentation

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.

Definition at line 77 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 200 of file NuclearTrackCorrector.cc.

200  {
201 }
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 311 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().

314  {
315 
316  TrajectoryStateOnSurface theInitialStateForRefitting;
317  //the starting state is the state closest to the first hit along seedDirection.
318 
319  //avoiding to use transientTrack, it should be faster;
320  TrajectoryStateOnSurface innerStateFromTrack=trajectoryStateTransform::innerStateOnSurface(*theT,*theG,theMF);
321  TrajectoryStateOnSurface outerStateFromTrack=trajectoryStateTransform::outerStateOnSurface(*theT,*theG,theMF);
322  TrajectoryStateOnSurface initialStateFromTrack =
323  ( (innerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() <
324  (outerStateFromTrack.globalPosition()-hits.front()->globalPosition()).mag2() ) ?
325  innerStateFromTrack: outerStateFromTrack;
326 
327  // error is rescaled, but correlation are kept.
328  initialStateFromTrack.rescaleError(100);
329  theInitialStateForRefitting = TrajectoryStateOnSurface(initialStateFromTrack.localParameters(),
330  initialStateFromTrack.localError(),
331  initialStateFromTrack.surface(),
332  theMF);
333  return theInitialStateForRefitting;
334 }
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 277 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().

277  {
278  Trajectory* theTraj = algoResults[0].first;
279  PropagationDirection seedDir = algoResults[0].second.second;
280 
281  TrajectoryStateOnSurface outertsos;
282  TrajectoryStateOnSurface innertsos;
283  unsigned int innerId, outerId;
284  if (theTraj->direction() == alongMomentum) {
285  outertsos = theTraj->lastMeasurement().updatedState();
286  innertsos = theTraj->firstMeasurement().updatedState();
287  outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
288  innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
289  } else {
290  outertsos = theTraj->firstMeasurement().updatedState();
291  innertsos = theTraj->lastMeasurement().updatedState();
292  outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
293  innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
294  }
295 
296  GlobalPoint v = outertsos.globalParameters().position();
297  GlobalVector p = outertsos.globalParameters().momentum();
298  math::XYZVector outmom( p.x(), p.y(), p.z() );
299  math::XYZPoint outpos( v.x(), v.y(), v.z() );
300  v = innertsos.globalParameters().position();
301  p = innertsos.globalParameters().momentum();
302  math::XYZVector inmom( p.x(), p.y(), p.z() );
303  math::XYZPoint inpos( v.x(), v.y(), v.z() );
304 
305  return reco::TrackExtra (outpos, outmom, true, inpos, inmom, true,
306  outertsos.curvilinearError(), outerId,
307  innertsos.curvilinearError(), innerId, seedDir);
308 
309 }
ConstRecHitPointer const & recHit() const
const CurvilinearTrajectoryError & curvilinearError() const
T y() const
Definition: PV3DBase.h:63
PropagationDirection
PropagationDirection const & direction() const
Definition: Trajectory.cc:118
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:181
T z() const
Definition: PV3DBase.h:64
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:194
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 237 of file NuclearTrackCorrector.cc.

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

Referenced by produce().

237  {
238 
239  const Trajectory* it = &newTraj;
240 
242  it->validRecHits( hits );
243 
244 
245  float ndof=0;
246  for(unsigned int h=0 ; h<hits.size() ; h++)
247  {
248  if( hits[h]->isValid() )
249  {
250  ndof = ndof + hits[h]->dimension() * hits[h]->weight();
251  }
252  else {
253  LogDebug("NuclearSeedGenerator") << " HIT IS INVALID ???";
254  }
255  }
256 
257 
258  ndof = ndof - 5;
259  reco::TrackRef theT = m_TrajToTrackCollection->operator[]( initialTrajRef );
260  LogDebug("NuclearSeedGenerator") << " TrackCorrector - number of valid hits" << hits.size() << "\n"
261  << " - number of hits from Track " << theT->recHitsSize() << "\n"
262  << " - number of valid hits from initial track " << theT->numberOfValidHits();
263 
264 
265  if( hits.size() > 1){
266 
267  TrajectoryStateOnSurface theInitialStateForRefitting = getInitialState(&(*theT),hits,theG.product(),theMF.product()
268 );
269 
270  reco::BeamSpot bs;
271  return theAlgo->buildTrack(theFitter.product(), thePropagator.product(), algoResults, hits, theInitialStateForRefitting ,it->seed(), ndof, bs, theT->seedRef());
272  }
273 
274  return false;
275 }
#define LogDebug(id)
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:275
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
TrackProducerAlgorithm< reco::Track > * theAlgo
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
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 203 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().

203  {
204 
205  bool needNewTraj=false;
206  reco::Vertex::Point vtx_pos = ni.vertex().position();
207  double vtx_pos_mag = sqrt (vtx_pos.X()*vtx_pos.X()+vtx_pos.Y()*vtx_pos.Y()+vtx_pos.Z()*vtx_pos.Z());
208  if(verbosity>=2) printf("Nuclear Interaction pos = %f\n",vtx_pos_mag );
209 
210 
211  newtrajectory = Trajectory(trajRef->seed(), alongMomentum);
212 
213  // Look all the Hits of the trajectory and keep only Hits before seeds
214  Trajectory::DataContainer Measurements = trajRef->measurements();
215  if(verbosity>=2)LogDebug("NuclearTrackCorrector")<<"Size of Measurements = "<<Measurements.size();
216 
217  for(unsigned int m=Measurements.size()-1 ;m!=(unsigned int)-1 ; m--){
218 
219  if(!Measurements[m].recHit()->isValid() )continue;
220  GlobalPoint hit_pos = theG->idToDet(Measurements[m].recHit()->geographicalId())->surface().toGlobal(Measurements[m].recHit()->localPosition());
221 
222  if(verbosity>=2)printf("Hit pos = %f",hit_pos.mag() );
223 
224  if(hit_pos.mag()>vtx_pos_mag){
225  if(verbosity>=2)printf(" X ");
226  needNewTraj=true;
227  }else{
228  newtrajectory.push(Measurements[m]);
229  }
230  if(verbosity>=2)printf("\n");
231  }
232 
233  return needNewTraj;
234 }
#define LogDebug(id)
const Point & position() const
position
Definition: Vertex.h:106
DataContainer const & measurements() const
Definition: Trajectory.h:203
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:39
const reco::Vertex & vertex() const
return the vertex
void push(const TrajectoryMeasurement &tm)
Definition: Trajectory.cc:30
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::second(), str_Input_NuclearInteraction, str_Input_Trajectory, AlCaHLTBitMon_QueryRunRegistry::string, swap_map(), theFitter, theG, theMF, thePropagator, 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  std::auto_ptr<TrackToTrackMap> Output_tracktrackmap ( new TrackToTrackMap );
75 
76 
77 
78 
79 
80  // Load Reccord
81  // --------------------------------------------------------------------------------------------------
82  std::string fitterName = conf_.getParameter<std::string>("Fitter");
83  iSetup.get<TrajectoryFitter::Record>().get(fitterName,theFitter);
84 
85  std::string propagatorName = conf_.getParameter<std::string>("Propagator");
86  iSetup.get<TrackingComponentsRecord>().get(propagatorName,thePropagator);
87 
88  iSetup.get<TrackerDigiGeometryRecord>().get(theG);
89 
91 
92  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
93 
94  // Load Inputs
95  // --------------------------------------------------------------------------------------------------
96  edm::Handle< TrajectoryCollection > temp_m_TrajectoryCollection;
97  iEvent.getByLabel( str_Input_Trajectory.c_str(), temp_m_TrajectoryCollection );
98  const TrajectoryCollection m_TrajectoryCollection = *(temp_m_TrajectoryCollection.product());
99 
100  edm::Handle< NuclearInteractionCollection > temp_m_NuclearInteractionCollection;
101  iEvent.getByLabel( str_Input_NuclearInteraction.c_str(), temp_m_NuclearInteractionCollection );
102  const NuclearInteractionCollection m_NuclearInteractionCollection = *(temp_m_NuclearInteractionCollection.product());
103 
104  edm::Handle< TrajTrackAssociationCollection > h_TrajToTrackCollection;
105  iEvent.getByLabel( str_Input_Trajectory.c_str(), h_TrajToTrackCollection );
106  m_TrajToTrackCollection = h_TrajToTrackCollection.product();
107 
108 
109  // Correct the trajectories (Remove trajectory's hits that are located after the nuclear interacion)
110  // --------------------------------------------------------------------------------------------------
111  if(verbosity>=1){
112  LogDebug("NuclearTrackCorrector")
113  <<"Number of trajectories = "<<m_TrajectoryCollection.size() <<std::endl
114  <<"Number of nuclear interactions = "<<m_NuclearInteractionCollection.size();
115  }
116 
117  std::map<reco::TrackRef,TrajectoryRef> m_TrackToTrajMap;
118  swap_map(temp_m_TrajectoryCollection, m_TrackToTrajMap);
119 
120  for(unsigned int i = 0 ; i < m_NuclearInteractionCollection.size() ; i++)
121  {
122  reco::NuclearInteraction ni = m_NuclearInteractionCollection[i];
123  if( ni.likelihood()<0.4) continue;
124 
125  reco::TrackRef primTrackRef = ni.primaryTrack().castTo<reco::TrackRef>();
126 
127  TrajectoryRef trajRef = m_TrackToTrajMap[primTrackRef];
128 
129  Trajectory newTraj;
130  if( newTrajNeeded(newTraj, trajRef, ni) ) {
131 
132  AlgoProductCollection algoResults;
133  bool isOK = getTrackFromTrajectory( newTraj , trajRef, algoResults);
134 
135  if( isOK ) {
136 
137  pair<unsigned int, unsigned int> tempory_pair;
138  tempory_pair.first = Output_track->size();
139  tempory_pair.second = i;
140  Indice_Map.push_back(tempory_pair);
141 
142  reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, i );
143  reco::TrackExtra newTrackExtra = getNewTrackExtra(algoResults);
144  (algoResults[0].second.first)->setExtra( teref );
145 
146  Output_track->push_back(*algoResults[0].second.first);
147  Output_trackextra->push_back( newTrackExtra );
148  Output_traj->push_back(newTraj);
149 
150  }
151  }
152  else {
154  Output_track->push_back(*primTrackRef);
155  Output_trackextra->push_back( *primTrackRef->extra() );
156  Output_traj->push_back(*trajRef);
157  }
158  }
159 
160  }
161  const edm::OrphanHandle<TrajectoryCollection> Handle_traj = iEvent.put(Output_traj);
162  const edm::OrphanHandle<reco::TrackCollection> Handle_tracks = iEvent.put(Output_track);
163  iEvent.put(Output_trackextra);
164 
165  // Make Maps between elements
166  // --------------------------------------------------------------------------------------------------
167  if(Handle_tracks->size() != Handle_traj->size() )
168  {
169  printf("ERROR Handle_tracks->size() != Handle_traj->size() \n");
170  return;
171  }
172 
173 
174 
175  for(unsigned int i = 0 ; i < Indice_Map.size() ; i++)
176  {
177  TrajectoryRef InTrajRef ( temp_m_TrajectoryCollection, Indice_Map[i].second );
178  TrajectoryRef OutTrajRef ( Handle_traj, Indice_Map[i].first );
179  reco::TrackRef TrackRef ( Handle_tracks, Indice_Map[i].first );
180 
181  Output_trajmap ->insert(OutTrajRef,InTrajRef);
182  Output_trackmap->insert(TrackRef,InTrajRef);
183 
184  try{
185  reco::TrackRef PrimaryTrackRef = m_TrajToTrackCollection->operator[]( InTrajRef );
186  Output_tracktrackmap->insert(TrackRef,PrimaryTrackRef);
187  }catch(edm::Exception event){}
188 
189  }
190  iEvent.put(Output_trajmap);
191  iEvent.put(Output_trackmap);
192  iEvent.put(Output_tracktrackmap);
193 
194 
195  if(verbosity>=3)printf("-----------------------\n");
196 }
#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:13
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)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
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
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:405
RefProd< PROD > getRefBeforePut()
Definition: Event.h:133
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
cast to a concrete type
Definition: RefToBase.h:241
const T & get() const
Definition: EventSetup.h:55
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:19
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 336 of file NuclearTrackCorrector.cc.

References i, m_TrajToTrackCollection, and query::result.

Referenced by produce().

336  {
337  for(unsigned int i = 0 ; i < trajColl->size() ; i++)
338  {
339  TrajectoryRef InTrajRef ( trajColl, i);
340  reco::TrackRef PrimaryTrackRef = m_TrajToTrackCollection->operator[]( InTrajRef );
341  result[ PrimaryTrackRef ] = InTrajRef;
342  }
343 }
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 126 of file NuclearTrackCorrector.h.

Referenced by produce().

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

Definition at line 119 of file NuclearTrackCorrector.h.

Referenced by produce().

int NuclearTrackCorrector::int_Input_Hit_Distance
private

Definition at line 114 of file NuclearTrackCorrector.h.

int NuclearTrackCorrector::KeepOnlyCorrectedTracks
private

Definition at line 117 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

const TrajTrackAssociationCollection* NuclearTrackCorrector::m_TrajToTrackCollection
private

Definition at line 130 of file NuclearTrackCorrector.h.

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

std::string NuclearTrackCorrector::str_Input_NuclearInteraction
private

Definition at line 113 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

std::string NuclearTrackCorrector::str_Input_Trajectory
private

Definition at line 112 of file NuclearTrackCorrector.h.

Referenced by NuclearTrackCorrector(), and produce().

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

Definition at line 129 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and NuclearTrackCorrector().

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

Definition at line 124 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

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

Definition at line 122 of file NuclearTrackCorrector.h.

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

TransientInitialStateEstimator* NuclearTrackCorrector::theInitialState
private

Definition at line 127 of file NuclearTrackCorrector.h.

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

Definition at line 123 of file NuclearTrackCorrector.h.

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

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

Definition at line 125 of file NuclearTrackCorrector.h.

Referenced by getTrackFromTrajectory(), and produce().

int NuclearTrackCorrector::verbosity
private

Definition at line 116 of file NuclearTrackCorrector.h.

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