BOSS 7.0.4
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtBody3 Class Reference

#include <EvtBody3.hh>

+ Inheritance diagram for EvtBody3:

Public Member Functions

 EvtBody3 ()
 
virtual ~EvtBody3 ()
 
void getName (std::string &name)
 
EvtDecayBaseclone ()
 
void initProbMax ()
 
void init ()
 
void decay (EvtParticle *p)
 
const char * setFileName ()
 
const char * setHpoint ()
 
const char * setDaugAng (int i)
 
int setDaugAngNo ()
 
int * setDaugPair ()
 
- Public Member Functions inherited from EvtDecayIncoherent
void makeDecay (EvtParticle *p)
 
virtual ~EvtDecayIncoherent ()
 
void setDaughterSpinDensity (int daughter)
 
int isDaughterSpinDensitySet (int daughter)
 
- Public Member Functions inherited from EvtDecayBase
virtual void getName (std::string &name)=0
 
virtual void decay (EvtParticle *p)=0
 
virtual void makeDecay (EvtParticle *p)=0
 
virtual EvtDecayBaseclone ()=0
 
virtual void init ()
 
virtual void initProbMax ()
 
virtual std::string commandName ()
 
virtual void command (std::string cmd)
 
double getProbMax (double prob)
 
double resetProbMax (double prob)
 
 EvtDecayBase ()
 
virtual ~EvtDecayBase ()
 
virtual bool matchingDecay (const EvtDecayBase &other) const
 
EvtId getParentId ()
 
double getBranchingFraction ()
 
void disableCheckQ ()
 
void checkQ ()
 
int getNDaug ()
 
EvtIdgetDaugs ()
 
EvtId getDaug (int i)
 
int getNArg ()
 
int getPHOTOS ()
 
void setPHOTOS ()
 
void setVerbose ()
 
void setSummary ()
 
double * getArgs ()
 
std::string * getArgsStr ()
 
double getArg (int j)
 
std::string getArgStr (int j)
 
std::string getModelName ()
 
int getDSum ()
 
int summary ()
 
int verbose ()
 
void saveDecayInfo (EvtId ipar, int ndaug, EvtId *daug, int narg, std::vector< std::string > &args, std::string name, double brfr)
 
void printSummary ()
 
void setProbMax (double prbmx)
 
void noProbMax ()
 
void checkNArg (int a1, int a2=-1, int a3=-1, int a4=-1)
 
void checkNDaug (int d1, int d2=-1)
 
void checkSpinParent (EvtSpinType::spintype sp)
 
void checkSpinDaughter (int d1, EvtSpinType::spintype sp)
 
virtual int nRealDaughters ()
 

Additional Inherited Members

- Static Public Member Functions inherited from EvtDecayBase
static void findMasses (EvtParticle *p, int ndaugs, EvtId daugs[10], double masses[10])
 
static void findMass (EvtParticle *p)
 
static double findMaxMass (EvtParticle *p)
 
- Protected Member Functions inherited from EvtDecayBase
bool daugsDecayedByParentModel ()
 
- Protected Attributes inherited from EvtDecayBase
bool _daugsDecayedByParentModel
 

Detailed Description

Definition at line 28 of file EvtBody3.hh.

Constructor & Destructor Documentation

◆ EvtBody3()

EvtBody3::EvtBody3 ( )
inline

Definition at line 32 of file EvtBody3.hh.

32{}

Referenced by clone().

◆ ~EvtBody3()

EvtBody3::~EvtBody3 ( )
virtual

Definition at line 57 of file EvtBody3.cc.

57{}

Member Function Documentation

◆ clone()

EvtDecayBase * EvtBody3::clone ( )
virtual

Implements EvtDecayBase.

Definition at line 65 of file EvtBody3.cc.

65 {
66
67 return new EvtBody3;
68
69}
EvtBody3()
Definition: EvtBody3.hh:32

◆ decay()

void EvtBody3::decay ( EvtParticle p)
virtual

Implements EvtDecayBase.

Definition at line 84 of file EvtBody3.cc.

84 {
85
86 const char* fl=setFileName();
87 const char* hp=setHpoint();
88 int* dp;
89 int nd1,nd2,ii,sn;
90
91
92 sn=setDaugAngNo();
93
94 if(sn==2){
95 nd1=0;
96 nd2=1;
97 }
98
99 if(sn==0){
100 nd1=1;
101 nd2=2;
102 }
103
104 if(sn==1){
105 nd1=0;
106 nd2=2;
107 }
108 const char* DA1=setDaugAng(nd1);
109 const char* DA2=setDaugAng(nd2);
110
111
112 dp=setDaugPair();
113 int d1=dp[0]; //m(d1,d2) pair at X axis
114 int d2=dp[1];
115 int d3=dp[2]; //m(d3,d4) pair at Y axis
116 int d4=dp[3];
117
118
119 TFile f(fl);
120 TH1F *did1= (TH1F*)f.Get(DA1);
121 TH1F *did2= (TH1F*)f.Get(DA2);
122 TH2F *hid = (TH2F*)f.Get(hp);
123
124 TAxis* d1x=did1->GetXaxis();
125 TAxis* d2x=did2->GetXaxis();
126 TAxis* xaxis = hid->GetXaxis();
127 TAxis* yaxis = hid->GetYaxis();
128
129 int BINSd1 =did1->GetXaxis()->GetLast();
130 int BINSd2 =did2->GetXaxis()->GetLast();
131 int BINSx =xaxis->GetLast();
132 int BINSy =yaxis->GetLast();
133 int BINS =BINSx*BINSy;
134
135 double av1,av2,avm1,avm2;
136 avm1=0.;
137 avm2=0.;
138 double yvalue,ymax=0.0;
139 int i,j,binxy;
140
141// look for maxmum for the first hist.1
142 for(i=1;i<BINSd1+1;i++){
143 av1=did1->GetBinContent(i);
144 if(av1>avm1) avm1=av1;
145 }
146
147// look for maxmum for the second hist.1
148 for(i=1;i<BINSd2+1;i++){
149 av2=did2->GetBinContent(i);
150 if(av2>avm2) avm2=av2;
151 }
152
153// look for maxmum for the Dalitz plot
154
155 for(i=1;i<BINSx+1;i++){
156 for(j=1;j<BINSy+1;j++){
157 binxy=hid->GetBin(i,j);
158 yvalue=hid->GetBinContent(binxy);
159// cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
160 if(yvalue>ymax) ymax=yvalue;
161 }
162 }
163
164loop:
166
167 EvtParticle *id1,*id2,*id3,*id4,*s1;
168 EvtVector4R pd1,pd2,pd3,pd4,ps;
169 EvtVector4R dp1,dp2;
170 double xmass2,ymass2;
171
172 id1 =p->getDaug(d1);
173 id2 =p->getDaug(d2);
174 id3 =p->getDaug(d3);
175 id4 =p->getDaug(d4);
176
177 pd1 =id1->getP4Lab();
178 pd2 =id2->getP4Lab();
179 pd3 =id3->getP4Lab();
180 pd4 =id4->getP4Lab();
181
182 dp1 =p->getDaug(nd1)->getP4Lab();
183 dp2 =p->getDaug(nd2)->getP4Lab();
184
185 xmass2=(pd1+pd2).mass2();
186 ymass2=(pd3+pd4).mass2();
187
188 int xbin = hid->GetXaxis()->FindBin(xmass2);
189 int ybin = hid->GetYaxis()->FindBin(ymass2);
190 int xybin= hid->GetBin(xbin,ybin);
191 double zvalue=hid->GetBinContent(xybin);
192 double xratio=zvalue/ymax;
193 if(xratio ==0 ) goto loop;
194 double rd1=EvtRandom::Flat(0.0, 1.0);
195 if(rd1>xratio) goto loop; //sigle out event by 2-D mass distribution histo.
196
197 double da1=dp1.get(3)/dp1.d3mag();
198 double da2=dp2.get(3)/dp2.d3mag();
199
200 int dbin1=did1->FindBin(da1);
201 int dbin2=did2->FindBin(da2);
202
203 double dr1=(did1->GetBinContent(dbin1))/avm1;
204 double dr2=(did2->GetBinContent(dbin2))/avm2;
205 if(dr1==0 || dr2==0) goto loop;
206 rd1=EvtRandom::Flat(0.0, 1.0);
207 if(rd1>dr1) goto loop; //sigle out event by the first angular distribution histo.
208
209 rd1=EvtRandom::Flat(0.0, 1.0);
210 if(rd1>dr2) goto loop; //sigle out event by the second angular distribution histo.
211
212 return ;
213}
const char * setFileName()
Definition: UserBody3.cc:11
int * setDaugPair()
Definition: UserBody3.cc:38
const char * setDaugAng(int i)
Definition: UserBody3.cc:23
int setDaugAngNo()
Definition: UserBody3.cc:33
const char * setHpoint()
Definition: UserBody3.cc:17
EvtId * getDaugs()
Definition: EvtDecayBase.hh:65
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:685
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 get(int i) const
Definition: EvtVector4R.hh:179
double d3mag() const
Definition: EvtVector4R.cc:186
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")

◆ getName()

void EvtBody3::getName ( std::string &  name)
virtual

Implements EvtDecayBase.

Definition at line 59 of file EvtBody3.cc.

59 {
60
61 model_name="Body3";
62
63}

◆ init()

void EvtBody3::init ( )
virtual

Reimplemented from EvtDecayBase.

Definition at line 72 of file EvtBody3.cc.

72 {
73
74 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
75 checkNArg(0);
77}
EvtId getParentId()
Definition: EvtDecayBase.hh:60
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
static EvtSpinType::spintype getSpinType(EvtId i)
Definition: EvtPDL.hh:61

◆ initProbMax()

void EvtBody3::initProbMax ( )
virtual

Reimplemented from EvtDecayBase.

Definition at line 78 of file EvtBody3.cc.

78 {
79
80 noProbMax();
81
82}
void noProbMax()

◆ setDaugAng()

const char * EvtBody3::setDaugAng ( int  i)

Definition at line 23 of file UserBody3.cc.

23 {
24
25const char* DaugAng0="costhetag"; //Daugher numbering defined at DEC file
26const char* DaugAng1="costhetav";
27const char* DaugAng2="xxx"; //not use, then set PN=2 in setDaugAngNo();
28 if(i==0) return DaugAng0;
29 if(i==1) return DaugAng1;
30 if(i==2) return DaugAng2;
31}

Referenced by decay().

◆ setDaugAngNo()

int EvtBody3::setDaugAngNo ( )

Definition at line 33 of file UserBody3.cc.

33 {
34 int PN=2; //e.g. DaugAng2="xxx", set PN=2; DaugAng0="xxx", set PN=0;
35 return PN;
36}

Referenced by decay().

◆ setDaugPair()

int * EvtBody3::setDaugPair ( )

Definition at line 38 of file UserBody3.cc.

38 {
39 static int DP[4];
40 DP[0]=0; // 0,1,2,... indexes for daughter particles
41 DP[1]=2;
42 DP[2]=1;
43 DP[3]=2;
44 return DP;
45}

Referenced by decay().

◆ setFileName()

const char * EvtBody3::setFileName ( )

Definition at line 11 of file UserBody3.cc.

11 {
12 const char* filename;
13 filename="./diy.root"; //specify the root histor. name
14 return filename;
15}

Referenced by decay().

◆ setHpoint()

const char * EvtBody3::setHpoint ( )

Definition at line 17 of file UserBody3.cc.

17 {
18 const char* hpoint;
19 hpoint="hdalitz"; //specify the histor. id
20 return hpoint;
21}

Referenced by decay().


The documentation for this class was generated from the following files: