BOSS 8.0.0
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 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 29 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 58 of file EvtBody3.cc.

58{}

Member Function Documentation

◆ clone()

EvtDecayBase * EvtBody3::clone ( )
virtual

Implements EvtDecayBase.

Definition at line 62 of file EvtBody3.cc.

62{ return new EvtBody3; }

◆ decay()

void EvtBody3::decay ( EvtParticle * p)
virtual

Implements EvtDecayBase.

Definition at line 72 of file EvtBody3.cc.

72 {
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")
const char * setFileName()
Definition UserBody3.cc:11
int * setDaugPair()
Definition UserBody3.cc:42
const char * setDaugAng(int i)
Definition UserBody3.cc:23
int setDaugAngNo()
Definition UserBody3.cc:37
const char * setHpoint()
Definition UserBody3.cc:17
EvtId * getDaugs()
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

◆ getName()

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

Implements EvtDecayBase.

Definition at line 60 of file EvtBody3.cc.

60{ model_name = "Body3"; }

◆ init()

void EvtBody3::init ( )
virtual

Reimplemented from EvtDecayBase.

Definition at line 64 of file EvtBody3.cc.

64 {
65
66 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
67 checkNArg( 0 );
69}
EvtId getParentId()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
static EvtSpinType::spintype getSpinType(EvtId i)
Definition EvtPDL.hh:66

◆ initProbMax()

void EvtBody3::initProbMax ( )
virtual

Reimplemented from EvtDecayBase.

Definition at line 70 of file EvtBody3.cc.

70{ noProbMax(); }

◆ setDaugAng()

const char * EvtBody3::setDaugAng ( int i)

Definition at line 23 of file UserBody3.cc.

23 {
24
25 const char* DaugAng0 = "costhetag"; // Daugher numbering defined at DEC file
26 const char* DaugAng1 = "costhetav";
27 const 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
32 cerr << __FILE__ << ":" << __LINE__ << ": "
33 << "Should not reach here!" << endl;
34 abort();
35}

Referenced by decay().

◆ setDaugAngNo()

int EvtBody3::setDaugAngNo ( )

Definition at line 37 of file UserBody3.cc.

37 {
38 int PN = 2; // e.g. DaugAng2="xxx", set PN=2; DaugAng0="xxx", set PN=0;
39 return PN;
40}

Referenced by decay().

◆ setDaugPair()

int * EvtBody3::setDaugPair ( )

Definition at line 42 of file UserBody3.cc.

42 {
43 static int DP[4];
44 DP[0] = 0; // 0,1,2,... indexes for daughter particles
45 DP[1] = 2;
46 DP[2] = 1;
47 DP[3] = 2;
48 return DP;
49}

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: