Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:43

0001 #include "Calibration/EcalTBTools/interface/TB06Reco.h"
0002 #include "Calibration/EcalTBTools/interface/TB06Tree.h"
0003 #include "TFile.h"
0004 #include "TTree.h"
0005 
0006 #include <iostream>
0007 
0008 TB06Tree::TB06Tree(const std::string &fileName, const std::string &treeName)
0009     : m_file(nullptr), m_tree(nullptr), m_data(nullptr), m_dataSize(0) {
0010   TDirectory *dir = gDirectory;
0011   m_file = new TFile(fileName.c_str(), "RECREATE");
0012   m_file->cd();
0013   m_tree = new TTree(treeName.c_str(), "Analysis tree");
0014   m_tree->SetAutoSave(10000000);
0015   dir->cd();
0016 
0017   //  m_tree->cd () ;
0018   m_data = new TClonesArray(TB06Reco::Class(), 1);
0019   m_data->ExpandCreateFast(1);
0020 
0021   //  m_tree->Branch ("EGCO", &m_data, 64000, 2) ;
0022   m_tree->Branch("TB06O", &m_data, 64000, 2);
0023   m_tree->Print();
0024 }
0025 
0026 // -------------------------------------------------------------------
0027 
0028 TB06Tree::~TB06Tree() {
0029   std::cout << "[TB06Tree][dtor] saving TTree " << m_tree->GetName() << " with " << m_tree->GetEntries() << " entries"
0030             << " on file: " << m_file->GetName() << std::endl;
0031 
0032   m_file->Write();
0033   delete m_tree;
0034   m_file->Close();
0035   delete m_file;
0036   delete m_data;
0037 }
0038 
0039 // -------------------------------------------------------------------
0040 
0041 //! to be called at each loop
0042 void TB06Tree::store(const int &tableIsMoving,
0043                      const int &run,
0044                      const int &event,
0045                      const int &S6adc,
0046                      const double &xhodo,
0047                      const double &yhodo,
0048                      const double &xslope,
0049                      const double &yslope,
0050                      const double &xquality,
0051                      const double &yquality,
0052                      const int &icMax,
0053                      const int &ietaMax,
0054                      const int &iphiMax,
0055                      const double &beamEnergy,
0056                      const double ampl[49]) {
0057   m_data->Clear();
0058   TB06Reco *entry = static_cast<TB06Reco *>(m_data->AddrAt(0));
0059 
0060   entry->reset();
0061   //  reset (entry->myCalibrationMap) ;
0062 
0063   entry->tableIsMoving = tableIsMoving;
0064   entry->run = run;
0065   entry->event = event;
0066   entry->S6ADC = S6adc;
0067 
0068   entry->MEXTLindex = icMax;
0069   entry->MEXTLeta = ietaMax;
0070   entry->MEXTLphi = iphiMax;
0071   entry->MEXTLenergy = ampl[24];
0072   entry->beamEnergy = beamEnergy;
0073 
0074   for (int eta = 0; eta < 7; ++eta)
0075     for (int phi = 0; phi < 7; ++phi) {
0076       // FIXME capire l'orientamento di phi!
0077       // FIXME capire se eta, phi iniziano da 1 o da 0
0078       entry->localMap[eta][phi] = ampl[eta * 7 + phi];
0079     }
0080 
0081   entry->xHodo = xhodo;
0082   entry->yHodo = yhodo;
0083   entry->xSlopeHodo = xslope;
0084   entry->ySlopeHodo = yslope;
0085   entry->xQualityHodo = xquality;
0086   entry->yQualityHodo = yquality;
0087 
0088   entry->convFactor = 0.;
0089 
0090   /*
0091   // loop over the 5x5 see (1)
0092   for (int xtal=0 ; xtal<25 ; ++xtal)
0093     {
0094       int ieta = xtal/5 + 3 ;
0095       int iphi = xtal%5 + 8 ;
0096       entry->myCalibrationMap[ieta][iphi] = ampl[xtal] ;
0097     } // loop over the 5x5
0098 
0099   entry->electron_Tr_Pmag_ = beamEnergy ;
0100 
0101         entry->centralCrystalEta_ = ietaMax ;
0102         entry->centralCrystalPhi_ = iphiMax ;
0103         entry->centralCrystalEnergy_ = ampl[12] ;
0104 
0105   // this is a trick
0106   entry->electron_Tr_Peta_ = xhodo ;
0107   entry->electron_Tr_Pphi_ = yhodo ;
0108   */
0109   m_tree->Fill();
0110 }
0111 
0112 // -------------------------------------------------------------------
0113 
0114 void TB06Tree::reset(float crystal[11][21]) {
0115   for (int eta = 0; eta < 11; ++eta) {
0116     for (int phi = 0; phi < 21; ++phi) {
0117       crystal[eta][phi] = -999.;
0118     }
0119   }
0120 }
0121 
0122 // -------------------------------------------------------------------
0123 
0124 void TB06Tree::check() {
0125   TB06Reco *entry = static_cast<TB06Reco *>(m_data->AddrAt(0));
0126 
0127   std::cout << "[TB06Tree][check]reading . . . \n";
0128   std::cout << "[TB06Tree][check] entry->run: " << entry->run << "\n";
0129   std::cout << "[TB06Tree][check] entry->event: " << entry->event << "\n";
0130   std::cout << "[TB06Tree][check] entry->tableIsMoving: " << entry->tableIsMoving << "\n";
0131   std::cout << "[TB06Tree][check] entry->MEXTLeta: " << entry->MEXTLeta << "\n";
0132   std::cout << "[TB06Tree][check] entry->MEXTLphi: " << entry->MEXTLphi << "\n";
0133   std::cout << "[TB06Tree][check] entry->MEXTLenergy: " << entry->MEXTLenergy << "\n";
0134 
0135   for (int eta = 0; eta < 7; ++eta)
0136     for (int phi = 0; phi < 7; ++phi)
0137       std::cout << "[TB06Tree][check]   entry->localMap[" << eta << "][" << phi << "]: " << entry->localMap[eta][phi]
0138                 << "\n";
0139 
0140   std::cout << "[TB06Tree][check] entry->xHodo: " << entry->xHodo << "\n";
0141   std::cout << "[TB06Tree][check] entry->yHodo: " << entry->yHodo << "\n";
0142   std::cout << "[TB06Tree][check] entry->xSlopeHodo: " << entry->xSlopeHodo << "\n";
0143   std::cout << "[TB06Tree][check] entry->ySlopeHodo: " << entry->ySlopeHodo << "\n";
0144   std::cout << "[TB06Tree][check] entry->xQualityHodo: " << entry->xQualityHodo << "\n";
0145   std::cout << "[TB06Tree][check] entry->yQualityHodo: " << entry->yQualityHodo << "\n";
0146   std::cout << "[TB06Tree][check] entry->convFactor: " << entry->convFactor << "\n";
0147 
0148   /* to be implemented with the right variables
0149   std::cout << "[TB06Tree][check] ------------------------" << std::endl ;
0150   std::cout << "[TB06Tree][check] " << entry->variable_name << std::endl ;
0151   */
0152 }
0153 
0154 /* (1) to fill the 25 crystals vector
0155 
0156    for (UInt_t icry=0 ; icry<25 ; ++icry)
0157      {
0158        UInt_t row = icry / 5 ;
0159        Int_t column = icry % 5 ;
0160        try
0161            {
0162              EBDetId tempo (maxHitId.ieta()+column-2,
0163                             maxHitId.iphi()+row-2,
0164                             EBDetId::ETAPHIMODE) ;
0165 
0166              Xtals5x5.push_back (tempo) ;
0167              amplitude [icry] = hits->find (Xtals5x5[icry])->energy () ;
0168 
0169            }
0170        catch ( std::runtime_error &e )
0171            {
0172              std::cout << "Cannot construct 5x5 matrix around EBDetId "
0173                      << maxHitId << std::endl ;
0174              return ;
0175            }
0176      } // loop over the 5x5 matrix
0177 
0178 
0179 */