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

#include <PndDpmDirect.h>

Inheritance diagram for PndDpmDirect:

Public Member Functions

 PndDpmDirect ()
 
 PndDpmDirect (Double_t Mom, Int_t Mode, Long_t Seed=-1)
 
 PndDpmDirect (Double_t Mom, Int_t Mode, Long_t Seed, Double_t ThtMin)
 
 PndDpmDirect (Double_t Mom, Int_t Mode, Double_t Rsigma, TF1 *DensityFunction, Long_t Seed=-1, Double_t ThtMin=0.001)
 
virtual ~PndDpmDirect ()
 
virtual Bool_t ReadEvent (FairPrimaryGenerator *primGen)
 
void SetUnstable (int pdg)
 
void SetStable (int pdg)
 

Private Member Functions

 ClassDef (PndDpmDirect, 1)
 

Private Attributes

double fMom
 
double fMode
 
double fSeed
 
int fGasmode
 
double fRsigma
 
double fThtMin
 
TF1 * fDensityFunction
 

Detailed Description

Definition at line 25 of file PndDpmDirect.h.

Constructor & Destructor Documentation

PndDpmDirect::PndDpmDirect ( )

Default constructor (should not be used)

Definition at line 45 of file PndDpmDirect.cxx.

45  {
46 
47 }
PndDpmDirect::PndDpmDirect ( Double_t  Mom,
Int_t  Mode,
Long_t  Seed = -1 
)

Standard constructor

Parameters
Momin GeV/C
Mode= 0. - No elastic scattering, only inelastic
Mode= 1. - Elastic and inelastic interactions
Mode= 2. - Only elastic scattering, no inelastic one

Definition at line 52 of file PndDpmDirect.cxx.

References Double_t, and CAMath::Log().

52  {
53 
54 //
55 // Calculate ThtMin first. For this we make a cut-off on the value of -t of 1e-2 GeV^2 (~100 MeV/c momentum)
56 // This estimated from a parametrization found in thesis of Thomas Wuerschig (figure 6.4, page 121):
57 // Roughly: 0.4 deg at 15 GeV/c and 4 deg at 1.5 GeV/c, lineair interpolation in double log-scale.
58 //
59 
60  Double_t logangle = TMath::Log(0.4)+(TMath::Log(15.)-TMath::Log(Mom))*(TMath::Log(4)-TMath::Log(0.4))/(TMath::Log(15)-TMath::Log(1.5));
61  Double_t ThtMin = TMath::Exp(logangle);
62 
63  PndDpmDirect(Mom, Mode, Seed, ThtMin) ;
64 }
Double_t
static T Log(const T &x)
Definition: PndCAMath.h:40
PndDpmDirect::PndDpmDirect ( Double_t  Mom,
Int_t  Mode,
Long_t  Seed,
Double_t  ThtMin 
)

Definition at line 66 of file PndDpmDirect.cxx.

References a, fDensityFunction, fSeed, and init1_().

66  {
67  fMom = Mom;
68  fMode = Mode;
69  fSeed = Seed;
70  fThtMin = ThtMin;
71 
72  if (fSeed < 0)
73  {
74  Long_t iSeed = gRandom->GetSeed();
75  int a = iSeed/100000;
76  fSeed=iSeed - a*100000 + a/100000.;
77  }
78 
79  fGasmode = 0;
80  fRsigma = 0.;
81  fDensityFunction = new TF1();
82 
83  cout << "<I> PndDpmDirect initialization" << endl;
84  cout << "<I> Momentum = " << fMom << endl;
85  cout << "<I> Seed = " << fSeed << endl;
86  cout << "<I> Mode = " << fMode << endl;
87  cout << "<I> Theta min = " << fThtMin<<endl;
88 
89  init1_(&fMom,&fSeed,&fMode, &fThtMin); // init the DPM generator
90 }
double fMode
Definition: PndDpmDirect.h:67
double fThtMin
Definition: PndDpmDirect.h:71
Int_t a
Definition: anaLmdDigi.C:126
TF1 * fDensityFunction
Definition: PndDpmDirect.h:73
int init1_(double *Plab, double *seed, double *Elastic, double *tetmin)
double fSeed
Definition: PndDpmDirect.h:68
double fRsigma
Definition: PndDpmDirect.h:70
PndDpmDirect::PndDpmDirect ( Double_t  Mom,
Int_t  Mode,
Double_t  Rsigma,
TF1 *  DensityFunction,
Long_t  Seed = -1,
Double_t  ThtMin = 0.001 
)

Definition at line 94 of file PndDpmDirect.cxx.

References a, fDensityFunction, fSeed, and init1_().

94  {
95  fMom = Mom;
96  fMode = Mode;
97  fSeed = Seed;
98  fThtMin = ThtMin;
99 
100  if (fSeed < 0)
101  {
102  Long_t iSeed = gRandom->GetSeed();
103  int a = iSeed/100000;
104  fSeed=iSeed - a*100000 + a/100000.;
105  }
106 
107  fGasmode = 1;
108  fRsigma = Rsigma;
109  fDensityFunction = DensityFunction;
110 
111  cout << "<I> PndDpmDirect initialization" << endl;
112  cout << "<I> Momentum = " << fMom << endl;
113  cout << "<I> Seed = " << fSeed << endl;
114  cout << "<I> Mode = " << fMode << endl;
115  cout << "<I> Gasmode = " << fGasmode << endl;
116  cout << "<I> Theta min = " << fThtMin<<endl;
117 
118  init1_(&fMom,&fSeed,&fMode, &fThtMin); // init the DPM generator
119 }
double fMode
Definition: PndDpmDirect.h:67
double fThtMin
Definition: PndDpmDirect.h:71
Int_t a
Definition: anaLmdDigi.C:126
TF1 * fDensityFunction
Definition: PndDpmDirect.h:73
int init1_(double *Plab, double *seed, double *Elastic, double *tetmin)
double fSeed
Definition: PndDpmDirect.h:68
double fRsigma
Definition: PndDpmDirect.h:70
PndDpmDirect::~PndDpmDirect ( )
virtual

Destructor

Definition at line 123 of file PndDpmDirect.cxx.

123  {
124 
125 }

Member Function Documentation

PndDpmDirect::ClassDef ( PndDpmDirect  ,
 
)
private
Bool_t PndDpmDirect::ReadEvent ( FairPrimaryGenerator *  primGen)
virtual

Generate one event using DPM

Parameters
primGenpointer to the FairPrimaryGenerator

Definition at line 143 of file PndDpmDirect.cxx.

References Double_t, dpm_gen_(), fDensityFunction, fSeed, fX, fY, fZ, and i.

143  {
144 
145  int npart, i;
146  Double_t fX, fY, fZ, radius;
147  double Px[1000],Py[1000],Pz[1000]; //,E[1000],Pm[1000],Wh[1000]; //[R.K. 01/2017] unused variables
148  int Id[1000];
149 
150  double Generator=0.; // Format in which events are produced (0=pythia, 1=pluto)
151 
152  //Double_t weight = 1.0; //[R.K. 01/2017] unused variable?
153  //Int_t activeCnt=0; //[R.K. 01/2017] unused variable?
154 
155  // run generator
156  dpm_gen_(&Generator, &fSeed);
157 
158  // Loop over all produced particles
159  npart = lujets_.n;
160 
161  for (i= 0; i< npart; ++i) {
162 
163  Id[i]=lujets_.k[i+1000];
164  Px[i]=lujets_.p[i];
165  Py[i]=lujets_.p[i+1000];
166  Pz[i]=lujets_.p[i+2000];
167  //Pm[i]=lujets_.p[i+4000];
168  //E[i]=lujets_.p[i+3000];
169  //Wh[i]=1.0;
170 
171  /* Check if fGasmode is set */
172  fX = 0.;
173  fY = 0.;
174  fZ = 0.;
175  if (fGasmode == 1) {
176 
177  // define position of track start
178  // Random 2D point in a circle of radius r (simple beamprofile)
179  radius = gRandom->Gaus(0,fRsigma);
180  gRandom->Circle(fX, fY, radius);
181 
182  // calculate fZ according to some (probability) density function of the gas
183  fZ=fDensityFunction->GetRandom();
184 
185  }
186 
187  // add track
188  //printf("- I -: new particle at: %f, %f, %f ...\n", fX, fY, fZ);
189  primGen->AddTrack(Id[i], Px[i], Py[i], Pz[i], fX, fY, fZ);
190 
191  }
192 
193  return kTRUE;
194 
195 }
Int_t i
Definition: run_full.C:25
Double_t fX
Definition: PndCaloDraw.cxx:34
Double_t fZ
Definition: PndCaloDraw.cxx:34
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
Double_t
TF1 * fDensityFunction
Definition: PndDpmDirect.h:73
Double_t fY
Definition: PndCaloDraw.cxx:34
double fSeed
Definition: PndDpmDirect.h:68
double fRsigma
Definition: PndDpmDirect.h:70
int dpm_gen_(double *Generator, double *seed)
void PndDpmDirect::SetStable ( int  pdg)

Definition at line 127 of file PndDpmDirect.cxx.

References chstatus_(), and status.

Referenced by prod_fsim().

128 {
129  int status=1;
130  chstatus_(&pdg, &status);
131 }
int chstatus_(int *iPDG, int *iStatus)
int status[10]
Definition: f_Init.h:28
void PndDpmDirect::SetUnstable ( int  pdg)

Definition at line 134 of file PndDpmDirect.cxx.

References chstatus_(), and status.

Referenced by prod_fsim(), quickfsimana(), simfast(), simfast_opt(), and tut_fastsim().

135 {
136  int status=0;
137  chstatus_(&pdg, &status);
138 }
int chstatus_(int *iPDG, int *iStatus)
int status[10]
Definition: f_Init.h:28

Member Data Documentation

TF1* PndDpmDirect::fDensityFunction
private

Definition at line 73 of file PndDpmDirect.h.

int PndDpmDirect::fGasmode
private

Definition at line 69 of file PndDpmDirect.h.

double PndDpmDirect::fMode
private

0. - No elastic scattering, only inelastic

  1. - Elastic and inelastic interactions
  2. - Only elastic scattering, no inelastic one

Definition at line 67 of file PndDpmDirect.h.

double PndDpmDirect::fMom
private

P_lab(GeV/c)

Definition at line 61 of file PndDpmDirect.h.

double PndDpmDirect::fRsigma
private

Definition at line 70 of file PndDpmDirect.h.

double PndDpmDirect::fSeed
private

Definition at line 68 of file PndDpmDirect.h.

double PndDpmDirect::fThtMin
private

Definition at line 71 of file PndDpmDirect.h.


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