Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryBase.h"
00017
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "FWCore/Framework/interface/EventSetup.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00023 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00024
00025 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00026
00027 #include "Alignment/ReferenceTrajectories/interface/DualKalmanTrajectory.h"
00028
00029 #include <algorithm>
00030
00031
00032 class DualKalmanFactory : public TrajectoryFactoryBase
00033 {
00034
00035 public:
00036
00037 DualKalmanFactory( const edm::ParameterSet & config );
00038 virtual ~DualKalmanFactory();
00039
00041 virtual const ReferenceTrajectoryCollection
00042 trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks,
00043 const reco::BeamSpot &beamSpot) const;
00044
00045 virtual const ReferenceTrajectoryCollection
00046 trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks,
00047 const ExternalPredictionCollection &external, const reco::BeamSpot &beamSpot) const;
00048
00049 virtual DualKalmanFactory* clone() const { return new DualKalmanFactory(*this); }
00050
00051 protected:
00052
00053 struct DualKalmanInput
00054 {
00055 TrajectoryStateOnSurface refTsos;
00056 Trajectory::DataContainer trajMeasurements;
00057 std::vector<unsigned int> fwdRecHitNums;
00058 std::vector<unsigned int> bwdRecHitNums;
00059 };
00060
00061 const DualKalmanInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const;
00062
00063
00064
00065
00066
00067 const double theMass;
00068 const int theResidMethod;
00069 };
00070
00071
00072
00073
00074
00075 DualKalmanFactory::DualKalmanFactory(const edm::ParameterSet &config)
00076 : TrajectoryFactoryBase(config), theMass(config.getParameter<double>("ParticleMass")),
00077 theResidMethod(config.getParameter<int>("ResidualMethod"))
00078 {
00079
00080 edm::LogInfo("Alignment") << "@SUB=DualKalmanFactory" << "Factory created.";
00081 }
00082
00083
00084
00085 DualKalmanFactory::~DualKalmanFactory() {}
00086
00087
00088
00089 const DualKalmanFactory::ReferenceTrajectoryCollection
00090 DualKalmanFactory::trajectories(const edm::EventSetup &setup,
00091 const ConstTrajTrackPairCollection &tracks,
00092 const reco::BeamSpot &beamSpot) const
00093 {
00094 ReferenceTrajectoryCollection trajectories;
00095
00096 edm::ESHandle<MagneticField> magneticField;
00097 setup.get<IdealMagneticFieldRecord>().get(magneticField);
00098
00099 ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00100
00101 while (itTracks != tracks.end()) {
00102 const DualKalmanInput input = this->referenceStateAndRecHits(*itTracks);
00103
00104 if (input.refTsos.isValid()) {
00105 ReferenceTrajectoryPtr ptr(new DualKalmanTrajectory(input.trajMeasurements,
00106 input.refTsos,
00107 input.fwdRecHitNums,
00108 input.bwdRecHitNums,
00109 magneticField.product(),
00110 this->materialEffects(),
00111 this->propagationDirection(),
00112 theMass, theUseBeamSpot, beamSpot,
00113 theResidMethod));
00114 trajectories.push_back(ptr);
00115 }
00116 ++itTracks;
00117 }
00118
00119 return trajectories;
00120 }
00121
00122
00123 const DualKalmanFactory::ReferenceTrajectoryCollection
00124 DualKalmanFactory::trajectories(const edm::EventSetup &setup,
00125 const ConstTrajTrackPairCollection &tracks,
00126 const ExternalPredictionCollection &external,
00127 const reco::BeamSpot &beamSpot) const
00128 {
00129 ReferenceTrajectoryCollection trajectories;
00130
00131 edm::LogError("Alignment") << "@SUB=DualKalmanFactory::trajectories"
00132 << "Not implemented with ExternalPrediction.";
00133 return trajectories;
00134
00135
00136 if (tracks.size() != external.size()) {
00137 edm::LogInfo("ReferenceTrajectories") << "@SUB=DualKalmanFactory::trajectories"
00138 << "Inconsistent input:\n"
00139 << "\tnumber of tracks = " << tracks.size()
00140 << "\tnumber of external predictions = " << external.size();
00141 return trajectories;
00142 }
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188 return trajectories;
00189 }
00190
00191
00192
00193 const DualKalmanFactory::DualKalmanInput
00194 DualKalmanFactory::referenceStateAndRecHits(const ConstTrajTrackPair &track) const
00195 {
00196
00197
00198
00199
00200
00201
00202
00203 DualKalmanInput input;
00204
00205
00206 input.trajMeasurements = this->orderedTrajectoryMeasurements(*track.first);
00207
00208
00209 std::vector<unsigned int> usedTrajMeasNums;
00210 for (unsigned int iM = 0; iM < input.trajMeasurements.size(); ++iM) {
00211 if (this->useRecHit(input.trajMeasurements[iM].recHit())) usedTrajMeasNums.push_back(iM);
00212 }
00213 unsigned int nRefStateMeas = usedTrajMeasNums.size()/2;
00214
00215
00216 for (unsigned int iMeas = 0; iMeas < usedTrajMeasNums.size(); ++iMeas) {
00217 if (iMeas < nRefStateMeas) {
00218 input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00219 } else if (iMeas > nRefStateMeas) {
00220 input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00221 } else {
00222 if (input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState().isValid()) {
00223 input.refTsos = input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState();
00224 input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00225 input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00226 } else {
00227
00228 ++nRefStateMeas;
00229 input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00230 }
00231 }
00232 }
00233
00234
00235 std::reverse(input.bwdRecHitNums.begin(), input.bwdRecHitNums.end());
00236
00237 return input;
00238 }
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, DualKalmanFactory, "DualKalmanFactory" );