Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-10 23:04:48

0001 #ifndef Alignment_MillePedeAlignmentAlgorithm_MillePedeAlignmentAlgorithm_h
0002 #define Alignment_MillePedeAlignmentAlgorithm_MillePedeAlignmentAlgorithm_h
0003 
0004 /// \class MillePedeAlignmentAlgorithm
0005 ///
0006 ///  CMSSW interface to pede: produces pede's binary input and steering file(s)
0007 ///
0008 ///  \author    : Gero Flucke
0009 ///  date       : October 2006
0010 ///  $Revision: 1.36 $
0011 ///  $Date: 2012/08/10 09:01:11 $
0012 ///  (last update by $Author: flucke $)
0013 
0014 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentAlgorithmPluginFactory.h"
0015 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentAlgorithmBase.h"
0016 
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/Framework/interface/ConsumesCollector.h"
0019 
0020 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0021 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0022 
0023 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectoryBase.h"
0024 
0025 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0026 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0027 
0028 #include "Geometry/TrackerGeometryBuilder/interface/PixelTopologyMap.h"
0029 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0030 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0031 
0032 #include "CondFormats/PCLConfig/interface/AlignPCLThresholdsHG.h"
0033 #include "CondFormats/DataRecord/interface/AlignPCLThresholdsHGRcd.h"
0034 
0035 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
0036 #include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
0037 
0038 #include <vector>
0039 #include <string>
0040 #include <memory>
0041 
0042 class Alignable;
0043 class AlignableTracker;
0044 class AlignableMuon;
0045 class AlignableExtras;
0046 
0047 class AlignmentParameters;
0048 class AlignableNavigator;
0049 class AlignableDetOrUnitPtr;
0050 class AlignmentUserVariables;
0051 
0052 class AlignmentParameterStore;
0053 
0054 class IntegratedCalibrationBase;
0055 
0056 class MillePedeMonitor;
0057 class PedeSteerer;
0058 class PedeLabelerBase;
0059 class Mille;
0060 class TrajectoryFactoryBase;
0061 
0062 // already from base class - and forward declaration does not work since typedef!
0063 /* class TkFittedLasBeamCollection; */
0064 /* class TsosVectorCollection; */
0065 
0066 class MillePedeAlignmentAlgorithm : public AlignmentAlgorithmBase {
0067 public:
0068   /// Constructor
0069   MillePedeAlignmentAlgorithm(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC);
0070 
0071   /// Destructor
0072   ~MillePedeAlignmentAlgorithm() override;
0073 
0074   /// Called at beginning of job
0075   void initialize(const edm::EventSetup &setup,
0076                   AlignableTracker *tracker,
0077                   AlignableMuon *muon,
0078                   AlignableExtras *extras,
0079                   AlignmentParameterStore *store) override;
0080 
0081   /// Returns whether MP supports calibrations
0082   bool supportsCalibrations() override;
0083   /// Pass integrated calibrations to Millepede (they are not owned by Millepede!)
0084   bool addCalibrations(const std::vector<IntegratedCalibrationBase *> &iCals) override;
0085 
0086   virtual bool storeThresholds(const int &nRecords,
0087                                const AlignPCLThresholdsHG::threshold_map &thresholdMap,
0088                                const AlignPCLThresholdsHG::param_map &floatMap);
0089 
0090   /// Called at end of job
0091   void terminate(const edm::EventSetup &iSetup) override;
0092   /// Called at end of job
0093   void terminate() override;
0094 
0095   /// Returns whether MP should process events in the current configuration
0096   bool processesEvents() override;
0097 
0098   /// Returns whether MP produced results to be stored
0099   bool storeAlignments() override;
0100 
0101   /// Run the algorithm on trajectories and tracks
0102   void run(const edm::EventSetup &setup, const EventInfo &eventInfo) override;
0103 
0104   /// called at begin of run
0105   void beginRun(const edm::Run &run, const edm::EventSetup &setup, bool changed) override;
0106 
0107   // TODO: This method does NOT match endRun() in base class! Nobody is
0108   //       calling this?
0109   /// Run on run products, e.g. TkLAS
0110   virtual void endRun(const EventInfo &, const EndRunInfo &,
0111                       const edm::EventSetup &);  //override;
0112 
0113   // This one will be called since it matches the interface of the base class
0114   void endRun(const EndRunInfo &runInfo, const edm::EventSetup &setup) override;
0115 
0116   /// called at begin of luminosity block (resets Mille binary in mille mode)
0117   void beginLuminosityBlock(const edm::EventSetup &) override;
0118 
0119   /// called at end of luminosity block
0120   void endLuminosityBlock(const edm::EventSetup &) override;
0121 
0122   /// Called in order to pass parameters to alignables for a specific run
0123   /// range in case the algorithm supports run range dependent alignment.
0124   bool setParametersForRunRange(const RunRange &runrange) override;
0125 
0126 private:
0127   enum MeasurementDirection { kLocalX = 0, kLocalY };
0128 
0129   /// fill mille for a trajectory, returning number of x/y hits ([0,0] if 'bad' trajectory)
0130   std::pair<unsigned int, unsigned int> addReferenceTrajectory(
0131       const edm::EventSetup &setup,
0132       const EventInfo &eventInfo,
0133       const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr);
0134 
0135   /// If hit is usable: callMille for x and (probably) y direction.
0136   /// If globalDerivatives fine: returns 2 if 2D-hit, 1 if 1D-hit, 0 if no Alignable for hit.
0137   /// Returns -1 if any problem (for params cf. globalDerivativesHierarchy)
0138   int addMeasurementData(const edm::EventSetup &setup,
0139                          const EventInfo &eventInfo,
0140                          const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0141                          unsigned int iHit,
0142                          AlignmentParameters *&params);
0143 
0144   /// Add global data (labels, derivatives) to GBL trajectory
0145   /// Returns -1 if any problem (for params cf. globalDerivativesHierarchy)
0146   int addGlobalData(const edm::EventSetup &setup,
0147                     const EventInfo &eventInfo,
0148                     const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0149                     unsigned int iHit,
0150                     gbl::GblPoint &gblPoint);
0151 
0152   /// Increase hit counting of MillePedeVariables behind each parVec[i]
0153   /// (and also for parameters higher in hierarchy),
0154   /// assuming 'parVec' and 'validHitVecY' to be parallel.
0155   /// Returns number of valid y-hits.
0156   unsigned int addHitCount(const std::vector<AlignmentParameters *> &parVec,
0157                            const std::vector<bool> &validHitVecY) const;
0158 
0159   /// adds data from reference trajectory from a specific Hit
0160   template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
0161   void addRefTrackData2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0162                          unsigned int iTrajHit,
0163                          Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
0164                          Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
0165                          Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM);
0166 
0167   /// adds data for virtual measurements from reference trajectory
0168   void addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas);
0169 
0170   /// adds data for a specific virtual measurement from reference trajectory
0171   template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
0172   void addRefTrackVirtualMeas1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0173                                 unsigned int iVirtualMeas,
0174                                 Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
0175                                 Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
0176                                 Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM);
0177 
0178   /// recursively adding derivatives and labels, false if problems
0179   bool globalDerivativesHierarchy(const EventInfo &eventInfo,
0180                                   const TrajectoryStateOnSurface &tsos,
0181                                   Alignable *ali,
0182                                   const AlignableDetOrUnitPtr &alidet,
0183                                   std::vector<float> &globalDerivativesX,
0184                                   std::vector<float> &globalDerivativesY,
0185                                   std::vector<int> &globalLabels,
0186                                   AlignmentParameters *&lowestParams) const;
0187 
0188   /// recursively adding derivatives (double) and labels, false if problems
0189   bool globalDerivativesHierarchy(const EventInfo &eventInfo,
0190                                   const TrajectoryStateOnSurface &tsos,
0191                                   Alignable *ali,
0192                                   const AlignableDetOrUnitPtr &alidet,
0193                                   std::vector<double> &globalDerivativesX,
0194                                   std::vector<double> &globalDerivativesY,
0195                                   std::vector<int> &globalLabels,
0196                                   AlignmentParameters *&lowestParams) const;
0197 
0198   /// adding derivatives from integrated calibrations
0199   void globalDerivativesCalibration(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
0200                                     const TrajectoryStateOnSurface &tsos,
0201                                     const edm::EventSetup &setup,
0202                                     const EventInfo &eventInfo,
0203                                     std::vector<float> &globalDerivativesX,
0204                                     std::vector<float> &globalDerivativesY,
0205                                     std::vector<int> &globalLabels) const;
0206 
0207   /// calls callMille1D or callMille2D
0208   int callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0209                 unsigned int iTrajHit,
0210                 const std::vector<int> &globalLabels,
0211                 const std::vector<float> &globalDerivativesX,
0212                 const std::vector<float> &globalDerivativesY);
0213 
0214   /// calls Mille for 1D hits
0215   int callMille1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0216                   unsigned int iTrajHit,
0217                   const std::vector<int> &globalLabels,
0218                   const std::vector<float> &globalDerivativesX);
0219 
0220   /// calls Mille for x and possibly y component of hit,
0221   /// y is skipped for non-real 2D (e.g. SiStripRecHit2D),
0222   /// for TID/TEC first diagonalises if correlation is larger than configurable
0223   int callMille2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0224                   unsigned int iTrajHit,
0225                   const std::vector<int> &globalLabels,
0226                   const std::vector<float> &globalDerivativesx,
0227                   const std::vector<float> &globalDerivativesy);
0228 
0229   template <typename CovarianceMatrix,
0230             typename LocalDerivativeMatrix,
0231             typename ResidualMatrix,
0232             typename GlobalDerivativeMatrix>
0233   void diagonalize(Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
0234                    Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM,
0235                    Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
0236                    Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM) const;
0237 
0238   // deals with the non matrix format of theFloatBufferX ...
0239   template <typename GlobalDerivativeMatrix>
0240   void makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
0241                            const std::vector<float> &globalDerivativesy,
0242                            Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM);
0243 
0244   //   void callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
0245   //               unsigned int iTrajHit, MeasurementDirection xOrY,
0246   //               const std::vector<float> &globalDerivatives, const std::vector<int> &globalLabels);
0247   /// true if hit belongs to 2D detector (currently tracker specific)
0248   bool is2D(const TransientTrackingRecHit::ConstRecHitPointer &recHit) const;
0249 
0250   /// read pede input defined by 'psetName', flag to create/not create MillePedeVariables
0251   bool readFromPede(const edm::ParameterSet &mprespset, bool setUserVars, const RunRange &runrange);
0252   bool areEmptyParams(const align::Alignables &alignables) const;
0253   unsigned int doIO(int loop) const;
0254   /// add MillePedeVariables for each AlignmentParameters (exception if no parameters...)
0255   void buildUserVariables(const align::Alignables &alignables) const;
0256 
0257   /// Generates list of files to read, given the list and dir from the configuration.
0258   /// This will automatically expand formatting directives, if they appear.
0259   std::vector<std::string> getExistingFormattedFiles(const std::vector<std::string> &plainFiles,
0260                                                      const std::string &theDir);
0261 
0262   void addLaserData(const EventInfo &eventInfo,
0263                     const TkFittedLasBeamCollection &tkLasBeams,
0264                     const TsosVectorCollection &tkLasBeamTsoses);
0265   void addLasBeam(const EventInfo &eventInfo,
0266                   const TkFittedLasBeam &lasBeam,
0267                   const std::vector<TrajectoryStateOnSurface> &tsoses);
0268 
0269   /// add measurement data from PXB survey
0270   void addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg);
0271 
0272   //
0273   bool areIOVsSpecified() const;
0274 
0275   //--------------------------------------------------------
0276   // Data members
0277   //--------------------------------------------------------
0278 
0279   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoToken_;
0280   const edm::ESGetToken<AlignPCLThresholdsHG, AlignPCLThresholdsHGRcd> aliThrToken_;
0281   const edm::ESGetToken<SiPixelQuality, SiPixelQualityFromDbRcd> siPixelQualityToken_;
0282   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0283 
0284   enum EModeBit { myMilleBit = 1 << 0, myPedeRunBit = 1 << 1, myPedeSteerBit = 1 << 2, myPedeReadBit = 1 << 3 };
0285   unsigned int decodeMode(const std::string &mode) const;
0286   bool isMode(unsigned int testMode) const { return (theMode & testMode); }
0287   bool addHitStatistics(int fromLoop, const std::string &outFile, const std::vector<std::string> &inFiles) const;
0288   bool addHits(const align::Alignables &alis, const std::vector<AlignmentUserVariables *> &mpVars) const;
0289 
0290   edm::ParameterSet theConfig;
0291   unsigned int theMode;
0292   std::string theDir;  /// directory for all kind of files
0293   AlignmentParameterStore *theAlignmentParameterStore;
0294   align::Alignables theAlignables;
0295   std::unique_ptr<AlignableNavigator> theAlignableNavigator;
0296   std::unique_ptr<MillePedeMonitor> theMonitor;
0297   std::unique_ptr<Mille> theMille;
0298   std::shared_ptr<PedeLabelerBase> thePedeLabels;
0299   std::unique_ptr<PedeSteerer> thePedeSteer;
0300   std::unique_ptr<TrajectoryFactoryBase> theTrajectoryFactory;
0301   std::vector<IntegratedCalibrationBase *> theCalibrations;
0302   std::shared_ptr<AlignPCLThresholdsHG> theThresholds;
0303   std::shared_ptr<PixelTopologyMap> pixelTopologyMap;
0304   std::shared_ptr<SiPixelQuality> pixelQuality;
0305   unsigned int theMinNumHits;
0306   double theMaximalCor2D;  /// maximal correlation allowed for 2D hit in TID/TEC.
0307                            /// If larger, the 2D measurement gets diagonalized!!!
0308   const align::RunNumber firstIOV_;
0309   const bool ignoreFirstIOVCheck_;
0310   const bool enableAlignableUpdates_;
0311   int theLastWrittenIov;  // keeping track for output trees...
0312   std::vector<float> theFloatBufferX;
0313   std::vector<float> theFloatBufferY;
0314   std::vector<int> theIntBuffer;
0315   bool theDoSurveyPixelBarrel;
0316   // CHK for GBL
0317   std::unique_ptr<gbl::MilleBinary> theBinary;
0318   bool theGblDoubleBinary;
0319 
0320   const bool runAtPCL_;
0321   const bool ignoreHitsWithoutGlobalDerivatives_;
0322   const bool skipGlobalPositionRcdCheck_;
0323 
0324   const align::RunRanges uniqueRunRanges_;
0325   const bool enforceSingleIOVInput_;
0326   std::vector<align::RunNumber> cachedRuns_;
0327   align::RunNumber lastProcessedRun_;
0328 };
0329 
0330 DEFINE_EDM_PLUGIN(AlignmentAlgorithmPluginFactory, MillePedeAlignmentAlgorithm, "MillePedeAlignmentAlgorithm");
0331 #endif