BOSS 8.0.0
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//
22#include "EvtMassH2.hh"
42#include <stdlib.h>
43#include <string>
44
45#include "TApplication.h"
46#include "TAxis.h"
47#include "TFile.h"
48#include "TH1.h"
49#include "TH2.h"
50#include "TROOT.h"
51// #include "CLHEP/config/CLHEP.h"
52// #include "CLHEP/config/TemplateFunctions.h"
53
54using std::endl;
55
57
58void EvtMassH2::getName( std::string& model_name ) { model_name = "MassH2"; }
59
61
63
64 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
65 checkNArg( 0 );
67}
69
71
72 const char* fl = setFileName();
73 const char* hp = setHpoint();
74 int* dp;
75 dp = setDaugPair();
76 int d1 = dp[0]; // m(d1,d2) pair at X axis
77 int d2 = dp[1];
78 int d3 = dp[2]; // m(d3,d4) pair at Y axis
79 int d4 = dp[3];
80
81 TFile f( fl );
82 TH2F* hid = (TH2F*)f.Get( hp );
83 TAxis* xaxis = hid->GetXaxis();
84 TAxis* yaxis = hid->GetYaxis();
85
86 int BINSx = xaxis->GetLast();
87 int BINSy = yaxis->GetLast();
88 int BINS = BINSx * BINSy;
89 double yvalue, ymax = 0.0;
90 int i, j, binxy;
91
92 for ( i = 1; i < BINSx + 1; i++ )
93 {
94 for ( j = 1; j < BINSy + 1; j++ )
95 {
96 binxy = hid->GetBin( i, j );
97 yvalue = hid->GetBinContent( binxy );
98 // cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
99 if ( yvalue > ymax ) ymax = yvalue;
100 }
101 }
102
103loop:
105
106 EvtParticle *id1, *id2, *id3, *id4, *s1;
107 EvtVector4R pd1, pd2, pd3, pd4, ps;
108 double xmass2, ymass2;
109
110 id1 = p->getDaug( d1 );
111 id2 = p->getDaug( d2 );
112 id3 = p->getDaug( d3 );
113 id4 = p->getDaug( d4 );
114
115 pd1 = id1->getP4Lab();
116 pd2 = id2->getP4Lab();
117 pd3 = id3->getP4Lab();
118 pd4 = id4->getP4Lab();
119
120 xmass2 = ( pd1 + pd2 ).mass2();
121 ymass2 = ( pd3 + pd4 ).mass2();
122
123 int xbin = hid->GetXaxis()->FindBin( xmass2 );
124 int ybin = hid->GetYaxis()->FindBin( ymass2 );
125 int xybin = hid->GetBin( xbin, ybin );
126 double zvalue = hid->GetBinContent( xybin );
127 double xratio = zvalue / ymax;
128 // cout<<"***************zvalue,ymax,xratio= "<<zvalue<<"; "<<ymax<<"; "<<xratio<<endl;
129 double rd1 = EvtRandom::Flat( 0.0, 1.0 );
130 if ( rd1 > xratio ) goto loop;
131 return;
132}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
EvtId getParentId()
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
virtual ~EvtMassH2()
Definition EvtMassH2.cc:56
void initProbMax()
Definition EvtMassH2.cc:68
EvtDecayBase * clone()
Definition EvtMassH2.cc:60
void init()
Definition EvtMassH2.cc:62
const char * setHpoint()
Definition UserMassH2.cc:16
int * setDaugPair()
Definition UserMassH2.cc:22
void decay(EvtParticle *p)
Definition EvtMassH2.cc:70
void getName(std::string &name)
Definition EvtMassH2.cc:58
const char * setFileName()
Definition UserMassH2.cc:10
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