BOSS 6.6.4.p01
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtMassH2.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 Dalitz
14// plot.
15//
16// Modification history:
17//
18// Ping R.-G. December, 2006 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 EvtMassH2::getName(std::string& model_name){
60
61 model_name="MassH2";
62
63}
64
66
67 return new EvtMassH2;
68
69}
70
71
73
74 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
75 checkNArg(0);
77}
79
80 noProbMax();
81
82}
83
85
86 const char* fl=setFileName();
87 const char* hp=setHpoint();
88 int* dp;
89 dp=setDaugPair();
90 int d1=dp[0]; //m(d1,d2) pair at X axis
91 int d2=dp[1];
92 int d3=dp[2]; //m(d3,d4) pair at Y axis
93 int d4=dp[3];
94
95 TFile f(fl);
96 TH2F *hid = (TH2F*)f.Get(hp);
97 TAxis* xaxis = hid->GetXaxis();
98 TAxis* yaxis = hid->GetYaxis();
99
100 int BINSx =xaxis->GetLast();
101 int BINSy =yaxis->GetLast();
102 int BINS =BINSx*BINSy;
103 double yvalue,ymax=0.0;
104 int i,j,binxy;
105
106 for(i=1;i<BINSx+1;i++){
107 for(j=1;j<BINSy+1;j++){
108 binxy=hid->GetBin(i,j);
109 yvalue=hid->GetBinContent(binxy);
110// cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
111 if(yvalue>ymax) ymax=yvalue;
112 }
113 }
114
115loop:
117
118 EvtParticle *id1,*id2,*id3,*id4,*s1;
119 EvtVector4R pd1,pd2,pd3,pd4,ps;
120 double xmass2,ymass2;
121
122 id1 =p->getDaug(d1);
123 id2 =p->getDaug(d2);
124 id3 =p->getDaug(d3);
125 id4 =p->getDaug(d4);
126
127 pd1 =id1->getP4Lab();
128 pd2 =id2->getP4Lab();
129 pd3 =id3->getP4Lab();
130 pd4 =id4->getP4Lab();
131
132 xmass2=(pd1+pd2).mass2();
133 ymass2=(pd3+pd4).mass2();
134
135 int xbin = hid->GetXaxis()->FindBin(xmass2);
136 int ybin = hid->GetYaxis()->FindBin(ymass2);
137 int xybin= hid->GetBin(xbin,ybin);
138 double zvalue=hid->GetBinContent(xybin);
139 double xratio=zvalue/ymax;
140// cout<<"***************zvalue,ymax,xratio= "<<zvalue<<"; "<<ymax<<"; "<<xratio<<endl;
141 double rd1=EvtRandom::Flat(0.0, 1.0);
142 if(rd1>xratio) goto loop;
143 return ;
144}
145
146
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)
virtual ~EvtMassH2()
Definition: EvtMassH2.cc:57
void initProbMax()
Definition: EvtMassH2.cc:78
EvtDecayBase * clone()
Definition: EvtMassH2.cc:65
void init()
Definition: EvtMassH2.cc:72
const char * setHpoint()
Definition: UserMassH2.cc:16
int * setDaugPair()
Definition: UserMassH2.cc:22
void decay(EvtParticle *p)
Definition: EvtMassH2.cc:84
void getName(std::string &name)
Definition: EvtMassH2.cc:59
const char * setFileName()
Definition: UserMassH2.cc:10
static EvtSpinType::spintype getSpinType(EvtId i)
Definition: EvtPDL.hh:61
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:683
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:73