CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
TrackerOfflineValidation Class Reference

#include <Alignment/Validator/src/TrackerOfflineValidation.cc>

Inheritance diagram for TrackerOfflineValidation:
edm::EDAnalyzer edm::EDConsumerBase

Classes

struct  DirectoryWrapper
 
struct  ModuleHistos
 
struct  SummaryContainer
 

Public Types

enum  HistogrammType {
  XResidual, NormXResidual, YResidual, XprimeResidual,
  NormXprimeResidual, YprimeResidual, NormYprimeResidual, XResidualProfile,
  YResidualProfile
}
 
- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 

Public Member Functions

 TrackerOfflineValidation (const edm::ParameterSet &)
 
 ~TrackerOfflineValidation ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookDirHists (DirectoryWrapper &tfd, const Alignable &ali, const TrackerTopology *tTopo)
 
void bookGlobalHists (DirectoryWrapper &tfd)
 
void bookHists (DirectoryWrapper &tfd, const Alignable &ali, const TrackerTopology *tTopo, align::StructureType type, int i)
 
TrackerOfflineValidation::SummaryContainer bookSummaryHists (DirectoryWrapper &tfd, const Alignable &ali, align::StructureType type, int i)
 
TH1 * bookTH1F (bool isTransient, DirectoryWrapper &tfd, const char *histName, const char *histTitle, int nBinsX, double lowX, double highX)
 
TProfile * bookTProfile (bool isTransient, DirectoryWrapper &tfd, const char *histName, const char *histTitle, int nBinsX, double lowX, double highX)
 
TProfile * bookTProfile (bool isTransient, DirectoryWrapper &tfd, const char *histName, const char *histTitle, int nBinsX, double lowX, double highX, double lowY, double highY)
 
virtual void checkBookHists (const edm::EventSetup &setup)
 
void collateSummaryHists (DirectoryWrapper &tfd, const Alignable &ali, int i, std::vector< TrackerOfflineValidation::SummaryContainer > &vLevelProfiles)
 
virtual void endJob () override
 
void fillTree (TTree &tree, const std::map< int, TrackerOfflineValidation::ModuleHistos > &moduleHist_, TkOffTreeVariables &treeMem, const TrackerGeometry &tkgeom, const TrackerTopology *tTopo)
 
std::pair< float, float > fitResiduals (TH1 *hist) const
 
float Fwhm (const TH1 *hist) const
 
void getBinning (uint32_t subDetId, TrackerOfflineValidation::HistogrammType residualtype, int &nBinsX, double &lowerBoundX, double &upperBoundX)
 
ModuleHistosgetHistStructFromMap (const DetId &detid)
 
template<class OBJECT_TYPE >
int GetIndex (const std::vector< OBJECT_TYPE * > &vec, const TString &name)
 
float getMedian (const TH1 *hist) const
 
bool isBarrel (uint32_t subDetId)
 
bool isDetOrDetUnit (align::StructureType type)
 
bool isEndCap (uint32_t subDetId)
 
bool isPixel (uint32_t subDetId)
 
void setSummaryBin (int bin, TH1 *targetHist, TH1 *sourceHist)
 
void summarizeBinInContainer (int bin, SummaryContainer &targetContainer, SummaryContainer &sourceContainer)
 
void summarizeBinInContainer (int bin, uint32_t subDetId, SummaryContainer &targetContainer, ModuleHistos &sourceContainer)
 

Private Attributes

TrackerValidationVariables avalidator_
 
const TrackerGeometrybareTkGeomPtr_
 
const bool dqmMode_
 
const edm::EventSetuplastSetup_
 
const bool lCoorHistOn_
 
const std::string moduleDirectory_
 
const bool moduleLevelHistsTransient_
 
const bool moduleLevelProfiles_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mPxbResiduals_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mPxeResiduals_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mTecResiduals_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mTibResiduals_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mTidResiduals_
 
std::map< int,
TrackerOfflineValidation::ModuleHistos
mTobResiduals_
 
const edm::ParameterSet parSet_
 
const bool stripYResiduals_
 
edm::ESHandle< TrackerGeometrytkGeom_
 
const bool useFit_
 
const bool useFwhm_
 
const bool useOverflowForRMS_
 
std::vector< TH1 * > vDeleteObjects_
 
std::vector< TH1 * > vTrack2DHistos_
 
std::vector< TH1 * > vTrackHistos_
 
std::vector< TH1 * > vTrackProfiles_
 

Additional Inherited Members

- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 78 of file TrackerOfflineValidation.cc.

Member Enumeration Documentation

Constructor & Destructor Documentation

TrackerOfflineValidation::TrackerOfflineValidation ( const edm::ParameterSet iConfig)
explicit

Definition at line 355 of file TrackerOfflineValidation.cc.

356  : parSet_(iConfig), bareTkGeomPtr_(0), lCoorHistOn_(parSet_.getParameter<bool>("localCoorHistosOn")),
357  moduleLevelHistsTransient_(parSet_.getParameter<bool>("moduleLevelHistsTransient")),
358  moduleLevelProfiles_(parSet_.getParameter<bool>("moduleLevelProfiles")),
359  stripYResiduals_(parSet_.getParameter<bool>("stripYResiduals")),
360  useFwhm_(parSet_.getParameter<bool>("useFwhm")),
361  useFit_(parSet_.getParameter<bool>("useFit")),
362  useOverflowForRMS_(parSet_.getParameter<bool>("useOverflowForRMS")),
363  dqmMode_(parSet_.getParameter<bool>("useInDqmMode")),
364  moduleDirectory_(parSet_.getParameter<std::string>("moduleDirectoryInOutput")),
365  lastSetup_(nullptr),
366  avalidator_(iConfig, consumesCollector())
367 {
368 }
const TrackerGeometry * bareTkGeomPtr_
T getParameter(std::string const &) const
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
TrackerValidationVariables avalidator_
const edm::ParameterSet parSet_
const edm::EventSetup * lastSetup_
TrackerOfflineValidation::~TrackerOfflineValidation ( )

Definition at line 371 of file TrackerOfflineValidation.cc.

References vDeleteObjects_.

372 {
373  // do anything here that needs to be done at desctruction time
374  // (e.g. close files, deallocate resources etc.)
375  for( std::vector<TH1*>::const_iterator it = vDeleteObjects_.begin(), itEnd = vDeleteObjects_.end();
376  it != itEnd;
377  ++it) delete *it;
378 }
std::vector< TH1 * > vDeleteObjects_

Member Function Documentation

void TrackerOfflineValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 973 of file TrackerOfflineValidation.cc.

References avalidator_, checkBookHists(), cond::rpcobgas::detid, TrackerValidationVariables::fillTrackQuantities(), getHistStructFromMap(), GetIndex(), isPixel(), lCoorHistOn_, TrackerOfflineValidation::ModuleHistos::LocalX, TrackerOfflineValidation::ModuleHistos::LocalY, moduleLevelProfiles_, TrackerOfflineValidation::ModuleHistos::NormResHisto, TrackerOfflineValidation::ModuleHistos::NormResXprimeHisto, TrackerOfflineValidation::ModuleHistos::NormResYprimeHisto, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, TrackerOfflineValidation::ModuleHistos::ResHisto, TrackerOfflineValidation::ModuleHistos::ResXprimeHisto, TrackerOfflineValidation::ModuleHistos::ResXvsXProfile, TrackerOfflineValidation::ModuleHistos::ResXvsYProfile, TrackerOfflineValidation::ModuleHistos::ResYHisto, TrackerOfflineValidation::ModuleHistos::ResYprimeHisto, TrackerOfflineValidation::ModuleHistos::ResYvsXProfile, TrackerOfflineValidation::ModuleHistos::ResYvsYProfile, stripYResiduals_, DetId::subdetId(), funct::tan(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, useOverflowForRMS_, vTrack2DHistos_, vTrackHistos_, and vTrackProfiles_.

974 {
975  if (useOverflowForRMS_)TH1::StatOverflows(kTRUE);
976  this->checkBookHists(iSetup); // check whether hists are booked and do so if not yet done
977 
978  std::vector<TrackerValidationVariables::AVTrackStruct> vTrackstruct;
979  avalidator_.fillTrackQuantities(iEvent, iSetup, vTrackstruct);
980 
981  for (std::vector<TrackerValidationVariables::AVTrackStruct>::const_iterator itT = vTrackstruct.begin();
982  itT != vTrackstruct.end();
983  ++itT) {
984 
985  // Fill 1D track histos
986  static const int etaindex = this->GetIndex(vTrackHistos_,"h_tracketa");
987  vTrackHistos_[etaindex]->Fill(itT->eta);
988  static const int phiindex = this->GetIndex(vTrackHistos_,"h_trackphi");
989  vTrackHistos_[phiindex]->Fill(itT->phi);
990  static const int numOfValidHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfValidHits");
991  vTrackHistos_[numOfValidHitsindex]->Fill(itT->numberOfValidHits);
992  static const int numOfLostHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfLostHits");
993  vTrackHistos_[numOfLostHitsindex]->Fill(itT->numberOfLostHits);
994  static const int kappaindex = this->GetIndex(vTrackHistos_,"h_curvature");
995  vTrackHistos_[kappaindex]->Fill(itT->kappa);
996  static const int kappaposindex = this->GetIndex(vTrackHistos_,"h_curvature_pos");
997  if (itT->charge > 0)
998  vTrackHistos_[kappaposindex]->Fill(fabs(itT->kappa));
999  static const int kappanegindex = this->GetIndex(vTrackHistos_,"h_curvature_neg");
1000  if (itT->charge < 0)
1001  vTrackHistos_[kappanegindex]->Fill(fabs(itT->kappa));
1002  static const int normchi2index = this->GetIndex(vTrackHistos_,"h_normchi2");
1003  vTrackHistos_[normchi2index]->Fill(itT->normchi2);
1004  static const int chi2index = this->GetIndex(vTrackHistos_,"h_chi2");
1005  vTrackHistos_[chi2index]->Fill(itT->chi2);
1006  static const int chi2Probindex = this->GetIndex(vTrackHistos_,"h_chi2Prob");
1007  vTrackHistos_[chi2Probindex]->Fill(itT->chi2Prob);
1008  static const int ptindex = this->GetIndex(vTrackHistos_,"h_pt");
1009  vTrackHistos_[ptindex]->Fill(itT->pt);
1010  if (itT->ptError != 0.) {
1011  static const int ptResolutionindex = this->GetIndex(vTrackHistos_,"h_ptResolution");
1012  vTrackHistos_[ptResolutionindex]->Fill(itT->ptError/itT->pt);
1013  }
1014  // Fill track profiles
1015  static const int d0phiindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_phi");
1016  vTrackProfiles_[d0phiindex]->Fill(itT->phi,itT->d0);
1017  static const int dzphiindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_phi");
1018  vTrackProfiles_[dzphiindex]->Fill(itT->phi,itT->dz);
1019  static const int d0etaindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_eta");
1020  vTrackProfiles_[d0etaindex]->Fill(itT->eta,itT->d0);
1021  static const int dzetaindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_eta");
1022  vTrackProfiles_[dzetaindex]->Fill(itT->eta,itT->dz);
1023  static const int chiphiindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_phi");
1024  vTrackProfiles_[chiphiindex]->Fill(itT->phi,itT->chi2);
1025  static const int chiProbphiindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_phi");
1026  vTrackProfiles_[chiProbphiindex]->Fill(itT->phi,itT->chi2Prob);
1027  static const int chiProbabsd0index = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_d0");
1028  vTrackProfiles_[chiProbabsd0index]->Fill(fabs(itT->d0),itT->chi2Prob);
1029  static const int normchiphiindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_phi");
1030  vTrackProfiles_[normchiphiindex]->Fill(itT->phi,itT->normchi2);
1031  static const int chietaindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_eta");
1032  vTrackProfiles_[chietaindex]->Fill(itT->eta,itT->chi2);
1033  static const int normchiptindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_pt");
1034  vTrackProfiles_[normchiptindex]->Fill(itT->pt,itT->normchi2);
1035  static const int normchipindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_p");
1036  vTrackProfiles_[normchipindex]->Fill(itT->p,itT->normchi2);
1037  static const int chiProbetaindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_eta");
1038  vTrackProfiles_[chiProbetaindex]->Fill(itT->eta,itT->chi2Prob);
1039  static const int normchietaindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_eta");
1040  vTrackProfiles_[normchietaindex]->Fill(itT->eta,itT->normchi2);
1041  static const int kappaphiindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_phi");
1042  vTrackProfiles_[kappaphiindex]->Fill(itT->phi,itT->kappa);
1043  static const int kappaetaindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_eta");
1044  vTrackProfiles_[kappaetaindex]->Fill(itT->eta,itT->kappa);
1045  static const int ptResphiindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_phi");
1046  vTrackProfiles_[ptResphiindex]->Fill(itT->phi,itT->ptError/itT->pt);
1047  static const int ptResetaindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_eta");
1048  vTrackProfiles_[ptResetaindex]->Fill(itT->eta,itT->ptError/itT->pt);
1049 
1050  // Fill 2D track histos
1051  static const int d0phiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_phi");
1052  vTrack2DHistos_[d0phiindex_2d]->Fill(itT->phi,itT->d0);
1053  static const int dzphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_phi");
1054  vTrack2DHistos_[dzphiindex_2d]->Fill(itT->phi,itT->dz);
1055  static const int d0etaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_eta");
1056  vTrack2DHistos_[d0etaindex_2d]->Fill(itT->eta,itT->d0);
1057  static const int dzetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_eta");
1058  vTrack2DHistos_[dzetaindex_2d]->Fill(itT->eta,itT->dz);
1059  static const int chiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_phi");
1060  vTrack2DHistos_[chiphiindex_2d]->Fill(itT->phi,itT->chi2);
1061  static const int chiProbphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_phi");
1062  vTrack2DHistos_[chiProbphiindex_2d]->Fill(itT->phi,itT->chi2Prob);
1063  static const int chiProbabsd0index_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_d0");
1064  vTrack2DHistos_[chiProbabsd0index_2d]->Fill(fabs(itT->d0),itT->chi2Prob);
1065  static const int normchiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_phi");
1066  vTrack2DHistos_[normchiphiindex_2d]->Fill(itT->phi,itT->normchi2);
1067  static const int chietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_eta");
1068  vTrack2DHistos_[chietaindex_2d]->Fill(itT->eta,itT->chi2);
1069  static const int chiProbetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_eta");
1070  vTrack2DHistos_[chiProbetaindex_2d]->Fill(itT->eta,itT->chi2Prob);
1071  static const int normchietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_eta");
1072  vTrack2DHistos_[normchietaindex_2d]->Fill(itT->eta,itT->normchi2);
1073  static const int kappaphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_phi");
1074  vTrack2DHistos_[kappaphiindex_2d]->Fill(itT->phi,itT->kappa);
1075  static const int kappaetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_eta");
1076  vTrack2DHistos_[kappaetaindex_2d]->Fill(itT->eta,itT->kappa);
1077  static const int normchi2kappa_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_kappa");
1078  vTrack2DHistos_[normchi2kappa_2d]->Fill(itT->normchi2,itT->kappa);
1079 
1080  // hit quantities: residuals, normalized residuals
1081  for (std::vector<TrackerValidationVariables::AVHitStruct>::const_iterator itH = itT->hits.begin();
1082  itH != itT->hits.end();
1083  ++itH) {
1084 
1085  DetId detid(itH->rawDetId);
1086  ModuleHistos &histStruct = this->getHistStructFromMap(detid);
1087 
1088  // fill histos in local coordinates if set in cf
1089  if (lCoorHistOn_) {
1090  histStruct.ResHisto->Fill(itH->resX);
1091  if(itH->resErrX != 0) histStruct.NormResHisto->Fill(itH->resX/itH->resErrX);
1092  if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
1093  histStruct.ResYHisto->Fill(itH->resY);
1094  // here add un-primed normalised y-residuals if wanted
1095  }
1096  }
1097  if (itH->resXprime != -999.) {
1098  histStruct.ResXprimeHisto->Fill(itH->resXprime);
1099 
1100  /******************************* Fill 2-D histo ResX vs momenta *****************************/
1101  if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
1102  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixB");
1103  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1104  }
1105  if (detid.subdetId() == PixelSubdetector::PixelEndcap) {
1106  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixE");
1107  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1108  }
1109  if (detid.subdetId() == StripSubdetector::TIB) {
1110  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TIB");
1111  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1112  }
1113  if (detid.subdetId() == StripSubdetector::TID) {
1114  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TID");
1115  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1116  }
1117  if (detid.subdetId() == StripSubdetector::TOB) {
1118  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TOB");
1119  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1120  }
1121  if (detid.subdetId() == StripSubdetector::TEC) {
1122  static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TEC");
1123  vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
1124  }
1125  /******************************************/
1126 
1127  if ( moduleLevelProfiles_ && itH->inside ) {
1128  float tgalpha = tan(itH->localAlpha);
1129  if ( fabs(tgalpha)!=0 ){
1130  histStruct.LocalX->Fill(itH->localXnorm, tgalpha*tgalpha);
1131  histStruct.LocalY->Fill(itH->localYnorm, tgalpha*tgalpha);
1132 /* if (this->isEndCap(detid.subdetId()) && !this->isPixel(detid.subdetId())) {
1133  if((itH->resX)*(itH->resXprime)>0){
1134  histStruct.ResXvsXProfile->Fill(itH->localXnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1135  histStruct.ResXvsYProfile->Fill(itH->localYnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1136  } else {
1137  histStruct.ResXvsXProfile->Fill(itH->localXnorm, (-1)*itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1138  histStruct.ResXvsYProfile->Fill(itH->localYnorm, (-1)*itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1139  }
1140 
1141  }else {
1142 */
1143  histStruct.ResXvsXProfile->Fill(itH->localXnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1144  histStruct.ResXvsYProfile->Fill(itH->localYnorm, itH->resXatTrkY/tgalpha, tgalpha*tgalpha);
1145 
1146 // }
1147 
1148  }
1149  }
1150 
1151  if(itH->resXprimeErr != 0 && itH->resXprimeErr != -999 ) {
1152  histStruct.NormResXprimeHisto->Fill(itH->resXprime/itH->resXprimeErr);
1153  }
1154  }
1155 
1156  if (itH->resYprime != -999.) {
1157  if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
1158  histStruct.ResYprimeHisto->Fill(itH->resYprime);
1159 
1160  /******************************* Fill 2-D histo ResY vs momenta *****************************/
1161  if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
1162  static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixB");
1163  vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
1164  }
1165  if (detid.subdetId() == PixelSubdetector::PixelEndcap) {
1166  static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixE");
1167  vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
1168  }
1169  /******************************************/
1170 
1171  if ( moduleLevelProfiles_ && itH->inside ) {
1172  float tgbeta = tan(itH->localBeta);
1173  if ( fabs(tgbeta)!=0 ){
1174 /* if (this->isEndCap(detid.subdetId()) && !this->isPixel(detid.subdetId())) {
1175 
1176  if((itH->resY)*(itH->resYprime)>0){
1177  histStruct.ResYvsXProfile->Fill(itH->localXnorm, itH->resYprime/tgbeta, tgbeta*tgbeta);
1178  histStruct.ResYvsYProfile->Fill(itH->localYnorm, itH->resYprime/tgbeta, tgbeta*tgbeta);
1179  } else {
1180  histStruct.ResYvsXProfile->Fill(itH->localXnorm, (-1)*itH->resYprime/tgbeta, tgbeta*tgbeta);
1181  histStruct.ResYvsYProfile->Fill(itH->localYnorm, (-1)*itH->resYprime/tgbeta, tgbeta*tgbeta);
1182  }
1183 
1184  }else {
1185 */
1186  histStruct.ResYvsXProfile->Fill(itH->localXnorm, itH->resY/tgbeta, tgbeta*tgbeta);
1187  histStruct.ResYvsYProfile->Fill(itH->localYnorm, itH->resY/tgbeta, tgbeta*tgbeta);
1188 // }
1189  }
1190  }
1191 
1192  if (itH->resYprimeErr != 0 && itH->resYprimeErr != -999. ) {
1193  histStruct.NormResYprimeHisto->Fill(itH->resYprime/itH->resYprimeErr);
1194  }
1195  }
1196  }
1197 
1198  } // finish loop over hit quantities
1199  } // finish loop over track quantities
1200 
1201  if (useOverflowForRMS_) TH1::StatOverflows(kFALSE);
1202 }
bool isPixel(uint32_t subDetId)
std::vector< TH1 * > vTrack2DHistos_
int GetIndex(const std::vector< OBJECT_TYPE * > &vec, const TString &name)
ModuleHistos & getHistStructFromMap(const DetId &detid)
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
virtual void checkBookHists(const edm::EventSetup &setup)
Definition: DetId.h:18
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
TrackerValidationVariables avalidator_
std::vector< TH1 * > vTrackProfiles_
void TrackerOfflineValidation::bookDirHists ( DirectoryWrapper tfd,
const Alignable ali,
const TrackerTopology tTopo 
)
private

Definition at line 599 of file TrackerOfflineValidation.cc.

References Alignable::alignableObjectId(), bookHists(), makeMuonMisalignmentScenario::components, Alignable::components(), HILowLumiHLTOfflineSource_cfi::dirname, dqmMode_, f, i, AlignableObjectId::idToString(), isDetOrDetUnit(), LogDebug, moduleDirectory_, findQualityFiles::size, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by checkBookHists().

600 {
601  std::vector<Alignable*> alivec(ali.components());
602  for(int i=0, iEnd = ali.components().size();i < iEnd; ++i) {
603  std::string structurename = AlignableObjectId::idToString((alivec)[i]->alignableObjectId());
604  LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
605  std::stringstream dirname;
606  dirname << structurename;
607  // add no suffix counter to Strip and Pixel, just aesthetics
608  if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << i+1;
609 
610  if (structurename.find("Endcap",0) != std::string::npos ) {
611  DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
612  bookHists(f, *(alivec)[i], tTopo, ali.alignableObjectId() , i);
613  bookDirHists( f, *(alivec)[i], tTopo);
614  } else if( !(this->isDetOrDetUnit( (alivec)[i]->alignableObjectId()) )
615  || alivec[i]->components().size() > 1) {
616  DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
617  bookHists(tfd, *(alivec)[i], tTopo, ali.alignableObjectId() , i);
618  bookDirHists( f, *(alivec)[i], tTopo);
619  } else {
620  bookHists(tfd, *(alivec)[i], tTopo, ali.alignableObjectId() , i);
621  }
622  }
623 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
virtual Alignables components() const =0
Return vector of all direct components.
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
double f[11][100]
bool isDetOrDetUnit(align::StructureType type)
void bookDirHists(DirectoryWrapper &tfd, const Alignable &ali, const TrackerTopology *tTopo)
static const char * idToString(align::StructureType type)
tuple size
Write out results.
void bookHists(DirectoryWrapper &tfd, const Alignable &ali, const TrackerTopology *tTopo, align::StructureType type, int i)
void TrackerOfflineValidation::bookGlobalHists ( DirectoryWrapper tfd)
private

Definition at line 427 of file TrackerOfflineValidation.cc.

References TrackerOfflineValidation::DirectoryWrapper::make(), vTrack2DHistos_, vTrackHistos_, and vTrackProfiles_.

Referenced by checkBookHists().

428 {
429 
430  vTrackHistos_.push_back(tfd.make<TH1F>("h_tracketa",
431  "Track #eta;#eta_{Track};Number of Tracks",
432  90,-3.,3.));
433  vTrackHistos_.push_back(tfd.make<TH1F>("h_trackphi",
434  "Track #phi;#phi_{Track};Number of Tracks",
435  90,-3.15,3.15));
436  vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfValidHits",
437  "Track # of valid hits;# of valid hits _{Track};Number of Tracks",
438  40,0.,40.));
439  vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfLostHits",
440  "Track # of lost hits;# of lost hits _{Track};Number of Tracks",
441  10,0.,10.));
442  vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature",
443  "Curvature #kappa;#kappa_{Track};Number of Tracks",
444  100,-.05,.05));
445  vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_pos",
446  "Curvature |#kappa| Positive Tracks;|#kappa_{pos Track}|;Number of Tracks",
447  100,.0,.05));
448  vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_neg",
449  "Curvature |#kappa| Negative Tracks;|#kappa_{neg Track}|;Number of Tracks",
450  100,.0,.05));
451  vTrackHistos_.push_back(tfd.make<TH1F>("h_diff_curvature",
452  "Curvature |#kappa| Tracks Difference;|#kappa_{Track}|;# Pos Tracks - # Neg Tracks",
453  100,.0,.05));
454  vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2",
455  "#chi^{2};#chi^{2}_{Track};Number of Tracks",
456  500,-0.01,500.));
457  vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2Prob",
458  "#chi^{2} probability;#chi^{2}prob_{Track};Number of Tracks",
459  100,0.0,1.));
460  vTrackHistos_.push_back(tfd.make<TH1F>("h_normchi2",
461  "#chi^{2}/ndof;#chi^{2}/ndof;Number of Tracks",
462  100,-0.01,10.));
463  vTrackHistos_.push_back(tfd.make<TH1F>("h_pt",
464  "p_{T}^{track};p_{T}^{track} [GeV];Number of Tracks",
465  250,0.,250));
466  vTrackHistos_.push_back(tfd.make<TH1F>("h_ptResolution",
467  "#delta_{p_{T}}/p_{T}^{track};#delta_{p_{T}}/p_{T}^{track};Number of Tracks",
468  100,0.,0.5));
469  vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_phi",
470  "Transverse Impact Parameter vs. #phi;#phi_{Track};#LT d_{0} #GT [cm]",
471  100,-3.15,3.15));
472  vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_phi",
473  "Longitudinal Impact Parameter vs. #phi;#phi_{Track};#LT d_{z} #GT [cm]",
474  100,-3.15,3.15));
475  vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_eta",
476  "Transverse Impact Parameter vs. #eta;#eta_{Track};#LT d_{0} #GT [cm]",
477  100,-3.15,3.15));
478  vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_eta",
479  "Longitudinal Impact Parameter vs. #eta;#eta_{Track};#LT d_{z} #GT [cm]",
480  100,-3.15,3.15));
481  vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_phi",
482  "#chi^{2} vs. #phi;#phi_{Track};#LT #chi^{2} #GT",
483  100,-3.15,3.15));
484  vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_phi",
485  "#chi^{2} probablility vs. #phi;#phi_{Track};#LT #chi^{2} probability#GT",
486  100,-3.15,3.15));
487  vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_d0",
488  "#chi^{2} probablility vs. |d_{0}|;|d_{0}|[cm];#LT #chi^{2} probability#GT",
489  100,0,80));
490  vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_phi",
491  "#chi^{2}/ndof vs. #phi;#phi_{Track};#LT #chi^{2}/ndof #GT",
492  100,-3.15,3.15));
493  vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_eta",
494  "#chi^{2} vs. #eta;#eta_{Track};#LT #chi^{2} #GT",
495  100,-3.15,3.15));
496  //variable binning for chi2/ndof vs. pT
497  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.};
498  vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_pt",
499  "norm #chi^{2} vs. p_{T}_{Track}; p_{T}_{Track};#LT #chi^{2}/ndof #GT",
500  18,xBins));
501 
502  vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_p",
503  "#chi^{2}/ndof vs. p_{Track};p_{Track};#LT #chi^{2}/ndof #GT",
504  18,xBins));
505  vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_eta",
506  "#chi^{2} probability vs. #eta;#eta_{Track};#LT #chi^{2} probability #GT",
507  100,-3.15,3.15));
508  vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_eta",
509  "#chi^{2}/ndof vs. #eta;#eta_{Track};#LT #chi^{2}/ndof #GT",
510  100,-3.15,3.15));
511  vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_phi",
512  "#kappa vs. #phi;#phi_{Track};#kappa",
513  100,-3.15,3.15));
514  vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_eta",
515  "#kappa vs. #eta;#eta_{Track};#kappa",
516  100,-3.15,3.15));
517  vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_phi",
518  "#delta_{p_{T}}/p_{T}^{track};#phi^{track};#delta_{p_{T}}/p_{T}^{track}",
519  100, -3.15,3.15));
520  vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_eta",
521  "#delta_{p_{T}}/p_{T}^{track};#eta^{track};#delta_{p_{T}}/p_{T}^{track}",
522  100, -3.15,3.15));
523 
524  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_phi",
525  "Transverse Impact Parameter vs. #phi;#phi_{Track};d_{0} [cm]",
526  100, -3.15, 3.15, 100,-1.,1.) );
527  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_phi",
528  "Longitudinal Impact Parameter vs. #phi;#phi_{Track};d_{z} [cm]",
529  100, -3.15, 3.15, 100,-100.,100.));
530  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_eta",
531  "Transverse Impact Parameter vs. #eta;#eta_{Track};d_{0} [cm]",
532  100, -3.15, 3.15, 100,-1.,1.));
533  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_eta",
534  "Longitudinal Impact Parameter vs. #eta;#eta_{Track};d_{z} [cm]",
535  100, -3.15, 3.15, 100,-100.,100.));
536  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_phi",
537  "#chi^{2} vs. #phi;#phi_{Track};#chi^{2}",
538  100, -3.15, 3.15, 500, 0., 500.));
539  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_phi",
540  "#chi^{2} probability vs. #phi;#phi_{Track};#chi^{2} probability",
541  100, -3.15, 3.15, 100, 0., 1.));
542  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_d0",
543  "#chi^{2} probability vs. |d_{0}|;|d_{0}| [cm];#chi^{2} probability",
544  100, 0, 80, 100, 0., 1.));
545  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_phi",
546  "#chi^{2}/ndof vs. #phi;#phi_{Track};#chi^{2}/ndof",
547  100, -3.15, 3.15, 100, 0., 10.));
548  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_eta",
549  "#chi^{2} vs. #eta;#eta_{Track};#chi^{2}",
550  100, -3.15, 3.15, 500, 0., 500.));
551  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_eta",
552  "#chi^{2} probaility vs. #eta;#eta_{Track};#chi^{2} probability",
553  100, -3.15, 3.15, 100, 0., 1.));
554  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_eta",
555  "#chi^{2}/ndof vs. #eta;#eta_{Track};#chi^{2}/ndof",
556  100,-3.15,3.15, 100, 0., 10.));
557  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_phi",
558  "#kappa vs. #phi;#phi_{Track};#kappa",
559  100,-3.15,3.15, 100, .0,.05));
560  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_eta",
561  "#kappa vs. #eta;#eta_{Track};#kappa",
562  100,-3.15,3.15, 100, .0,.05));
563  vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_kappa",
564  "#kappa vs. #chi^{2}/ndof;#chi^{2}/ndof;#kappa",
565  100,0.,10, 100,-.03,.03));
566 
567  /****************** Definition of 2-D Histos of ResX vs momenta ****************************/
568  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixB",
569  "res_{x'} vs momentum in BPix;p [GeV]; res_{x'}",
570  15,0.,15., 200, -0.1,0.1));
571  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixE",
572  "res_{x'} vs momentum in FPix;p [GeV]; res_{x'}",
573  15,0.,15., 200, -0.1,0.1));
574  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TIB",
575  "res_{x'} vs momentum in TIB;p [GeV]; res_{x'}",
576  15,0.,15., 200, -0.1,0.1));
577  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TID",
578  "res_{x'} vs momentum in TID;p [GeV]; res_{x'}",
579  15,0.,15., 200, -0.1,0.1));
580  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TOB",
581  "res_{x'} vs momentum in TOB;p [GeV]; res_{x'}",
582  15,0.,15., 200, -0.1,0.1));
583  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TEC",
584  "res_{x'} vs momentum in TEC;p [GeV]; res_{x'}",
585  15,0.,15., 200, -0.1,0.1));
586 
587  /****************** Definition of 2-D Histos of ResY vs momenta ****************************/
588  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixB",
589  "res_{y'} vs momentum in BPix;p [GeV]; res_{y'}",
590  15,0.,15., 200, -0.1,0.1));
591  vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixE",
592  "res_{y'} vs momentum in FPix;p [GeV]; res_{y'}",
593  15,0.,15., 200, -0.1,0.1));
594 
595 }
std::vector< TH1 * > vTrack2DHistos_
std::vector< TH1 * > vTrackProfiles_
void TrackerOfflineValidation::bookHists ( DirectoryWrapper tfd,
const Alignable ali,
const TrackerTopology tTopo,
align::StructureType  type,
int  i 
)
private

Definition at line 627 of file TrackerOfflineValidation.cc.

References align::AlignableDetUnit, Alignable::alignableObjectId(), bookTH1F(), bookTProfile(), dqmMode_, getBinning(), getHistStructFromMap(), Alignable::id(), align::invalid, isBarrel(), isDetOrDetUnit(), isEndCap(), isPixel(), lCoorHistOn_, TrackerOfflineValidation::ModuleHistos::LocalX, TrackerOfflineValidation::ModuleHistos::LocalY, moduleLevelHistsTransient_, moduleLevelProfiles_, pileupCalc::nbins, TrackerOfflineValidation::ModuleHistos::NormResHisto, TrackerOfflineValidation::ModuleHistos::NormResXprimeHisto, TrackerOfflineValidation::ModuleHistos::NormResYprimeHisto, NormXprimeResidual, NormXResidual, NormYprimeResidual, TrackerOfflineValidation::ModuleHistos::ResHisto, TrackerOfflineValidation::ModuleHistos::ResXprimeHisto, TrackerOfflineValidation::ModuleHistos::ResXvsXProfile, TrackerOfflineValidation::ModuleHistos::ResXvsYProfile, TrackerOfflineValidation::ModuleHistos::ResYHisto, TrackerOfflineValidation::ModuleHistos::ResYprimeHisto, TrackerOfflineValidation::ModuleHistos::ResYvsXProfile, TrackerOfflineValidation::ModuleHistos::ResYvsYProfile, AlCaHLTBitMon_QueryRunRegistry::string, stripYResiduals_, TrackerAlignableId::typeAndLayerFromDetId(), SiStripMonitorClusterAlca_cfi::xmax, SiStripMonitorClusterAlca_cfi::xmin, XprimeResidual, XResidual, XResidualProfile, SiStripMonitorClusterAlca_cfi::ymax, SiStripMonitorClusterAlca_cfi::ymin, YprimeResidual, YResidual, and YResidualProfile.

Referenced by bookDirHists().

628 {
629  TrackerAlignableId aliid;
630  const DetId id = ali.id();
631 
632  // comparing subdetandlayer to subdetIds gives a warning at compile time
633  // -> subdetandlayer could also be pair<uint,uint> but this has to be adapted
634  // in AlignableObjId
635  std::pair<int,int> subdetandlayer = aliid.typeAndLayerFromDetId(id, tTopo);
636 
638 
639  // are we on or just above det, detunit level respectively?
640  if (type == align::AlignableDetUnit )subtype = type;
641  else if( this->isDetOrDetUnit(ali.alignableObjectId()) ) subtype = ali.alignableObjectId();
642 
643  // construct histogramm title and name
644  std::stringstream histoname, histotitle, normhistoname, normhistotitle,
645  yhistoname, yhistotitle,
646  xprimehistoname, xprimehistotitle, normxprimehistoname, normxprimehistotitle,
647  yprimehistoname, yprimehistotitle, normyprimehistoname, normyprimehistotitle,
648  localxname, localxtitle, localyname, localytitle,
649  resxvsxprofilename, resxvsxprofiletitle, resyvsxprofilename, resyvsxprofiletitle,
650  resxvsyprofilename, resxvsyprofiletitle, resyvsyprofilename, resyvsyprofiletitle;
651 
652  std::string wheel_or_layer;
653 
654  if( this->isEndCap(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_wheel_";
655  else if ( this->isBarrel(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_layer_";
656  else edm::LogWarning("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookHists"
657  << "Unknown subdetid: " << subdetandlayer.first;
658 
659  histoname << "h_residuals_subdet_" << subdetandlayer.first
660  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
661  yhistoname << "h_y_residuals_subdet_" << subdetandlayer.first
662  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
663  xprimehistoname << "h_xprime_residuals_subdet_" << subdetandlayer.first
664  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
665  yprimehistoname << "h_yprime_residuals_subdet_" << subdetandlayer.first
666  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
667  normhistoname << "h_normresiduals_subdet_" << subdetandlayer.first
668  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
669  normxprimehistoname << "h_normxprimeresiduals_subdet_" << subdetandlayer.first
670  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
671  normyprimehistoname << "h_normyprimeresiduals_subdet_" << subdetandlayer.first
672  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
673  histotitle << "X Residual for module " << id.rawId() << ";x_{tr} - x_{hit} [cm]";
674  yhistotitle << "Y Residual for module " << id.rawId() << ";y_{tr} - y_{hit} [cm]";
675  normhistotitle << "Normalized Residual for module " << id.rawId() << ";x_{tr} - x_{hit}/#sigma";
676  xprimehistotitle << "X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})' [cm]";
677  normxprimehistotitle << "Normalized X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})'/#sigma";
678  yprimehistotitle << "Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})' [cm]";
679  normyprimehistotitle << "Normalized Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})'/#sigma";
680 
681  if ( moduleLevelProfiles_ ) {
682  localxname << "h_localx_subdet_" << subdetandlayer.first
683  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
684  localyname << "h_localy_subdet_" << subdetandlayer.first
685  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
686  localxtitle << "u local for module " << id.rawId() << "; u_{tr,r}";
687  localytitle << "v local for module " << id.rawId() << "; v_{tr,r}";
688 
689  resxvsxprofilename << "p_residuals_x_vs_x_subdet_" << subdetandlayer.first
690  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
691  resyvsxprofilename << "p_residuals_y_vs_x_subdet_" << subdetandlayer.first
692  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
693  resxvsyprofilename << "p_residuals_x_vs_y_subdet_" << subdetandlayer.first
694  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
695  resyvsyprofilename << "p_residuals_y_vs_y_subdet_" << subdetandlayer.first
696  << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
697  resxvsxprofiletitle << "U Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
698  resyvsxprofiletitle << "V Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(v_{tr} - v_{hit})/tan#beta [cm]";
699  resxvsyprofiletitle << "U Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
700  resyvsyprofiletitle << "V Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(v_{tr} - v_{hit})/tan#beta [cm]";
701  }
702 
703  if( this->isDetOrDetUnit( subtype ) ) {
704  ModuleHistos &histStruct = this->getHistStructFromMap(id);
705  int nbins = 0;
706  double xmin = 0., xmax = 0.;
707  double ymin = -0.1, ymax = 0.1;
708 
709  // do not allow transient hists in DQM mode
710  bool moduleLevelHistsTransient(moduleLevelHistsTransient_);
711  if (dqmMode_) moduleLevelHistsTransient = false;
712 
713  // decide via cfg if hists in local coordinates should be booked
714  if(lCoorHistOn_) {
715  this->getBinning(id.subdetId(), XResidual, nbins, xmin, xmax);
716  histStruct.ResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
717  histoname.str().c_str(),histotitle.str().c_str(),
718  nbins, xmin, xmax);
719  this->getBinning(id.subdetId(), NormXResidual, nbins, xmin, xmax);
720  histStruct.NormResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
721  normhistoname.str().c_str(),normhistotitle.str().c_str(),
722  nbins, xmin, xmax);
723  }
724  this->getBinning(id.subdetId(), XprimeResidual, nbins, xmin, xmax);
725  histStruct.ResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
726  xprimehistoname.str().c_str(),xprimehistotitle.str().c_str(),
727  nbins, xmin, xmax);
728  this->getBinning(id.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
729  histStruct.NormResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
730  normxprimehistoname.str().c_str(),normxprimehistotitle.str().c_str(),
731  nbins, xmin, xmax);
732 
733  if ( moduleLevelProfiles_ ) {
734  this->getBinning(id.subdetId(), XResidualProfile, nbins, xmin, xmax);
735 
736  histStruct.LocalX = this->bookTH1F(moduleLevelHistsTransient, tfd,
737  localxname.str().c_str(),localxtitle.str().c_str(),
738  nbins, xmin, xmax);
739  histStruct.LocalY = this->bookTH1F(moduleLevelHistsTransient, tfd,
740  localyname.str().c_str(),localytitle.str().c_str(),
741  nbins, xmin, xmax);
742  histStruct.ResXvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
743  resxvsxprofilename.str().c_str(),resxvsxprofiletitle.str().c_str(),
744  nbins, xmin, xmax, ymin, ymax);
745  histStruct.ResXvsXProfile->Sumw2(); // to be filled with weights, so uncertainties need sum of square of weights
746  histStruct.ResXvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
747  resxvsyprofilename.str().c_str(),resxvsyprofiletitle.str().c_str(),
748  nbins, xmin, xmax, ymin, ymax);
749  histStruct.ResXvsYProfile->Sumw2(); // to be filled with weights, so uncertainties need sum of square of weights
750  }
751 
752  if( this->isPixel(subdetandlayer.first) || stripYResiduals_ ) {
753  this->getBinning(id.subdetId(), YprimeResidual, nbins, xmin, xmax);
754  histStruct.ResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
755  yprimehistoname.str().c_str(),yprimehistotitle.str().c_str(),
756  nbins, xmin, xmax);
757  if (lCoorHistOn_) { // un-primed y-residual
758  this->getBinning(id.subdetId(), YResidual, nbins, xmin, xmax);
759  histStruct.ResYHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
760  yhistoname.str().c_str(), yhistotitle.str().c_str(),
761  nbins, xmin, xmax);
762  }
763  this->getBinning(id.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
764  histStruct.NormResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
765  normyprimehistoname.str().c_str(),normyprimehistotitle.str().c_str(),
766  nbins, xmin, xmax);
767  // Here we could add un-primed normalised y-residuals if(lCoorHistOn_)...
768  if ( moduleLevelProfiles_ ) {
769  this->getBinning(id.subdetId(), YResidualProfile, nbins, xmin, xmax);
770 
771  histStruct.ResYvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
772  resyvsxprofilename.str().c_str(),resyvsxprofiletitle.str().c_str(),
773  nbins, xmin, xmax, ymin, ymax);
774  histStruct.ResYvsXProfile->Sumw2(); // to be filled with weights, so uncertainties need sum of square of weights
775  histStruct.ResYvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
776  resyvsyprofilename.str().c_str(),resyvsyprofiletitle.str().c_str(),
777  nbins, xmin, xmax, ymin, ymax);
778  histStruct.ResYvsYProfile->Sumw2(); // to be filled with weights, so uncertainties need sum of square of weights
779  }
780  }
781  }
782 }
type
Definition: HCALResponse.h:21
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:180
std::pair< int, int > typeAndLayerFromDetId(const DetId &detId, const TrackerTopology *tTopo) const
bool isPixel(uint32_t subDetId)
TH1 * bookTH1F(bool isTransient, DirectoryWrapper &tfd, const char *histName, const char *histTitle, int nBinsX, double lowX, double highX)
ModuleHistos & getHistStructFromMap(const DetId &detid)
void getBinning(uint32_t subDetId, TrackerOfflineValidation::HistogrammType residualtype, int &nBinsX, double &lowerBoundX, double &upperBoundX)
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
TProfile * bookTProfile(bool isTransient, DirectoryWrapper &tfd, const char *histName, const char *histTitle, int nBinsX, double lowX, double highX)
bool isDetOrDetUnit(align::StructureType type)
Definition: DetId.h:18
bool isEndCap(uint32_t subDetId)
bool isBarrel(uint32_t subDetId)
TrackerOfflineValidation::SummaryContainer TrackerOfflineValidation::bookSummaryHists ( DirectoryWrapper tfd,
const Alignable ali,
align::StructureType  type,
int  i 
)
private

Definition at line 1300 of file TrackerOfflineValidation.cc.

References align::AlignableDet, align::AlignableDetUnit, Alignable::alignableObjectId(), Alignable::components(), cond::rpcobgas::detid, getBinning(), getHistStructFromMap(), Alignable::id(), AlignableObjectId::idToString(), isPixel(), j, relval_steps::k, TrackerOfflineValidation::DirectoryWrapper::make(), pileupCalc::nbins, TrackerOfflineValidation::ModuleHistos::NormResXprimeHisto, TrackerOfflineValidation::ModuleHistos::NormResYprimeHisto, NormXprimeResidual, NormYprimeResidual, TrackerOfflineValidation::ModuleHistos::ResXprimeHisto, TrackerOfflineValidation::ModuleHistos::ResYprimeHisto, stripYResiduals_, DetId::subdetId(), summarizeBinInContainer(), indexGen::title, SiStripMonitorClusterAlca_cfi::xmax, SiStripMonitorClusterAlca_cfi::xmin, XprimeResidual, and YprimeResidual.

Referenced by collateSummaryHists().

1302 {
1303  const uint aliSize = ali.components().size();
1304  const align::StructureType alitype = ali.alignableObjectId();
1305  const align::StructureType subtype = ali.components()[0]->alignableObjectId();
1306  const char *aliTypeName = AlignableObjectId::idToString(alitype); // lifetime of char* OK
1307  const char *aliSubtypeName = AlignableObjectId::idToString(subtype);
1308  const char *typeName = AlignableObjectId::idToString(type);
1309 
1310  const DetId aliDetId = ali.id();
1311  // y residuals only if pixel or specially requested for strip:
1312  const bool bookResidY = this->isPixel(aliDetId.subdetId()) || stripYResiduals_;
1313 
1314  SummaryContainer sumContainer;
1315 
1316  // Book summary hists with one bin per component,
1317  // but special case for Det with two DetUnit that we want to summarize one level up
1318  // (e.g. in TOBRods with 12 bins for 6 stereo and 6 rphi DetUnit.)
1319  // component of ali is not Det or Det with just one components
1320  const uint subcompSize = ali.components()[0]->components().size();
1321  if (subtype != align::AlignableDet || subcompSize == 1) { // Det with 1 comp. should not exist anymore...
1322  const TString title(Form("Summary for substructures in %s %d;%s;",aliTypeName,i,aliSubtypeName));
1323 
1324  sumContainer.summaryXResiduals_ = tfd.make<TH1F>(Form("h_summaryX%s_%d",aliTypeName,i),
1325  title + "#LT #Delta x' #GT",
1326  aliSize, 0.5, aliSize+0.5);
1327  sumContainer.summaryNormXResiduals_ = tfd.make<TH1F>(Form("h_summaryNormX%s_%d",aliTypeName,i),
1328  title + "#LT #Delta x'/#sigma #GT",
1329  aliSize,0.5,aliSize+0.5);
1330 
1331  if (bookResidY) {
1332  sumContainer.summaryYResiduals_ = tfd.make<TH1F>(Form("h_summaryY%s_%d",aliTypeName,i),
1333  title + "#LT #Delta y' #GT",
1334  aliSize, 0.5, aliSize+0.5);
1335  sumContainer.summaryNormYResiduals_ = tfd.make<TH1F>(Form("h_summaryNormY%s_%d",aliTypeName,i),
1336  title + "#LT #Delta y'/#sigma #GT",
1337  aliSize,0.5,aliSize+0.5);
1338  }
1339 
1340  } else if (subtype == align::AlignableDet && subcompSize > 1) { // fixed: was aliSize before
1341  if (subcompSize != 2) { // strange... expect only 2 DetUnits in DS layers
1342  // this 2 is hardcoded factor 2 in binning below and also assumed later on
1343  edm::LogError("Alignment") << "@SUB=bookSummaryHists"
1344  << "Det with " << subcompSize << " components";
1345  }
1346  // title contains x-title
1347  const TString title(Form("Summary for substructures in %s %d;%s;", aliTypeName, i,
1348  AlignableObjectId::idToString(ali.components()[0]->components()[0]->alignableObjectId())));
1349 
1350  sumContainer.summaryXResiduals_
1351  = tfd.make<TH1F>(Form("h_summaryX%s_%d", aliTypeName, i),
1352  title + "#LT #Delta x' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
1353  sumContainer.summaryNormXResiduals_
1354  = tfd.make<TH1F>(Form("h_summaryNormX%s_%d", aliTypeName, i),
1355  title + "#LT #Delta x'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
1356 
1357  if (bookResidY) {
1358  sumContainer.summaryYResiduals_
1359  = tfd.make<TH1F>(Form("h_summaryY%s_%d", aliTypeName, i),
1360  title + "#LT #Delta y' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
1361  sumContainer.summaryNormYResiduals_
1362  = tfd.make<TH1F>(Form("h_summaryNormY%s_%d", aliTypeName, i),
1363  title + "#LT #Delta y'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
1364  }
1365 
1366  } else {
1367  edm::LogError("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookSummaryHists"
1368  << "No summary histogramm for hierarchy level "
1369  << aliTypeName << " in subdet " << aliDetId.subdetId();
1370  }
1371 
1372  // Now book hists that just sum up the residual histograms from lower levels.
1373  // Axis title is copied from lowest level module of structure.
1374  // Should be safe that y-hists are only touched if non-null pointers...
1375  int nbins = 0;
1376  double xmin = 0., xmax = 0.;
1377  const TString sumTitle(Form("Residual for %s %d in %s;", aliTypeName, i, typeName));
1378  const ModuleHistos &xTitHists = this->getHistStructFromMap(aliDetId); // for x-axis titles
1379  this->getBinning(aliDetId.subdetId(), XprimeResidual, nbins, xmin, xmax);
1380 
1381  sumContainer.sumXResiduals_ = tfd.make<TH1F>(Form("h_Xprime_%s_%d", aliTypeName, i),
1382  sumTitle + xTitHists.ResXprimeHisto->GetXaxis()->GetTitle(),
1383  nbins, xmin, xmax);
1384 
1385  this->getBinning(aliDetId.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
1386  sumContainer.sumNormXResiduals_ = tfd.make<TH1F>(Form("h_NormXprime_%s_%d",aliTypeName,i),
1387  sumTitle + xTitHists.NormResXprimeHisto->GetXaxis()->GetTitle(),
1388  nbins, xmin, xmax);
1389  if (bookResidY) {
1390  this->getBinning(aliDetId.subdetId(), YprimeResidual, nbins, xmin, xmax);
1391  sumContainer.sumYResiduals_ = tfd.make<TH1F>(Form("h_Yprime_%s_%d",aliTypeName,i),
1392  sumTitle + xTitHists.ResYprimeHisto->GetXaxis()->GetTitle(),
1393  nbins, xmin, xmax);
1394 
1395  this->getBinning(aliDetId.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
1396  sumContainer.sumNormYResiduals_ = tfd.make<TH1F>(Form("h_NormYprime_%s_%d",aliTypeName,i),
1397  sumTitle + xTitHists.NormResYprimeHisto->GetXaxis()->GetTitle(),
1398  nbins, xmin, xmax);
1399  }
1400 
1401  // If we are at the lowest level, we already sum up and fill the summary.
1402 
1403  // special case I: For DetUnits and Dets with only one subcomponent start filling summary histos
1404  if( ( subtype == align::AlignableDet && subcompSize == 1) || subtype == align::AlignableDetUnit ) {
1405  for(uint k = 0; k < aliSize; ++k) {
1406  DetId detid = ali.components()[k]->id();
1407  ModuleHistos &histStruct = this->getHistStructFromMap(detid);
1408  this->summarizeBinInContainer(k+1, detid.subdetId() ,sumContainer, histStruct );
1409  sumContainer.sumXResiduals_->Add(histStruct.ResXprimeHisto);
1410  sumContainer.sumNormXResiduals_->Add(histStruct.NormResXprimeHisto);
1411  if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
1412  sumContainer.sumYResiduals_->Add(histStruct.ResYprimeHisto);
1413  sumContainer.sumNormYResiduals_->Add(histStruct.NormResYprimeHisto);
1414  }
1415  }
1416  } else if( subtype == align::AlignableDet && subcompSize > 1) { // fixed: was aliSize before
1417  // special case II: Fill summary histos for dets with two detunits
1418  for(uint k = 0; k < aliSize; ++k) {
1419  for(uint j = 0; j < subcompSize; ++j) { // assumes all have same size (as binning does)
1420  DetId detid = ali.components()[k]->components()[j]->id();
1421  ModuleHistos &histStruct = this->getHistStructFromMap(detid);
1422  this->summarizeBinInContainer(2*k+j+1, detid.subdetId() ,sumContainer, histStruct );
1423  sumContainer.sumXResiduals_->Add( histStruct.ResXprimeHisto);
1424  sumContainer.sumNormXResiduals_->Add( histStruct.NormResXprimeHisto);
1425  if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
1426  sumContainer.sumYResiduals_->Add( histStruct.ResYprimeHisto);
1427  sumContainer.sumNormYResiduals_->Add( histStruct.NormResYprimeHisto);
1428  }
1429  }
1430  }
1431  }
1432  return sumContainer;
1433 }
type
Definition: HCALResponse.h:21
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:180
int i
Definition: DBlmapReader.cc:9
bool isPixel(uint32_t subDetId)
ModuleHistos & getHistStructFromMap(const DetId &detid)
virtual Alignables components() const =0
Return vector of all direct components.
void getBinning(uint32_t subDetId, TrackerOfflineValidation::HistogrammType residualtype, int &nBinsX, double &lowerBoundX, double &upperBoundX)
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
int j
Definition: DBlmapReader.cc:9
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
Definition: DetId.h:18
void summarizeBinInContainer(int bin, SummaryContainer &targetContainer, SummaryContainer &sourceContainer)
static const char * idToString(align::StructureType type)
TH1 * TrackerOfflineValidation::bookTH1F ( bool  isTransient,
DirectoryWrapper tfd,
const char *  histName,
const char *  histTitle,
int  nBinsX,
double  lowX,
double  highX 
)
private

Definition at line 785 of file TrackerOfflineValidation.cc.

References TrackerOfflineValidation::DirectoryWrapper::make(), and vDeleteObjects_.

Referenced by bookHists().

787 {
788  if (isTransient) {
789  vDeleteObjects_.push_back(new TH1F(histName, histTitle, nBinsX, lowX, highX));
790  return vDeleteObjects_.back(); // return last element of vector
791  }
792  else
793  return tfd.make<TH1F>(histName, histTitle, nBinsX, lowX, highX);
794 }
std::vector< TH1 * > vDeleteObjects_
TProfile * TrackerOfflineValidation::bookTProfile ( bool  isTransient,
DirectoryWrapper tfd,
const char *  histName,
const char *  histTitle,
int  nBinsX,
double  lowX,
double  highX 
)
private

Definition at line 796 of file TrackerOfflineValidation.cc.

References TrackerOfflineValidation::DirectoryWrapper::make(), and vDeleteObjects_.

Referenced by bookHists().

798 {
799  if (isTransient) {
800  TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX);
801  vDeleteObjects_.push_back(profile);
802  return profile;
803  }
804  else
805  return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX);
806 }
std::vector< TH1 * > vDeleteObjects_
TProfile * TrackerOfflineValidation::bookTProfile ( bool  isTransient,
DirectoryWrapper tfd,
const char *  histName,
const char *  histTitle,
int  nBinsX,
double  lowX,
double  highX,
double  lowY,
double  highY 
)
private

Definition at line 809 of file TrackerOfflineValidation.cc.

References TrackerOfflineValidation::DirectoryWrapper::make(), and vDeleteObjects_.

811 {
812  if (isTransient) {
813  TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
814  vDeleteObjects_.push_back(profile);
815  return profile;
816  }
817  else
818  return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
819 }
std::vector< TH1 * > vDeleteObjects_
void TrackerOfflineValidation::checkBookHists ( const edm::EventSetup setup)
privatevirtual

Definition at line 388 of file TrackerOfflineValidation.cc.

References bareTkGeomPtr_, bookDirHists(), bookGlobalHists(), TrackerGeometry::detIds(), TrackerGeometry::detUnitIds(), dqmMode_, edm::EventSetup::get(), moduleDirectory_, edm::ESHandle< class >::product(), AlCaHLTBitMon_QueryRunRegistry::string, and tkGeom_.

Referenced by analyze().

389 {
390  es.get<TrackerDigiGeometryRecord>().get( tkGeom_ );
391  const TrackerGeometry *newBareTkGeomPtr = &(*tkGeom_);
392  if (newBareTkGeomPtr == bareTkGeomPtr_) return; // already booked hists, nothing changed
393 
394  if (!bareTkGeomPtr_) { // pointer not yet set: called the first time => book hists
395 
396  //Retrieve tracker topology from geometry
397  edm::ESHandle<TrackerTopology> tTopoHandle;
398  es.get<TrackerTopologyRcd>().get(tTopoHandle);
399  const TrackerTopology* const tTopo = tTopoHandle.product();
400 
401  // construct alignable tracker to get access to alignable hierarchy
402  AlignableTracker aliTracker(&(*tkGeom_), tTopo);
403 
404  edm::LogInfo("TrackerOfflineValidation") << "There are " << newBareTkGeomPtr->detIds().size()
405  << " dets in the Geometry record.\n"
406  << "Out of these "<<newBareTkGeomPtr->detUnitIds().size()
407  <<" are detUnits";
408 
409  // Book Histogramms for global track quantities
410  std::string globDir("GlobalTrackVariables");
411  DirectoryWrapper trackglobal(globDir,moduleDirectory_,dqmMode_);
412  this->bookGlobalHists(trackglobal);
413 
414  // recursively book histogramms on lowest level
415  DirectoryWrapper tfdw("",moduleDirectory_,dqmMode_);
416  this->bookDirHists(tfdw, aliTracker, tTopo);
417  }
418  else { // histograms booked, but changed TrackerGeometry?
419  edm::LogWarning("GeometryChange") << "@SUB=checkBookHists"
420  << "TrackerGeometry changed, but will not re-book hists!";
421  }
422  bareTkGeomPtr_ = newBareTkGeomPtr;
423 }
const TrackerGeometry * bareTkGeomPtr_
edm::ESHandle< TrackerGeometry > tkGeom_
virtual const DetIdContainer & detIds() const
Returm a vector of all GeomDet DetIds (including those of GeomDetUnits)
T const * product() const
Definition: ESHandle.h:86
void bookGlobalHists(DirectoryWrapper &tfd)
virtual const DetIdContainer & detUnitIds() const
Returm a vector of all GeomDetUnit DetIds.
void bookDirHists(DirectoryWrapper &tfd, const Alignable &ali, const TrackerTopology *tTopo)
void TrackerOfflineValidation::collateSummaryHists ( DirectoryWrapper tfd,
const Alignable ali,
int  i,
std::vector< TrackerOfflineValidation::SummaryContainer > &  vLevelProfiles 
)
private

Definition at line 1254 of file TrackerOfflineValidation.cc.

References Alignable::alignableObjectId(), bookSummaryHists(), makeMuonMisalignmentScenario::components, Alignable::components(), HILowLumiHLTOfflineSource_cfi::dirname, dqmMode_, f, fitResiduals(), AlignableObjectId::idToString(), isDetOrDetUnit(), LogDebug, moduleDirectory_, gen::n, findQualityFiles::size, AlCaHLTBitMon_QueryRunRegistry::string, and summarizeBinInContainer().

Referenced by endJob().

1256 {
1257  std::vector<Alignable*> alivec(ali.components());
1258  if( this->isDetOrDetUnit((alivec)[0]->alignableObjectId()) ) return;
1259 
1260  for(int iComp=0, iCompEnd = ali.components().size();iComp < iCompEnd; ++iComp) {
1261  std::vector< TrackerOfflineValidation::SummaryContainer > vProfiles;
1262  std::string structurename = AlignableObjectId::idToString((alivec)[iComp]->alignableObjectId());
1263 
1264  LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
1265  std::stringstream dirname;
1266  dirname << structurename;
1267 
1268  // add no suffix counter to strip and pixel -> just aesthetics
1269  if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << iComp+1;
1270 
1271  if( !(this->isDetOrDetUnit( (alivec)[iComp]->alignableObjectId()) )
1272  || (alivec)[0]->components().size() > 1 ) {
1273  DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
1274  this->collateSummaryHists( f, *(alivec)[iComp], i, vProfiles);
1275  vLevelProfiles.push_back(this->bookSummaryHists(tfd, *(alivec[iComp]), ali.alignableObjectId(), iComp+1));
1276  TH1 *hY = vLevelProfiles[iComp].sumYResiduals_;
1277  TH1 *hNormY = vLevelProfiles[iComp].sumNormYResiduals_;
1278  for(uint n = 0; n < vProfiles.size(); ++n) {
1279  this->summarizeBinInContainer(n+1, vLevelProfiles[iComp], vProfiles[n] );
1280  vLevelProfiles[iComp].sumXResiduals_->Add(vProfiles[n].sumXResiduals_);
1281  vLevelProfiles[iComp].sumNormXResiduals_->Add(vProfiles[n].sumNormXResiduals_);
1282  if (hY) hY->Add(vProfiles[n].sumYResiduals_); // only if existing
1283  if (hNormY) hNormY->Add(vProfiles[n].sumNormYResiduals_); // dito (pxl, stripYResiduals_)
1284  }
1285  if(dqmMode_)continue; // No fits in dqmMode
1286  //add fit values to stat box
1287  this->fitResiduals(vLevelProfiles[iComp].sumXResiduals_);
1288  this->fitResiduals(vLevelProfiles[iComp].sumNormXResiduals_);
1289  if (hY) this->fitResiduals(hY); // only if existing (pixel or stripYResiduals_)
1290  if (hNormY) this->fitResiduals(hNormY); // dito
1291  } else {
1292  // nothing to be done for det or detunits
1293  continue;
1294  }
1295  }
1296 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
void collateSummaryHists(DirectoryWrapper &tfd, const Alignable &ali, int i, std::vector< TrackerOfflineValidation::SummaryContainer > &vLevelProfiles)
virtual Alignables components() const =0
Return vector of all direct components.
TrackerOfflineValidation::SummaryContainer bookSummaryHists(DirectoryWrapper &tfd, const Alignable &ali, align::StructureType type, int i)
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
double f[11][100]
bool isDetOrDetUnit(align::StructureType type)
std::pair< float, float > fitResiduals(TH1 *hist) const
void summarizeBinInContainer(int bin, SummaryContainer &targetContainer, SummaryContainer &sourceContainer)
static const char * idToString(align::StructureType type)
tuple size
Write out results.
void TrackerOfflineValidation::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 1207 of file TrackerOfflineValidation.cc.

References collateSummaryHists(), dqmMode_, f, fillTree(), edm::EventSetup::get(), GetIndex(), lastSetup_, TFileService::make(), moduleDirectory_, mPxbResiduals_, mPxeResiduals_, mTecResiduals_, mTibResiduals_, mTidResiduals_, mTobResiduals_, edm::ESHandle< class >::product(), tkGeom_, MainPageGenerator::tree, and vTrackHistos_.

1208 {
1209 
1210  if (!tkGeom_.product()) return;
1211 
1212  //Retrieve tracker topology from geometry
1213  edm::ESHandle<TrackerTopology> tTopoHandle;
1214  lastSetup_->get<TrackerTopologyRcd>().get(tTopoHandle);
1215  const TrackerTopology* const tTopo = tTopoHandle.product();
1216 
1217  AlignableTracker aliTracker(&(*tkGeom_), tTopo);
1218 
1219  static const int kappadiffindex = this->GetIndex(vTrackHistos_,"h_diff_curvature");
1220  vTrackHistos_[kappadiffindex]->Add(vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_neg")],
1221  vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_pos")],-1,1);
1222 
1223  // Collate Information for Subdetectors
1224  // create summary histogramms recursively
1225  std::vector<TrackerOfflineValidation::SummaryContainer> vTrackerprofiles;
1226  DirectoryWrapper f("",moduleDirectory_,dqmMode_);
1227  this->collateSummaryHists(f,(aliTracker), 0, vTrackerprofiles);
1228 
1229  if (dqmMode_) return;
1230  // Should be excluded in dqmMode, since TTree is not usable
1231  // In dqmMode tree operations are are sourced out to the additional module TrackerOfflineValidationSummary
1232 
1234  TTree *tree = fs->make<TTree>("TkOffVal","TkOffVal");
1235 
1236  TkOffTreeVariables *treeMemPtr = new TkOffTreeVariables;
1237  // We create branches for all members of 'TkOffTreeVariables' (even if not needed).
1238  // This works because we have a dictionary for 'TkOffTreeVariables'
1239  // (see src/classes_def.xml and src/classes.h):
1240  tree->Branch("TkOffTreeVariables", &treeMemPtr); // address of pointer!
1241 
1242  this->fillTree(*tree, mPxbResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1243  this->fillTree(*tree, mPxeResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1244  this->fillTree(*tree, mTibResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1245  this->fillTree(*tree, mTidResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1246  this->fillTree(*tree, mTobResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1247  this->fillTree(*tree, mTecResiduals_, *treeMemPtr, *tkGeom_, tTopo);
1248 
1249  delete treeMemPtr; treeMemPtr = 0;
1250 }
std::map< int, TrackerOfflineValidation::ModuleHistos > mPxeResiduals_
void collateSummaryHists(DirectoryWrapper &tfd, const Alignable &ali, int i, std::vector< TrackerOfflineValidation::SummaryContainer > &vLevelProfiles)
std::map< int, TrackerOfflineValidation::ModuleHistos > mTecResiduals_
container to hold data to be written into TTree
edm::ESHandle< TrackerGeometry > tkGeom_
int GetIndex(const std::vector< OBJECT_TYPE * > &vec, const TString &name)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void fillTree(TTree &tree, const std::map< int, TrackerOfflineValidation::ModuleHistos > &moduleHist_, TkOffTreeVariables &treeMem, const TrackerGeometry &tkgeom, const TrackerTopology *tTopo)
std::map< int, TrackerOfflineValidation::ModuleHistos > mTibResiduals_
std::map< int, TrackerOfflineValidation::ModuleHistos > mTidResiduals_
double f[11][100]
std::map< int, TrackerOfflineValidation::ModuleHistos > mTobResiduals_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
std::map< int, TrackerOfflineValidation::ModuleHistos > mPxbResiduals_
const edm::EventSetup * lastSetup_
void TrackerOfflineValidation::fillTree ( TTree &  tree,
const std::map< int, TrackerOfflineValidation::ModuleHistos > &  moduleHist_,
TkOffTreeVariables treeMem,
const TrackerGeometry tkgeom,
const TrackerTopology tTopo 
)
private

Definition at line 1467 of file TrackerOfflineValidation.cc.

References TkOffTreeVariables::blade, TkOffTreeVariables::chi2PerDofX, TkOffTreeVariables::chi2PerDofY, TkOffTreeVariables::clear(), SiPixelRawToDigiRegional_cfi::deltaPhi, dPhi(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, TkOffTreeVariables::entries, PV3DBase< T, PVType, FrameType >::eta(), TkOffTreeVariables::fitMeanNormX, TkOffTreeVariables::fitMeanNormY, TkOffTreeVariables::fitMeanX, TkOffTreeVariables::fitMeanY, fitResiduals(), TkOffTreeVariables::fitSigmaNormX, TkOffTreeVariables::fitSigmaNormY, TkOffTreeVariables::fitSigmaX, TkOffTreeVariables::fitSigmaY, Fwhm(), getMedian(), h, TkOffTreeVariables::half, TkOffTreeVariables::histNameLocalX, TkOffTreeVariables::histNameLocalY, TkOffTreeVariables::histNameNormLocalX, TkOffTreeVariables::histNameNormX, TkOffTreeVariables::histNameNormY, TkOffTreeVariables::histNameX, TkOffTreeVariables::histNameY, TrackerGeometry::idToDet(), TkOffTreeVariables::isDoubleSide, TkOffTreeVariables::isStereo, TkOffTreeVariables::layer, lCoorHistOn_, TkOffTreeVariables::meanLocalX, TkOffTreeVariables::meanNormLocalX, TkOffTreeVariables::meanNormX, TkOffTreeVariables::meanNormY, TkOffTreeVariables::meanResXvsX, TkOffTreeVariables::meanResXvsY, TkOffTreeVariables::meanResYvsX, TkOffTreeVariables::meanResYvsY, TkOffTreeVariables::meanX, TkOffTreeVariables::meanY, TkOffTreeVariables::medianX, TkOffTreeVariables::medianY, TkOffTreeVariables::module, TkOffTreeVariables::moduleId, moduleLevelProfiles_, TkOffTreeVariables::numberOfOutliers, TkOffTreeVariables::numberOfOverflows, TkOffTreeVariables::numberOfUnderflows, TkOffTreeVariables::outerInner, TkOffTreeVariables::panel, PV3DBase< T, PVType, FrameType >::perp(), TkOffTreeVariables::petal, PV3DBase< T, PVType, FrameType >::phi(), TkOffTreeVariables::phiDirection, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, TkOffTreeVariables::posEta, GeomDet::position(), TkOffTreeVariables::posPhi, TkOffTreeVariables::posR, TkOffTreeVariables::posX, TkOffTreeVariables::posY, TkOffTreeVariables::posZ, TkOffTreeVariables::profileNameResXvsX, TkOffTreeVariables::profileNameResXvsY, TkOffTreeVariables::profileNameResYvsX, TkOffTreeVariables::profileNameResYvsY, TrackerTopology::pxbLadder(), TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), TrackerTopology::pxfBlade(), TrackerTopology::pxfDisk(), TrackerTopology::pxfModule(), TrackerTopology::pxfPanel(), TrackerTopology::pxfSide(), DetId::rawId(), TkOffTreeVariables::rDirection, TkOffTreeVariables::ring, TkOffTreeVariables::rmsLocalX, TkOffTreeVariables::rmsNormLocalX, TkOffTreeVariables::rmsNormX, TkOffTreeVariables::rmsNormY, TkOffTreeVariables::rmsResXvsX, TkOffTreeVariables::rmsResXvsY, TkOffTreeVariables::rmsResYvsX, TkOffTreeVariables::rmsResYvsY, TkOffTreeVariables::rmsX, TkOffTreeVariables::rmsY, TkOffTreeVariables::rod, TkOffTreeVariables::rOrZDirection, TkOffTreeVariables::side, TkOffTreeVariables::sigmaNormX, DetId::subdetId(), TkOffTreeVariables::subDetId, GeomDet::surface(), StripSubdetector::TEC, TrackerTopology::tecIsDoubleSide(), TrackerTopology::tecModule(), TrackerTopology::tecPetalInfo(), TrackerTopology::tecRing(), TrackerTopology::tecSide(), TrackerTopology::tecStereo(), TrackerTopology::tecWheel(), StripSubdetector::TIB, TrackerTopology::tibIsDoubleSide(), TrackerTopology::tibLayer(), TrackerTopology::tibModule(), TrackerTopology::tibStereo(), TrackerTopology::tibStringInfo(), StripSubdetector::TID, TrackerTopology::tidIsDoubleSide(), TrackerTopology::tidModuleInfo(), TrackerTopology::tidRing(), TrackerTopology::tidSide(), TrackerTopology::tidStereo(), TrackerTopology::tidWheel(), StripSubdetector::TOB, TrackerTopology::tobIsDoubleSide(), TrackerTopology::tobLayer(), TrackerTopology::tobModule(), TrackerTopology::tobRodInfo(), TrackerTopology::tobStereo(), Surface::toGlobal(), useFit_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and TkOffTreeVariables::zDirection.

Referenced by endJob(), and core.AutoFillTreeProducer.AutoFillTreeProducer::process().

1470 {
1471 
1472  for(std::map<int, TrackerOfflineValidation::ModuleHistos>::const_iterator it = moduleHist_.begin(),
1473  itEnd= moduleHist_.end(); it != itEnd;++it ) {
1474  treeMem.clear(); // make empty/default
1475 
1476  //variables concerning the tracker components/hierarchy levels
1477  DetId detId_ = it->first;
1478  treeMem.moduleId = detId_;
1479  treeMem.subDetId = detId_.subdetId();
1480  treeMem.isDoubleSide =0;
1481 
1482  if(treeMem.subDetId == PixelSubdetector::PixelBarrel){
1483  unsigned int whichHalfBarrel(1), rawId(detId_.rawId()); //DetId does not know about halfBarrels is PXB ...
1484  if( (rawId>=302056964 && rawId<302059300) || (rawId>=302123268 && rawId<302127140) ||
1485  (rawId>=302189572 && rawId<302194980) ) whichHalfBarrel=2;
1486  treeMem.layer = tTopo->pxbLayer(detId_);
1487  treeMem.half = whichHalfBarrel;
1488  treeMem.rod = tTopo->pxbLadder(detId_); // ... so, ladder is not per halfBarrel-Layer, but per barrel-layer!
1489  treeMem.module = tTopo->pxbModule(detId_);
1490  } else if(treeMem.subDetId == PixelSubdetector::PixelEndcap){
1491  unsigned int whichHalfCylinder(1), rawId(detId_.rawId()); //DetId does not kmow about halfCylinders in PXF
1492  if( (rawId>=352394500 && rawId<352406032) || (rawId>=352460036 && rawId<352471568) ||
1493  (rawId>=344005892 && rawId<344017424) || (rawId>=344071428 && rawId<344082960) ) whichHalfCylinder=2;
1494  treeMem.layer = tTopo->pxfDisk(detId_);
1495  treeMem.side = tTopo->pxfSide(detId_);
1496  treeMem.half = whichHalfCylinder;
1497  treeMem.blade = tTopo->pxfBlade(detId_);
1498  treeMem.panel = tTopo->pxfPanel(detId_);
1499  treeMem.module = tTopo->pxfModule(detId_);
1500  } else if(treeMem.subDetId == StripSubdetector::TIB){
1501  unsigned int whichHalfShell(1), rawId(detId_.rawId()); //DetId does not kmow about halfShells in TIB
1502  if ( (rawId>=369120484 && rawId<369120688) || (rawId>=369121540 && rawId<369121776) ||
1503  (rawId>=369136932 && rawId<369137200) || (rawId>=369137988 && rawId<369138288) ||
1504  (rawId>=369153396 && rawId<369153744) || (rawId>=369154436 && rawId<369154800) ||
1505  (rawId>=369169844 && rawId<369170256) || (rawId>=369170900 && rawId<369171344) ||
1506  (rawId>=369124580 && rawId<369124784) || (rawId>=369125636 && rawId<369125872) ||
1507  (rawId>=369141028 && rawId<369141296) || (rawId>=369142084 && rawId<369142384) ||
1508  (rawId>=369157492 && rawId<369157840) || (rawId>=369158532 && rawId<369158896) ||
1509  (rawId>=369173940 && rawId<369174352) || (rawId>=369174996 && rawId<369175440) ) whichHalfShell=2;
1510  treeMem.layer = tTopo->tibLayer(detId_);
1511  treeMem.side = tTopo->tibStringInfo(detId_)[0];
1512  treeMem.half = whichHalfShell;
1513  treeMem.rod = tTopo->tibStringInfo(detId_)[2];
1514  treeMem.outerInner = tTopo->tibStringInfo(detId_)[1];
1515  treeMem.module = tTopo->tibModule(detId_);
1516  treeMem.isStereo = tTopo->tibStereo(detId_);
1517  treeMem.isDoubleSide = tTopo->tibIsDoubleSide(detId_);
1518  } else if(treeMem.subDetId == StripSubdetector::TID){
1519  treeMem.layer = tTopo->tidWheel(detId_);
1520  treeMem.side = tTopo->tidSide(detId_);
1521  treeMem.ring = tTopo->tidRing(detId_);
1522  treeMem.outerInner = tTopo->tidModuleInfo(detId_)[0];
1523  treeMem.module = tTopo->tidModuleInfo(detId_)[1];
1524  treeMem.isStereo = tTopo->tidStereo(detId_);
1525  treeMem.isDoubleSide = tTopo->tidIsDoubleSide(detId_);
1526  } else if(treeMem.subDetId == StripSubdetector::TOB){
1527  treeMem.layer = tTopo->tobLayer(detId_);
1528  treeMem.side = tTopo->tobRodInfo(detId_)[0];
1529  treeMem.rod = tTopo->tobRodInfo(detId_)[1];
1530  treeMem.module = tTopo->tobModule(detId_);
1531  treeMem.isStereo = tTopo->tobStereo(detId_);
1532  treeMem.isDoubleSide = tTopo->tobIsDoubleSide(detId_);
1533  } else if(treeMem.subDetId == StripSubdetector::TEC) {
1534  treeMem.layer = tTopo->tecWheel(detId_);
1535  treeMem.side = tTopo->tecSide(detId_);
1536  treeMem.ring = tTopo->tecRing(detId_);
1537  treeMem.petal = tTopo->tecPetalInfo(detId_)[1];
1538  treeMem.outerInner = tTopo->tecPetalInfo(detId_)[0];
1539  treeMem.module = tTopo->tecModule(detId_);
1540  treeMem.isStereo = tTopo->tecStereo(detId_);
1541  treeMem.isDoubleSide = tTopo->tecIsDoubleSide(detId_);
1542  }
1543 
1544  //variables concerning the tracker geometry
1545  const Surface::PositionType &gPModule = tkgeom.idToDet(detId_)->position();
1546  treeMem.posPhi = gPModule.phi();
1547  treeMem.posEta = gPModule.eta();
1548  treeMem.posR = gPModule.perp();
1549  treeMem.posX = gPModule.x();
1550  treeMem.posY = gPModule.y();
1551  treeMem.posZ = gPModule.z();
1552 
1553  const Surface& surface = tkgeom.idToDet(detId_)->surface();
1554 
1555  //global Orientation of local coordinate system of dets/detUnits
1556  LocalPoint lUDirection(1.,0.,0.), lVDirection(0.,1.,0.), lWDirection(0.,0.,1.);
1557  GlobalPoint gUDirection = surface.toGlobal(lUDirection),
1558  gVDirection = surface.toGlobal(lVDirection),
1559  gWDirection = surface.toGlobal(lWDirection);
1560  double dR(999.), dPhi(999.), dZ(999.);
1562  dR = gWDirection.perp() - gPModule.perp();
1563  dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
1564  dZ = gVDirection.z() - gPModule.z();
1565  if(dZ>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
1566  }else if(treeMem.subDetId==PixelSubdetector::PixelEndcap){
1567  dR = gUDirection.perp() - gPModule.perp();
1568  dPhi = deltaPhi(gVDirection.phi(),gPModule.phi());
1569  dZ = gWDirection.z() - gPModule.z();
1570  if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
1571  }else if(treeMem.subDetId==StripSubdetector::TID || treeMem.subDetId==StripSubdetector::TEC){
1572  dR = gVDirection.perp() - gPModule.perp();
1573  dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
1574  dZ = gWDirection.z() - gPModule.z();
1575  if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
1576  }
1577  if(dR>=0.)treeMem.rDirection = 1; else treeMem.rDirection = -1;
1578  if(dPhi>=0.)treeMem.phiDirection = 1; else treeMem.phiDirection = -1;
1579  if(dZ>=0.)treeMem.zDirection = 1; else treeMem.zDirection = -1;
1580 
1581  //mean and RMS values (extracted from histograms(Xprime on module level)
1582  treeMem.entries = static_cast<UInt_t>(it->second.ResXprimeHisto->GetEntries());
1583  treeMem.meanX = it->second.ResXprimeHisto->GetMean();
1584  treeMem.rmsX = it->second.ResXprimeHisto->GetRMS();
1585  //treeMem.sigmaX = Fwhm(it->second.ResXprimeHisto)/2.355;
1586 
1587  if (useFit_) {
1588  //call fit function which returns mean and sigma from the fit
1589  //for absolute residuals
1590  std::pair<float,float> fitResult1 = this->fitResiduals(it->second.ResXprimeHisto);
1591  treeMem.fitMeanX = fitResult1.first;
1592  treeMem.fitSigmaX = fitResult1.second;
1593  //for normalized residuals
1594  std::pair<float,float> fitResult2 = this->fitResiduals(it->second.NormResXprimeHisto);
1595  treeMem.fitMeanNormX = fitResult2.first;
1596  treeMem.fitSigmaNormX = fitResult2.second;
1597  }
1598 
1599  //get median for absolute residuals
1600  treeMem.medianX = this->getMedian(it->second.ResXprimeHisto);
1601 
1602  int numberOfBins=it->second.ResXprimeHisto->GetNbinsX();
1603  treeMem.numberOfUnderflows = it->second.ResXprimeHisto->GetBinContent(0);
1604  treeMem.numberOfOverflows = it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
1605  treeMem.numberOfOutliers = it->second.ResXprimeHisto->GetBinContent(0)+it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
1606 
1607  //mean and RMS values (extracted from histograms(normalized Xprime on module level)
1608  treeMem.meanNormX = it->second.NormResXprimeHisto->GetMean();
1609  treeMem.rmsNormX = it->second.NormResXprimeHisto->GetRMS();
1610 
1611  double stats[20];
1612  it->second.NormResXprimeHisto->GetStats(stats);
1613  // GF treeMem.chi2PerDofX = stats[3]/(stats[0]-1);
1614  if (stats[0]) treeMem.chi2PerDofX = stats[3]/stats[0];
1615 
1616  treeMem.sigmaNormX = Fwhm(it->second.NormResXprimeHisto)/2.355;
1617  treeMem.histNameX = it->second.ResXprimeHisto->GetName();
1618  treeMem.histNameNormX = it->second.NormResXprimeHisto->GetName();
1619 
1620  // fill tree variables in local coordinates if set in cfg
1621  if(lCoorHistOn_) {
1622  treeMem.meanLocalX = it->second.ResHisto->GetMean();
1623  treeMem.rmsLocalX = it->second.ResHisto->GetRMS();
1624  treeMem.meanNormLocalX = it->second.NormResHisto->GetMean();
1625  treeMem.rmsNormLocalX = it->second.NormResHisto->GetRMS();
1626 
1627  treeMem.histNameLocalX = it->second.ResHisto->GetName();
1628  treeMem.histNameNormLocalX = it->second.NormResHisto->GetName();
1629  if (it->second.ResYHisto) treeMem.histNameLocalY = it->second.ResYHisto->GetName();
1630  }
1631 
1632  // mean and RMS values in local y (extracted from histograms(normalized Yprime on module level)
1633  // might exist in pixel only
1634  if (it->second.ResYprimeHisto) {//(stripYResiduals_)
1635  TH1 *h = it->second.ResYprimeHisto;
1636  treeMem.meanY = h->GetMean();
1637  treeMem.rmsY = h->GetRMS();
1638 
1639  if (useFit_) { // fit function which returns mean and sigma from the fit
1640  std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
1641  treeMem.fitMeanY = fitMeanSigma.first;
1642  treeMem.fitSigmaY = fitMeanSigma.second;
1643  }
1644 
1645  //get median for absolute residuals
1646  treeMem.medianY = this->getMedian(h);
1647 
1648  treeMem.histNameY = h->GetName();
1649  }
1650  if (it->second.NormResYprimeHisto) {
1651  TH1 *h = it->second.NormResYprimeHisto;
1652  treeMem.meanNormY = h->GetMean();
1653  treeMem.rmsNormY = h->GetRMS();
1654  h->GetStats(stats); // stats buffer defined above
1655  if (stats[0]) treeMem.chi2PerDofY = stats[3]/stats[0];
1656 
1657  if (useFit_) { // fit function which returns mean and sigma from the fit
1658  std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
1659  treeMem.fitMeanNormY = fitMeanSigma.first;
1660  treeMem.fitSigmaNormY = fitMeanSigma.second;
1661  }
1662  treeMem.histNameNormY = h->GetName();
1663  }
1664 
1665  if (moduleLevelProfiles_) {
1666  if (it->second.ResXvsXProfile) {
1667  TH1 *h = it->second.ResXvsXProfile;
1668  treeMem.meanResXvsX = h->GetMean();
1669  treeMem.rmsResXvsX = h->GetRMS();
1670  treeMem.profileNameResXvsX = h->GetName();
1671  }
1672  if (it->second.ResXvsYProfile) {
1673  TH1 *h = it->second.ResXvsYProfile;
1674  treeMem.meanResXvsY = h->GetMean();
1675  treeMem.rmsResXvsY = h->GetRMS();
1676  treeMem.profileNameResXvsY = h->GetName();
1677  }
1678  if (it->second.ResYvsXProfile) {
1679  TH1 *h = it->second.ResYvsXProfile;
1680  treeMem.meanResYvsX = h->GetMean();
1681  treeMem.rmsResYvsX = h->GetRMS();
1682  treeMem.profileNameResYvsX = h->GetName();
1683  }
1684  if (it->second.ResYvsYProfile) {
1685  TH1 *h = it->second.ResYvsYProfile;
1686  treeMem.meanResYvsY = h->GetMean();
1687  treeMem.rmsResYvsY = h->GetRMS();
1688  treeMem.profileNameResYvsY = h->GetName();
1689  }
1690  }
1691 
1692  tree.Fill();
1693  }
1694 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
bool tecIsDoubleSide(const DetId &id) const
bool tobIsDoubleSide(const DetId &id) const
bool tibIsDoubleSide(const DetId &id) const
T perp() const
Definition: PV3DBase.h:72
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
std::vector< unsigned int > tidModuleInfo(const DetId &id) const
unsigned int pxfDisk(const DetId &id) const
std::string profileNameResXvsX
unsigned int tecRing(const DetId &id) const
ring id
uint32_t tobStereo(const DetId &id) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
unsigned int pxbLadder(const DetId &id) const
T y() const
Definition: PV3DBase.h:63
unsigned int tidWheel(const DetId &id) const
unsigned int pxbModule(const DetId &id) const
float getMedian(const TH1 *hist) const
std::vector< unsigned int > tibStringInfo(const DetId &id) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
std::string profileNameResYvsX
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
std::vector< unsigned int > tecPetalInfo(const DetId &id) const
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:46
unsigned int tidSide(const DetId &id) const
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
uint32_t tidStereo(const DetId &id) const
std::vector< unsigned int > tobRodInfo(const DetId &id) const
T z() const
Definition: PV3DBase.h:64
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::string histNameLocalY
std::string histNameNormLocalX
unsigned int tibModule(const DetId &id) const
unsigned int pxfModule(const DetId &id) const
unsigned int pxbLayer(const DetId &id) const
unsigned int tecModule(const DetId &id) const
Definition: DetId.h:18
std::pair< float, float > fitResiduals(TH1 *hist) const
std::string histNameLocalX
bool tidIsDoubleSide(const DetId &id) const
unsigned int tobModule(const DetId &id) const
std::string profileNameResXvsY
T eta() const
Definition: PV3DBase.h:76
void clear()
set to empty values
uint32_t tecStereo(const DetId &id) const
unsigned int pxfSide(const DetId &id) const
uint32_t tibStereo(const DetId &id) const
T x() const
Definition: PV3DBase.h:62
float Fwhm(const TH1 *hist) const
unsigned int tecWheel(const DetId &id) const
unsigned int pxfPanel(const DetId &id) const
unsigned int pxfBlade(const DetId &id) const
std::string profileNameResYvsY
unsigned int tobLayer(const DetId &id) const
unsigned int tecSide(const DetId &id) const
virtual const TrackerGeomDet * idToDet(DetId) const
std::pair< float, float > TrackerOfflineValidation::fitResiduals ( TH1 *  hist) const
private

Definition at line 1698 of file TrackerOfflineValidation.cc.

References alignCSCRings::e, timingPdfMaker::mean, and cms::Exception::what().

Referenced by collateSummaryHists(), and fillTree().

1699 {
1700  std::pair<float,float> fitResult(9999., 9999.);
1701  if (!hist || hist->GetEntries() < 20) return fitResult;
1702 
1703  float mean = hist->GetMean();
1704  float sigma = hist->GetRMS();
1705 
1706  try { // for < CMSSW_2_2_0 since ROOT warnings from fit are converted to exceptions
1707  // Remove the try/catch for more recent CMSSW!
1708  // first fit: two RMS around mean
1709  TF1 func("tmp", "gaus", mean - 2.*sigma, mean + 2.*sigma);
1710  if (0 == hist->Fit(&func,"QNR")) { // N: do not blow up file by storing fit!
1711  mean = func.GetParameter(1);
1712  sigma = func.GetParameter(2);
1713  // second fit: three sigma of first fit around mean of first fit
1714  func.SetRange(mean - 3.*sigma, mean + 3.*sigma);
1715  // I: integral gives more correct results if binning is too wide
1716  // L: Likelihood can treat empty bins correctly (if hist not weighted...)
1717  if (0 == hist->Fit(&func, "Q0LR")) {
1718  if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
1719  hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
1720  }
1721  fitResult.first = func.GetParameter(1);
1722  fitResult.second = func.GetParameter(2);
1723  }
1724  }
1725  } catch (cms::Exception const & e) {
1726  edm::LogWarning("Alignment") << "@SUB=TrackerOfflineValidation::fitResiduals"
1727  << "Caught this exception during ROOT fit: "
1728  << e.what();
1729  }
1730  return fitResult;
1731 }
virtual char const * what() const
Definition: Exception.cc:141
float TrackerOfflineValidation::Fwhm ( const TH1 *  hist) const
private

Definition at line 1437 of file TrackerOfflineValidation.cc.

References i, and bookConverter::max.

Referenced by fillTree(), and setSummaryBin().

1438 {
1439  float max = hist->GetMaximum();
1440  int left = -1, right = -1;
1441  for(unsigned int i = 1, iEnd = hist->GetNbinsX(); i <= iEnd; ++i) {
1442  if(hist->GetBinContent(i) < max/2. && hist->GetBinContent(i+1) > max/2. && left == -1) {
1443  if(max/2. - hist->GetBinContent(i) < hist->GetBinContent(i+1) - max/2.) {
1444  left = i;
1445  ++i;
1446  } else {
1447  left = i+1;
1448  ++i;
1449  }
1450  }
1451  if(left != -1 && right == -1) {
1452  if(hist->GetBinContent(i) > max/2. && hist->GetBinContent(i+1) < max/2.) {
1453  if( hist->GetBinContent(i) - max/2. < max/2. - hist->GetBinContent(i+1)) {
1454  right = i;
1455  } else {
1456  right = i+1;
1457  }
1458 
1459  }
1460  }
1461  }
1462  return hist->GetXaxis()->GetBinCenter(right) - hist->GetXaxis()->GetBinCenter(left);
1463 }
int i
Definition: DBlmapReader.cc:9
void TrackerOfflineValidation::getBinning ( uint32_t  subDetId,
TrackerOfflineValidation::HistogrammType  residualtype,
int &  nBinsX,
double &  lowerBoundX,
double &  upperBoundX 
)
private

Definition at line 852 of file TrackerOfflineValidation.cc.

References edm::ParameterSet::getParameter(), isPixel(), NormXprimeResidual, NormXResidual, NormYprimeResidual, parSet_, XprimeResidual, XResidual, XResidualProfile, YprimeResidual, YResidual, and YResidualProfile.

Referenced by bookHists(), and bookSummaryHists().

855 {
856  // determine if
857  const bool isPixel = this->isPixel(subDetId);
858 
859  edm::ParameterSet binningPSet;
860 
861  switch(residualType)
862  {
863  case XResidual :
864  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XResPixelModules");
865  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XResStripModules");
866  break;
867  case NormXResidual :
868  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXResPixelModules");
869  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXResStripModules");
870  break;
871  case XprimeResidual :
872  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResPixelModules");
873  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResStripModules");
874  break;
875  case NormXprimeResidual :
876  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResPixelModules");
877  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResStripModules");
878  break;
879  case YResidual : // borrow y-residual binning from yprime
880  case YprimeResidual :
881  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1YResPixelModules");
882  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1YResStripModules");
883  break;
884  /* case NormYResidual :*/
885  case NormYprimeResidual :
886  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormYResPixelModules");
887  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormYResStripModules");
888  break;
889  case XResidualProfile :
890  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileXResPixelModules");
891  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileXResStripModules");
892  break;
893  case YResidualProfile :
894  if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileYResPixelModules");
895  else binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileYResStripModules");
896  break;
897  }
898  nBinsX = binningPSet.getParameter<int32_t>("Nbinx");
899  lowerBoundX = binningPSet.getParameter<double>("xmin");
900  upperBoundX = binningPSet.getParameter<double>("xmax");
901 }
T getParameter(std::string const &) const
bool isPixel(uint32_t subDetId)
unsigned int subDetId[18]
const edm::ParameterSet parSet_
TrackerOfflineValidation::ModuleHistos & TrackerOfflineValidation::getHistStructFromMap ( const DetId detid)
private

Definition at line 944 of file TrackerOfflineValidation.cc.

References DetId::det(), edm::hlt::Exception, mPxbResiduals_, mPxeResiduals_, mTecResiduals_, mTibResiduals_, mTidResiduals_, mTobResiduals_, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::rawId(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.

Referenced by analyze(), bookHists(), and bookSummaryHists().

945 {
946  // get a struct with histograms from the respective map
947  // if no object exist, the reference is automatically created by the map
948  // throw exception if non-tracker id is passed
949  uint subdetid = detid.subdetId();
950  if(subdetid == PixelSubdetector::PixelBarrel ) {
951  return mPxbResiduals_[detid.rawId()];
952  } else if (subdetid == PixelSubdetector::PixelEndcap) {
953  return mPxeResiduals_[detid.rawId()];
954  } else if(subdetid == StripSubdetector::TIB) {
955  return mTibResiduals_[detid.rawId()];
956  } else if(subdetid == StripSubdetector::TID) {
957  return mTidResiduals_[detid.rawId()];
958  } else if(subdetid == StripSubdetector::TOB) {
959  return mTobResiduals_[detid.rawId()];
960  } else if(subdetid == StripSubdetector::TEC) {
961  return mTecResiduals_[detid.rawId()];
962  } else {
963  throw cms::Exception("Geometry Error")
964  << "[TrackerOfflineValidation] Error, tried to get reference for non-tracker subdet " << subdetid
965  << " from detector " << detid.det();
966  return mPxbResiduals_[0];
967  }
968 }
std::map< int, TrackerOfflineValidation::ModuleHistos > mPxeResiduals_
std::map< int, TrackerOfflineValidation::ModuleHistos > mTecResiduals_
std::map< int, TrackerOfflineValidation::ModuleHistos > mTibResiduals_
std::map< int, TrackerOfflineValidation::ModuleHistos > mTidResiduals_
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::map< int, TrackerOfflineValidation::ModuleHistos > mTobResiduals_
std::map< int, TrackerOfflineValidation::ModuleHistos > mPxbResiduals_
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
template<class OBJECT_TYPE >
int TrackerOfflineValidation::GetIndex ( const std::vector< OBJECT_TYPE * > &  vec,
const TString &  name 
)
private

Definition at line 296 of file TrackerOfflineValidation.cc.

References mergeVDriftHistosByStation::name, and query::result.

Referenced by analyze(), and endJob().

297 {
298  int result = 0;
299  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
300  iter != iterEnd; ++iter, ++result) {
301  if (*iter && (*iter)->GetName() == name) return result;
302  }
303  edm::LogError("Alignment") << "@SUB=TrackerOfflineValidation::GetIndex" << " could not find " << name;
304  return -1;
305 }
tuple result
Definition: query.py:137
float TrackerOfflineValidation::getMedian ( const TH1 *  hist) const
private

Definition at line 1735 of file TrackerOfflineValidation.cc.

References j, pileupCalc::nbins, x, and y.

Referenced by fillTree().

1736 {
1737  float median = 999;
1738  int nbins = histo->GetNbinsX();
1739 
1740  //extract median from histogram
1741  double *x = new double[nbins];
1742  double *y = new double[nbins];
1743  for (int j = 0; j < nbins; j++) {
1744  x[j] = histo->GetBinCenter(j+1);
1745  y[j] = histo->GetBinContent(j+1);
1746  }
1747  median = TMath::Median(nbins, x, y);
1748 
1749  delete[] x; x = 0;
1750  delete [] y; y = 0;
1751 
1752  return median;
1753 
1754 }
int j
Definition: DBlmapReader.cc:9
bool TrackerOfflineValidation::isBarrel ( uint32_t  subDetId)
private
bool TrackerOfflineValidation::isDetOrDetUnit ( align::StructureType  type)
private
bool TrackerOfflineValidation::isEndCap ( uint32_t  subDetId)
private
bool TrackerOfflineValidation::isPixel ( uint32_t  subDetId)
private
void TrackerOfflineValidation::setSummaryBin ( int  bin,
TH1 *  targetHist,
TH1 *  sourceHist 
)
private

Definition at line 905 of file TrackerOfflineValidation.cc.

References Fwhm(), and useFwhm_.

Referenced by summarizeBinInContainer().

906 {
907  if(targetHist && sourceHist) {
908  targetHist->SetBinContent(bin, sourceHist->GetMean(1));
909  if(useFwhm_) targetHist->SetBinError(bin, Fwhm(sourceHist)/2.);
910  else targetHist->SetBinError(bin, sourceHist->GetRMS(1) );
911  }
912  else return;
913 }
float Fwhm(const TH1 *hist) const
void TrackerOfflineValidation::summarizeBinInContainer ( int  bin,
SummaryContainer targetContainer,
SummaryContainer sourceContainer 
)
private

Definition at line 917 of file TrackerOfflineValidation.cc.

References setSummaryBin(), TrackerOfflineValidation::SummaryContainer::summaryNormXResiduals_, TrackerOfflineValidation::SummaryContainer::summaryNormYResiduals_, TrackerOfflineValidation::SummaryContainer::summaryXResiduals_, TrackerOfflineValidation::SummaryContainer::summaryYResiduals_, TrackerOfflineValidation::SummaryContainer::sumNormXResiduals_, TrackerOfflineValidation::SummaryContainer::sumNormYResiduals_, TrackerOfflineValidation::SummaryContainer::sumXResiduals_, and TrackerOfflineValidation::SummaryContainer::sumYResiduals_.

Referenced by bookSummaryHists(), and collateSummaryHists().

919 {
920  this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.sumXResiduals_);
921  this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.sumNormXResiduals_);
922  // If no y-residual hists, just returns:
923  this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.sumYResiduals_);
924  this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.sumNormYResiduals_);
925 }
void setSummaryBin(int bin, TH1 *targetHist, TH1 *sourceHist)
void TrackerOfflineValidation::summarizeBinInContainer ( int  bin,
uint32_t  subDetId,
SummaryContainer targetContainer,
ModuleHistos sourceContainer 
)
private

Definition at line 929 of file TrackerOfflineValidation.cc.

References isPixel(), TrackerOfflineValidation::ModuleHistos::NormResXprimeHisto, TrackerOfflineValidation::ModuleHistos::NormResYprimeHisto, TrackerOfflineValidation::ModuleHistos::ResXprimeHisto, TrackerOfflineValidation::ModuleHistos::ResYprimeHisto, setSummaryBin(), stripYResiduals_, TrackerOfflineValidation::SummaryContainer::summaryNormXResiduals_, TrackerOfflineValidation::SummaryContainer::summaryNormYResiduals_, TrackerOfflineValidation::SummaryContainer::summaryXResiduals_, and TrackerOfflineValidation::SummaryContainer::summaryYResiduals_.

932 {
933  // takes two summary Containers and sets summaryBins for all histograms
934  this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.ResXprimeHisto);
935  this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.NormResXprimeHisto);
936  if( this->isPixel(subDetId) || stripYResiduals_ ) {
937  this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.ResYprimeHisto);
938  this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.NormResYprimeHisto);
939  }
940 }
bool isPixel(uint32_t subDetId)
void setSummaryBin(int bin, TH1 *targetHist, TH1 *sourceHist)
unsigned int subDetId[18]

Member Data Documentation

TrackerValidationVariables TrackerOfflineValidation::avalidator_
private

Definition at line 283 of file TrackerOfflineValidation.cc.

Referenced by analyze().

const TrackerGeometry* TrackerOfflineValidation::bareTkGeomPtr_
private

Definition at line 254 of file TrackerOfflineValidation.cc.

Referenced by checkBookHists().

const bool TrackerOfflineValidation::dqmMode_
private
const edm::EventSetup* TrackerOfflineValidation::lastSetup_
private

Definition at line 281 of file TrackerOfflineValidation.cc.

Referenced by endJob().

const bool TrackerOfflineValidation::lCoorHistOn_
private

Definition at line 257 of file TrackerOfflineValidation.cc.

Referenced by analyze(), bookHists(), and fillTree().

const std::string TrackerOfflineValidation::moduleDirectory_
private
const bool TrackerOfflineValidation::moduleLevelHistsTransient_
private

Definition at line 258 of file TrackerOfflineValidation.cc.

Referenced by bookHists().

const bool TrackerOfflineValidation::moduleLevelProfiles_
private

Definition at line 259 of file TrackerOfflineValidation.cc.

Referenced by analyze(), bookHists(), and fillTree().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mPxbResiduals_
private

Definition at line 274 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mPxeResiduals_
private

Definition at line 275 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mTecResiduals_
private

Definition at line 279 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mTibResiduals_
private

Definition at line 276 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mTidResiduals_
private

Definition at line 277 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

std::map<int,TrackerOfflineValidation::ModuleHistos> TrackerOfflineValidation::mTobResiduals_
private

Definition at line 278 of file TrackerOfflineValidation.cc.

Referenced by endJob(), and getHistStructFromMap().

const edm::ParameterSet TrackerOfflineValidation::parSet_
private

Definition at line 252 of file TrackerOfflineValidation.cc.

Referenced by getBinning().

const bool TrackerOfflineValidation::stripYResiduals_
private
edm::ESHandle<TrackerGeometry> TrackerOfflineValidation::tkGeom_
private

Definition at line 253 of file TrackerOfflineValidation.cc.

Referenced by checkBookHists(), and endJob().

const bool TrackerOfflineValidation::useFit_
private

Definition at line 262 of file TrackerOfflineValidation.cc.

Referenced by fillTree().

const bool TrackerOfflineValidation::useFwhm_
private

Definition at line 261 of file TrackerOfflineValidation.cc.

Referenced by setSummaryBin().

const bool TrackerOfflineValidation::useOverflowForRMS_
private

Definition at line 263 of file TrackerOfflineValidation.cc.

Referenced by analyze().

std::vector<TH1*> TrackerOfflineValidation::vDeleteObjects_
private

Definition at line 268 of file TrackerOfflineValidation.cc.

Referenced by bookTH1F(), bookTProfile(), and ~TrackerOfflineValidation().

std::vector<TH1*> TrackerOfflineValidation::vTrack2DHistos_
private

Definition at line 272 of file TrackerOfflineValidation.cc.

Referenced by analyze(), and bookGlobalHists().

std::vector<TH1*> TrackerOfflineValidation::vTrackHistos_
private

Definition at line 270 of file TrackerOfflineValidation.cc.

Referenced by analyze(), bookGlobalHists(), and endJob().

std::vector<TH1*> TrackerOfflineValidation::vTrackProfiles_
private

Definition at line 271 of file TrackerOfflineValidation.cc.

Referenced by analyze(), and bookGlobalHists().