Usage of Pluto build-in templates: beam smearing


{
    //This macro demonstrates how the beam smearing model
    //of Pluto can be used.
    //
    //In addition, one can play around what happens to the 
    //reconstructed eta mass if the unkown beam momentum is used


    PBeamSmearing *smear = new PBeamSmearing("beam_smear", "Beam smearing");
    smear->SetReaction("p + p");

    //Now one has to set the beam parameters
    //There are 3 options:
    //   * standard set (from the old Pluto)
    //     only angular smearing
    //   * A distribution function for angular
    //     smearing (the rho in polar coordinates)
    //   * A momentum smearing function
    //     (stretching, e.g. ==1 means no change)
    // 
    
    //Set "standard" parameters:
    //smear->SetBeamParameters(1., 0. , 0.5);
    //Parameters 
    // 1.) Tilt theta in degree
    // 2.) Tilt phi in degree
    // 3.) Sigma (of gaus) in degree

    //An angular spot
    //Range of theta in unites of degree
    //N.B. the radial distance probability is not included, so
    //if you want to have a flat distribution
    //one has to take this into account by multiplying with x
    //This produces a flat beam spot around around +/- 1deg
    smear->SetAngularSmearing(new TF1("delme","1 *x", 0, 1.));

    //momentum smearing +/- 10%
    //smear->SetMomentumSmearing(new TF1("delme", "1", 0.9, 1.1));

    makeDistributionManager()->Add(smear);

    PReaction my_reaction("2.2", "p", "p", "p p eta [dilepton [e+ e-] g]", "eta_dalitz", 1, 0, 0, 0);

    //This histogram shows the beam profile:
    TH2F *histo1 = new TH2F ("histo1", "Px vs. Py of beam", 100, -.1, .1, 100, -.1, .1);
    my_reaction.Do(histo1, "_x = [p + p]->Px(); _y  = [p + p]->Py();");
    
    //This histogram shows how a wrong assumption about the beam
    //momentum can influence the reconstruction in an exclusive
    //reaction
    TH1F *histo2 = new TH1F ("histo2", "Reconstructed eta mass", 100, 0.3, 0.7);
    my_reaction.Do("wrong_cm = P3E(0.000000,0.000000,2.994728,4.076545);");
    my_reaction.Do(histo2, "_x = (wrong_cm - ([p,1] + [p,2]))->M();");

    my_reaction.Print();   //The "Print()" statement is optional
    my_reaction.Loop(10000);
    
    histo1->Draw();
    
    new TCanvas();

    histo2->Draw();

}
 beamsmearing_example.C:1
 beamsmearing_example.C:2
 beamsmearing_example.C:3
 beamsmearing_example.C:4
 beamsmearing_example.C:5
 beamsmearing_example.C:6
 beamsmearing_example.C:7
 beamsmearing_example.C:8
 beamsmearing_example.C:9
 beamsmearing_example.C:10
 beamsmearing_example.C:11
 beamsmearing_example.C:12
 beamsmearing_example.C:13
 beamsmearing_example.C:14
 beamsmearing_example.C:15
 beamsmearing_example.C:16
 beamsmearing_example.C:17
 beamsmearing_example.C:18
 beamsmearing_example.C:19
 beamsmearing_example.C:20
 beamsmearing_example.C:21
 beamsmearing_example.C:22
 beamsmearing_example.C:23
 beamsmearing_example.C:24
 beamsmearing_example.C:25
 beamsmearing_example.C:26
 beamsmearing_example.C:27
 beamsmearing_example.C:28
 beamsmearing_example.C:29
 beamsmearing_example.C:30
 beamsmearing_example.C:31
 beamsmearing_example.C:32
 beamsmearing_example.C:33
 beamsmearing_example.C:34
 beamsmearing_example.C:35
 beamsmearing_example.C:36
 beamsmearing_example.C:37
 beamsmearing_example.C:38
 beamsmearing_example.C:39
 beamsmearing_example.C:40
 beamsmearing_example.C:41
 beamsmearing_example.C:42
 beamsmearing_example.C:43
 beamsmearing_example.C:44
 beamsmearing_example.C:45
 beamsmearing_example.C:46
 beamsmearing_example.C:47
 beamsmearing_example.C:48
 beamsmearing_example.C:49
 beamsmearing_example.C:50
 beamsmearing_example.C:51
 beamsmearing_example.C:52
 beamsmearing_example.C:53
 beamsmearing_example.C:54
 beamsmearing_example.C:55
 beamsmearing_example.C:56
 beamsmearing_example.C:57
 beamsmearing_example.C:58
 beamsmearing_example.C:59
 beamsmearing_example.C:60
 beamsmearing_example.C:61
 beamsmearing_example.C:62
 beamsmearing_example.C:63
 beamsmearing_example.C:64
 beamsmearing_example.C:65
 beamsmearing_example.C:66