X-Git-Url: http://git.uio.no/git/?a=blobdiff_plain;f=EVGEN%2FAliGenSlowNucleons.cxx;h=56033125638a1e2a0c883f6d9b2b4bace59f8a0d;hb=b733a8b0eafbb2f4dd2bf2c66280b6bd8c1bf69d;hp=7b9a52b78507c476cb17a748120e798e2491e1e6;hpb=1ffc7bfa00d00db9de5ebc83dc14c76972a52637;p=u%2Fmrichter%2FAliRoot.git diff --git a/EVGEN/AliGenSlowNucleons.cxx b/EVGEN/AliGenSlowNucleons.cxx index 7b9a52b7850..56033125638 100644 --- a/EVGEN/AliGenSlowNucleons.cxx +++ b/EVGEN/AliGenSlowNucleons.cxx @@ -30,8 +30,13 @@ #include #include #include +#include +#include "AliConst.h" #include "AliCollisionGeometry.h" +#include "AliStack.h" +#include "AliRun.h" +#include "AliMC.h" #include "AliGenSlowNucleons.h" #include "AliSlowNucleonModel.h" @@ -60,6 +65,9 @@ AliGenSlowNucleons::AliGenSlowNucleons() fThetaDistribution(), fCosThetaGrayHist(), fCosTheta(), + fBeamCrossingAngle(0.), + fBeamDivergence(0.), + fBeamDivEvent(0.), fSlowNucleonModel(0) { // Default constructor @@ -88,7 +96,11 @@ AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart) fThetaDistribution(), fCosThetaGrayHist(), fCosTheta(), + fBeamCrossingAngle(0.), + fBeamDivergence(0.), + fBeamDivEvent(0.), fSlowNucleonModel(new AliSlowNucleonModel()) + { // Constructor fName = "SlowNucleons"; @@ -130,6 +142,8 @@ void AliGenSlowNucleons::Init() "(2./3.14159265358979312)/(exp(2./3.14159265358979312)-exp(-2./3.14159265358979312))*exp(2.*x/3.14159265358979312)", -1., 1.); } + + printf("\n AliGenSlowNucleons: applying crossing angle %f mrad to slow nucleons\n",fBeamCrossingAngle*1000.); } void AliGenSlowNucleons::FinishRun() @@ -140,7 +154,7 @@ void AliGenSlowNucleons::FinishRun() TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700); c->Divide(2,1); c->cd(1); - fDebugHist1->Draw(); + fDebugHist1->Draw("colz"); c->cd(2); fDebugHist2->Draw(); c->cd(3); @@ -164,12 +178,14 @@ void AliGenSlowNucleons::Generate() // Int_t nnw = fCollisionGeometry->NNw(); // Int_t nwnw = fCollisionGeometry->NwNw(); - fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn); + //fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn); + fSlowNucleonModel->GetNumberOfSlowNucleons2(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn); if (fDebug) { //printf("Collision Geometry %f %d %d %d %d\n", b, nn, nwn, nnw, nwnw); printf("Slow nucleons: %d grayp %d grayn %d blackp %d blackn \n", fNgp, fNgn, fNbp, fNbn); - fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NwN(), 1.); + fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NNw(), 1.); fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.); + } } @@ -180,7 +196,12 @@ void AliGenSlowNucleons::Generate() Float_t polar [3] = {0., 0., 0.}; Int_t nt, i, j; Int_t kf; - + + // Extracting 1 value per event for the divergence angle + Double_t rvec = gRandom->Gaus(0.0, 1.0); + fBeamDivEvent = fBeamDivergence * TMath::Abs(rvec); + printf("\n AliGenSlowNucleons: applying beam divergence %f mrad to slow nucleons\n",fBeamDivEvent*1000.); + if(fVertexSmear == kPerEvent) { Vertex(); for (j=0; j < 3; j++) origin[j] = fVertex[j]; @@ -195,8 +216,9 @@ void AliGenSlowNucleons::Generate() GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta); if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta)); PushTrack(fTrackIt, -1, kf, p, origin, polar, - time, kPNoProcess, nt, 1.,1); + time, kPNoProcess, nt, 1.,-2); KeepTrack(nt); + SetProcessID(nt,kGrayProcess); } // // Gray neutrons @@ -207,8 +229,9 @@ void AliGenSlowNucleons::Generate() GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta); if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta)); PushTrack(fTrackIt, -1, kf, p, origin, polar, - time, kPNoProcess, nt, 1.,1); + time, kPNoProcess, nt, 1.,-2); KeepTrack(nt); + SetProcessID(nt,kGrayProcess); } // // Black protons @@ -218,8 +241,9 @@ void AliGenSlowNucleons::Generate() for(i = 0; i < fNbp; i++) { GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta); PushTrack(fTrackIt, -1, kf, p, origin, polar, - time, kPNoProcess, nt, 1.,1); + time, kPNoProcess, nt, 1.,-1); KeepTrack(nt); + SetProcessID(nt,kBlackProcess); } // // Black neutrons @@ -229,8 +253,9 @@ void AliGenSlowNucleons::Generate() for(i = 0; i < fNbn; i++) { GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta); PushTrack(fTrackIt, -1, kf, p, origin, polar, - time, kPNoProcess, nt, 1.,1); + time, kPNoProcess, nt, 1.,-1); KeepTrack(nt); + SetProcessID(nt,kBlackProcess); } } @@ -288,16 +313,22 @@ void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, q[0] = p * TMath::Sin(theta) * TMath::Cos(phi); q[1] = p * TMath::Sin(theta) * TMath::Sin(phi); q[2] = p * TMath::Cos(theta); + //if(fDebug==1) printf("\n Momentum in RS of the moving source: p = (%f, %f, %f)\n",q[0],q[1],q[2]); /* Transform to system of the target nucleus */ /* beta is passed as negative, because the gray nucleons are slowed down */ Lorentz(m, -beta, q); - + //if(fDebug==1) printf(" Momentum in RS of the target nucleus: p = (%f, %f, %f)\n",q[0],q[1],q[2]); + /* Transform to laboratory system */ Lorentz(m, fBeta, q); q[2] *= fProtonDirection; - + if(fDebug==1)printf("\n Momentum after LHC boost: p = (%f, %f, %f)\n",q[0],q[1],q[2]); + + if(fBeamCrossingAngle>0.) BeamCrossDivergence(1, q); // applying crossing angle + if(fBeamDivergence>0.) BeamCrossDivergence(2, q); // applying divergence + } Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T) @@ -309,6 +340,7 @@ Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T) } +//_____________________________________________________________________________ void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q) { /* Lorentz transform in the direction of q[2] */ @@ -319,3 +351,93 @@ void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q) //printf(" \t beta %1.10f gamma %f energy %f -> p_z = %f\n",beta, gamma, energy,q[2]); } +//_____________________________________________________________________________ +void AliGenSlowNucleons::BeamCrossDivergence(Int_t iwhat, Float_t *pLab) +{ + // Applying beam divergence and crossing angle + // + Double_t pmod = TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]+pLab[2]*pLab[2]); + + Double_t tetdiv = 0.; + Double_t fidiv = 0.; + if(iwhat==1){ + tetdiv = fBeamCrossingAngle; + fidiv = k2PI/4.; + } + else if(iwhat==2){ + tetdiv = fBeamDivEvent; + fidiv = (gRandom->Rndm())*k2PI; + } + + Double_t tetpart = TMath::ATan2(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]), pLab[2]); + Double_t fipart=0.; + if(TMath::Abs(pLab[1])>0. || TMath::Abs(pLab[0])>0.) fipart = TMath::ATan2(pLab[1], pLab[0]); + if(fipart<0.) {fipart = fipart+k2PI;} + tetdiv = tetdiv*kRaddeg; + fidiv = fidiv*kRaddeg; + tetpart = tetpart*kRaddeg; + fipart = fipart*kRaddeg; + + Double_t angleSum[2]={0., 0.}; + AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum); + + Double_t tetsum = angleSum[0]; + Double_t fisum = angleSum[1]; + //printf("tetpart %f fipart %f tetdiv %f fidiv %f angleSum %f %f\n",tetpart,fipart,tetdiv,fidiv,angleSum[0],angleSum[1]); + tetsum = tetsum*kDegrad; + fisum = fisum*kDegrad; + + pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum); + pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum); + pLab[2] = pmod*TMath::Cos(tetsum); + if(fDebug==1){ + if(iwhat==1) printf(" Beam crossing angle %f mrad ", fBeamCrossingAngle*1000.); + else if(iwhat==2) printf(" Beam divergence %f mrad ", fBeamDivEvent*1000.); + printf(" p = (%f, %f, %f)\n",pLab[0],pLab[1],pLab[2]); + } +} + +//_____________________________________________________________________________ +void AliGenSlowNucleons::AddAngle(Double_t theta1, Double_t phi1, + Double_t theta2, Double_t phi2, Double_t *angleSum) +{ + // Calculating the sum of 2 angles + Double_t temp = -1.; + Double_t conv = 180./TMath::ACos(temp); + + Double_t ct1 = TMath::Cos(theta1/conv); + Double_t st1 = TMath::Sin(theta1/conv); + Double_t cp1 = TMath::Cos(phi1/conv); + Double_t sp1 = TMath::Sin(phi1/conv); + Double_t ct2 = TMath::Cos(theta2/conv); + Double_t st2 = TMath::Sin(theta2/conv); + Double_t cp2 = TMath::Cos(phi2/conv); + Double_t sp2 = TMath::Sin(phi2/conv); + Double_t cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2; + Double_t cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2; + Double_t cz = ct1*ct2-st1*st2*cp2; + + Double_t rtetsum = TMath::ACos(cz); + Double_t tetsum = conv*rtetsum; + if(TMath::Abs(tetsum)<1e-4 || tetsum==180.) return; + + temp = cx/TMath::Sin(rtetsum); + if(temp>1.) temp=1.; + if(temp<-1.) temp=-1.; + Double_t fisum = conv*TMath::ACos(temp); + if(cy<0) {fisum = 360.-fisum;} + + angleSum[0] = tetsum; + angleSum[1] = fisum; +} + +//_____________________________________________________________________________ +void AliGenSlowNucleons::SetProcessID(Int_t nt, UInt_t process) +{ + // Tag the particle as + // gray or black + if (fStack) + fStack->Particle(nt)->SetUniqueID(process); + else + gAlice->GetMCApp()->Particle(nt)->SetUniqueID(process); +}