BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtD0toKpietaPlot.cc
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: EvtD0toKpietaPlot.cc
12//
13// Modification history:
14//
15// Liaoyuan Dong August, 2022 Module created
16//
17//------------------------------------------------------------------------
21#include "EvtGenBase/EvtPDL.hh"
25#include <stdlib.h>
26#include <string>
27
28#include "TFile.h"
29#include "TApplication.h"
30#include "TROOT.h"
31#include "TH1.h"
32#include "TH2.h"
33#include "TAxis.h"
34using std::endl;
35
37
38void EvtD0toKpietaPlot::getName(std::string& model_name){
39 model_name="D0toKpietaPlot";
40}
41
43 return new EvtD0toKpietaPlot;
44}
45
47
48 checkNArg(0);
49
50 bool idN = getDaugs()[0]==EvtPDL::getId(std::string("K-"))||getDaugs()[0]==EvtPDL::getId(std::string("K+"));
51 bool idKs = getDaugs()[1]==EvtPDL::getId(std::string("pi+"))||getDaugs()[1]==EvtPDL::getId(std::string("pi-")) ;
52 bool idPi = getDaugs()[2]==EvtPDL::getId(std::string("eta"));
53 if(!(idN && idKs && idPi ) ){std::cout<<"EvtD0toKpietaPlot: the daughter sequence should be K- pi+ eta"<<std::endl;abort();}
55
56Xmin = 0.401;
57Xmax = 1.734;
58Xwid = 0.0888667;
59Ymin = 0.473;
60Ymax = 1.88;
61Ywid = 0.0938;
62avm1 = 4.0;
63double HisPDFtmp[17][17] =
64{
65{0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0},
66{0, 1e-07, 0.026749, 0.00605681, 0.0634184, 0.0892256, 0.202502, 0.176861, 0.14297, 0.0871013, 0.0597803, 0.0816327, 1e-07, 1e-07, 1e-07, 1e-07, 0},
67{0, 0.0333333, 0.0170782, 0.0151869, 0.0481322, 0.0805231, 0.205394, 0.116217, 0.0875592, 0.129291, 0.0702469, 0.0752794, 0.100899, 0.0771605, 0.0869565, 1e-07, 0},
68{0, 0.0280093, 0.0390357, 0.0253937, 0.0261579, 0.071521, 0.141498, 0.104404, 0.0741497, 0.0714744, 0.0822128, 0.0975891, 0.155818, 0.159259, 0.185316, 0.148148, 0},
69{0, 0.113668, 0.119763, 0.0787132, 0.102704, 0.112065, 0.135281, 0.0820707, 0.0419041, 0.079242, 0.0900836, 0.150758, 0.264155, 0.316919, 0.317745, 0.431149, 0},
70{0, 0.589684, 0.371107, 0.293253, 0.168061, 0.120212, 0.187226, 0.103501, 0.0238857, 0.0322682, 0.0430699, 0.131763, 0.254585, 0.314903, 0.422812, 0.531974, -0.0925926},
71{0, 0.171536, 0.092114, 0.0382082, 0.0237244, 1e-07, 0.0313214, 0.025242, 0.0140927, 0.019018, 0.0117647, 0.00504009, 0.00752062, 0.0128317, 0.0120773, 0.00329218, 0},
72{0, 0.0398736, 0.0122582, 0.000112007, 0.0110455, 0.0132114, 0.0126157, 0.00567562, 0.00792236, 0.0138889, 0.0132507, 0.0179872, 0.00873761, 0.00735294, 1e-07, 1e-07, 0},
73{0, 0.0320175, 0.0250501, 0.0182082, 0.0198841, 0.0312431, 0.0178492, 0.00932836, 0.00484048, 0.0196806, 0.0125213, 0.0167646, 0.0242894, 0.0175224, 0.00443081, 0.012037, 0},
74{0, 0.0281924, 0.0348485, 0.055883, 0.025371, 0.0317179, 0.0222222, 0.0132327, 1e-07, 0.00288889, 0.00841751, 0.0165652, 0.0143026, 0.013714, 0.000664011, 1e-07, 0},
75{0, 0.0814623, 0.0519035, 0.0857828, 0.050823, 0.0422381, 0.038486, 0.00943396, 0.00109127, 0.003367, 0.00705829, 0.00633333, 0.00933246, 0.00857237, 0.0114638, 1e-07, 0},
76{0, 0.0945946, 0.103036, 0.0933035, 0.0920969, 0.0668599, 0.0459987, 0.0155599, 0.00218477, 1e-07, 0.00224791, 0.00819475, 0.00309358, 0.00826446, 1e-07, 1e-07, 0},
77{0, 1e-07, 0.0715244, 0.084398, 0.0704767, 0.0615041, 0.0658281, 0.00491735, 0.00484436, 0.00573099, 0.00257898, 0.003927, 0.0204082, 1e-07, 1e-07, 1e-07, 0},
78{0, 1e-07, 0.0763572, 0.0943127, 0.0825968, 0.0734052, 0.0765793, 0.0305351, 0.0121083, 0.00391802, 0.00722106, 1e-07, 1e-07, 1e-07, 1e-07, 1e-07, 0},
79{0, 1e-07, 1e-07, 0.111518, 0.0960134, 0.128118, 0.0622588, 0.0360153, 0.0130845, 0.0106838, 1e-07, 1e-07, 1e-07, 1e-07, 1e-07, 1e-07, 0},
80{0, 1e-07, 1e-07, 0.173913, 0.120973, 0.121418, 0.121403, 0.0548894, 0.0436284, 0.142857, 1e-07, 1e-07, 1e-07, 1e-07, 1e-07, 1e-07, 0},
81{0, 0, 0, 0, 0, 0.232639, 0.0972222, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}
82};
83
84for(int i=0;i<17;i++){
85 for (int j=0; j<17; j++) {
86 HisPDF[i][j] = HisPDFtmp[i][j];
87 }
88}
89
90}
91
93 noProbMax();
94}
95
97
98loop:
100
101 EvtParticle *id1,*id2,*id3;
102 EvtVector4R pd1,pd2,pd3;
103 double xmass13,xmass12, xmass23;
104
105 id1 =p->getDaug(0);
106 id2 =p->getDaug(1);
107 id3 =p->getDaug(2);
108
109 pd1 =id1->getP4Lab();
110 pd2 =id2->getP4Lab();
111 pd3 =id3->getP4Lab();
112
113 xmass12=(pd1+pd2).mass()*(pd1+pd2).mass(); // M_ksopi
114// xmass13=(pd1+pd3).mass()*(pd1+pd3).mass(); // M_ksow
115 xmass23=(pd2+pd3).mass()*(pd2+pd3).mass(); // M_piw
116
117 int xbin = FindXBin(xmass12);
118 int ybin = FindYBin(xmass23);
119 double xratio12=HisPDF[xbin][ybin]/avm1;
120
121 if(xratio12 <=0) goto loop;
122
123 double rd12=EvtRandom::Flat(0.0, 1.0);
124 if(rd12>xratio12) goto loop;
125
126 return ;
127}
128
130 if (mass2 < Xmin) { return 0; }
131 else if (mass2>=Xmax) { return 16; }
132 else { return int((mass2-Xmin)/Xwid)+1; }
133}
134
136 if (mass2 < Ymin) { return 0; }
137 else if (mass2>=Ymax) { return 16; }
138 else { return int((mass2-Ymin)/Ywid)+1; }
139}
double mass
int FindYBin(double mass2)
void getName(std::string &name)
int FindXBin(double mass2)
void decay(EvtParticle *p)
EvtDecayBase * clone()
virtual ~EvtD0toKpietaPlot()
void noProbMax()
EvtId getParentId()
Definition: EvtDecayBase.hh:60
EvtId * getDaugs()
Definition: EvtDecayBase.hh:65
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
static EvtSpinType::spintype getSpinType(EvtId i)
Definition: EvtPDL.hh:61
static EvtId getId(const std::string &name)
Definition: EvtPDL.cc:287
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:684
EvtParticle * getDaug(int i)
Definition: EvtParticle.cc:84
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
static double Flat()
Definition: EvtRandom.cc:74