1 //-------------------------------------------------------------------------
2 // Implementation of Class AliCollisionNormalization
4 // This class is used to store the vertex ditributions in the data
5 // and in Monte Carlo, needed to compute the real number of
6 // collisions a given sample is corresponding to.
7 // The strategy matches what described in CERN-THESIS-2009-033 p 119.
9 // Author: Michele Floris, CERN
10 //-------------------------------------------------------------------------
12 #include "AliCollisionNormalization.h"
13 #include "AliPhysicsSelection.h"
17 #include "AliGenPythiaEventHeader.h"
18 #include "AliGenDPMjetEventHeader.h"
19 #include "AliGenEventHeader.h"
20 #include "AliMCEvent.h"
22 ClassImp(AliCollisionNormalization)
24 const char * AliCollisionNormalization::fProcLabel[] = {"SD","DD","ND", "Unknown"};
26 AliCollisionNormalization::AliCollisionNormalization() :
48 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
49 fHistVzMCGen[iproc] = 0;
50 fHistVzMCRec[iproc] = 0;
51 fHistVzMCTrg[iproc] = 0;
58 AliCollisionNormalization::AliCollisionNormalization(Int_t nbinz, Float_t minz, Float_t maxz):
75 // ctor, allows setting binning
80 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
81 fHistVzMCGen[iproc] = 0;
82 fHistVzMCRec[iproc] = 0;
83 fHistVzMCTrg[iproc] = 0;
89 AliCollisionNormalization::AliCollisionNormalization(const char * dataFile, const char * dataListName,
90 const char * mcFile, const char * mcListName,
91 const char * eventStatFile) :
108 // ctor, loads histograms from file
109 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
110 fHistVzMCGen[iproc] = 0;
111 fHistVzMCRec[iproc] = 0;
112 fHistVzMCTrg[iproc] = 0;
116 TFile * fdata = new TFile (dataFile);
117 TFile * fmc = new TFile (mcFile );
118 TFile * fstat = new TFile(eventStatFile);
120 if (!fdata->IsOpen() || !fmc->IsOpen() || !fstat->IsOpen()) {
121 AliFatal("Cannot open input file(s)");
124 TList * ldata = (TList*) fdata->Get(dataListName);
125 TList * lmc = (TList*) fmc ->Get(mcListName );
127 AliCollisionNormalization * cndata = (AliCollisionNormalization*) ldata->FindObject("AliCollisionNormalization");
128 AliCollisionNormalization * cnmc = (AliCollisionNormalization*) lmc ->FindObject("AliCollisionNormalization");
131 // Assign or book all histos
132 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
133 fHistVzMCGen[iproc]= cnmc->GetVzMCGen(iproc);
134 fHistVzMCRec[iproc]= cnmc->GetVzMCRec(iproc);
135 fHistVzMCTrg[iproc]= cnmc->GetVzMCTrg(iproc);
137 fHistVzData = cndata->GetVzData();
138 fHistProcTypes = cnmc->GetHistProcTypes();
140 fHistStatBin0 = (TH1F*) fstat->Get("fHistStatistics_Bin0");
141 fHistStat = (TH1F*) fstat->Get("fHistStatistics");
146 AliCollisionNormalization::~AliCollisionNormalization(){
149 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
150 if(fHistVzMCGen[iproc]) { delete fHistVzMCGen[iproc] ; fHistVzMCGen[iproc] =0;}
151 if(fHistVzMCRec[iproc]) { delete fHistVzMCRec[iproc] ; fHistVzMCRec[iproc] =0;}
152 if(fHistVzMCTrg[iproc]) { delete fHistVzMCTrg[iproc] ; fHistVzMCTrg[iproc] =0;}
155 if(fHistVzData ) { delete fHistVzData ; fHistVzData =0;}
156 if(fHistStatBin0 ) { delete fHistStatBin0 ; fHistStatBin0 =0;}
157 if(fHistStat ) { delete fHistStat ; fHistStat =0;}
158 if(fHistProcTypes ) { delete fHistProcTypes ; fHistProcTypes =0;}
162 void AliCollisionNormalization::BookAllHistos(){
164 // Book histos of vz distributions vs multiplicity
165 // if vzOnly == kTRUE, it returns a 1d histo with vz dist only
167 // Do not attach histos to the current directory
168 Bool_t oldStatus = TH1::AddDirectoryStatus();
169 TH1::AddDirectory(kFALSE);
171 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
172 fHistVzMCGen [iproc] = (TH2F*) BookVzHisto(TString("fHistVzMCGen")+ fProcLabel[iproc] ,"Vz distribution of generated events vs rec multiplicity ");
173 fHistVzMCRec [iproc] = (TH2F*) BookVzHisto(TString("fHistVzMCRec")+ fProcLabel[iproc] ,"Vz distribution of reconstructed events vs rec multiplicity");
174 fHistVzMCTrg [iproc] = (TH2F*) BookVzHisto(TString("fHistVzMCTrg")+ fProcLabel[iproc] ,"Vz distribution of triggered events vs rec multiplicity ");
176 fHistVzData = (TH2F*) BookVzHisto("fHistVzData" ,"Vz distribution of triggered events vs rec multiplicity ");
177 fHistProcTypes = new TH1F ("fHistProcTypes", "Number of events in the different process classes", kNProcs, -0.5 , kNProcs-0.5);
179 fHistProcTypes->GetXaxis()->SetBinLabel(kProcSD+1,"SD");
180 fHistProcTypes->GetXaxis()->SetBinLabel(kProcND+1,"ND");
181 fHistProcTypes->GetXaxis()->SetBinLabel(kProcDD+1,"DD");
182 fHistProcTypes->GetXaxis()->SetBinLabel(kProcUnknown+1,"Unknown");
184 TH1::AddDirectory(oldStatus);
188 TH1 * AliCollisionNormalization::BookVzHisto(const char * name , const char * title, Bool_t vzOnly) {
190 // Book histos of vz distributions vs multiplicity
191 // if vzOnly == kTRUE, it returns a 1d histo with vz dist only
194 // Do not attach histos to the current directory
195 Bool_t oldStatus = TH1::AddDirectoryStatus();
196 TH1::AddDirectory(kFALSE);
199 Double_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-7,-5.5,-4,-3,-2,-1,0,1,2,3,4,5.5,7,10,15,20,25,30};
201 h = new TH1F(name,title,22,binLimitsVtx);
203 h = new TH2F(name,title,22,binLimitsVtx,100,Double_t(-0.5),Double_t(99.5));
206 h->SetXTitle("V_{z} (cm)");
207 h->SetYTitle("n_{trk}");
210 TH1::AddDirectory(oldStatus);
216 TH2F * AliCollisionNormalization::GetVzMCGen (Int_t procType) {
218 // returns MC gen histo. If proc type < 0 sums over all proc types, reweighting the XS
220 if(procType>=0) return fHistVzMCGen[procType] ;
222 TH2F * sum = (TH2F*) fHistVzMCGen[kProcSD]->Clone();
225 for(Int_t iproc = 0; iproc < kProcUnknown; iproc++){
226 sum->Add(fHistVzMCGen[iproc],GetProcessWeight(iproc));
231 TH2F * AliCollisionNormalization::GetVzMCRec (Int_t procType) {
233 // returns MC rec histo. If proc type < 0 sums over all proc types, reweighting the XS
235 if(procType>=0) return fHistVzMCRec[procType] ;
237 TH2F * sum = (TH2F*) fHistVzMCRec[kProcSD]->Clone();
240 for(Int_t iproc = 0; iproc < kProcUnknown; iproc++){
241 sum->Add(fHistVzMCRec[iproc],GetProcessWeight(iproc));
249 TH2F * AliCollisionNormalization::GetVzMCTrg (Int_t procType) {
251 // returns MC trg histo. If proc type < 0 sums over all proc types, reweighting the XS
253 if(procType>=0) return fHistVzMCTrg[procType] ;
255 TH2F * sum = (TH2F*) fHistVzMCTrg[kProcSD]->Clone();
258 for(Int_t iproc = 0; iproc < kProcUnknown; iproc++){
259 sum->Add(fHistVzMCTrg[iproc],GetProcessWeight(iproc));
266 Double_t AliCollisionNormalization::ComputeNint() {
268 // Compute the number of collisions applying all corrections
269 // TODO: check error propagation
271 TH1 * hVzData = fHistVzData->ProjectionX("_px",2,-1,"E"); // Skip zero bin
272 Int_t allEventsWithVertex = (Int_t) fHistVzData->Integral(0, fHistVzData->GetNbinsX()+1,
273 2, fHistVzData->GetNbinsY()+1); // include under/overflow!, skip 0 bin
275 // Assign histos reweighted with measured XS
276 TH2F * histVzMCGen = GetVzMCGen(-1);
277 TH2F * histVzMCTrg = GetVzMCTrg(-1);
278 TH2F * histVzMCRec = GetVzMCRec(-1);
280 // Before we start: print number of input events to the physics
281 // selection: this allows to crosscheck that all runs were
282 // successfully processed (useful if running on grid: you may have a
283 // crash without noticing it).
284 AliInfo(Form("Input Events (No cuts: %d, After Phys. Sel.:%d)",
285 Int_t(fHistStat->GetBinContent(1,1)),
286 Int_t(fHistStat->GetBinContent(fHistStat->GetNbinsX(),1)))); // Fixme: will this change with a different trigger scheme?
288 // Get or compute BG. This assumes the CINT1B suite
289 Double_t triggeredEventsWith0MultWithBG = fHistVzData->Integral(0, fHistVzData->GetNbinsX()+1,1, 1);
290 Double_t bg = 0; // This will include beam gas + accidentals
291 if (fHistStatBin0->GetNbinsY() > 4) { // FIXME: we need a better criterion to decide...
292 AliInfo("Using BG computed by Physics Selection");
294 // CHECK IF THE COMPUTATION OF BG OFFSET IS STILL OK, IN CASE OF CHANGES TO THE PHYSICS SELECTION
296 Int_t nbiny = fHistStatBin0->GetNbinsY();
297 for(Int_t ibiny =1; ibiny <= nbiny; ibiny++){
299 printf("%d, %s\n", bgOffset, fHistStatBin0->GetYaxis()->GetBinLabel(ibiny) );
301 if(!strncmp("All B", fHistStatBin0->GetYaxis()->GetBinLabel(ibiny),5)) break;
302 if((ibiny+1)==nbiny) AliFatal("Cannot compute bg offset");
305 if(fVerbose > 2) AliInfo(Form("BG Offset: %d",bgOffset));
308 bg = fHistStatBin0->GetBinContent(fHistStatBin0->GetNbinsX(), bgOffset+AliPhysicsSelection::kStatRowBG);
309 bg += fHistStatBin0->GetBinContent(fHistStatBin0->GetNbinsX(),bgOffset+AliPhysicsSelection::kStatRowAcc);
310 Int_t cint1B = (Int_t) fHistStatBin0->GetBinContent(fHistStatBin0->GetNbinsX(),1);
311 if (cint1B != Int_t(triggeredEventsWith0MultWithBG)) {
312 AliWarning(Form("Events in bin0 from physics selection and local counter not consistent: %d - %d", cint1B, Int_t(triggeredEventsWith0MultWithBG)));
315 AliInfo("Computing BG using CINT1A/B/C/E, ignoring intensities");
316 AliError("This will only work for early runs!!! USE BG FROM PHYSICS SELECTION INSTEAD");
317 Int_t icol = fHistStatBin0->GetNbinsX();
318 Int_t cint1B = (Int_t) fHistStatBin0->GetBinContent(icol,1);
319 Int_t cint1A = (Int_t) fHistStatBin0->GetBinContent(icol,2);
320 Int_t cint1C = (Int_t) fHistStatBin0->GetBinContent(icol,3);
321 Int_t cint1E = (Int_t) fHistStatBin0->GetBinContent(icol,4);
322 bg = cint1A + cint1C-2*cint1E ; // FIXME: to be changed to take into account ratios of events
323 if (cint1B != triggeredEventsWith0MultWithBG) {
324 AliWarning(Form("Events in bin0 from physics selection and local counter not consistent: %d - %d", cint1B, (Int_t)triggeredEventsWith0MultWithBG));
328 Double_t triggeredEventsWith0Mult = triggeredEventsWith0MultWithBG - bg;
329 if(fVerbose > 0) AliInfo(Form("Measured events with vertex: %d",allEventsWithVertex));
330 Double_t bin0 = fHistVzData->Integral(0, fHistVzData->GetNbinsX()+1,
332 if(fVerbose > 0) AliInfo(Form("Zero Bin, Meas: %2.2f, BG: %2.2f, Meas - BG: %2.2f", bin0, bg, triggeredEventsWith0Mult));
334 // This pointers are here so that I use the same names used in Jan Fiete's code (allows for faster/easier comparison)
336 TH2* eTrig = histVzMCTrg;
337 TH2* eTrigVtx = histVzMCRec;
338 TH1* eTrigVtx_projx = eTrigVtx->ProjectionX("eTrigVtx_projx", 2, eTrigVtx->GetNbinsX()+1);
340 // compute trigger and vertex efficiency
341 TH2 * effVtxTrig = (TH2*) histVzMCRec->Clone("effVtxTrig");
343 effVtxTrig->Divide(histVzMCRec,histVzMCGen,1,1,"B");
344 // Apply correction to data to get corrected events
345 TH2 * correctedEvents = (TH2*) fHistVzData->Clone("correctedEvents");
346 correctedEvents->Divide(effVtxTrig);
348 // TH1 * correctedEvents = fHistVzData->ProjectionX("eTrigVtx_projx", 2, eTrigVtx->GetNbinsX()+1);
350 // loop over vertex bins
351 for (Int_t i = 1; i <= fHistVzData->GetNbinsX(); i++)
353 Double_t alpha = (Double_t) hVzData->GetBinContent(i) / allEventsWithVertex;
354 Double_t events = alpha * triggeredEventsWith0Mult;
356 if (histVzMCRec->GetBinContent(i,1) == 0)
359 Double_t fZ = eTrigVtx_projx->Integral(0, eTrigVtx_projx->GetNbinsX()+1) / eTrigVtx_projx->GetBinContent(i) *
360 eTrig->GetBinContent(i, 1) / eTrig->Integral(0, eTrig->GetNbinsX()+1, 1, 1);
364 // multiply with trigger correction
365 Double_t correction = histVzMCGen->GetBinContent(i,1)/histVzMCTrg->GetBinContent(i,1);
366 events *= correction;
368 if (fVerbose > 1) Printf(" Bin %d, alpha is %.2f%%, fZ is %.3f, number of events with 0 mult.: %.2f (MC comparison: %.2f)", i, alpha * 100., fZ, events,
369 histVzMCRec->GetBinContent(i,1));
370 correctedEvents->SetBinContent(i, 1, events);
375 // Integrate correctedEvents over full z range
376 Double_t allEvents = correctedEvents->Integral(0, correctedEvents->GetNbinsX()+1,0, correctedEvents->GetNbinsY()+1);
377 // Integrate correctedEvents over needed z range
378 Double_t allEventsZrange = correctedEvents->Integral(correctedEvents->GetXaxis()->FindBin(-fZRange), correctedEvents->GetXaxis()->FindBin(fZRange),
379 0, correctedEvents->GetNbinsY()+1);
381 if(fVerbose > 1) AliInfo(Form("Results in |Vz| < %3.3f",fZRange));
382 if(fVerbose > 1) AliInfo(Form(" Events in Bin0: %2.2f, With > 1 track: %2.2f, All corrected: %2.2f",
383 correctedEvents->Integral(correctedEvents->GetXaxis()->FindBin(-fZRange), correctedEvents->GetXaxis()->FindBin(fZRange),1,1),
384 correctedEvents->Integral(correctedEvents->GetXaxis()->FindBin(-fZRange), correctedEvents->GetXaxis()->FindBin(fZRange),
385 2,correctedEvents->GetNbinsX()+1),
388 Int_t nbin = histVzMCRec->GetNbinsX();
389 AliInfo(Form("Efficiency in the zero bin: %3.3f", histVzMCRec->Integral(0,nbin+1,1,1)/histVzMCGen->Integral(0,nbin+1,1,1) ));
393 AliInfo(Form("Number of collisions in full phase space: %2.2f", allEvents));
394 // effVtxTrig->Draw();
396 // correctedEvents->Draw(); // FIXME: debug
401 Long64_t AliCollisionNormalization::Merge(TCollection* list)
403 // Merge a list of AliCollisionNormalization objects with this
404 // (needed for PROOF).
405 // Returns the number of merged objects (including this).
413 TIterator* iter = list->MakeIterator();
416 // collections of all histograms
417 const Int_t nHists = kNProcs*3+5;
418 TList collections[nHists];
421 while ((obj = iter->Next())) {
423 AliCollisionNormalization* entry = dynamic_cast<AliCollisionNormalization*> (obj);
428 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
429 if (entry->fHistVzMCGen[iproc] ) collections[++ihist].Add(entry->fHistVzMCGen[iproc] );
430 if (entry->fHistVzMCRec[iproc] ) collections[++ihist].Add(entry->fHistVzMCRec[iproc] );
431 if (entry->fHistVzMCTrg[iproc] ) collections[++ihist].Add(entry->fHistVzMCTrg[iproc] );
433 if (entry->fHistVzData ) collections[++ihist].Add(entry->fHistVzData );
434 if (entry->fHistProcTypes ) collections[++ihist].Add(entry->fHistProcTypes );
435 if (entry->fHistStatBin0 ) collections[++ihist].Add(entry->fHistStatBin0 );
436 if (entry->fHistStat ) collections[++ihist].Add(entry->fHistStat );
442 for(Int_t iproc = 0; iproc < kNProcs; iproc++){
443 if (fHistVzMCGen[iproc] ) fHistVzMCGen[iproc] ->Merge(&collections[++ihist]);
444 if (fHistVzMCRec[iproc] ) fHistVzMCRec[iproc] ->Merge(&collections[++ihist]);
445 if (fHistVzMCTrg[iproc] ) fHistVzMCTrg[iproc] ->Merge(&collections[++ihist]);
448 if (fHistVzData ) fHistVzData ->Merge(&collections[++ihist]);
449 if (fHistProcTypes ) fHistProcTypes ->Merge(&collections[++ihist]);
450 if (fHistStatBin0 ) fHistStatBin0 ->Merge(&collections[++ihist]);
451 if (fHistStat ) fHistStat ->Merge(&collections[++ihist]);
459 void AliCollisionNormalization::FillVzMCGen(Float_t vz, Int_t ntrk, AliMCEvent * mcEvt) {
461 // Fill MC gen histo and the process types statistics
463 Int_t evtype = GetProcessType(mcEvt);
464 // When I fill gen histos, I also fill statistics of process types (used to detemine ratios afterwards).
465 fHistProcTypes->Fill(evtype);
466 fHistVzMCGen[evtype]->Fill(vz,ntrk);
469 void AliCollisionNormalization::FillVzMCRec(Float_t vz, Int_t ntrk, AliMCEvent * mcEvt){
472 Int_t evtype = GetProcessType(mcEvt);
473 fHistVzMCRec[evtype]->Fill(vz,ntrk);
476 void AliCollisionNormalization::FillVzMCTrg(Float_t vz, Int_t ntrk, AliMCEvent * mcEvt) {
479 Int_t evtype = GetProcessType(mcEvt);
480 fHistVzMCTrg[evtype]->Fill(vz,ntrk);
484 Int_t AliCollisionNormalization::GetProcessType(AliMCEvent * mcEvt) {
486 // Determine if the event was generated with pythia or phojet and return the process type
488 // Check if mcEvt is fine
490 AliFatal("NULL mc event");
493 // Determine if it was a pythia or phojet header, and return the correct process type
494 AliGenPythiaEventHeader * headPy = 0;
495 AliGenDPMjetEventHeader * headPho = 0;
496 AliGenEventHeader * htmp = mcEvt->GenEventHeader();
498 AliFatal("Cannot Get MC Header!!");
500 if( TString(htmp->IsA()->GetName()) == "AliGenPythiaEventHeader") {
501 headPy = (AliGenPythiaEventHeader*) htmp;
502 } else if (TString(htmp->IsA()->GetName()) == "AliGenDPMjetEventHeader") {
503 headPho = (AliGenDPMjetEventHeader*) htmp;
505 AliError("Unknown header");
508 // Determine process type
510 if(headPy->ProcessType() == 92 || headPy->ProcessType() == 93) {
513 } else if (headPy->ProcessType() == 94) {
514 // double diffractive
517 else if(headPy->ProcessType() != 92 && headPy->ProcessType() != 93 && headPy->ProcessType() != 94) {
521 } else if (headPho) {
522 if(headPho->ProcessType() == 5 || headPho->ProcessType() == 6 ) {
525 } else if (headPho->ProcessType() == 7) {
526 // double diffractive
528 } else if(headPho->ProcessType() != 5 && headPho->ProcessType() != 6 && headPho->ProcessType() != 7 ) {
535 // no process type found?
536 AliError(Form("Unknown header: %s", htmp->IsA()->GetName()));
540 Double_t AliCollisionNormalization::GetProcessWeight(Int_t proc) {
542 // Return a weight to be used when combining the MC histos to
543 // compute efficiency, defined as the ratio XS in generator / XS
546 Float_t ref_SD, ref_DD, ref_ND, error_SD, error_DD, error_ND;
547 GetRelativeFractions(fReferenceXS,ref_SD, ref_DD, ref_ND, error_SD, error_DD, error_ND);
549 static Double_t total = fHistProcTypes->Integral();
550 if (fHistProcTypes->GetBinContent(fHistProcTypes->FindBin(kProcUnknown)) > 0) {
551 AliError("There were unknown processes!!!");
553 static Double_t SD = fHistProcTypes->GetBinContent(fHistProcTypes->FindBin(kProcSD))/total;
554 static Double_t DD = fHistProcTypes->GetBinContent(fHistProcTypes->FindBin(kProcDD))/total;
555 static Double_t ND = fHistProcTypes->GetBinContent(fHistProcTypes->FindBin(kProcND))/total;
558 AliInfo(Form("Total MC evts: %f",total));
559 AliInfo(Form(" Frac SD %4.4f", SD));
560 AliInfo(Form(" Frac DD %4.4f", DD));
561 AliInfo(Form(" Frac ND %4.4f", ND));
575 AliError("Unknown process");
582 void AliCollisionNormalization::GetRelativeFractions(Int_t origin, Float_t& ref_SD, Float_t& ref_DD, Float_t& ref_ND, Float_t& error_SD, Float_t& error_DD, Float_t& error_ND)
584 // Returns fraction of XS (SD, ND, DD) and corresponding error
585 // Stolen from Jan Fiete's drawSystematics macro
588 // -1 = Pythia (test)
595 Double_t epsilon = 0.0001;
596 if(TMath::Abs(fEnergy-900)<epsilon) {
600 case -10: // Pythia default at 7 GeV, 50% error
601 AliInfo("PYTHIA x-sections");
602 ref_SD = 0.192637; error_SD = ref_SD * 0.5;
603 ref_DD = 0.129877; error_DD = ref_DD * 0.5;
604 ref_ND = 0.677486; error_ND = 0;
607 case -1: // Pythia default at 900 GeV, as test
608 AliInfo("PYTHIA x-sections");
615 AliInfo("UA5 x-sections a la first paper");
616 ref_SD = 0.153; error_SD = 0.05;
617 ref_DD = 0.080; error_DD = 0.05;
618 ref_ND = 0.767; error_ND = 0;
622 AliInfo("UA5 x-sections hadron level definition for Pythia");
623 // Fractions in Pythia with UA5 cuts selection for SD
625 // SD: 0.188588 --> this should be 15.3
627 ref_SD = 0.224 * 0.153 / 0.189; error_SD = 0.023 * 0.224 / 0.189;
628 ref_DD = 0.095; error_DD = 0.06;
629 ref_ND = 1.0 - ref_SD - ref_DD; error_ND = 0;
633 AliInfo("UA5 x-sections hadron level definition for Phojet");
634 // Fractions in Phojet with UA5 cuts selection for SD
636 // SD: 0.151601 --> this should be 15.3
638 ref_SD = 0.191 * 0.153 / 0.152; error_SD = 0.023 * 0.191 / 0.152;
639 ref_DD = 0.095; error_DD = 0.06;
640 ref_ND = 1.0 - ref_SD - ref_DD; error_ND = 0;
642 case 2: // tel-aviv model
643 AliInfo("Tel-aviv model x-sections");
646 ref_ND = 1 - ref_SD - ref_DD;
649 case 3: // durham model
650 AliInfo("Durham model x-sections");
653 ref_ND = 1 - ref_SD - ref_DD;
656 AliFatal(Form("Unknown origin %d, Energy %f", origin, fEnergy));
659 else if(TMath::Abs(fEnergy-1800)<epsilon) {
662 case 20: // E710, 1.8 TeV
663 AliInfo("E710 x-sections hadron level definition for Pythia");
665 // SD: 0.166590 --> this should be 15.9
667 ref_SD = 0.217 * 0.159 / 0.167; error_SD = 0.024 * 0.217 / 0.167;
668 ref_DD = 0.075 * 1.43; error_DD = 0.02 * 1.43;
669 ref_ND = 1.0 - ref_SD - ref_DD; error_ND = 0;
672 case 21: // E710, 1.8 TeV
673 AliInfo("E710 x-sections hadron level definition for Phojet");
675 // SD: 0.125506 --> this should be 15.9
677 ref_SD = 0.161 * 0.159 / 0.126; error_SD = 0.024 * 0.161 / 0.126;
678 ref_DD = 0.075 * 1.43; error_DD = 0.02 * 1.43;
679 ref_ND = 1.0 - ref_SD - ref_DD; error_ND = 0;
682 case 1: // data 1.8 TeV
683 AliInfo("??? x-sections");
686 ref_ND = 1 - ref_SD - ref_DD;
689 AliFatal(Form("Unknown origin %d, Energy %f", origin, fEnergy));
693 AliFatal(Form("Unknown energy %f", fEnergy));