Demonstrate how batch scripting can be used to define filters


{

    //First define our histograms
    TH1F * histo1     = new TH1F ("histo1","p p missing mass",100,0.1,1.0);
    TH1F * histo1_cut = new TH1F ("histo1_cut","p p missing mass",100,0.1,1.0);

    TH1F * histo2      = new TH1F ("histo2","p p missing momentum",100,0.1,1.0);
    TH1F * histo2_cut  = new TH1F ("histo2_cut","p p missing momentum",100,0.1,1.0);
    TH1F * histo2_cut2 = new TH1F ("histo2_cut2","p p missing momentum",100,0.1,1.0);

    TH1F * histo3 = new TH1F ("histo3","Theta of the first proton",100,0,180);

    //Define the reaction
    PReaction my_reaction("3.5","p","p","p p w [pi+ pi- pi0]");
    
    //The batch commands
    my_reaction.Do("pp_miss = [p + p] - ([p,1] + [p,2]); total_miss = [p + p] - ([p,1] + [p,2] + [pi+] + [pi-])");

    //my_reaction.Do("[p,1]->Print();[p,2]->Print(); echo *************");
    //A simple missing mass/momentum histogram
    my_reaction.Do(histo1,"_x= pp_miss->M2()");
    my_reaction.Do(histo2,"_x= pp_miss->P()");

    //Mathematical operations can be used:
    my_reaction.Do(histo3,"_x= ([p,1]->Theta() * 180.)/TMath::Pi()");

    //*******Conditional histogram
    //Only omegas having a certain momentum:
    my_reaction.Do(histo1_cut,"if pp_miss->P() > 0.3; _x= pp_miss->M2()");
    //Control (to be on the save side!):
    my_reaction.Do(histo2_cut,"if pp_miss->P() > 0.3; _x= pp_miss->P()");
    
    
    //*******Conditional Reaction: filters (in our case connected with a control histo)
    //Boolean operations:
    my_reaction.Do(histo2_cut2,"#mom = 0.; if pp_miss->P() > 0.3 && pp_miss->P() < 0.6 ; #mom = 1.");

    my_reaction.Print();

    //Start event loop:
    cout << my_reaction.Loop(50000) << " events recorded" << endl;

    //Plot our on-line histograms
    //******* Starting from here all commands are standard ROOT ***********
    TCanvas *c1 = new TCanvas("c1","Canvas");
    c1->SetLogy(1);
    c1->Divide(2);
    c1->cd(1);
    histo1->Draw("");
    histo1_cut->Draw("same");
    histo1_cut->SetLineColor(2);

    c1->cd(2);
    histo2->Draw("");
    histo2_cut->Draw("same");
    histo2_cut->SetLineColor(2);
    histo2_cut2->Draw("same");
    histo2_cut2->SetLineColor(3);
}
 batch_filters.C:1
 batch_filters.C:2
 batch_filters.C:3
 batch_filters.C:4
 batch_filters.C:5
 batch_filters.C:6
 batch_filters.C:7
 batch_filters.C:8
 batch_filters.C:9
 batch_filters.C:10
 batch_filters.C:11
 batch_filters.C:12
 batch_filters.C:13
 batch_filters.C:14
 batch_filters.C:15
 batch_filters.C:16
 batch_filters.C:17
 batch_filters.C:18
 batch_filters.C:19
 batch_filters.C:20
 batch_filters.C:21
 batch_filters.C:22
 batch_filters.C:23
 batch_filters.C:24
 batch_filters.C:25
 batch_filters.C:26
 batch_filters.C:27
 batch_filters.C:28
 batch_filters.C:29
 batch_filters.C:30
 batch_filters.C:31
 batch_filters.C:32
 batch_filters.C:33
 batch_filters.C:34
 batch_filters.C:35
 batch_filters.C:36
 batch_filters.C:37
 batch_filters.C:38
 batch_filters.C:39
 batch_filters.C:40
 batch_filters.C:41
 batch_filters.C:42
 batch_filters.C:43
 batch_filters.C:44
 batch_filters.C:45
 batch_filters.C:46
 batch_filters.C:47
 batch_filters.C:48
 batch_filters.C:49
 batch_filters.C:50
 batch_filters.C:51
 batch_filters.C:52
 batch_filters.C:53
 batch_filters.C:54
 batch_filters.C:55
 batch_filters.C:56
 batch_filters.C:57
 batch_filters.C:58
 batch_filters.C:59
 batch_filters.C:60
 batch_filters.C:61