]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TMEVSIM/TMevSim.cxx
6832949a89e074b2007a94e51feba17f6a27bf8e
[u/mrichter/AliRoot.git] / TMEVSIM / TMevSim.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.1  2001/03/25 10:15:23  morsch
19 Root interface to MevSim code as TGenerator realisation (Sylwester Radomski et al.)
20
21 */
22
23 ////////////////////////////////////////////////////////////////////////////
24 // 
25 // TMevSim
26 //
27 // TMevSim is an interface class between the event generator MEVSIM and 
28 // the ROOT system. The current implementation is based on the 6.11.2000
29 // version provided by Lanny Ray on /afs/cern.ch/user/y/yiota/ray/mult_gen.
30 // 
31 // Authors of MEVSIM:
32 // For The STAR Collaboration
33 //
34 //  Lanny Ray      
35 //     Dept. of Physics
36 //     The University of Texas at Austin
37 //     Austin, Texas 78712
38 //     (512) 471-6107
39 //     ray@physics.utexas.edu
40 //
41 //  Ron Longacre   email:
42 //
43 //
44 ////////////////////////////////////////////////////////////////////////////
45 //
46 // I.    OVERVIEW
47 //
48 //    This code is intended to provide a quick means of producing
49 //    uncorrelated simulated events for event-by-event studies,
50 //    detector acceptance and efficiency studies, etc.  The
51 //    user selects the number of events, the one-particle distribution
52 //    model, the Geant particles to include, the ranges in transverse
53 //    momentum, pseudorapidity and azimuthal angle, the mean
54 //    multiplicity for each particle type for the event run, the
55 //    mean temperature, Rapidity width, etc., and the standard deviations
56 //    for the event-to-event variation in the model parameters.
57 //    Note that these events are produced in the c.m. frame only.
58 //    
59 //    Anisotropic flow may also be simulated by introducing explicit
60 //    phi-dependence (azimuthal angle) in the particle distributions.  
61 //    The assumed model is taken from Poskanzer and Voloshin, Phys. Rev.
62 //    C58, 1671 (1998), Eq.(1), where we use,
63 //
64 //         E d^3N/dp^3 = (1/2*pi*pt)*[d^2N/dpt*dy]
65 //            * [1 + SUM(n=1,nflowterms){2*Vn*cos[n(phi-PSIr)]}]
66 //
67 //    with up to 'nflowterms' (currently set to 6, see file
68 //    Parameter_values.inc) Fourier components allowed.  Vn are
69 //    coefficients and PSIr is the reaction plane angle.
70 //    The algebraic signs of the Vn terms for n=odd are reversed
71 //    from their input values for particles with rapidity (y) < 0
72 //    as suggested in Poskanzer and Voloshin.
73 //    The flow parameters can depend on pt and rapidity (y) according
74 //    to the model suggested by Art Poskanzer (Feb. 2000) and as
75 //    defined in the Function Vn_pt_y.
76 //
77 //    The user may select either to have the same multiplicity per
78 //    particle type for each event or to let the multiplicity vary
79 //    randomly according to a Poisson distribution. In addition, an
80 //    overall multiplicative scale factor can be applied to each
81 //    particle ID's multiplicity (same factor applied to each PID).
82 //    This scaling can vary randomly according to a Gaussian from
83 //    event-to-event.  This is to simulate trigger acceptance
84 //    fluctuations.  Similarly the
85 //    parameters of the one-particle distribution models may either
86 //    be fixed to the same value for each event or allowed to randomly
87 //    vary about a specified mean with a specified standard deviation
88 //    and assuming a Gaussian distribution.
89 //
90 //    With respect to the reaction plane and anisotropic flow simulation,
91 //    the user may select among four options:
92 //       (1) ignore reaction plane and anisotropic flow effects; all
93 //           distributions will be azimuthally invariant, on average.
94 //       (2) assume a fixed reaction plane angle, PSIr, for all events
95 //           in the run.
96 //       (3) assume a Gaussian distribution with specified mean and
97 //           standard deviation for the reaction plane angles for the
98 //           events in the run.  PSIr is randomly determined for each
99 //           event.
100 //       (4) assume uniformly distributed, random values for the reaction  
101 //           plane angles from 0 to 360 deg., for each event in the run.
102 //
103 //    The user may also select the anisotropic flow parameters, Vn,
104 //    to either be fixed for each event, or to randomly vary from event
105 //    to event according to a Gaussian distribution where the user must
106 //    specify the mean and std. dev.  For both cases the input file must
107 //    list the 'nflowterms' (e.g. 6) values of the mean and Std.dev. for
108 //    the Vn parameters for all particle ID types included in the run.
109 //
110 //    The available list of particles has been increased to permit a
111 //    number of meson and baryon resonances.  For those with broad widths
112 //    the code samples the mass distribution for the resonance and outputs
113 //    the resonance mass for each instance in a special kinematic file
114 //    (see file unit=9, filename = 'mult_gen.kin').  The resonance shapes
115 //    are approximately Breit-Wigner and are specific for each resonance 
116 //    case.  The additional particle/resonances include: rho(+,-,0),
117 //    omega(0), eta', phi, J/Psi, Delta(-,0,+,++) and K*(+,-,0).  Masses
118 //    are sampled for the rho, omega, phi, Deltas and D*s. 
119 //    Refer to SUBR: Particle_prop and Particle_mass for the explicit
120 //    parameters, resonance shape models, and sampling ranges.
121 //
122 //    The input is from a file, named 'mult_gen.in'.  The output is
123 //    loaded into a file named 'mult_gen.out' which includes run
124 //    header information, event header information and the EVENT: and
125 //    TRACK: formats as in the new STAR TEXT Format for event generator
126 //    input to GSTAR.  A log file, 'mult_gen.log' is also written which
127 //    may contain error messages.  Normally this file should be empty
128 //    after a successful run.  These filenames can easily be changed
129 //    to more suitable names by the script that runs the program or
130 //    by hand.
131 //
132 //
133 // II.   ALGORITHM
134 //
135 //
136 //
137 //    The method for generating random multiplicities and model parameter
138 //    values involves the following steps:
139 //       (1) The Poisson or Gaussian distributions are computed and
140 //           loaded into function f().
141 //       (2) The distribution f(x') is integrated from xmin to x
142 //           and saved from x = xmin to x = xmax.  The range and mesh
143 //           spaces are specified by the user.
144 //       (3) The integral of f is normalized to unity where 
145 //           integral[f(x')](at x = xmin) = 0.0
146 //           integral[f(x')](at x = xmax) = 1.0
147 //       (4) A random number generator is called which delivers values
148 //           between 0.0 and 1.0.  
149 //       (5) We consider the coordinate x (from xmin to xmax) to be
150 //           dependent on the integral[f].  Using the random number
151 //           for the selected value of integral[f] the value of x
152 //           is obtained by interpolation.
153 //
154 //    An interpolation subroutine from Rubin Landau, Oregon State Univ.,
155 //    is used to do this interpolation; it involves uneven mesh point 
156 //    spacing.
157 //
158 //    The method for generating the particle momenta uses the
159 //    standard random elimination method and involves the following
160 //    steps:
161 //
162 //    For model_type = 1,2,3,4 which are functions of pt,y (see following):
163 //       (1) The y range is computed using the pseudorapidity (eta)
164 //           range and includes ample cushioning around the sides
165 //           along the eta acceptance edges.
166 //       (2) The transverse momentum (pt) and rapidity (y) are
167 //           randomly chosen within the specified ranges.
168 //       (3) The pseudorapidity is computed for this (pt,y) value
169 //           (and the mass for each pid) and checked against the
170 //           pseudorapidity acceptance range.
171 //       (4) If the pseudorapidity is within range then the one-particle
172 //           model distribution is calculated at this point and its ratio
173 //           to the maximum value throughout (pt,eta) acceptance region
174 //           is calculated.
175 //       (5) Another random number is called and if less than the ratio
176 //           from step#4 the particle momentum is used; if not, then 
177 //           another trial value of (pt,y) is obtained.
178 //       (6) This continues until the required multiplicity for the
179 //           specific event and particle type has been satisfied.
180 //       (7) This process is repeated for the requested number of particle
181 //           types and events.
182 //
183 //    For model_type = 5,6 (see following) which are input bin-by-bin
184 //    in pt,eta:
185 //       (1) The transverse momentum (pt) and pseudorapidity (eta) are 
186 //           randomly chosen within the specified ranges.
187 //       (2) The one-particle model distribution is calculated at this
188 //           point and its ratio to the maximum value throughout the
189 //           (pt,eta) region is calculated.
190 //       (3) Another random number is called and if less than the ratio
191 //           from step(2) the particle momentum is used; if not then
192 //           another trial value of (pt,eta) is obtained.
193 //       (4) This continues until the required multiplicity for the 
194 //           specific event and particle type has been satisfied.
195 //       (5) This process is repeated for the requested number of particle
196 //           types and events. 
197 //
198 //    Problematic parameter values are tested, bad input values are checked
199 //    and in some cases may be changed so that the program will not crash.
200 //    In some cases the code execution is stopped.
201 //    Some distributions and/or unusual model parameter values may cause the
202 //    code to hang up due to the poor performance of the "elimination"
203 //    method for very strongly peaked distributions.  These are tested for
204 //    certain problematic values and if necessary these events are aborted.
205 //    A message, "*** Event No.    2903 ABORTED:" for example is printed
206 //    in the 'mult_gen.out' file.  Temperatures .le. 0.01 GeV and rapidity
207 //    width parameters .le. 0.01 will cause the event to abort.
208 //
209 //
210 //
211 // III.  DESCRIPTION OF THE INPUT:
212 //
213 //
214 //    The input is described below in the 'read' statements and also in
215 //    the sample input file.  Some additional comments are as follows:
216 //
217 //    (1) n_events - Selected number of events in run. Can be anything
218 //                   .ge. 1.
219 //    (2) n_pid_type - Number of particle ID types to include in the
220 //                     particle list. e.g. pi(+) and pi(-) are counted
221 //                     separately.  The limit is set by parameter npid
222 //                     in the accompanying include file 'Parameter_values.inc'
223 //                     and is presently set at 20.
224 //    (3) model_type - equals 1,2,3,4,5 or 6 so far.  See comments in
225 //                     Function dNdpty to see what is calculated.
226 //                     The models included are:
227 //                   = 1, Factorized mt exponential, Gaussian rapidity model
228 //                   = 2, Pratt non-expanding, spherical thermal source model
229 //                   = 3, Bertsch non-expanding spherical thermal source model
230 //                   = 4, Pratt spherically expanding, thermally equilibrated
231 //                        source model.
232 //                   = 5, Factorized pt and eta distributions input bin-by-bin.
233 //                   = 6, Fully 2D pt,eta distributions input bin-by-bin.
234 //                        NOTE: model_type = 1-4 are functions of (pt,y)
235 //                              model_type = 5,6 are functions of (pt,eta)
236 //    (4) reac_plane_cntrl - Can be either 1,2,3 or 4 where:
237 //                         = 1 to ignore reaction plane and anisotropic flow,
238 //                             all distributions will be azimuthally symm.
239 //                         = 2 to use a fixed reaction plane angle for all
240 //                             events in the run.
241 //                         = 3 to assume a randomly varying reaction plane
242 //                             angle for each event as determined by a
243 //                             Gaussian distribution.
244 //                         = 4 to assume a randomly varying reaction plane
245 //                             for each event in the run as determined by
246 //                             a uniform distribution from 0 to 360 deg.
247 //    (5) PSIr_mean, PSIr_stdev - Reaction plane angle mean and Gaussian
248 //                                std.dev. (both are in degrees) for cases
249 //                                with reac_plane_cntrl = 2 (use mean value)
250 //                                and 3.  Note: these are read in regardless
251 //                                of the value of reac_plane_cntrl.
252 //    (6) MultFac_mean, MultFac_stdev - Overall multiplicity scaling factor
253 //                                      for all PID types; mean and std.dev.;
254 //                                      for trigger fluctuations event-to-evt.
255 //    (7) pt_cut_min,pt_cut_max - Range of transverse momentum in GeV/c.
256 //    (8) eta_cut_min,eta_cut_max - Pseudorapidity range
257 //    (9) phi_cut_min,phi_cut_max - Azimuthal angular range in degrees.
258 //   (10) n_stdev_mult - Number of standard deviations about the mean value
259 //                       of multiplicity to include in the random event-to-
260 //                       event selection process.  The maximum number of
261 //                       steps that can be covered is determined by
262 //                       parameter n_mult_max_steps in the accompanying
263 //                       include file 'Parameter_values.inc' which is
264 //                       presently set at 1000, but the true upper limit for
265 //                       this is n_mult_max_steps - 1 = 999.
266 //   (11) n_stdev_temp - Same, except for the "Temperature" parameter.
267 //   (12) n_stdev_sigma- Same, except for the rapidity width parameter.
268 //   (13) n_stdev_expvel - Same, except for the expansion velocity parameter.
269 //   (14) n_stdev_PSIr   - Same, except for the reaction plane angle
270 //   (15) n_stdev_Vn     - Same, except for the anisotropy coefficients, Vn.
271 //   (16) n_stdev_MultFac - Same, except for the multiplicity scaling factor.
272 //   (17) n_integ_pts - Number of mesh points to use in the random model
273 //                      parameter selection process.  The upper limit is
274 //                      set by parameter nmax_integ in the accompanying
275 //                      include file 'Parameter_values.inc' which is presently
276 //                      set at 100, but the true upper limit for n_integ_pts
277 //                      is nmax_integ - 1 = 99. 
278 //   (18) n_scan_pts  - Number of mesh points to use to scan the (pt,y)
279 //                      dependence of the model distributions looking for
280 //                      the maximum value.  The 2-D grid has
281 //                      n_scan_pts * n_scan_pts points; no limit to size of
282 //                      n_scan_pts.
283 //   (19) irand       - Starting random number seed.
284 //
285 //**************************************************************************
286 //   FOR MODEL_TYPE = 1,2,3 or 4:
287 //   Input the following 7 lines for each particle type; repeat these
288 //   set of lines n_pid_type times:
289 //
290 //        (a) gpid - Geant Particle ID code number
291 //        (b) mult_mean,mult_variance_control - Mean multiplicity and
292 //                                              variance control where:
293 //            mult_variance_control = 0 for no variance in multiplicity 
294 //            mult_variance_control = 1 to allow Poisson distribution for
295 //                                      particle multiplicities for all events.
296 //            Note that a hard limit exists for the maximum possible
297 //            multiplicity for a given particle type per event.  This is
298 //            determined by parameter factorial_max in accompanying include
299 //            file 'common_facfac.inc' and is presently set at 10000.
300 //        (c) Temp_mean, Temp_stdev - Temperature parameter mean (in GeV)
301 //            and standard deviation (Gaussian distribution assumed).
302 //        (d) sigma_mean, sigma_stdev - Rapidity distribution width (sigma)
303 //            parameter mean and standard deviation (Gaussian distribution
304 //            assumed).
305 //        (e) expvel_mean, expvel_stdev - S. Pratt expansion velocity
306 //            (in units of c) mean and standard deviation (Gaussian 
307 //            distribution assumed).
308 //        (f) Vn_mean(k);  k=1,4  - Anisotropic flow parameters, mean values
309 //                                  for Fourier component n=1.
310 //        (g) Vn_stdev(k); k=1,4  - Anisotropic flow parameters, std.dev.
311 //                                  values for Fourier component n=1.
312 //
313 //            Repeat the last two lines of input for remaining Fourier
314 //            components n=2,3...6.  Include all 6 sets of parameters
315 //            even if these are not used by the model for Vn(pt,y) (set
316 //            unused parameter means and std.dev. to 0.0).  List 4 values
317 //            on every line, even though for n=even the 4th quantity is
318 //            not used.
319 //
320 //**************************************************************************
321 //   FOR MODEL_TYPE = 5 input the following set of lines for each particle
322 //                      type; repeat these n_pid_type times.
323 //
324 //        (a) gpid - Geant Particle ID code number
325 //        (b) mult_mean,mult_variance_control - Mean multiplicity and
326 //                                              variance control where:
327 //            mult_variance_control = 0 for no variance in multiplicity
328 //            mult_variance_control = 1 to allow Poisson distribution for
329 //                                      particle multiplicities for all events.
330 //        (c) pt_start, eta_start - minimum starting values for pt, eta 
331 //                                  input for the bin-by-bin distributions.
332 //        (d) n_pt_bins, n_eta_bins - # input pt and eta bins.
333 //        (e) delta_pt, pt_bin - pt bin size and function value, repeat for
334 //                               each pt bin.
335 //        (f) delta_eta, eta_bin - eta bin size and function value, repeat
336 //                                 for each eta bin.
337 //        (g) Vn_mean(k);  k=1,4  - Anisotropic flow parameters, mean values
338 //                                  for Fourier component n=1.
339 //        (h) Vn_stdev(k); k=1,4  - Anisotropic flow parameters, std.dev.
340 //                                  values for Fourier component n=1.
341 //
342 //            Repeat the last two lines of input for remaining Fourier
343 //            components n=2,3...6.  Include all 6 sets of parameters
344 //            even if these are not used by the model for Vn(pt,y) (set
345 //            unused parameter means and std.dev. to 0.0).  List 4 values
346 //            on every line, even though for n=even the 4th quantity is
347 //            not used.
348 //
349 //        NOTE: The pt, eta ranges must fully include the requested ranges
350 //              in input #4 and 5 above; else the code execution will stop.
351 //
352 //        Also, variable bin sizes are permitted for the input distributions.
353 //
354 //        Also, this input distribution is used for all events in the run;
355 //        no fluctuations in this "parent" distribution are allowed from 
356 //        event-to-event.
357 //
358 //**************************************************************************
359 //   FOR MODEL_TYPE = 6 input the following set of lines for each particle
360 //                      type; repeat these n_pid_type times.
361 //
362 //        (a) gpid - Geant Particle ID code number
363 //        (b) mult_mean,mult_variance_control - Mean multiplicity and
364 //                                              variance control where:
365 //            mult_variance_control = 0 for no variance in multiplicity
366 //            mult_variance_control = 1 to allow Poisson distribution for
367 //                                      particle multiplicities for all events.
368 //        (c) pt_start, eta_start - minimum starting values for pt, eta
369 //                                  input for the bin-by-bin distributions.
370 //        (d) n_pt_bins, n_eta_bins - # input pt and eta bins.
371 //        (e) delta_pt - pt bin size, repeat for each pt bin. 
372 //        (f) delta_eta - eta bin size, repeat for each eta bin.
373 //        (g) i,j,pt_eta_bin(i,j) - read pt (index = i) and eta (index = j)
374 //                                  bin numbers and bin value for full 2D space
375 //        (h) Vn_mean(k);  k=1,4  - Anisotropic flow parameters, mean values
376 //                                  for Fourier component n=1.
377 //        (i) Vn_stdev(k); k=1,4  - Anisotropic flow parameters, std.dev.
378 //                                  values for Fourier component n=1.
379 //
380 //            Repeat the last two lines of input for remaining Fourier
381 //            components n=2,3...6.  Include all 6 sets of parameters
382 //            even if these are not used by the model for Vn(pt,y) (set
383 //            unused parameter means and std.dev. to 0.0).  List 4 values
384 //            on every line, even though for n=even the 4th quantity is
385 //            not used.
386 //
387 //        NOTE: The pt, eta ranges must fully include the requested ranges
388 //              in input #4 and 5 above; else the code execution will stop.
389 //
390 //        Also, variable bin sizes are permitted for the input distributions.
391 //
392 //        Also, this input distribution is used for all events in the run;
393 //        no fluctuations in this "parent" distribution are allowed from
394 //        event-to-event.
395 //
396 ///////////////////////////////////////////////////////////////////////////////
397
398
399
400
401 #include <fstream.h>
402 #include <iomanip.h>
403 #include "TMevSim.h"
404
405 #include "MevSimCommon.h"
406 #include "TParticle.h"
407 #include "TFile.h"
408
409 #ifndef WIN32
410 # define multgen multgen_
411 # define type_of_call
412 #else
413 # define multgen MULTGEN
414 # define type_of_call _stdcall
415 #endif
416
417
418 ClassImp(TMevSim)
419
420
421 extern "C" void type_of_call multgen();
422
423 //______________________________________________________________________________
424 TMevSim::TMevSim(Int_t nEvents, Int_t modelType, Int_t reacPlaneCntrl,
425            Float_t psiRMean, Float_t psiRStDev, Float_t multFacMean, Float_t multFacStDev,
426            Float_t ptCutMin, Float_t ptCutMax, Float_t etaCutMin, Float_t etaCutMax, 
427            Float_t phiCutMin, Float_t phiCutMax, Int_t irand) : TGenerator("MevSim", "MevSim")
428 {
429 // TMevSim constructor: initializes all the event-wide variables of MevSim with
430 // user supplied values, or with the default ones (declared in the header file).
431 // It also allocates space for the array which will store parameters specific to
432 // each particle species.
433 // Caution: Setting nEvents > 1 will have no effect, since only the last generated 
434 // event will be stored in POUT COMMON, and therefore only one event can be 
435 // accessible at a time. 
436  
437    fNEvents = nEvents;
438    fModelType = modelType;
439    fReacPlaneCntrl = reacPlaneCntrl;
440    fPsiRMean = psiRMean;
441    fPsiRStDev = psiRStDev;
442    fMultFacMean = multFacMean;
443    fMultFacStDev = multFacStDev;
444    fPtCutMin = ptCutMin;
445    fPtCutMax = ptCutMax;
446    fEtaCutMin = etaCutMin;
447    fEtaCutMax = etaCutMax;
448    fPhiCutMin = phiCutMin;
449    fPhiCutMax = phiCutMax;
450    fNStDevMult = fNStDevTemp = fNStDevSigma = fNStDevExpVel = fNStdDevPSIr = fNStDevVn = fNStDevMultFac = 3.0;
451    fNIntegPts = 100;
452    fNScanPts = 100;
453    firand = irand;
454    fParticleTypeParameters = new TClonesArray("TMevSimPartTypeParams",10);
455    DefineParticles();
456 }
457 //______________________________________________________________________________
458 TMevSim::~TMevSim()
459 {
460 // TMevSim destructor: destroys the object and all the particle information stored
461 // in the list.
462    
463   if (fParticleTypeParameters) {
464     fParticleTypeParameters->Clear();
465     delete fParticleTypeParameters;
466     fParticleTypeParameters = 0;
467   }
468 }
469 //______________________________________________________________________________
470 TMevSim::TMevSim(TMevSim& mevsim) {
471 // The copy constructor
472    
473    *this = mevsim;
474 }
475 //______________________________________________________________________________
476 TMevSim& TMevSim::operator=(TMevSim& mevsim) {
477 // An assignment operator: initializes all the event-wide variables of MevSim with
478 // the ones from a copied object. It also copies the parameters specific to
479 // each particle species.
480    
481    fNEvents = mevsim.GetNEvents();
482    fModelType = mevsim.GetModelType();
483    fReacPlaneCntrl = mevsim.GetReacPlaneCntrl();
484    fPsiRMean = mevsim.GetPsiRMean();
485    fPsiRStDev = mevsim.GetPsiRStDev();
486    fMultFacMean = mevsim.GetMultFacMean();
487    fMultFacStDev = mevsim.GetMultFacStDev();
488    fPtCutMin = mevsim.GetPtCutMin();
489    fPtCutMax = mevsim.GetPtCutMax();
490    fEtaCutMin = mevsim.GetEtaCutMin();
491    fEtaCutMax = mevsim.GetEtaCutMax();
492    fPhiCutMin = mevsim.GetPhiCutMin();
493    fPhiCutMax = mevsim.GetPhiCutMax();
494    fNStDevMult = mevsim.GetNStDevMult();
495    fNStDevTemp = mevsim.GetNStDevTemp();
496    fNStDevSigma =GetNStDevSigma();
497    fNStDevExpVel = mevsim.GetNStDevExpVel();
498    fNStdDevPSIr = mevsim.GetNStDevPSIr(); 
499    fNStDevVn =  mevsim.GetNStDevVn();
500    fNStDevMultFac =  mevsim.GetNStDevMultFac();
501    fNIntegPts = mevsim.GetNintegPts();
502    fNScanPts = mevsim.GetNScanPts();
503    firand = mevsim.firand;
504    fParticleTypeParameters = new TClonesArray("TMevSimPartTypeParams",mevsim.GetNPidTypes());
505    for (int i=0; i< mevsim.GetNPidTypes(); i++) 
506      {  
507         TMevSimPartTypeParams *temp = 0;
508         mevsim.GetPartTypeParamsByIndex(i,temp);
509         fParticleTypeParameters->AddLast(temp);
510      }
511    DefineParticles();
512    return (*this);
513 }
514 //______________________________________________________________________________
515 void        TMevSim::Initialize() {
516 // TMevSim initialization: creates an input file for the FORTRAN
517 // program MevSim. Converts all the event-wide information and particle
518 // specific information to the format readable by MevSim and writes it
519 // to disk in current directory.
520 // Caution: At least one TMevSimPartTypeParams object must be created and 
521 // added to the collection before event generation can start.
522    
523   TMevSimPartTypeParams * params = 0;
524   
525
526    ofstream *file = new ofstream("mult_gen.in",ios::out | ios::trunc);
527    // Write out the parameters to the pramameter file
528    *file << "  " << fNEvents << "  ! Number of Events \n";
529    *file << "  " << GetNPidTypes() << " \n";
530    *file << "  " << fModelType << " \n";
531    *file << "  " << fReacPlaneCntrl << " \n";
532    file->setf(ios::showpoint);
533    *file << "  " << fPsiRMean << "   " << fPsiRStDev << " \n";
534    *file << "  " << fMultFacMean << "   " << fMultFacStDev << " \n";
535    *file << "  " << fPtCutMin << "   " << fPtCutMax << " \n";
536    *file << "  " << fEtaCutMin << "   " << fEtaCutMax << " \n";
537    *file << "  " << fPhiCutMin << "   " << fPhiCutMax << " \n";
538    *file << "  " << fNStDevMult << " \n";
539    *file << "  " << fNStDevTemp << " \n";
540    *file << "  " << fNStDevSigma << " \n";
541    *file << "  " << fNStDevExpVel << " \n";
542    *file << "  " << fNStdDevPSIr << " \n";
543    *file << "  " << fNStDevVn << " \n";
544    *file << "  " << fNStDevMultFac << " \n";
545    *file << "  " << fNIntegPts << " \n";
546    *file << "  " << fNScanPts << " \n";
547    *file << "  " << firand << " \n";
548    // Write out particle specific information
549    for (Int_t i=0; i< (fParticleTypeParameters->GetLast() + 1); i++) {
550
551      params = (TMevSimPartTypeParams *) ((*fParticleTypeParameters)[i]);
552      
553      *file << "  " << params->GetGPid() << "       ! Particle GEANT Pid \n";
554      *file << "  " << params->GetMultMean() << "  " << params->GetMultVarianceControl() << "  \n";
555      *file << "  " << params->GetTempMean() << "  " << params->GetTempStDev() << "  \n";
556      *file << "  " << params->GetSigmaMean() << "  " << params->GetSigmaStDev() << "  \n";
557      *file << "  " << params->GetExpVelMean() << "  " << params->GetExpVelStDev() << "  \n";
558
559      for (Int_t cnt1 = 0; cnt1 < NFLOWTERMS; cnt1++) {
560        *file << "  ";
561        Int_t cnt2;
562        for (cnt2 = 0; cnt2 < 4; cnt2++) *file << params->GetVnMeanComponent(cnt1, cnt2) << "  ";
563        *file << "  \n  ";
564        for (cnt2 = 0; cnt2 < 4; cnt2++) *file << params->GetVnStDevComponent(cnt1, cnt2) << "  ";
565        *file << "  \n";
566      }
567    }
568    file->close();
569
570 }
571 //______________________________________________________________________________
572 void        TMevSim::GenerateEvent() {
573 // Generates one MevSim event. TMevSim::Initialize() must be called prior
574 // to calling this function.
575    
576    cout << "Calling FORTRAN multgen()" << endl;
577    multgen();
578 }
579
580 //______________________________________________________________________________
581 Int_t TMevSim::ImportParticles(TClonesArray *particles, Option_t *option)
582 {
583 // Read in particles created by MevSim into the TClonesArray(). The Initialize()
584 // and GenrateEvent() functions must be called prior to calling this funtion.
585 // The particles are read from the COMMON POUT. Right now the only provided 
586 // information is Geant PID, 3 momentum components and the energy of the particle.
587    
588    if (particles == 0) return 0;
589    TClonesArray &Particles = *particles;
590    Particles.Clear();
591    
592    Int_t totpart = 0;
593    for (Int_t nrpidtype=0; nrpidtype < (fParticleTypeParameters->GetLast() + 1); nrpidtype++) {
594       Int_t nrpart = 0;
595       Int_t pidcode = ((TMevSimPartTypeParams *) (*fParticleTypeParameters)[nrpidtype])->GetGPid();
596       while ((TRACK.pout[(4*nrpart+3)*NPID+nrpidtype] > 0.0) || (TRACK.pout[(4*nrpart)*NPID+nrpidtype] != 0.0)) {
597          int poffset = 4*nrpart*NPID+nrpidtype;
598          Float_t px = TRACK.pout[poffset];
599          poffset += NPID;
600          Float_t py = TRACK.pout[poffset];
601          poffset += NPID;
602          Float_t pz = TRACK.pout[poffset];
603          poffset += NPID;
604          Float_t mass = TRACK.pout[poffset];
605          new(Particles[totpart+nrpart]) TParticle(
606                                           PDGFromId(pidcode),  // Get the PDG ID from GEANT ID
607                                           0,
608                                           0,
609                                           0,
610                                           0,
611                                           0,
612                                           px,
613                                           py,
614                                           pz,
615                                           sqrt(mass*mass+px*px+py*py+pz*pz),                               
616                                           0,
617                                           0,
618                                           0,
619                                           0);
620          nrpart++;
621       }
622       totpart += nrpart;
623    }
624    return totpart;
625 }
626 //______________________________________________________________________________
627 void        TMevSim::SetNEvents(Int_t nEvents ) { 
628 // Sets the number of events to be generated by MevSim.
629 // Caution: Setting nEvents > 1 will have no effect, since only the last generated 
630 // event will be stored in POUT COMMON, and therefore only one event can be 
631 // accessible at a time. 
632
633    fNEvents = nEvents;
634 }
635 //______________________________________________________________________________
636 Int_t       TMevSim::GetNEvents() const {
637    return fNEvents;
638 }
639 //______________________________________________________________________________
640 Int_t       TMevSim::GetNPidTypes() const {
641    return fParticleTypeParameters->GetLast()+1;
642 }
643 //______________________________________________________________________________
644 void        TMevSim::SetModelType(Int_t modelType) {
645    fModelType  = modelType;
646 }
647 //______________________________________________________________________________
648 Int_t       TMevSim::GetModelType() const {
649    return fModelType;
650 }
651 //______________________________________________________________________________
652 void        TMevSim::SetReacPlaneCntrl(Int_t reacPlaneCntrl) {
653    fReacPlaneCntrl = reacPlaneCntrl; 
654 }
655 //______________________________________________________________________________
656 Int_t       TMevSim::GetReacPlaneCntrl() const {
657    return fReacPlaneCntrl;
658 }
659 //______________________________________________________________________________
660 void        TMevSim::SetPsiRParams(Float_t psiRMean,  Float_t psiRStDev) {
661    fPsiRMean = psiRMean;
662    fPsiRStDev = psiRStDev;
663 }
664 //______________________________________________________________________________
665 Float_t     TMevSim::GetPsiRMean() const { 
666    return fPsiRMean;
667 }
668 //______________________________________________________________________________
669 Float_t     TMevSim::GetPsiRStDev() const { 
670    return fPsiRStDev;
671 }
672 //______________________________________________________________________________
673 void        TMevSim::SetMultFacParams(Float_t multFacMean,  Float_t multFacStDev) {
674    fMultFacMean = multFacMean;
675    fMultFacStDev = multFacStDev;
676 }
677 //______________________________________________________________________________
678 Float_t     TMevSim::GetMultFacMean() const { 
679    return fMultFacMean;     
680 }
681 //______________________________________________________________________________
682 Float_t     TMevSim::GetMultFacStDev() const { 
683    return fMultFacStDev;
684 }
685 //______________________________________________________________________________
686 void        TMevSim::SetPtCutRange(Float_t ptCutMin,  Float_t ptCutMax) { 
687    fPtCutMin = ptCutMin;
688    fPtCutMax = ptCutMax;
689 }
690 //______________________________________________________________________________
691 Float_t     TMevSim::GetPtCutMin() const { 
692    return fPtCutMin;
693 }
694 //______________________________________________________________________________
695 Float_t     TMevSim::GetPtCutMax() const { 
696    return fPtCutMax;
697 }
698 //______________________________________________________________________________
699 void        TMevSim::SetEtaCutRange(Float_t etaCutMin,  Float_t etaCutMax) {     fEtaCutMin = etaCutMin;
700    fEtaCutMax = etaCutMax;
701 }
702
703 //______________________________________________________________________________
704    Float_t     TMevSim::GetEtaCutMin() const { 
705      return fEtaCutMin;
706 }
707 //______________________________________________________________________________
708    Float_t     TMevSim::GetEtaCutMax() const {
709      return fEtaCutMax;
710 }
711 //______________________________________________________________________________
712 void        TMevSim::SetPhiCutRange(Float_t phiCutMin,  Float_t phiCutMax) {
713    fPhiCutMin = phiCutMin;
714    fPhiCutMax = phiCutMax;
715 }
716 //______________________________________________________________________________
717 Float_t     TMevSim::GetPhiCutMin() const { 
718    return fPhiCutMin;
719 }
720 //______________________________________________________________________________
721 Float_t     TMevSim::GetPhiCutMax() const { 
722    return fPhiCutMax;
723 }
724 //______________________________________________________________________________
725 void        TMevSim::SetNStDevMult(Float_t nStDevMult) { 
726    fNStDevMult = nStDevMult;
727 }
728 //______________________________________________________________________________
729 Float_t       TMevSim::GetNStDevMult() const { 
730    return  fNStDevMult;
731 }
732 //______________________________________________________________________________
733 void        TMevSim::SetNStDevTemp(Float_t nStDevTemp) { 
734    fNStDevTemp = nStDevTemp;
735 }
736 //______________________________________________________________________________
737 Float_t       TMevSim::GetNStDevTemp() const { 
738    return fNStDevTemp;
739 }
740 //______________________________________________________________________________
741 void        TMevSim::SetNStDevSigma(Float_t nStDevSigma) { 
742    fNStDevSigma = nStDevSigma;
743 }
744 //______________________________________________________________________________
745 Float_t       TMevSim::GetNStDevSigma() const { 
746    return fNStDevSigma;  
747 }
748 //______________________________________________________________________________
749 void        TMevSim::SetNStDevExpVel(Float_t nStDevExpVel) { 
750    fNStDevExpVel = nStDevExpVel;
751 }
752 //______________________________________________________________________________
753 Float_t       TMevSim::GetNStDevExpVel() const { 
754    return fNStDevExpVel;
755 }
756 //______________________________________________________________________________
757 void        TMevSim::SetNStDevPSIr(Float_t nStDevPSIr) { 
758    fNStdDevPSIr = nStDevPSIr;
759 }
760 //______________________________________________________________________________
761 Float_t       TMevSim::GetNStDevPSIr() const { 
762    return fNStdDevPSIr;
763 }
764 //______________________________________________________________________________
765 void        TMevSim::SetNStDevVn(Float_t nStDevVn) { 
766    fNStDevVn = nStDevVn;
767 }
768 //______________________________________________________________________________
769 Float_t       TMevSim::GetNStDevVn() const { 
770    return fNStDevVn; 
771 }
772 //______________________________________________________________________________
773 void        TMevSim::SetNStDevMultFac(Float_t nStDevMultFac) { 
774    fNStDevMultFac = nStDevMultFac;
775 }
776 //______________________________________________________________________________
777 Float_t       TMevSim::GetNStDevMultFac() const { 
778    return fNStDevMultFac;  
779 }
780 //______________________________________________________________________________
781 void        TMevSim::SetNIntegPts(Int_t nIntegPts) { 
782    fNIntegPts = nIntegPts;
783 }
784 //______________________________________________________________________________
785 Int_t       TMevSim::GetNintegPts() const { 
786    return fNIntegPts;
787 }
788 //______________________________________________________________________________
789 void        TMevSim::SetNScanPts(Int_t nScanPts) { 
790    fNScanPts = nScanPts;
791 }
792 //______________________________________________________________________________
793 Int_t       TMevSim::GetNScanPts() const { 
794    return fNScanPts;  
795 }
796 //______________________________________________________________________________
797 void        TMevSim::AddPartTypeParams(TMevSimPartTypeParams *params) {
798 // Add the particle specied parameters and the end of the list.
799    
800   //cout << params << " " <<  fParticleTypeParameters <<  endl;
801
802   //fParticleTypeParameters->Dump(); 
803   params->Dump();
804   
805   Int_t last = fParticleTypeParameters->GetLast();
806   new ((*fParticleTypeParameters)[last+1]) TMevSimPartTypeParams(*params);
807 }
808 //______________________________________________________________________________
809 void        TMevSim::SetPartTypeParams(Int_t index, TMevSimPartTypeParams *params) 
810 {
811 // Create the new copy particle species parameters provided by params, and store
812 // them in the position designated by index. 
813    
814    *((TMevSimPartTypeParams *) ((*fParticleTypeParameters)[index])) = *params;
815 }
816 //______________________________________________________________________________
817 void        TMevSim::GetPartTypeParamsByIndex(Int_t index, TMevSimPartTypeParams *params) 
818 {
819 // Return the particle parameters stored in the list at the postion index.   
820 // Returns NULL if index is out of bounds.
821    
822    if ((index < fParticleTypeParameters->GetLast()) && (index >= 0))
823         params = (TMevSimPartTypeParams *) (*fParticleTypeParameters)[index];
824    else
825      params = NULL;
826 }
827 //______________________________________________________________________________
828 void        TMevSim::GetPartTypeParamsByGPid(Int_t gpid, TMevSimPartTypeParams *params) 
829 {
830 // Return the particle parameters for the particle with Geant PID gpid.
831 // Returns NULL if the parameters for such particle do not exist in the list.   
832    
833    Int_t i = -1;
834    
835    while (++i <= fParticleTypeParameters->GetLast())
836      {
837         if (((TMevSimPartTypeParams *) (*fParticleTypeParameters)[i])->GetGPid() == gpid)
838           {
839              params = (TMevSimPartTypeParams *) (*fParticleTypeParameters)[i];
840              return;
841           }
842      }
843    params = NULL;
844    return;
845 }
846 //_____________________________________________________________________________
847 Int_t TMevSim::PDGFromId(Int_t id) const 
848 {
849   //
850   // Return PDG code and pseudo ENDF code from Geant3 code
851   //
852   if(id>0 && id<fNPDGCodes) return fPDGCode[id];
853   else return -1;
854 }
855 //_____________________________________________________________________________
856 void TMevSim::DefineParticles() 
857 {
858   //
859   // Load standard numbers for GEANT particles and PDG conversion
860   fPDGCode[fNPDGCodes++]=-99;   //  0 = unused location
861   fPDGCode[fNPDGCodes++]=22;    //  1 = photon
862   fPDGCode[fNPDGCodes++]=-11;   //  2 = positron
863   fPDGCode[fNPDGCodes++]=11;    //  3 = electron
864   fPDGCode[fNPDGCodes++]=12;    //  4 = neutrino e
865   fPDGCode[fNPDGCodes++]=-13;   //  5 = muon +
866   fPDGCode[fNPDGCodes++]=13;    //  6 = muon -
867   fPDGCode[fNPDGCodes++]=111;   //  7 = pi0
868   fPDGCode[fNPDGCodes++]=211;   //  8 = pi+
869   fPDGCode[fNPDGCodes++]=-211;  //  9 = pi-
870   fPDGCode[fNPDGCodes++]=130;   // 10 = Kaon Long
871   fPDGCode[fNPDGCodes++]=321;   // 11 = Kaon +
872   fPDGCode[fNPDGCodes++]=-321;  // 12 = Kaon -
873   fPDGCode[fNPDGCodes++]=2112;  // 13 = Neutron
874   fPDGCode[fNPDGCodes++]=2212;  // 14 = Proton
875   fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
876   fPDGCode[fNPDGCodes++]=310;   // 16 = Kaon Short
877   fPDGCode[fNPDGCodes++]=221;   // 17 = Eta
878   fPDGCode[fNPDGCodes++]=3122;  // 18 = Lambda
879   fPDGCode[fNPDGCodes++]=3222;  // 19 = Sigma +
880   fPDGCode[fNPDGCodes++]=3212;  // 20 = Sigma 0
881   fPDGCode[fNPDGCodes++]=3112;  // 21 = Sigma -
882   fPDGCode[fNPDGCodes++]=3322;  // 22 = Xi0
883   fPDGCode[fNPDGCodes++]=3312;  // 23 = Xi-
884   fPDGCode[fNPDGCodes++]=3334;  // 24 = Omega-
885   fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
886   fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
887   fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
888   fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
889   fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
890   fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
891   fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
892   fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
893 }
894