CMS 3D CMS Logo

MillePedeAlignmentAlgorithm.h
Go to the documentation of this file.
1 #ifndef Alignment_MillePedeAlignmentAlgorithm_MillePedeAlignmentAlgorithm_h
2 #define Alignment_MillePedeAlignmentAlgorithm_MillePedeAlignmentAlgorithm_h
3 
13 
16 
18 
21 
24 
25 #include <vector>
26 #include <string>
27 #include <memory>
28 
29 class Alignable;
30 class AlignableTracker;
31 class AlignableMuon;
32 class AlignableExtras;
33 
35 class AlignableNavigator;
38 
40 
42 
43 class MillePedeMonitor;
44 class PedeSteerer;
45 class PedeLabelerBase;
46 class Mille;
48 
49 // already from base class - and forward declaration does not work since typedef!
50 /* class TkFittedLasBeamCollection; */
51 /* class TsosVectorCollection; */
52 
54 public:
57 
60 
62  void initialize(const edm::EventSetup &setup,
65  AlignableExtras *extras,
66  AlignmentParameterStore *store) override;
67 
69  bool supportsCalibrations() override;
71  bool addCalibrations(const std::vector<IntegratedCalibrationBase *> &iCals) override;
72 
73  virtual bool storeThresholds(const int &nRecords, const AlignPCLThresholds::threshold_map &thresholdMap);
74 
76  void terminate(const edm::EventSetup &iSetup) override;
78  void terminate() override;
79 
81  bool processesEvents() override;
82 
84  bool storeAlignments() override;
85 
87  void run(const edm::EventSetup &setup, const EventInfo &eventInfo) override;
88 
90  void beginRun(const edm::Run &run, const edm::EventSetup &setup, bool changed) override;
91 
92  // TODO: This method does NOT match endRun() in base class! Nobody is
93  // calling this?
95  virtual void endRun(const EventInfo &, const EndRunInfo &,
96  const edm::EventSetup &); //override;
97 
98  // This one will be called since it matches the interface of the base class
99  void endRun(const EndRunInfo &runInfo, const edm::EventSetup &setup) override;
100 
102  void beginLuminosityBlock(const edm::EventSetup &) override;
103 
105  void endLuminosityBlock(const edm::EventSetup &) override;
106 
107  /* virtual void beginLuminosityBlock(const edm::EventSetup &setup) {} */
108  /* virtual void endLuminosityBlock(const edm::EventSetup &setup) {} */
109 
112  bool setParametersForRunRange(const RunRange &runrange) override;
113 
114 private:
116 
118  std::pair<unsigned int, unsigned int> addReferenceTrajectory(
119  const edm::EventSetup &setup,
120  const EventInfo &eventInfo,
122 
126  int addMeasurementData(const edm::EventSetup &setup,
127  const EventInfo &eventInfo,
129  unsigned int iHit,
131 
134  int addGlobalData(const edm::EventSetup &setup,
135  const EventInfo &eventInfo,
137  unsigned int iHit,
138  gbl::GblPoint &gblPoint);
139 
144  unsigned int addHitCount(const std::vector<AlignmentParameters *> &parVec,
145  const std::vector<bool> &validHitVecY) const;
146 
148  template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
150  unsigned int iTrajHit,
151  Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
152  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
153  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM);
154 
156  void addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas);
157 
159  template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
161  unsigned int iVirtualMeas,
162  Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
163  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
164  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM);
165 
167  bool globalDerivativesHierarchy(const EventInfo &eventInfo,
168  const TrajectoryStateOnSurface &tsos,
169  Alignable *ali,
170  const AlignableDetOrUnitPtr &alidet,
171  std::vector<float> &globalDerivativesX,
172  std::vector<float> &globalDerivativesY,
173  std::vector<int> &globalLabels,
174  AlignmentParameters *&lowestParams) const;
175 
177  bool globalDerivativesHierarchy(const EventInfo &eventInfo,
178  const TrajectoryStateOnSurface &tsos,
179  Alignable *ali,
180  const AlignableDetOrUnitPtr &alidet,
181  std::vector<double> &globalDerivativesX,
182  std::vector<double> &globalDerivativesY,
183  std::vector<int> &globalLabels,
184  AlignmentParameters *&lowestParams) const;
185 
188  const TrajectoryStateOnSurface &tsos,
189  const edm::EventSetup &setup,
190  const EventInfo &eventInfo,
191  std::vector<float> &globalDerivativesX,
192  std::vector<float> &globalDerivativesY,
193  std::vector<int> &globalLabels) const;
194 
197  unsigned int iTrajHit,
198  const std::vector<int> &globalLabels,
199  const std::vector<float> &globalDerivativesX,
200  const std::vector<float> &globalDerivativesY);
201 
204  unsigned int iTrajHit,
205  const std::vector<int> &globalLabels,
206  const std::vector<float> &globalDerivativesX);
207 
212  unsigned int iTrajHit,
213  const std::vector<int> &globalLabels,
214  const std::vector<float> &globalDerivativesx,
215  const std::vector<float> &globalDerivativesy);
216 
217  template <typename CovarianceMatrix,
218  typename LocalDerivativeMatrix,
219  typename ResidualMatrix,
220  typename GlobalDerivativeMatrix>
221  void diagonalize(Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
222  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM,
223  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
224  Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM) const;
225 
226  // deals with the non matrix format of theFloatBufferX ...
227  template <typename GlobalDerivativeMatrix>
228  void makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
229  const std::vector<float> &globalDerivativesy,
230  Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM);
231 
232  // void callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
233  // unsigned int iTrajHit, MeasurementDirection xOrY,
234  // const std::vector<float> &globalDerivatives, const std::vector<int> &globalLabels);
236  bool is2D(const TransientTrackingRecHit::ConstRecHitPointer &recHit) const;
237 
239  bool readFromPede(const edm::ParameterSet &mprespset, bool setUserVars, const RunRange &runrange);
240  bool areEmptyParams(const align::Alignables &alignables) const;
241  unsigned int doIO(int loop) const;
243  void buildUserVariables(const align::Alignables &alignables) const;
244 
247  std::vector<std::string> getExistingFormattedFiles(const std::vector<std::string> &plainFiles,
248  const std::string &theDir);
249 
250  void addLaserData(const EventInfo &eventInfo,
251  const TkFittedLasBeamCollection &tkLasBeams,
252  const TsosVectorCollection &tkLasBeamTsoses);
253  void addLasBeam(const EventInfo &eventInfo,
254  const TkFittedLasBeam &lasBeam,
255  const std::vector<TrajectoryStateOnSurface> &tsoses);
256 
258  void addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg);
259 
260  //
261  bool areIOVsSpecified() const;
262 
263  //--------------------------------------------------------
264  // Data members
265  //--------------------------------------------------------
266  enum EModeBit { myMilleBit = 1 << 0, myPedeRunBit = 1 << 1, myPedeSteerBit = 1 << 2, myPedeReadBit = 1 << 3 };
267  unsigned int decodeMode(const std::string &mode) const;
268  bool isMode(unsigned int testMode) const { return (theMode & testMode); }
269  bool addHitStatistics(int fromLoop, const std::string &outFile, const std::vector<std::string> &inFiles) const;
270  bool addHits(const align::Alignables &alis, const std::vector<AlignmentUserVariables *> &mpVars) const;
271 
273  unsigned int theMode;
277  std::unique_ptr<AlignableNavigator> theAlignableNavigator;
278  std::unique_ptr<MillePedeMonitor> theMonitor;
279  std::unique_ptr<Mille> theMille;
280  std::shared_ptr<PedeLabelerBase> thePedeLabels;
281  std::unique_ptr<PedeSteerer> thePedeSteer;
282  std::unique_ptr<TrajectoryFactoryBase> theTrajectoryFactory;
283  std::vector<IntegratedCalibrationBase *> theCalibrations;
284  std::shared_ptr<AlignPCLThresholds> theThresholds;
285  unsigned int theMinNumHits;
291  int theLastWrittenIov; // keeping track for output trees...
292  std::vector<float> theFloatBufferX;
293  std::vector<float> theFloatBufferY;
294  std::vector<int> theIntBuffer;
296  // CHK for GBL
297  std::unique_ptr<gbl::MilleBinary> theBinary;
299 
300  const bool runAtPCL_;
303 
306  std::vector<align::RunNumber> cachedRuns_;
308 };
309 
311 #endif
bool processesEvents() override
Returns whether MP should process events in the current configuration.
void globalDerivativesCalibration(const TransientTrackingRecHit::ConstRecHitPointer &recHit, const TrajectoryStateOnSurface &tsos, const edm::EventSetup &setup, const EventInfo &eventInfo, std::vector< float > &globalDerivativesX, std::vector< float > &globalDerivativesY, std::vector< int > &globalLabels) const
adding derivatives from integrated calibrations
bool supportsCalibrations() override
Returns whether MP supports calibrations.
std::unique_ptr< MillePedeMonitor > theMonitor
bool setParametersForRunRange(const RunRange &runrange) override
virtual void endRun(const EventInfo &, const EndRunInfo &, const edm::EventSetup &)
Run on run products, e.g. TkLAS.
std::shared_ptr< AlignPCLThresholds > theThresholds
void diagonalize(Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< GlobalDerivativeMatrix > &aGlobalDerivativesM) const
std::map< std::string, AlignPCLThreshold > threshold_map
bool areEmptyParams(const align::Alignables &alignables) const
void endLuminosityBlock(const edm::EventSetup &) override
called at end of luminosity block
void buildUserVariables(const align::Alignables &alignables) const
add MillePedeVariables for each AlignmentParameters (exception if no parameters...)
std::unique_ptr< TrajectoryFactoryBase > theTrajectoryFactory
std::vector< RunRange > RunRanges
Definition: Utilities.h:39
std::unique_ptr< gbl::MilleBinary > theBinary
int addGlobalData(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iHit, gbl::GblPoint &gblPoint)
bool globalDerivativesHierarchy(const EventInfo &eventInfo, const TrajectoryStateOnSurface &tsos, Alignable *ali, const AlignableDetOrUnitPtr &alidet, std::vector< float > &globalDerivativesX, std::vector< float > &globalDerivativesY, std::vector< int > &globalLabels, AlignmentParameters *&lowestParams) const
recursively adding derivatives and labels, false if problems
std::shared_ptr< PedeLabelerBase > thePedeLabels
std::vector< std::string > getExistingFormattedFiles(const std::vector< std::string > &plainFiles, const std::string &theDir)
unsigned int addHitCount(const std::vector< AlignmentParameters * > &parVec, const std::vector< bool > &validHitVecY) const
bool addCalibrations(const std::vector< IntegratedCalibrationBase * > &iCals) override
Pass integrated calibrations to Millepede (they are not owned by Millepede!)
int callMille2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesx, const std::vector< float > &globalDerivativesy)
math::Error< 5 >::type CovarianceMatrix
Interface/Base class for alignment algorithms, each alignment algorithm has to be derived from this c...
void addLasBeam(const EventInfo &eventInfo, const TkFittedLasBeam &lasBeam, const std::vector< TrajectoryStateOnSurface > &tsoses)
define event information passed to algorithms
unsigned int doIO(int loop) const
int callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesX, const std::vector< float > &globalDerivativesY)
calls callMille1D or callMille2D
void beginLuminosityBlock(const edm::EventSetup &) override
called at begin of luminosity block (resets Mille binary in mille mode)
std::vector< align::RunNumber > cachedRuns_
int callMille1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesX)
calls Mille for 1D hits
void addLaserData(const EventInfo &eventInfo, const TkFittedLasBeamCollection &tkLasBeams, const TsosVectorCollection &tkLasBeamTsoses)
int addMeasurementData(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iHit, AlignmentParameters *&params)
bool isMode(unsigned int testMode) const
void beginRun(const edm::Run &run, const edm::EventSetup &setup, bool changed) override
called at begin of run
bool addHits(const align::Alignables &alis, const std::vector< AlignmentUserVariables * > &mpVars) const
void makeGlobDerivMatrix(const std::vector< float > &globalDerivativesx, const std::vector< float > &globalDerivativesy, Eigen::MatrixBase< GlobalDerivativeMatrix > &aGlobalDerivativesM)
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
AlignmentParameterStore * theAlignmentParameterStore
directory for all kind of files
bool is2D(const TransientTrackingRecHit::ConstRecHitPointer &recHit) const
true if hit belongs to 2D detector (currently tracker specific)
void addRefTrackVirtualMeas1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas, Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM)
adds data for a specific virtual measurement from reference trajectory
void addRefTrackData2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM)
adds data from reference trajectory from a specific Hit
(Abstract) Base class for alignment algorithm user variables
std::unique_ptr< PedeSteerer > thePedeSteer
bool storeAlignments() override
Returns whether MP produced results to be stored.
void terminate() override
Called at end of job.
void initialize(const edm::EventSetup &setup, AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras, AlignmentParameterStore *store) override
Called at beginning of job.
bool addHitStatistics(int fromLoop, const std::string &outFile, const std::vector< std::string > &inFiles) const
void addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas)
adds data for virtual measurements from reference trajectory
bool readFromPede(const edm::ParameterSet &mprespset, bool setUserVars, const RunRange &runrange)
read pede input defined by &#39;psetName&#39;, flag to create/not create MillePedeVariables ...
std::pair< unsigned int, unsigned int > addReferenceTrajectory(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr)
fill mille for a trajectory, returning number of x/y hits ([0,0] if &#39;bad&#39; trajectory) ...
virtual bool storeThresholds(const int &nRecords, const AlignPCLThresholds::threshold_map &thresholdMap)
std::vector< IntegratedCalibrationBase * > theCalibrations
std::vector< Alignable * > Alignables
Definition: Utilities.h:31
void run(const edm::EventSetup &setup, const EventInfo &eventInfo) override
Run the algorithm on trajectories and tracks.
MillePedeAlignmentAlgorithm(const edm::ParameterSet &cfg)
Constructor.
~MillePedeAlignmentAlgorithm() override
Destructor.
std::vector< TkFittedLasBeam > TkFittedLasBeamCollection
Definition: Mille.h:26
#define DEFINE_EDM_PLUGIN(factory, type, name)
void addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg)
add measurement data from PXB survey
unsigned int decodeMode(const std::string &mode) const
std::unique_ptr< AlignableNavigator > theAlignableNavigator
std::vector< std::vector< TrajectoryStateOnSurface > > TsosVectorCollection
define run information passed to algorithms (in endRun)
Constructor of the full muon geometry.
Definition: AlignableMuon.h:33
Definition: Run.h:45
cond::RealTimeType< cond::runnumber >::type RunNumber
Definition: Utilities.h:37