[352] | 1 | /***********************************************************************
|
---|
| 2 | ** **
|
---|
| 3 | ** /----------------------------------------------\ **
|
---|
| 4 | ** | Delphes, a framework for the fast simulation | **
|
---|
| 5 | ** | of a generic collider experiment | **
|
---|
[443] | 6 | ** \------------- arXiv:0903.2225v1 ------------/ **
|
---|
[352] | 7 | ** **
|
---|
| 8 | ** **
|
---|
| 9 | ** This package uses: **
|
---|
| 10 | ** ------------------ **
|
---|
[443] | 11 | ** ROOT: Nucl. Inst. & Meth. in Phys. Res. A389 (1997) 81-86 **
|
---|
| 12 | ** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
|
---|
| 13 | ** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
|
---|
[352] | 14 | ** FROG: [hep-ex/0901.2718v1] **
|
---|
[443] | 15 | ** HepMC: Comput. Phys. Commun.134 (2001) 41 **
|
---|
[352] | 16 | ** **
|
---|
| 17 | ** ------------------------------------------------------------------ **
|
---|
| 18 | ** **
|
---|
| 19 | ** Main authors: **
|
---|
| 20 | ** ------------- **
|
---|
| 21 | ** **
|
---|
[443] | 22 | ** Severine Ovyn Xavier Rouby **
|
---|
| 23 | ** severine.ovyn@uclouvain.be xavier.rouby@cern **
|
---|
[352] | 24 | ** **
|
---|
[443] | 25 | ** Center for Particle Physics and Phenomenology (CP3) **
|
---|
| 26 | ** Universite catholique de Louvain (UCL) **
|
---|
| 27 | ** Louvain-la-Neuve, Belgium **
|
---|
| 28 | ** **
|
---|
[352] | 29 | ** Copyright (C) 2008-2009, **
|
---|
[443] | 30 | ** All rights reserved. **
|
---|
[352] | 31 | ** **
|
---|
| 32 | ***********************************************************************/
|
---|
| 33 |
|
---|
| 34 | #include <iostream>
|
---|
| 35 | #include <utility>
|
---|
| 36 | #include <deque>
|
---|
| 37 | #include <sstream>
|
---|
| 38 | #include <fstream>
|
---|
[419] | 39 | #include <cmath>
|
---|
[352] | 40 |
|
---|
| 41 | #include "TROOT.h"
|
---|
| 42 | #include "TTree.h"
|
---|
| 43 | #include "TBranch.h"
|
---|
| 44 | #include "TLeaf.h"
|
---|
| 45 | #include "TString.h"
|
---|
| 46 | #include "TLorentzVector.h"
|
---|
| 47 | #include "TChain.h"
|
---|
| 48 |
|
---|
| 49 | #include "BlockClasses.h"
|
---|
| 50 | #include "ExRootTreeWriter.h"
|
---|
| 51 | #include "ExRootTreeBranch.h"
|
---|
| 52 |
|
---|
| 53 | #include "HEPEVTConverter.h"
|
---|
| 54 |
|
---|
| 55 | using namespace std;
|
---|
| 56 |
|
---|
| 57 | //------------------------------------------------------------------------------
|
---|
| 58 |
|
---|
| 59 | //Nevents not yet implemented! 08.03.2009
|
---|
| 60 |
|
---|
| 61 | HEPTreeReader::HEPTreeReader(TTree *tree, HEPEvent *event) : fChain(tree), fCurrentTree(-1), fEvent(event)
|
---|
| 62 | {
|
---|
| 63 | if(!fChain) return;
|
---|
| 64 |
|
---|
| 65 | TBranch *branch;
|
---|
| 66 | TLeaf *leaf;
|
---|
| 67 | TString name;
|
---|
| 68 | Int_t i;
|
---|
| 69 |
|
---|
| 70 | name = "Nhep";
|
---|
| 71 | branch = fChain->GetBranch(name);
|
---|
| 72 | branch->SetAddress(&(event->Nhep));
|
---|
| 73 | fBranches.push_back(make_pair(name, branch));
|
---|
| 74 |
|
---|
| 75 | TString intNames[3] = {"Idhep", "Jsmhep", "Jsdhep"};
|
---|
| 76 | Int_t **intData[3] = {&event->Idhep, &event->Jsmhep, &event->Jsdhep};
|
---|
| 77 |
|
---|
| 78 | for(i = 0; i < 3; ++i)
|
---|
| 79 | {
|
---|
| 80 | name = intNames[i];
|
---|
| 81 | branch = fChain->GetBranch(name);
|
---|
| 82 | leaf = branch->GetLeaf(name);
|
---|
| 83 | *intData[i] = new Int_t[leaf->GetNdata()];
|
---|
| 84 | branch->SetAddress(*intData[i]);
|
---|
| 85 | fBranches.push_back(make_pair(name, branch));
|
---|
| 86 | }
|
---|
| 87 |
|
---|
| 88 | TString floatNames[2] = {"Phep", "Vhep"};
|
---|
| 89 | Float_t **floatData[2] = {&event->Phep, &event->Vhep};
|
---|
| 90 |
|
---|
| 91 | for(i = 0; i < 2; ++i)
|
---|
| 92 | {
|
---|
| 93 | name = floatNames[i];
|
---|
| 94 | branch = fChain->GetBranch(name);
|
---|
| 95 | leaf = branch->GetLeaf(name);
|
---|
| 96 | *floatData[i] = new Float_t[leaf->GetNdata()];
|
---|
| 97 | branch->SetAddress(*floatData[i]);
|
---|
| 98 | fBranches.push_back(make_pair(name, branch));
|
---|
| 99 | }
|
---|
| 100 | }
|
---|
| 101 |
|
---|
| 102 | //------------------------------------------------------------------------------
|
---|
| 103 |
|
---|
| 104 | HEPTreeReader::~HEPTreeReader()
|
---|
| 105 | {
|
---|
| 106 | if(!fChain) return;
|
---|
| 107 |
|
---|
| 108 | delete[] fEvent->Idhep;
|
---|
| 109 | delete[] fEvent->Jsmhep;
|
---|
| 110 | delete[] fEvent->Jsdhep;
|
---|
| 111 | delete[] fEvent->Phep;
|
---|
| 112 | delete[] fEvent->Vhep;
|
---|
| 113 | }
|
---|
| 114 |
|
---|
| 115 | //------------------------------------------------------------------------------
|
---|
| 116 |
|
---|
| 117 | Bool_t HEPTreeReader::ReadEntry(Long64_t element)
|
---|
| 118 | {
|
---|
| 119 | if(!fChain) return kFALSE;
|
---|
| 120 |
|
---|
| 121 | Int_t treeEntry = fChain->LoadTree(element);
|
---|
| 122 | if(treeEntry < 0) return kFALSE;
|
---|
| 123 |
|
---|
| 124 | if(fChain->IsA() == TChain::Class())
|
---|
| 125 | {
|
---|
| 126 | TChain *chain = static_cast<TChain*>(fChain);
|
---|
| 127 | if(chain->GetTreeNumber() != fCurrentTree)
|
---|
| 128 | {
|
---|
| 129 | fCurrentTree = chain->GetTreeNumber();
|
---|
| 130 | Notify();
|
---|
| 131 | }
|
---|
| 132 | }
|
---|
| 133 |
|
---|
| 134 | deque< pair<TString, TBranch*> >::iterator it_deque;
|
---|
| 135 | TBranch *branch;
|
---|
| 136 |
|
---|
| 137 | for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
|
---|
| 138 | {
|
---|
| 139 | branch = it_deque->second;
|
---|
| 140 | if(branch)
|
---|
| 141 | {
|
---|
| 142 | branch->GetEntry(treeEntry);
|
---|
| 143 | }
|
---|
| 144 | }
|
---|
| 145 |
|
---|
| 146 | return kTRUE;
|
---|
| 147 | }
|
---|
| 148 |
|
---|
| 149 | //------------------------------------------------------------------------------
|
---|
| 150 |
|
---|
| 151 | void HEPTreeReader::Notify()
|
---|
| 152 | {
|
---|
| 153 | // Called when loading a new file.
|
---|
| 154 | // Get branch pointers.
|
---|
| 155 | if(!fChain) return;
|
---|
| 156 |
|
---|
| 157 | deque< pair<TString, TBranch*> >::iterator it_deque;
|
---|
| 158 |
|
---|
| 159 | for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
|
---|
| 160 | {
|
---|
| 161 | it_deque->second = fChain->GetBranch(it_deque->first);
|
---|
| 162 | if(!it_deque->second)
|
---|
| 163 | {
|
---|
| 164 | cerr << left << setw(40) <<"** WARNING: cannot get branch: "<<""
|
---|
| 165 | << left << setw(20) << it_deque->first <<""
|
---|
| 166 | << right << setw(9) <<" **"<<""<<endl;
|
---|
| 167 | }
|
---|
| 168 | }
|
---|
| 169 | }
|
---|
| 170 |
|
---|
| 171 |
|
---|
| 172 |
|
---|
[380] | 173 | HEPEVTConverter::HEPEVTConverter(const string& inputFileList, const string& outputFileName, const PdgTable& pdg, const int& Nevents) : DataConverter(pdg,Nevents)
|
---|
[352] | 174 | {
|
---|
| 175 | string buffer;
|
---|
| 176 |
|
---|
| 177 | gROOT->SetBatch();
|
---|
| 178 |
|
---|
| 179 | ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
|
---|
| 180 | ExRootTreeBranch *branchGen = treeWriter->NewBranch("Particle", TRootC::GenParticle::Class());
|
---|
| 181 |
|
---|
| 182 | ifstream infile(inputFileList.c_str());
|
---|
| 183 | if(!infile.is_open())
|
---|
| 184 | {
|
---|
| 185 | cerr << left << setw(30) <<"** ERROR: Can't open "<<""
|
---|
| 186 | << left << setw(20) << inputFileList <<""
|
---|
| 187 | << right << setw(19) <<"for input **"<<""<<endl;
|
---|
| 188 | exit(1);
|
---|
| 189 | }
|
---|
| 190 | while(1)
|
---|
| 191 | {
|
---|
| 192 |
|
---|
| 193 | TChain *chain = new TChain("h101");
|
---|
| 194 | infile >> buffer;
|
---|
| 195 | ifstream checking_the_file(buffer.c_str());
|
---|
| 196 | if(!checking_the_file.good())
|
---|
| 197 | {
|
---|
| 198 | cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
|
---|
| 199 | << left << setw(20) << buffer <<""
|
---|
| 200 | << right << setw(19) <<"for input **"<<""<<endl;
|
---|
| 201 | continue;
|
---|
| 202 | }
|
---|
| 203 | else checking_the_file.close();
|
---|
| 204 |
|
---|
| 205 | if(!infile.good()) break;
|
---|
| 206 | chain->Add(buffer.c_str());
|
---|
| 207 |
|
---|
| 208 | HEPEvent event;
|
---|
| 209 | HEPTreeReader *treeReader = new HEPTreeReader(chain, &event);
|
---|
| 210 |
|
---|
| 211 | Long64_t entry, allEntries = treeReader->GetEntries();
|
---|
| 212 | Int_t address, particle;
|
---|
| 213 | Float_t signEta;
|
---|
| 214 |
|
---|
| 215 | TRootC::GenParticle *element;
|
---|
[419] | 216 | allEntries = (Nevt<0)?allEntries:min((int)allEntries,Nevt); // do not miss the "+1"
|
---|
[352] | 217 | // Loop over all events
|
---|
| 218 | for(entry = 0; entry < allEntries; ++entry)
|
---|
| 219 | {
|
---|
| 220 | // Load selected branches with data from specified event
|
---|
| 221 | treeReader->ReadEntry(entry);
|
---|
| 222 | treeWriter->Clear();
|
---|
| 223 |
|
---|
| 224 | for(particle = 0; particle < event.Nhep; ++particle)
|
---|
| 225 | {
|
---|
| 226 | element = (TRootC::GenParticle*) branchGen->NewEntry();
|
---|
| 227 |
|
---|
| 228 | element->PID = event.Idhep[particle];
|
---|
| 229 | element->Status = event.Jsmhep[particle]/16000000
|
---|
| 230 | + event.Jsdhep[particle]/16000000*100;
|
---|
| 231 | element->M1 = (event.Jsmhep[particle]%16000000)/4000 -1;
|
---|
| 232 | element->M2 = event.Jsmhep[particle]%4000 -1;
|
---|
| 233 | element->D1 = (event.Jsdhep[particle]%16000000)/4000 -1;
|
---|
| 234 | element->D2 = event.Jsdhep[particle]%4000 -1;
|
---|
| 235 | address = particle*5;
|
---|
| 236 | element->E = event.Phep[address + 3];
|
---|
| 237 | element->Px = event.Phep[address + 0];
|
---|
| 238 | element->Py = event.Phep[address + 1];
|
---|
| 239 | element->Pz = event.Phep[address + 2];
|
---|
| 240 | element->M = event.Phep[address + 4];
|
---|
| 241 |
|
---|
| 242 | TLorentzVector vector(element->Px, element->Py, element->Pz, element->E);
|
---|
| 243 | element->PT = vector.Perp();
|
---|
| 244 | signEta = (element->Pz >= 0.0) ? 1.0 : -1.0;
|
---|
| 245 | element->Eta = vector.CosTheta()*vector.CosTheta() == 1.0 ? signEta*999.9 : vector.Eta();
|
---|
| 246 | element->Phi = vector.Phi();
|
---|
| 247 |
|
---|
| 248 | address = particle*4;
|
---|
| 249 | element->T = event.Vhep[address + 3];
|
---|
| 250 | element->X = event.Vhep[address + 0];
|
---|
| 251 | element->Y = event.Vhep[address + 1];
|
---|
| 252 | element->Z = event.Vhep[address + 2];
|
---|
| 253 | }
|
---|
| 254 | treeWriter->Fill();
|
---|
| 255 | }
|
---|
| 256 | delete chain;
|
---|
| 257 | delete treeReader;
|
---|
| 258 | }
|
---|
| 259 | treeWriter->Write();
|
---|
| 260 |
|
---|
| 261 | delete treeWriter;
|
---|
| 262 |
|
---|
| 263 | }
|
---|
| 264 |
|
---|
| 265 |
|
---|
| 266 |
|
---|
| 267 |
|
---|