BOSS
7.0.3
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtmH2.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: EvtMassH1.cc
12
//
13
// Description: Routine to decay a particle using a scatter
14
// plot forn n-body decays (n>3).
15
//
16
// Modification history:
17
//
18
// Ping R.-G. Oct. 2011 Module created
19
//
20
//------------------------------------------------------------------------
21
//
22
#include "
EvtGenBase/EvtPatches.hh
"
23
#include <stdlib.h>
24
#include "
EvtGenBase/EvtParticle.hh
"
25
#include "
EvtGenBase/EvtGenKine.hh
"
26
#include "
EvtGenBase/EvtPDL.hh
"
27
#include "
EvtGenBase/EvtVector4C.hh
"
28
#include "
EvtGenBase/EvtVector4R.hh
"
29
#include "
EvtGenBase/EvtTensor4C.hh
"
30
#include "
EvtGenBase/EvtDiracParticle.hh
"
31
#include "
EvtGenBase/EvtScalarParticle.hh
"
32
#include "
EvtGenBase/EvtVectorParticle.hh
"
33
#include "
EvtGenBase/EvtTensorParticle.hh
"
34
#include "
EvtGenBase/EvtPhotonParticle.hh
"
35
#include "
EvtGenBase/EvtNeutrinoParticle.hh
"
36
#include "
EvtGenBase/EvtStringParticle.hh
"
37
#include "
EvtGenBase/EvtRaritaSchwingerParticle.hh
"
38
#include "
EvtGenBase/EvtHighSpinParticle.hh
"
39
#include "
EvtGenBase/EvtReport.hh
"
40
#include "
EvtGenBase/EvtHelSys.hh
"
41
#include "
EvtGenModels/EvtmH2.hh
"
42
#include "
EvtGenBase/EvtRandom.hh
"
43
#include <string>
44
45
#include "TH1.h"
46
#include "TAxis.h"
47
#include "TH2.h"
48
#include "TFile.h"
49
#include "TApplication.h"
50
#include "TROOT.h"
51
//#include "CLHEP/config/CLHEP.h"
52
//#include "CLHEP/config/TemplateFunctions.h"
53
54
55
using
std::endl;
56
57
EvtmH2::~EvtmH2
() {}
58
59
void
EvtmH2::getName
(std::string& model_name){
60
61
model_name=
"mH2"
;
62
63
}
64
65
EvtDecayBase
*
EvtmH2::clone
(){
66
67
return
new
EvtmH2
;
68
69
}
70
71
72
void
EvtmH2::init
(){
73
74
// check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
75
checkNArg
(2);
76
EvtSpinType::spintype
parenttype =
EvtPDL::getSpinType
(
getParentId
());
77
nbx =
getArg
(0);
78
nby =
getArg
(1);
79
}
80
void
EvtmH2::initProbMax
(){
81
82
noProbMax
();
83
84
}
85
86
void
EvtmH2::decay
(
EvtParticle
*p ){
87
88
const
char
* fl=
setFileName
();
89
const
char
* hp=
setHpoint
();
90
91
TFile f(fl);
92
TH2F *hid = (TH2F*)f.Get(
"mH2"
);
93
TAxis* xaxis = hid->GetXaxis();
94
TAxis* yaxis = hid->GetYaxis();
95
96
int
BINSx =xaxis->GetLast();
97
int
BINSy =yaxis->GetLast();
98
int
BINS =BINSx*BINSy;
99
double
yvalue,ymax=0.0;
100
int
i,j,binxy;
101
102
for
(i=1;i<BINSx+1;i++){
103
for
(j=1;j<BINSy+1;j++){
104
binxy=hid->GetBin(i,j);
105
yvalue=hid->GetBinContent(binxy);
106
// cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
107
if
(yvalue>ymax) ymax=yvalue;
108
}
109
}
110
111
loop:
112
p->
initializePhaseSpace
(
getNDaug
(),
getDaugs
());
113
114
if
(p->
getNDaug
()!= nbx+nby) {std::cout<<
"The number of specified particles is not equal the number of decay daughters "
<<endl;::abort();}
115
116
EvtVector4R
pt,pt2;
117
double
xmass
,ymass;
118
119
pt=p->
getDaug
(0)->
getP4Lab
();
120
for
(
int
ii=1;ii<nbx;ii++){
121
pt=pt+p->
getDaug
(ii)->
getP4Lab
();
122
}
123
xmass
=pt.
mass
();
124
125
pt2=p->
getDaug
(nbx)->
getP4Lab
();
126
for
(
int
jj=nbx+1;jj<nbx+nby;jj++) pt2=pt2+p->
getDaug
(jj)->
getP4Lab
();
127
ymass=pt2.
mass
();
128
129
130
int
xbin = hid->GetXaxis()->FindBin(
xmass
);
131
int
ybin = hid->GetYaxis()->FindBin(ymass);
132
int
xybin= hid->GetBin(xbin,ybin);
133
double
zvalue=hid->GetBinContent(xybin);
134
double
xratio=zvalue/ymax;
135
// cout<<"***************zvalue,ymax,xratio= "<<zvalue<<"; "<<ymax<<"; "<<xratio<<endl;
136
double
rd1=
EvtRandom::Flat
(0.0, 1.0);
137
if
(rd1>xratio)
goto
loop;
138
return ;
139
}
140
141
EvtDiracParticle.hh
EvtGenKine.hh
EvtHelSys.hh
EvtHighSpinParticle.hh
EvtNeutrinoParticle.hh
EvtPDL.hh
EvtParticle.hh
EvtPatches.hh
EvtPhotonParticle.hh
EvtRandom.hh
EvtRaritaSchwingerParticle.hh
EvtReport.hh
EvtScalarParticle.hh
EvtStringParticle.hh
EvtTensor4C.hh
EvtTensorParticle.hh
EvtVector4C.hh
EvtVector4R.hh
EvtVectorParticle.hh
EvtmH2.hh
xmass
const double xmass[5]
Definition:
Gam4pikp.cxx:50
EvtDecayBase
Definition:
EvtDecayBase.hh:33
EvtDecayBase::getArg
double getArg(int j)
Definition:
EvtDecayBase.cc:564
EvtDecayBase::noProbMax
void noProbMax()
Definition:
EvtDecayBase.cc:304
EvtDecayBase::getParentId
EvtId getParentId()
Definition:
EvtDecayBase.hh:60
EvtDecayBase::getNDaug
int getNDaug()
Definition:
EvtDecayBase.hh:64
EvtDecayBase::getDaugs
EvtId * getDaugs()
Definition:
EvtDecayBase.hh:65
EvtDecayBase::checkNArg
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
Definition:
EvtDecayBase.cc:482
EvtPDL::getSpinType
static EvtSpinType::spintype getSpinType(EvtId i)
Definition:
EvtPDL.hh:61
EvtParticle
Definition:
EvtParticle.hh:42
EvtParticle::getP4Lab
EvtVector4R getP4Lab()
Definition:
EvtParticle.cc:685
EvtParticle::getNDaug
int getNDaug() const
Definition:
EvtParticle.cc:125
EvtParticle::getDaug
EvtParticle * getDaug(int i)
Definition:
EvtParticle.cc:85
EvtParticle::initializePhaseSpace
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
Definition:
EvtParticle.cc:1071
EvtRandom::Flat
static double Flat()
Definition:
EvtRandom.cc:73
EvtSpinType::spintype
spintype
Definition:
EvtSpinType.hh:31
EvtVector4R
Definition:
EvtVector4R.hh:29
EvtVector4R::mass
double mass() const
Definition:
EvtVector4R.cc:39
EvtmH2::initProbMax
void initProbMax()
Definition:
EvtmH2.cc:80
EvtmH2::init
void init()
Definition:
EvtmH2.cc:72
EvtmH2::clone
EvtDecayBase * clone()
Definition:
EvtmH2.cc:65
EvtmH2::setFileName
const char * setFileName()
Definition:
UsermH2.cc:10
EvtmH2::getName
void getName(std::string &name)
Definition:
EvtmH2.cc:59
EvtmH2::setHpoint
const char * setHpoint()
Definition:
UsermH2.cc:16
EvtmH2::EvtmH2
EvtmH2()
Definition:
EvtmH2.hh:33
EvtmH2::decay
void decay(EvtParticle *p)
Definition:
EvtmH2.cc:86
EvtmH2::~EvtmH2
virtual ~EvtmH2()
Definition:
EvtmH2.cc:57
source
Generator
BesEvtGen
BesEvtGen-00-03-90
src
EvtGen
EvtGenModels
EvtmH2.cc
Generated by
1.9.6