CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HIPAlignmentAlgorithm.h
Go to the documentation of this file.
1 #ifndef Alignment_HIPAlignmentAlgorithm_HIPAlignmentAlgorithm_h
2 #define Alignment_HIPAlignmentAlgorithm_HIPAlignmentAlgorithm_h
3 
4 #include <vector>
13 #include "Riostream.h"
14 
21 
24 
28 
31 #include "TFormula.h"
32 
33 class TFile;
34 class TTree;
35 
37 public:
40 
42  ~HIPAlignmentAlgorithm() override{};
43 
45  void initialize(const edm::EventSetup& setup,
46  AlignableTracker* tracker,
48  AlignableExtras* extras,
49  AlignmentParameterStore* store) override;
50 
52  void terminate(const edm::EventSetup& setup) override;
53 
55  void startNewLoop(void) override;
56 
58  void run(const edm::EventSetup& setup, const EventInfo& eventInfo) override;
59 
60 private:
61  bool processHit1D(const AlignableDetOrUnitPtr& alidet,
62  const Alignable* ali,
63  const HIPAlignableSpecificParameters* alispecifics,
64  const TrajectoryStateOnSurface& tsos,
65  const TrackingRecHit* hit,
66  double hitwt);
67 
68  bool processHit2D(const AlignableDetOrUnitPtr& alidet,
69  const Alignable* ali,
70  const HIPAlignableSpecificParameters* alispecifics,
71  const TrajectoryStateOnSurface& tsos,
72  const TrackingRecHit* hit,
73  double hitwt);
74 
76  void writeIterationFile(std::string filename, int iter);
77  void setAlignmentPositionError(void);
78  double calcAPE(double* par, int iter, int function);
79  void bookRoot(void);
81  bool calcParameters(Alignable* ali, int setDet, double start, double step);
82  void collector(void);
83  void collectMonitorTrees(const std::vector<std::string>& filenames);
84 
86 
87  // private data members
90 
91  std::unique_ptr<AlignableObjectId> alignableObjectId_;
94  std::unique_ptr<AlignableNavigator> theAlignableDetAccessor;
95 
97  int ioerr;
99 
100  // steering parameters
101 
102  // verbosity flag
103  const bool verbose;
104  // Monitor configuration
107  // names of IO root files
110 
112  std::vector<unsigned> theIOVrangeSet;
113 
114  // alignment position error parameters
116  std::vector<edm::ParameterSet> theAPEParameterSet;
117  std::vector<std::pair<align::Alignables, std::vector<double> > > theAPEParameters;
118 
119  // Default alignment specifications
120  // - min number of hits on alignable to calc parameters
121  // - max allowed rel error on parameter (else not used)
122  // - max allowed pull (residual / uncertainty) on a hit used in alignment
124 
126  std::vector<edm::ParameterSet> theCutsPerComponent;
127  std::vector<HIPAlignableSpecificParameters> theAlignableSpecifics;
128 
129  // collector mode (parallel processing)
133  int theDataGroup; // The data type specified in the cfg
137  std::vector<double> SetScanDet;
138 
139  std::unique_ptr<TFormula> theEtaFormula;
140 
141  const std::vector<std::string> surveyResiduals_;
142  std::vector<align::StructureType> theLevels; // for survey residuals
143 
144  // root tree variables
146  TTree* theTrackMonitorTree; // event-wise tree
147  TTree* theHitMonitorTree; // hit-wise tree
149  TTree* theAlignablesMonitorTree; // alignable-wise tree
151  TTree* theSurveyTree; // survey tree
152 
153  // common variables for monitor trees
155 
156  // variables for alignable-wise tree
163  unsigned int m2_nsurfdef;
164  std::vector<float> m2_surfDef;
165 
166  // variables for survey tree
169  float m3_par[6];
170 };
171 
172 #endif
std::vector< HIPAlignableSpecificParameters > theAlignableSpecifics
HIPAlignableSpecificParameters defaultAlignableSpecs
tuple cfg
Definition: looper.py:296
AlignmentParameterStore * theAlignmentParameterStore
HIPAlignableSpecificParameters * findAlignableSpecs(const Alignable *ali)
void run(const edm::EventSetup &setup, const EventInfo &eventInfo) override
Run the algorithm.
void initialize(const edm::EventSetup &setup, AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras, AlignmentParameterStore *store) override
Call at beginning of job.
uint32_t ID
Definition: Definitions.h:24
align::StructureType m2_ObjId
void writeIterationFile(std::string filename, int iter)
~HIPAlignmentAlgorithm() override
Destructor.
Interface/Base class for alignment algorithms, each alignment algorithm has to be derived from this c...
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken2_
std::unique_ptr< AlignableNavigator > theAlignableDetAccessor
define event information passed to algorithms
std::vector< unsigned > theIOVrangeSet
bool calcParameters(Alignable *ali, int setDet, double start, double step)
std::vector< edm::ParameterSet > theCutsPerComponent
int readIterationFile(std::string filename)
std::vector< double > SetScanDet
void collectMonitorTrees(const std::vector< std::string > &filenames)
double calcAPE(double *par, int iter, int function)
std::unique_ptr< AlignableObjectId > alignableObjectId_
HIPAlignmentAlgorithm(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
Constructor.
std::vector< edm::ParameterSet > theAPEParameterSet
align::StructureType m3_ObjId
std::unique_ptr< TFormula > theEtaFormula
std::vector< Alignable * > Alignables
Definition: Utilities.h:31
const edm::ESGetToken< TrackerTopology, IdealGeometryRecord > topoToken_
align::Alignables theAlignables
virtual void terminate()
Called at end of job (must be implemented in derived class)
std::vector< float > m2_surfDef
HIPMonitorConfig theMonitorConfig
void fillAlignablesMonitor(const edm::EventSetup &setup)
tuple filename
Definition: lut2db_cfg.py:20
void startNewLoop(void) override
Called at start of new loop.
const std::vector< std::string > surveyResiduals_
step
Definition: StallMonitor.cc:98
bool processHit2D(const AlignableDetOrUnitPtr &alidet, const Alignable *ali, const HIPAlignableSpecificParameters *alispecifics, const TrajectoryStateOnSurface &tsos, const TrackingRecHit *hit, double hitwt)
bool processHit1D(const AlignableDetOrUnitPtr &alidet, const Alignable *ali, const HIPAlignableSpecificParameters *alispecifics, const TrajectoryStateOnSurface &tsos, const TrackingRecHit *hit, double hitwt)
std::vector< std::pair< align::Alignables, std::vector< double > > > theAPEParameters
Constructor of the full muon geometry.
Definition: AlignableMuon.h:38
std::vector< align::StructureType > theLevels
SurfaceDeformationFactory::Type m2_dtype