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