BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtD0toKSpi0etaPlot.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: EvtD0toKspi0etaplot.cc
12//
13// Modification history:
14//
15// Liaoyuan Dong August, 2022 Module created
16//
17//------------------------------------------------------------------------
25#include <stdlib.h>
26#include <string>
27
28#include "TApplication.h"
29#include "TAxis.h"
30#include "TFile.h"
31#include "TH1.h"
32#include "TH2.h"
33#include "TROOT.h"
34using std::endl;
35
37
38void EvtD0toKSpi0etaPlot::getName( std::string& model_name ) {
39 model_name = "D0toKSpi0etaPlot";
40}
41
43
45
46 checkNArg( 0 );
47
48 bool idN = getDaugs()[0] == EvtPDL::getId( std::string( "K_S0" ) ) ||
49 getDaugs()[0] == EvtPDL::getId( std::string( "K_L0" ) );
50 bool idKs = getDaugs()[1] == EvtPDL::getId( std::string( "pi0" ) );
51 bool idPi = getDaugs()[2] == EvtPDL::getId( std::string( "eta" ) );
52 if ( !( idN && idKs && idPi ) )
53 {
54 std::cout << "EvtD0toKSpi0etaPlot: the daughter sequence should be Ks/Kl pi0 eta"
55 << std::endl;
56 abort();
57 }
58
60
61 Xmin = 0.4;
62 Xmax = 1.734;
63 Xwid = 0.1334;
64 Ymin = 0.466;
65 Ymax = 1.869;
66 Ywid = 0.1403;
67 avm1 = 0.125;
68 double HisPDFtmp[12][12] = {
69 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
70 { 0, 0.00527169, 0.0113872, 0.0301807, 0.0804598, 0.0617421, 0.0427827, 0.0540283,
71 0.0165289, 0, 0, 0 },
72 { 0, 0.020155, 0.0159653, 0.0296171, 0.0808365, 0.0416933, 0.0265935, 0.0193752,
73 0.0441176, 0.0259259, 0.0287037, 0 },
74 { 0, 0.0417141, 0.0513983, 0.0479749, 0.0659032, 0.0323741, 0.0175809, 0.0304969,
75 0.0254219, 0.054984, 0.0801812, 0.0493827 },
76 { 0, 0.0388568, 0.00614035, 0.0154096, 0.028749, 0.0281987, 0.0182064, 0.00844951,
77 0.00682977, 0.0027885, 0.00239316, -0.0173611 },
78 { 0, 0.0100519, 0.0169912, 0.0200723, 0.0272782, 0.0163995, 0.0231309, 0.0170477,
79 0.013863, 0.000892459, 0.00420168, 0 },
80 { 0, 0.0361953, 0.0251366, 0.0173922, 0.0399432, 0.00877736, 0.0105188, 0.0110497,
81 0.0107359, 0.00134903, 0.00448934, 0 },
82 { 0, 0.0078219, 0.0234282, 0.0324155, 0.0451149, 0.011002, 0.0149195, 0.00172265,
83 0.00821524, 0, 0, 0 },
84 { 0, 0.021721, 0.0280722, 0.0341469, 0.0507545, 0.0136752, 0.011148, 0.0100806,
85 0.00383142, 0, 0, 0 },
86 { 0, 0, 0.0237679, 0.0316915, 0.0674883, 0.018785, 0.0192308, 0.0037594, 0.047619, 0, 0,
87 0 },
88 { 0, 0, 0.0307692, 0.0461138, 0.0755956, 0.0250696, 0, 0, 0, 0, 0, 0 },
89 { 0, 0, 0, 0, 0.0938697, 0.125, 0, 0, 0, 0, 0, 0 },
90 };
91
92 for ( int i = 0; i < 12; i++ )
93 {
94 for ( int j = 0; j < 12; j++ ) { HisPDF[i][j] = HisPDFtmp[i][j]; }
95 }
96}
97
99
101
102loop:
104
105 EvtParticle *id1, *id2, *id3;
106 EvtVector4R pd1, pd2, pd3;
107 double xmass13, xmass12, xmass23;
108
109 id1 = p->getDaug( 0 );
110 id2 = p->getDaug( 1 );
111 id3 = p->getDaug( 2 );
112
113 pd1 = id1->getP4Lab();
114 pd2 = id2->getP4Lab();
115 pd3 = id3->getP4Lab();
116
117 xmass12 = ( pd1 + pd2 ).mass() * ( pd1 + pd2 ).mass(); // M_ksopi
118 // xmass13=(pd1+pd3).mass()*(pd1+pd3).mass(); // M_ksow
119 xmass23 = ( pd2 + pd3 ).mass() * ( pd2 + pd3 ).mass(); // M_piw
120
121 int xbin = FindXBin( xmass12 );
122 int ybin = FindYBin( xmass23 );
123 double xratio12 = HisPDF[xbin][ybin] / avm1;
124
125 if ( xratio12 <= 0 ) goto loop;
126
127 double rd12 = EvtRandom::Flat( 0.0, 1.0 );
128 if ( rd12 > xratio12 ) goto loop;
129
130 return;
131}
132
133int EvtD0toKSpi0etaPlot::FindXBin( double mass2 ) {
134 if ( mass2 < Xmin ) { return 0; }
135 else if ( mass2 >= Xmax ) { return 11; }
136 else { return int( ( mass2 - Xmin ) / Xwid ) + 1; }
137}
138
139int EvtD0toKSpi0etaPlot::FindYBin( double mass2 ) {
140 if ( mass2 < Ymin ) { return 0; }
141 else if ( mass2 >= Ymax ) { return 11; }
142 else { return int( ( mass2 - Ymin ) / Ywid ) + 1; }
143}
double mass
void getName(std::string &name)
int FindYBin(double mass2)
int FindXBin(double mass2)
void decay(EvtParticle *p)
EvtId getParentId()
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
static EvtSpinType::spintype getSpinType(EvtId i)
Definition EvtPDL.hh:66
static EvtId getId(const std::string &name)
Definition EvtPDL.cc:272
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:69