FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairBoxGenerator.h
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 // -------------------------------------------------------------------------
9 // ----- FairBoxGenerator header file -----
10 // ----- Created 09/09/04 by Yu.Kharlov -----
11 // -------------------------------------------------------------------------
12 
21 /* $Id: FairBoxGenerator.h,v 1.3 2006/07/14 11:23:57 kharlov Exp $ */
22 
23 /* History of cvs commits:
24  *
25  * $Log: FairBoxGenerator.h,v $
26  * Revision 1.3 2006/07/14 11:23:57 kharlov
27  * Add protection for simultaneously set ranges; split vertex and kinematics ranges
28  *
29  * Revision 1.2 2006/03/29 16:25:50 kharlov
30  * New functionality added
31  *
32  */
33 
34 #ifndef FAIR_BOXGENERATOR_H
35 #define FAIR_BOXGENERATOR_H
36 
37 #include "FairBaseMCGenerator.h"
38 #include "FairGenerator.h" // for FairGenerator
39 
40 #include <Rtypes.h> // for Double32_t, Bool_t, kTRUE, etc
41 
43 
45 {
46  public:
49 
54  FairBoxGenerator(Int_t pdgid, Int_t mult = 1);
55 
57  virtual ~FairBoxGenerator() {}
58 
61  void SetPRange(Double32_t pmin = 0, Double32_t pmax = 10)
62  {
63  fPMin = pmin;
64  fPMax = pmax;
65  fPRangeIsSet = kTRUE;
66  }
67 
68  void SetPtRange(Double32_t ptmin = 0, Double32_t ptmax = 10)
69  {
70  fPtMin = ptmin;
71  fPtMax = ptmax;
72  fPtRangeIsSet = kTRUE;
73  }
74 
75  void SetEkinRange(Double32_t kmin = 0, Double32_t kmax = 10)
76  {
77  fEkinMin = kmin;
78  fEkinMax = kmax;
79  fEkinRangeIsSet = kTRUE;
80  }
81 
82  void SetPhiRange(Double32_t phimin = 0, Double32_t phimax = 360)
83  {
84  fPhiMin = phimin;
85  fPhiMax = phimax;
86  }
87 
88  void SetEtaRange(Double32_t etamin = -5, Double32_t etamax = 7)
89  {
90  fEtaMin = etamin;
91  fEtaMax = etamax;
92  fEtaRangeIsSet = kTRUE;
93  }
94 
95  void SetYRange(Double32_t ymin = -5, Double32_t ymax = 7)
96  {
97  fYMin = ymin;
98  fYMax = ymax;
99  fYRangeIsSet = kTRUE;
100  }
101 
102  void SetThetaRange(Double32_t thetamin = 0, Double32_t thetamax = 90)
103  {
104  fThetaMin = thetamin;
105  fThetaMax = thetamax;
106  fThetaRangeIsSet = kTRUE;
107  }
108 
109  void SetCosTheta() { fCosThetaIsSet = kTRUE; }
110 
111  void SetXYZ(Double32_t x = 0, Double32_t y = 0, Double32_t z = 0);
112 
113  void SetBoxXYZ(Double32_t x1 = 0, Double32_t y1 = 0, Double32_t x2 = 0, Double32_t y2 = 0, Double32_t z = 0);
114 
119  void SetDebug(Bool_t /*debug*/) {}
121  Bool_t Init();
122 
126  virtual Bool_t ReadEvent(FairPrimaryGenerator* primGen);
127 
129  virtual FairGenerator* CloneGenerator() const;
130 
131  protected:
133  FairBoxGenerator(const FairBoxGenerator&) = default;
134  FairBoxGenerator& operator=(const FairBoxGenerator&) = default;
135 
136  private:
137  Double32_t fPtMin, fPtMax; // Transverse momentum range [GeV]
138  Double32_t fPhiMin, fPhiMax; // Azimuth angle range [degree]
139  Double32_t fEtaMin, fEtaMax; // Pseudorapidity range in lab system
140  Double32_t fYMin, fYMax; // Rapidity range in lab system
141  Double32_t fPMin, fPMax; // Momentum range in lab system
142  Double32_t fThetaMin, fThetaMax; // Polar angle range in lab system [degree]
143  Double32_t fEkinMin, fEkinMax; // Kinetic Energy range in lab system [GeV]
144 
145  Bool_t fEtaRangeIsSet; // True if eta range is set
146  Bool_t fYRangeIsSet; // True if rapidity range is set
147  Bool_t fThetaRangeIsSet; // True if theta range is set
148  Bool_t fCosThetaIsSet; // True if uniform distribution in
149  // cos(theta) is set (default -> not set)
150  Bool_t fPtRangeIsSet; // True if transverse momentum range is set
151  Bool_t fPRangeIsSet; // True if abs.momentum range is set
152  Bool_t fEkinRangeIsSet; // True if kinetic energy range is set
153 
154  ClassDef(FairBoxGenerator, 5);
155 };
156 
157 #endif
void SetDebug(Bool_t)
void SetPtRange(Double32_t ptmin=0, Double32_t ptmax=10)
void SetEkinRange(Double32_t kmin=0, Double32_t kmax=10)
void SetYRange(Double32_t ymin=-5, Double32_t ymax=7)
void SetPRange(Double32_t pmin=0, Double32_t pmax=10)
void SetEtaRange(Double32_t etamin=-5, Double32_t etamax=7)
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
void SetXYZ(Double32_t x=0, Double32_t y=0, Double32_t z=0)
virtual FairGenerator * CloneGenerator() const
virtual ~FairBoxGenerator()
void SetBoxXYZ(Double32_t x1=0, Double32_t y1=0, Double32_t x2=0, Double32_t y2=0, Double32_t z=0)
void SetThetaRange(Double32_t thetamin=0, Double32_t thetamax=90)
FairBoxGenerator & operator=(const FairBoxGenerator &)=default
void SetPhiRange(Double32_t phimin=0, Double32_t phimax=360)