Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
PhaseSpace.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_physics
3/// \notebook -js
4/// Example of use of TGenPhaseSpace
5///
6/// \macro_image
7/// \macro_code
8///
9/// \author Valerio Filippini
10
11void PhaseSpace() {
12
13 TLorentzVector target(0.0, 0.0, 0.0, 0.938);
14 TLorentzVector beam(0.0, 0.0, .65, .65);
15 TLorentzVector W = beam + target;
16
17 //(Momentum, Energy units are Gev/C, GeV)
18 Double_t masses[3] = { 0.938, 0.139, 0.139} ;
19
21 event.SetDecay(W, 3, masses);
22
23 TH2F *h2 = new TH2F("h2","h2", 50,1.1,1.8, 50,1.1,1.8);
24
25 for (Int_t n=0;n<100000;n++) {
26 Double_t weight = event.Generate();
27
28 TLorentzVector *pProton = event.GetDecay(0);
29
30 TLorentzVector *pPip = event.GetDecay(1);
31 TLorentzVector *pPim = event.GetDecay(2);
32
33 TLorentzVector pPPip = *pProton + *pPip;
34 TLorentzVector pPPim = *pProton + *pPim;
35
36 h2->Fill(pPPip.M2() ,pPPim.M2() ,weight);
37 }
38 h2->Draw();
39}
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
Utility class to generate n-body event, with constant cross-section (default) or with Fermi energy de...
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3074
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:251
Int_t Fill(Double_t)
Invalid Fill method.
Definition TH2.cxx:358
Double_t M2() const
const Int_t n
Definition legend1.C:16