00001
00002
00003
00004
00005
00006
00007 #include "Validation/TrackerRecHits/interface/SiStripRecHitsValid.h"
00008 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00009
00010
00011 #include "DataFormats/DetId/interface/DetId.h"
00012 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00013 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00014 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00015 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00016 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00017 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00018 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00019
00020
00021
00022 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00023 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
00024 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00025 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00026 #include "DataFormats/Common/interface/OwnVector.h"
00027 #include "DQMServices/Core/interface/DQMStore.h"
00028
00029 using namespace std;
00030 using namespace edm;
00031
00032 namespace helper {
00033 struct GetDetId {
00034 template<typename X>
00035 DetId operator()(const X &x) { return DetId(x.detId()); }
00036 };
00037
00038 template<typename T>
00039 std::pair<typename T::DetSet::const_iterator, typename T::DetSet::const_iterator>
00040 getRange(const T &detset, const DetId &id) {
00041 typedef std::pair<typename T::DetSet::const_iterator, typename T::DetSet::const_iterator> return_type;
00042 typename T::const_iterator match = detset.find(id);
00043 if (match == detset.end()) return return_type();
00044 typename T::DetSet hits = *match;
00045 return return_type(hits.begin(), hits.end());
00046 }
00047 }
00048
00049
00050
00051
00052 SiStripRecHitsValid::SiStripRecHitsValid(const ParameterSet& ps) :
00053 dbe_(0),
00054 conf_(ps),
00055 matchedRecHits_( ps.getParameter<edm::InputTag>("matchedRecHits") ),
00056 rphiRecHits_( ps.getParameter<edm::InputTag>("rphiRecHits") ),
00057 stereoRecHits_( ps.getParameter<edm::InputTag>("stereoRecHits") ) {
00058
00059 outputFile_ = ps.getUntrackedParameter<string>("outputFile", "sistriprechitshisto.root");
00060 dbe_ = Service<DQMStore>().operator->();
00061
00062 dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Strip/SISTRIP");
00063
00064 meNumTotRphi = dbe_->book1D("NumTotRphi","Num of RecHits rphi",100, 0, 10000);
00065 meNumTotSas = dbe_->book1D("NumTotSas","Num of RecHits sas",100, 0, 10000);
00066 meNumTotMatched = dbe_->book1D("NumTotMatched","Num of RecHits rmatched",100, 0, 10000);
00067
00068 dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Strip/TIB");
00069 meNumRphiTIB = dbe_->book1D("NumRphiTIB","Num of RecHits rphi", 100, 0, 1000.);
00070 meNumSasTIB = dbe_->book1D("NumSasTIB","Num of RecHits sas", 100, 0, 1000.);
00071 meNumMatchedTIB = dbe_->book1D("NumMatchedTIB","Num of RecHits matched", 100, 0, 1000.);
00072
00073
00074 for(int i = 0 ;i<4 ; i++) {
00075 Char_t histo[200];
00076 sprintf(histo,"Nstp_rphi_layer%dtib",i+1);
00077 meNstpRphiTIB[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00078 sprintf(histo,"Adc_rphi_layer%dtib",i+1);
00079 meAdcRphiTIB[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00080 sprintf(histo,"Posx_rphi_layer%dtib",i+1);
00081 mePosxRphiTIB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00082 sprintf(histo,"Errx_rphi_layer%dtib",i+1);
00083 meErrxRphiTIB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0,0.01);
00084 sprintf(histo,"Res_rphi_layer%dtib",i+1);
00085 meResRphiTIB[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.02,+0.02);
00086 sprintf(histo,"Pull_LF_rphi_layer%dtib",i+1);
00087 mePullLFRphiTIB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00088 sprintf(histo,"Pull_MF_rphi_layer%dtib",i+1);
00089 mePullMFRphiTIB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00090 sprintf(histo,"Chi2_rphi_layer%dtib",i+1);
00091 meChi2RphiTIB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00092 }
00093
00094
00095 for(int i = 0 ;i<2 ; i++) {
00096 Char_t histo[200];
00097 sprintf(histo,"Nstp_sas_layer%dtib",i+1);
00098 meNstpSasTIB[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00099 sprintf(histo,"Adc_sas_layer%dtib",i+1);
00100 meAdcSasTIB[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00101 sprintf(histo,"Posx_sas_layer%dtib",i+1);
00102 mePosxSasTIB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00103 sprintf(histo,"Errx_sas_layer%dtib",i+1);
00104 meErrxSasTIB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0.,0.01);
00105 sprintf(histo,"Res_sas_layer%dtib",i+1);
00106 meResSasTIB[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.02,+0.02);
00107 sprintf(histo,"Pull_LF_sas_layer%dtib",i+1);
00108 mePullLFSasTIB[i] = dbe_->book1D(histo,"Pull",100,-4.,4.);
00109 sprintf(histo,"Pull_MF_sas_layer%dtib",i+1);
00110 mePullMFSasTIB[i] = dbe_->book1D(histo,"Pull",100,-4.,4.);
00111 sprintf(histo,"Chi2_sas_layer%dtib",i+1);
00112 meChi2SasTIB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00113
00114 sprintf(histo,"Posx_matched_layer%dtib",i+1);
00115 mePosxMatchedTIB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0, +6.0);
00116 sprintf(histo,"Posy_matched_layer%dtib",i+1);
00117 mePosyMatchedTIB[i] = dbe_->book1D(histo,"RecHit y coord.",100,-6.0, +6.0);
00118 sprintf(histo,"Errx_matched_layer%dtib",i+1);
00119 meErrxMatchedTIB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0., 0.01);
00120 sprintf(histo,"Erry_matched_layer%dtib",i+1);
00121 meErryMatchedTIB[i] = dbe_->book1D(histo,"RecHit err(y) coord.",100,0., 0.05);
00122 sprintf(histo,"Resx_matched_layer%dtib",i+1);
00123 meResxMatchedTIB[i] = dbe_->book1D(histo,"RecHit Res(x) coord.",100,-0.02, +0.02);
00124 sprintf(histo,"Resy_matched_layer%dtib",i+1);
00125 meResyMatchedTIB[i] = dbe_->book1D(histo,"RecHit Res(y) coord.",100,-1., +1.);
00126 sprintf(histo,"Chi2_matched_layer%dtib",i+1);
00127 meChi2MatchedTIB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00128 }
00129
00130 dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Strip/TOB");
00131 meNumRphiTOB = dbe_->book1D("NumRphiTOB","Num of RecHits rphi", 100, 0, 1000.);
00132 meNumSasTOB = dbe_->book1D("NumSasTOB","Num of RecHits sas", 100, 0, 1000.);
00133 meNumMatchedTOB = dbe_->book1D("NumMatchedTOB","Num of RecHits matched", 100, 0, 1000.);
00134
00135 for(int i = 0 ;i<6 ; i++) {
00136 Char_t histo[200];
00137 sprintf(histo,"Nstp_rphi_layer%dtob",i+1);
00138 meNstpRphiTOB[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00139 sprintf(histo,"Adc_rphi_layer%dtob",i+1);
00140 meAdcRphiTOB[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00141 sprintf(histo,"Posx_rphi_layer%dtob",i+1);
00142 mePosxRphiTOB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00143 sprintf(histo,"Errx_rphi_layer%dtob",i+1);
00144 meErrxRphiTOB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0,0.01);
00145 sprintf(histo,"Res_rphi_layer%dtob",i+1);
00146 meResRphiTOB[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.02,+0.02);
00147 sprintf(histo,"Pull_LF_rphi_layer%dtob",i+1);
00148 mePullLFRphiTOB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00149 sprintf(histo,"Pull_MF_rphi_layer%dtob",i+1);
00150 mePullMFRphiTOB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00151 sprintf(histo,"Chi2_rphi_layer%dtob",i+1);
00152 meChi2RphiTOB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00153 }
00154
00155
00156 for(int i = 0 ;i<2 ; i++) {
00157 Char_t histo[200];
00158 sprintf(histo,"Nstp_sas_layer%dtob",i+1);
00159 meNstpSasTOB[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00160 sprintf(histo,"Adc_sas_layer%dtob",i+1);
00161 meAdcSasTOB[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00162 sprintf(histo,"Posx_sas_layer%dtob",i+1);
00163 mePosxSasTOB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00164 sprintf(histo,"Errx_sas_layer%dtob",i+1);
00165 meErrxSasTOB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0.,0.01);
00166 sprintf(histo,"Res_sas_layer%dtob",i+1);
00167 meResSasTOB[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.02,+0.02);
00168 sprintf(histo,"Pull_LF_sas_layer%dtob",i+1);
00169 mePullLFSasTOB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00170 sprintf(histo,"Pull_MF_sas_layer%dtob",i+1);
00171 mePullMFSasTOB[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00172 sprintf(histo,"Chi2_sas_layer%dtob",i+1);
00173 meChi2SasTOB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00174
00175 sprintf(histo,"Posx_matched_layer%dtob",i+1);
00176 mePosxMatchedTOB[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0, +6.0);
00177 sprintf(histo,"Posy_matched_layer%dtob",i+1);
00178 mePosyMatchedTOB[i] = dbe_->book1D(histo,"RecHit y coord.",100,-10.0, +10.0);
00179 sprintf(histo,"Errx_matched_layer%dtob",i+1);
00180 meErrxMatchedTOB[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0., 0.01);
00181 sprintf(histo,"Erry_matched_layer%dtob",i+1);
00182 meErryMatchedTOB[i] = dbe_->book1D(histo,"RecHit err(y) coord.",100,0., 0.1);
00183 sprintf(histo,"Resx_matched_layer%dtob",i+1);
00184 meResxMatchedTOB[i] = dbe_->book1D(histo,"RecHit Res(x) coord.",100,-0.02, +0.02);
00185 sprintf(histo,"Resy_matched_layer%dtob",i+1);
00186 meResyMatchedTOB[i] = dbe_->book1D(histo,"RecHit Res(y) coord.",100,-1., +1.);
00187 sprintf(histo,"Chi2_matched_layer%dtob",i+1);
00188 meChi2MatchedTOB[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00189 }
00190
00191 dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Strip/TID");
00192 meNumRphiTID = dbe_->book1D("NumRphiTID","Num of RecHits rphi", 100, 0, 1000.);
00193 meNumSasTID = dbe_->book1D("NumSasTID","Num of RecHits sas", 100, 0, 1000.);
00194 meNumMatchedTID = dbe_->book1D("NumMatchedTID","Num of RecHits matched", 100, 0, 1000.);
00195
00196
00197 for(int i = 0 ;i<3 ; i++) {
00198 Char_t histo[200];
00199 sprintf(histo,"Nstp_rphi_layer%dtid",i+1);
00200 meNstpRphiTID[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00201 sprintf(histo,"Adc_rphi_layer%dtid",i+1);
00202 meAdcRphiTID[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00203 sprintf(histo,"Posx_rphi_layer%dtid",i+1);
00204 mePosxRphiTID[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00205 sprintf(histo,"Errx_rphi_layer%dtid",i+1);
00206 meErrxRphiTID[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0,0.5);
00207 sprintf(histo,"Res_rphi_layer%dtid",i+1);
00208 meResRphiTID[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.5,+0.5);
00209 sprintf(histo,"Pull_LF_rphi_layer%dtid",i+1);
00210 mePullLFRphiTID[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00211 sprintf(histo,"Pull_MF_rphi_layer%dtid",i+1);
00212 mePullMFRphiTID[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00213 sprintf(histo,"Chi2_rphi_layer%dtid",i+1);
00214 meChi2RphiTID[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00215 }
00216
00217
00218 for(int i = 0 ;i<2 ; i++) {
00219 Char_t histo[200];
00220 sprintf(histo,"Nstp_sas_layer%dtid",i+1);
00221 meNstpSasTID[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00222 sprintf(histo,"Adc_sas_layer%dtid",i+1);
00223 meAdcSasTID[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00224 sprintf(histo,"Posx_sas_layer%dtid",i+1);
00225 mePosxSasTID[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00226 sprintf(histo,"Errx_sas_layer%dtid",i+1);
00227 meErrxSasTID[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0.,0.5);
00228 sprintf(histo,"Res_sas_layer%dtid",i+1);
00229 meResSasTID[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.5,+0.5);
00230 sprintf(histo,"Pull_LF_sas_layer%dtid",i+1);
00231 mePullLFSasTID[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00232 sprintf(histo,"Pull_MF_sas_layer%dtid",i+1);
00233 mePullMFSasTID[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00234 sprintf(histo,"Chi2_sas_layer%dtid",i+1);
00235 meChi2SasTID[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00236
00237 sprintf(histo,"Posx_matched_layer%dtid",i+1);
00238 mePosxMatchedTID[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0, +6.0);
00239 sprintf(histo,"Posy_matched_layer%dtid",i+1);
00240 mePosyMatchedTID[i] = dbe_->book1D(histo,"RecHit y coord.",100,-6.0, +6.0);
00241 sprintf(histo,"Errx_matched_layer%dtid",i+1);
00242 meErrxMatchedTID[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0., 0.05);
00243 sprintf(histo,"Erry_matched_layer%dtid",i+1);
00244 meErryMatchedTID[i] = dbe_->book1D(histo,"RecHit err(y) coord.",100,0., 0.1);
00245 sprintf(histo,"Resx_matched_layer%dtid",i+1);
00246 meResxMatchedTID[i] = dbe_->book1D(histo,"RecHit Res(x) coord.",100,-0.2, +0.2);
00247 sprintf(histo,"Resy_matched_layer%dtid",i+1);
00248 meResyMatchedTID[i] = dbe_->book1D(histo,"RecHit Res(y) coord.",100,-1., +1.);
00249 sprintf(histo,"Chi2_matched_layer%dtid",i+1);
00250 meChi2MatchedTID[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00251 }
00252
00253 dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Strip/TEC");
00254 meNumRphiTEC = dbe_->book1D("NumRphiTEC","Num of RecHits rphi", 100, 0, 1000.);
00255 meNumSasTEC = dbe_->book1D("NumSasTEC","Num of RecHits sas", 100, 0, 1000.);
00256 meNumMatchedTEC = dbe_->book1D("NumMatchedTEC","Num of RecHits matched", 100, 0, 1000.);
00257
00258
00259 for(int i = 0 ;i<7 ; i++) {
00260 Char_t histo[200];
00261 sprintf(histo,"Nstp_rphi_layer%dtec",i+1);
00262 meNstpRphiTEC[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00263 sprintf(histo,"Adc_rphi_layer%dtec",i+1);
00264 meAdcRphiTEC[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00265 sprintf(histo,"Posx_rphi_layer%dtec",i+1);
00266 mePosxRphiTEC[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00267 sprintf(histo,"Errx_rphi_layer%dtec",i+1);
00268 meErrxRphiTEC[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0,0.5);
00269 sprintf(histo,"Res_rphi_layer%dtec",i+1);
00270 meResRphiTEC[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.5,+0.5);
00271 sprintf(histo,"Pull_LF_rphi_layer%dtec",i+1);
00272 mePullLFRphiTEC[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00273 sprintf(histo,"Pull_MF_rphi_layer%dtec",i+1);
00274 mePullMFRphiTEC[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00275 sprintf(histo,"Chi2_rphi_layer%dtec",i+1);
00276 meChi2RphiTEC[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00277 }
00278
00279
00280 for(int i = 0 ;i<5 ; i++) {
00281 if(i == 0 || i == 1 || i == 4) {
00282 Char_t histo[200];
00283 sprintf(histo,"Nstp_sas_layer%dtec",i+1);
00284 meNstpSasTEC[i] = dbe_->book1D(histo,"RecHit Cluster Size",10,0.5,10.5);
00285 sprintf(histo,"Adc_sas_layer%dtec",i+1);
00286 meAdcSasTEC[i] = dbe_->book1D(histo,"RecHit Cluster Charge",100,0.,300.);
00287 sprintf(histo,"Posx_sas_layer%dtec",i+1);
00288 mePosxSasTEC[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0,+6.0);
00289 sprintf(histo,"Errx_sas_layer%dtec",i+1);
00290 meErrxSasTEC[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0.,0.5);
00291 sprintf(histo,"Res_sas_layer%dtec",i+1);
00292 meResSasTEC[i] = dbe_->book1D(histo,"RecHit Residual",100,-0.5,+0.5);
00293 sprintf(histo,"Pull_LF_sas_layer%dtec",i+1);
00294 mePullLFSasTEC[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00295 sprintf(histo,"Pull_MF_sas_layer%dtec",i+1);
00296 mePullMFSasTEC[i] = dbe_->book1D(histo,"Pull",100,-5.,5.);
00297 sprintf(histo,"Chi2_sas_layer%dtec",i+1);
00298 meChi2SasTEC[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00299
00300 sprintf(histo,"Posx_matched_layer%dtec",i+1);
00301 mePosxMatchedTEC[i] = dbe_->book1D(histo,"RecHit x coord.",100,-6.0, +6.0);
00302 sprintf(histo,"Posy_matched_layer%dtec",i+1);
00303 mePosyMatchedTEC[i] = dbe_->book1D(histo,"RecHit y coord.",100,-8.0, +8.0);
00304 sprintf(histo,"Errx_matched_layer%dtec",i+1);
00305 meErrxMatchedTEC[i] = dbe_->book1D(histo,"RecHit err(x) coord.",100,0., 0.05);
00306 sprintf(histo,"Erry_matched_layer%dtec",i+1);
00307 meErryMatchedTEC[i] = dbe_->book1D(histo,"RecHit err(y) coord.",100,0., 0.1);
00308 sprintf(histo,"Resx_matched_layer%dtec",i+1);
00309 meResxMatchedTEC[i] = dbe_->book1D(histo,"RecHit Res(x) coord.",100,-0.2, +0.2);
00310 sprintf(histo,"Resy_matched_layer%dtec",i+1);
00311 meResyMatchedTEC[i] = dbe_->book1D(histo,"RecHit Res(y) coord.",100,-1., +1.);
00312 sprintf(histo,"Chi2_matched_layer%dtec",i+1);
00313 meChi2MatchedTEC[i] = dbe_->book1D(histo,"RecHit Chi2 test",100,0., 50);
00314 }
00315 }
00316 }
00317
00318
00319 SiStripRecHitsValid::~SiStripRecHitsValid(){
00320
00321 }
00322
00323 void SiStripRecHitsValid::beginJob(){
00324
00325 }
00326
00327 void SiStripRecHitsValid::endJob() {
00328 if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00329 }
00330
00331
00332 void SiStripRecHitsValid::analyze(const edm::Event& e, const edm::EventSetup& es) {
00333
00334 LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
00335
00336
00337
00338
00339 std::string rechitProducer = "SiStripRecHits2D";
00340
00341
00342 edm::Handle<SiStripMatchedRecHit2DCollection> rechitsmatched;
00343 edm::Handle<SiStripRecHit2DCollection> rechitsrphi;
00344 edm::Handle<SiStripRecHit2DCollection> rechitsstereo;
00345 e.getByLabel(matchedRecHits_, rechitsmatched);
00346 e.getByLabel(rphiRecHits_, rechitsrphi);
00347 e.getByLabel(stereoRecHits_, rechitsstereo);
00348
00349 int numrechitrphi =0;
00350 int numrechitsas =0;
00351 int numrechitmatched=0;
00352
00353 int totTibnumrechitrphi=0;
00354 int totTibnumrechitsas=0;
00355 int totTibnumrechitmatched=0;
00356 int totTobnumrechitrphi=0;
00357 int totTobnumrechitsas=0;
00358 int totTobnumrechitmatched=0;
00359 int totTidnumrechitrphi=0;
00360 int totTidnumrechitsas=0;
00361 int totTidnumrechitmatched=0;
00362 int totTecnumrechitrphi=0;
00363 int totTecnumrechitsas=0;
00364 int totTecnumrechitmatched=0;
00365 int totrechitrphi =0;
00366 int totrechitsas =0;
00367 int totrechitmatched =0;
00368
00369 TrackerHitAssociator associate(e, conf_);
00370
00371 edm::ESHandle<TrackerGeometry> pDD;
00372 es.get<TrackerDigiGeometryRecord> ().get (pDD);
00373 const TrackerGeometry &tracker(*pDD);
00374
00375 std::vector<DetId> IDs;
00376 IDs.reserve(rechitsrphi->size() + rechitsmatched->size() + rechitsstereo->size());
00377 std::transform(rechitsrphi->begin(), rechitsrphi->end(), std::back_inserter(IDs), helper::GetDetId() );
00378 std::transform(rechitsstereo->begin(), rechitsstereo->end(), std::back_inserter(IDs), helper::GetDetId() );
00379 std::transform(rechitsmatched->begin(), rechitsmatched->end(), std::back_inserter(IDs), helper::GetDetId() );
00380
00381
00382 for(std::vector<DetId>::const_iterator it = IDs.begin(); it != IDs.end(); ++it ){
00383 uint32_t myid=((*it).rawId());
00384 DetId detid = ((*it));
00385
00386
00387 for(int i=0; i<MAXHIT; i++){
00388 rechitrphix[i] =0;
00389 rechitrphierrx[i] =0;
00390 rechitrphiy[i] =0;
00391 rechitrphiz[i] =0;
00392 rechitsasx[i] =0;
00393 rechitsaserrx[i] =0;
00394 rechitsasy[i] =0;
00395 rechitsasz[i] =0;
00396 clusizrphi[i] =0;
00397 clusizsas[i] =0;
00398 cluchgrphi[i] =0;
00399 cluchgsas[i] =0;
00400 rechitrphires[i]=-999.;
00401 rechitsasres[i]=-999.;
00402 rechitrphipullMF[i]=-999.;
00403 rechitsaspullMF[i]=-999.;
00404 chi2rphi[i] =0;
00405 chi2sas[i]=0;
00406 rechitmatchedx[i] =0;
00407 rechitmatchedy[i] =0;
00408 rechitmatchedz[i] =0;
00409 rechitmatchederrxx[i] =0;
00410 rechitmatchederrxy[i] =0;
00411 rechitmatchederryy[i] =0;
00412 rechitmatchedresx[i]=-999;
00413 rechitmatchedresy[i]=-999;
00414 chi2matched[i]=0;
00415 }
00416
00417 numrechitrphi =0;
00418
00419 std::pair<SiStripRecHit2DCollection::DetSet::const_iterator,SiStripRecHit2DCollection::DetSet::const_iterator> rechitrphiRange = helper::getRange(*rechitsrphi, detid);
00420 SiStripRecHit2DCollection::DetSet::const_iterator rechitrphiRangeIteratorBegin = rechitrphiRange.first;
00421 SiStripRecHit2DCollection::DetSet::const_iterator rechitrphiRangeIteratorEnd = rechitrphiRange.second;
00422 SiStripRecHit2DCollection::DetSet::const_iterator iterrphi=rechitrphiRangeIteratorBegin;
00423
00424 numrechitrphi = rechitrphiRangeIteratorEnd - rechitrphiRangeIteratorBegin;
00425
00426
00427 int i=0;
00428 int i2=0;
00429
00430 if(numrechitrphi > 0 ){
00431 totrechitrphi+=numrechitrphi;
00432 for(iterrphi=rechitrphiRangeIteratorBegin; iterrphi!=rechitrphiRangeIteratorEnd;++iterrphi){
00433 const GeomDetUnit * det = tracker.idToDetUnit(detid);
00434 const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)(det);
00435 const StripTopology &topol=(StripTopology&)stripdet->topology();
00436 SiStripRecHit2D const rechit=*iterrphi;
00437 LocalPoint position=rechit.localPosition();
00438 LocalError error=rechit.localPositionError();
00439 MeasurementPoint Mposition;
00440 MeasurementError Merror;
00441 Mposition = topol.measurementPosition(position);
00442 Merror = topol.measurementError(position,error);
00443 SiStripRecHit2D::ClusterRef clust=rechit.cluster();
00444 int clusiz=0;
00445 int totcharge=0;
00446 clusiz = clust->amplitudes().size();
00447 const std::vector<uint8_t> amplitudes=clust->amplitudes();
00448 for(size_t ia=0; ia<amplitudes.size();ia++){
00449 totcharge+=amplitudes[ia];
00450 }
00451 rechitrphix[i] = position.x();
00452 rechitrphiy[i] = position.y();
00453 rechitrphiz[i] = position.z();
00454 rechitrphierrx[i] = error.xx();
00455 clusizrphi[i] = clusiz;
00456 cluchgrphi[i] = totcharge;
00457
00458 matched.clear();
00459 matched = associate.associateHit(rechit);
00460 float mindist = 999999;
00461 float dist = 999999;
00462 PSimHit closest;
00463 if(!matched.empty()){
00464 for(vector<PSimHit>::const_iterator m=matched.begin(); m<matched.end(); m++){
00465 dist = fabs(rechitrphix[i] - (*m).localPosition().x());
00466 if(dist<mindist){
00467 mindist = dist;
00468 closest = (*m);
00469 }
00470 rechitrphires[i] = rechitrphix[i] - closest.localPosition().x();
00471 }
00472 rechitrphipullMF[i] = (Mposition.x() - (topol.measurementPosition(closest.localPosition())).x())/sqrt(Merror.uu());
00473
00474
00475
00476 AlgebraicVector rhparameters = rechit.parameters();
00477 AlgebraicVector shparameters(2);
00478 shparameters[0] = closest.localPosition().x();
00479 shparameters[1] = closest.localPosition().y();
00480 AlgebraicVector r(rhparameters - shparameters);
00481 AlgebraicSymMatrix R = rechit.parametersError();
00482 int ierr;
00483 R.invert(ierr);
00484 double est = R.similarity(r);
00485
00486
00487
00488
00489
00490
00491 chi2rphi[i2] = est;
00492 i2++;
00493 }
00494 i++;
00495 }
00496 }
00497
00498
00499 int j=0;
00500 int j2=0;
00501 numrechitsas=0;
00502 std::pair<SiStripRecHit2DCollection::DetSet::const_iterator,SiStripRecHit2DCollection::DetSet::const_iterator> rechitsasRange = helper::getRange(*rechitsstereo, detid);
00503 SiStripRecHit2DCollection::DetSet::const_iterator rechitsasRangeIteratorBegin = rechitsasRange.first;
00504 SiStripRecHit2DCollection::DetSet::const_iterator rechitsasRangeIteratorEnd = rechitsasRange.second;
00505 SiStripRecHit2DCollection::DetSet::const_iterator itersas=rechitsasRangeIteratorBegin;
00506 numrechitsas = rechitsasRangeIteratorEnd - rechitsasRangeIteratorBegin;
00507 if(numrechitsas > 0){
00508 totrechitsas+=numrechitsas;
00509 for(itersas=rechitsasRangeIteratorBegin; itersas!=rechitsasRangeIteratorEnd;++itersas){
00510 const GeomDetUnit * det = tracker.idToDetUnit(detid);
00511 const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)(det);
00512 const StripTopology &topol=(StripTopology&)stripdet->topology();
00513 SiStripRecHit2D const rechit=*itersas;
00514 LocalPoint position=rechit.localPosition();
00515 LocalError error=rechit.localPositionError();
00516 MeasurementPoint Mposition;
00517 MeasurementError Merror;
00518 Mposition = topol.measurementPosition(position);
00519 Merror = topol.measurementError(position,error);
00520 SiStripRecHit2D::ClusterRef clust=rechit.cluster(); int clusiz=0;
00521 int totcharge=0;
00522 clusiz = clust->amplitudes().size();
00523 const std::vector<uint8_t> amplitudes=clust->amplitudes();
00524 for(size_t ia=0; ia<amplitudes.size();ia++){
00525 totcharge+=amplitudes[ia];
00526 }
00527
00528 rechitsasx[j] = position.x();
00529 rechitsasy[j] = position.y();
00530 rechitsasz[j] = position.z();
00531 rechitsaserrx[j] = error.xx();
00532 clusizsas[j] = clusiz;
00533 cluchgsas[j] = totcharge;
00534
00535 float mindist = 999999;
00536 float dist = 999999;
00537 PSimHit closest;
00538 matched.clear();
00539 matched = associate.associateHit(rechit);
00540 if(!matched.empty()){
00541 for(vector<PSimHit>::const_iterator m=matched.begin(); m<matched.end(); m++){
00542 dist = fabs(rechitsasx[j] - (*m).localPosition().x());
00543 if(dist<mindist){
00544 mindist = dist;
00545 closest = (*m);
00546 }
00547 rechitsasres[j] = rechitsasx[j] - closest.localPosition().x();
00548 }
00549 rechitsaspullMF[j] = (Mposition.x() - (topol.measurementPosition(closest.localPosition())).x())/sqrt(Merror.uu());
00550
00551
00552 AlgebraicVector rhparameters = rechit.parameters();
00553 AlgebraicVector shparameters(2);
00554 shparameters[0] = closest.localPosition().x();
00555 shparameters[1] = closest.localPosition().y();
00556 AlgebraicVector r(rhparameters - shparameters);
00557 AlgebraicSymMatrix R = rechit.parametersError();
00558 int ierr;
00559 R.invert(ierr);
00560 double est = R.similarity(r);
00561
00562
00563
00564
00565
00566
00567 chi2sas[j2] = est;
00568 j2++;
00569 }
00570 j++;
00571 }
00572 }
00573
00574
00575
00576 int k=0;
00577 int k2=0;
00578
00579
00580 numrechitmatched=0;
00581 std::pair<SiStripMatchedRecHit2DCollection::DetSet::const_iterator,SiStripMatchedRecHit2DCollection::DetSet::const_iterator> rechitmatchedRange = helper::getRange(*rechitsmatched, detid);
00582 SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorBegin = rechitmatchedRange.first;
00583 SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorEnd = rechitmatchedRange.second;
00584 SiStripMatchedRecHit2DCollection::DetSet::const_iterator itermatched=rechitmatchedRangeIteratorBegin;
00585 numrechitmatched = rechitmatchedRangeIteratorEnd - rechitmatchedRangeIteratorBegin;
00586 if(numrechitmatched > 0){
00587 totrechitmatched +=numrechitmatched;
00588
00589 for(itermatched=rechitmatchedRangeIteratorBegin; itermatched!=rechitmatchedRangeIteratorEnd;++itermatched){
00590 SiStripMatchedRecHit2D const rechit=*itermatched;
00591 LocalPoint position=rechit.localPosition();
00592 LocalError error=rechit.localPositionError();
00593
00594 float mindist = 999999;
00595 float distx = 999999;
00596 float disty = 999999;
00597 float dist = 999999;
00598 std::pair<LocalPoint,LocalVector> closestPair;
00599 matched.clear();
00600
00601
00602
00603
00604
00605 rechitmatchedx[k] = position.x();
00606 rechitmatchedy[k] = position.y();
00607 rechitmatchedz[k] = position.z();
00608 rechitmatchederrxx[k] = error.xx();
00609 rechitmatchederrxy[k] = error.xy();
00610 rechitmatchederryy[k] = error.yy();
00611
00612
00613 matched = associate.associateHit(rechit);
00614
00615
00616 if(!matched.empty()){
00617
00618 const GluedGeomDet* gluedDet = (const GluedGeomDet*)tracker.idToDet(rechit.geographicalId());
00619 const StripGeomDetUnit* partnerstripdet =(StripGeomDetUnit*) gluedDet->stereoDet();
00620 std::pair<LocalPoint,LocalVector> hitPair;
00621
00622
00623
00624 for(vector<PSimHit>::const_iterator m=matched.begin(); m<matched.end(); m++){
00625
00626 hitPair= projectHit((*m),partnerstripdet,gluedDet->surface());
00627 distx = fabs(rechitmatchedx[k] - hitPair.first.x());
00628 disty = fabs(rechitmatchedy[k] - hitPair.first.y());
00629 dist = sqrt(distx*distx+disty*disty);
00630
00631
00632 if(dist<mindist){
00633 mindist = dist;
00634 closestPair = hitPair;
00635 }
00636 }
00637
00638
00639 rechitmatchedresx[k] = rechitmatchedx[k] - closestPair.first.x();
00640 rechitmatchedresy[k] = rechitmatchedy[k] - closestPair.first.y();
00641
00642
00643
00644
00645 AlgebraicVector rhparameters = rechit.parameters();
00646 LocalPoint sh = closestPair.first;
00647 AlgebraicVector shparameters(2);
00648 shparameters[0] = sh.x();
00649 shparameters[1] = sh.y();
00650 AlgebraicVector r(rhparameters - shparameters);
00651 AlgebraicSymMatrix R = rechit.parametersError();
00652 int ierr;
00653 R.invert(ierr);
00654 double est = R.similarity(r);
00655
00656
00657
00658
00659
00660 chi2matched[k2] = est;
00661 k2++;
00662 }
00663
00664 k++;
00665 }
00666 }
00667
00668
00669 if(numrechitrphi>0 || numrechitsas>0 || numrechitmatched){
00670 if (detid.subdetId() == int(StripSubdetector::TIB)){
00671 TIBDetId tibid(myid);
00672 int Tibnumrechitrphi = numrechitrphi;
00673 int Tibnumrechitsas = numrechitsas;
00674 int Tibnumrechitmatched = numrechitmatched;
00675
00676 int Tibnumrechitrphi2 = i2;
00677 int Tibnumrechitsas2 = j2;
00678 int Tibnumrechitmatched2 = k2;
00679
00680 totTibnumrechitrphi +=numrechitrphi;
00681 totTibnumrechitsas +=numrechitsas;
00682 totTibnumrechitmatched +=numrechitmatched;
00683
00684 int ilay = tibid.layer() - 1;
00685
00686 if(tibid.stereo()==0){
00687 for(int k = 0; k<Tibnumrechitrphi; k++){
00688 meNstpRphiTIB[ilay]->Fill(clusizrphi[k]);
00689 meAdcRphiTIB[ilay]->Fill(cluchgrphi[k]);
00690 mePosxRphiTIB[ilay]->Fill(rechitrphix[k]);
00691 meErrxRphiTIB[ilay]->Fill(sqrt(rechitrphierrx[k]));
00692 meResRphiTIB[ilay]->Fill(rechitrphires[k]);
00693 mePullLFRphiTIB[ilay]->Fill(rechitrphires[k]/sqrt(rechitrphierrx[k]));
00694 mePullMFRphiTIB[ilay]->Fill(rechitrphipullMF[k]);
00695
00696 }
00697
00698 for(int k = 0; k<Tibnumrechitrphi2; k++){
00699 meChi2RphiTIB[ilay]->Fill(chi2rphi[k]);
00700 }
00701
00702 } else if(tibid.stereo()==1){
00703 for(int kk = 0; kk < Tibnumrechitsas; kk++)
00704 {
00705 meNstpSasTIB[ilay]->Fill(clusizsas[kk]);
00706 meAdcSasTIB[ilay]->Fill(cluchgsas[kk]);
00707 mePosxSasTIB[ilay]->Fill(rechitsasx[kk]);
00708 meErrxSasTIB[ilay]->Fill(sqrt(rechitsaserrx[kk]));
00709 meResSasTIB[ilay]->Fill(rechitsasres[kk]);
00710 mePullLFSasTIB[ilay]->Fill(rechitsasres[kk]/sqrt(rechitsaserrx[kk]));
00711 mePullMFSasTIB[ilay]->Fill(rechitsaspullMF[kk]);
00712 }
00713 for(int l = 0; l<Tibnumrechitsas2; l++){
00714 meChi2SasTIB[ilay]->Fill(chi2sas[l]);
00715 }
00716
00717 }
00718 if(Tibnumrechitmatched>0){
00719 for(int kkk = 0; kkk<Tibnumrechitmatched; kkk++)
00720 {
00721 mePosxMatchedTIB[ilay]->Fill(rechitmatchedx[kkk]);
00722 mePosyMatchedTIB[ilay]->Fill(rechitmatchedy[kkk]);
00723 meErrxMatchedTIB[ilay]->Fill(sqrt(rechitmatchederrxx[kkk]));
00724 meErryMatchedTIB[ilay]->Fill(sqrt(rechitmatchederryy[kkk]));
00725 meResxMatchedTIB[ilay]->Fill(rechitmatchedresx[kkk]);
00726 meResyMatchedTIB[ilay]->Fill(rechitmatchedresy[kkk]);
00727 }
00728 for(int l = 0; l<Tibnumrechitmatched2; l++){
00729 meChi2MatchedTIB[ilay]->Fill(chi2matched[l]);
00730 }
00731
00732 }
00733 }
00734
00735
00736 if (detid.subdetId() == int(StripSubdetector::TOB)){
00737 TOBDetId tobid(myid);
00738 int Tobnumrechitrphi = numrechitrphi;
00739 int Tobnumrechitsas = numrechitsas;
00740 int Tobnumrechitmatched = numrechitmatched;
00741 totTobnumrechitrphi +=numrechitrphi;
00742 totTobnumrechitsas +=numrechitsas;
00743 totTobnumrechitmatched +=numrechitmatched;
00744
00745 int Tobnumrechitrphi2 = i2;
00746 int Tobnumrechitsas2 = j2;
00747 int Tobnumrechitmatched2 = k2;
00748
00749 int ilay = tobid.layer() - 1;
00750
00751 if(tobid.stereo()==0){
00752 for(int k = 0; k<Tobnumrechitrphi; k++){
00753 meNstpRphiTOB[ilay]->Fill(clusizrphi[k]);
00754 meAdcRphiTOB[ilay]->Fill(cluchgrphi[k]);
00755 mePosxRphiTOB[ilay]->Fill(rechitrphix[k]);
00756 meErrxRphiTOB[ilay]->Fill(sqrt(rechitrphierrx[k]));
00757 meResRphiTOB[ilay]->Fill(rechitrphires[k]);
00758 mePullLFRphiTOB[ilay]->Fill(rechitrphires[k]/sqrt(rechitrphierrx[k]));
00759 mePullMFRphiTOB[ilay]->Fill(rechitrphipullMF[k]);
00760 }
00761 for(int l = 0; l<Tobnumrechitrphi2; l++){
00762 meChi2RphiTOB[ilay]->Fill(chi2rphi[l]);
00763 }
00764 } else if(tobid.stereo()==1){
00765 for(int kk = 0; kk < Tobnumrechitsas; kk++)
00766 {
00767 meNstpSasTOB[ilay]->Fill(clusizsas[kk]);
00768 meAdcSasTOB[ilay]->Fill(cluchgsas[kk]);
00769 mePosxSasTOB[ilay]->Fill(rechitsasx[kk]);
00770 meErrxSasTOB[ilay]->Fill(sqrt(rechitsaserrx[kk]));
00771 meResSasTOB[ilay]->Fill(rechitsasres[kk]);
00772 mePullLFSasTOB[ilay]->Fill(rechitsasres[kk]/sqrt(rechitsaserrx[kk]));
00773 mePullMFSasTOB[ilay]->Fill(rechitsaspullMF[kk]);
00774 }
00775 for(int l = 0; l<Tobnumrechitsas2; l++){
00776 meChi2SasTOB[ilay]->Fill(chi2sas[l]);
00777 }
00778
00779 }
00780 if(Tobnumrechitmatched>0){
00781 for(int kkk = 0; kkk<Tobnumrechitmatched; kkk++)
00782 {
00783 mePosxMatchedTOB[ilay]->Fill(rechitmatchedx[kkk]);
00784 mePosyMatchedTOB[ilay]->Fill(rechitmatchedy[kkk]);
00785 meErrxMatchedTOB[ilay]->Fill(sqrt(rechitmatchederrxx[kkk]));
00786 meErryMatchedTOB[ilay]->Fill(sqrt(rechitmatchederryy[kkk]));
00787 meResxMatchedTOB[ilay]->Fill(rechitmatchedresx[kkk]);
00788 meResyMatchedTOB[ilay]->Fill(rechitmatchedresy[kkk]);
00789 }
00790 for(int l = 0; l<Tobnumrechitmatched2; l++){
00791 meChi2MatchedTOB[ilay]->Fill(chi2matched[l]);
00792 }
00793
00794 }
00795 }
00796 if (detid.subdetId() == int(StripSubdetector::TID)){
00797 TIDDetId tidid(myid);
00798 int Tidnumrechitrphi = numrechitrphi;
00799 int Tidnumrechitsas = numrechitsas;
00800 int Tidnumrechitmatched = numrechitmatched;
00801 totTidnumrechitrphi +=numrechitrphi;
00802 totTidnumrechitsas +=numrechitsas;
00803 totTidnumrechitmatched +=numrechitmatched;
00804
00805 int Tidnumrechitrphi2 = i2;
00806 int Tidnumrechitsas2 = j2;
00807 int Tidnumrechitmatched2 = k2;
00808
00809 int ilay = tidid.ring() - 1;
00810
00811 if(tidid.stereo()==0){
00812 for(int k = 0; k<Tidnumrechitrphi; k++){
00813 meNstpRphiTID[ilay]->Fill(clusizrphi[k]);
00814 meAdcRphiTID[ilay]->Fill(cluchgrphi[k]);
00815 mePosxRphiTID[ilay]->Fill(rechitrphix[k]);
00816 meErrxRphiTID[ilay]->Fill(sqrt(rechitrphierrx[k]));
00817 meResRphiTID[ilay]->Fill(rechitrphires[k]);
00818 mePullLFRphiTID[ilay]->Fill(rechitrphires[k]/sqrt(rechitrphierrx[k]));
00819 mePullMFRphiTID[ilay]->Fill(rechitrphipullMF[k]);
00820 }
00821 for(int l = 0; l<Tidnumrechitrphi2; l++){
00822 meChi2RphiTID[ilay]->Fill(chi2rphi[l]);
00823 }
00824
00825 } else if(tidid.stereo()==1){
00826 for(int kk = 0; kk < Tidnumrechitsas; kk++)
00827 {
00828 meNstpSasTID[ilay]->Fill(clusizsas[kk]);
00829 meAdcSasTID[ilay]->Fill(cluchgsas[kk]);
00830 mePosxSasTID[ilay]->Fill(rechitsasx[kk]);
00831 meErrxSasTID[ilay]->Fill(sqrt(rechitsaserrx[kk]));
00832 meResSasTID[ilay]->Fill(rechitsasres[kk]);
00833 mePullLFSasTID[ilay]->Fill(rechitsasres[kk]/sqrt(rechitsaserrx[kk]));
00834 mePullMFSasTID[ilay]->Fill(rechitsaspullMF[kk]);
00835 }
00836 for(int l = 0; l<Tidnumrechitsas2; l++){
00837 meChi2SasTID[ilay]->Fill(chi2sas[l]);
00838 }
00839
00840 }
00841 if(Tidnumrechitmatched>0){
00842 for(int kkk = 0; kkk<Tidnumrechitmatched; kkk++)
00843 {
00844 mePosxMatchedTID[ilay]->Fill(rechitmatchedx[kkk]);
00845 mePosyMatchedTID[ilay]->Fill(rechitmatchedy[kkk]);
00846 meErrxMatchedTID[ilay]->Fill(sqrt(rechitmatchederrxx[kkk]));
00847 meErryMatchedTID[ilay]->Fill(sqrt(rechitmatchederryy[kkk]));
00848 meResxMatchedTID[ilay]->Fill(rechitmatchedresx[kkk]);
00849 meResyMatchedTID[ilay]->Fill(rechitmatchedresy[kkk]);
00850 }
00851 for(int l = 0; l<Tidnumrechitmatched2; l++){
00852 meChi2MatchedTID[ilay]->Fill(chi2matched[l]);
00853 }
00854
00855 }
00856 }
00857 if (detid.subdetId() == int(StripSubdetector::TEC)){
00858 TECDetId tecid(myid);
00859 int Tecnumrechitrphi = numrechitrphi;
00860 int Tecnumrechitsas = numrechitsas;
00861 int Tecnumrechitmatched = numrechitmatched;
00862 totTecnumrechitrphi +=numrechitrphi;
00863 totTecnumrechitsas +=numrechitsas;
00864 totTecnumrechitmatched +=numrechitmatched;
00865
00866 int Tecnumrechitrphi2 = i2;
00867 int Tecnumrechitsas2 = j2;
00868 int Tecnumrechitmatched2 = k2;
00869
00870 int ilay = tecid.ring() - 1;
00871
00872 if(tecid.stereo()==0){
00873 for(int k = 0; k<Tecnumrechitrphi; k++){
00874 meNstpRphiTEC[ilay]->Fill(clusizrphi[k]);
00875 meAdcRphiTEC[ilay]->Fill(cluchgrphi[k]);
00876 mePosxRphiTEC[ilay]->Fill(rechitrphix[k]);
00877 meErrxRphiTEC[ilay]->Fill(sqrt(rechitrphierrx[k]));
00878 meResRphiTEC[ilay]->Fill(rechitrphires[k]);
00879 mePullLFRphiTEC[ilay]->Fill(rechitrphires[k]/sqrt(rechitrphierrx[k]));
00880 mePullMFRphiTEC[ilay]->Fill(rechitrphipullMF[k]);
00881 }
00882 for(int l = 0; l<Tecnumrechitrphi2; l++){
00883 meChi2RphiTEC[ilay]->Fill(chi2rphi[l]);
00884 }
00885
00886 } else if(tecid.stereo()==1){
00887 for(int kk = 0; kk < Tecnumrechitsas; kk++)
00888 {
00889 meNstpSasTEC[ilay]->Fill(clusizsas[kk]);
00890 meAdcSasTEC[ilay]->Fill(cluchgsas[kk]);
00891 mePosxSasTEC[ilay]->Fill(rechitsasx[kk]);
00892 meErrxSasTEC[ilay]->Fill(sqrt(rechitsaserrx[kk]));
00893 meResSasTEC[ilay]->Fill(rechitsasres[kk]);
00894 mePullLFSasTEC[ilay]->Fill(rechitsasres[kk]/sqrt(rechitsaserrx[kk]));
00895 mePullMFSasTEC[ilay]->Fill(rechitsaspullMF[kk]);
00896 }
00897 for(int l = 0; l<Tecnumrechitsas2; l++){
00898 meChi2SasTEC[ilay]->Fill(chi2sas[l]);
00899 }
00900
00901 }
00902 if(Tecnumrechitmatched>0){
00903 for(int kkk = 0; kkk<Tecnumrechitmatched; kkk++)
00904 {
00905 mePosxMatchedTEC[ilay]->Fill(rechitmatchedx[kkk]);
00906 mePosyMatchedTEC[ilay]->Fill(rechitmatchedy[kkk]);
00907 meErrxMatchedTEC[ilay]->Fill(sqrt(rechitmatchederrxx[kkk]));
00908 meErryMatchedTEC[ilay]->Fill(sqrt(rechitmatchederryy[kkk]));
00909 meResxMatchedTEC[ilay]->Fill(rechitmatchedresx[kkk]);
00910 meResyMatchedTEC[ilay]->Fill(rechitmatchedresy[kkk]);
00911 }
00912 for(int l = 0; l<Tecnumrechitmatched2; l++){
00913 meChi2MatchedTEC[ilay]->Fill(chi2matched[l]);
00914 }
00915
00916 }
00917 }
00918
00919 }
00920 }
00921
00922
00923
00924 meNumRphiTIB->Fill(totTibnumrechitrphi);
00925 meNumSasTIB->Fill(totTibnumrechitsas);
00926 meNumMatchedTIB->Fill(totTibnumrechitmatched);
00927 meNumRphiTOB->Fill(totTobnumrechitrphi);
00928 meNumSasTOB->Fill(totTobnumrechitsas);
00929 meNumMatchedTOB->Fill(totTobnumrechitmatched);
00930 meNumRphiTID->Fill(totTidnumrechitrphi);
00931 meNumSasTID->Fill(totTidnumrechitsas);
00932 meNumMatchedTID->Fill(totTidnumrechitmatched);
00933 meNumRphiTEC->Fill(totTecnumrechitrphi);
00934 meNumSasTEC->Fill(totTecnumrechitsas);
00935 meNumMatchedTEC->Fill(totTecnumrechitmatched);
00936
00937 meNumTotRphi->Fill(totrechitrphi);
00938 meNumTotSas->Fill(totrechitsas);
00939 meNumTotMatched->Fill(totrechitmatched);
00940
00941
00942 }
00943
00944
00945
00946 std::pair<LocalPoint,LocalVector> SiStripRecHitsValid::projectHit( const PSimHit& hit, const StripGeomDetUnit* stripDet,
00947 const BoundPlane& plane)
00948 {
00949
00950
00951
00952 const StripTopology& topol = stripDet->specificTopology();
00953 GlobalPoint globalpos= stripDet->surface().toGlobal(hit.localPosition());
00954 LocalPoint localHit = plane.toLocal(globalpos);
00955
00956 LocalVector locdir=hit.localDirection();
00957
00958
00959 GlobalVector globaldir= stripDet->surface().toGlobal(locdir);
00960 LocalVector dir=plane.toLocal(globaldir);
00961 float scale = -localHit.z() / dir.z();
00962
00963 LocalPoint projectedPos = localHit + scale*dir;
00964
00965
00966
00967 float selfAngle = topol.stripAngle( topol.strip( hit.localPosition()));
00968
00969 LocalVector stripDir( sin(selfAngle), cos(selfAngle), 0);
00970
00971 LocalVector localStripDir( plane.toLocal(stripDet->surface().toGlobal( stripDir)));
00972
00973 return std::pair<LocalPoint,LocalVector>( projectedPos, localStripDir);
00974 }
00975
00976