BOSS 7.1.1
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDtoKSpietaPlot.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: EvtDtoKSpietaPlot.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 "TAxis.h"
33#include "TH2.h"
34using std::endl;
35
37
38void EvtDtoKSpietaPlot::getName(std::string& model_name){
39 model_name="DtoKSpietaPlot";
40}
41
45
46
48
49 checkNArg(0);
50
51 bool idN = getDaugs()[0]==EvtPDL::getId(std::string("K_S0"))||getDaugs()[0]==EvtPDL::getId(std::string("K_L0"));
52 bool idKs = getDaugs()[1]==EvtPDL::getId(std::string("pi+"))||getDaugs()[1]==EvtPDL::getId(std::string("pi-")) ;
53 bool idPi = getDaugs()[2]==EvtPDL::getId(std::string("eta"));
54 if(!(idN && idKs && idPi ) ){std::cout<<"EvtDtoKSpietaPlot: the daughter sequence should be Ks/Kl pi+ eta"<<std::endl;abort();}
56
57Xmin = 0.406;
58Xmax = 1.747;
59Xwid = 0.1341;
60Ymin = 0.473;
61Ymax = 1.882;
62Ywid = 0.1409;
63avm1 = 0.240741;
64double HisPDFtmp[12][12] =
65{
66{0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0},
67{0, 1e-07, 0.0119599, 0.0296515, 0.0444187, 0.0196774, 0.00679957, 0.00950145, 0.00704225, 0.0833333, 1e-07, 0},
68{0, 0.00265756, 0.00594999, 0.04354, 0.0266106, 0.00888889, 0.0104167, 0.00903452, 0.0110731, 0.00418803, 1e-07, 0},
69{0, 0.00629828, 0.0145145, 0.0280041, 0.0236976, 0.0160962, 0.00604595, 0.0225225, 0.00838338, 0.00828402, 0.0102041, 0},
70{0, 0.0058566, 0.0214311, 0.0238095, 0.0421477, 0.0116178, 0.0188645, 0.0262201, 0.00917431, 0.00817757, 0.00310856, 0},
71{0, 0.00679852, 0.0101833, 0.0106561, 0.0188282, 0.0231757, 0.0246205, 0.0146294, 0.0124968, 0.0152171, 0.00802826, 0},
72{0, 0.01116, 0.00961875, 0.0153356, 0.0245832, 0.0207115, 0.0180312, 0.01606, 0.0126437, 0.00146566, 0.00326797, 0},
73{0, 0.005919, 0.00744745, 0.0114328, 0.0330096, 0.0306159, 0.0280269, 0.0182468, 0.0111238, 0.00387888, 1e-07, 0},
74{0, 0.00185185, 0.00144012, 0.0189673, 0.0478907, 0.0233768, 0.0250272, 0.0130943, 0.0166364, 1e-07, 1e-07, 0},
75{0, 1e-07, 0.00462388, 0.0170154, 0.0512206, 0.0407905, 0.0289623, 0.0307394, 1e-07, 1e-07, 1e-07, 0},
76{0, 1e-07, 0.00657895, 0.0152134, 0.0479368, 0.0486111, 0.0318021, 0.166667, 1e-07, 1e-07, 1e-07, 0},
77{0, 0, 0, 0, 0.0689655, 0.240741, 0, 0, 0, 0, 0, 0}
78};
79
80for(int i=0;i<12;i++){
81 for (int j=0; j<12; j++) {
82 HisPDF[i][j] = HisPDFtmp[i][j];
83 }
84}
85
86}
90
92
93loop:
95
96 EvtParticle *id1,*id2,*id3;
97 EvtVector4R pd1,pd2,pd3;
98 double xmass13,xmass12, xmass23;
99
100 id1 =p->getDaug(0);
101 id2 =p->getDaug(1);
102 id3 =p->getDaug(2);
103
104 pd1 =id1->getP4Lab();
105 pd2 =id2->getP4Lab();
106 pd3 =id3->getP4Lab();
107
108 xmass12=(pd1+pd2).mass()*(pd1+pd2).mass(); // M_ksopi
109// xmass13=(pd1+pd3).mass()*(pd1+pd3).mass(); // M_ksow
110 xmass23=(pd2+pd3).mass()*(pd2+pd3).mass(); // M_piw
111
112 int xbin = FindXBin(xmass12);
113 int ybin = FindYBin(xmass23);
114 double xratio12=HisPDF[xbin][ybin]/avm1;
115
116 if(xratio12 <=0) goto loop;
117
118 double rd12=EvtRandom::Flat(0.0, 1.0);
119 if(rd12>xratio12) goto loop;
120
121 return ;
122}
123
125 if (mass2 < Xmin) { return 0; }
126 else if (mass2>=Xmax) { return 11; }
127 else { return int((mass2-Xmin)/Xwid)+1; }
128}
129
131 if (mass2 < Ymin) { return 0; }
132 else if (mass2>=Ymax) { return 11; }
133 else { return int((mass2-Ymin)/Ywid)+1; }
134}
double mass
EvtId getParentId()
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
void getName(std::string &name)
EvtDecayBase * clone()
void decay(EvtParticle *p)
int FindYBin(double mass2)
int FindXBin(double mass2)
static EvtSpinType::spintype getSpinType(EvtId i)
Definition EvtPDL.hh:61
static EvtId getId(const std::string &name)
Definition EvtPDL.cc:287
EvtVector4R getP4Lab()
EvtParticle * getDaug(int i)
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
static double Flat()
Definition EvtRandom.cc:74