FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndCorrDistGenerator Class Reference

#include <PndCorrDistGenerator.h>

Inheritance diagram for PndCorrDistGenerator:

Public Member Functions

 PndCorrDistGenerator ()
 
 PndCorrDistGenerator (const Char_t *fileName)
 
virtual ~PndCorrDistGenerator ()
 
virtual Bool_t ReadEvent (FairPrimaryGenerator *primGen)
 
void SetParam ()
 
void SetThetaRange (Double_t thetLow=0., Double_t thetHigh=0.)
 

Private Member Functions

void CloseInput ()
 
Double_t MaxBoltDistP (Double_t MeanP)
 
Double_t MeanMomentum (Double_t thet)
 
 ClassDef (PndCorrDistGenerator, 3)
 

Private Attributes

Int_t iEvent
 
const Char_t * fFileName
 Event number. More...
 
TFile * fInputFile
 Input file name. More...
 
TH2F * fInputHist
 Pointer to input file. More...
 
TClonesArray * fParticles
 Pointer to input histogramm 2D. More...
 
Int_t fPdgType
 Particle array from PLUTO. More...
 
Bool_t fParam
 
Double_t fTheLow
 
Double_t fTheHigh
 

Detailed Description

Definition at line 31 of file PndCorrDistGenerator.h.

Constructor & Destructor Documentation

PndCorrDistGenerator::PndCorrDistGenerator ( )

Default constructor (should not be used)

Definition at line 24 of file PndCorrDistGenerator.cxx.

24  {
25  iEvent = 0;
26  fInputFile = NULL;
27  fInputHist = NULL;
28  fPdgType = 0;
29  fParam = kFALSE;
30  fTheLow =0.;
31  fTheHigh =0.;
32 }
TH2F * fInputHist
Pointer to input file.
Int_t fPdgType
Particle array from PLUTO.
TFile * fInputFile
Input file name.
PndCorrDistGenerator::PndCorrDistGenerator ( const Char_t *  fileName)

Standard constructor

Parameters
fileNameThe input (PLUTO) file name

Definition at line 38 of file PndCorrDistGenerator.cxx.

38  {
39  iEvent = 0;
40  fFileName = fileName;
41  fInputFile = new TFile(fFileName);
42  if(fParam){cout<<"choosing parametrization modell"<<endl;}
43  else{
44  cout<<" Random events gerated from a TH2D P(th) dist-->GiBUU Modell"<<endl;
45  fInputHist = (TH2F*) fInputFile->Get("source");
46  }
47 
48 }
TH2F * fInputHist
Pointer to input file.
const Char_t * fFileName
Event number.
TFile * fInputFile
Input file name.
PndCorrDistGenerator::~PndCorrDistGenerator ( )
virtual

Destructor

Definition at line 53 of file PndCorrDistGenerator.cxx.

53  {
54  CloseInput();
55 }

Member Function Documentation

PndCorrDistGenerator::ClassDef ( PndCorrDistGenerator  ,
 
)
private
void PndCorrDistGenerator::CloseInput ( )
private

Private method CloseInput. Just for convenience. Closes the input file properly. Called from destructor and from ReadEvent.

Definition at line 146 of file PndCorrDistGenerator.cxx.

146  {
147  if ( fInputFile ) {
148  cout << "-I PndCorrDistGenerator: Closing input file " << fFileName
149  << endl;
150  fInputFile->Close();
151  delete fInputFile;
152  }
153  fInputFile = NULL;
154 }
const Char_t * fFileName
Event number.
TFile * fInputFile
Input file name.
Double_t PndCorrDistGenerator::MaxBoltDistP ( Double_t  MeanP)
private

Definition at line 166 of file PndCorrDistGenerator.cxx.

References Double_t, sigma, CAMath::Sqrt(), and v.

166  {
167 
168  Double_t Et=0.,EK=0.,v=0.,gamma=0.,zet=0.,sigma=0.;
169 
170  Et = TMath::Sqrt((MeanP*MeanP)+ (1.321*1.321)); // rel. energy, speed of light c = 1.
171  v = MeanP/Et; // v = beta*c/Etotal
172  gamma = 1./TMath::Sqrt(1.-(v*v)); // gamma factor
173  EK = gamma*1.321 - 1.321; //rel. kinetic energy, mass of Xi_minus = 1.321 GeV/c2
174  zet = EK/1.321;
175 
176 
177  sigma = TMath::Sqrt(3*zet);
178  Double_t limLow,limHigh;
179  limLow = MeanP*0.882;
180  limHigh = MeanP*1.085;
181  //Double_t beta =
182 
183  TF1* MBol = new TF1("mb","TMath::Gaus(x,[0],[1])",limLow,limHigh);
184  // TF1* MBol = new TF1("mb","TMath::Exp(-TMath::Sqrt(1.+ (TMath::Power(x/[0],2))))/(4*TMath::Power([0],3)*TMath::Pi()*[1]*TMath::BesselK(2,1/[1]))",-2.,3);
185  //TF1* MBol = new TF1("mb","x*x*[1]*TMath::Exp(-x/[0])/([0]*TMath::BesselK(2,1/[0]))",0.,3);
186  //MBol->SetParameter(0,1.321);
187  //MBol->SetParameter(1,zet);
188  // MBol->SetParameter(0,zet);
189  //MBol->SetParameter(1,v);
190  MBol->SetParameter(0,MeanP);
191  MBol->SetParameter(1,sigma);
192 
193  cout << "-I PndCorrDistGenerator: Calculating Momentum distribution via MB <v >" << endl;
194  Double_t fitval = MBol->GetRandom();
195 
196  return fitval;
197 
198 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t sigma[nsteps]
Definition: dedx_bands.C:65
__m128 v
Definition: P4_F32vec4.h:4
Double_t
Double_t PndCorrDistGenerator::MeanMomentum ( Double_t  thet)
private

Definition at line 157 of file PndCorrDistGenerator.cxx.

References Double_t.

157  {
158 
159  cout << "-I PndCorrDistGenerator: Calculating Mean Momentum as function of theta" << endl;
160  //f(x) = exp(p0+p1*x)
161  Double_t fitval = TMath::Exp((4.93964e-01)-(thet*1.41427e-02));
162  return fitval;
163 
164 }
Double_t
Bool_t PndCorrDistGenerator::ReadEvent ( FairPrimaryGenerator *  primGen)
virtual

Reads on event from the input file and pushes the tracks onto the stack. Abstract method in base class.

Parameters
primGenpointer to the FairPrimaryGenerator

Definition at line 60 of file PndCorrDistGenerator.cxx.

References CAMath::Cos(), Double_t, exp(), p, phi, pz, CAMath::Sin(), and theta.

60  {
61 
62  // Check for input file
63  if ( ! fInputFile ) {
64  cout << "-E PndCorrDistGenerator: Input file nor open!" << endl;
65  return kFALSE;
66  }
67 
68  //defining kinematics variables, angle is expressed in degrees.
69 
70  Double_t phi, p=0.0,theta=0.0;
71  Double_t pz=0.,py=0.,px=0.;
72 
73 
74  phi = gRandom->Uniform(0,360)* TMath::DegToRad();
75 
76  TF1 *lan;//Landau
77  TF1 *exp; // exponential
78  TF1 *MBol;
79  // relativistic Maxwell Boltzmann dist <v> (Maxwell-Juettner dist), mean velocity
80  // For more INFO see, http://en.wikipedia.org/wiki/Maxwell–Jüttner_distribution
81 
82  Double_t mean_p=0.;
83 
84  if(fParam){
85 
86  //--> Landau, dsigma/dtheta dist
87 
88  //Constant 1.03980e+03 //Peak height (kFALSE-->0), no integral,
89  //MPV 1.10093e+01
90  //Sigma 4.65284e+00
91 
92  lan = new TF1("lan","1.03980e+03*TMath::Landau(x,[0],[1],0)",fTheLow,fTheHigh);
93  lan->SetParameters(1.10093e+01,4.65284e+00);
94 
95  // getting theta
96  theta = lan->GetRandom();
97 
98  // <p> --- > Max. Boltzmann dist --> p
99 
100  //more or less -->exponential <P>(theta)
101  //Constant 4.93964e-01
102  // Slope -1.41427e-02
103 
104 
105  mean_p = MeanMomentum(theta);// getting <p>
106  p = MaxBoltDistP(mean_p);
107 
108  pz = p*TMath::Cos(theta* TMath::DegToRad());
109  py = p*TMath::Sin(theta* TMath::DegToRad())*TMath::Sin(phi);
110  px = p*TMath::Sin(theta* TMath::DegToRad())*TMath::Cos(phi);
111 
112  }else{
113  //getting p and theta randomly from a TH2F object
114  //both variables are correlated to each other
115  //phi is independtly generated over the full angular range.
116 
117  fInputHist->GetRandom2(theta,p);
118 
119 
120 
121  pz = p*TMath::Cos(theta* TMath::DegToRad());
122  py = p*TMath::Sin(theta* TMath::DegToRad())*TMath::Sin(phi);
123  px = p*TMath::Sin(theta* TMath::DegToRad())*TMath::Cos(phi);
124  }
125  // vertex coordinates given during generation via prim generator
126 
127  // Give track to PrimaryGenerator
128  // multiplocity equal to 1
129  // pdg is -3312, corresponding to a Xi minus
130 
131  fPdgType==3312;
132 
133  //cout<<" generating Xi minus "<<px<<" "<<py<<" "<<pz<<endl;
134 
135  primGen->AddTrack(3312, px, py, pz, 0., 0., 0.);
136 
137 
138 
139  return kTRUE;
140 
141 }
Double_t p
Definition: anasim.C:58
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
Double_t MeanMomentum(Double_t thet)
static T Sin(const T &x)
Definition: PndCAMath.h:42
TH2F * fInputHist
Pointer to input file.
static T Cos(const T &x)
Definition: PndCAMath.h:43
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
Int_t fPdgType
Particle array from PLUTO.
Double_t
TFile * fInputFile
Input file name.
Double_t MaxBoltDistP(Double_t MeanP)
double pz[39]
Definition: pipisigmas.h:14
void PndCorrDistGenerator::SetParam ( )
inline

Definition at line 55 of file PndCorrDistGenerator.h.

References fParam.

Referenced by run(), run_test(), runSimHF_GiB_DC(), and sim_pid_runSimHF_GiB_DC().

55 {fParam=kTRUE;};
void PndCorrDistGenerator::SetThetaRange ( Double_t  thetLow = 0.,
Double_t  thetHigh = 0. 
)
inline

Definition at line 56 of file PndCorrDistGenerator.h.

Referenced by run(), run_test(), runSimHF_GiB_DC(), and sim_pid_runSimHF_GiB_DC().

57  {fTheLow = thetLow; fTheHigh = thetHigh;};

Member Data Documentation

const Char_t* PndCorrDistGenerator::fFileName
private

Event number.

Definition at line 63 of file PndCorrDistGenerator.h.

TFile* PndCorrDistGenerator::fInputFile
private

Input file name.

Definition at line 64 of file PndCorrDistGenerator.h.

TH2F* PndCorrDistGenerator::fInputHist
private

Pointer to input file.

Definition at line 65 of file PndCorrDistGenerator.h.

Bool_t PndCorrDistGenerator::fParam
private

Definition at line 68 of file PndCorrDistGenerator.h.

Referenced by SetParam().

TClonesArray* PndCorrDistGenerator::fParticles
private

Pointer to input histogramm 2D.

Definition at line 66 of file PndCorrDistGenerator.h.

Int_t PndCorrDistGenerator::fPdgType
private

Particle array from PLUTO.

Definition at line 67 of file PndCorrDistGenerator.h.

Double_t PndCorrDistGenerator::fTheHigh
private

Definition at line 69 of file PndCorrDistGenerator.h.

Double_t PndCorrDistGenerator::fTheLow
private

Definition at line 69 of file PndCorrDistGenerator.h.

Int_t PndCorrDistGenerator::iEvent
private

Definition at line 57 of file PndCorrDistGenerator.h.


The documentation for this class was generated from the following files: