CMS 3D CMS Logo

ReferenceTrajectoryFactory.cc
Go to the documentation of this file.
11 
13 
15 
19 
21 {
22 public:
24  ~ReferenceTrajectoryFactory() override;
25 
29  const reco::BeamSpot &beamSpot) const override;
30 
32  const ConstTrajTrackPairCollection &tracks,
34  const reco::BeamSpot &beamSpot) const override;
35 
36  ReferenceTrajectoryFactory* clone() const override { return new ReferenceTrajectoryFactory(*this); }
37 
38 protected:
40  const TrajectoryFactoryBase *bzeroFactory() const;
41 
42  double theMass;
45 };
46 
50 
52  TrajectoryFactoryBase( config ),
53  theMass(config.getParameter<double>("ParticleMass")),
54  theUseBzeroIfFieldOff(config.getParameter<bool>("UseBzeroIfFieldOff")),
56 {
57  edm::LogInfo("Alignment") << "@SUB=ReferenceTrajectoryFactory"
58  << "mass: " << theMass
59  << "\nusing Bzero if |B| = 0: "
60  << (theUseBzeroIfFieldOff ? "yes" : "no");
61 }
62 
64  TrajectoryFactoryBase(other),
65  theMass(other.theMass),
67  theBzeroFactory(nullptr) // copy data members, but no double pointing to same Bzero factory...
68 {
69 }
70 
72 {
73  delete theBzeroFactory;
74 }
75 
76 
80  const reco::BeamSpot &beamSpot) const
81 {
83  setup.get< IdealMagneticFieldRecord >().get( magneticField );
84  if (theUseBzeroIfFieldOff && magneticField->inTesla(GlobalPoint(0.,0.,0.)).mag2() < 1.e-6) {
85  return this->bzeroFactory()->trajectories(setup, tracks, beamSpot);
86  }
87 
89 
90  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
91 
92  while ( itTracks != tracks.end() )
93  {
94  TrajectoryInput input = this->innermostStateAndRecHits( *itTracks );
95 
96  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
97  if ( input.first.isValid() )
98  {
100  config.useBeamSpot = useBeamSpot_;
101  config.includeAPEs = includeAPEs_;
103  // set the flag for reversing the RecHits to false, since they are already in the correct order.
104  config.hitsAreReverse = false;
105  trajectories.push_back(ReferenceTrajectoryPtr(new ReferenceTrajectory(input.first, input.second,
106  magneticField.product(),
107  beamSpot, config)));
108  }
109 
110  ++itTracks;
111  }
112 
113  return trajectories;
114 }
115 
116 
121  const reco::BeamSpot &beamSpot) const
122 {
124 
125  if ( tracks.size() != external.size() )
126  {
127  edm::LogInfo("ReferenceTrajectories") << "@SUB=ReferenceTrajectoryFactory::trajectories"
128  << "Inconsistent input:\n"
129  << "\tnumber of tracks = " << tracks.size()
130  << "\tnumber of external predictions = " << external.size();
131  return trajectories;
132  }
133 
135  setup.get< IdealMagneticFieldRecord >().get( magneticField );
136  if (theUseBzeroIfFieldOff && magneticField->inTesla(GlobalPoint(0.,0.,0.)).mag2() < 1.e-6) {
137  return this->bzeroFactory()->trajectories(setup, tracks, external, beamSpot);
138  }
139 
140  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
141  ExternalPredictionCollection::const_iterator itExternal = external.begin();
142 
143  while ( itTracks != tracks.end() )
144  {
146  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
147  if ( input.first.isValid() )
148  {
149  if ( (*itExternal).isValid() && sameSurface( (*itExternal).surface(), input.first.surface() ) )
150  {
152  config.useBeamSpot = useBeamSpot_;
153  config.includeAPEs = includeAPEs_;
155  // set the flag for reversing the RecHits to false, since they are already in the correct order.
156  config.hitsAreReverse = false;
157  ReferenceTrajectoryPtr refTraj(new ReferenceTrajectory(*itExternal, input.second,
158  magneticField.product(),
159  beamSpot, config));
160 
161  AlgebraicSymMatrix externalParamErrors( asHepMatrix<5>( (*itExternal).localError().matrix() ) );
162  refTraj->setParameterErrors( externalParamErrors );
163  trajectories.push_back( refTraj );
164  }
165  else
166  {
168  config.useBeamSpot = useBeamSpot_;
169  config.includeAPEs = includeAPEs_;
171  config.hitsAreReverse = false;
172  trajectories.push_back(ReferenceTrajectoryPtr(new ReferenceTrajectory(input.first, input.second,
173  magneticField.product(),
174  beamSpot, config)));
175  }
176  }
177 
178  ++itTracks;
179  ++itExternal;
180  }
181 
182  return trajectories;
183 }
184 
186 {
187  if (!theBzeroFactory) {
188  const edm::ParameterSet &myPset = this->configuration();
189  edm::LogInfo("Alignment") << "@SUB=ReferenceTrajectoryFactory::bzeroFactory"
190  << "Using BzeroReferenceTrajectoryFactory for some (all?) events.";
191  // We take the config of this factory, copy it, replace its name and add
192  // the momentum parameter as expected by BzeroReferenceTrajectoryFactory and create it:
194  pset.copyForModify(myPset);
195  // next two lines not needed, but may help to better understand log file:
196  pset.eraseSimpleParameter("TrajectoryFactoryName");
197  pset.addParameter("TrajectoryFactoryName", std::string("BzeroReferenceTrajectoryFactory"));
198  pset.addParameter("MomentumEstimate", myPset.getParameter<double>("MomentumEstimateFieldOff"));
200  }
201  return theBzeroFactory;
202 }
203 
T getParameter(std::string const &) const
T mag2() const
Definition: PV3DBase.h:66
const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const override
Produce the reference trajectories.
MaterialEffects materialEffects(void) const
ReferenceTrajectoryFactory * clone() const override
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
#define nullptr
static std::string const input
Definition: EdmProvDump.cc:44
config
Definition: looper.py:287
bool sameSurface(const Surface &s1, const Surface &s2) const
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:144
BzeroReferenceTrajectoryFactory
BzeroReferenceTrajectoryFactory.
const edm::ParameterSet & configuration() const
void eraseSimpleParameter(std::string const &name)
void copyForModify(ParameterSet const &other)
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const =0
const TrajectoryFactoryBase * bzeroFactory() const
virtual const TrajectoryInput innermostStateAndRecHits(const ConstTrajTrackPair &track) const
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
const T & get() const
Definition: EventSetup.h:58
const TrajectoryFactoryBase * theBzeroFactory
ReferenceTrajectoryFactory(const edm::ParameterSet &config)
PropagationDirection propagationDirection(void) const
CLHEP::HepSymMatrix AlgebraicSymMatrix
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::pair< TrajectoryStateOnSurface, TransientTrackingRecHit::ConstRecHitContainer > TrajectoryInput
std::vector< TrajectoryStateOnSurface > ExternalPredictionCollection
ReferenceTrajectoryBase::ReferenceTrajectoryPtr ReferenceTrajectoryPtr
T const * product() const
Definition: ESHandle.h:86