BOSS 7.0.4
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtAngH2.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:
14// This model allows user to specify scatter plot id, x-axis is the angular distribution(in Lab system) for daugher(0)
15//, y-axix is the the angular distribution(in Lab system) for daugher(1). They are corrected by detection efficiency.
16// The scatter plot is filled with cos(theta_0) vs. cos(theta_1),where subscript 0,1 denote the daughter index
17// if son_0 and son_1 are identitical particles, to distinguish them by E_0<E_1
18// Modification history:
19//
20// Ping R.-G. December, 2006 Module created
21//
22//------------------------------------------------------------------------
23//
25#include <stdlib.h>
28#include "EvtGenBase/EvtPDL.hh"
45#include <string>
46
47#include "TH1.h"
48#include "TAxis.h"
49#include "TH2.h"
50#include "TFile.h"
51#include "TApplication.h"
52#include "TROOT.h"
53//#include "CLHEP/config/CLHEP.h"
54//#include "CLHEP/config/TemplateFunctions.h"
55
56
57using std::endl;
58
60
61void EvtAngH2::getName(std::string& model_name){
62
63 model_name="AngH2";
64
65}
66
68
69 return new EvtAngH2;
70
71}
72
73
75
76 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
77 checkNArg(0);
79}
81
82 noProbMax();
83
84}
85
87
88 const char* fl=setFileName();
89 const char* hp=setHpoint();
90 int* dp;
91
92 TFile f(fl);
93 TH2F *hid = (TH2F*)f.Get(hp);
94 TAxis* xaxis = hid->GetXaxis();
95 TAxis* yaxis = hid->GetYaxis();
96
97 int BINSx =xaxis->GetLast();
98 int BINSy =yaxis->GetLast();
99 int BINS =BINSx*BINSy;
100 double yvalue,ymax=0.0;
101 int i,j,binxy;
102
103 for(i=1;i<BINSx+1;i++){
104 for(j=1;j<BINSy+1;j++){
105 binxy=hid->GetBin(i,j);
106 yvalue=hid->GetBinContent(binxy);
107// cout<<"binxy,yvalue = "<<binxy<<"; "<<yvalue<<endl;
108 if(yvalue>ymax) ymax=yvalue;
109 }
110 }
111
112loop:
114
115 EvtParticle *id1,*id2,*id3,*id4,*s1;
116 EvtVector4R pd1,pd2,pd3,pd4,ps;
117 double xcostheta,ycostheta;
118
119 id1 =p->getDaug(0);
120 id2 =p->getDaug(1);
121 id3 =p->getDaug(2);
122
123
124 pd1 =id1->getP4Lab();
125 pd2 =id2->getP4Lab();
126 pd3 =id3->getP4Lab();
127
128 xcostheta=pd1.get(3)/pd1.d3mag();
129 ycostheta=pd2.get(3)/pd2.d3mag();
130 if( EvtPDL::name(p->getDaug(0)->getId()) == EvtPDL::name(p->getDaug(1)->getId()) ){
131 if(pd1.get(0)>pd2.get(0)){
132 xcostheta=pd2.get(3)/pd2.d3mag();
133 ycostheta=pd1.get(3)/pd1.d3mag();
134 }
135 }
136 int xbin = hid->GetXaxis()->FindBin(xcostheta);
137 int ybin = hid->GetYaxis()->FindBin(ycostheta);
138 int xybin= hid->GetBin(xbin,ybin);
139 double zvalue=hid->GetBinContent(xybin);
140 double xratio=zvalue/ymax;
141// cout<<"***************zvalue,ymax,xratio= "<<zvalue<<"; "<<ymax<<"; "<<xratio<<endl;
142 double rd1=EvtRandom::Flat(0.0, 1.0);
143 if(rd1>xratio) goto loop;
144 return ;
145}
146
147
virtual ~EvtAngH2()
Definition: EvtAngH2.cc:59
void init()
Definition: EvtAngH2.cc:74
void decay(EvtParticle *p)
Definition: EvtAngH2.cc:86
const char * setHpoint()
Definition: UserAngH2.cc:16
void getName(std::string &name)
Definition: EvtAngH2.cc:61
EvtAngH2()
Definition: EvtAngH2.hh:34
EvtDecayBase * clone()
Definition: EvtAngH2.cc:67
const char * setFileName()
Definition: UserAngH2.cc:10
void initProbMax()
Definition: EvtAngH2.cc:80
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)
static std::string name(EvtId i)
Definition: EvtPDL.hh:64
static EvtSpinType::spintype getSpinType(EvtId i)
Definition: EvtPDL.hh:61
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:685
EvtId getId() const
Definition: EvtParticle.cc:113
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")