Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:48

0001 #include "TFile.h"
0002 #include "TTree.h"
0003 #include "TTreeReader.h"
0004 #include "TTreeReaderArray.h"
0005 #include <vector>
0006 
0007 // takes "new" form of HF shower library and makes it a v3 version
0008 
0009 void convertHFShowerLibrary() {
0010   TFile *nF=new TFile("david.root","RECREATE","hiThere",1);
0011   TTree *nT=new TTree("HFSimHits","HF simple tree");
0012   std::vector<float> *parts=new std::vector<float>();
0013   std::vector<float> *partsHad=new std::vector<float>();
0014   nT->Branch("emParticles", &parts);
0015   nT->GetBranch("emParticles")->SetBasketSize(1);
0016   nT->Branch("hadParticles", &partsHad);
0017   nT->GetBranch("hadParticles")->SetBasketSize(1);
0018 
0019   TDirectory *target=gDirectory;
0020   TFile *oF=TFile::Open("HFShowerLibrary_npmt_noatt_eta4_16en_v3_orig.root");
0021   TTree *t=(TTree*)oF->Get("HFSimHits");
0022   TTreeReader fReader(t);
0023   TTreeReaderArray<Float_t> b1(fReader, "emParticles.position_.fCoordinates.fX");
0024   TTreeReaderArray<Float_t> b2(fReader, "emParticles.position_.fCoordinates.fY");
0025   TTreeReaderArray<Float_t> b3(fReader, "emParticles.position_.fCoordinates.fZ");
0026   TTreeReaderArray<Float_t> b4(fReader, "emParticles.lambda_");
0027   TTreeReaderArray<Float_t> b5(fReader, "emParticles.time_");
0028 
0029   TTreeReaderArray<Float_t> h1(fReader, "hadParticles.position_.fCoordinates.fX");
0030   TTreeReaderArray<Float_t> h2(fReader, "hadParticles.position_.fCoordinates.fY");
0031   TTreeReaderArray<Float_t> h3(fReader, "hadParticles.position_.fCoordinates.fZ");
0032   TTreeReaderArray<Float_t> h4(fReader, "hadParticles.lambda_");
0033   TTreeReaderArray<Float_t> h5(fReader, "hadParticles.time_");
0034 
0035 
0036   target->cd();
0037   while ( fReader.Next() ) {
0038     parts->clear();
0039     unsigned int s=b1.GetSize();
0040     parts->resize(5*s);
0041     for ( unsigned int i=0; i<b1.GetSize(); i++) {
0042       (*parts)[i]=(b1[i]);
0043       (*parts)[i+1*s]=(b2[i]);
0044       (*parts)[i+2*s]=(b3[i]);
0045       (*parts)[i+3*s]=(b4[i]);
0046       (*parts)[i+4*s]=(b5[i]);
0047     }  
0048 
0049     partsHad->clear();
0050     s=h1.GetSize();
0051     partsHad->resize(5*s);
0052     for ( unsigned int i=0; i<h1.GetSize(); i++) {
0053       (*partsHad)[i]=(h1[i]);
0054       (*partsHad)[i+1*s]=(h2[i]);
0055       (*partsHad)[i+2*s]=(h3[i]);
0056       (*partsHad)[i+3*s]=(h4[i]);
0057       (*partsHad)[i+4*s]=(h5[i]);
0058     }  
0059 
0060     nT->Fill();
0061   }
0062 
0063   nT->Write();
0064   nF->Close();
0065 }