Total cross sections: The dp cocktail at 1.25AGeV

//
//
{

//for explanations what is going on here see pp_cocktail_1.25.C
    
    makeDistributionManager()->Exec("elementary");
    makeDistributionManager()->Exec("brems: kaptari; weighting");

//makeDistributionManager()->Exec("brems: delta");
    makeDistributionManager()->Exec("brems: sum");  
//makeDistributionManager()->Exec("brems: elastic");
    
//makeDistributionManager()->Exec("brems: fsi");  
    
    PTCrossWeight * my_cross = 
	new PTCrossWeight("n + p_to_n_p_pi0_pi0/tcross",
			  "Cross section for double pi0 production",-1);
    my_cross->SetCrossSection(0.1 * 0.001 * 2); //0.1mb times 2 (see below) 
    makeDistributionManager()->Add(my_cross);
    
    makeDistributionManager()->Exec("dalitz_mod: krivoruchenko");
    makeDistributionManager()->Exec("dalitz_mod: static_br_thresh=0.100 ; flat_generator");
    
//Once we are done, print the collection just to see if the plugins worked:    
    
    makeDistributionManager()->Print("root");

//Here we have to add all isospin channels

    PReaction *my_reaction = 
	new PReaction("2.50","d","p","n D+ [p pi0 [g dilepton [e+ e-]]] (p)","dp_1.25",1,0,0,0);
    my_reaction->AddReaction("p D0 [n pi0 [g dilepton [e+ e-]]] (p)");
    my_reaction->AddReaction("p n pi0 pi0 [g dilepton [e+ e-]] (p)");

    //Add Delta, if you want that:
    //my_reaction->AddReaction("n D+ [p dilepton [e+ e-]] (p)");
    //my_reaction->AddReaction("p D0 [n dilepton [e+ e-]] (p)");
    
    //Add Bremsstrahlung, if you want that:
    my_reaction->AddReaction("n p dilepton [e+ e-] (p)");
    //keep in mind double counting, if you enabled brems with "delta" or "sum"

//Attach a control histogram
    TH1F * pn_sum = 
	new TH1F ("pn_sum","pn DiLepton mass",100,0.,0.6);
    pn_sum->Sumw2();

    my_reaction->Do(pn_sum,"_x=[dilepton]->M()");

    my_reaction->Print();

    my_reaction->Preheating(100);
    my_reaction->Loop(10000);

    //Subtreshold eta production
    //make an extra reaction to have an independent histo
    PReaction *my_reaction2 = 
	new PReaction("2.50","d","p","n p eta [g dilepton [e+ e-]] (p)",NULL,1,0,0,0);
    my_reaction2->AddReaction("d eta [g dilepton [e+ e-]] (p)");

    my_reaction2->Print();

//Only eta
    TH1F * eta_sum = 
	new TH1F ("eta_sum","pn DiLepton mass (coherent sum)",100,0.,0.6);

    eta_sum->Sumw2();

    my_reaction2->Do(eta_sum,"_x=[dilepton]->M()");

    my_reaction2->Preheating(100);

    my_reaction2->Loop(10000);

    
//It is very crucial to correct the histo for binsize
    pn_sum->Add(eta_sum);
    PUtils::correct(pn_sum);
    PUtils::correct(eta_sum);


    TCanvas *c1 = new TCanvas("ee_mass", "ee invmass",800,800);
    c1->SetLogy(1);
    pn_sum->Draw();
    eta_sum->Draw("same");

}
 dp_cocktail_1.25.C:1
 dp_cocktail_1.25.C:2
 dp_cocktail_1.25.C:3
 dp_cocktail_1.25.C:4
 dp_cocktail_1.25.C:5
 dp_cocktail_1.25.C:6
 dp_cocktail_1.25.C:7
 dp_cocktail_1.25.C:8
 dp_cocktail_1.25.C:9
 dp_cocktail_1.25.C:10
 dp_cocktail_1.25.C:11
 dp_cocktail_1.25.C:12
 dp_cocktail_1.25.C:13
 dp_cocktail_1.25.C:14
 dp_cocktail_1.25.C:15
 dp_cocktail_1.25.C:16
 dp_cocktail_1.25.C:17
 dp_cocktail_1.25.C:18
 dp_cocktail_1.25.C:19
 dp_cocktail_1.25.C:20
 dp_cocktail_1.25.C:21
 dp_cocktail_1.25.C:22
 dp_cocktail_1.25.C:23
 dp_cocktail_1.25.C:24
 dp_cocktail_1.25.C:25
 dp_cocktail_1.25.C:26
 dp_cocktail_1.25.C:27
 dp_cocktail_1.25.C:28
 dp_cocktail_1.25.C:29
 dp_cocktail_1.25.C:30
 dp_cocktail_1.25.C:31
 dp_cocktail_1.25.C:32
 dp_cocktail_1.25.C:33
 dp_cocktail_1.25.C:34
 dp_cocktail_1.25.C:35
 dp_cocktail_1.25.C:36
 dp_cocktail_1.25.C:37
 dp_cocktail_1.25.C:38
 dp_cocktail_1.25.C:39
 dp_cocktail_1.25.C:40
 dp_cocktail_1.25.C:41
 dp_cocktail_1.25.C:42
 dp_cocktail_1.25.C:43
 dp_cocktail_1.25.C:44
 dp_cocktail_1.25.C:45
 dp_cocktail_1.25.C:46
 dp_cocktail_1.25.C:47
 dp_cocktail_1.25.C:48
 dp_cocktail_1.25.C:49
 dp_cocktail_1.25.C:50
 dp_cocktail_1.25.C:51
 dp_cocktail_1.25.C:52
 dp_cocktail_1.25.C:53
 dp_cocktail_1.25.C:54
 dp_cocktail_1.25.C:55
 dp_cocktail_1.25.C:56
 dp_cocktail_1.25.C:57
 dp_cocktail_1.25.C:58
 dp_cocktail_1.25.C:59
 dp_cocktail_1.25.C:60
 dp_cocktail_1.25.C:61
 dp_cocktail_1.25.C:62
 dp_cocktail_1.25.C:63
 dp_cocktail_1.25.C:64
 dp_cocktail_1.25.C:65
 dp_cocktail_1.25.C:66
 dp_cocktail_1.25.C:67
 dp_cocktail_1.25.C:68
 dp_cocktail_1.25.C:69
 dp_cocktail_1.25.C:70
 dp_cocktail_1.25.C:71
 dp_cocktail_1.25.C:72
 dp_cocktail_1.25.C:73
 dp_cocktail_1.25.C:74
 dp_cocktail_1.25.C:75
 dp_cocktail_1.25.C:76
 dp_cocktail_1.25.C:77
 dp_cocktail_1.25.C:78
 dp_cocktail_1.25.C:79
 dp_cocktail_1.25.C:80
 dp_cocktail_1.25.C:81
 dp_cocktail_1.25.C:82
 dp_cocktail_1.25.C:83
 dp_cocktail_1.25.C:84
 dp_cocktail_1.25.C:85
 dp_cocktail_1.25.C:86
 dp_cocktail_1.25.C:87
 dp_cocktail_1.25.C:88
 dp_cocktail_1.25.C:89