CMS 3D CMS Logo

AlignmentAlgorithmBase.h
Go to the documentation of this file.
1 
2 #ifndef Alignment_CommonAlignmentAlgorithm_AlignmentAlgorithmBase_h
3 #define Alignment_CommonAlignmentAlgorithm_AlignmentAlgorithmBase_h
4 
19 #include <vector>
20 #include <utility>
21 
22 class AlignableTracker;
23 class AlignableMuon;
24 class AlignableExtras;
27 class Trajectory;
28 // These data formats cannot be forward declared since they are typedef's,
29 // so include the headers that define the typedef's
30 // (no need to include in dependencies in BuildFile):
31 // class TsosVectorCollection;
32 // class TkFittedLasBeamCollection;
33 // class AliClusterValueMap;
39 
40 namespace edm {
41  class EventSetup;
42  class ParameterSet;
43 } // namespace edm
44 namespace reco {
45  class Track;
46  class BeamSpot;
47 } // namespace reco
48 
49 /*** Global typedefs part I (see EOF for part II) ***/
50 typedef std::pair<const Trajectory *, const reco::Track *> ConstTrajTrackPair;
51 typedef std::vector<ConstTrajTrackPair> ConstTrajTrackPairs;
52 
53 typedef std::vector<IntegratedCalibrationBase *> Calibrations;
54 
56 public:
57  // TODO: DEPRECATED: For not breaking the interface, used in serveral files.
58  // If possible use the global typedefs above.
59  // With global typedefs one does not have to typedef again like
60  // 'typedef AlignmentAlgorithmBase::ConstTrajTrackPair ConstTrajTrackPair;'
61  // in other files.
62  typedef std::pair<const Trajectory *, const reco::Track *> ConstTrajTrackPair;
63  typedef std::vector<ConstTrajTrackPair> ConstTrajTrackPairCollection;
66 
68  class EventInfo {
69  public:
70  EventInfo(const edm::EventID &theEventId,
71  const ConstTrajTrackPairCollection &theTrajTrackPairs,
72  const reco::BeamSpot &theBeamSpot,
73  const AliClusterValueMap *theClusterValueMap)
74  : eventId_(theEventId),
75  trajTrackPairs_(theTrajTrackPairs),
76  beamSpot_(theBeamSpot),
77  clusterValueMap_(theClusterValueMap) {}
78 
79  const edm::EventID eventId() const { return eventId_; }
80  const ConstTrajTrackPairCollection &trajTrackPairs() const { return trajTrackPairs_; }
81  const reco::BeamSpot &beamSpot() const { return beamSpot_; }
82  const AliClusterValueMap *clusterValueMap() const { return clusterValueMap_; }
83 
84  private:
86  const ConstTrajTrackPairCollection &trajTrackPairs_;
89  };
90 
92  class EndRunInfo {
93  public:
94  EndRunInfo(const edm::RunID &theRunId,
95  const TkFittedLasBeamCollection *theTkLasBeams,
96  const TsosVectorCollection *theTkLasBeamTsoses)
97  : runId_(theRunId), tkLasBeams_(theTkLasBeams), tkLasBeamTsoses_(theTkLasBeamTsoses) {}
98 
99  const edm::RunID runId() const { return runId_; }
100  const TkFittedLasBeamCollection *tkLasBeams() const { return tkLasBeams_; }
101  const TsosVectorCollection *tkLasBeamTsoses() const { return tkLasBeamTsoses_; }
102 
103  private:
107  };
108 
111 
114 
116  virtual void initialize(const edm::EventSetup &setup,
119  AlignableExtras *extras,
120  AlignmentParameterStore *store) = 0;
121 
124  virtual bool supportsCalibrations() { return false; }
127  virtual bool addCalibrations(const Calibrations &) { return false; }
128 
130  virtual bool processesEvents() { return true; }
131 
133  virtual bool storeAlignments() { return true; }
134 
135  // TODO: DEPRECATED: Actually, there are no iterative algorithms, use
136  // initialze() and terminate()
139  virtual void startNewLoop() {}
140 
142  virtual void terminate(const edm::EventSetup &iSetup) = 0;
144  virtual void terminate() {}
145 
147  virtual void run(const edm::EventSetup &setup, const EventInfo &eventInfo) = 0;
148 
150  virtual void beginRun(const edm::Run &, const edm::EventSetup &, bool changed){};
151 
153  virtual void endRun(const EndRunInfo &runInfo, const edm::EventSetup &setup){};
154 
156  virtual void beginLuminosityBlock(const edm::EventSetup &setup){};
157 
159  virtual void endLuminosityBlock(const edm::EventSetup &setup){};
160 
163  virtual bool setParametersForRunRange(const RunRange &rr) { return false; };
164 };
165 
166 /*** Global typedefs part II ***/
169 
170 #endif
virtual bool processesEvents()
Returns whether algorithm proccesses events in current configuration.
virtual void beginLuminosityBlock(const edm::EventSetup &setup)
called at begin of luminosity block (no lumi block info passed yet)
static AlgebraicMatrix initialize()
std::vector< ConstTrajTrackPair > ConstTrajTrackPairs
const TsosVectorCollection * tkLasBeamTsoses_
might be null!
virtual void endLuminosityBlock(const edm::EventSetup &setup)
called at end of luminosity block (no lumi block info passed yet)
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
EndRunInfo(const edm::RunID &theRunId, const TkFittedLasBeamCollection *theTkLasBeams, const TsosVectorCollection *theTkLasBeamTsoses)
const edm::EventID eventId() const
std::pair< const Trajectory *, const reco::Track * > ConstTrajTrackPair
const ConstTrajTrackPairCollection & trajTrackPairs() const
define event information passed to algorithms
virtual void endRun(const EndRunInfo &runInfo, const edm::EventSetup &setup)
called at end of run - order of arguments like in EDProducer etc.
virtual bool addCalibrations(const Calibrations &)
virtual bool setParametersForRunRange(const RunRange &rr)
virtual ~AlignmentAlgorithmBase()
Destructor.
std::vector< IntegratedCalibrationBase * > Calibrations
virtual void beginRun(const edm::Run &, const edm::EventSetup &, bool changed)
called at begin of run
AlignmentAlgorithmBase(const edm::ParameterSet &)
Constructor.
const TkFittedLasBeamCollection * tkLasBeams() const
AlignmentAlgorithmBase::EndRunInfo EndRunInfo
const TkFittedLasBeamCollection * tkLasBeams_
const edm::RunID runId_
might be null!
const edm::EventID eventId_
might be null!
std::pair< RunNumber, RunNumber > RunRange
Definition: Utilities.h:39
const AliClusterValueMap * clusterValueMap_
AlignmentAlgorithmBase::EventInfo EventInfo
fixed size matrix
virtual void terminate()
Called at end of job (must be implemented in derived class)
HLT enums.
const AliClusterValueMap * clusterValueMap() const
std::vector< TkFittedLasBeam > TkFittedLasBeamCollection
const reco::BeamSpot & beamSpot() const
EventInfo(const edm::EventID &theEventId, const ConstTrajTrackPairCollection &theTrajTrackPairs, const reco::BeamSpot &theBeamSpot, const AliClusterValueMap *theClusterValueMap)
virtual bool storeAlignments()
Returns whether algorithm produced results to be stored.
std::pair< const Trajectory *, const reco::Track * > ConstTrajTrackPair
eventInfo
add run, event number and lumi section
std::vector< std::vector< TrajectoryStateOnSurface > > TsosVectorCollection
define run information passed to algorithms (in endRun)
const TsosVectorCollection * tkLasBeamTsoses() const
might be null!
Constructor of the full muon geometry.
Definition: AlignableMuon.h:37
const ConstTrajTrackPairCollection & trajTrackPairs_
Definition: Run.h:45
cond::RealTimeType< cond::runnumber >::type RunNumber
Definition: Utilities.h:38
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection