BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtHAngSam3.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 the invariant mass distribution of histogram
14// dimension one.
15//
16// Modification history:
17//
18// Ping R.-G. December, 2006 Module created
19//
20//------------------------------------------------------------------------
21//
22#include "EvtHAngSam3.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 EvtHAngSam3::getName( std::string& model_name ) { model_name = "HAngSam3"; }
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
76 TFile f( fl );
77 TH1F* hid = (TH1F*)f.Get( hp );
78 TAxis* xaxis = hid->GetXaxis();
79
80 double BLE = xaxis->GetBinLowEdge( 1 );
81 int BINS = xaxis->GetLast();
82 double yvalue[20000];
83 static double ymax = 0.0;
84 int i;
85 double Ntotal = 0, yc[20000];
86 static int icount = 0;
87
88 if ( icount == 0 )
89 {
90 for ( i = 1; i < BINS + 1; i++ )
91 {
92 yvalue[i] = hid->GetBinContent( i );
93 if ( yvalue[i] > ymax ) ymax = yvalue[i];
94 }
95 icount++;
96 ymax = ymax * 1.2;
97 }
98
99loop:
101
102 EvtParticle *d1, *d2, *d3;
103 EvtVector4R pd1, pd2, pd3, pp, nor;
104 double costheta, cosphi;
105
106 d1 = p->getDaug( 0 );
107 d2 = p->getDaug( 1 );
108 pd1 = d1->getP4();
109 pd2 = d2->getP4();
110 pp = p->getP4();
111 double ppmag = pp.d3mag();
112 nor = pd1.cross( pd2 );
113
114 if ( ppmag < 0.0001 ) { costheta = nor.get( 3 ) / nor.d3mag(); }
115 else { costheta = nor.dot( pp ) / nor.d3mag() / pp.d3mag(); }
116
117 int xbin = hid->FindBin( costheta );
118 double xratio = ( hid->GetBinContent( xbin ) ) / ymax;
119
120 double rd1 = EvtRandom::Flat( 0.0, 1.0 );
121 if ( rd1 > xratio ) goto loop;
122 return;
123}
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)
EvtDecayBase * clone()
const char * setFileName()
void decay(EvtParticle *p)
void initProbMax()
virtual ~EvtHAngSam3()
void getName(std::string &name)
const char * setHpoint()
static EvtSpinType::spintype getSpinType(EvtId i)
Definition EvtPDL.hh:66
const EvtVector4R & getP4() 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 dot(const EvtVector4R &v2) const
double get(int i) const
EvtVector4R cross(const EvtVector4R &v2)
double d3mag() const