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