BOSS 7.0.2
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(double El, double Eh,int mode);
111 double gamHXSection_er(double El,double Eh);
112
113 void findMaxXS(EvtParticle *p);
114 double difgamXs(EvtParticle* p); //differential cross section for gamma + hadrons
115 double difgamXs(double mhds,double sintheta);
116 double Mhad_sampling(double *x,double *y);
117 double ISR_ang_integrate(double x,double theta);
118 double ISR_ang_sampling(double x);
120 void SetP4(EvtParticle *part, double mhdr,double xeng, double theta); //set the gamma and gamma* momentum according sampled results
121 void SetP4Rvalue(EvtParticle *part, double mhdr,double xeng, double theta); //set the gamma and gamma* momentum according sampled results
122 bool gam_sampling(EvtParticle *p);
123 bool xs_sampling(double xs);
124 bool xs_sampling(double xs,double xs1);
125 bool baryon_sampling(EvtVector4R pcm, EvtVector4R pi);//baryon angular distri. 1+cos^2\theta
126 bool meson_sampling(EvtVector4R pcm, EvtVector4R pi); //meson angular distri. sin^2\theta
127 bool VP_sampling(EvtVector4R pcm, EvtVector4R pi); //VP ang. dist. 1+cos^2theta
128 bool angularSampling(EvtParticle* part);
129 bool photonSampling(EvtParticle* part);
130 double baryonAng(double mx);
131 double Rad1(double s, double x);
132 double Rad2(double s, double x);
133 double Rad2difXs(EvtParticle *p);
134 double Rad2difXs(double s, double x);
135 double Rad1difXs(EvtParticle *p);
136 double Ros_xs(double mx, double bree,EvtId pid);
138 static double _cms; //energy of CMS of ee beam
139 static double XS_max;// maxium of cross section in experiment
140 static double SetMthr;
141
142 static int getMode;
143 double Li2(double x);
144 double SoftPhoton_xs(double s,double b);
145 double lgr(double *x,double *y,int n,double t);
146 bool islgr(double *x,double *y,int n,double t);
147 double LLr(double *x,double *y,int n,double t);
148 int selectMode(std::vector<int> vmod, double mhds);
149 void findMaxXS(double mhds );
150 bool gam_sampling(double mhds,double sintheta);
151 void resetResMass();
152 void getResMass();
153 bool checkdecay(EvtParticle* p);
154 double sumExc(double mx);
155 void showResMass();
156 int getNdaugs(){return _ndaugs;}
157 EvtId* getDaugId(){return daugs;}
158 int getSelectedMode(){return _selectedMode;}
159 static int conexcmode;
160 double narrowRXS(double mxL,double mxH);
161 double selectMass();
162 double addNarrowRXS(double mhi, double binwidth);
163 void ReadVP();
164 double getVP(double cms);
165 void mk_VXS(double Esig,double Egamcut,double EgamH,int midx); //make a set of observed cross section for exclusive processes,including narrow resonance
166 int get_mode_index(int mode);
167 double getObsXsection(double mhds,int mode);
168 double Egam2Mhds(double Egam);
169 std::vector<EvtId> get_mode(int mode);
170 void writeDecayTabel();
171 void checkEvtRatio();
172 int getModeIndex(int m);
173private:
174
175 int _mode,_ndaugs,radflag,testflag;
176 EvtId daugs[10],gamId;//daugs[0]~dagus[_ndaugs-1] are hadrons, daugs[_ndaugs] is ISR gamma
177 static double _xs0,_xs1; //cross section for 0 and 1-photon processes
178 static double _er0,_er1; //cross section for 0 and 1-photon processes
179 static int _nevt;
180 std::vector<double> ISRXS,ISRM;
181 std::vector<bool> ISRFLAG;
182 EvtParticle* gamH;
183 double maxXS;//maximum of diffrential cross section respective to cos\theta and mhds for ee->gamma hadrons
184 double differ,differ2,Rad2Xs;
185 std::string _unit;
186 std::vector<double> BR_ee;
187 std::vector<EvtId > ResId,ISRID;
188
189 double Egamcut;
190 //-- for deguggint to make root file
191 TFile *myfile;
192 Double_t pgam[4],phds[4],ph1[4],ph2[4],mhds,sumxs;
193 Double_t mass1,mass2,costheta,selectmode;
194 Double_t cosp,cosk;
195 Int_t imode;
196 TTree *xs;
197 bool mydbg; //handler w/o debugging
198 TGraphErrors *mygr;
199 TH1F* myth,*Xobs,*Xsum;
200 // for calulate the correction factor
201 int pdgcode;
202 std::string file;
203 EvtId son[10],pid;
204 int nson;
205 double vph; //for vaccuam polarization calculation
206 double AF[600],AA[600],MH[600],RadXS[600],EgamH;
207
208 double mjsi,mpsip,mpsipp,mphi,momega,mrho0,mrho3s,momega2s;
209 double wjsi,wpsip,wpsipp,wphi,womega,wrho0,wrho3s,womega2s;
210 double _mhdL;
211 double cmsspread;
212 int _selectedMode;
213 std::vector<int> _modeFlag;
214 bool VISR;
215 std::vector<int > vmode;
216 static std::vector<std::vector <double> > VXS;
217 std::vector<double> vpx,vpr,vpi;
218 double Mthr;
219 EvtParticle *theparent;
220};
221
222#endif
223
const Int_t n
XmlRpcServer s
Definition: HelloServer.cpp:11
bool checkdecay(EvtParticle *p)
Definition: EvtConExc.cc:2728
double Li2(double x)
Definition: EvtConExc.cc:2391
int getNdaugs()
Definition: EvtConExc.hh:156
void findMaxXS(EvtParticle *p)
Definition: EvtConExc.cc:2021
static int getMode
Definition: EvtConExc.hh:142
static int conexcmode
Definition: EvtConExc.hh:159
double addNarrowRXS(double mhi, double binwidth)
Definition: EvtConExc.cc:2889
void init()
Definition: EvtConExc.cc:189
int getSelectedMode()
Definition: EvtConExc.hh:158
double narrowRXS(double mxL, double mxH)
Definition: EvtConExc.cc:2855
bool VP_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:2147
void init_Br_ee()
Definition: EvtConExc.cc:2265
double gamHXSection_er(double El, double Eh)
Definition: EvtConExc.cc:2003
bool islgr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:2412
void showResMass()
Definition: EvtConExc.cc:2717
double lgr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:2397
void ReadVP()
Definition: EvtConExc.cc:2944
double baryonAng(double mx)
Definition: EvtConExc.cc:2814
double Ros_xs(double mx, double bree, EvtId pid)
Definition: EvtConExc.cc:2297
double Rad1(double s, double x)
Definition: EvtConExc.cc:2158
static EvtXsection * myxsection
Definition: EvtConExc.hh:137
void SetP4Rvalue(EvtParticle *part, double mhdr, double xeng, double theta)
Definition: EvtConExc.cc:2525
int getModeIndex(int m)
Definition: EvtConExc.cc:3128
bool meson_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:2136
bool photonSampling(EvtParticle *part)
Definition: EvtConExc.cc:2824
double ISR_ang_integrate(double x, double theta)
Definition: EvtConExc.cc:2466
static double SetMthr
Definition: EvtConExc.hh:140
bool xs_sampling(double xs)
Definition: EvtConExc.cc:2107
int selectMode(std::vector< int > vmod, double mhds)
Definition: EvtConExc.cc:2569
double gamHXSection(EvtParticle *p, double El, double Eh, int nmc=100000)
Definition: EvtConExc.cc:1876
double Rad1difXs(EvtParticle *p)
Definition: EvtConExc.cc:2237
double sumExc(double mx)
Definition: EvtConExc.cc:2742
void init_mode(int mode)
Definition: EvtConExc.cc:624
bool baryon_sampling(EvtVector4R pcm, EvtVector4R pi)
Definition: EvtConExc.cc:2119
int get_mode_index(int mode)
Definition: EvtConExc.cc:2997
double LLr(double *x, double *y, int n, double t)
Definition: EvtConExc.cc:2427
double Rad2difXs(EvtParticle *p)
Definition: EvtConExc.cc:2194
double SoftPhoton_xs(double s, double b)
Definition: EvtConExc.cc:2366
bool angularSampling(EvtParticle *part)
Definition: EvtConExc.cc:2778
EvtId * getDaugId()
Definition: EvtConExc.hh:157
double difgamXs(EvtParticle *p)
Definition: EvtConExc.cc:2066
double ISR_ang_sampling(double x)
Definition: EvtConExc.cc:2484
void checkEvtRatio()
Definition: EvtConExc.cc:3100
void initProbMax()
Definition: EvtConExc.cc:1544
double Mhad_sampling(double *x, double *y)
Definition: EvtConExc.cc:2443
static double _cms
Definition: EvtConExc.hh:138
void writeDecayTabel()
Definition: EvtConExc.cc:3025
std::vector< EvtId > get_mode(int mode)
Definition: EvtConExc.cc:1080
void resetResMass()
Definition: EvtConExc.cc:2646
double getObsXsection(double mhds, int mode)
Definition: EvtConExc.cc:3006
void mk_VXS(double Esig, double Egamcut, double EgamH, int midx)
Definition: EvtConExc.cc:2967
bool hadron_angle_sampling(EvtVector4R ppi, EvtVector4R pcm)
Definition: EvtConExc.cc:1858
virtual ~EvtConExc()
Definition: EvtConExc.cc:166
void getResMass()
Definition: EvtConExc.cc:2681
bool gam_sampling(EvtParticle *p)
Definition: EvtConExc.cc:2093
double Rad2(double s, double x)
Definition: EvtConExc.cc:2170
double selectMass()
Definition: EvtConExc.cc:2898
double getVP(double cms)
Definition: EvtConExc.cc:2923
void getName(std::string &name)
Definition: EvtConExc.cc:176
void decay(EvtParticle *p)
Definition: EvtConExc.cc:1550
static double XS_max
Definition: EvtConExc.hh:139
double Egam2Mhds(double Egam)
Definition: EvtConExc.cc:3019
EvtDecayBase * clone()
Definition: EvtConExc.cc:182
void SetP4(EvtParticle *part, double mhdr, double xeng, double theta)
Definition: EvtConExc.cc:2508
Definition: EvtId.hh:27
int t()
Definition: t.c:1