]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TEPEMGEN/TEpEmGen.cxx
Including header files without dirs
[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 "TEpEmGen.h"
34 #include "TClonesArray.h"
35 #include "TParticle.h"
36 #include "TEcommon.h"
37
38 #ifndef WIN32
39 # define ee_init  ee_init_
40 # define ee_event ee_event_
41 #else
42 # define ee_init  EE_INIT
43 # define ee_event EE_EVENT
44 #endif
45
46 extern "C" {
47   void ee_init  (Double_t &ymin, Double_t &ymax, Double_t &xmin, Double_t &xmax);
48   void ee_event (Double_t &ymin, Double_t &ymax, Double_t &xmin, Double_t &xmax,
49                  Double_t &yE,   Double_t &yP,   Double_t &xE,   Double_t &xP, 
50                  Double_t &phi,  Double_t &w);
51 }
52
53 ClassImp(TEpEmGen)
54
55 //------------------------------------------------------------------------------
56 TEpEmGen::TEpEmGen() : TGenerator("TEpEmGen","TEpEmGen")
57 {
58 // TEpEmGen constructor: creates a TClonesArray in which it will store all
59 // particles. Note that there may be only one functional TEpEmGen object
60 // at a time, so it's not use to create more than one instance of it.
61
62   delete fParticles; // was allocated as TObjArray in TGenerator
63   fParticles = new TClonesArray("TParticle",2);
64
65 }
66
67 //------------------------------------------------------------------------------
68 TEpEmGen::~TEpEmGen()
69 {
70   // Destroys the object, deletes and disposes all TParticles currently on list.
71   if (fParticles) {
72     fParticles->Delete();
73     delete fParticles;
74     fParticles = 0;
75   }
76 }
77
78 //______________________________________________________________________________
79 void TEpEmGen::GenerateEvent(Double_t ymin, Double_t ymax, Double_t ptmin, Double_t ptmax,
80                              Double_t &yElectron, Double_t &yPositron,
81                              Double_t &xElectron, Double_t &xPositron,
82                              Double_t &phi12,     Double_t &weight)
83 {
84   //produce one event
85   ee_event(ymin,ymax,ptmin,ptmax,
86            yElectron,yPositron,xElectron,xPositron,
87            phi12,weight);
88 }
89
90 //______________________________________________________________________________
91 void TEpEmGen::Initialize(Double_t ymin, Double_t ymax, Double_t ptmin, Double_t ptmax)
92 {
93   // Initialize EpEmGen
94   Double_t ptminMeV = ptmin*1000;
95   Double_t ptmaxMeV = ptmax*1000;
96   ee_init(ymin,ymax,ptminMeV,ptmaxMeV);
97 }
98
99 //______________________________________________________________________________
100 Double_t TEpEmGen::GetXsection()
101 {
102   // Return cross section accumulated so far
103   return EEVENT.Xsecttot;
104 }
105
106 //______________________________________________________________________________
107 Double_t TEpEmGen::GetDsection()
108 {
109   // Return cross section error accumulated so far
110   return EEVENT.Dsecttot;
111 }