CMS 3D CMS Logo

CrossSection.h
Go to the documentation of this file.
1 //-*-C++-*-
2 //-*-CrossSection.h-*-
3 // Written by James Monk and Andrew Pilkington
5 
6 #ifndef CROSSSECTION_HH
7 #define CROSSSECTION_HH
8 
9 #include <cmath>
10 #include <complex>
11 #include <cstdlib>
12 #include <iostream>
13 #include <map>
14 #include <string>
15 #include <utility>
16 #include <vector>
17 
18 //#include "CLHEP/config/CLHEP.h"
19 #include "CLHEP/Vector/LorentzVector.h"
20 #include "CLHEP/Vector/ThreeVector.h"
21 
23 //#include "GeneratorInterface/ExhumeInterface/interface/PythiaRecord.h"
25 
27 
28 //#include "CLHEP/HepMC/include/PythiaWrapper6_2.h"
29 //#include "CLHEP/HepMC/ConvertHEPEVT.h"
30 //#include "CLHEP/HepMC/CBhepevt.h"
31 
32 namespace CLHEP {
33 class HepRandomEngine;
34 }
35 
37 namespace Exhume{
38 
39  typedef std::pair<double,double> fEntry;
40  typedef std::pair<const char*,char*> PCharPair;
41  typedef std::pair<const char*,const void*> PConstVoidPair;
42 
43  class CrossSection{
44 
45  public:
46 
48  virtual ~CrossSection();
49  virtual void MaximiseSubParameters()=0;
50  virtual void SetPartons()=0;
51  virtual void SetSubParameters()=0;
52  virtual double SubParameterRange()=0;
53  virtual double SubParameterWeight()=0;
54 
55  double AlphaS(const double&);
56 
57  inline void SetRandomEngine(CLHEP::HepRandomEngine* engine){randomEngine = engine;}
58 
59  inline double GetRg(const double &x_, const double &Qt){
60 
61  return(Rg_(x_, Qt) );
62  };
63 
64  //.....
65  inline double Differential(){
66 
67  if(x1>1.0 || x2 > 1.0){
68  return(0.0);
69  }
70 
71  //return( SubProcess() );
72 
73  //Factor of 2/sqrt(sHat) because we integrate dM not dln(M^2)
74  return( 2.0 * InvSqrtsHat * Lumi() * SubProcess() * exp(B*(t1 + t2)));
75  };
76  //.....
77 
78  inline double GetB(){
79  return(B);
80  };
81 
82  inline double GetRoot_s(){
83  return(root_s);
84  };
85 
86  inline std::map<double,double> Getfg2Map(){
87  return(fg2Map);
88  };
89 
90  inline CLHEP::HepLorentzVector GetProton1(){
91  return(Proton1);
92  };
93 
94  inline CLHEP::HepLorentzVector GetProton2(){
95  return(Proton2);
96  };
97 
98  inline std::vector<Particle> GetPartons(){
99  return(Partons);
100  };
101 
102  inline double Getx1(){
103  return(x1);
104  };
105 
106  inline double Getx2(){
107  return(x2);
108  };
109 
110  inline double Gett1(){
111  return(t1);
112  };
113 
114  inline double Gett2(){
115  return(t2);
116  };
117 
118  inline double GetsHat(){
119  return(sHat);
120  };
121 
122  inline double GetSqrtsHat(){
123  return(SqrtsHat);
124  };
125 
126  inline double GetEta(){
127  return(y);
128  };
129 
130  inline double GetPhi1(){
131  return(Phi1);
132  };
133 
134  inline double GetPhi2(){
135  return(Phi2);
136  };
137 
138  void Hadronise();
139  void SetKinematics(const double&, const double&,
140  const double&, const double&,
141  const double&, const double&);
142 
144  return(Name);
145  };
146 
147  private:
148 
149  void AlphaSInit();
150 
151  double Fg1Fg2(const double&, const double&, const double&);
152  double Fg1Fg2(const int&, const double&, const double&);
153 
154  inline double Fg_Qt2(const double &Qt2_,
155  const double &x_, const double &xp_){
156 
157  double grad = 5.0 * (Txg(1.1 * Qt2_, x_) - Txg(0.9 * Qt2_, x_) ) / Qt2_;
158 
159  //protect from Qt2 == 0
160 
161  if(grad!=grad){
162 
163  return(0.0);
164  }
165  return(grad);
166  };
167 
168  template <typename T_>
169  inline void insert(const std::string _name_, const T_ _x_){
170  PMap.insert
171  (std::pair<std::string,PConstVoidPair>
172  (_name_,PConstVoidPair(typeid(_x_).name(),_x_)));
173  }
174 
175  virtual double SubProcess() = 0;
176 
177  void LumiInit();
178  double Lumi();
179  double Lumi_();
180  void NoMem();
181  void ReadCard(const std::string&);
182  double Splitting(const double&);
183  double T(const double&);
184  double TFast(const double&);
185 
186  bool TBegin, TInterpolate;
187  std::map<double, std::map<double, double> > TMap2d;
188  std::map<double, std::map<double, double> >::iterator TMjuHigh, TMjuLow;
189 
190  double Txg(const double&, const double&);
191 
192  double Rg1Rg2(const double&);
193  double Rg_(const double&, double);
194 
195  std::map<double, std::map<double, double> > RgMap2d;
196  std::map<double, std::map<double, double> >::iterator RgHigh[2], RgLow[2];
197  bool RgInterpolate[2], RgBegin;
198 
199  double B,
200  Freeze,
201  LambdaQCD,
202  Rg,
203  Survive;
204 
205  double PDF;
206 
207  double ASFreeze, ASConst;
208  //.............................
209  //Used for Lumi function:
210 
211  double MinQt2, MidQt2, InvMidQt2, InvMidQt4, InvMaxQt2,
212  LumSimpsIncr;
213  unsigned int LumAccuracy, LumNStart;
214  int LumNSimps, LumNSimps_1;
215  double *LumSimpsFunc, *_Qt2, *_Qt, *_KtLow, *_KtHigh,
216  *_AlphaS, *_CfAs_2PIRg, *_NcAs_2PI;
217  double LumConst;
218  double Inv2PI, Nc_2PI, Cf_2PIRg;
219 
220  std::map<std::string,PConstVoidPair> PMap;
221  std::map<const char*,char*> TypeMap;
222  std::map<double, double> fg2Map;
223 
224  //...............................
225  //Used for T:
226  double TConst;
227  int Tn, Tn_1;
228  double *TFunc;
229 
230  //...............................
231  //Used for Splitting:
232  double Inv3;
233 
234  int Proton1Id, Proton2Id;
235 
236  protected:
237 
239 
240  std::complex<double> F0(const double&);
241  std::complex<double> f(const double&);
242  std::complex<double> Fsf(const double&);
243 
244  //PPhi is azimuthal angle between protons.
245  //InvSqrtsHat = 1/sHat
246  //Others are obvious.
247  double x1, x1p, x2, x2p, t1, t2, sHat, SqrtsHat,
248  sHat2, InvsHat, InvsHat2, InvSqrtsHat, y, PPhi, Phi1, Phi2,
249  Mju2, Mju, LnMju2,
250  Pt1, Pt2, Pt1DotPt2,
251  x1x2, ey;
252 
253  double AlphaEw,
254  gw,
255  HiggsVev,
256  LambdaW;
257 
258  double BottomMass,
259  CharmMass,
260  StrangeMass,
261  TopMass;
262 
263  double MuonMass,
264  TauMass;
265 
266  double HiggsMass,
267  WMass,
268  ZMass;
269 
271  double root_s,
272  s, Invs;
273 
274  CLHEP::HepLorentzVector CentralVector;
275  CLHEP::HepLorentzVector Proton1, Proton2, P1In, P2In;
276  /*
277  struct Particle {
278  CLHEP::HepLorentzVector p;
279  int id;
280  };
281  */
282 
283  std::vector<Particle> Partons;
284  double Invsx1x2, InvV1MinusV2;
285 
286  double Gev2fb;
287 
289 
290  CLHEP::HepRandomEngine* randomEngine;
291  };
292 }
293 
294 #endif
295 
CLHEP::HepRandomEngine * randomEngine
Definition: Dummies.cc:7
CLHEP::HepLorentzVector CentralVector
Definition: CrossSection.h:274
void insert(const std::string _name_, const T_ _x_)
Definition: CrossSection.h:169
std::map< std::string, PConstVoidPair > PMap
Definition: CrossSection.h:220
std::string lhapdfSetPath_
Definition: CrossSection.h:288
std::string GetName()
Definition: CrossSection.h:143
std::vector< Particle > Partons
Definition: CrossSection.h:283
CLHEP::HepLorentzVector Proton2
Definition: CrossSection.h:275
std::map< double, double > fg2Map
Definition: CrossSection.h:222
CLHEP::HepLorentzVector GetProton1()
Definition: CrossSection.h:90
std::vector< Particle > GetPartons()
Definition: CrossSection.h:98
std::map< double, std::map< double, double > > TMap2d
Definition: CrossSection.h:187
std::map< const char *, char * > TypeMap
Definition: CrossSection.h:221
std::map< double, std::map< double, double > >::iterator TMjuLow
Definition: CrossSection.h:188
std::pair< const char *, const void * > PConstVoidPair
Definition: CrossSection.h:41
std::pair< const char *, char * > PCharPair
Definition: CrossSection.h:40
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
double f[11][100]
static const std::string B
double Fg_Qt2(const double &Qt2_, const double &x_, const double &xp_)
Definition: CrossSection.h:154
CLHEP::HepLorentzVector GetProton2()
Definition: CrossSection.h:94
std::pair< double, double > fEntry
Definition: CrossSection.h:39
CLHEP::HepRandomEngine * randomEngine
Definition: CrossSection.h:290
const double MuonMass
void SetRandomEngine(CLHEP::HepRandomEngine *engine)
Definition: CrossSection.h:57
std::map< double, std::map< double, double > > RgMap2d
Definition: CrossSection.h:195
unsigned int LumNStart
Definition: CrossSection.h:213
std::map< double, double > Getfg2Map()
Definition: CrossSection.h:86
long double T
double GetRg(const double &x_, const double &Qt)
Definition: CrossSection.h:59