// ------------------------------------------------------------------------- // ----- PndEvtGenDirect source file ----- // ----- ----- // ------------------------------------------------------------------------- #include #include "TClonesArray.h" #include "TFile.h" #include "TLorentzVector.h" #include "TTree.h" #include "TVector3.h" #include "TParticle.h" #include "FairPrimaryGenerator.h" #include "PndEvtGenDirect.h" #include "TRandom.h" #include "EvtGenBase/EvtPatches.hh" #include "EvtGenBase/EvtPatches.hh" #include "EvtGenBase/EvtParticleFactory.hh" #include "EvtGenBase/EvtStdHep.hh" #include "EvtGenBase/EvtParticle.hh" #include "EvtGenBase/EvtRandom.hh" #include "EvtGenBase/EvtRandomEngine.hh" //#include "EvtGenBase/EvtReport.hh" #include #include #include using namespace std; using std::endl; using std::ofstream; using std::cout; #include "TRandom3.h" //Define random number fcn used by Jetset extern "C" { extern float rlu_(); extern float begran_(int *); } float rlu_(){ return EvtRandom::Flat(); } float begran_(int *){ return EvtRandom::Flat(); } //define class for generating random nubers class EvtRootRandomEngine:public EvtRandomEngine{ public: EvtRootRandomEngine(int s=0) {seed=s;} double random(); int seed; }; double EvtRootRandomEngine::random(){ static TRandom3 randengine(seed); return randengine.Rndm(); } // ----- Default constructor ------------------------------------------ PndEvtGenDirect::PndEvtGenDirect() { SetName("PndEvtGenDirect"); fStoreTree=false; verbose=0; } // ------------------------------------------------------------------------ // ----- Standard constructor ----------------------------------------- PndEvtGenDirect::PndEvtGenDirect(TString particle,TString decfile,Double_t Mom, Long_t Seed,TString defaultDECAY,TString defaultPDL) { PndEvtGenDirect(); cout << " PndEvtGenDirect"< Particle: "< decfile: "<0) cout << " pbar-Momentum: "< Momentum: "< CMS energy: "< Rnd Seed: "<GetSeed(); cout << " Rnd Seed changed to " << Seed << endl; } myRandomEngine=new EvtRootRandomEngine(Seed); myGenerator=new EvtGen(defaultDECAY,defaultPDL,myRandomEngine); //If I wanted a user decay file, I would read it in now. if(decfile!="") myGenerator->readUDecay(decfile.Data()); PART=EvtPDL::getId(std::string(particle.Data())); if ( (particle.Contains("pbarp") || particle.Contains("pbard")) && Mom==0) { cerr <<"\033[5m\033[31m -E ****** FATAL ERROR: is '" << particle.Data() << "'; MUST give pbar momentum or cms energy!\033[0m"<0){ fMomentum = val; if ( particle.Contains("pbarpSystem") ) fEnergy = mp+sqrt(fMomentum*fMomentum+mp*mp); if ( particle.Contains("pbardSystem") ) fEnergy = md+sqrt(fMomentum*fMomentum+mp*mp); } else //val is -E_cm { val=-val; fEnergy = val*val/(2*mp); fMomentum = sqrt(fEnergy*fEnergy-val*val); } cout <<"\n############# Generating with following conditions:\n\n"; cout <<"incident 4-mom : ("<setDiagonalSpinDensity(); // Generate the event myGenerator->generateDecay(parent); // Write out the results evtstdhep.init(); parent->makeStdHep(evtstdhep); Bool_t plotflag; plotflag=false; //print out some status info if (verbose>1 ||(verbose==1 && (evtnr<10 || ((evtnr+1)%100)==0))){ cout << "PndEvtGenDirect::ReadEvent "<printParticle(); // report(INFO,"EvtGen") << "event Number\t"<< evtnr << evtstdhep << endl; cout << "event Number\t"<< evtnr << evtstdhep << endl; cout << evtnr << "\t" << evtstdhep.getNPart(); cout < cm conversion fY=vxyz.get(2)/10.; // mm -> cm conversion fZ=vxyz.get(3)/10.; // mm -> cm conversion fE=pxyz.get(0); Px=pxyz.get(1); Py=pxyz.get(2); Pz=pxyz.get(3); if(plotflag) printf("- I -: new particle at: %f, %f, %f (%f)-> %f %f %f (%f) ID %d ##Daughters %d %d Mothers %d %d\n", fX, fY, fZ, fT,Px, Py, Pz, fE, Id, nFD, nLD,evtstdhep.getFirstMother(i),evtstdhep.getLastMother(i)); if(fStoreTree){ primGen->AddTrack(Id, Px, Py, Pz, fX, fY, fZ, evtstdhep.getFirstMother(i),(nFD==-1 && nLD==-1),fE); }else{ primGen->AddTrack(Id, Px, Py, Pz, fX, fY, fZ);// default -1, true } } } if(plotflag) cout <<"==== compare end ==="<deleteTree(); evtnr++; return kTRUE; } // ------------------------------------------------------------------------ ClassImp(PndEvtGenDirect)