]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TEPEMGEN/TEpEmGen.cxx
Completely reengineered version of CMake build system (Johny)
[u/mrichter/AliRoot.git] / TEPEMGEN / TEpEmGen.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2002, 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  * Copyright(c) 1997, 1998, 2002, Adrian Alscher and Kai Hencken          *
17  * See $ALICE_ROOT/EpEmGen/diffcross.f for full Copyright notice          *
18  *                                                                        *
19  *                                                                        *
20  * Copyright(c) 2002 Kai Hencken, Yuri Kharlov, Serguei Sadovsky          *
21  * See $ALICE_ROOT/EpEmGen/epemgen.f for full Copyright notice            *
22  *                                                                        *
23  **************************************************************************/
24 //------------------------------------------------------------------------
25 // TEpEmGen is an interface class to fortran event generator of
26 // single e+e- pair production in ultraperipheral PbPb collisions
27 // at 5.5 GeV/c
28 //%
29 // Yuri.Kharlov@cern.ch
30 // 9 October 2002
31 //------------------------------------------------------------------------
32
33 #include "TRandom.h"
34
35 #include "TEpEmGen.h"
36 #include "TClonesArray.h"
37 #include "TParticle.h"
38 #include "TEcommon.h"
39
40 #ifndef WIN32
41 # define ee_init  ee_init_
42 # define ee_event ee_event_
43 # define eernd    eernd_
44 #else
45 # define ee_init  EE_INIT
46 # define ee_event EE_EVENT
47 # define eernd    EERND
48 #endif
49
50 extern "C" {
51   void ee_init  (Double_t &ymin, Double_t &ymax, Double_t &xmin, Double_t &xmax);
52   void ee_event (Double_t &ymin, Double_t &ymax, Double_t &xmin, Double_t &xmax,
53                  Double_t &yE,   Double_t &yP,   Double_t &xE,   Double_t &xP, 
54                  Double_t &phi,  Double_t &w);
55   Double_t eernd(Int_t*) {
56     Double_t r;
57     do r=gRandom->Rndm(); while(0 >= r || r >= 1);
58     return r;
59   }
60 }
61
62 ClassImp(TEpEmGen)
63
64 //------------------------------------------------------------------------------
65 TEpEmGen::TEpEmGen() : TGenerator("TEpEmGen","TEpEmGen")
66 {
67 // TEpEmGen constructor: creates a TClonesArray in which it will store all
68 // particles. Note that there may be only one functional TEpEmGen object
69 // at a time, so it's not use to create more than one instance of it.
70
71   delete fParticles; // was allocated as TObjArray in TGenerator
72   fParticles = new TClonesArray("TParticle",2);
73
74 }
75
76 //------------------------------------------------------------------------------
77 TEpEmGen::~TEpEmGen()
78 {
79   // Destroys the object, deletes and disposes all TParticles currently on list.
80   if (fParticles) {
81     fParticles->Delete();
82     delete fParticles;
83     fParticles = 0;
84   }
85 }
86
87 //______________________________________________________________________________
88 void TEpEmGen::GenerateEvent(Double_t ymin, Double_t ymax, Double_t ptmin, Double_t ptmax,
89                              Double_t &yElectron, Double_t &yPositron,
90                              Double_t &xElectron, Double_t &xPositron,
91                              Double_t &phi12,     Double_t &weight)
92 {
93   //produce one event
94   ee_event(ymin,ymax,ptmin,ptmax,
95            yElectron,yPositron,xElectron,xPositron,
96            phi12,weight);
97 }
98
99 //______________________________________________________________________________
100 void TEpEmGen::Initialize(Double_t ymin, Double_t ymax, Double_t ptmin, Double_t ptmax)
101 {
102   // Initialize EpEmGen
103   Double_t ptminMeV = ptmin*1000;
104   Double_t ptmaxMeV = ptmax*1000;
105   ee_init(ymin,ymax,ptminMeV,ptmaxMeV);
106 }
107
108 //______________________________________________________________________________
109 Double_t TEpEmGen::GetXsection()
110 {
111   // Return cross section accumulated so far
112   return EEVENT.Xsecttot;
113 }
114
115 //______________________________________________________________________________
116 Double_t TEpEmGen::GetDsection()
117 {
118   // Return cross section error accumulated so far
119   return EEVENT.Dsecttot;
120 }