CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/CalibMuon/DTCalibration/interface/vDriftHistos.h

Go to the documentation of this file.
00001 #ifndef vDriftHistos_H
00002 #define vDriftHistos_H
00003 
00004 #include "TH1.h"
00005 #include "TString.h"
00006 #include "TFile.h"
00007 #include "DTTMax.h"
00008 #include <string>
00009 
00010 // A set of histograms on chamber angle and position
00011 class h4DSegm{
00012  public:
00013   h4DSegm(std::string name_){
00014     TString N = name_.c_str();
00015     name=name_.c_str();
00016     h4DSegmXPosInCham     = new TH1F(N+"_h4DSegmXPosInCham", 
00017                                     "4D Segment x position (cm) in Chamber RF", 200, -200, 200); 
00018     h4DSegmYPosInCham     = new TH1F(N+"_h4DSegmYPosInCham", 
00019                                     "4D Segment y position (cm) in Chamber RF", 200, -200, 200); 
00020     h4DSegmPhiAngleInCham   = new TH1F(N+"_h4DSegmPhiAngleInCham",  
00021                                     "4D Segment phi angle (rad) in Chamber RF", 180, -180, 180); 
00022     h4DSegmThetaAngleInCham   = new TH1F(N+"_h4DSegmThetaAngleInCham",  
00023                                     "4D Segment theta angle (rad) in Chamber RF", 180, -180, 180); 
00024     h4DSegmImpactAngleInCham   = new TH1F(N+"_h4DSegmImpactAngleInCham",  
00025                                     "4D Segment impact angle (rad) in Chamber RF", 180, -180, 180); 
00026   }
00027  h4DSegm(TString name_, TFile* file){
00028     name=name_;
00029     h4DSegmXPosInCham  = (TH1F *) file->Get(name+"_h4DSegmXPosInCham"); 
00030     h4DSegmYPosInCham  = (TH1F *) file->Get(name+"_h4DSegmYPosInCham"); 
00031     h4DSegmPhiAngleInCham  = (TH1F *) file->Get(name+"_h4DSegmPhiAngleInCham"); 
00032     h4DSegmThetaAngleInCham  = (TH1F *) file->Get(name+"_h4DSegmThetaAngleInCham"); 
00033     h4DSegmImpactAngleInCham  = (TH1F *) file->Get(name+"_h4DSegmImpactAngleInCham"); 
00034  }
00035  ~h4DSegm(){
00036     delete h4DSegmXPosInCham;     
00037     delete h4DSegmYPosInCham;     
00038     delete h4DSegmPhiAngleInCham;   
00039     delete h4DSegmThetaAngleInCham;   
00040     delete h4DSegmImpactAngleInCham;   
00041  } 
00042 void Fill(float x, float y, float phi, float theta, float impact) {
00043     h4DSegmXPosInCham->Fill(x); 
00044     h4DSegmYPosInCham->Fill(y); 
00045     h4DSegmPhiAngleInCham->Fill(phi);   
00046     h4DSegmThetaAngleInCham->Fill(theta);   
00047     h4DSegmImpactAngleInCham->Fill(impact);   
00048 } 
00049 void Fill(float x, float phi) {
00050     h4DSegmXPosInCham->Fill(x); 
00051     h4DSegmPhiAngleInCham->Fill(phi);   
00052 } 
00053  void Write() {
00054     h4DSegmXPosInCham->Write();     
00055     h4DSegmYPosInCham->Write();     
00056     h4DSegmPhiAngleInCham->Write();   
00057     h4DSegmThetaAngleInCham->Write();   
00058     h4DSegmImpactAngleInCham->Write();   
00059   }
00060  public:
00061 
00062   TH1F *h4DSegmXPosInCham;     
00063   TH1F *h4DSegmYPosInCham;     
00064   TH1F *h4DSegmPhiAngleInCham;   
00065   TH1F *h4DSegmThetaAngleInCham;   
00066   TH1F *h4DSegmImpactAngleInCham;   
00067 
00068   TString name;
00069 };
00070 
00071 
00072 // A set of histograms on SL angle and position
00073 class h2DSegm{
00074  public:
00075   h2DSegm(std::string name_){
00076     TString N = name_.c_str();
00077     name=name_.c_str();
00078     h2DSegmPosInCham     = new TH1F(N+"_h2DSegmPosInCham", 
00079                                     "2D Segment position (cm) in Chamber RF", 200, -200, 200); 
00080     h2DSegmAngleInCham   = new TH1F(N+"_h2DSegmAngleInCham",  
00081                                     "2D Segment angle (rad) in Chamber RF", 200, -2, 2); 
00082     h2DSegmCosAngleInCham   = new TH1F(N+"_h2DSegmCosAngleInCham",  
00083                                        "2D Segment cos(angle) in Chamber RF", 200, -2, 2); 
00084   }
00085   h2DSegm(TString name_, TFile* file){
00086     name=name_;
00087 
00088     h2DSegmPosInCham  = (TH1F *) file->Get(name+"_h2DSegmPosInCham"); 
00089     h2DSegmAngleInCham  = (TH1F *) file->Get(name+"_h2DSegmAngleInCham"); 
00090     h2DSegmCosAngleInCham  = (TH1F *) file->Get(name+"_h2DSegmCosAngleInCham"); 
00091   }
00092   ~h2DSegm(){
00093     delete h2DSegmPosInCham;     
00094     delete h2DSegmAngleInCham;   
00095     delete h2DSegmCosAngleInCham;   
00096   }
00097   void Fill(float pos, float localAngle) {
00098 
00099     h2DSegmPosInCham->Fill(pos); 
00100     h2DSegmAngleInCham->Fill(atan(localAngle));   
00101     h2DSegmCosAngleInCham->Fill(cos(atan(localAngle)));   
00102   }
00103   void Write() {
00104     h2DSegmPosInCham->Write();     
00105     h2DSegmAngleInCham->Write();   
00106     h2DSegmCosAngleInCham->Write();   
00107   }
00108  public:
00109 
00110   TH1F *h2DSegmPosInCham;     
00111   TH1F *h2DSegmAngleInCham;   
00112   TH1F *h2DSegmCosAngleInCham;   
00113 
00114   TString name;
00115 };
00116 
00117 // A set of histograms on SL Tmax
00118 class hTMaxCell{
00119  public:
00120   hTMaxCell(TString name_){
00121     name = name_;
00122 
00123     // book TMax histograms 
00124     hTmax123      = new TH1F (name+"_Tmax123", "Tmax123 value", 2000, -1000., 1000.);
00125     hTmax124s72   = new TH1F (name+"_Tmax124_s72", "Tmax124 sigma=sqrt(7/2) value", 2000, -1000., 1000.);
00126     hTmax124s78   = new TH1F (name+"_Tmax124_s78", "Tmax124 sigma=sqrt(7/8) value", 2000, -1000., 1000.);
00127     hTmax134s72   = new TH1F (name+"_Tmax134_s72", "Tmax134 sigma=sqrt(7/2) value", 2000, -1000., 1000.);
00128     hTmax134s78   = new TH1F (name+"_Tmax134_s78", "Tmax134 sigma=sqrt(7/8) value", 2000, -1000., 1000.);
00129     hTmax234      = new TH1F (name+"_Tmax234", "Tmax234 value", 2000, -1000., 1000.);
00130     hTmax_3t0  = new TH1F (name+"_3t0", "Tmax+3*Delta(t0)", 2000, -1000., 1000.); 
00131     hTmax_3t0_0  = new TH1F (name+"_3t0_0", "Tmax+3*Delta(t0); 3 hits", 2000, -1000., 1000.); 
00132     hTmax_3t0_1  = new TH1F (name+"_3t0_1", "Tmax+3*Delta(t0); one t<5ns", 2000, -1000., 1000.); 
00133     hTmax_3t0_2  = new TH1F (name+"_3t0_2", "Tmax+3*Delta(t0); one t<10ns", 2000, -1000., 1000.); 
00134     hTmax_3t0_3  = new TH1F (name+"_3t0_3", "Tmax+3*Delta(t0); one t<20ns", 2000, -1000., 1000.); 
00135     hTmax_3t0_4  = new TH1F (name+"_3t0_4", "Tmax+3*Delta(t0); one t<50ns", 2000, -1000., 1000.); 
00136     hTmax_3t0_5  = new TH1F (name+"_3t0_5", "Tmax+3*Delta(t0); all t>50ns", 2000, -1000., 1000.); 
00137     hTmax_2t0  = new TH1F (name+"_2t0", "Tmax+2*Delta(t0)", 2000, -1000., 1000.);
00138     hTmax_2t0_0  = new TH1F (name+"_2t0_0", "Tmax+2*Delta(t0); 3 hits", 2000, -1000., 1000.); 
00139     hTmax_2t0_1  = new TH1F (name+"_2t0_1", "Tmax+2*Delta(t0); one t<5ns", 2000, -1000., 1000.); 
00140     hTmax_2t0_2  = new TH1F (name+"_2t0_2", "Tmax+2*Delta(t0); one t<10ns", 2000, -1000., 1000.); 
00141     hTmax_2t0_3  = new TH1F (name+"_2t0_3", "Tmax+2*Delta(t0); one t<20ns", 2000, -1000., 1000.); 
00142     hTmax_2t0_4  = new TH1F (name+"_2t0_4", "Tmax+2*Delta(t0); one t<50ns", 2000, -1000., 1000.); 
00143     hTmax_2t0_5  = new TH1F (name+"_2t0_5", "Tmax+2*Delta(t0); all t>50ns", 2000, -1000., 1000.); 
00144     hTmax_t0   = new TH1F (name+"_t0", "Tmax+Delta(t0)", 2000, -1000., 1000.);
00145     hTmax_t0_0  = new TH1F (name+"_t0_0", "Tmax+Delta(t0); 3 hits", 2000, -1000., 1000.); 
00146     hTmax_t0_1  = new TH1F (name+"_t0_1", "Tmax+Delta(t0); one t<5ns", 2000, -1000., 1000.); 
00147     hTmax_t0_2  = new TH1F (name+"_t0_2", "Tmax+Delta(t0); one t<10ns", 2000, -1000., 1000.); 
00148     hTmax_t0_3  = new TH1F (name+"_t0_3", "Tmax+Delta(t0); one t<20ns", 2000, -1000., 1000.); 
00149     hTmax_t0_4  = new TH1F (name+"_t0_4", "Tmax+Delta(t0); one t<50ns", 2000, -1000., 1000.); 
00150     hTmax_t0_5  = new TH1F (name+"_t0_5", "Tmax+Delta(t0); all t>50ns", 2000, -1000., 1000.); 
00151     hTmax_0    = new TH1F (name+"_0", "Tmax", 2000, -1000., 1000.);
00152   }
00153 
00154 
00155   hTMaxCell (TString name_, TFile* file){
00156     name=name_;
00157     hTmax123      = (TH1F *) file->Get(name+"_Tmax123");
00158     hTmax124s72   = (TH1F *) file->Get(name+"_Tmax124_s72");
00159     hTmax124s78   = (TH1F *) file->Get(name+"_Tmax124_s78");
00160     hTmax134s72   = (TH1F *) file->Get(name+"_Tmax134_s72");
00161     hTmax134s78   = (TH1F *) file->Get(name+"_Tmax134_s78");
00162     hTmax234      = (TH1F *) file->Get(name+"_Tmax234");
00163     hTmax_3t0  = (TH1F *) file->Get(name+"_3t0");
00164     hTmax_3t0_0  = (TH1F *) file->Get(name+"_3t0_0");
00165     hTmax_3t0_1  = (TH1F *) file->Get(name+"_3t0_1");
00166     hTmax_3t0_2  = (TH1F *) file->Get(name+"_3t0_2");
00167     hTmax_3t0_3  = (TH1F *) file->Get(name+"_3t0_3");
00168     hTmax_3t0_4  = (TH1F *) file->Get(name+"_3t0_4");
00169     hTmax_3t0_5  = (TH1F *) file->Get(name+"_3t0_5");
00170     hTmax_2t0  = (TH1F *) file->Get(name+"_2t0");
00171     hTmax_2t0_0  = (TH1F *) file->Get(name+"_2t0_0");
00172     hTmax_2t0_1  = (TH1F *) file->Get(name+"_2t0_1");
00173     hTmax_2t0_2  = (TH1F *) file->Get(name+"_2t0_2");
00174     hTmax_2t0_3  = (TH1F *) file->Get(name+"_2t0_3");
00175     hTmax_2t0_4  = (TH1F *) file->Get(name+"_2t0_4");
00176     hTmax_2t0_5  = (TH1F *) file->Get(name+"_2t0_5");
00177     hTmax_t0  = (TH1F *) file->Get(name+"_t0");
00178     hTmax_t0_1  = (TH1F *) file->Get(name+"_t0_1");
00179     hTmax_t0_2  = (TH1F *) file->Get(name+"_t0_2");
00180     hTmax_t0_3  = (TH1F *) file->Get(name+"_t0_3");
00181     hTmax_t0_4  = (TH1F *) file->Get(name+"_t0_4");
00182     hTmax_t0_5  = (TH1F *) file->Get(name+"_t0_5");
00183     hTmax_0  = (TH1F *) file->Get(name+"_0");
00184 
00185   }
00186 
00187  
00188   ~hTMaxCell(){
00189     delete hTmax123;
00190     delete hTmax124s72;
00191     delete hTmax124s78;
00192     delete hTmax134s72;
00193     delete hTmax134s78;
00194     delete hTmax234;
00195     delete hTmax_3t0;
00196     delete hTmax_3t0_0;
00197     delete hTmax_3t0_1;
00198     delete hTmax_3t0_2;
00199     delete hTmax_3t0_3;
00200     delete hTmax_3t0_4;
00201     delete hTmax_3t0_5;
00202     delete hTmax_2t0;
00203     delete hTmax_2t0_0;
00204     delete hTmax_2t0_1;
00205     delete hTmax_2t0_2;
00206     delete hTmax_2t0_3;
00207     delete hTmax_2t0_4;
00208     delete hTmax_2t0_5;
00209     delete hTmax_t0;
00210     delete hTmax_t0_0;
00211     delete hTmax_t0_1;
00212     delete hTmax_t0_2;
00213     delete hTmax_t0_3;
00214     delete hTmax_t0_4;
00215     delete hTmax_t0_5;
00216     delete hTmax_0;
00217 
00218   }
00219 
00220   void Fill(float tmax123, float tmax124, float tmax134, float tmax234,
00221             dttmaxenums::SigmaFactor s124,// Give the factor relating the width of the TMax distribution
00222             dttmaxenums::SigmaFactor s134,// and the cell resolution (for tmax123 and tmax234 is always= sqrt(3/2)).
00223             unsigned t0_123, // Give the "quantity" of Delta(t0) included in the tmax 
00224             unsigned t0_124, // formula used for Tmax123 or Tma124 or Tma134 or Tma234
00225             unsigned t0_134,
00226             unsigned t0_234,
00227             unsigned hSubGroup//different t0 hists(at least one hit within a given distance from the wire)
00228             ) {
00229  
00230     if(tmax123 > 0.) {
00231       hTmax123->Fill(tmax123);
00232       if(t0_123==1) {
00233         hTmax_t0->Fill(tmax123);
00234         switch(hSubGroup) {
00235         case 0: hTmax_t0_0->Fill(tmax123); break;
00236         case 1: hTmax_t0_1->Fill(tmax123); break;
00237         case 2: hTmax_t0_2->Fill(tmax123); break;
00238         case 3: hTmax_t0_3->Fill(tmax123); break;
00239         case 4: hTmax_t0_4->Fill(tmax123); break;
00240         case 99: hTmax_t0_5->Fill(tmax123); break;
00241         }
00242       }
00243       else {
00244         hTmax_2t0->Fill(tmax123);
00245         switch(hSubGroup) {
00246         case 0: hTmax_2t0_0->Fill(tmax123); break;
00247         case 1: hTmax_2t0_1->Fill(tmax123); break;
00248         case 2: hTmax_2t0_2->Fill(tmax123); break;
00249         case 3: hTmax_2t0_3->Fill(tmax123); break;
00250         case 4: hTmax_2t0_4->Fill(tmax123); break;
00251         case 99: hTmax_2t0_5->Fill(tmax123); break;
00252         }
00253       }
00254     }
00255     if(tmax124 > 0.) {
00256       (s124==dttmaxenums::r72)? hTmax124s72->Fill(tmax124):hTmax124s78->Fill(tmax124);
00257       if(t0_124==0)  
00258         hTmax_0->Fill(tmax124);
00259       else if(t0_124==1) {
00260         hTmax_t0->Fill(tmax124); 
00261         switch(hSubGroup) {
00262         case 0: hTmax_t0_0->Fill(tmax124); break;
00263         case 1: hTmax_t0_1->Fill(tmax124); break;
00264         case 2: hTmax_t0_2->Fill(tmax124); break;
00265         case 3: hTmax_t0_3->Fill(tmax124); break;
00266         case 4: hTmax_t0_4->Fill(tmax124); break;
00267         case 99: hTmax_t0_5->Fill(tmax124); break;
00268         }
00269       }
00270       else if(t0_124== 2) {
00271         hTmax_2t0->Fill(tmax124);
00272         switch(hSubGroup) {
00273         case 0: hTmax_2t0_0->Fill(tmax124); break;
00274         case 1: hTmax_2t0_1->Fill(tmax124); break;
00275         case 2: hTmax_2t0_2->Fill(tmax124); break;
00276         case 3: hTmax_2t0_3->Fill(tmax124); break;
00277         case 4: hTmax_2t0_4->Fill(tmax124); break;
00278         case 99: hTmax_2t0_5->Fill(tmax124); break;
00279         }
00280       }
00281       else if(t0_124==3) {
00282         hTmax_3t0->Fill(tmax124); 
00283         switch(hSubGroup) {
00284         case 0: hTmax_3t0_0->Fill(tmax124); break;
00285         case 1: hTmax_3t0_1->Fill(tmax124); break;
00286         case 2: hTmax_3t0_2->Fill(tmax124); break;
00287         case 3: hTmax_3t0_3->Fill(tmax124); break;
00288         case 4: hTmax_3t0_4->Fill(tmax124); break;
00289         case 99: hTmax_3t0_5->Fill(tmax124);break;
00290         }      
00291       }
00292     }
00293     if(tmax134 > 0.) {
00294       (s134==dttmaxenums::r72)? hTmax134s72->Fill(tmax134):hTmax134s78->Fill(tmax134);
00295       if(t0_134==0)  
00296         hTmax_0->Fill(tmax134);
00297       else if(t0_134==1) {
00298         hTmax_t0->Fill(tmax134); 
00299         switch(hSubGroup) {
00300         case 0: hTmax_t0_0->Fill(tmax134); break;
00301         case 1: hTmax_t0_1->Fill(tmax134); break;
00302         case 2: hTmax_t0_2->Fill(tmax134); break;
00303         case 3: hTmax_t0_3->Fill(tmax134); break;
00304         case 4: hTmax_t0_4->Fill(tmax134); break;
00305         case 99: hTmax_t0_5->Fill(tmax134); break;
00306         }
00307       }
00308       else if(t0_134== 2) {
00309         hTmax_2t0->Fill(tmax134);
00310         switch(hSubGroup) {
00311         case 0: hTmax_2t0_0->Fill(tmax134); break;
00312         case 1: hTmax_2t0_1->Fill(tmax134); break;
00313         case 2: hTmax_2t0_2->Fill(tmax134); break;
00314         case 3: hTmax_2t0_3->Fill(tmax134); break;
00315         case 4: hTmax_2t0_4->Fill(tmax134); break;
00316         case 99: hTmax_2t0_5->Fill(tmax134); break;
00317         }
00318       }
00319       else if(t0_134==3) {
00320         hTmax_3t0->Fill(tmax134); 
00321         switch(hSubGroup) {
00322         case 0: hTmax_3t0_0->Fill(tmax134); break;
00323         case 1: hTmax_3t0_1->Fill(tmax134); break;
00324         case 2: hTmax_3t0_2->Fill(tmax134); break;
00325         case 3: hTmax_3t0_3->Fill(tmax134); break;
00326         case 4: hTmax_3t0_4->Fill(tmax134); break;
00327         case 99: hTmax_3t0_5->Fill(tmax134); break;
00328         }      
00329       }
00330     }
00331     if(tmax234 > 0.) {
00332       hTmax234->Fill(tmax234);
00333       if(t0_234==1) {
00334         hTmax_t0->Fill(tmax234);
00335         switch(hSubGroup) {
00336         case 0: hTmax_t0_0->Fill(tmax234); break;
00337         case 1: hTmax_t0_1->Fill(tmax234); break;
00338         case 2: hTmax_t0_2->Fill(tmax234); break;
00339         case 3: hTmax_t0_3->Fill(tmax234); break;
00340         case 4: hTmax_t0_4->Fill(tmax234); break;
00341         case 99: hTmax_t0_5->Fill(tmax234);break;
00342         }
00343       }
00344       else {
00345         hTmax_2t0->Fill(tmax234);
00346         switch(hSubGroup) {
00347         case 0: hTmax_2t0_0->Fill(tmax234); break;
00348         case 1: hTmax_2t0_1->Fill(tmax234); break;
00349         case 2: hTmax_2t0_2->Fill(tmax234); break;
00350         case 3: hTmax_2t0_3->Fill(tmax234); break;
00351         case 4: hTmax_2t0_4->Fill(tmax234); break;
00352         case 99: hTmax_2t0_5->Fill(tmax234); break;
00353         }
00354       }
00355     }      
00356   }
00357 
00358   void Write() { 
00359     // write the Tmax histograms 
00360     hTmax123->Write();
00361     hTmax124s72->Write();
00362     hTmax124s78->Write();
00363     hTmax134s72->Write();
00364     hTmax134s78->Write();
00365     hTmax234->Write(); 
00366     hTmax_3t0->Write();
00367     hTmax_3t0_0->Write();
00368     hTmax_3t0_1->Write();
00369     hTmax_3t0_2->Write();
00370     hTmax_3t0_3->Write();
00371     hTmax_3t0_4->Write();
00372     hTmax_3t0_5->Write();
00373     hTmax_2t0->Write();
00374     hTmax_2t0_0->Write();
00375     hTmax_2t0_1->Write();
00376     hTmax_2t0_2->Write();
00377     hTmax_2t0_3->Write();
00378     hTmax_2t0_4->Write();
00379     hTmax_2t0_5->Write();
00380     hTmax_t0->Write();
00381     hTmax_t0_0->Write();
00382     hTmax_t0_1->Write();
00383     hTmax_t0_2->Write();
00384     hTmax_t0_3->Write();
00385     hTmax_t0_4->Write();
00386     hTmax_t0_5->Write();
00387     hTmax_0->Write();
00388 
00389   }
00390 
00391   int GetT0Factor(TH1F* hist) {
00392     unsigned t0 = 999;
00393 
00394     if (hist == hTmax_3t0)
00395       t0 = 3;
00396     else if(hist == hTmax_2t0)
00397       t0 = 2;
00398     else if(hist ==  hTmax_t0)
00399       t0 = 1;    
00400     else if (hist == hTmax_0)
00401       t0 = 0;
00402        
00403     return t0;
00404   }
00405  
00406   TH1F* hTmax123;
00407   TH1F* hTmax124s72;
00408   TH1F* hTmax124s78;
00409   TH1F* hTmax134s72;
00410   TH1F* hTmax134s78;
00411   TH1F* hTmax234;
00412   TH1F* hTmax_3t0;
00413   TH1F* hTmax_3t0_0;
00414   TH1F* hTmax_3t0_1;
00415   TH1F* hTmax_3t0_2;
00416   TH1F* hTmax_3t0_3;
00417   TH1F* hTmax_3t0_4;
00418   TH1F* hTmax_3t0_5;
00419   TH1F* hTmax_2t0;
00420   TH1F* hTmax_2t0_0;
00421   TH1F* hTmax_2t0_1;
00422   TH1F* hTmax_2t0_2;
00423   TH1F* hTmax_2t0_3;
00424   TH1F* hTmax_2t0_4;
00425   TH1F* hTmax_2t0_5;
00426   TH1F* hTmax_t0;
00427   TH1F* hTmax_t0_0;
00428   TH1F* hTmax_t0_1;
00429   TH1F* hTmax_t0_2;
00430   TH1F* hTmax_t0_3;
00431   TH1F* hTmax_t0_4;
00432   TH1F* hTmax_t0_5;
00433   TH1F* hTmax_0;
00434 
00435   TString name;  
00436 };
00437 
00438 #endif