// ------------------------------------------------------------------------- // ----- TpcRootGiBUUGenerator source file ----- // ----- Created 03/04/12 by S. Doerheim ----- // ------------------------------------------------------------------------- #include "TpcRootGiBUUGenerator.h" #include "FairPrimaryGenerator.h" #include "TDatabasePDG.h" #include "TRandom3.h" #include "TMath.h" #include #include using std::cout; using std::endl; // ----- Default constructor ------------------------------------------ TpcRootGiBUUGenerator::TpcRootGiBUUGenerator() :FairGenerator(), fInputFile(NULL), fInputTree(NULL), filter_on_PDG(false), fFileName(""), pdgToFilter(0), evNb(0), evNbFile(0), targetThickness(5), vetoR(2) { } // ------------------------------------------------------------------------ // ----- Standard constructor ----------------------------------------- TpcRootGiBUUGenerator::TpcRootGiBUUGenerator(const char* fileName) :FairGenerator(), fInputFile(0), fFileName(fileName), filter_on_PDG(false), pdgToFilter(0), targetThickness(5), vetoR(2) { cout << "-I TpcRootGiBUUGenerator: Opening input file " << fileName << endl; fInputFile = TFile::Open(fFileName,"read"); if ( fInputFile->IsZombie() ) { Fatal("TpcRootGiBUUGenerator","Cannot open input file."); } fInputTree=(TTree*)fInputFile->Get("h47"); if ( fInputTree==0 || fInputTree->IsZombie() ) { Fatal("TpcRootGiBUUGenerator","Cannot open input tree."); } fInputTree->SetBranchAddress("Evt",&Evt); fInputTree->SetBranchAddress("Mul",&Mul); fInputTree->SetBranchAddress("Pid",&Pid); fInputTree->SetBranchAddress("Px",&Px); fInputTree->SetBranchAddress("Py",&Py); fInputTree->SetBranchAddress("Pz",&Pz); fInputTree->SetBranchAddress("E",&E); fInputTree->SetBranchAddress("M",&M); evNb=0; evNbFile=0; fPDG=TDatabasePDG::Instance(); } // ------------------------------------------------------------------------ // ----- Destructor --------------------------------------------------- TpcRootGiBUUGenerator::~TpcRootGiBUUGenerator() { CloseInput(); } // ------------------------------------------------------------------------ // ----- Public method ReadEvent -------------------------------------- Bool_t TpcRootGiBUUGenerator::ReadEvent(FairPrimaryGenerator* primGen) { std::vector mom; std::vector pdgs; std::vector energies; Bool_t foundParticle=false; //std::cout<<"Searching for pdgID:"<AddTrack(pdgs[itrack], mom[itrack].X(), mom[itrack].Y(),mom[itrack].Z(), vx, vy, vz,-1,true,energies[itrack]); } evNb++; return kTRUE; } }else{ return kFALSE; } } return kFALSE; } // ------------------------------------------------------------------------ Bool_t TpcRootGiBUUGenerator::ReadFileEvent(std::vector &mom, std::vector &pdgs, std::vector &energies, Bool_t &found){ mom.clear(); pdgs.clear(); energies.clear(); // Check for input file if ( ! fInputFile->IsOpen()) { cout << "-E TpcRootGiBUUGenerator: Input file not open!" << endl; return kFALSE; } // If end of input file is reached : close it and abort run if ( evNbFile>=fInputTree->GetEntriesFast() ) { cout << "-I TpcRootGiBUUGenerator: End of input file reached " << endl; CloseInput(); return kFALSE; } // Define event variable to be read from file Int_t ntracks = 0, dummyInt = 0; Double_t dummyDouble = 0., weight=0.; // Define track variables to be read from file Int_t pdgID = 0; fInputTree->GetEvent(evNbFile); ntracks=Mul; realEvNbFile=Evt; for (Int_t itrack=0; itrackGetEvent(evNbFile+itrack); if(Pid==60221){ continue; } double pdgMass= fPDG->GetParticle(Pid)->Mass(); if(Pid==pdgToFilter){ found=kTRUE; } TVector3 p(Px,Py,Pz); E=TMath::Sqrt(M*M+p.Mag2()); mom.push_back(p); pdgs.push_back(Pid); energies.push_back(E); } evNbFile+=ntracks; return kTRUE; } // ----- Private method CloseInput ------------------------------------ void TpcRootGiBUUGenerator::CloseInput() { if ( fInputFile ) { if ( fInputFile->IsOpen() ) { cout << "-I TpcRootGiBUUGenerator: Closing input file " << fFileName << endl; fInputFile->Close(); } delete fInputFile; fInputFile = NULL; } } // ------------------------------------------------------------------------ void TpcRootGiBUUGenerator::GenerateVertex(double& vx,double& vy,double& vz){ //TODO configurable target double v_phi=gRandom->Uniform(0,2*TMath::Pi()); double v_r=gRandom->Uniform(0,vetoR); TVector2 a; a.SetMagPhi(v_r,v_phi); vx=a.X(); vy=a.Y(); vz=gRandom->Uniform(-targetThickness/2,targetThickness/2); } bool TpcRootGiBUUGenerator::SetStartEvent(unsigned int eventNbr){ std::cout<<"********************************** TpcRootGiBUUGenerator::SetStartEvent( "<GetEntriesFast()){ // Define event variable to be read from file // Read event header line from input file fInputTree->GetEvent(evNbFile); std::string buffer=""; int mulCopy=Mul; evNbFile=evNbFile+Mul; realEvNbFile=Evt; evNb=Evt; }else{ std::cout<<"end of file reached while searching for start event"<