BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDtoKSpiomegaPlot.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: EvtDtoKSpiomegaPlot.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 EvtDtoKSpiomegaPlot::getName( std::string& model_name ) {
39 model_name = "DtoKSpiomegaPlot";
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( "pi+" ) ) ||
51 getDaugs()[1] == EvtPDL::getId( std::string( "pi-" ) );
52 bool idPi = getDaugs()[2] == EvtPDL::getId( std::string( "omega" ) );
53 if ( !( idN && idKs && idPi ) )
54 {
55 std::cout << "EvtDtoKSpiomegaPlot: the daughter sequence should be K_S0/K_L0 pi+ omega"
56 << std::endl;
57 abort();
58 }
60
61 Xmin = 0.405769;
62 Xmax = 1.18157;
63 Xwid = 0.096975;
64 Ymin = 0.850084;
65 Ymax = 1.88238;
66 Ywid = 0.129038;
67 avm1 = 0.235;
68 double HisPDFtmp[10][10] = {
69 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
70 { 0, 1e-08, 0.120721, 0.052381, 0.0146277, 1e-08, 0.0325342, 0.235, 1e-08, 0 },
71 { 0, 0.048494, 0.0426205, 0.057598, 0.0610795, 0.0782258, 0.0776099, 0.0327206,
72 0.00078125, 0 },
73 { 0.111111, 0.00310078, 0.0329032, 0.08, 0.0677989, 0.139752, 0.134032, 0.0708333,
74 0.0407328, 0 },
75 { -0.108333, 0.0421875, 0.0802147, 0.075, 0.0754011, 0.187206, 0.12803, 0.0151351,
76 0.0235849, -0.0928571 },
77 { 0, 0.0457386, 0.033237, 0.0667614, 0.126403, 0.173194, 0.108482, 0.0508333, 0.0395349,
78 -0.325 },
79 { 0, 0.05, 0.0637931, 0.0892663, 0.141294, 0.125867, 0.111432, 0.0390288, 0.0578571, 0 },
80 { 0, 1e-08, 1e-08, 0.0670455, 0.0672059, 0.0942187, 0.0969551, 0.0503425, 1e-08, 0 },
81 { 0, 1e-08, 1e-08, 0.0186937, 0.0445402, 0.027, 0.00136364, 1e-08, 1e-08, 0 },
82 { 0, 0, 0, -0.08125, 0.0769231, -0.0464286, 0, 0, 0, 0 } };
83
84 for ( int i = 0; i < 10; i++ )
85 {
86 for ( int j = 0; j < 10; j++ ) { HisPDF[i][j] = HisPDFtmp[i][j]; }
87 }
88}
89
91
93
94loop:
96
97 EvtParticle *id1, *id2, *id3;
98 EvtVector4R pd1, pd2, pd3;
99 double xmass13, xmass12, xmass23;
100
101 id1 = p->getDaug( 0 );
102 id2 = p->getDaug( 1 );
103 id3 = p->getDaug( 2 );
104
105 pd1 = id1->getP4Lab();
106 pd2 = id2->getP4Lab();
107 pd3 = id3->getP4Lab();
108
109 xmass12 = ( pd1 + pd2 ).mass() * ( pd1 + pd2 ).mass(); // M_ksopi
110 // xmass13=(pd1+pd3).mass()*(pd1+pd3).mass(); // M_ksow
111 xmass23 = ( pd2 + pd3 ).mass() * ( pd2 + pd3 ).mass(); // M_piw
112
113 int xbin = FindXBin( xmass12 );
114 int ybin = FindYBin( xmass23 );
115 double xratio12 = HisPDF[xbin][ybin] / avm1;
116
117 if ( xratio12 <= 0 ) goto loop;
118
119 double rd12 = EvtRandom::Flat( 0.0, 1.0 );
120 if ( rd12 > xratio12 ) goto loop;
121
122 return;
123}
124
125int EvtDtoKSpiomegaPlot::FindXBin( double mass2 ) {
126 if ( mass2 < Xmin ) { return 0; }
127 else if ( mass2 >= Xmax ) { return 9; }
128 else { return int( ( mass2 - Xmin ) / Xwid ) + 1; }
129}
130
131int EvtDtoKSpiomegaPlot::FindYBin( double mass2 ) {
132 if ( mass2 < Ymin ) { return 0; }
133 else if ( mass2 >= Ymax ) { return 9; }
134 else { return int( ( mass2 - Ymin ) / Ywid ) + 1; }
135}
double mass
EvtId getParentId()
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
int FindYBin(double mass2)
void getName(std::string &name)
void decay(EvtParticle *p)
int FindXBin(double mass2)
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