CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Alignment/MuonAlignmentAlgorithms/interface/MuonChamberResidual.h

Go to the documentation of this file.
00001 #ifndef Alignment_MuonAlignmentAlgorithms_MuonChamberResidual_H
00002 #define Alignment_MuonAlignmentAlgorithms_MuonChamberResidual_H
00003 
00010 #include "FWCore/Framework/interface/ESHandle.h"
00011 #include "Alignment/CommonAlignment/interface/AlignableNavigator.h"
00012 #include "Alignment/CommonAlignment/interface/Alignable.h"
00013 #include "Alignment/CommonAlignment/interface/AlignmentParameters.h"
00014 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00015 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00016 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00017 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00018 #include "DataFormats/DetId/interface/DetId.h"
00019 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
00020 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
00021 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
00022 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00023 
00024 class MuonChamberResidual {
00025 public:
00026   MuonChamberResidual(edm::ESHandle<GlobalTrackingGeometry> globalGeometry, AlignableNavigator *navigator, DetId chamberId, AlignableDetOrUnitPtr chamberAlignable)
00027     : m_globalGeometry(globalGeometry)
00028     , m_navigator(navigator)
00029     , m_chamberId(chamberId)
00030     , m_chamberAlignable(chamberAlignable)
00031     , m_numHits(0)
00032     , m_residual_1(0.)
00033     , m_residual_x(0.)
00034     , m_residual_y(0.)
00035     , m_residual_xx(0.)
00036     , m_residual_xy(0.)
00037     , m_trackx_1(0.)
00038     , m_trackx_x(0.)
00039     , m_trackx_y(0.)
00040     , m_trackx_xx(0.)
00041     , m_trackx_xy(0.)
00042     , m_tracky_1(0.)
00043     , m_tracky_x(0.)
00044     , m_tracky_y(0.)
00045     , m_tracky_xx(0.)
00046     , m_tracky_xy(0.)
00047   {};
00048 
00049   virtual ~MuonChamberResidual() {};
00050 
00051   enum {
00052     kDT13,
00053     kDT2,
00054     kCSC
00055   };
00056 
00057   virtual void addResidual(const TrajectoryStateOnSurface *tsos, const TransientTrackingRecHit *hit) = 0;
00058   virtual double signConvention(const unsigned int rawId=0) const = 0;
00059 
00060   DetId chamberId() const { return m_chamberId; };
00061   AlignableDetOrUnitPtr chamberAlignable() const { return m_chamberAlignable; };
00062   virtual int type() const = 0;
00063 
00064   int numHits() const { return m_numHits; };
00065 
00066   double residual() const {
00067     assert(m_numHits > 1);
00068     double delta = m_residual_1*m_residual_xx - m_residual_x*m_residual_x;
00069     return (m_residual_xx*m_residual_y - m_residual_x*m_residual_xy) / delta;
00070   };
00071 
00072   double residual_error() const {
00073     assert(m_numHits > 1);
00074     double delta = m_residual_1*m_residual_xx - m_residual_x*m_residual_x;
00075     return sqrt(m_residual_xx / delta);
00076   };
00077 
00078   double resslope() const {
00079     assert(m_numHits > 1);
00080     double delta = m_residual_1*m_residual_xx - m_residual_x*m_residual_x;
00081     return (m_residual_1*m_residual_xy - m_residual_x*m_residual_y) / delta;
00082   };
00083 
00084   double resslope_error() const {
00085     assert(m_numHits > 1);
00086     double delta = m_residual_1*m_residual_xx - m_residual_x*m_residual_x;
00087     return sqrt(m_residual_1 / delta);
00088   };
00089 
00090   double chi2() const {
00091     double output = 0.;
00092     double a = residual();
00093     double b = resslope();
00094 
00095     std::vector<double>::const_iterator x = m_individual_x.begin();
00096     std::vector<double>::const_iterator y = m_individual_y.begin();
00097     std::vector<double>::const_iterator w = m_individual_weight.begin();
00098     for (;  x != m_individual_x.end();  ++x, ++y, ++w) {
00099       output += pow((*y) - a - b*(*x), 2) * (*w);
00100     }
00101     return output;
00102   };
00103 
00104   int ndof() const {
00105     return m_individual_x.size() - 2;
00106   };
00107 
00108   double trackdxdz() const {
00109     assert(m_numHits > 0);
00110     double delta = m_trackx_1*m_trackx_xx - m_trackx_x*m_trackx_x;
00111     return (m_trackx_1*m_trackx_xy - m_trackx_x*m_trackx_y) / delta;
00112   };
00113 
00114   double trackdydz() const {
00115     assert(m_numHits > 0);
00116     double delta = m_tracky_1*m_tracky_xx - m_tracky_x*m_tracky_x;
00117     return (m_tracky_1*m_tracky_xy - m_tracky_x*m_tracky_y) / delta;
00118   };
00119 
00120   double trackx() const {
00121     assert(m_numHits > 0);
00122     double delta = m_trackx_1*m_trackx_xx - m_trackx_x*m_trackx_x;
00123     return (m_trackx_xx*m_trackx_y - m_trackx_x*m_trackx_xy) / delta;
00124   };
00125 
00126   double tracky() const {
00127     assert(m_numHits > 0);
00128     double delta = m_tracky_1*m_tracky_xx - m_tracky_x*m_tracky_x;
00129     return (m_tracky_xx*m_tracky_y - m_tracky_x*m_tracky_xy) / delta;
00130   };
00131 
00132   GlobalPoint global_trackpos() {
00133     return chamberAlignable()->surface().toGlobal(LocalPoint(trackx(), tracky(), 0.));
00134   };
00135 
00136   double hitresid(int i) const {
00137     assert(0 <= i  &&  i < int(m_localIDs.size()));
00138     return m_localResids[i];
00139   }
00140 
00141   double global_residual() const {
00142     return residual() * signConvention();
00143   };
00144 
00145   double global_resslope() const {
00146     return resslope() * signConvention();
00147   };
00148 
00149   double global_hitresid(int i) const {
00150     return hitresid(i) * signConvention(m_localIDs[i].rawId());
00151   };
00152 
00153   int hitlayer(int i) const {  // only difference between DTs and CSCs is the DetId subclass
00154     assert(0 <= i  &&  i < int(m_localIDs.size()));
00155     if (m_chamberId.subdetId() == MuonSubdetId::DT) {
00156       DTLayerId layerId(m_localIDs[i].rawId());
00157       return 4*(layerId.superlayer() - 1) + layerId.layer();
00158     }
00159     else if (m_chamberId.subdetId() == MuonSubdetId::CSC) {
00160       CSCDetId layerId(m_localIDs[i].rawId());
00161       return layerId.layer();
00162     }
00163     else assert(false);
00164   };
00165 
00166   double hitposition(int i) const {
00167     assert(0 <= i  &&  i < int(m_localIDs.size()));
00168     if (m_chamberId.subdetId() == MuonSubdetId::DT) {
00169       GlobalPoint pos = m_globalGeometry->idToDet(m_localIDs[i])->position();
00170       return sqrt(pow(pos.x(), 2) + pow(pos.y(), 2));                   // R for DTs
00171     }
00172     else if (m_chamberId.subdetId() == MuonSubdetId::CSC) {
00173       return m_globalGeometry->idToDet(m_localIDs[i])->position().z();  // Z for CSCs
00174     }
00175     else assert(false);
00176   };
00177 
00178   DetId localid(int i) const {
00179     return m_localIDs[i];
00180   };
00181 
00182 protected:
00183   edm::ESHandle<GlobalTrackingGeometry> m_globalGeometry;
00184   AlignableNavigator *m_navigator;
00185   DetId m_chamberId;
00186   AlignableDetOrUnitPtr m_chamberAlignable;
00187 
00188   int m_numHits;
00189   double m_residual_1;
00190   double m_residual_x;
00191   double m_residual_y;
00192   double m_residual_xx;
00193   double m_residual_xy;
00194   double m_trackx_1;
00195   double m_trackx_x;
00196   double m_trackx_y;
00197   double m_trackx_xx;
00198   double m_trackx_xy;
00199   double m_tracky_1;
00200   double m_tracky_x;
00201   double m_tracky_y;
00202   double m_tracky_xx;
00203   double m_tracky_xy;
00204   std::vector<DetId> m_localIDs;
00205   std::vector<double> m_localResids;
00206   std::vector<double> m_individual_x;
00207   std::vector<double> m_individual_y;
00208   std::vector<double> m_individual_weight;
00209 };
00210 
00211 #endif // Alignment_MuonAlignmentAlgorithms_MuonChamberResidual_H