BOSS 8.0.0
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
15// Lab system) for daugher(0)
16//, y-axix is the the angular distribution(in Lab system) for daugher(1). They are corrected
17// by detection efficiency.
18// The scatter plot is filled with cos(theta_0) vs. cos(theta_1),where subscript 0,1 denote
19// the daughter index if son_0 and son_1 are identitical particles, to distinguish them by
20// E_0<E_1
21// Modification history:
22//
23// Ping R.-G. December, 2006 Module created
24//
25//------------------------------------------------------------------------
26//
27#include "EvtAngH2.hh"
47#include <stdlib.h>
48#include <string>
49
50#include "TApplication.h"
51#include "TAxis.h"
52#include "TFile.h"
53#include "TH1.h"
54#include "TH2.h"
55#include "TROOT.h"
56// #include "CLHEP/config/CLHEP.h"
57// #include "CLHEP/config/TemplateFunctions.h"
58
59using std::endl;
60
62
63void EvtAngH2::getName( std::string& model_name ) { model_name = "AngH2"; }
64
66
68
69 // check that there are 4 arguments: Invariant mass part. Index: i,j, histor. file name, Hid
70 checkNArg( 0 );
72}
74
76
77 const char* fl = setFileName();
78 const char* hp = setHpoint();
79 int* dp;
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 xcostheta, ycostheta;
109
110 id1 = p->getDaug( 0 );
111 id2 = p->getDaug( 1 );
112 id3 = p->getDaug( 2 );
113
114 pd1 = id1->getP4Lab();
115 pd2 = id2->getP4Lab();
116 pd3 = id3->getP4Lab();
117
118 xcostheta = pd1.get( 3 ) / pd1.d3mag();
119 ycostheta = pd2.get( 3 ) / pd2.d3mag();
120 if ( EvtPDL::name( p->getDaug( 0 )->getId() ) == EvtPDL::name( p->getDaug( 1 )->getId() ) )
121 {
122 if ( pd1.get( 0 ) > pd2.get( 0 ) )
123 {
124 xcostheta = pd2.get( 3 ) / pd2.d3mag();
125 ycostheta = pd1.get( 3 ) / pd1.d3mag();
126 }
127 }
128 int xbin = hid->GetXaxis()->FindBin( xcostheta );
129 int ybin = hid->GetYaxis()->FindBin( ycostheta );
130 int xybin = hid->GetBin( xbin, ybin );
131 double zvalue = hid->GetBinContent( xybin );
132 double xratio = zvalue / ymax;
133 // cout<<"***************zvalue,ymax,xratio= "<<zvalue<<"; "<<ymax<<"; "<<xratio<<endl;
134 double rd1 = EvtRandom::Flat( 0.0, 1.0 );
135 if ( rd1 > xratio ) goto loop;
136 return;
137}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
virtual ~EvtAngH2()
Definition EvtAngH2.cc:61
void init()
Definition EvtAngH2.cc:67
void decay(EvtParticle *p)
Definition EvtAngH2.cc:75
const char * setHpoint()
Definition UserAngH2.cc:16
void getName(std::string &name)
Definition EvtAngH2.cc:63
EvtDecayBase * clone()
Definition EvtAngH2.cc:65
const char * setFileName()
Definition UserAngH2.cc:10
void initProbMax()
Definition EvtAngH2.cc:73
EvtId getParentId()
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
static std::string name(EvtId i)
Definition EvtPDL.hh:70
static EvtSpinType::spintype getSpinType(EvtId i)
Definition EvtPDL.hh:66
EvtVector4R getP4Lab()
EvtId getId() const
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