BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtBody3.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: EvtBody3.cc
12//
13// Description: Routine to decay a particle into three bodies using the Dalitz plots and two
14// particle
15// angular distributions.
16//
17// Modification history:
18//
19// Ping R.-G. December, 2006 Module created
20//
21//------------------------------------------------------------------------
22//
23#include "EvtBody3.hh"
43#include <stdlib.h>
44#include <string>
45
46#include "TApplication.h"
47#include "TAxis.h"
48#include "TFile.h"
49#include "TH1.h"
50#include "TH2.h"
51#include "TROOT.h"
52// #include "CLHEP/config/CLHEP.h"
53// #include "CLHEP/config/TemplateFunctions.h"
54
55// using namespace CLHEP;
56using std::endl;
57
59
60void EvtBody3::getName( std::string& model_name ) { model_name = "Body3"; }
61
63
65
66 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
67 checkNArg( 0 );
69}
71
73
74 const char* fl = setFileName();
75 const char* hp = setHpoint();
76 int* dp;
77 int nd1, nd2, ii, sn;
78
79 sn = setDaugAngNo();
80
81 if ( sn == 2 )
82 {
83 nd1 = 0;
84 nd2 = 1;
85 }
86
87 if ( sn == 0 )
88 {
89 nd1 = 1;
90 nd2 = 2;
91 }
92
93 if ( sn == 1 )
94 {
95 nd1 = 0;
96 nd2 = 2;
97 }
98 const char* DA1 = setDaugAng( nd1 );
99 const char* DA2 = setDaugAng( nd2 );
100
101 dp = setDaugPair();
102 int d1 = dp[0]; // m(d1,d2) pair at X axis
103 int d2 = dp[1];
104 int d3 = dp[2]; // m(d3,d4) pair at Y axis
105 int d4 = dp[3];
106
107 TFile f( fl );
108 TH1F* did1 = (TH1F*)f.Get( DA1 );
109 TH1F* did2 = (TH1F*)f.Get( DA2 );
110 TH2F* hid = (TH2F*)f.Get( hp );
111
112 TAxis* d1x = did1->GetXaxis();
113 TAxis* d2x = did2->GetXaxis();
114 TAxis* xaxis = hid->GetXaxis();
115 TAxis* yaxis = hid->GetYaxis();
116
117 int BINSd1 = did1->GetXaxis()->GetLast();
118 int BINSd2 = did2->GetXaxis()->GetLast();
119 int BINSx = xaxis->GetLast();
120 int BINSy = yaxis->GetLast();
121 int BINS = BINSx * BINSy;
122
123 double av1, av2, avm1, avm2;
124 avm1 = 0.;
125 avm2 = 0.;
126 double yvalue, ymax = 0.0;
127 int i, j, binxy;
128
129 // look for maxmum for the first hist.1
130 for ( i = 1; i < BINSd1 + 1; i++ )
131 {
132 av1 = did1->GetBinContent( i );
133 if ( av1 > avm1 ) avm1 = av1;
134 }
135
136 // look for maxmum for the second hist.1
137 for ( i = 1; i < BINSd2 + 1; i++ )
138 {
139 av2 = did2->GetBinContent( i );
140 if ( av2 > avm2 ) avm2 = av2;
141 }
142
143 // look for maxmum for the Dalitz plot
144
145 for ( i = 1; i < BINSx + 1; i++ )
146 {
147 for ( j = 1; j < BINSy + 1; j++ )
148 {
149 binxy = hid->GetBin( i, j );
150 yvalue = hid->GetBinContent( binxy );
151 // cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
152 if ( yvalue > ymax ) ymax = yvalue;
153 }
154 }
155
156loop:
158
159 EvtParticle *id1, *id2, *id3, *id4, *s1;
160 EvtVector4R pd1, pd2, pd3, pd4, ps;
161 EvtVector4R dp1, dp2;
162 double xmass2, ymass2;
163
164 id1 = p->getDaug( d1 );
165 id2 = p->getDaug( d2 );
166 id3 = p->getDaug( d3 );
167 id4 = p->getDaug( d4 );
168
169 pd1 = id1->getP4Lab();
170 pd2 = id2->getP4Lab();
171 pd3 = id3->getP4Lab();
172 pd4 = id4->getP4Lab();
173
174 dp1 = p->getDaug( nd1 )->getP4Lab();
175 dp2 = p->getDaug( nd2 )->getP4Lab();
176
177 xmass2 = ( pd1 + pd2 ).mass2();
178 ymass2 = ( pd3 + pd4 ).mass2();
179
180 int xbin = hid->GetXaxis()->FindBin( xmass2 );
181 int ybin = hid->GetYaxis()->FindBin( ymass2 );
182 int xybin = hid->GetBin( xbin, ybin );
183 double zvalue = hid->GetBinContent( xybin );
184 double xratio = zvalue / ymax;
185 if ( xratio == 0 ) goto loop;
186 double rd1 = EvtRandom::Flat( 0.0, 1.0 );
187 if ( rd1 > xratio ) goto loop; // sigle out event by 2-D mass distribution histo.
188
189 double da1 = dp1.get( 3 ) / dp1.d3mag();
190 double da2 = dp2.get( 3 ) / dp2.d3mag();
191
192 int dbin1 = did1->FindBin( da1 );
193 int dbin2 = did2->FindBin( da2 );
194
195 double dr1 = ( did1->GetBinContent( dbin1 ) ) / avm1;
196 double dr2 = ( did2->GetBinContent( dbin2 ) ) / avm2;
197 if ( dr1 == 0 || dr2 == 0 ) goto loop;
198 rd1 = EvtRandom::Flat( 0.0, 1.0 );
199 if ( rd1 > dr1 ) goto loop; // sigle out event by the first angular distribution histo.
200
201 rd1 = EvtRandom::Flat( 0.0, 1.0 );
202 if ( rd1 > dr2 ) goto loop; // sigle out event by the second angular distribution histo.
203
204 return;
205}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
EvtDecayBase * clone()
Definition EvtBody3.cc:62
const char * setFileName()
Definition UserBody3.cc:11
void getName(std::string &name)
Definition EvtBody3.cc:60
void init()
Definition EvtBody3.cc:64
void decay(EvtParticle *p)
Definition EvtBody3.cc:72
int * setDaugPair()
Definition UserBody3.cc:42
const char * setDaugAng(int i)
Definition UserBody3.cc:23
virtual ~EvtBody3()
Definition EvtBody3.cc:58
int setDaugAngNo()
Definition UserBody3.cc:37
void initProbMax()
Definition EvtBody3.cc:70
const char * setHpoint()
Definition UserBody3.cc:17
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
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
double get(int i) const
double d3mag() const