BOSS 7.0.1
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtConExc.hh
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2//
3// Environment:
4// This software is part of models developed at BES collaboration
5// based on the EvtGen framework. If you use all or part
6// of it, please give an appropriate acknowledgement.
7//
8// Copyright Information: See EvtGen/BesCopyright
9// Copyright (A) 2006 Ping Rong-Gang @IHEP
10//
11// Module: EvtConExc.hh
12//
13// Description: To define cross section for the continuum exclusive process
14// Experimental cross section taken from PRD73,012005, PRD76,092006, also
15// see a review: Rev. Mod. Phys. 83,1545
16 /*******************--- mode definition: also see EvtXsection.cc
17 0: ppbar
18 1: nnbar
19 2: Lambda0 anti-Lambda0
20 3: Sigma0 anti-Sigma0
21 4: Lambda0 anti-Sigma0
22 5: Sigma0 anti-Lambda0
23 6: pi+ pi-
24 7: pi+ pi- pi0
25 8: K+K- pi0
26 9: KsK+pi-
27 10: KsK-pi+
28 11: K+K-eta
29 12: 2(pi+pi-)
30 13: pi+pi-2pi0
31 14: K+K-pi+pi-
32 15: K+K-2pi0
33 16: 2(K+K-)
34 17: 2(pi+pi-)pi0
35 18: 2(pi+pi-)eta
36 19: K+K-pi+pi-pi0
37 20: K+K-pi+pi-eta
38 21: 3(pi+pi-)
39 22: 2(pi+pi-pi0)
40 23: phi eta
41 24: phi pi0
42 25: K+K*-
43 26: K-K*+
44 27: K_SK*0-bar
45 28: K*0(892)K+pi-
46 29: K*0(892)K-pi+
47 30: K*+K-pi0
48 31: K*-K+pi0
49 32: K_2*(1430)0 K+pi-
50 33: K_2*(1430)0 K-pi+
51 34: K+K-rho
52 35: phi pi+pi-
53 36: phi f0(980)
54 37: eta pi+pi-
55 38: omega pi+ pi-
56 39: omega f0(980)
57 40: eta' pi+ pi-
58 41: f_1(1285)pi+pi-
59 42: omega K+K-
60 43: omega pi+pi-pi0
61 *************************************/
62// Modification history:
63//
64// Ping R.-G. Nov., 2012 Module created
65//
66//------------------------------------------------------------------------
67//
68#ifndef EVTCONEXC_HH
69#define EVTCONEXC_HH
70
71//#include <string.h>
72#include "EvtGenBase/EvtId.hh"
75//------ including root
76#include "TH1.h"
77#include "TH2.h"
78#include "TH3.h"
79#include "TFile.h"
80#include "TROOT.h"
81#include "TTree.h"
82#include "TGraphErrors.h"
83#include "TDirectory.h"
84
85class EvtParticle;
86
88
89public:
90
92 } //constructor
93 //---
94 virtual ~EvtConExc();
95
96 void getName(std::string& name);
97
99
100 void initProbMax();
101
102 void init();
103 void init_Br_ee();
104
105 void decay(EvtParticle *p);
106 void init_mode(int mode);
107 double gamHXSection(EvtParticle* p, double El, double Eh, int nmc=100000);
108 double gamHXSection(double s, double El, double Eh, int nmc=100000);
109 double gamHXSection(double El, double Eh);
110 double gamHXSection_er(double El,double Eh);
111
112 void findMaxXS(EvtParticle *p);
113 double difgamXs(EvtParticle* p); //differential cross section for gamma + hadrons
114 double difgamXs(double mhds,double sintheta);
115 double Mhad_sampling(double *x,double *y);
116 double ISR_ang_integrate(double x,double theta);
117 double ISR_ang_sampling(double x);
119 void SetP4(EvtParticle *part, double mhdr,double xeng, double theta); //set the gamma and gamma* momentum according sampled results
120 void SetP4Rvalue(EvtParticle *part, double mhdr,double xeng, double theta); //set the gamma and gamma* momentum according sampled results
121 bool gam_sampling(EvtParticle *p);
122 bool xs_sampling(double xs);
123 bool xs_sampling(double xs,double xs1);
124 bool baryon_sampling(EvtVector4R pcm, EvtVector4R pi);//baryon angular distri. 1+cos^2\theta
125 bool meson_sampling(EvtVector4R pcm, EvtVector4R pi); //meson angular distri. sin^2\theta
126 bool VP_sampling(EvtVector4R pcm, EvtVector4R pi); //VP ang. dist. 1+cos^2theta
127 bool angularSampling(EvtParticle* part);
128 bool photonSampling(EvtParticle* part);
129 double baryonAng(double mx);
130 double Rad1(double s, double x);
131 double Rad2(double s, double x);
132 double Rad2difXs(EvtParticle *p);
133 double Rad2difXs(double s, double x);
134 double Rad1difXs(EvtParticle *p);
135 double Ros_xs(double mx, double bree,EvtId pid);
137 static double _cms; //energy of CMS of ee beam
138 static double XS_max;// maxium of cross section in experiment
139
140 static int getMode;
141 double Li2(double x);
142 double SoftPhoton_xs(double s,double b);
143 double lgr(double *x,double *y,int n,double t);
144 bool islgr(double *x,double *y,int n,double t);
145 double LLr(double *x,double *y,int n,double t);
146 int selectMode(std::vector<int> vmod, double mhds);
147 void findMaxXS(double mhds );
148 bool gam_sampling(double mhds,double sintheta);
149 void resetResMass();
150 void getResMass();
151 bool checkdecay(EvtParticle* p);
152 double sumExc(double mx);
153 void showResMass();
154 int getNdaugs(){return _ndaugs;}
155 EvtId* getDaugId(){return daugs;}
156 int getSelectedMode(){return _selectedMode;}
157 static int conexcmode;
158 double narrowRXS(double mxL,double mxH);
159 double selectMass();
160 double addNarrowRXS(double mhi, double binwidth);
161 void ReadVP();
162 double getVP(double cms);
163private:
164
165 int _mode,_ndaugs,radflag,testflag;
166 EvtId daugs[10],gamId;//daugs[0]~dagus[_ndaugs-1] are hadrons, daugs[_ndaugs] is ISR gamma
167 static double _xs0,_xs1; //cross section for 0 and 1-photon processes
168 static double _er0,_er1; //cross section for 0 and 1-photon processes
169 static int _nevt;
170 std::vector<double> ISRXS,ISRM;
171 std::vector<bool> ISRFLAG;
172 EvtParticle* gamH;
173 double maxXS;//maximum of diffrential cross section respective to cos\theta and mhds for ee->gamma hadrons
174 double differ,differ2,Rad2Xs;
175 std::string _unit;
176 std::vector<double> BR_ee;
177 std::vector<EvtId > ResId,ISRID;
178
179 double Egamcut;
180 //-- for deguggint to make root file
181 TFile *myfile;
182 Double_t pgam[4],phds[4],ph1[4],ph2[4],mhds,sumxs;
183 Double_t mass1,mass2,costheta,selectmode;
184 Double_t cosp,cosk;
185 Int_t imode;
186 TTree *xs;
187 bool mydbg; //handler w/o debugging
188 TGraphErrors *mygr;
189 TH1F* myth,*Xobs,*Xsum;
190 // for calulate the correction factor
191 int pdgcode;
192 std::string file;
193 EvtId son[10],pid;
194 int nson;
195 double vph; //for vaccuam polarization calculation
196 double AF[600],AA[600],MH[600],RadXS[600],EgamH;
197
198 double mjsi,mpsip,mpsipp,mphi,momega,mrho0,mrho3s,momega2s;
199 double wjsi,wpsip,wpsipp,wphi,womega,wrho0,wrho3s,womega2s;
200 double _mhdL;
201 double cmsspread;
202 int _selectedMode;
203 std::vector<int> _modeFlag;
204 bool VISR;
205
206 std::vector<double> vpx,vpr,vpi;
207};
208
209#endif
210
const Int_t n
XmlRpcServer s
Definition: HelloServer.cpp:11
bool checkdecay(EvtParticle *p)
Definition: EvtConExc.cc:2143
double Li2(double x)
Definition: EvtConExc.cc:1836
int getNdaugs()
Definition: EvtConExc.hh:154
void findMaxXS(EvtParticle *p)
Definition: EvtConExc.cc:1466
static int getMode
Definition: EvtConExc.hh:140
static int conexcmode
Definition: EvtConExc.hh:157
double addNarrowRXS(double mhi, double binwidth)
Definition: EvtConExc.cc:2300
void init()
Definition: EvtConExc.cc:187
int getSelectedMode()
Definition: EvtConExc.hh:156
double narrowRXS(double mxL, double mxH)
Definition: EvtConExc.cc:2266
bool VP_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:1592
void init_Br_ee()
Definition: EvtConExc.cc:1710
double gamHXSection_er(double El, double Eh)
Definition: EvtConExc.cc:1448
bool islgr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:1857
void showResMass()
Definition: EvtConExc.cc:2132
double lgr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:1842
void ReadVP()
Definition: EvtConExc.cc:2355
double baryonAng(double mx)
Definition: EvtConExc.cc:2225
double Ros_xs(double mx, double bree, EvtId pid)
Definition: EvtConExc.cc:1742
double Rad1(double s, double x)
Definition: EvtConExc.cc:1603
static EvtXsection * myxsection
Definition: EvtConExc.hh:136
void SetP4Rvalue(EvtParticle *part, double mhdr, double xeng, double theta)
Definition: EvtConExc.cc:1970
bool meson_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:1581
bool photonSampling(EvtParticle *part)
Definition: EvtConExc.cc:2235
double ISR_ang_integrate(double x, double theta)
Definition: EvtConExc.cc:1911
bool xs_sampling(double xs)
Definition: EvtConExc.cc:1552
int selectMode(std::vector< int > vmod, double mhds)
Definition: EvtConExc.cc:2014
double gamHXSection(EvtParticle *p, double El, double Eh, int nmc=100000)
Definition: EvtConExc.cc:1340
double Rad1difXs(EvtParticle *p)
Definition: EvtConExc.cc:1682
double sumExc(double mx)
Definition: EvtConExc.cc:2157
void init_mode(int mode)
Definition: EvtConExc.cc:576
bool baryon_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:1564
double LLr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:1872
double Rad2difXs(EvtParticle *p)
Definition: EvtConExc.cc:1639
double SoftPhoton_xs(double s, double b)
Definition: EvtConExc.cc:1811
bool angularSampling(EvtParticle *part)
Definition: EvtConExc.cc:2189
EvtId * getDaugId()
Definition: EvtConExc.hh:155
double difgamXs(EvtParticle *p)
Definition: EvtConExc.cc:1511
double ISR_ang_sampling(double x)
Definition: EvtConExc.cc:1929
void initProbMax()
Definition: EvtConExc.cc:1024
double Mhad_sampling(double *x, double *y)
Definition: EvtConExc.cc:1888
static double _cms
Definition: EvtConExc.hh:137
void resetResMass()
Definition: EvtConExc.cc:2061
bool hadron_angle_sampling(EvtVector4R ppi, EvtVector4R pcm)
Definition: EvtConExc.cc:1322
virtual ~EvtConExc()
Definition: EvtConExc.cc:165
void getResMass()
Definition: EvtConExc.cc:2096
bool gam_sampling(EvtParticle *p)
Definition: EvtConExc.cc:1538
double Rad2(double s, double x)
Definition: EvtConExc.cc:1615
double selectMass()
Definition: EvtConExc.cc:2309
double getVP(double cms)
Definition: EvtConExc.cc:2334
void getName(std::string &name)
Definition: EvtConExc.cc:174
void decay(EvtParticle *p)
Definition: EvtConExc.cc:1030
static double XS_max
Definition: EvtConExc.hh:138
EvtDecayBase * clone()
Definition: EvtConExc.cc:180
void SetP4(EvtParticle *part, double mhdr, double xeng, double theta)
Definition: EvtConExc.cc:1953
Definition: EvtId.hh:27
int t()
Definition: t.c:1