Choosing the rho0 -> e+e- partial width model


{

    TCanvas *c1 = new TCanvas();
    c1->SetLogy(1);

    //Remove VDM M**3 scaling
    //Not part of standard pluto -> Add this by hand
    PFixedDecay *pmodel = new PFixedDecay("rho0_novdm_e-_e+",
					  "Rho decay without VDM", -1);
    pmodel->Add("rho0, parent");
    pmodel->Add("e+, daughter");
    pmodel->Add("e-, daughter");
    makeDistributionManager()->Add(pmodel);

    //Print out what we have:
    makeDistributionManager()->Print("decay_models");

    //Choosing the model using the distribution manager
    makeDistributionManager()->Enable("rho_picutoff_e-_e+");

    //Write the primary models to the data base:
    makeDistributionManager()->LinkDB();

    //Read the primary model for the rho meson:
    PChannelModel *rho = makeDynamicData()->GetParticleModel(makeStaticData()->GetParticleID("rho0"));

    //Print rho0 properties:
    makeStaticData()->PrintParticle("rho0");

    //Set the decay index of the primary model:
    rho->SetParameter(1, 85);
    rho->SetDidx(85);
    
    //Save rho model as ROOT will not evaluate the GetRandom again
    PChannelModel *rho2 = (PChannelModel *)rho->Clone();
    PChannelModel *rho3 = (PChannelModel *)rho->Clone();

    TH1F *h1 = new TH1F("rho0", "rho0", 100, 0., 1.2);
    rho->GetRandom();
    for (int i=0; i<10000; i++) 
	h1->Fill(rho->GetRandom());
    h1->Draw();
    
    //Change it!
    makeDistributionManager()->Enable("rho0_ee_e-_e+");
    makeDistributionManager()->LinkDB();

    TH1F *h2 = new TH1F("rho0_2", "rho0_2", 100, 0., 1.2);
    for (int i=0; i<15000; i++) 
	h2->Fill(rho2->GetRandom());
    h2->SetLineStyle(7);
    h2->Draw("same");

    makeDistributionManager()->Enable("rho0_novdm_e-_e+");
    makeDistributionManager()->LinkDB();

    TH1F *h3 = new TH1F("rho0_3", "rho0_3", 100, 0., 1.2);
    for (int i=0; i<10000; i++) 
	h3->Fill(rho3->GetRandom());
    h3->SetLineStyle(2);
    h3->Draw("same");

    //scale all histograms that they match in the limit of higher masses
    h2->Scale(h1->GetBinContent(90) / h2->GetBinContent(90));
    h3->Scale(h1->GetBinContent(90) / h3->GetBinContent(90));
}
 rho_mass_shape.C:1
 rho_mass_shape.C:2
 rho_mass_shape.C:3
 rho_mass_shape.C:4
 rho_mass_shape.C:5
 rho_mass_shape.C:6
 rho_mass_shape.C:7
 rho_mass_shape.C:8
 rho_mass_shape.C:9
 rho_mass_shape.C:10
 rho_mass_shape.C:11
 rho_mass_shape.C:12
 rho_mass_shape.C:13
 rho_mass_shape.C:14
 rho_mass_shape.C:15
 rho_mass_shape.C:16
 rho_mass_shape.C:17
 rho_mass_shape.C:18
 rho_mass_shape.C:19
 rho_mass_shape.C:20
 rho_mass_shape.C:21
 rho_mass_shape.C:22
 rho_mass_shape.C:23
 rho_mass_shape.C:24
 rho_mass_shape.C:25
 rho_mass_shape.C:26
 rho_mass_shape.C:27
 rho_mass_shape.C:28
 rho_mass_shape.C:29
 rho_mass_shape.C:30
 rho_mass_shape.C:31
 rho_mass_shape.C:32
 rho_mass_shape.C:33
 rho_mass_shape.C:34
 rho_mass_shape.C:35
 rho_mass_shape.C:36
 rho_mass_shape.C:37
 rho_mass_shape.C:38
 rho_mass_shape.C:39
 rho_mass_shape.C:40
 rho_mass_shape.C:41
 rho_mass_shape.C:42
 rho_mass_shape.C:43
 rho_mass_shape.C:44
 rho_mass_shape.C:45
 rho_mass_shape.C:46
 rho_mass_shape.C:47
 rho_mass_shape.C:48
 rho_mass_shape.C:49
 rho_mass_shape.C:50
 rho_mass_shape.C:51
 rho_mass_shape.C:52
 rho_mass_shape.C:53
 rho_mass_shape.C:54
 rho_mass_shape.C:55
 rho_mass_shape.C:56
 rho_mass_shape.C:57
 rho_mass_shape.C:58
 rho_mass_shape.C:59
 rho_mass_shape.C:60
 rho_mass_shape.C:61
 rho_mass_shape.C:62
 rho_mass_shape.C:63
 rho_mass_shape.C:64
 rho_mass_shape.C:65
 rho_mass_shape.C:66
 rho_mass_shape.C:67
 rho_mass_shape.C:68