CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/Calibration/EcalTBTools/src/TB06Tree.cc

Go to the documentation of this file.
00001 #include "Calibration/EcalTBTools/interface/TB06Tree.h"
00002 #include "TFile.h" 
00003 #include "TTree.h" 
00004 #include "Calibration/EcalTBTools/interface/TB06Reco.h"
00005 
00006 #include <iostream>
00007 
00008 
00009 TB06Tree::TB06Tree (const std::string & fileName, 
00010                     const std::string & treeName):
00011   m_file (0), m_tree (0), m_data (0), m_dataSize (0)  
00012 {
00013   TDirectory *dir = gDirectory ;
00014   m_file = new TFile (fileName.c_str (),"RECREATE") ;
00015   m_file->cd () ;
00016   m_tree = new TTree (treeName.c_str(),"Analysis tree") ;
00017   m_tree->SetAutoSave (10000000) ;
00018   dir->cd () ;
00019 
00020   //  m_tree->cd () ;
00021   m_data = new TClonesArray (TB06Reco::Class (), 1) ;
00022   m_data->ExpandCreateFast (1) ;
00023 
00024   //  m_tree->Branch ("EGCO", &m_data, 64000, 2) ;
00025   m_tree->Branch ("TB06O", &m_data, 64000, 2) ;
00026   m_tree->Print () ;
00027  
00028 }
00029 
00030 
00031 // -------------------------------------------------------------------      
00032 
00033 TB06Tree::~TB06Tree () 
00034 {
00035   std::cout << "[TB06Tree][dtor] saving TTree " << m_tree->GetName ()
00036             << " with " << m_tree->GetEntries () << " entries"
00037             << " on file: " << m_file->GetName () << std::endl ;
00038 
00039   m_file->Write () ;
00040   delete m_tree ;
00041   m_file->Close () ;
00042   delete m_file ;
00043   delete m_data ;
00044 }
00045 
00046 
00047 // -------------------------------------------------------------------   
00048 
00050   void TB06Tree::store (const int & tableIsMoving,
00051                       const int & run, const int & event,
00052                       const int & S6adc ,
00053                       const double & xhodo, const double & yhodo, 
00054                       const double & xslope, const double & yslope, 
00055                       const double & xquality, const double & yquality,
00056                       const int & icMax,
00057                       const int & ietaMax, const int & iphiMax,
00058                       const double & beamEnergy, 
00059                       const double ampl[49]) 
00060 {
00061 
00062   m_data->Clear () ;
00063   TB06Reco * entry = static_cast<TB06Reco*> (m_data->AddrAt (0)) ;
00064 
00065   entry->reset () ;
00066   //  reset (entry->myCalibrationMap) ;
00067 
00068   entry->tableIsMoving = tableIsMoving ;
00069   entry->run = run ;
00070   entry->event = event ;
00071   entry->S6ADC = S6adc ;
00072 
00073   entry->MEXTLindex = icMax ;
00074   entry->MEXTLeta = ietaMax ;
00075   entry->MEXTLphi = iphiMax ;
00076   entry->MEXTLenergy = ampl[24] ;
00077   entry->beamEnergy = beamEnergy ;
00078 
00079   for (int eta = 0 ; eta<7 ; ++eta)
00080     for (int phi = 0 ; phi<7 ; ++phi)
00081       {
00082         // FIXME capire l'orientamento di phi!
00083         // FIXME capire se eta, phi iniziano da 1 o da 0
00084         entry->localMap[eta][phi] = ampl[eta*7+phi] ;
00085       }
00086 
00087   entry->xHodo = xhodo ;
00088   entry->yHodo = yhodo ;
00089   entry->xSlopeHodo = xslope ;
00090   entry->ySlopeHodo = yslope ;
00091   entry->xQualityHodo = xquality ;
00092   entry->yQualityHodo = yquality ;
00093 
00094   entry->convFactor = 0. ;
00095 
00096   /*
00097   // loop over the 5x5 see (1)
00098   for (int xtal=0 ; xtal<25 ; ++xtal) 
00099     {
00100       int ieta = xtal/5 + 3 ;
00101       int iphi = xtal%5 + 8 ;
00102       entry->myCalibrationMap[ieta][iphi] = ampl[xtal] ;
00103     } // loop over the 5x5
00104 
00105   entry->electron_Tr_Pmag_ = beamEnergy ;
00106   
00107         entry->centralCrystalEta_ = ietaMax ;
00108         entry->centralCrystalPhi_ = iphiMax ;
00109         entry->centralCrystalEnergy_ = ampl[12] ; 
00110 
00111   // this is a trick
00112   entry->electron_Tr_Peta_ = xhodo ;
00113   entry->electron_Tr_Pphi_ = yhodo ;
00114   */
00115   m_tree->Fill () ;
00116 }
00117 
00118 
00119 // -------------------------------------------------------------------
00120 
00121 
00122 void TB06Tree::reset (float crystal[11][21])
00123 {
00124   for (int eta =0 ; eta<11 ; ++eta)    
00125     {
00126       for (int phi =0 ; phi<21 ; ++phi) 
00127         {
00128           crystal[eta][phi] = -999. ;  
00129         }   
00130     }
00131 }
00132 
00133 
00134 // -------------------------------------------------------------------
00135 
00136 
00137 void TB06Tree::check ()
00138 {
00139   TB06Reco * entry = static_cast<TB06Reco*> (m_data->AddrAt (0)) ;
00140 
00141   std::cout << "[TB06Tree][check]reading . . . \n" ;
00142   std::cout << "[TB06Tree][check] entry->run: " << entry->run << "\n" ;
00143   std::cout << "[TB06Tree][check] entry->event: " << entry->event << "\n" ;
00144   std::cout << "[TB06Tree][check] entry->tableIsMoving: " << entry->tableIsMoving << "\n" ;
00145   std::cout << "[TB06Tree][check] entry->MEXTLeta: " << entry->MEXTLeta << "\n" ;
00146   std::cout << "[TB06Tree][check] entry->MEXTLphi: " << entry->MEXTLphi << "\n" ;
00147   std::cout << "[TB06Tree][check] entry->MEXTLenergy: " << entry->MEXTLenergy << "\n" ;
00148 
00149   for (int eta = 0 ; eta<7 ; ++eta)
00150       for (int phi = 0 ; phi<7 ; ++phi)
00151         std::cout << "[TB06Tree][check]   entry->localMap[" << eta
00152                   << "][" << phi << "]: "
00153                   << entry->localMap[eta][phi] << "\n" ;
00154 
00155   std::cout << "[TB06Tree][check] entry->xHodo: " << entry->xHodo << "\n" ;
00156   std::cout << "[TB06Tree][check] entry->yHodo: " << entry->yHodo << "\n" ;
00157   std::cout << "[TB06Tree][check] entry->xSlopeHodo: " << entry->xSlopeHodo << "\n" ;
00158   std::cout << "[TB06Tree][check] entry->ySlopeHodo: " << entry->ySlopeHodo << "\n" ;
00159   std::cout << "[TB06Tree][check] entry->xQualityHodo: " << entry->xQualityHodo << "\n" ;
00160   std::cout << "[TB06Tree][check] entry->yQualityHodo: " << entry->yQualityHodo << "\n" ;
00161   std::cout << "[TB06Tree][check] entry->convFactor: " << entry->convFactor << "\n" ;
00162 
00163   /* to be implemented with the right variables
00164   std::cout << "[TB06Tree][check] ------------------------" << std::endl ;
00165   std::cout << "[TB06Tree][check] " << entry->variable_name << std::endl ;
00166   */
00167 }
00168 
00169 
00170 
00171 /* (1) to fill the 25 crystals vector
00172 
00173    for (UInt_t icry=0 ; icry<25 ; ++icry)
00174      {
00175        UInt_t row = icry / 5 ;
00176        Int_t column = icry % 5 ;
00177        try
00178            {
00179              EBDetId tempo (maxHitId.ieta()+column-2, 
00180                             maxHitId.iphi()+row-2, 
00181                             EBDetId::ETAPHIMODE) ;
00182 
00183              Xtals5x5.push_back (tempo) ;
00184              amplitude [icry] = hits->find (Xtals5x5[icry])->energy () ;
00185              
00186            }
00187        catch ( std::runtime_error &e )
00188            {
00189              std::cout << "Cannot construct 5x5 matrix around EBDetId " 
00190                      << maxHitId << std::endl ;
00191              return ;
00192            }
00193      } // loop over the 5x5 matrix
00194 
00195 
00196 */