00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include <memory>
00022 #include <map>
00023 #include <sstream>
00024 #include <math.h>
00025 #include <utility>
00026 #include <vector>
00027 #include <iostream>
00028
00029
00030 #include "TH1.h"
00031 #include "TH2.h"
00032 #include "TProfile.h"
00033 #include "TFile.h"
00034 #include "TTree.h"
00035 #include "TF1.h"
00036 #include "TMath.h"
00037
00038
00039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00040 #include "FWCore/Framework/interface/Frameworkfwd.h"
00041 #include "FWCore/Framework/interface/EDAnalyzer.h"
00042 #include "FWCore/Framework/interface/EventSetup.h"
00043 #include "FWCore/Framework/interface/Event.h"
00044 #include "FWCore/Framework/interface/MakerMacros.h"
00045 #include "FWCore/Framework/interface/ESHandle.h"
00046 #include "FWCore/ServiceRegistry/interface/Service.h"
00047
00048 #include "DataFormats/DetId/interface/DetId.h"
00049 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00050 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00051 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00052 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00053 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00054 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00055 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00056 #include "DataFormats/Math/interface/deltaPhi.h"
00057
00058 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00059 #include "CommonTools/Utils/interface/TFileDirectory.h"
00060
00061 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00062 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00063 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00064
00065 #include "DQMServices/Core/interface/DQMStore.h"
00066 #include "DQMServices/Core/interface/MonitorElement.h"
00067
00068 #include "Alignment/CommonAlignment/interface/Alignable.h"
00069 #include "Alignment/CommonAlignment/interface/Utilities.h"
00070 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"
00071 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00072 #include "Alignment/TrackerAlignment/interface/TrackerAlignableId.h"
00073
00074 #include "Alignment/OfflineValidation/interface/TrackerValidationVariables.h"
00075 #include "Alignment/OfflineValidation/interface/TkOffTreeVariables.h"
00076
00077
00078
00079
00080 class TrackerOfflineValidation : public edm::EDAnalyzer {
00081 public:
00082 explicit TrackerOfflineValidation(const edm::ParameterSet&);
00083 ~TrackerOfflineValidation();
00084
00085 enum HistogrammType { XResidual, NormXResidual,
00086 YResidual,
00087 XprimeResidual, NormXprimeResidual,
00088 YprimeResidual, NormYprimeResidual,
00089 XResidualProfile, YResidualProfile };
00090
00091 private:
00092
00093 struct ModuleHistos{
00094 ModuleHistos() : ResHisto(), NormResHisto(), ResYHisto(),
00095 ResXprimeHisto(), NormResXprimeHisto(),
00096 ResYprimeHisto(), NormResYprimeHisto(),
00097 ResXvsXProfile(), ResXvsYProfile(),
00098 ResYvsXProfile(), ResYvsYProfile(),
00099 LocalX(), LocalY() {}
00100 TH1* ResHisto;
00101 TH1* NormResHisto;
00102 TH1* ResYHisto;
00103
00104 TH1* ResXprimeHisto;
00105 TH1* NormResXprimeHisto;
00106 TH1* ResYprimeHisto;
00107 TH1* NormResYprimeHisto;
00108
00109 TProfile* ResXvsXProfile;
00110 TProfile* ResXvsYProfile;
00111 TProfile* ResYvsXProfile;
00112 TProfile* ResYvsYProfile;
00113
00114 TH1* LocalX;
00115 TH1* LocalY;
00116 };
00117
00118
00119 struct SummaryContainer{
00120 SummaryContainer() : sumXResiduals_(), summaryXResiduals_(),
00121 sumNormXResiduals_(), summaryNormXResiduals_(),
00122 sumYResiduals_(), summaryYResiduals_() ,
00123 sumNormYResiduals_(), summaryNormYResiduals_() {}
00124
00125 TH1* sumXResiduals_;
00126 TH1* summaryXResiduals_;
00127 TH1* sumNormXResiduals_;
00128 TH1* summaryNormXResiduals_;
00129 TH1* sumYResiduals_;
00130 TH1* summaryYResiduals_;
00131 TH1* sumNormYResiduals_;
00132 TH1* summaryNormYResiduals_;
00133 };
00134
00135
00136 struct DirectoryWrapper{
00137 DirectoryWrapper(const DirectoryWrapper& upDir,const std::string& newDir,
00138 const std::string& basedir,bool useDqmMode)
00139 : tfd(0),
00140 dqmMode(useDqmMode),
00141 theDbe(0) {
00142 if (newDir.length()!=0){
00143 if(upDir.directoryString.length()!=0)directoryString=upDir.directoryString+"/"+newDir;
00144 else directoryString = newDir;
00145 }
00146 else
00147 directoryString=upDir.directoryString;
00148
00149 if (!dqmMode){
00150 if (newDir.length()==0) tfd.reset(&(*upDir.tfd));
00151 else
00152 tfd.reset(new TFileDirectory(upDir.tfd->mkdir(newDir)));
00153 }
00154 else {
00155 theDbe=edm::Service<DQMStore>().operator->();
00156 }
00157 }
00158
00159 DirectoryWrapper(const std::string& newDir,const std::string& basedir,bool useDqmMode)
00160 : tfd(0),
00161 dqmMode(useDqmMode),
00162 theDbe(0) {
00163 if (!dqmMode){
00164 edm::Service<TFileService> fs;
00165 if (newDir.length()==0){
00166 tfd.reset(new TFileDirectory(static_cast<TFileDirectory&>(*fs)));
00167 }
00168 else {
00169 tfd.reset(new TFileDirectory(fs->mkdir(newDir)));
00170 directoryString=newDir;
00171 }
00172 }
00173 else {
00174 if (newDir.length()!=0){
00175 if(basedir.length()!=0)directoryString=basedir+"/"+newDir;
00176 else directoryString = newDir;
00177 }
00178 else directoryString=basedir;
00179 theDbe=edm::Service<DQMStore>().operator->();
00180 }
00181 }
00182
00183 template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX);
00184 template <typename T> TH1* make(const char* name,const char* title,int nBinX,double *xBins);
00185 template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,int nBinY,double minBinY,double maxBinY);
00186 template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,double minBinY,double maxBinY);
00187
00188 std::auto_ptr<TFileDirectory> tfd;
00189 std::string directoryString;
00190 const bool dqmMode;
00191 DQMStore* theDbe;
00192 };
00193
00194
00195
00196
00197
00198 virtual void analyze(const edm::Event&, const edm::EventSetup&);
00199 virtual void endJob();
00200
00201 virtual void checkBookHists(const edm::EventSetup& setup);
00202
00203 void bookGlobalHists(DirectoryWrapper& tfd);
00204 void bookDirHists(DirectoryWrapper& tfd, const Alignable& ali);
00205 void bookHists(DirectoryWrapper& tfd, const Alignable& ali, align::StructureType type, int i);
00206
00207 void collateSummaryHists( DirectoryWrapper& tfd, const Alignable& ali, int i,
00208 std::vector<TrackerOfflineValidation::SummaryContainer>& vLevelProfiles);
00209
00210 void fillTree(TTree& tree, const std::map<int, TrackerOfflineValidation::ModuleHistos>& moduleHist_,
00211 TkOffTreeVariables& treeMem, const TrackerGeometry& tkgeom);
00212
00213 TrackerOfflineValidation::SummaryContainer bookSummaryHists(DirectoryWrapper& tfd,
00214 const Alignable& ali,
00215 align::StructureType type, int i);
00216
00217 ModuleHistos& getHistStructFromMap(const DetId& detid);
00218
00219 bool isBarrel(uint32_t subDetId);
00220 bool isEndCap(uint32_t subDetId);
00221 bool isPixel(uint32_t subDetId);
00222 bool isDetOrDetUnit(align::StructureType type);
00223
00224 TH1* bookTH1F(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00225 int nBinsX, double lowX, double highX);
00226
00227 TProfile* bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00228 int nBinsX, double lowX, double highX);
00229
00230 TProfile* bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00231 int nBinsX, double lowX, double highX, double lowY, double highY);
00232
00233 void getBinning(uint32_t subDetId, TrackerOfflineValidation::HistogrammType residualtype,
00234 int& nBinsX, double& lowerBoundX, double& upperBoundX);
00235
00236 void summarizeBinInContainer(int bin, SummaryContainer& targetContainer,
00237 SummaryContainer& sourceContainer);
00238
00239 void summarizeBinInContainer(int bin, uint32_t subDetId, SummaryContainer& targetContainer,
00240 ModuleHistos& sourceContainer);
00241
00242 void setSummaryBin(int bin, TH1* targetHist, TH1* sourceHist);
00243
00244 float Fwhm(const TH1* hist) const;
00245 std::pair<float,float> fitResiduals(TH1* hist) const;
00246 float getMedian( const TH1* hist) const;
00247
00248
00249 template <class OBJECT_TYPE> int GetIndex(const std::vector<OBJECT_TYPE*>& vec, const TString& name);
00250
00251
00252
00253
00254 const edm::ParameterSet parSet_;
00255 edm::ESHandle<TrackerGeometry> tkGeom_;
00256 const TrackerGeometry *bareTkGeomPtr_;
00257
00258
00259 const bool lCoorHistOn_;
00260 const bool moduleLevelHistsTransient_;
00261 const bool moduleLevelProfiles_;
00262 const bool stripYResiduals_;
00263 const bool useFwhm_;
00264 const bool useFit_;
00265 const bool useOverflowForRMS_;
00266 const bool dqmMode_;
00267 const std::string moduleDirectory_;
00268
00269
00270 std::vector<TH1*> vDeleteObjects_;
00271
00272 std::vector<TH1*> vTrackHistos_;
00273 std::vector<TH1*> vTrackProfiles_;
00274 std::vector<TH1*> vTrack2DHistos_;
00275
00276 std::map<int,TrackerOfflineValidation::ModuleHistos> mPxbResiduals_;
00277 std::map<int,TrackerOfflineValidation::ModuleHistos> mPxeResiduals_;
00278 std::map<int,TrackerOfflineValidation::ModuleHistos> mTibResiduals_;
00279 std::map<int,TrackerOfflineValidation::ModuleHistos> mTidResiduals_;
00280 std::map<int,TrackerOfflineValidation::ModuleHistos> mTobResiduals_;
00281 std::map<int,TrackerOfflineValidation::ModuleHistos> mTecResiduals_;
00282 };
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293 template <class OBJECT_TYPE>
00294 int TrackerOfflineValidation::GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name)
00295 {
00296 int result = 0;
00297 for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
00298 iter != iterEnd; ++iter, ++result) {
00299 if (*iter && (*iter)->GetName() == name) return result;
00300 }
00301 edm::LogError("Alignment") << "@SUB=TrackerOfflineValidation::GetIndex" << " could not find " << name;
00302 return -1;
00303 }
00304
00305
00306 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TH1F>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX){
00307 if(dqmMode){theDbe->setCurrentFolder(directoryString); return theDbe->book1D(name,title,nBinX,minBinX,maxBinX)->getTH1();}
00308 else{return tfd->make<TH1F>(name,title,nBinX,minBinX,maxBinX);}
00309 }
00310
00311 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name,const char* title,int nBinX,double *xBins){
00312 if(dqmMode){
00313 theDbe->setCurrentFolder(directoryString);
00314
00315 TProfile *tmpProfile=new TProfile(name,title,nBinX,xBins);
00316 tmpProfile->SetDirectory(0);
00317 return theDbe->bookProfile(name,tmpProfile)->getTH1();
00318 }
00319 else{return tfd->make<TProfile>(name,title,nBinX,xBins);}
00320 }
00321
00322 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX){
00323 if(dqmMode){
00324 theDbe->setCurrentFolder(directoryString);
00325
00326 TProfile *tmpProfile=new TProfile(name,title,nBinX,minBinX,maxBinX);
00327 tmpProfile->SetDirectory(0);
00328 return theDbe->bookProfile(name,tmpProfile)->getTH1();
00329 }
00330 else{return tfd->make<TProfile>(name,title,nBinX,minBinX,maxBinX);}
00331 }
00332
00333 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name ,const char* title,int nbinX,double minX ,double maxX,double minY,double maxY){
00334 if(dqmMode){
00335 theDbe->setCurrentFolder(directoryString);
00336 int dummy(0);
00337 return (theDbe->bookProfile(name,title,nbinX,minX,maxX,dummy,minY,maxY)->getTH1());
00338 }
00339 else{
00340 return tfd->make<TProfile>(name,title,nbinX,minX,maxX,minY,maxY);
00341 }
00342 }
00343
00344 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TH2F>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,int nBinY,double minBinY,double maxBinY){
00345 if(dqmMode){theDbe->setCurrentFolder(directoryString); return theDbe->book2D(name,title,nBinX,minBinX,maxBinX,nBinY,minBinY,maxBinY)->getTH1();}
00346 else{return tfd->make<TH2F>(name,title,nBinX,minBinX,maxBinX,nBinY,minBinY,maxBinY);}
00347 }
00348
00349
00350
00351
00352
00353 TrackerOfflineValidation::TrackerOfflineValidation(const edm::ParameterSet& iConfig)
00354 : parSet_(iConfig), bareTkGeomPtr_(0), lCoorHistOn_(parSet_.getParameter<bool>("localCoorHistosOn")),
00355 moduleLevelHistsTransient_(parSet_.getParameter<bool>("moduleLevelHistsTransient")),
00356 moduleLevelProfiles_(parSet_.getParameter<bool>("moduleLevelProfiles")),
00357 stripYResiduals_(parSet_.getParameter<bool>("stripYResiduals")),
00358 useFwhm_(parSet_.getParameter<bool>("useFwhm")),
00359 useFit_(parSet_.getParameter<bool>("useFit")),
00360 useOverflowForRMS_(parSet_.getParameter<bool>("useOverflowForRMS")),
00361 dqmMode_(parSet_.getParameter<bool>("useInDqmMode")),
00362 moduleDirectory_(parSet_.getParameter<std::string>("moduleDirectoryInOutput"))
00363 {
00364 }
00365
00366
00367 TrackerOfflineValidation::~TrackerOfflineValidation()
00368 {
00369
00370
00371 for( std::vector<TH1*>::const_iterator it = vDeleteObjects_.begin(), itEnd = vDeleteObjects_.end();
00372 it != itEnd;
00373 ++it) delete *it;
00374 }
00375
00376
00377
00378
00379
00380
00381
00382
00383 void
00384 TrackerOfflineValidation::checkBookHists(const edm::EventSetup& es)
00385 {
00386 es.get<TrackerDigiGeometryRecord>().get( tkGeom_ );
00387 const TrackerGeometry *newBareTkGeomPtr = &(*tkGeom_);
00388 if (newBareTkGeomPtr == bareTkGeomPtr_) return;
00389
00390 if (!bareTkGeomPtr_) {
00391
00392 AlignableTracker aliTracker(&(*tkGeom_));
00393
00394 edm::LogInfo("TrackerOfflineValidation") << "There are " << newBareTkGeomPtr->detIds().size()
00395 << " dets in the Geometry record.\n"
00396 << "Out of these "<<newBareTkGeomPtr->detUnitIds().size()
00397 <<" are detUnits";
00398
00399
00400 std::string globDir("GlobalTrackVariables");
00401 DirectoryWrapper trackglobal(globDir,moduleDirectory_,dqmMode_);
00402 this->bookGlobalHists(trackglobal);
00403
00404
00405 DirectoryWrapper tfdw("",moduleDirectory_,dqmMode_);
00406 this->bookDirHists(tfdw, aliTracker);
00407 }
00408 else {
00409 edm::LogWarning("GeometryChange") << "@SUB=checkBookHists"
00410 << "TrackerGeometry changed, but will not re-book hists!";
00411 }
00412 bareTkGeomPtr_ = newBareTkGeomPtr;
00413 }
00414
00415
00416 void
00417 TrackerOfflineValidation::bookGlobalHists(DirectoryWrapper& tfd )
00418 {
00419
00420 vTrackHistos_.push_back(tfd.make<TH1F>("h_tracketa",
00421 "Track #eta;#eta_{Track};Number of Tracks",
00422 90,-3.,3.));
00423 vTrackHistos_.push_back(tfd.make<TH1F>("h_trackphi",
00424 "Track #phi;#phi_{Track};Number of Tracks",
00425 90,-3.15,3.15));
00426 vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfValidHits",
00427 "Track # of valid hits;# of valid hits _{Track};Number of Tracks",
00428 40,0.,40.));
00429 vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfLostHits",
00430 "Track # of lost hits;# of lost hits _{Track};Number of Tracks",
00431 10,0.,10.));
00432 vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature",
00433 "Curvature #kappa;#kappa_{Track};Number of Tracks",
00434 100,-.05,.05));
00435 vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_pos",
00436 "Curvature |#kappa| Positive Tracks;|#kappa_{pos Track}|;Number of Tracks",
00437 100,.0,.05));
00438 vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_neg",
00439 "Curvature |#kappa| Negative Tracks;|#kappa_{neg Track}|;Number of Tracks",
00440 100,.0,.05));
00441 vTrackHistos_.push_back(tfd.make<TH1F>("h_diff_curvature",
00442 "Curvature |#kappa| Tracks Difference;|#kappa_{Track}|;# Pos Tracks - # Neg Tracks",
00443 100,.0,.05));
00444 vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2",
00445 "#chi^{2};#chi^{2}_{Track};Number of Tracks",
00446 500,-0.01,500.));
00447 vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2Prob",
00448 "#chi^{2} probability;#chi^{2}prob_{Track};Number of Tracks",
00449 100,0.0,1.));
00450 vTrackHistos_.push_back(tfd.make<TH1F>("h_normchi2",
00451 "#chi^{2}/ndof;#chi^{2}/ndof;Number of Tracks",
00452 100,-0.01,10.));
00453 vTrackHistos_.push_back(tfd.make<TH1F>("h_pt",
00454 "p_{T}^{track};p_{T}^{track} [GeV];Number of Tracks",
00455 250,0.,250));
00456 vTrackHistos_.push_back(tfd.make<TH1F>("h_ptResolution",
00457 "#delta{p_{T}/p_{T}^{track}};#delta_{p_{T}/p_{T}^{track}};Number of Tracks",
00458 100,0.,0.5));
00459
00460 vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_phi",
00461 "Transverse Impact Parameter vs. #phi;#phi_{Track};#LT d_{0} #GT [cm]",
00462 100,-3.15,3.15));
00463 vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_phi",
00464 "Longitudinal Impact Parameter vs. #phi;#phi_{Track};#LT d_{z} #GT [cm]",
00465 100,-3.15,3.15));
00466 vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_eta",
00467 "Transverse Impact Parameter vs. #eta;#eta_{Track};#LT d_{0} #GT [cm]",
00468 100,-3.15,3.15));
00469 vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_eta",
00470 "Longitudinal Impact Parameter vs. #eta;#eta_{Track};#LT d_{z} #GT [cm]",
00471 100,-3.15,3.15));
00472 vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_phi",
00473 "#chi^{2} vs. #phi;#phi_{Track};#LT #chi^{2} #GT",
00474 100,-3.15,3.15));
00475 vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_phi",
00476 "#chi^{2} probablility vs. #phi;#phi_{Track};#LT #chi^{2} probability#GT",
00477 100,-3.15,3.15));
00478 vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_d0",
00479 "#chi^{2} probablility vs. |d_{0}|;|d_{0}|[cm];#LT #chi^{2} probability#GT",
00480 100,0,80));
00481 vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_phi",
00482 "#chi^{2}/ndof vs. #phi;#phi_{Track};#LT #chi^{2}/ndof #GT",
00483 100,-3.15,3.15));
00484 vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_eta",
00485 "#chi^{2} vs. #eta;#eta_{Track};#LT #chi^{2} #GT",
00486 100,-3.15,3.15));
00487
00488 double xBins[19]={0.,0.15,0.5,1.,1.5,2.,2.5,3.,3.5,4.,4.5,5.,7.,10.,15.,25.,40.,100.,200.};
00489 vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_pt",
00490 "norm #chi^{2} vs. p_{T}_{Track}; p_{T}_{Track};#LT #chi^{2}/ndof #GT",
00491 18,xBins));
00492
00493 vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_p",
00494 "#chi^{2}/ndof vs. p_{Track};p_{Track};#LT #chi^{2}/ndof #GT",
00495 18,xBins));
00496 vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_eta",
00497 "#chi^{2} probability vs. #eta;#eta_{Track};#LT #chi^{2} probability #GT",
00498 100,-3.15,3.15));
00499 vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_eta",
00500 "#chi^{2}/ndof vs. #eta;#eta_{Track};#LT #chi^{2}/ndof #GT",
00501 100,-3.15,3.15));
00502 vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_phi",
00503 "#kappa vs. #phi;#phi_{Track};#kappa",
00504 100,-3.15,3.15));
00505 vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_eta",
00506 "#kappa vs. #eta;#eta_{Track};#kappa",
00507 100,-3.15,3.15));
00508 vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_phi",
00509 "#delta_{p_{T}}/p_{T}^{track};#phi^{track};#delta_{p_{T}}/p_{T}^{track}",
00510 100, -3.15,3.15));
00511 vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_eta",
00512 "#delta_{p_{T}}/p_{T}^{track};#eta^{track};#delta_{p_{T}}/p_{T}^{track}",
00513 100, -3.15,3.15));
00514
00515 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_phi",
00516 "Transverse Impact Parameter vs. #phi;#phi_{Track};d_{0} [cm]",
00517 100, -3.15, 3.15, 100,-1.,1.) );
00518 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_phi",
00519 "Longitudinal Impact Parameter vs. #phi;#phi_{Track};d_{z} [cm]",
00520 100, -3.15, 3.15, 100,-100.,100.));
00521 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_eta",
00522 "Transverse Impact Parameter vs. #eta;#eta_{Track};d_{0} [cm]",
00523 100, -3.15, 3.15, 100,-1.,1.));
00524 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_eta",
00525 "Longitudinal Impact Parameter vs. #eta;#eta_{Track};d_{z} [cm]",
00526 100, -3.15, 3.15, 100,-100.,100.));
00527 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_phi",
00528 "#chi^{2} vs. #phi;#phi_{Track};#chi^{2}",
00529 100, -3.15, 3.15, 500, 0., 500.));
00530 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_phi",
00531 "#chi^{2} probability vs. #phi;#phi_{Track};#chi^{2} probability",
00532 100, -3.15, 3.15, 100, 0., 1.));
00533 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_d0",
00534 "#chi^{2} probability vs. |d_{0}|;|d_{0}| [cm];#chi^{2} probability",
00535 100, 0, 80, 100, 0., 1.));
00536 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_phi",
00537 "#chi^{2}/ndof vs. #phi;#phi_{Track};#chi^{2}/ndof",
00538 100, -3.15, 3.15, 100, 0., 10.));
00539 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_eta",
00540 "#chi^{2} vs. #eta;#eta_{Track};#chi^{2}",
00541 100, -3.15, 3.15, 500, 0., 500.));
00542 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_eta",
00543 "#chi^{2} probaility vs. #eta;#eta_{Track};#chi^{2} probability",
00544 100, -3.15, 3.15, 100, 0., 1.));
00545 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_eta",
00546 "#chi^{2}/ndof vs. #eta;#eta_{Track};#chi^{2}/ndof",
00547 100,-3.15,3.15, 100, 0., 10.));
00548 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_phi",
00549 "#kappa vs. #phi;#phi_{Track};#kappa",
00550 100,-3.15,3.15, 100, .0,.05));
00551 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_eta",
00552 "#kappa vs. #eta;#eta_{Track};#kappa",
00553 100,-3.15,3.15, 100, .0,.05));
00554 vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_kappa",
00555 "#kappa vs. #chi^{2}/ndof;#chi^{2}/ndof;#kappa",
00556 100,0.,10, 100,-.03,.03));
00557
00558
00559 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixB",
00560 "#momentum vs. #resX in pixB;#momentum;#resX",
00561 15,0.,15., 200, -0.1,0.1));
00562 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixE",
00563 "#momentum vs. #resX in pixE;#momentum;#resX",
00564 15,0.,15., 200, -0.1,0.1));
00565 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TIB",
00566 "#momentum vs. #resX in TIB;#momentum;#resX",
00567 15,0.,15., 200, -0.1,0.1));
00568 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TID",
00569 "#momentum vs. #resX in TID;#momentum;#resX",
00570 15,0.,15., 200, -0.1,0.1));
00571 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TOB",
00572 "#momentum vs. #resX in TOB;#momentum;#resX",
00573 15,0.,15., 200, -0.1,0.1));
00574 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TEC",
00575 "#momentum vs. #resX in TEC;#momentum;#resX",
00576 15,0.,15., 200, -0.1,0.1));
00577
00578
00579 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixB",
00580 "#momentum vs. #resY in pixB;#momentum;#resY",
00581 15,0.,15., 200, -0.1,0.1));
00582 vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixE",
00583 "#momentum vs. #resY in pixE;#momentum;#resY",
00584 15,0.,15., 200, -0.1,0.1));
00585
00586 }
00587
00588
00589 void
00590 TrackerOfflineValidation::bookDirHists(DirectoryWrapper& tfd, const Alignable& ali)
00591 {
00592 std::vector<Alignable*> alivec(ali.components());
00593 for(int i=0, iEnd = ali.components().size();i < iEnd; ++i) {
00594 std::string structurename = AlignableObjectId::idToString((alivec)[i]->alignableObjectId());
00595 LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
00596 std::stringstream dirname;
00597 dirname << structurename;
00598
00599 if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << i+1;
00600
00601 if (structurename.find("Endcap",0) != std::string::npos ) {
00602 DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
00603 bookHists(f, *(alivec)[i], ali.alignableObjectId() , i);
00604 bookDirHists( f, *(alivec)[i]);
00605 } else if( !(this->isDetOrDetUnit( (alivec)[i]->alignableObjectId()) )
00606 || alivec[i]->components().size() > 1) {
00607 DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
00608 bookHists(tfd, *(alivec)[i], ali.alignableObjectId() , i);
00609 bookDirHists( f, *(alivec)[i]);
00610 } else {
00611 bookHists(tfd, *(alivec)[i], ali.alignableObjectId() , i);
00612 }
00613 }
00614 }
00615
00616
00617 void
00618 TrackerOfflineValidation::bookHists(DirectoryWrapper& tfd, const Alignable& ali, align::StructureType type, int i)
00619 {
00620 TrackerAlignableId aliid;
00621 const DetId id = ali.id();
00622
00623
00624
00625
00626 std::pair<int,int> subdetandlayer = aliid.typeAndLayerFromDetId(id);
00627
00628 align::StructureType subtype = align::invalid;
00629
00630
00631 if (type == align::AlignableDetUnit )subtype = type;
00632 else if( this->isDetOrDetUnit(ali.alignableObjectId()) ) subtype = ali.alignableObjectId();
00633
00634
00635 std::stringstream histoname, histotitle, normhistoname, normhistotitle,
00636 yhistoname, yhistotitle,
00637 xprimehistoname, xprimehistotitle, normxprimehistoname, normxprimehistotitle,
00638 yprimehistoname, yprimehistotitle, normyprimehistoname, normyprimehistotitle,
00639 localxname, localxtitle, localyname, localytitle,
00640 resxvsxprofilename, resxvsxprofiletitle, resyvsxprofilename, resyvsxprofiletitle,
00641 resxvsyprofilename, resxvsyprofiletitle, resyvsyprofilename, resyvsyprofiletitle;
00642
00643 std::string wheel_or_layer;
00644
00645 if( this->isEndCap(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_wheel_";
00646 else if ( this->isBarrel(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_layer_";
00647 else edm::LogWarning("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookHists"
00648 << "Unknown subdetid: " << subdetandlayer.first;
00649
00650 histoname << "h_residuals_subdet_" << subdetandlayer.first
00651 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00652 yhistoname << "h_y_residuals_subdet_" << subdetandlayer.first
00653 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00654 xprimehistoname << "h_xprime_residuals_subdet_" << subdetandlayer.first
00655 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00656 yprimehistoname << "h_yprime_residuals_subdet_" << subdetandlayer.first
00657 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00658 normhistoname << "h_normresiduals_subdet_" << subdetandlayer.first
00659 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00660 normxprimehistoname << "h_normxprimeresiduals_subdet_" << subdetandlayer.first
00661 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00662 normyprimehistoname << "h_normyprimeresiduals_subdet_" << subdetandlayer.first
00663 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00664 histotitle << "X Residual for module " << id.rawId() << ";x_{tr} - x_{hit} [cm]";
00665 yhistotitle << "Y Residual for module " << id.rawId() << ";y_{tr} - y_{hit} [cm]";
00666 normhistotitle << "Normalized Residual for module " << id.rawId() << ";x_{tr} - x_{hit}/#sigma";
00667 xprimehistotitle << "X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})' [cm]";
00668 normxprimehistotitle << "Normalized X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})'/#sigma";
00669 yprimehistotitle << "Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})' [cm]";
00670 normyprimehistotitle << "Normalized Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})'/#sigma";
00671
00672 if ( moduleLevelProfiles_ ) {
00673 localxname << "h_localx_subdet_" << subdetandlayer.first
00674 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00675 localyname << "h_localy_subdet_" << subdetandlayer.first
00676 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00677 localxtitle << "u local for module " << id.rawId() << "; u_{tr,r}";
00678 localytitle << "v local for module " << id.rawId() << "; v_{tr,r}";
00679
00680 resxvsxprofilename << "p_residuals_x_vs_x_subdet_" << subdetandlayer.first
00681 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00682 resyvsxprofilename << "p_residuals_y_vs_x_subdet_" << subdetandlayer.first
00683 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00684 resxvsyprofilename << "p_residuals_x_vs_y_subdet_" << subdetandlayer.first
00685 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00686 resyvsyprofilename << "p_residuals_y_vs_y_subdet_" << subdetandlayer.first
00687 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00688 resxvsxprofiletitle << "U Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
00689 resyvsxprofiletitle << "V Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(v_{tr} - v_{hit})/tan#beta [cm]";
00690 resxvsyprofiletitle << "U Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
00691 resyvsyprofiletitle << "V Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(v_{tr} - v_{hit})/tan#beta [cm]";
00692 }
00693
00694 if( this->isDetOrDetUnit( subtype ) ) {
00695 ModuleHistos &histStruct = this->getHistStructFromMap(id);
00696 int nbins = 0;
00697 double xmin = 0., xmax = 0.;
00698 double ymin = -0.1, ymax = 0.1;
00699
00700
00701 bool moduleLevelHistsTransient(moduleLevelHistsTransient_);
00702 if (dqmMode_) moduleLevelHistsTransient = false;
00703
00704
00705 if(lCoorHistOn_) {
00706 this->getBinning(id.subdetId(), XResidual, nbins, xmin, xmax);
00707 histStruct.ResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00708 histoname.str().c_str(),histotitle.str().c_str(),
00709 nbins, xmin, xmax);
00710 this->getBinning(id.subdetId(), NormXResidual, nbins, xmin, xmax);
00711 histStruct.NormResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00712 normhistoname.str().c_str(),normhistotitle.str().c_str(),
00713 nbins, xmin, xmax);
00714 }
00715 this->getBinning(id.subdetId(), XprimeResidual, nbins, xmin, xmax);
00716 histStruct.ResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00717 xprimehistoname.str().c_str(),xprimehistotitle.str().c_str(),
00718 nbins, xmin, xmax);
00719 this->getBinning(id.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
00720 histStruct.NormResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00721 normxprimehistoname.str().c_str(),normxprimehistotitle.str().c_str(),
00722 nbins, xmin, xmax);
00723
00724 if ( moduleLevelProfiles_ ) {
00725 this->getBinning(id.subdetId(), XResidualProfile, nbins, xmin, xmax);
00726
00727 histStruct.LocalX = this->bookTH1F(moduleLevelHistsTransient, tfd,
00728 localxname.str().c_str(),localxtitle.str().c_str(),
00729 nbins, xmin, xmax);
00730 histStruct.LocalY = this->bookTH1F(moduleLevelHistsTransient, tfd,
00731 localyname.str().c_str(),localytitle.str().c_str(),
00732 nbins, xmin, xmax);
00733 histStruct.ResXvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00734 resxvsxprofilename.str().c_str(),resxvsxprofiletitle.str().c_str(),
00735 nbins, xmin, xmax, ymin, ymax);
00736 histStruct.ResXvsXProfile->Sumw2();
00737 histStruct.ResXvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00738 resxvsyprofilename.str().c_str(),resxvsyprofiletitle.str().c_str(),
00739 nbins, xmin, xmax, ymin, ymax);
00740 histStruct.ResXvsYProfile->Sumw2();
00741 }
00742
00743 if( this->isPixel(subdetandlayer.first) || stripYResiduals_ ) {
00744 this->getBinning(id.subdetId(), YprimeResidual, nbins, xmin, xmax);
00745 histStruct.ResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00746 yprimehistoname.str().c_str(),yprimehistotitle.str().c_str(),
00747 nbins, xmin, xmax);
00748 if (lCoorHistOn_) {
00749 this->getBinning(id.subdetId(), YResidual, nbins, xmin, xmax);
00750 histStruct.ResYHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00751 yhistoname.str().c_str(), yhistotitle.str().c_str(),
00752 nbins, xmin, xmax);
00753 }
00754 this->getBinning(id.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
00755 histStruct.NormResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00756 normyprimehistoname.str().c_str(),normyprimehistotitle.str().c_str(),
00757 nbins, xmin, xmax);
00758
00759 if ( moduleLevelProfiles_ ) {
00760 this->getBinning(id.subdetId(), YResidualProfile, nbins, xmin, xmax);
00761
00762 histStruct.ResYvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00763 resyvsxprofilename.str().c_str(),resyvsxprofiletitle.str().c_str(),
00764 nbins, xmin, xmax, ymin, ymax);
00765 histStruct.ResYvsXProfile->Sumw2();
00766 histStruct.ResYvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00767 resyvsyprofilename.str().c_str(),resyvsyprofiletitle.str().c_str(),
00768 nbins, xmin, xmax, ymin, ymax);
00769 histStruct.ResYvsYProfile->Sumw2();
00770 }
00771 }
00772 }
00773 }
00774
00775
00776 TH1* TrackerOfflineValidation::bookTH1F(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00777 int nBinsX, double lowX, double highX)
00778 {
00779 if (isTransient) {
00780 vDeleteObjects_.push_back(new TH1F(histName, histTitle, nBinsX, lowX, highX));
00781 return vDeleteObjects_.back();
00782 }
00783 else
00784 return tfd.make<TH1F>(histName, histTitle, nBinsX, lowX, highX);
00785 }
00786
00787 TProfile* TrackerOfflineValidation::bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00788 int nBinsX, double lowX, double highX)
00789 {
00790 if (isTransient) {
00791 TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX);
00792 vDeleteObjects_.push_back(profile);
00793 return profile;
00794 }
00795 else
00796 return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX);
00797 }
00798
00799
00800 TProfile* TrackerOfflineValidation::bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle,
00801 int nBinsX, double lowX, double highX, double lowY, double highY)
00802 {
00803 if (isTransient) {
00804 TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
00805 vDeleteObjects_.push_back(profile);
00806 return profile;
00807 }
00808 else
00809 return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
00810 }
00811
00812 bool TrackerOfflineValidation::isBarrel(uint32_t subDetId)
00813 {
00814 return (subDetId == StripSubdetector::TIB ||
00815 subDetId == StripSubdetector::TOB ||
00816 subDetId == PixelSubdetector::PixelBarrel );
00817 }
00818
00819
00820 bool TrackerOfflineValidation::isEndCap(uint32_t subDetId)
00821 {
00822 return ( subDetId == StripSubdetector::TID ||
00823 subDetId == StripSubdetector::TEC ||
00824 subDetId == PixelSubdetector::PixelEndcap );
00825 }
00826
00827
00828 bool TrackerOfflineValidation::isPixel(uint32_t subDetId)
00829 {
00830 return (subDetId == PixelSubdetector::PixelBarrel ||
00831 subDetId == PixelSubdetector::PixelEndcap );
00832 }
00833
00834
00835 bool TrackerOfflineValidation::isDetOrDetUnit(align::StructureType type)
00836 {
00837 return ( type == align::AlignableDet ||
00838 type == align::AlignableDetUnit );
00839 }
00840
00841
00842 void
00843 TrackerOfflineValidation::getBinning(uint32_t subDetId,
00844 TrackerOfflineValidation::HistogrammType residualType,
00845 int& nBinsX, double& lowerBoundX, double& upperBoundX)
00846 {
00847
00848 const bool isPixel = this->isPixel(subDetId);
00849
00850 edm::ParameterSet binningPSet;
00851
00852 switch(residualType)
00853 {
00854 case XResidual :
00855 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XResPixelModules");
00856 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XResStripModules");
00857 break;
00858 case NormXResidual :
00859 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXResPixelModules");
00860 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXResStripModules");
00861 break;
00862 case XprimeResidual :
00863 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResPixelModules");
00864 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResStripModules");
00865 break;
00866 case NormXprimeResidual :
00867 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResPixelModules");
00868 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResStripModules");
00869 break;
00870 case YResidual :
00871 case YprimeResidual :
00872 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1YResPixelModules");
00873 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1YResStripModules");
00874 break;
00875
00876 case NormYprimeResidual :
00877 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormYResPixelModules");
00878 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormYResStripModules");
00879 break;
00880 case XResidualProfile :
00881 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileXResPixelModules");
00882 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileXResStripModules");
00883 break;
00884 case YResidualProfile :
00885 if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileYResPixelModules");
00886 else binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileYResStripModules");
00887 break;
00888 }
00889 nBinsX = binningPSet.getParameter<int32_t>("Nbinx");
00890 lowerBoundX = binningPSet.getParameter<double>("xmin");
00891 upperBoundX = binningPSet.getParameter<double>("xmax");
00892 }
00893
00894
00895 void
00896 TrackerOfflineValidation::setSummaryBin(int bin, TH1* targetHist, TH1* sourceHist)
00897 {
00898 if(targetHist && sourceHist) {
00899 targetHist->SetBinContent(bin, sourceHist->GetMean(1));
00900 if(useFwhm_) targetHist->SetBinError(bin, Fwhm(sourceHist)/2.);
00901 else targetHist->SetBinError(bin, sourceHist->GetRMS(1) );
00902 }
00903 else return;
00904 }
00905
00906
00907 void
00908 TrackerOfflineValidation::summarizeBinInContainer( int bin, SummaryContainer& targetContainer,
00909 SummaryContainer& sourceContainer)
00910 {
00911 this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.sumXResiduals_);
00912 this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.sumNormXResiduals_);
00913
00914 this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.sumYResiduals_);
00915 this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.sumNormYResiduals_);
00916 }
00917
00918
00919 void
00920 TrackerOfflineValidation::summarizeBinInContainer( int bin, uint32_t subDetId,
00921 SummaryContainer& targetContainer,
00922 ModuleHistos& sourceContainer)
00923 {
00924
00925 this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.ResXprimeHisto);
00926 this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.NormResXprimeHisto);
00927 if( this->isPixel(subDetId) || stripYResiduals_ ) {
00928 this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.ResYprimeHisto);
00929 this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.NormResYprimeHisto);
00930 }
00931 }
00932
00933
00934 TrackerOfflineValidation::ModuleHistos&
00935 TrackerOfflineValidation::getHistStructFromMap(const DetId& detid)
00936 {
00937
00938
00939
00940 uint subdetid = detid.subdetId();
00941 if(subdetid == PixelSubdetector::PixelBarrel ) {
00942 return mPxbResiduals_[detid.rawId()];
00943 } else if (subdetid == PixelSubdetector::PixelEndcap) {
00944 return mPxeResiduals_[detid.rawId()];
00945 } else if(subdetid == StripSubdetector::TIB) {
00946 return mTibResiduals_[detid.rawId()];
00947 } else if(subdetid == StripSubdetector::TID) {
00948 return mTidResiduals_[detid.rawId()];
00949 } else if(subdetid == StripSubdetector::TOB) {
00950 return mTobResiduals_[detid.rawId()];
00951 } else if(subdetid == StripSubdetector::TEC) {
00952 return mTecResiduals_[detid.rawId()];
00953 } else {
00954 throw cms::Exception("Geometry Error")
00955 << "[TrackerOfflineValidation] Error, tried to get reference for non-tracker subdet " << subdetid
00956 << " from detector " << detid.det();
00957 return mPxbResiduals_[0];
00958 }
00959 }
00960
00961
00962
00963 void
00964 TrackerOfflineValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00965 {
00966 if (useOverflowForRMS_)TH1::StatOverflows(kTRUE);
00967 this->checkBookHists(iSetup);
00968
00969 TrackerValidationVariables avalidator_(iSetup,parSet_);
00970
00971 std::vector<TrackerValidationVariables::AVTrackStruct> vTrackstruct;
00972 avalidator_.fillTrackQuantities(iEvent, vTrackstruct);
00973
00974 for (std::vector<TrackerValidationVariables::AVTrackStruct>::const_iterator itT = vTrackstruct.begin();
00975 itT != vTrackstruct.end();
00976 ++itT) {
00977
00978
00979 static const int etaindex = this->GetIndex(vTrackHistos_,"h_tracketa");
00980 vTrackHistos_[etaindex]->Fill(itT->eta);
00981 static const int phiindex = this->GetIndex(vTrackHistos_,"h_trackphi");
00982 vTrackHistos_[phiindex]->Fill(itT->phi);
00983 static const int numOfValidHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfValidHits");
00984 vTrackHistos_[numOfValidHitsindex]->Fill(itT->numberOfValidHits);
00985 static const int numOfLostHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfLostHits");
00986 vTrackHistos_[numOfLostHitsindex]->Fill(itT->numberOfLostHits);
00987 static const int kappaindex = this->GetIndex(vTrackHistos_,"h_curvature");
00988 vTrackHistos_[kappaindex]->Fill(itT->kappa);
00989 static const int kappaposindex = this->GetIndex(vTrackHistos_,"h_curvature_pos");
00990 if (itT->charge > 0)
00991 vTrackHistos_[kappaposindex]->Fill(fabs(itT->kappa));
00992 static const int kappanegindex = this->GetIndex(vTrackHistos_,"h_curvature_neg");
00993 if (itT->charge < 0)
00994 vTrackHistos_[kappanegindex]->Fill(fabs(itT->kappa));
00995 static const int normchi2index = this->GetIndex(vTrackHistos_,"h_normchi2");
00996 vTrackHistos_[normchi2index]->Fill(itT->normchi2);
00997 static const int chi2index = this->GetIndex(vTrackHistos_,"h_chi2");
00998 vTrackHistos_[chi2index]->Fill(itT->chi2);
00999 static const int chi2Probindex = this->GetIndex(vTrackHistos_,"h_chi2Prob");
01000 vTrackHistos_[chi2Probindex]->Fill(itT->chi2Prob);
01001 static const int ptindex = this->GetIndex(vTrackHistos_,"h_pt");
01002 vTrackHistos_[ptindex]->Fill(itT->pt);
01003 if (itT->ptError != 0.) {
01004 static const int ptResolutionindex = this->GetIndex(vTrackHistos_,"h_ptResolution");
01005 vTrackHistos_[ptResolutionindex]->Fill(itT->ptError/itT->pt);
01006 }
01007
01008 static const int d0phiindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_phi");
01009 vTrackProfiles_[d0phiindex]->Fill(itT->phi,itT->d0);
01010 static const int dzphiindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_phi");
01011 vTrackProfiles_[dzphiindex]->Fill(itT->phi,itT->dz);
01012 static const int d0etaindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_eta");
01013 vTrackProfiles_[d0etaindex]->Fill(itT->eta,itT->d0);
01014 static const int dzetaindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_eta");
01015 vTrackProfiles_[dzetaindex]->Fill(itT->eta,itT->dz);
01016 static const int chiphiindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_phi");
01017 vTrackProfiles_[chiphiindex]->Fill(itT->phi,itT->chi2);
01018 static const int chiProbphiindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_phi");
01019 vTrackProfiles_[chiProbphiindex]->Fill(itT->phi,itT->chi2Prob);
01020 static const int chiProbabsd0index = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_d0");
01021 vTrackProfiles_[chiProbabsd0index]->Fill(fabs(itT->d0),itT->chi2Prob);
01022 static const int normchiphiindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_phi");
01023 vTrackProfiles_[normchiphiindex]->Fill(itT->phi,itT->normchi2);
01024 static const int chietaindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_eta");
01025 vTrackProfiles_[chietaindex]->Fill(itT->eta,itT->chi2);
01026 static const int normchiptindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_pt");
01027 vTrackProfiles_[normchiptindex]->Fill(itT->pt,itT->normchi2);
01028 static const int normchipindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_p");
01029 vTrackProfiles_[normchipindex]->Fill(itT->p,itT->normchi2);
01030 static const int chiProbetaindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_eta");
01031 vTrackProfiles_[chiProbetaindex]->Fill(itT->eta,itT->chi2Prob);
01032 static const int normchietaindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_eta");
01033 vTrackProfiles_[normchietaindex]->Fill(itT->eta,itT->normchi2);
01034 static const int kappaphiindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_phi");
01035 vTrackProfiles_[kappaphiindex]->Fill(itT->phi,itT->kappa);
01036 static const int kappaetaindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_eta");
01037 vTrackProfiles_[kappaetaindex]->Fill(itT->eta,itT->kappa);
01038 static const int ptResphiindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_phi");
01039 vTrackProfiles_[ptResphiindex]->Fill(itT->phi,itT->ptError/itT->pt);
01040 static const int ptResetaindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_eta");
01041 vTrackProfiles_[ptResetaindex]->Fill(itT->eta,itT->ptError/itT->pt);
01042
01043
01044 static const int d0phiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_phi");
01045 vTrack2DHistos_[d0phiindex_2d]->Fill(itT->phi,itT->d0);
01046 static const int dzphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_phi");
01047 vTrack2DHistos_[dzphiindex_2d]->Fill(itT->phi,itT->dz);
01048 static const int d0etaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_eta");
01049 vTrack2DHistos_[d0etaindex_2d]->Fill(itT->eta,itT->d0);
01050 static const int dzetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_eta");
01051 vTrack2DHistos_[dzetaindex_2d]->Fill(itT->eta,itT->dz);
01052 static const int chiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_phi");
01053 vTrack2DHistos_[chiphiindex_2d]->Fill(itT->phi,itT->chi2);
01054 static const int chiProbphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_phi");
01055 vTrack2DHistos_[chiProbphiindex_2d]->Fill(itT->phi,itT->chi2Prob);
01056 static const int chiProbabsd0index_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_d0");
01057 vTrack2DHistos_[chiProbabsd0index_2d]->Fill(fabs(itT->d0),itT->chi2Prob);
01058 static const int normchiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_phi");
01059 vTrack2DHistos_[normchiphiindex_2d]->Fill(itT->phi,itT->normchi2);
01060 static const int chietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_eta");
01061 vTrack2DHistos_[chietaindex_2d]->Fill(itT->eta,itT->chi2);
01062 static const int chiProbetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_eta");
01063 vTrack2DHistos_[chiProbetaindex_2d]->Fill(itT->eta,itT->chi2Prob);
01064 static const int normchietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_eta");
01065 vTrack2DHistos_[normchietaindex_2d]->Fill(itT->eta,itT->normchi2);
01066 static const int kappaphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_phi");
01067 vTrack2DHistos_[kappaphiindex_2d]->Fill(itT->phi,itT->kappa);
01068 static const int kappaetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_eta");
01069 vTrack2DHistos_[kappaetaindex_2d]->Fill(itT->eta,itT->kappa);
01070 static const int normchi2kappa_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_kappa");
01071 vTrack2DHistos_[normchi2kappa_2d]->Fill(itT->normchi2,itT->kappa);
01072
01073
01074 for (std::vector<TrackerValidationVariables::AVHitStruct>::const_iterator itH = itT->hits.begin();
01075 itH != itT->hits.end();
01076 ++itH) {
01077
01078 DetId detid(itH->rawDetId);
01079 ModuleHistos &histStruct = this->getHistStructFromMap(detid);
01080
01081
01082 if (lCoorHistOn_) {
01083 histStruct.ResHisto->Fill(itH->resX);
01084 if(itH->resErrX != 0) histStruct.NormResHisto->Fill(itH->resX/itH->resErrX);
01085 if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01086 histStruct.ResYHisto->Fill(itH->resY);
01087
01088 }
01089 }
01090 if (itH->resXprime != -999.) {
01091 histStruct.ResXprimeHisto->Fill(itH->resXprime);
01092
01093
01094 if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
01095 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixB");
01096 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01097 }
01098 if (detid.subdetId() == PixelSubdetector::PixelEndcap) {
01099 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixE");
01100 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01101 }
01102 if (detid.subdetId() == StripSubdetector::TIB) {
01103 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TIB");
01104 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01105 }
01106 if (detid.subdetId() == StripSubdetector::TID) {
01107 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TID");
01108 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01109 }
01110 if (detid.subdetId() == StripSubdetector::TOB) {
01111 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TOB");
01112 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01113 }
01114 if (detid.subdetId() == StripSubdetector::TEC) {
01115 static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TEC");
01116 vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01117 }
01118
01119
01120 if ( moduleLevelProfiles_ && itH->inside ) {
01121 float tgalpha = tan(itH->localAlpha);
01122 if ( fabs(tgalpha)!=0 ){
01123 histStruct.LocalX->Fill(itH->localXnorm, tgalpha*tgalpha);
01124 histStruct.LocalY->Fill(itH->localYnorm, tgalpha*tgalpha);
01125
01126
01127
01128
01129
01130
01131
01132
01133
01134
01135
01136 histStruct.ResXvsXProfile->Fill(itH->localXnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
01137 histStruct.ResXvsYProfile->Fill(itH->localYnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
01138
01139
01140
01141 }
01142 }
01143
01144 if(itH->resXprimeErr != 0 && itH->resXprimeErr != -999 ) {
01145 histStruct.NormResXprimeHisto->Fill(itH->resXprime/itH->resXprimeErr);
01146 }
01147 }
01148
01149 if (itH->resYprime != -999.) {
01150 if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01151 histStruct.ResYprimeHisto->Fill(itH->resYprime);
01152
01153
01154 if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
01155 static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixB");
01156 vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
01157 }
01158 if (detid.subdetId() == PixelSubdetector::PixelEndcap) {
01159 static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixE");
01160 vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
01161 }
01162
01163
01164 if ( moduleLevelProfiles_ && itH->inside ) {
01165 float tgbeta = tan(itH->localBeta);
01166 if ( fabs(tgbeta)!=0 ){
01167
01168
01169
01170
01171
01172
01173
01174
01175
01176
01177
01178
01179 histStruct.ResYvsXProfile->Fill(itH->localXnorm, itH->resY/tgbeta, tgbeta*tgbeta);
01180 histStruct.ResYvsYProfile->Fill(itH->localYnorm, itH->resY/tgbeta, tgbeta*tgbeta);
01181
01182 }
01183 }
01184
01185 if (itH->resYprimeErr != 0 && itH->resYprimeErr != -999. ) {
01186 histStruct.NormResYprimeHisto->Fill(itH->resYprime/itH->resYprimeErr);
01187 }
01188 }
01189 }
01190
01191 }
01192 }
01193
01194 if (useOverflowForRMS_) TH1::StatOverflows(kFALSE);
01195 }
01196
01197
01198
01199 void
01200 TrackerOfflineValidation::endJob()
01201 {
01202
01203 if (!tkGeom_.product()) return;
01204
01205 AlignableTracker aliTracker(&(*tkGeom_));
01206
01207 static const int kappadiffindex = this->GetIndex(vTrackHistos_,"h_diff_curvature");
01208 vTrackHistos_[kappadiffindex]->Add(vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_neg")],
01209 vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_pos")],-1,1);
01210
01211
01212
01213 std::vector<TrackerOfflineValidation::SummaryContainer> vTrackerprofiles;
01214 DirectoryWrapper f("",moduleDirectory_,dqmMode_);
01215 this->collateSummaryHists(f,(aliTracker), 0, vTrackerprofiles);
01216
01217 if (dqmMode_) return;
01218
01219
01220
01221 edm::Service<TFileService> fs;
01222 TTree *tree = fs->make<TTree>("TkOffVal","TkOffVal");
01223
01224 TkOffTreeVariables *treeMemPtr = new TkOffTreeVariables;
01225
01226
01227
01228 tree->Branch("TkOffTreeVariables", &treeMemPtr);
01229
01230 this->fillTree(*tree, mPxbResiduals_, *treeMemPtr, *tkGeom_);
01231 this->fillTree(*tree, mPxeResiduals_, *treeMemPtr, *tkGeom_);
01232 this->fillTree(*tree, mTibResiduals_, *treeMemPtr, *tkGeom_);
01233 this->fillTree(*tree, mTidResiduals_, *treeMemPtr, *tkGeom_);
01234 this->fillTree(*tree, mTobResiduals_, *treeMemPtr, *tkGeom_);
01235 this->fillTree(*tree, mTecResiduals_, *treeMemPtr, *tkGeom_);
01236
01237 delete treeMemPtr; treeMemPtr = 0;
01238 }
01239
01240
01241 void
01242 TrackerOfflineValidation::collateSummaryHists( DirectoryWrapper& tfd, const Alignable& ali, int i,
01243 std::vector<TrackerOfflineValidation::SummaryContainer>& vLevelProfiles)
01244 {
01245 std::vector<Alignable*> alivec(ali.components());
01246 if( this->isDetOrDetUnit((alivec)[0]->alignableObjectId()) ) return;
01247
01248 for(int iComp=0, iCompEnd = ali.components().size();iComp < iCompEnd; ++iComp) {
01249 std::vector< TrackerOfflineValidation::SummaryContainer > vProfiles;
01250 std::string structurename = AlignableObjectId::idToString((alivec)[iComp]->alignableObjectId());
01251
01252 LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
01253 std::stringstream dirname;
01254 dirname << structurename;
01255
01256
01257 if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << iComp+1;
01258
01259 if( !(this->isDetOrDetUnit( (alivec)[iComp]->alignableObjectId()) )
01260 || (alivec)[0]->components().size() > 1 ) {
01261 DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
01262 this->collateSummaryHists( f, *(alivec)[iComp], i, vProfiles);
01263 vLevelProfiles.push_back(this->bookSummaryHists(tfd, *(alivec[iComp]), ali.alignableObjectId(), iComp+1));
01264 TH1 *hY = vLevelProfiles[iComp].sumYResiduals_;
01265 TH1 *hNormY = vLevelProfiles[iComp].sumNormYResiduals_;
01266 for(uint n = 0; n < vProfiles.size(); ++n) {
01267 this->summarizeBinInContainer(n+1, vLevelProfiles[iComp], vProfiles[n] );
01268 vLevelProfiles[iComp].sumXResiduals_->Add(vProfiles[n].sumXResiduals_);
01269 vLevelProfiles[iComp].sumNormXResiduals_->Add(vProfiles[n].sumNormXResiduals_);
01270 if (hY) hY->Add(vProfiles[n].sumYResiduals_);
01271 if (hNormY) hNormY->Add(vProfiles[n].sumNormYResiduals_);
01272 }
01273 if(dqmMode_)continue;
01274
01275 this->fitResiduals(vLevelProfiles[iComp].sumXResiduals_);
01276 this->fitResiduals(vLevelProfiles[iComp].sumNormXResiduals_);
01277 if (hY) this->fitResiduals(hY);
01278 if (hNormY) this->fitResiduals(hNormY);
01279 } else {
01280
01281 continue;
01282 }
01283 }
01284 }
01285
01286
01287 TrackerOfflineValidation::SummaryContainer
01288 TrackerOfflineValidation::bookSummaryHists(DirectoryWrapper& tfd, const Alignable& ali,
01289 align::StructureType type, int i)
01290 {
01291 const uint aliSize = ali.components().size();
01292 const align::StructureType alitype = ali.alignableObjectId();
01293 const align::StructureType subtype = ali.components()[0]->alignableObjectId();
01294 const char *aliTypeName = AlignableObjectId::idToString(alitype);
01295 const char *aliSubtypeName = AlignableObjectId::idToString(subtype);
01296 const char *typeName = AlignableObjectId::idToString(type);
01297
01298 const DetId aliDetId = ali.id();
01299
01300 const bool bookResidY = this->isPixel(aliDetId.subdetId()) || stripYResiduals_;
01301
01302 SummaryContainer sumContainer;
01303
01304
01305
01306
01307
01308 const uint subcompSize = ali.components()[0]->components().size();
01309 if (subtype != align::AlignableDet || subcompSize == 1) {
01310 const TString title(Form("Summary for substructures in %s %d;%s;",aliTypeName,i,aliSubtypeName));
01311
01312 sumContainer.summaryXResiduals_ = tfd.make<TH1F>(Form("h_summaryX%s_%d",aliTypeName,i),
01313 title + "#LT #Delta x' #GT",
01314 aliSize, 0.5, aliSize+0.5);
01315 sumContainer.summaryNormXResiduals_ = tfd.make<TH1F>(Form("h_summaryNormX%s_%d",aliTypeName,i),
01316 title + "#LT #Delta x'/#sigma #GT",
01317 aliSize,0.5,aliSize+0.5);
01318
01319 if (bookResidY) {
01320 sumContainer.summaryYResiduals_ = tfd.make<TH1F>(Form("h_summaryY%s_%d",aliTypeName,i),
01321 title + "#LT #Delta y' #GT",
01322 aliSize, 0.5, aliSize+0.5);
01323 sumContainer.summaryNormYResiduals_ = tfd.make<TH1F>(Form("h_summaryNormY%s_%d",aliTypeName,i),
01324 title + "#LT #Delta y'/#sigma #GT",
01325 aliSize,0.5,aliSize+0.5);
01326 }
01327
01328 } else if (subtype == align::AlignableDet && subcompSize > 1) {
01329 if (subcompSize != 2) {
01330
01331 edm::LogError("Alignment") << "@SUB=bookSummaryHists"
01332 << "Det with " << subcompSize << " components";
01333 }
01334
01335 const TString title(Form("Summary for substructures in %s %d;%s;", aliTypeName, i,
01336 AlignableObjectId::idToString(ali.components()[0]->components()[0]->alignableObjectId())));
01337
01338 sumContainer.summaryXResiduals_
01339 = tfd.make<TH1F>(Form("h_summaryX%s_%d", aliTypeName, i),
01340 title + "#LT #Delta x' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01341 sumContainer.summaryNormXResiduals_
01342 = tfd.make<TH1F>(Form("h_summaryNormX%s_%d", aliTypeName, i),
01343 title + "#LT #Delta x'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01344
01345 if (bookResidY) {
01346 sumContainer.summaryYResiduals_
01347 = tfd.make<TH1F>(Form("h_summaryY%s_%d", aliTypeName, i),
01348 title + "#LT #Delta y' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01349 sumContainer.summaryNormYResiduals_
01350 = tfd.make<TH1F>(Form("h_summaryNormY%s_%d", aliTypeName, i),
01351 title + "#LT #Delta y'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01352 }
01353
01354 } else {
01355 edm::LogError("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookSummaryHists"
01356 << "No summary histogramm for hierarchy level "
01357 << aliTypeName << " in subdet " << aliDetId.subdetId();
01358 }
01359
01360
01361
01362
01363 int nbins = 0;
01364 double xmin = 0., xmax = 0.;
01365 const TString sumTitle(Form("Residual for %s %d in %s;", aliTypeName, i, typeName));
01366 const ModuleHistos &xTitHists = this->getHistStructFromMap(aliDetId);
01367 this->getBinning(aliDetId.subdetId(), XprimeResidual, nbins, xmin, xmax);
01368
01369 sumContainer.sumXResiduals_ = tfd.make<TH1F>(Form("h_Xprime_%s_%d", aliTypeName, i),
01370 sumTitle + xTitHists.ResXprimeHisto->GetXaxis()->GetTitle(),
01371 nbins, xmin, xmax);
01372
01373 this->getBinning(aliDetId.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
01374 sumContainer.sumNormXResiduals_ = tfd.make<TH1F>(Form("h_NormXprime_%s_%d",aliTypeName,i),
01375 sumTitle + xTitHists.NormResXprimeHisto->GetXaxis()->GetTitle(),
01376 nbins, xmin, xmax);
01377 if (bookResidY) {
01378 this->getBinning(aliDetId.subdetId(), YprimeResidual, nbins, xmin, xmax);
01379 sumContainer.sumYResiduals_ = tfd.make<TH1F>(Form("h_Yprime_%s_%d",aliTypeName,i),
01380 sumTitle + xTitHists.ResYprimeHisto->GetXaxis()->GetTitle(),
01381 nbins, xmin, xmax);
01382
01383 this->getBinning(aliDetId.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
01384 sumContainer.sumNormYResiduals_ = tfd.make<TH1F>(Form("h_NormYprime_%s_%d",aliTypeName,i),
01385 sumTitle + xTitHists.NormResYprimeHisto->GetXaxis()->GetTitle(),
01386 nbins, xmin, xmax);
01387 }
01388
01389
01390
01391
01392 if( ( subtype == align::AlignableDet && subcompSize == 1) || subtype == align::AlignableDetUnit ) {
01393 for(uint k = 0; k < aliSize; ++k) {
01394 DetId detid = ali.components()[k]->id();
01395 ModuleHistos &histStruct = this->getHistStructFromMap(detid);
01396 this->summarizeBinInContainer(k+1, detid.subdetId() ,sumContainer, histStruct );
01397 sumContainer.sumXResiduals_->Add(histStruct.ResXprimeHisto);
01398 sumContainer.sumNormXResiduals_->Add(histStruct.NormResXprimeHisto);
01399 if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01400 sumContainer.sumYResiduals_->Add(histStruct.ResYprimeHisto);
01401 sumContainer.sumNormYResiduals_->Add(histStruct.NormResYprimeHisto);
01402 }
01403 }
01404 } else if( subtype == align::AlignableDet && subcompSize > 1) {
01405
01406 for(uint k = 0; k < aliSize; ++k) {
01407 for(uint j = 0; j < subcompSize; ++j) {
01408 DetId detid = ali.components()[k]->components()[j]->id();
01409 ModuleHistos &histStruct = this->getHistStructFromMap(detid);
01410 this->summarizeBinInContainer(2*k+j+1, detid.subdetId() ,sumContainer, histStruct );
01411 sumContainer.sumXResiduals_->Add( histStruct.ResXprimeHisto);
01412 sumContainer.sumNormXResiduals_->Add( histStruct.NormResXprimeHisto);
01413 if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01414 sumContainer.sumYResiduals_->Add( histStruct.ResYprimeHisto);
01415 sumContainer.sumNormYResiduals_->Add( histStruct.NormResYprimeHisto);
01416 }
01417 }
01418 }
01419 }
01420 return sumContainer;
01421 }
01422
01423
01424 float
01425 TrackerOfflineValidation::Fwhm (const TH1* hist) const
01426 {
01427 float max = hist->GetMaximum();
01428 int left = -1, right = -1;
01429 for(unsigned int i = 1, iEnd = hist->GetNbinsX(); i <= iEnd; ++i) {
01430 if(hist->GetBinContent(i) < max/2. && hist->GetBinContent(i+1) > max/2. && left == -1) {
01431 if(max/2. - hist->GetBinContent(i) < hist->GetBinContent(i+1) - max/2.) {
01432 left = i;
01433 ++i;
01434 } else {
01435 left = i+1;
01436 ++i;
01437 }
01438 }
01439 if(left != -1 && right == -1) {
01440 if(hist->GetBinContent(i) > max/2. && hist->GetBinContent(i+1) < max/2.) {
01441 if( hist->GetBinContent(i) - max/2. < max/2. - hist->GetBinContent(i+1)) {
01442 right = i;
01443 } else {
01444 right = i+1;
01445 }
01446
01447 }
01448 }
01449 }
01450 return hist->GetXaxis()->GetBinCenter(right) - hist->GetXaxis()->GetBinCenter(left);
01451 }
01452
01453
01454 void
01455 TrackerOfflineValidation::fillTree(TTree& tree,
01456 const std::map<int, TrackerOfflineValidation::ModuleHistos>& moduleHist_,
01457 TkOffTreeVariables &treeMem, const TrackerGeometry& tkgeom)
01458 {
01459
01460 for(std::map<int, TrackerOfflineValidation::ModuleHistos>::const_iterator it = moduleHist_.begin(),
01461 itEnd= moduleHist_.end(); it != itEnd;++it ) {
01462 treeMem.clear();
01463
01464
01465 DetId detId_ = it->first;
01466 treeMem.moduleId = detId_;
01467 treeMem.subDetId = detId_.subdetId();
01468 treeMem.isDoubleSide =0;
01469
01470 if(treeMem.subDetId == PixelSubdetector::PixelBarrel){
01471 PXBDetId pxbId(detId_);
01472 unsigned int whichHalfBarrel(1), rawId(detId_.rawId());
01473 if( (rawId>=302056964 && rawId<302059300) || (rawId>=302123268 && rawId<302127140) ||
01474 (rawId>=302189572 && rawId<302194980) ) whichHalfBarrel=2;
01475 treeMem.layer = pxbId.layer();
01476 treeMem.half = whichHalfBarrel;
01477 treeMem.rod = pxbId.ladder();
01478 treeMem.module = pxbId.module();
01479 } else if(treeMem.subDetId == PixelSubdetector::PixelEndcap){
01480 PXFDetId pxfId(detId_);
01481 unsigned int whichHalfCylinder(1), rawId(detId_.rawId());
01482 if( (rawId>=352394500 && rawId<352406032) || (rawId>=352460036 && rawId<352471568) ||
01483 (rawId>=344005892 && rawId<344017424) || (rawId>=344071428 && rawId<344082960) ) whichHalfCylinder=2;
01484 treeMem.layer = pxfId.disk();
01485 treeMem.side = pxfId.side();
01486 treeMem.half = whichHalfCylinder;
01487 treeMem.blade = pxfId.blade();
01488 treeMem.panel = pxfId.panel();
01489 treeMem.module = pxfId.module();
01490 } else if(treeMem.subDetId == StripSubdetector::TIB){
01491 TIBDetId tibId(detId_);
01492 unsigned int whichHalfShell(1), rawId(detId_.rawId());
01493 if ( (rawId>=369120484 && rawId<369120688) || (rawId>=369121540 && rawId<369121776) ||
01494 (rawId>=369136932 && rawId<369137200) || (rawId>=369137988 && rawId<369138288) ||
01495 (rawId>=369153396 && rawId<369153744) || (rawId>=369154436 && rawId<369154800) ||
01496 (rawId>=369169844 && rawId<369170256) || (rawId>=369170900 && rawId<369171344) ||
01497 (rawId>=369124580 && rawId<369124784) || (rawId>=369125636 && rawId<369125872) ||
01498 (rawId>=369141028 && rawId<369141296) || (rawId>=369142084 && rawId<369142384) ||
01499 (rawId>=369157492 && rawId<369157840) || (rawId>=369158532 && rawId<369158896) ||
01500 (rawId>=369173940 && rawId<369174352) || (rawId>=369174996 && rawId<369175440) ) whichHalfShell=2;
01501 treeMem.layer = tibId.layer();
01502 treeMem.side = tibId.string()[0];
01503 treeMem.half = whichHalfShell;
01504 treeMem.rod = tibId.string()[2];
01505 treeMem.outerInner = tibId.string()[1];
01506 treeMem.module = tibId.module();
01507 treeMem.isStereo = tibId.stereo();
01508 treeMem.isDoubleSide = tibId.isDoubleSide();
01509 } else if(treeMem.subDetId == StripSubdetector::TID){
01510 TIDDetId tidId(detId_);
01511 treeMem.layer = tidId.wheel();
01512 treeMem.side = tidId.side();
01513 treeMem.ring = tidId.ring();
01514 treeMem.outerInner = tidId.module()[0];
01515 treeMem.module = tidId.module()[1];
01516 treeMem.isStereo = tidId.stereo();
01517 treeMem.isDoubleSide = tidId.isDoubleSide();
01518 } else if(treeMem.subDetId == StripSubdetector::TOB){
01519 TOBDetId tobId(detId_);
01520 treeMem.layer = tobId.layer();
01521 treeMem.side = tobId.rod()[0];
01522 treeMem.rod = tobId.rod()[1];
01523 treeMem.module = tobId.module();
01524 treeMem.isStereo = tobId.stereo();
01525 treeMem.isDoubleSide = tobId.isDoubleSide();
01526 } else if(treeMem.subDetId == StripSubdetector::TEC) {
01527 TECDetId tecId(detId_);
01528 treeMem.layer = tecId.wheel();
01529 treeMem.side = tecId.side();
01530 treeMem.ring = tecId.ring();
01531 treeMem.petal = tecId.petal()[1];
01532 treeMem.outerInner = tecId.petal()[0];
01533 treeMem.module = tecId.module();
01534 treeMem.isStereo = tecId.stereo();
01535 treeMem.isDoubleSide = tecId.isDoubleSide();
01536 }
01537
01538
01539 const Surface::PositionType &gPModule = tkgeom.idToDet(detId_)->position();
01540 treeMem.posPhi = gPModule.phi();
01541 treeMem.posEta = gPModule.eta();
01542 treeMem.posR = gPModule.perp();
01543 treeMem.posX = gPModule.x();
01544 treeMem.posY = gPModule.y();
01545 treeMem.posZ = gPModule.z();
01546
01547 const Surface& surface = tkgeom.idToDet(detId_)->surface();
01548
01549
01550 LocalPoint lUDirection(1.,0.,0.), lVDirection(0.,1.,0.), lWDirection(0.,0.,1.);
01551 GlobalPoint gUDirection = surface.toGlobal(lUDirection),
01552 gVDirection = surface.toGlobal(lVDirection),
01553 gWDirection = surface.toGlobal(lWDirection);
01554 double dR(999.), dPhi(999.), dZ(999.);
01555 if(treeMem.subDetId==PixelSubdetector::PixelBarrel || treeMem.subDetId==StripSubdetector::TIB || treeMem.subDetId==StripSubdetector::TOB){
01556 dR = gWDirection.perp() - gPModule.perp();
01557 dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
01558 dZ = gVDirection.z() - gPModule.z();
01559 if(dZ>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01560 }else if(treeMem.subDetId==PixelSubdetector::PixelEndcap){
01561 dR = gUDirection.perp() - gPModule.perp();
01562 dPhi = deltaPhi(gVDirection.phi(),gPModule.phi());
01563 dZ = gWDirection.z() - gPModule.z();
01564 if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01565 }else if(treeMem.subDetId==StripSubdetector::TID || treeMem.subDetId==StripSubdetector::TEC){
01566 dR = gVDirection.perp() - gPModule.perp();
01567 dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
01568 dZ = gWDirection.z() - gPModule.z();
01569 if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01570 }
01571 if(dR>=0.)treeMem.rDirection = 1; else treeMem.rDirection = -1;
01572 if(dPhi>=0.)treeMem.phiDirection = 1; else treeMem.phiDirection = -1;
01573 if(dZ>=0.)treeMem.zDirection = 1; else treeMem.zDirection = -1;
01574
01575
01576 treeMem.entries = static_cast<UInt_t>(it->second.ResXprimeHisto->GetEntries());
01577 treeMem.meanX = it->second.ResXprimeHisto->GetMean();
01578 treeMem.rmsX = it->second.ResXprimeHisto->GetRMS();
01579
01580
01581 if (useFit_) {
01582
01583
01584 std::pair<float,float> fitResult1 = this->fitResiduals(it->second.ResXprimeHisto);
01585 treeMem.fitMeanX = fitResult1.first;
01586 treeMem.fitSigmaX = fitResult1.second;
01587
01588 std::pair<float,float> fitResult2 = this->fitResiduals(it->second.NormResXprimeHisto);
01589 treeMem.fitMeanNormX = fitResult2.first;
01590 treeMem.fitSigmaNormX = fitResult2.second;
01591 }
01592
01593
01594 treeMem.medianX = this->getMedian(it->second.ResXprimeHisto);
01595
01596 int numberOfBins=it->second.ResXprimeHisto->GetNbinsX();
01597 treeMem.numberOfUnderflows = it->second.ResXprimeHisto->GetBinContent(0);
01598 treeMem.numberOfOverflows = it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
01599 treeMem.numberOfOutliers = it->second.ResXprimeHisto->GetBinContent(0)+it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
01600
01601
01602 treeMem.meanNormX = it->second.NormResXprimeHisto->GetMean();
01603 treeMem.rmsNormX = it->second.NormResXprimeHisto->GetRMS();
01604
01605 double stats[20];
01606 it->second.NormResXprimeHisto->GetStats(stats);
01607
01608 if (stats[0]) treeMem.chi2PerDofX = stats[3]/stats[0];
01609
01610 treeMem.sigmaNormX = Fwhm(it->second.NormResXprimeHisto)/2.355;
01611 treeMem.histNameX = it->second.ResXprimeHisto->GetName();
01612 treeMem.histNameNormX = it->second.NormResXprimeHisto->GetName();
01613
01614
01615 if(lCoorHistOn_) {
01616 treeMem.meanLocalX = it->second.ResHisto->GetMean();
01617 treeMem.rmsLocalX = it->second.ResHisto->GetRMS();
01618 treeMem.meanNormLocalX = it->second.NormResHisto->GetMean();
01619 treeMem.rmsNormLocalX = it->second.NormResHisto->GetRMS();
01620
01621 treeMem.histNameLocalX = it->second.ResHisto->GetName();
01622 treeMem.histNameNormLocalX = it->second.NormResHisto->GetName();
01623 if (it->second.ResYHisto) treeMem.histNameLocalY = it->second.ResYHisto->GetName();
01624 }
01625
01626
01627
01628 if (it->second.ResYprimeHisto) {
01629 TH1 *h = it->second.ResYprimeHisto;
01630 treeMem.meanY = h->GetMean();
01631 treeMem.rmsY = h->GetRMS();
01632
01633 if (useFit_) {
01634 std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
01635 treeMem.fitMeanY = fitMeanSigma.first;
01636 treeMem.fitSigmaY = fitMeanSigma.second;
01637 }
01638
01639
01640 treeMem.medianY = this->getMedian(h);
01641
01642 treeMem.histNameY = h->GetName();
01643 }
01644 if (it->second.NormResYprimeHisto) {
01645 TH1 *h = it->second.NormResYprimeHisto;
01646 treeMem.meanNormY = h->GetMean();
01647 treeMem.rmsNormY = h->GetRMS();
01648 h->GetStats(stats);
01649 if (stats[0]) treeMem.chi2PerDofY = stats[3]/stats[0];
01650
01651 if (useFit_) {
01652 std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
01653 treeMem.fitMeanNormY = fitMeanSigma.first;
01654 treeMem.fitSigmaNormY = fitMeanSigma.second;
01655 }
01656 treeMem.histNameNormY = h->GetName();
01657 }
01658
01659 if (moduleLevelProfiles_) {
01660 if (it->second.ResXvsXProfile) {
01661 TH1 *h = it->second.ResXvsXProfile;
01662 treeMem.meanResXvsX = h->GetMean();
01663 treeMem.rmsResXvsX = h->GetRMS();
01664 treeMem.profileNameResXvsX = h->GetName();
01665 }
01666 if (it->second.ResXvsYProfile) {
01667 TH1 *h = it->second.ResXvsYProfile;
01668 treeMem.meanResXvsY = h->GetMean();
01669 treeMem.rmsResXvsY = h->GetRMS();
01670 treeMem.profileNameResXvsY = h->GetName();
01671 }
01672 if (it->second.ResYvsXProfile) {
01673 TH1 *h = it->second.ResYvsXProfile;
01674 treeMem.meanResYvsX = h->GetMean();
01675 treeMem.rmsResYvsX = h->GetRMS();
01676 treeMem.profileNameResYvsX = h->GetName();
01677 }
01678 if (it->second.ResYvsYProfile) {
01679 TH1 *h = it->second.ResYvsYProfile;
01680 treeMem.meanResYvsY = h->GetMean();
01681 treeMem.rmsResYvsY = h->GetRMS();
01682 treeMem.profileNameResYvsY = h->GetName();
01683 }
01684 }
01685
01686 tree.Fill();
01687 }
01688 }
01689
01690
01691 std::pair<float,float>
01692 TrackerOfflineValidation::fitResiduals(TH1* hist) const
01693 {
01694 std::pair<float,float> fitResult(9999., 9999.);
01695 if (!hist || hist->GetEntries() < 20) return fitResult;
01696
01697 float mean = hist->GetMean();
01698 float sigma = hist->GetRMS();
01699
01700 try {
01701
01702
01703 TF1 func("tmp", "gaus", mean - 2.*sigma, mean + 2.*sigma);
01704 if (0 == hist->Fit(&func,"QNR")) {
01705 mean = func.GetParameter(1);
01706 sigma = func.GetParameter(2);
01707
01708 func.SetRange(mean - 3.*sigma, mean + 3.*sigma);
01709
01710
01711 if (0 == hist->Fit(&func, "Q0LR")) {
01712 if (hist->GetFunction(func.GetName())) {
01713 hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
01714 }
01715 fitResult.first = func.GetParameter(1);
01716 fitResult.second = func.GetParameter(2);
01717 }
01718 }
01719 } catch (cms::Exception const & e) {
01720 edm::LogWarning("Alignment") << "@SUB=TrackerOfflineValidation::fitResiduals"
01721 << "Caught this exception during ROOT fit: "
01722 << e.what();
01723 }
01724 return fitResult;
01725 }
01726
01727
01728 float
01729 TrackerOfflineValidation::getMedian(const TH1* histo) const
01730 {
01731 float median = 999;
01732 int nbins = histo->GetNbinsX();
01733
01734
01735 double *x = new double[nbins];
01736 double *y = new double[nbins];
01737 for (int j = 0; j < nbins; j++) {
01738 x[j] = histo->GetBinCenter(j+1);
01739 y[j] = histo->GetBinContent(j+1);
01740 }
01741 median = TMath::Median(nbins, x, y);
01742
01743 delete[] x; x = 0;
01744 delete [] y; y = 0;
01745
01746 return median;
01747
01748 }
01749
01750 DEFINE_FWK_MODULE(TrackerOfflineValidation);