BOSS 7.0.6
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//
23#include <stdlib.h>
26#include "EvtGenBase/EvtPDL.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
55using std::endl;
56
58
59void EvtmH2::getName(std::string& model_name){
60
61 model_name="mH2";
62
63}
64
66
67 return new EvtmH2;
68
69}
70
71
73
74 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
75 checkNArg(2);
77 nbx = getArg(0);
78 nby = getArg(1);
79}
81
82 noProbMax();
83
84}
85
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
111loop:
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
const double xmass[5]
Definition: Gam4pikp.cxx:50
double getArg(int j)
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
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:685
int getNDaug() const
Definition: EvtParticle.cc:125
EvtParticle * getDaug(int i)
Definition: EvtParticle.cc:85
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
static double Flat()
Definition: EvtRandom.cc:74
double mass() const
Definition: EvtVector4R.cc:39
void initProbMax()
Definition: EvtmH2.cc:80
void init()
Definition: EvtmH2.cc:72
EvtDecayBase * clone()
Definition: EvtmH2.cc:65
const char * setFileName()
Definition: UsermH2.cc:10
void getName(std::string &name)
Definition: EvtmH2.cc:59
const char * setHpoint()
Definition: UsermH2.cc:16
EvtmH2()
Definition: EvtmH2.hh:33
void decay(EvtParticle *p)
Definition: EvtmH2.cc:86
virtual ~EvtmH2()
Definition: EvtmH2.cc:57
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")