CMS 3D CMS Logo

MillePedeMonitor Class Reference

monitoring of MillePedeAlignmentAlgorithm and its input tracks More...

#include <Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeMonitor.h>

List of all members.

Public Member Functions

void fillDerivatives (const TransientTrackingRecHit::ConstRecHitPointer &recHit, const float *localDerivs, unsigned int nLocal, const float *globalDerivs, unsigned int nGlobal)
void fillFrameToFrame (const AlignableDetOrUnitPtr &aliDet, const Alignable *ali)
void fillRefTrajectory (const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr)
void fillResiduals (const TransientTrackingRecHit::ConstRecHitPointer &recHit, const TrajectoryStateOnSurface &tsos, unsigned int nHit, float residuum, float sigma, bool isY)
void fillTrack (const reco::Track *track)
void fillUsedTrack (const reco::Track *track, unsigned int nHitX, unsigned int nHitY)
 MillePedeMonitor (TDirectory *rootDir)
 MillePedeMonitor (const char *rootFile="trackMonitor.root")
 ~MillePedeMonitor ()

Private Member Functions

template<class OBJECT_TYPE>
void addToDirectory (const std::vector< OBJECT_TYPE * > &objs, TDirectory *dir) const
template<class OBJECT_TYPE>
std::vector< OBJECT_TYPE * > cloneHists (const std::vector< OBJECT_TYPE * > &orgs, const TString &namAd, const TString &titAd) const
bool equidistLogBins (double *bins, int nBins, double first, double last) const
void fillResidualHists (const std::vector< TH1 * > &hists, float phiSensToNorm, float residuum, float sigma)
void fillResidualHitHists (const std::vector< TH1 * > &hists, float angle, float residuum, float sigma, unsigned int nHit)
void fillTrack (const reco::Track *track, std::vector< TH1 * > &trackHists1D, std::vector< TH2 * > &trackHists2D)
template<class OBJECT_TYPE>
int GetIndex (const std::vector< OBJECT_TYPE * > &vec, const TString &name)
bool init (TDirectory *directory)

Private Attributes

bool myDeleteDir
std::vector< TH2 * > myDerivHists2D
std::vector< TH2 * > myFrame2FrameHists2D
std::vector< TH2 * > myResidHists2D
std::vector< std::vector< TH1 * > > myResidHistsVec1DX
std::vector< std::vector< TH1 * > > myResidHistsVec1DY
 [0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
std::vector< TH1 * > myResidHitHists1DX
 [0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
std::vector< TH1 * > myResidHitHists1DY
TDirectory * myRootDir
std::vector< TH1 * > myTrackHists1D
std::vector< TH2 * > myTrackHists2D
std::vector< TH1 * > myTrajectoryHists1D
std::vector< TH2 * > myTrajectoryHists2D
std::vector< TH1 * > myUsedTrackHists1D
std::vector< TH2 * > myUsedTrackHists2D


Detailed Description

monitoring of MillePedeAlignmentAlgorithm and its input tracks

Author:
: Gero Flucke date : October 2006
Revision
1.10
Date
2008/03/27 17:49:16
(last update by
Author
flucke
)

Definition at line 42 of file MillePedeMonitor.h.


Constructor & Destructor Documentation

MillePedeMonitor::MillePedeMonitor ( const char *  rootFile = "trackMonitor.root"  ) 

Definition at line 34 of file MillePedeMonitor.cc.

References init(), myDeleteDir, and myRootDir.

00035   : myRootDir(0), myDeleteDir(false)
00036 {
00037   myRootDir = TFile::Open(rootFileName, "recreate");
00038   myDeleteDir = true;
00039 
00040   this->init(myRootDir);
00041 }

MillePedeMonitor::MillePedeMonitor ( TDirectory *  rootDir  ) 

Definition at line 44 of file MillePedeMonitor.cc.

References init(), myDeleteDir, and myRootDir.

00045   : myRootDir(0), myDeleteDir(false)
00046 {
00047   //  cout << "MillePedeMonitor using input TDirectory" << endl;
00048 
00049   myRootDir = rootDir;
00050   myDeleteDir = false;
00051 
00052   this->init(myRootDir);
00053 }

MillePedeMonitor::~MillePedeMonitor (  ) 

Definition at line 56 of file MillePedeMonitor.cc.

References myDeleteDir, and myRootDir.

00057 {
00058 
00059   myRootDir->Write();
00060   if (myDeleteDir) delete myRootDir; //hists are deleted with their directory
00061 }


Member Function Documentation

template<class OBJECT_TYPE>
void MillePedeMonitor::addToDirectory ( const std::vector< OBJECT_TYPE * > &  objs,
TDirectory *  dir 
) const [inline, private]

Definition at line 132 of file MillePedeMonitor.h.

References iter.

Referenced by init().

00134 {
00135   // OBJECT_TYPE is required to have method SetDirectory(TDirectory *dir)
00136   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = obs.begin(), iterEnd = obs.end();
00137        iter != iterEnd; ++iter) {
00138     if (*iter) (*iter)->SetDirectory(dir);
00139   }
00140 }

template<class OBJECT_TYPE>
std::vector< OBJECT_TYPE * > MillePedeMonitor::cloneHists ( const std::vector< OBJECT_TYPE * > &  orgs,
const TString &  namAd,
const TString &  titAd 
) const [inline, private]

Definition at line 112 of file MillePedeMonitor.h.

References iter, and HLT_VtxMuL3::result.

Referenced by init().

00115 {
00116   // OBJECT_TYPE required to have methods Clone(const char*), GetName(), SetTitle(const char*) and GetTitle()
00117   std::vector<OBJECT_TYPE*> result;
00118   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = orgs.begin(), iterEnd = orgs.end();
00119        iter != iterEnd; ++iter) {
00120     if (!(*iter)) continue;
00121     result.push_back(static_cast<OBJECT_TYPE*>((*iter)->Clone(namAd + (*iter)->GetName())));
00122     if (result.back()) result.back()->SetTitle((*iter)->GetTitle() + titAd);
00123     else edm::LogError("Alignment") <<"@SUB=MillePedeMonitor::cloneHists"
00124                                     << "out of memory?";
00125   }
00126   
00127   return result;
00128 }

bool MillePedeMonitor::equidistLogBins ( double *  bins,
int  nBins,
double  first,
double  last 
) const [private]

Definition at line 445 of file MillePedeMonitor.cc.

References i.

Referenced by init().

00447 {
00448   // Filling 'bins' with borders of 'nBins' bins between 'first' and 'last'
00449   // that are equidistant when viewed in log scale,
00450   // so 'bins' must have length nBins+1;
00451   // If 'first', 'last' or 'nBins' are not positive, failure is reported.
00452 
00453   if (nBins < 1 || first <= 0. || last <= 0.) return false;
00454 
00455   bins[0] = first;
00456   bins[nBins] = last;
00457   const double firstLog = TMath::Log10(bins[0]);
00458   const double lastLog  = TMath::Log10(bins[nBins]);
00459   for (int i = 1; i < nBins; ++i) {
00460     bins[i] = TMath::Power(10., firstLog + i*(lastLog-firstLog)/(nBins));
00461   }
00462 
00463   return true;
00464 }

void MillePedeMonitor::fillDerivatives ( const TransientTrackingRecHit::ConstRecHitPointer recHit,
const float *  localDerivs,
unsigned int  nLocal,
const float *  globalDerivs,
unsigned int  nGlobal 
)

Definition at line 682 of file MillePedeMonitor.cc.

References GetIndex(), i, myDerivHists2D, and phi.

Referenced by MillePedeAlignmentAlgorithm::callMille2D().

00685 {
00686   const double phi = recHit->det()->position().phi();
00687 
00688   static const int iLocPar = this->GetIndex(myDerivHists2D, "localDerivsPar");
00689   static const int iLocPhi = this->GetIndex(myDerivHists2D, "localDerivsPhi");
00690   static const int iLocParLog = this->GetIndex(myDerivHists2D, "localDerivsParLog");
00691   static const int iLocPhiLog = this->GetIndex(myDerivHists2D, "localDerivsPhiLog");
00692   for (unsigned int i = 0; i < nLocal; ++i) {
00693     myDerivHists2D[iLocPar]->Fill(i, localDerivs[i]);
00694     myDerivHists2D[iLocPhi]->Fill(phi, localDerivs[i]);
00695     if (localDerivs[i]) {
00696       myDerivHists2D[iLocParLog]->Fill(i, TMath::Abs(localDerivs[i]));
00697       myDerivHists2D[iLocPhiLog]->Fill(phi, TMath::Abs(localDerivs[i]));
00698     }
00699   }
00700 
00701   static const int iGlobPar = this->GetIndex(myDerivHists2D, "globalDerivsPar");
00702   static const int iGlobPhi = this->GetIndex(myDerivHists2D, "globalDerivsPhi");
00703   static const int iGlobParLog = this->GetIndex(myDerivHists2D, "globalDerivsParLog");
00704   static const int iGlobPhiLog = this->GetIndex(myDerivHists2D, "globalDerivsPhiLog");
00705 //   static const int iGlobPhiLog2 = this->GetIndex(myDerivHists2D, "globalDerivsPhiLog2");
00706   for (unsigned int i = 0; i < nGlobal; ++i) {
00707     myDerivHists2D[iGlobPar]->Fill(i, globalDerivs[i]);
00708     myDerivHists2D[iGlobPhi]->Fill(phi, globalDerivs[i]);
00709     if (globalDerivs[i]) {
00710       myDerivHists2D[iGlobParLog]->Fill(i, TMath::Abs(globalDerivs[i]));
00711       myDerivHists2D[iGlobPhiLog]->Fill(phi, TMath::Abs(globalDerivs[i]));
00712 //       myDerivHists2D[iGlobPhiLog2]->Fill(phi, TMath::Abs(globalDerivs[i]));
00713     }
00714   }
00715 
00716 }

void MillePedeMonitor::fillFrameToFrame ( const AlignableDetOrUnitPtr aliDet,
const Alignable ali 
)

Definition at line 867 of file MillePedeMonitor.cc.

References first, FrameToFrameDerivative::frameToFrameDerivative(), GetIndex(), i, j, myFrame2FrameHists2D, phi, and r.

Referenced by MillePedeAlignmentAlgorithm::globalDerivativesHierarchy().

00868 {
00869   // get derivative of higher level structure w.r.t. det
00870   FrameToFrameDerivative ftfd;
00871   const AlgebraicMatrix frameDeriv = ftfd.frameToFrameDerivative(aliDet, ali);
00872   //const AlgebraicMatrix frameDeriv = ftfd.frameToFrameDerivativeAtOrgRot(aliDet, ali);
00873 
00874   static const int iF2f = this->GetIndex(myFrame2FrameHists2D, "frame2frame");
00875   static const int iF2fAbs = this->GetIndex(myFrame2FrameHists2D, "frame2frameAbs");
00876   static int iF2fIjPhi[6][6], iF2fIjPhiLog[6][6], iF2fIjR[6][6], iF2fIjRLog[6][6];
00877   static bool first = true;
00878   if (first) {
00879     for (unsigned int i = 0; i < 6; ++i) {
00880       for (unsigned int j = 0; j < 6; ++j) {
00881         iF2fIjPhi[i][j] = this->GetIndex(myFrame2FrameHists2D, Form("frame2framePhi%d%d", i, j));
00882         iF2fIjPhiLog[i][j]=this->GetIndex(myFrame2FrameHists2D, Form("frame2framePhiLog%d%d",i,j));
00883         iF2fIjR[i][j] = this->GetIndex(myFrame2FrameHists2D, Form("frame2frameR%d%d", i, j));
00884         iF2fIjRLog[i][j]=this->GetIndex(myFrame2FrameHists2D, Form("frame2frameRLog%d%d",i,j));
00885       }
00886     }
00887     first = false;
00888   }
00889   
00890   const double phi = aliDet->globalPosition().phi(); // after misalignment...
00891   const double r = aliDet->globalPosition().perp(); // after misalignment...
00892   for (unsigned int i = 0; i < 6; ++i) {
00893     for (unsigned int j = 0; j < 6; ++j) {
00894       myFrame2FrameHists2D[iF2f]->Fill(i, j, frameDeriv[i][j]);
00895       myFrame2FrameHists2D[iF2fIjPhi[i][j]]->Fill(phi, frameDeriv[i][j]);
00896       myFrame2FrameHists2D[iF2fIjR[i][j]]->Fill(r, frameDeriv[i][j]);
00897       if (frameDeriv[i][j]) {
00898         myFrame2FrameHists2D[iF2fAbs]->Fill(i, j, TMath::Abs(frameDeriv[i][j]));
00899         myFrame2FrameHists2D[iF2fIjPhiLog[i][j]]->Fill(phi, TMath::Abs(frameDeriv[i][j]));
00900         myFrame2FrameHists2D[iF2fIjRLog[i][j]]->Fill(r, TMath::Abs(frameDeriv[i][j]));
00901       }
00902     }
00903   }
00904 }

void MillePedeMonitor::fillRefTrajectory ( const ReferenceTrajectoryBase::ReferenceTrajectoryPtr refTrajPtr  ) 

Definition at line 559 of file MillePedeMonitor.cc.

References GetIndex(), myTrajectoryHists1D, myTrajectoryHists2D, phi, and rho.

Referenced by MillePedeAlignmentAlgorithm::run().

00560 {
00561 
00562   static const int iValid = this->GetIndex(myTrajectoryHists1D, "validRefTraj");
00563   if (refTrajPtr->isValid()) {
00564     myTrajectoryHists1D[iValid]->Fill(1.);
00565   } else {
00566     myTrajectoryHists1D[iValid]->Fill(0.);
00567     return;
00568   }
00569 
00570   const AlgebraicSymMatrix &covMeasLoc = refTrajPtr->measurementErrors();
00571   const AlgebraicVector &measurementsLoc = refTrajPtr->measurements();
00572   const AlgebraicVector &trajectoryLoc = refTrajPtr->trajectoryPositions();
00573   const TransientTrackingRecHit::ConstRecHitContainer &recHits = refTrajPtr->recHits();
00574   const AlgebraicMatrix &derivatives = refTrajPtr->derivatives();
00575 
00576   for (int iRow = 0; iRow < covMeasLoc.num_row(); ++iRow) {
00577     const double residuum = measurementsLoc[iRow] - trajectoryLoc[iRow];
00578     const double covMeasLocRow = covMeasLoc[iRow][iRow];
00579     const bool is2DhitRow = (!recHits[iRow/2]->detUnit() // FIXME: as in MillePedeAlignmentAlgorithm::is2D()
00580                              || recHits[iRow/2]->detUnit()->type().isTrackerPixel());
00581     //the GeomDetUnit* is zero for composite hits (matched hits in the tracker,...). 
00582     if (TMath::Even(iRow)) { // local x
00583       static const int iMeasLocX = this->GetIndex(myTrajectoryHists1D, "measLocX");
00584       myTrajectoryHists1D[iMeasLocX]->Fill(measurementsLoc[iRow]);
00585       static const int iTrajLocX = this->GetIndex(myTrajectoryHists1D, "trajLocX");
00586       myTrajectoryHists1D[iTrajLocX]->Fill(trajectoryLoc[iRow]);
00587       static const int iResidLocX = this->GetIndex(myTrajectoryHists1D, "residLocX");
00588       myTrajectoryHists1D[iResidLocX]->Fill(residuum);
00589 
00590       static const int iMeasLocXvsHit = this->GetIndex(myTrajectoryHists2D, "measLocXvsHit");
00591       myTrajectoryHists2D[iMeasLocXvsHit]->Fill(iRow/2, measurementsLoc[iRow]);
00592       static const int iTrajLocXvsHit = this->GetIndex(myTrajectoryHists2D, "trajLocXvsHit");
00593       myTrajectoryHists2D[iTrajLocXvsHit]->Fill(iRow/2, trajectoryLoc[iRow]);
00594       static const int iResidLocXvsHit = this->GetIndex(myTrajectoryHists2D, "residLocXvsHit");
00595       myTrajectoryHists2D[iResidLocXvsHit]->Fill(iRow/2, residuum);
00596 
00597       if (covMeasLocRow > 0.) { 
00598         static const int iReduResidLocX = this->GetIndex(myTrajectoryHists1D, "reduResidLocX");
00599         myTrajectoryHists1D[iReduResidLocX]->Fill(residuum / TMath::Sqrt(covMeasLocRow));
00600         static const int iReduResidLocXvsHit = 
00601           this->GetIndex(myTrajectoryHists2D, "reduResidLocXvsHit");
00602         myTrajectoryHists2D[iReduResidLocXvsHit]->Fill(iRow/2, residuum/TMath::Sqrt(covMeasLocRow));
00603       }
00604     } else if (is2DhitRow) { // local y, 2D detectors only
00605 
00606       static const int iMeasLocY = this->GetIndex(myTrajectoryHists1D, "measLocY");
00607       myTrajectoryHists1D[iMeasLocY]->Fill(measurementsLoc[iRow]);
00608       static const int iTrajLocY = this->GetIndex(myTrajectoryHists1D, "trajLocY");
00609       myTrajectoryHists1D[iTrajLocY]->Fill(trajectoryLoc[iRow]);
00610       static const int iResidLocY = this->GetIndex(myTrajectoryHists1D, "residLocY");
00611       myTrajectoryHists1D[iResidLocY]->Fill(residuum);
00612 
00613       static const int iMeasLocYvsHit = this->GetIndex(myTrajectoryHists2D, "measLocYvsHit");
00614       myTrajectoryHists2D[iMeasLocYvsHit]->Fill(iRow/2, measurementsLoc[iRow]);
00615       static const int iTrajLocYvsHit = this->GetIndex(myTrajectoryHists2D, "trajLocYvsHit");
00616       myTrajectoryHists2D[iTrajLocYvsHit]->Fill(iRow/2, trajectoryLoc[iRow]);
00617       static const int iResidLocYvsHit = this->GetIndex(myTrajectoryHists2D, "residLocYvsHit");
00618       myTrajectoryHists2D[iResidLocYvsHit]->Fill(iRow/2, residuum);
00619 
00620       if (covMeasLocRow > 0.) {
00621         static const int iReduResidLocY = this->GetIndex(myTrajectoryHists1D, "reduResidLocY");
00622         myTrajectoryHists1D[iReduResidLocY]->Fill(residuum / TMath::Sqrt(covMeasLocRow));
00623         static const int iReduResidLocYvsHit = this->GetIndex(myTrajectoryHists2D,
00624                                                               "reduResidLocYvsHit");
00625         myTrajectoryHists2D[iReduResidLocYvsHit]->Fill(iRow/2, residuum/TMath::Sqrt(covMeasLocRow));
00626       }
00627     }
00628 
00629     float nHitRow = iRow/2; // '/2' not '/2.'!
00630     if (TMath::Odd(iRow)) nHitRow += 0.5; // y-hit gets 0.5
00631     // correlations
00632     for (int iCol = iRow+1; iCol < covMeasLoc.num_col(); ++iCol) {
00633       double rho = TMath::Sqrt(covMeasLocRow + covMeasLoc[iCol][iCol]);
00634       rho = (0. == rho ? -2 : covMeasLoc[iRow][iCol] / rho);
00635       float nHitCol = iCol/2; //cf. comment nHitRow
00636       if (TMath::Odd(iCol)) nHitCol += 0.5; // dito
00637       //      myProfileCorr->Fill(nHitRow, nHitCol, TMath::Abs(rho));
00638       if (0. == rho) continue;
00639       static const int iProfileCorr = this->GetIndex(myTrajectoryHists2D, "profCorr");
00640       myTrajectoryHists2D[iProfileCorr]->Fill(nHitRow, nHitCol, TMath::Abs(rho));
00641       if (iRow+1 == iCol && TMath::Even(iRow)) { // i.e. if iRow is x and iCol the same hit's y
00642 //      static const int iProfileCorrXy = this->GetIndex(myTrajectoryHists2D,"profCorrOffXy");
00643 //      myTrajectoryHists2D[iProfileCorrOffXy]->Fill(iRow/2, rho);
00644         static const int iHitCorrXy = this->GetIndex(myTrajectoryHists2D, "hitCorrXy");
00645         myTrajectoryHists2D[iHitCorrXy]->Fill(iRow/2, rho);
00646         static const int iHitCorrXyLog = this->GetIndex(myTrajectoryHists2D, "hitCorrXyLog");
00647         myTrajectoryHists2D[iHitCorrXyLog]->Fill(iRow/2, TMath::Abs(rho));
00648         if (is2DhitRow) {
00649           static const int iHitCorrXyValid = this->GetIndex(myTrajectoryHists2D, "hitCorrXyValid");
00650           myTrajectoryHists2D[iHitCorrXyValid]->Fill(iRow/2, rho); // nhitRow??
00651           static const int iHitCorrXyLogValid = this->GetIndex(myTrajectoryHists2D,
00652                                                                "hitCorrXyLogValid");
00653           myTrajectoryHists2D[iHitCorrXyLogValid]->Fill(iRow/2, TMath::Abs(rho)); // nhitRow??
00654         }
00655       } else {
00656         static const int iProfCorrOffXy = this->GetIndex(myTrajectoryHists2D, "profCorrOffXy");
00657         myTrajectoryHists2D[iProfCorrOffXy]->Fill(nHitRow, nHitCol, TMath::Abs(rho));
00658         static const int iHitCorrOffBlock = this->GetIndex(myTrajectoryHists2D, "hitCorrOffBlock");
00659         myTrajectoryHists2D[iHitCorrOffBlock]->Fill(nHitRow, rho);
00660         static const int iHitCorOffBlkLg = this->GetIndex(myTrajectoryHists2D,"hitCorrOffBlockLog");
00661         myTrajectoryHists2D[iHitCorOffBlkLg]->Fill(nHitRow, TMath::Abs(rho));
00662       }
00663     } // end loop on columns of covariance
00664     
00665     // derivatives
00666     for (int iCol = 0; iCol < derivatives.num_col(); ++iCol) {
00667       static const int iProfDerivatives = this->GetIndex(myTrajectoryHists2D, "profDerivatives");
00668       myTrajectoryHists2D[iProfDerivatives]->Fill(nHitRow, iCol, derivatives[iRow][iCol]);
00669       static const int iDerivatives = this->GetIndex(myTrajectoryHists2D, "derivatives");
00670       myTrajectoryHists2D[iDerivatives]->Fill(iCol, derivatives[iRow][iCol]);
00671       static const int iDerivativesLog = this->GetIndex(myTrajectoryHists2D, "derivativesLog");
00672       myTrajectoryHists2D[iDerivativesLog]->Fill(iCol, TMath::Abs(derivatives[iRow][iCol]));
00673       static const int iDerivativesPhi = this->GetIndex(myTrajectoryHists2D, "derivativesVsPhi");
00674       myTrajectoryHists2D[iDerivativesPhi]->Fill(recHits[iRow/2]->det()->position().phi(),
00675                                                  derivatives[iRow][iCol]);
00676     }
00677   } // end loop on rows of covarianvce
00678 
00679 }

void MillePedeMonitor::fillResidualHists ( const std::vector< TH1 * > &  hists,
float  phiSensToNorm,
float  residuum,
float  sigma 
) [private]

Definition at line 794 of file MillePedeMonitor.cc.

References GetIndex().

Referenced by fillResiduals().

00796 {
00797   // Histogram indices are calculated at first call, so make sure the vector of hists at the first
00798   // call has the correct names (i.e. without subdet extension) and all later calls have the
00799   // same order of hists...
00800   
00801   static const int iRes = this->GetIndex(hists, "resid");
00802   hists[iRes]->Fill(residuum);
00803   static const int iSigma = this->GetIndex(hists, "sigma");
00804   hists[iSigma]->Fill(sigma);
00805   static const int iSigmaVsAngle = this->GetIndex(hists, "sigmaVsAngle");
00806   hists[iSigmaVsAngle]->Fill(phiSensToNorm, sigma);
00807   static const int iResidVsAngle = this->GetIndex(hists, "residVsAngle");
00808   hists[iResidVsAngle]->Fill(phiSensToNorm, residuum);
00809   static const int iReduRes = this->GetIndex(hists, "reduResid");
00810   static const int iReduResidVsAngle = this->GetIndex(hists, "reduResidVsAngle");
00811   if (sigma) {
00812     hists[iReduRes]->Fill(residuum/sigma);
00813     hists[iReduResidVsAngle]->Fill(phiSensToNorm, residuum/sigma);
00814   }
00815   static const int iAngle = this->GetIndex(hists, "angle");
00816   hists[iAngle]->Fill(phiSensToNorm);
00817   
00818   if (phiSensToNorm > TMath::DegToRad()*45.) {
00819     static const int iResGt45 = this->GetIndex(hists, "residGt45");
00820     hists[iResGt45]->Fill(residuum);
00821     static const int iSigmaGt45 = this->GetIndex(hists, "sigmaGt45");
00822     hists[iSigmaGt45]->Fill(sigma);
00823     static const int iReduResGt45 = this->GetIndex(hists, "reduResidGt45");
00824     if (sigma) hists[iReduResGt45]->Fill(residuum/sigma);
00825   } else {
00826     static const int iResLt45 = this->GetIndex(hists, "residLt45");
00827     hists[iResLt45]->Fill(residuum);
00828     static const int iSigmaLt45 = this->GetIndex(hists, "sigmaLt45");
00829     hists[iSigmaLt45]->Fill(sigma);
00830     static const int iReduResLt45 = this->GetIndex(hists, "reduResidLt45");
00831     if (sigma) hists[iReduResLt45]->Fill(residuum/sigma);
00832   }
00833 }

void MillePedeMonitor::fillResidualHitHists ( const std::vector< TH1 * > &  hists,
float  angle,
float  residuum,
float  sigma,
unsigned int  nHit 
) [private]

Definition at line 836 of file MillePedeMonitor.cc.

References GetIndex(), and i.

Referenced by fillResiduals().

00838 {
00839   // Histogram indices are calculated at first call, so make sure the vector of hists at the first
00840   // call has the correct names (i.e. without subdet extension) and all later calls have the
00841   // same order of hists...
00842   
00843   static const unsigned int maxNhit = 29;  // 0...29 foreseen in initialisation...
00844   static int iResHit[maxNhit+1] = {-1};
00845   static int iSigmaHit[maxNhit+1] = {-1};
00846   static int iReduResHit[maxNhit+1] = {-1};
00847   static int iAngleHit[maxNhit+1] = {-1};
00848   if (iResHit[0] == -1) { // first event only...
00849     for (unsigned int i = 0; i <= maxNhit; ++i) {
00850       iResHit[i] = this->GetIndex(hists, Form("resid_%d", i));
00851       iSigmaHit[i] = this->GetIndex(hists, Form("sigma_%d", i));
00852       iReduResHit[i] = this->GetIndex(hists, Form("reduResid_%d", i));
00853       iAngleHit[i] = this->GetIndex(hists, Form("angle_%d", i));
00854     }
00855   }
00856   if (nHit > maxNhit) nHit = maxNhit; // limit of hists
00857 
00858   hists[iResHit[nHit]]->Fill(residuum);
00859   hists[iSigmaHit[nHit]]->Fill(sigma);
00860   if (sigma) {
00861     hists[iReduResHit[nHit]]->Fill(residuum/sigma);
00862   }
00863   hists[iAngleHit[nHit]]->Fill(angle);
00864 }

void MillePedeMonitor::fillResiduals ( const TransientTrackingRecHit::ConstRecHitPointer recHit,
const TrajectoryStateOnSurface tsos,
unsigned int  nHit,
float  residuum,
float  sigma,
bool  isY 
)

Definition at line 720 of file MillePedeMonitor.cc.

References fillResidualHists(), fillResidualHitHists(), GetIndex(), TrajectoryStateOnSurface::localDirection(), myResidHists2D, myResidHistsVec1DX, myResidHistsVec1DY, myResidHitHists1DX, myResidHitHists1DY, phi, and DetId::Tracker.

Referenced by MillePedeAlignmentAlgorithm::callMille2D().

00723 {
00724   // isY == false: x-measurements
00725   // isY == true:  y-measurements
00726   const GlobalPoint detPos(recHit->det()->position());
00727   const double phi = detPos.phi();
00728 //  const double rSigned = (detPos.y() > 0. ? detPos.perp() : -detPos.perp());
00729 
00730   static const int iResPhi = this->GetIndex(myResidHists2D, "residPhi");
00731   myResidHists2D[iResPhi]->Fill(phi, residuum);
00732   static const int iSigPhi = this->GetIndex(myResidHists2D, "sigmaPhi");
00733   myResidHists2D[iSigPhi]->Fill(phi, sigma);
00734   if (sigma) {
00735     static const int iReduResPhi = this->GetIndex(myResidHists2D, "reduResidPhi");
00736     myResidHists2D[iReduResPhi]->Fill(phi, residuum/sigma);
00737   }
00738 
00739 //   if (isY) {
00740 //     static const int iResYXy = this->GetIndex(myResidHists2D, "residYProfXy");
00741 //     myResidHists2D[iResYXy]->Fill(detPos.x(), detPos.y(), TMath::Abs(residuum));
00742 //     static const int iResYZr = this->GetIndex(myResidHists2D, "residYProfZr");
00743 //     myResidHists2D[iResYZr]->Fill(detPos.z(), rSigned, TMath::Abs(residuum));
00744 //     static const int iSigmaYXy = this->GetIndex(myResidHists2D, "sigmaYProfXy");
00745 //     myResidHists2D[iSigmaYXy]->Fill(detPos.x(), detPos.y(), sigma);
00746 //     static const int iSigmaYZr = this->GetIndex(myResidHists2D, "sigmaYProfZr");
00747 //     myResidHists2D[iSigmaYZr]->Fill(detPos.z(), rSigned, sigma);
00748 //     if (sigma) {
00749 //       static const int iReduResYXy = this->GetIndex(myResidHists2D, "reduResidYProfXy");
00750 //       myResidHists2D[iReduResYXy]->Fill(detPos.x(), detPos.y(), TMath::Abs(residuum/sigma));
00751 //       static const int iReduResYZr = this->GetIndex(myResidHists2D, "reduResidYProfZr");
00752 //       myResidHists2D[iReduResYZr]->Fill(detPos.z(), rSigned, TMath::Abs(residuum/sigma));
00753 //     }
00754 //   } else {
00755 //     static const int iResXXy = this->GetIndex(myResidHists2D, "residXProfXy");
00756 //     myResidHists2D[iResXXy]->Fill(detPos.x(), detPos.y(), TMath::Abs(residuum));
00757 //     static const int iResXZr = this->GetIndex(myResidHists2D, "residXProfZr");
00758 //     myResidHists2D[iResXZr]->Fill(detPos.z(), rSigned, TMath::Abs(residuum));
00759 //     static const int iSigmaXXy = this->GetIndex(myResidHists2D, "sigmaXProfXy");
00760 //     myResidHists2D[iSigmaXXy]->Fill(detPos.x(), detPos.y(), sigma);
00761 //     static const int iSigmaXZr = this->GetIndex(myResidHists2D, "sigmaXProfZr");
00762 //     myResidHists2D[iSigmaXZr]->Fill(detPos.z(), rSigned, sigma);
00763 //     if (sigma) {
00764 //       static const int iReduResXXy = this->GetIndex(myResidHists2D, "reduResidXProfXy");
00765 //       myResidHists2D[iReduResXXy]->Fill(detPos.x(), detPos.y(), TMath::Abs(residuum/sigma));
00766 //       static const int iReduResXZr = this->GetIndex(myResidHists2D, "reduResidXProfZr");
00767 //       myResidHists2D[iReduResXZr]->Fill(detPos.z(), rSigned, TMath::Abs(residuum/sigma));
00768 //     }
00769 //   }
00770 
00771   const LocalVector mom(tsos.localDirection()); // mom.z()==0. impossible for TSOS:
00772   const float phiSensToNorm = TMath::ATan(TMath::Abs((isY ? mom.y() : mom.x())/mom.z()));
00773 
00774   std::vector<std::vector<TH1*> > &histArrayVec = (isY ? myResidHistsVec1DY : myResidHistsVec1DX);
00775   std::vector<TH1*> &hitHists =                   (isY ? myResidHitHists1DY : myResidHitHists1DX);
00776 
00777   // call with histArrayVec[0] first: 'static' inside is referring to "subdet-less" names (X/Y irrelevant)
00778   this->fillResidualHists(histArrayVec[0], phiSensToNorm, residuum, sigma);
00779   this->fillResidualHitHists(hitHists, phiSensToNorm, residuum, sigma, nHit);
00780 
00781   if (recHit->det()->geographicalId().det() == DetId::Tracker) {
00782     //   const GeomDet::SubDetector subDetNum = recHit->det()->subDetector();
00783     unsigned int subDetNum = recHit->det()->geographicalId().subdetId();
00784     if (subDetNum < histArrayVec.size() && subDetNum > 0) {
00785       this->fillResidualHists(histArrayVec[subDetNum], phiSensToNorm, residuum, sigma);
00786     } else {
00787       edm::LogWarning("Alignment") << "@SUB=MillePedeMonitor::fillResiduals"
00788                                    << "Expect subDetNum from 1 to 6, got " << subDetNum;
00789     }
00790   }
00791 }

void MillePedeMonitor::fillTrack ( const reco::Track track,
std::vector< TH1 * > &  trackHists1D,
std::vector< TH2 * > &  trackHists2D 
) [private]

Definition at line 488 of file MillePedeMonitor.cc.

References reco::TrackBase::d0(), reco::TrackBase::d0Error(), reco::TrackBase::dz(), reco::TrackBase::dzError(), GetIndex(), reco::Track::innerOk(), reco::Track::innerPosition(), reco::TrackBase::momentum(), reco::TrackBase::normalizedChi2(), reco::TrackBase::numberOfLostHits(), reco::TrackBase::numberOfValidHits(), reco::Track::outerOk(), reco::Track::outerPosition(), p, and reco::TrackBase::pt().

00490 {
00491   if (!track) return;
00492 
00493   const reco::TrackBase::Vector p(track->momentum());
00494 
00495   static const int iPtLog = this->GetIndex(trackHists1D, "ptTrackLogBins");
00496   trackHists1D[iPtLog]->Fill(track->pt());
00497   static const int iPt = this->GetIndex(trackHists1D, "ptTrack");
00498   trackHists1D[iPt]->Fill(track->pt());
00499   static const int iP = this->GetIndex(trackHists1D, "pTrack");
00500   trackHists1D[iP]->Fill(p.R());
00501   static const int iEta = this->GetIndex(trackHists1D, "etaTrack");
00502   trackHists1D[iEta]->Fill(p.Eta());
00503   static const int iPhi = this->GetIndex(trackHists1D, "phiTrack");
00504   trackHists1D[iPhi]->Fill(p.Phi());
00505 
00506   static const int iNhit = this->GetIndex(trackHists1D, "nHitTrack");
00507   trackHists1D[iNhit]->Fill(track->numberOfValidHits());
00508   static const int iNhitInvalid = this->GetIndex(trackHists1D, "nHitInvalidTrack");
00509   trackHists1D[iNhitInvalid]->Fill(track->numberOfLostHits());
00510   if (track->innerOk()) {
00511     const reco::TrackBase::Point firstPoint(track->innerPosition());
00512     static const int iR1 = this->GetIndex(trackHists1D, "r1Track");
00513     trackHists1D[iR1]->Fill(firstPoint.Rho());
00514     const double rSigned1 = (firstPoint.y() > 0 ? firstPoint.Rho() : -firstPoint.Rho());
00515     static const int iR1Signed = this->GetIndex(trackHists1D, "r1TrackSigned");
00516     trackHists1D[iR1Signed]->Fill(rSigned1);
00517     static const int iZ1 = this->GetIndex(trackHists1D, "z1Track");
00518     trackHists1D[iZ1]->Fill(firstPoint.Z());
00519     static const int iZ1Full = this->GetIndex(trackHists1D, "z1TrackFull");
00520     trackHists1D[iZ1Full]->Fill(firstPoint.Z());
00521     static const int iY1 = this->GetIndex(trackHists1D, "y1Track");
00522     trackHists1D[iY1]->Fill(firstPoint.Y());
00523     static const int iPhi1 = this->GetIndex(trackHists1D, "phi1Track");
00524     trackHists1D[iPhi1]->Fill(firstPoint.phi());
00525     static const int iRz1Full = this->GetIndex(trackHists2D, "rz1TrackFull");
00526     trackHists2D[iRz1Full]->Fill(firstPoint.Z(), rSigned1);
00527     static const int iXy1 = this->GetIndex(trackHists2D, "xy1Track");
00528     trackHists2D[iXy1]->Fill(firstPoint.X(), firstPoint.Y());
00529   }
00530   if (track->outerOk()) {
00531     const reco::TrackBase::Point lastPoint(track->outerPosition());
00532     static const int iRlast = this->GetIndex(trackHists1D, "rLastTrack");
00533     trackHists1D[iRlast]->Fill(lastPoint.Rho());
00534     static const int iZlast = this->GetIndex(trackHists1D, "zLastTrack");
00535     trackHists1D[iZlast]->Fill(lastPoint.Z());
00536     static const int iYlast = this->GetIndex(trackHists1D, "yLastTrack");
00537     trackHists1D[iYlast]->Fill(lastPoint.Y());
00538     static const int iPhiLast = this->GetIndex(trackHists1D, "phiLastTrack");
00539     trackHists1D[iPhiLast]->Fill(lastPoint.phi());
00540   }
00541 
00542   static const int iChi2Ndf = this->GetIndex(trackHists1D, "chi2PerNdf");
00543   trackHists1D[iChi2Ndf]->Fill(track->normalizedChi2());
00544 
00545   static const int iImpZ = this->GetIndex(trackHists1D, "impParZ");
00546   trackHists1D[iImpZ]->Fill(track->dz());
00547   static const int iImpZerr = this->GetIndex(trackHists1D, "impParErrZ");
00548   trackHists1D[iImpZerr]->Fill(track->dzError());
00549 
00550 
00551   static const int iImpRphi = this->GetIndex(trackHists1D, "impParRphi");
00552   trackHists1D[iImpRphi]->Fill(track->d0());
00553   static const int iImpRphiErr = this->GetIndex(trackHists1D, "impParErrRphi");
00554   trackHists1D[iImpRphiErr]->Fill(track->d0Error());
00555 
00556 }

void MillePedeMonitor::fillTrack ( const reco::Track track  ) 

Definition at line 467 of file MillePedeMonitor.cc.

References myTrackHists1D, and myTrackHists2D.

Referenced by fillUsedTrack(), and MillePedeAlignmentAlgorithm::run().

00468 {
00469   this->fillTrack(track, myTrackHists1D, myTrackHists2D);
00470 }

void MillePedeMonitor::fillUsedTrack ( const reco::Track track,
unsigned int  nHitX,
unsigned int  nHitY 
)

Definition at line 473 of file MillePedeMonitor.cc.

References fillTrack(), GetIndex(), myUsedTrackHists1D, and myUsedTrackHists2D.

Referenced by MillePedeAlignmentAlgorithm::run().

00475 {
00476   if (!track) return;
00477   this->fillTrack(track, myUsedTrackHists1D, myUsedTrackHists2D);
00478 
00479   // these hist exist only for 'used track' hists:
00480   static const int iUsedX = this->GetIndex(myUsedTrackHists1D, "usedHitsX");
00481   myUsedTrackHists1D[iUsedX]->Fill(nHitX);
00482   static const int iUsedY = this->GetIndex(myUsedTrackHists1D, "usedHitsY");
00483   myUsedTrackHists1D[iUsedY]->Fill(nHitY);
00484 
00485 }

template<class OBJECT_TYPE>
int MillePedeMonitor::GetIndex ( const std::vector< OBJECT_TYPE * > &  vec,
const TString &  name 
) [inline, private]

Definition at line 100 of file MillePedeMonitor.h.

References iter, and HLT_VtxMuL3::result.

Referenced by fillDerivatives(), fillFrameToFrame(), fillRefTrajectory(), fillResidualHists(), fillResidualHitHists(), fillResiduals(), fillTrack(), and fillUsedTrack().

00101 {
00102   int result = 0;
00103   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
00104        iter != iterEnd; ++iter, ++result) {
00105     if (*iter && (*iter)->GetName() == name) return result;
00106   }
00107   edm::LogError("Alignment") << "@SUB=MillePedeMonitor::GetIndex" << " could not find " << name;
00108   return -1;
00109 }

bool MillePedeMonitor::init ( TDirectory *  directory  )  [private]

Definition at line 64 of file MillePedeMonitor.cc.

References addToDirectory(), cloneHists(), equidistLogBins(), h, i, j, myDerivHists2D, myFrame2FrameHists2D, myResidHists2D, myResidHistsVec1DX, myResidHistsVec1DY, myResidHitHists1DX, myResidHitHists1DY, myTrackHists1D, myTrackHists2D, myTrajectoryHists1D, myTrajectoryHists2D, myUsedTrackHists1D, myUsedTrackHists2D, and Pi.

Referenced by MillePedeMonitor().

00065 {
00066   if (!directory) return false;
00067   TDirectory *oldDir = gDirectory;
00068 
00069   const int kNumBins = 20;
00070   double binsPt[kNumBins+1] = {0.}; // fully initialised with 0.
00071   if (!this->equidistLogBins(binsPt, kNumBins, 0.8, 100.)) {
00072 //     cerr << "MillePedeMonitor::init: problem with log bins" << endl;
00073   }
00074   const int nHits = 25;
00075 
00076   myTrackHists1D.push_back(new TH1F("ptTrackLogBins",  "p_{t}(track);p_{t} [GeV]",
00077                                     kNumBins, binsPt));
00078 
00079   myTrackHists1D.push_back(new TH1F("ptTrack",  "p_{t}(track);p_{t} [GeV]",
00080                                     kNumBins, binsPt[0], binsPt[kNumBins]));
00081   myTrackHists1D.push_back(new TH1F("pTrack",  "p(track);p [GeV]",
00082                                     kNumBins, binsPt[0], 1.3*binsPt[kNumBins]));
00083   myTrackHists1D.push_back(new TH1F("etaTrack", "#eta(track);#eta", 26, -2.6, 2.6));
00084   myTrackHists1D.push_back(new TH1F("phiTrack", "#phi(track);#phi", 15, -TMath::Pi(), TMath::Pi()));
00085 
00086   myTrackHists1D.push_back(new TH1F("nHitTrack", "N_{hit}(track);N_{hit}", 30, 5., 35.));
00087   myTrackHists1D.push_back(new TH1F("nHitInvalidTrack", "N_{hit, invalid}(track)", 5, 0., 5.));
00088   myTrackHists1D.push_back(new TH1F("r1Track", "r(1st hit);r [cm]", 20, 0., 20.));
00089   myTrackHists1D.push_back(new TH1F("phi1Track", "#phi(1st hit);#phi",
00090                                     30, -TMath::Pi(), TMath::Pi()));
00091   myTrackHists1D.push_back(new TH1F("y1Track", "y(1st hit);y [cm]", 40, -120., +120.));
00092   myTrackHists1D.push_back(new TH1F("z1Track", "z(1st hit);z [cm]", 20, -50, +50));
00093   myTrackHists1D.push_back(new TH1F("r1TrackSigned", "r(1st hit);r_{#pm} [cm]",
00094                                     40, -120., 120.));
00095   myTrackHists1D.push_back(new TH1F("z1TrackFull", "z(1st hit);z [cm]", 20, -300., +300.));
00096   myTrackHists1D.push_back(new TH1F("rLastTrack", "r(last hit);r [cm]", 20, 20., 120.));
00097   myTrackHists1D.push_back(new TH1F("phiLastTrack", "#phi(last hit);#phi",
00098                                     30, -TMath::Pi(), TMath::Pi()));
00099   myTrackHists1D.push_back(new TH1F("yLastTrack", "y(last hit);y [cm]", 40, -120., +120.));
00100   myTrackHists1D.push_back(new TH1F("zLastTrack", "z(last hit);z [cm]", 30, -300., +300.));
00101   myTrackHists1D.push_back(new TH1F("chi2PerNdf", "#chi^{2}/ndf;#chi^{2}/ndf", 500, 0., 50.));
00102   myTrackHists1D.push_back(new TH1F("impParZ", "impact parameter in z", 20, -20., 20.));
00103   myTrackHists1D.push_back(new TH1F("impParErrZ", "error of impact parameter in z",
00104                                     40, 0., 0.06));  
00105   myTrackHists1D.push_back(new TH1F("impParRphi", "impact parameter in r#phi", 51, -0.05, .05));
00106   myTrackHists1D.push_back(new TH1F("impParErrRphi", "error of impact parameter in r#phi",
00107                                     50, 0., 0.01));
00108 
00109   myTrackHists2D.push_back(new TH2F("rz1TrackFull", "rz(1st hit);z [cm]; r_{#pm} [cm]",
00110                                     40, -282., +282., 40, -115., +115.));
00111   myTrackHists2D.push_back(new TH2F("xy1Track", "xy(1st hit);x [cm]; y [cm]",
00112                                     40, -115., +115., 40, -115., +115.));
00113   
00114   TDirectory *dirTracks = directory->mkdir("trackHists", "input tracks");
00115   this->addToDirectory(myTrackHists1D, dirTracks);
00116   this->addToDirectory(myTrackHists2D, dirTracks);
00117 
00118 // used track 
00119   myUsedTrackHists1D = this->cloneHists(myTrackHists1D, "used", " (used tracks)");
00120   myUsedTrackHists2D = this->cloneHists(myTrackHists2D, "used", " (used tracks)");
00121   // must be after clone: index in vector!
00122   myUsedTrackHists1D.push_back(new TH1F("usedHitsX", "n(x-hits) transferred to pede;n(x-hits)", nHits, 0., nHits));
00123   myUsedTrackHists1D.push_back(new TH1F("usedHitsY", "n(y-hits) transferred to pede;n(y-hits)", nHits-10, 0., nHits-10));
00124 
00125   TDirectory *dirUsedTracks = directory->mkdir("usedTrackHists", "used tracks");
00126   this->addToDirectory(myUsedTrackHists1D, dirUsedTracks);
00127   this->addToDirectory(myUsedTrackHists2D, dirUsedTracks);
00128 
00129 // ReferenceTrajectory
00130   myTrajectoryHists1D.push_back(new TH1F("validRefTraj", "validity of ReferenceTrajectory",
00131                                          2, 0., 2.));
00132 
00133   myTrajectoryHists2D.push_back(new TProfile2D("profCorr",
00134                                                "mean of |#rho|, #rho#neq0;hit x/y;hit x/y;",
00135                                                2*nHits, 0., nHits, 2*nHits, 0., nHits));
00136   myTrajectoryHists2D.push_back
00137     (new TProfile2D("profCorrOffXy", "mean of |#rho|, #rho#neq0, no xy_{hit};hit x/y;hit x/y;",
00138                     2*nHits, 0., nHits, 2*nHits, 0., nHits));
00139 
00140   myTrajectoryHists2D.push_back(new TH2F("hitCorrOffBlock",
00141                                          "hit correlations: off-block-diagonals;N(hit);#rho",
00142                                          2*nHits, 0., nHits, 81, -.06, .06));
00143   TArrayD logBins(102);
00144   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-11, .1);
00145   myTrajectoryHists2D.push_back(new TH2F("hitCorrOffBlockLog",
00146                                          "hit correlations: off-block-diagonals;N(hit);|#rho|",
00147                                          2*nHits, 0., nHits, logBins.GetSize()-1, logBins.GetArray()));
00148 
00149   myTrajectoryHists2D.push_back(new TH2F("hitCorrXy", "hit correlations: xy;N(hit);#rho",
00150                                          nHits, 0., nHits, 81, -.5, .5));
00151   myTrajectoryHists2D.push_back
00152     (new TH2F("hitCorrXyValid", "hit correlations: xy, 2D-det.;N(hit);#rho",
00153               nHits, 0., nHits, 81, -.02, .02));
00154   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-10, .5);
00155   myTrajectoryHists2D.push_back(new TH2F("hitCorrXyLog", "hit correlations: xy;N(hit);|#rho|",
00156                                          nHits, 0., nHits, logBins.GetSize()-1, logBins.GetArray()));
00157   myTrajectoryHists2D.push_back
00158     (new TH2F("hitCorrXyLogValid", "hit correlations: xy, 2D-det.;N(hit);|#rho|",
00159               nHits, 0., nHits, logBins.GetSize()-1, logBins.GetArray()));
00160 
00161 
00162   myTrajectoryHists1D.push_back(new TH1F("measLocX", "local x measurements;x", 101, -6., 6.));
00163   myTrajectoryHists1D.push_back(new TH1F("measLocY", "local y measurements, 2D-det.;y",
00164                                          101, -10., 10.));
00165   myTrajectoryHists1D.push_back(new TH1F("trajLocX", "local x trajectory;x", 101, -6., 6.));
00166   myTrajectoryHists1D.push_back(new TH1F("trajLocY", "local y trajectory, 2D-det.;y",
00167                                          101, -10., 10.));
00168 
00169   myTrajectoryHists1D.push_back(new TH1F("residLocX", "local x residual;#Deltax", 101, -.75, .75));
00170   myTrajectoryHists1D.push_back(new TH1F("residLocY", "local y residual, 2D-det.;#Deltay",
00171                                          101, -2., 2.));
00172   myTrajectoryHists1D.push_back(new TH1F("reduResidLocX", "local x reduced residual;#Deltax/#sigma",
00173                                          101, -20., 20.));
00174   myTrajectoryHists1D.push_back
00175     (new TH1F("reduResidLocY", "local y reduced residual, 2D-det.;#Deltay/#sigma", 101, -20., 20.));
00176 
00177   // 2D vs. hit
00178   myTrajectoryHists2D.push_back(new TH2F("measLocXvsHit", "local x measurements;hit;x", 
00179                                          nHits, 0., nHits, 101, -6., 6.));
00180   myTrajectoryHists2D.push_back(new TH2F("measLocYvsHit", "local y measurements, 2D-det.;hit;y",
00181                                          nHits, 0., nHits, 101, -10., 10.));
00182   myTrajectoryHists2D.push_back(new TH2F("trajLocXvsHit", "local x trajectory;hit;x",
00183                                          nHits, 0., nHits, 101, -6., 6.));
00184   myTrajectoryHists2D.push_back(new TH2F("trajLocYvsHit", "local y trajectory, 2D-det.;hit;y",
00185                                          nHits, 0., nHits,  101, -10., 10.));
00186 
00187   myTrajectoryHists2D.push_back(new TH2F("residLocXvsHit", "local x residual;hit;#Deltax",
00188                                          nHits, 0., nHits, 101, -.75, .75));
00189   myTrajectoryHists2D.push_back(new TH2F("residLocYvsHit", "local y residual, 2D-det.;hit;#Deltay",
00190                                          nHits, 0., nHits, 101, -1., 1.));
00191   myTrajectoryHists2D.push_back
00192     (new TH2F("reduResidLocXvsHit", "local x reduced residual;hit;#Deltax/#sigma",
00193               nHits, 0., nHits, 101, -20., 20.));
00194   myTrajectoryHists2D.push_back
00195     (new TH2F("reduResidLocYvsHit", "local y reduced residual, 2D-det.;hit;#Deltay/#sigma",
00196               nHits, 0., nHits, 101, -20., 20.));
00197 
00198 
00199   myTrajectoryHists2D.push_back(new TProfile2D("profDerivatives",
00200                                                "mean derivatives;hit x/y;parameter;",
00201                                                2*nHits, 0., nHits, 10, 0., 10.));
00202 
00203   myTrajectoryHists2D.push_back
00204     (new TH2F("derivatives", "derivative;parameter;#partial(x/y)_{local}/#partial(param)",
00205               10, 0., 10., 101, -20., 20.));
00206   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-12, 100.);
00207   myTrajectoryHists2D.push_back
00208     (new TH2F("derivativesLog", "|derivative|;parameter;|#partial(x/y)_{local}/#partial(param)|",
00209               10, 0., 10., logBins.GetSize()-1, logBins.GetArray()));
00210   myTrajectoryHists2D.push_back
00211     (new TH2F("derivativesVsPhi", 
00212               "derivatives vs. #phi;#phi(geomDet);#partial(x/y)_{local}/#partial(param)",
00213               50, -TMath::Pi(), TMath::Pi(), 101, -300., 300.));
00214   //  myTrajectoryHists2D.back()->SetBit(TH1::kCanRebin);
00215 
00216   TDirectory *dirTraject = directory->mkdir("refTrajectoryHists", "ReferenceTrajectory's");
00217   this->addToDirectory(myTrajectoryHists2D, dirTraject);
00218   this->addToDirectory(myTrajectoryHists1D, dirTraject);
00219 
00220 // derivatives hists
00221   myDerivHists2D.push_back
00222     (new TH2F("localDerivsPar","local derivatives vs. paramater;parameter;#partial/#partial(param)",
00223               6, 0., 6., 101, -200., 200.));
00224   myDerivHists2D.push_back
00225     (new TH2F("localDerivsPhi","local derivatives vs. #phi(det);#phi(det);#partial/#partial(param)",
00226               51, -TMath::Pi(), TMath::Pi(), 101, -150., 150.));
00227   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-13, 150.);
00228   myDerivHists2D.push_back
00229     (new TH2F("localDerivsParLog",
00230               "local derivatives (#neq 0) vs. parameter;parameter;|#partial/#partial(param)|",
00231               6, 0., 6., logBins.GetSize()-1, logBins.GetArray()));
00232   myDerivHists2D.push_back
00233     (new TH2F("localDerivsPhiLog",
00234               "local derivatives (#neq 0) vs. #phi(det);#phi(det);|#partial/#partial(param)|",
00235               51, -TMath::Pi(), TMath::Pi(), logBins.GetSize()-1, logBins.GetArray()));
00236   myDerivHists2D.push_back
00237     (new TH2F("globalDerivsPar",
00238               "global derivatives vs. paramater;parameter;#partial/#partial(param)",
00239               6, 0., 6., 101, -.01, .01));
00240   myDerivHists2D.push_back
00241     (new TH2F("globalDerivsPhi",
00242               "global derivatives vs. #phi(det);#phi(det);#partial/#partial(param)",
00243               51, -TMath::Pi(), TMath::Pi(), 101, -0.01, 0.01));
00244   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-36, 0.04);
00245   myDerivHists2D.push_back
00246     (new TH2F("globalDerivsParLog",
00247               "global derivatives (#neq 0) vs. parameter;parameter;|#partial/#partial(param)|",
00248               6, 0., 6., logBins.GetSize()-1, logBins.GetArray()));
00249   myDerivHists2D.push_back
00250     (new TH2F("globalDerivsPhiLog",
00251               "global derivatives (#neq 0) vs. #phi(det);#phi(det);|#partial/#partial(param)|",
00252               51, -TMath::Pi(), TMath::Pi(), logBins.GetSize()-1, logBins.GetArray()));
00253 //   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.e-40, 1.E-35);
00254 //   myDerivHists2D.push_back
00255 //     (new TH2F("globalDerivsPhiLog2",
00256 //               "global derivatives (#neq 0) vs. #phi(det);#phi(det);|#partial/#partial(param)|",
00257 //               51, -TMath::Pi(), TMath::Pi(), logBins.GetSize()-1, logBins.GetArray()));
00258 
00259   TDirectory *dirDerivs = directory->mkdir("derivatives", "derivatives etc.");
00260   this->addToDirectory(myDerivHists2D, dirDerivs);
00261 
00262 // residual hists
00263   myResidHists2D.push_back(new TH2F("residPhi","residuum vs. #phi(det);#phi(det);residuum[cm]",
00264                                     51, -TMath::Pi(), TMath::Pi(), 101, -.5, .5));
00265   myResidHists2D.push_back(new TH2F("sigmaPhi","#sigma vs. #phi(det);#phi(det);#sigma[cm]",
00266                                     51, -TMath::Pi(), TMath::Pi(), 101, .0, .1));
00267   myResidHists2D.push_back(new TH2F("reduResidPhi",
00268                                     "residuum/#sigma vs. #phi(det);#phi(det);residuum/#sigma",
00269                                     51, -TMath::Pi(), TMath::Pi(), 101, -14., 14.));
00270   
00271 //   myResidHists2D.push_back(new TProfile2D("residXProfXy",
00272 //                                        "mean |residuum| (u);x [cm];y [cm];#LT|residuum|#GT",
00273 //                                        25, -110., 110., 25, -110., 110.));
00274 //   myResidHists2D.push_back(new TProfile2D("residXProfZr",
00275 //                                        "mean |residuum| (u);z [cm];r_{#pm} [cm];#LT|residuum|#GT",
00276 //                                        25, -275., 275., 25, -110., 110.));
00277 //   myResidHists2D.push_back(new TProfile2D("residYProfXy",
00278 //                                        "mean |residuum| (v);x [cm];y [cm];#LT|residuum|#GT",
00279 //                                        25, -110., 110., 25, -110., 110.));
00280 //   myResidHists2D.push_back(new TProfile2D("residYProfZr",
00281 //                                        "mean |residuum| (v);z [cm];r_{#pm} [cm];#LT|residuum|#GT",
00282 //                                        25, -275., 275., 25, -110., 110.));
00283 //   myResidHists2D.push_back(new TProfile2D("reduResidXProfXy",
00284 //                                        "mean |residuum/#sigma| (u);x [cm];y [cm];#LT|res./#sigma|#GT",
00285 //                                        25, -110., 110., 25, -110., 110.));
00286 //   myResidHists2D.push_back(new TProfile2D("reduResidXProfZr",
00287 //                                        "mean |residuum/#sigma| (u);z [cm];r_{#pm} [cm];#LT|res./#sigma|#GT",
00288 //                                        25, -275., 275., 25, -110., 110.));
00289 //   myResidHists2D.push_back(new TProfile2D("reduResidYProfXy",
00290 //                                        "mean |residuum/#sigma| (v);x [cm];y [cm];#LT|res./#sigma|#GT",
00291 //                                        25, -110., 110., 25, -110., 110.));
00292 //   myResidHists2D.push_back(new TProfile2D("reduResidYProfZr",
00293 //                                        "mean |residuum/#sigma| (v);z [cm];r_{#pm} [cm];#LT|res./#sigma|#GT",
00294 //                                        25, -275., 275., 25, -110., 110.));
00295 //   myResidHists2D.push_back(new TProfile2D("sigmaXProfXy",
00296 //                                        "mean sigma (u);x [cm];y [cm];#LT#sigma#GT",
00297 //                                        25, -110., 110., 25, -110., 110.));
00298 //   myResidHists2D.push_back(new TProfile2D("sigmaXProfZr",
00299 //                                        "mean sigma (u);z [cm];r_{#pm} [cm];#LT#sigma#GT",
00300 //                                        25, -275., 275., 25, -110., 110.));
00301 //   myResidHists2D.push_back(new TProfile2D("sigmaYProfXy",
00302 //                                        "mean sigma (v);x [cm];y [cm];#LT#sigma#GT",
00303 //                                        25, -110., 110., 25, -110., 110.));
00304 //   myResidHists2D.push_back(new TProfile2D("sigmaYProfZr",
00305 //                                        "mean sigma (v);z [cm];r_{#pm} [cm];#LT#sigma#GT",
00306 //                                        25, -275., 275., 25, -110., 110.));
00307   
00308   TDirectory *dirResid = directory->mkdir("residuals", "hit residuals, sigma,...");
00309   this->addToDirectory(myResidHists2D, dirResid);
00310 
00311 // Residuum, hit sigma and res./sigma for all sensor/track angles and separated for large and
00312 // small angles with respect to the sensor normal in sensitive direction. 
00313   // Here for x-measurements:
00314   std::vector<TH1*> allResidHistsX;
00315   allResidHistsX.push_back(new TH1F("resid", "hit residuals;residuum [cm]", 51, -.05, .05));
00316   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00317   allResidHistsX.push_back(new TH1F("sigma", "hit uncertainties;#sigma [cm]", 50, 0., .02));
00318   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00319   allResidHistsX.push_back(new TH1F("reduResid", "reduced hit residuals;res./#sigma",
00320                                     51, -3., 3.));
00321   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00322   allResidHistsX.push_back(new TH1F("angle", "#phi_{tr} wrt normal (sens. plane);#phi_{n}^{sens}",
00323                                     50, 0., TMath::PiOver2()));
00324   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00325   allResidHistsX.push_back(new TH2F("residVsAngle",
00326                                     "residuum vs. #phi_{tr} wrt normal (sens. plane);#phi_{n}^{sens};residuum [cm]",
00327                                     50, 0., TMath::PiOver2(), 51, -1., 1.));
00328   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-6, 100.);
00329   allResidHistsX.push_back(new TH2F("sigmaVsAngle",
00330                                     "#sigma vs. #phi_{tr} wrt normal (sens. plane);#phi_{n}^{sens};#sigma [cm]",
00331                                     50, 0., TMath::PiOver2(), logBins.GetSize()-1, logBins.GetArray()));
00332   allResidHistsX.push_back(new TH2F("reduResidVsAngle",
00333                                     "reduced residuum vs. #phi_{tr} wrt normal (sens. plane);#phi_{n}^{sens};res./#sigma",
00334                                     50, 0., TMath::PiOver2(), 51, -15., 15.));
00335 
00336   allResidHistsX.push_back(new TH1F("residGt45",
00337                                     "hit residuals (#phi_{n}^{sens}>45#circ);residuum [cm]",
00338                                     51, -.05, .05));
00339   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00340   allResidHistsX.push_back(new TH1F("sigmaGt45",
00341                                     "hit uncertainties(#phi_{n}^{sens}>45#circ);#sigma [cm]",
00342                                     50, 0., .02));
00343   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00344   allResidHistsX.push_back(new TH1F("reduResidGt45",
00345                                     "reduced hit residuals(#phi_{n}^{sens}>45#circ);res./#sigma",
00346                                     51,-3.,3.));
00347   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00348   allResidHistsX.push_back(new TH1F("residLt45",
00349                                     "hit residuals (#phi_{n}^{sens}<45#circ);residuum [cm]",
00350                                     51, -.15, .15));
00351   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00352   allResidHistsX.push_back(new TH1F("sigmaLt45",
00353                                     "hit uncertainties(#phi_{n}^{sens}<45#circ);#sigma [cm]",
00354                                     50, 0., .01));
00355   allResidHistsX.back()->SetBit(TH1::kCanRebin);
00356   allResidHistsX.push_back(new TH1F("reduResidLt45",
00357                                     "reduced hit residuals(#phi_{n}^{sens}<45#circ);res./#sigma",
00358                                     51,-3.,3.));
00359   myResidHistsVec1DX.push_back(allResidHistsX); // at [0] for all subdets together...
00360   // ... then separately - indices/order like DetId.subdetId() in tracker:
00361   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TPB", " x-coord. in pixel barrel"));
00362   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TPE", " x-coord. in pixel discs"));
00363   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TIB", " x-coord. in TIB"));
00364   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TID", " x-coord. in TID"));
00365   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TOB", " x-coord. in TOB"));
00366   myResidHistsVec1DX.push_back(this->cloneHists(allResidHistsX, "TEC", " x-coord. in TEC"));
00367   // finally, differential in hit number (but subdet independent)
00368   for (unsigned int iHit = 0; iHit < 30; ++iHit) { // 4: for each hit fill only angle independent plots
00369     for (unsigned int iHist = 0; iHist < 4 && iHist < allResidHistsX.size(); ++iHist) {
00370       TH1 *h = allResidHistsX[iHist];
00371       myResidHitHists1DX.push_back(static_cast<TH1*>(h->Clone(Form("%s_%d", h->GetName(), iHit))));
00372       myResidHitHists1DX.back()->SetTitle(Form("%s, hit %d", h->GetTitle(), iHit));
00373     }
00374   }
00375 
00376   TDirectory *dirResidX =
00377     (dirResid ? dirResid : directory)->mkdir("X", "hit residuals etc. for x measurements");
00378   this->addToDirectory(myResidHitHists1DX, dirResidX);
00379   for (std::vector<std::vector<TH1*> >::iterator vecIter = myResidHistsVec1DX.begin(),
00380          vecIterEnd = myResidHistsVec1DX.end(); vecIter != vecIterEnd; ++vecIter) {
00381     this->addToDirectory(*vecIter, dirResidX);
00382   }
00383 
00384   // Now clone the same as above for y-ccordinate:
00385   myResidHistsVec1DY.push_back(this->cloneHists(allResidHistsX, "", " y-coord."));// all subdets
00386   std::vector<TH1*> &yHists1D = allResidHistsX;//myResidHistsVec1DY.back(); crashes? ROOT?
00387   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TPB", " y-coord. in pixel barrel"));
00388   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TPE", " y-coord. in pixel discs"));
00389   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TIB", " y-coord. in TIB"));
00390   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TID", " y-coord. in TID"));
00391   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TOB", " y-coord. in TOB"));
00392   myResidHistsVec1DY.push_back(this->cloneHists(yHists1D, "TEC", " y-coord. in TEC"));
00393   myResidHitHists1DY = this->cloneHists(myResidHitHists1DX, "", " y-coord.");// diff. in nHit
00394 
00395   TDirectory *dirResidY = 
00396     (dirResid ? dirResid : directory)->mkdir("Y", "hit residuals etc. for y measurements");
00397   this->addToDirectory(myResidHitHists1DY, dirResidY);
00398   for (std::vector<std::vector<TH1*> >::iterator vecIter = myResidHistsVec1DY.begin(),
00399          vecIterEnd = myResidHistsVec1DY.end(); vecIter != vecIterEnd; ++vecIter) {
00400     this->addToDirectory(*vecIter, dirResidY);
00401   }
00402 
00403   // farme-to-frame derivatives
00404   myFrame2FrameHists2D.push_back(new TProfile2D("frame2frame",
00405                                                 "mean frame to frame derivatives;col;row",
00406                                                 6, 0., 6., 6, 0., 6.));
00407   myFrame2FrameHists2D.push_back(new TProfile2D("frame2frameAbs",
00408                                                 "mean |frame to frame derivatives|, #neq0;col;row",
00409                                                 6, 0., 6., 6, 0., 6.));
00410 
00411   this->equidistLogBins(logBins.GetArray(), logBins.GetSize()-1, 1.E-36, 100.);
00412   for (unsigned int i = 0; i < 6; ++i) {
00413     for (unsigned int j = 0; j < 6; ++j) {
00414       myFrame2FrameHists2D.push_back
00415         (new TH2F(Form("frame2framePhi%d%d", i, j),
00416                   Form("frame to frame derivatives, %d%d;#phi(aliDet);deriv",i,j),
00417                   51, -TMath::Pi(), TMath::Pi(), 10, 0., 1.));
00418       myFrame2FrameHists2D.back()->SetBit(TH1::kCanRebin);
00419       myFrame2FrameHists2D.push_back
00420         (new TH2F(Form("frame2frameR%d%d", i, j),
00421                   Form("frame to frame derivatives, %d%d;r(aliDet);deriv",i,j),
00422                   51, 0., 110., 10, 0., 1.));
00423       myFrame2FrameHists2D.back()->SetBit(TH1::kCanRebin);
00424 
00425       myFrame2FrameHists2D.push_back
00426         (new TH2F(Form("frame2framePhiLog%d%d", i, j),
00427                   Form("frame to frame |derivatives|, %d%d, #neq0;#phi(aliDet);deriv",i,j),
00428                   51, -TMath::Pi(), TMath::Pi(), logBins.GetSize()-1, logBins.GetArray()));
00429       myFrame2FrameHists2D.push_back
00430         (new TH2F(Form("frame2frameRLog%d%d", i, j),
00431                   Form("frame to frame |derivatives|, %d%d, #neq0;r(aliDet);deriv",i,j),
00432                   51, 0., 110., logBins.GetSize()-1, logBins.GetArray()));
00433     }
00434   }
00435 
00436   TDirectory *dirF2f = directory->mkdir("frame2FrameHists", "derivatives etc.");
00437   this->addToDirectory(myFrame2FrameHists2D, dirF2f);
00438 
00439   oldDir->cd();
00440   return true;
00441 }


Member Data Documentation

bool MillePedeMonitor::myDeleteDir [private]

Definition at line 81 of file MillePedeMonitor.h.

Referenced by MillePedeMonitor(), and ~MillePedeMonitor().

std::vector<TH2*> MillePedeMonitor::myDerivHists2D [private]

Definition at line 89 of file MillePedeMonitor.h.

Referenced by fillDerivatives(), and init().

std::vector<TH2*> MillePedeMonitor::myFrame2FrameHists2D [private]

Definition at line 95 of file MillePedeMonitor.h.

Referenced by fillFrameToFrame(), and init().

std::vector<TH2*> MillePedeMonitor::myResidHists2D [private]

Definition at line 90 of file MillePedeMonitor.h.

Referenced by fillResiduals(), and init().

std::vector<std::vector<TH1*> > MillePedeMonitor::myResidHistsVec1DX [private]

Definition at line 91 of file MillePedeMonitor.h.

Referenced by fillResiduals(), and init().

std::vector<std::vector<TH1*> > MillePedeMonitor::myResidHistsVec1DY [private]

[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC

Definition at line 92 of file MillePedeMonitor.h.

Referenced by fillResiduals(), and init().

std::vector<TH1*> MillePedeMonitor::myResidHitHists1DX [private]

[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC

Definition at line 93 of file MillePedeMonitor.h.

Referenced by fillResiduals(), and init().

std::vector<TH1*> MillePedeMonitor::myResidHitHists1DY [private]

Definition at line 94 of file MillePedeMonitor.h.

Referenced by fillResiduals(), and init().

TDirectory* MillePedeMonitor::myRootDir [private]

Definition at line 80 of file MillePedeMonitor.h.

Referenced by MillePedeMonitor(), and ~MillePedeMonitor().

std::vector<TH1*> MillePedeMonitor::myTrackHists1D [private]

Definition at line 83 of file MillePedeMonitor.h.

Referenced by fillTrack(), and init().

std::vector<TH2*> MillePedeMonitor::myTrackHists2D [private]

Definition at line 84 of file MillePedeMonitor.h.

Referenced by fillTrack(), and init().

std::vector<TH1*> MillePedeMonitor::myTrajectoryHists1D [private]

Definition at line 87 of file MillePedeMonitor.h.

Referenced by fillRefTrajectory(), and init().

std::vector<TH2*> MillePedeMonitor::myTrajectoryHists2D [private]

Definition at line 88 of file MillePedeMonitor.h.

Referenced by fillRefTrajectory(), and init().

std::vector<TH1*> MillePedeMonitor::myUsedTrackHists1D [private]

Definition at line 85 of file MillePedeMonitor.h.

Referenced by fillUsedTrack(), and init().

std::vector<TH2*> MillePedeMonitor::myUsedTrackHists2D [private]

Definition at line 86 of file MillePedeMonitor.h.

Referenced by fillUsedTrack(), and init().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:28:16 2009 for CMSSW by  doxygen 1.5.4