]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliForwardMCCorrectionsTask.cxx
Set branches
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardMCCorrectionsTask.cxx
1 // 
2 // Calculate the corrections in the forward regions
3 // 
4 // Inputs: 
5 //   - AliESDEvent 
6 //
7 // Outputs: 
8 //   - AliAODForwardMult 
9 // 
10 // Histograms 
11 //   
12 // Corrections used 
13 // 
14 #include "AliForwardMCCorrectionsTask.h"
15 #include "AliTriggerAnalysis.h"
16 #include "AliPhysicsSelection.h"
17 #include "AliLog.h"
18 #include "AliHeader.h"
19 #include "AliGenEventHeader.h"
20 #include "AliESDEvent.h"
21 #include "AliAODHandler.h"
22 #include "AliMultiplicity.h"
23 #include "AliInputEventHandler.h"
24 #include "AliStack.h"
25 #include "AliMCEvent.h"
26 #include "AliAODForwardMult.h"
27 #include "AliFMDStripIndex.h"
28 #include "AliFMDCorrSecondaryMap.h"
29 #include <TH1.h>
30 #include <TH2D.h>
31 #include <TDirectory.h>
32 #include <TTree.h>
33 #include <TList.h>
34 #include <TROOT.h>
35 #include <iostream>
36
37 //====================================================================
38 namespace {
39   const char* GetEventName(Bool_t tr, Bool_t vtx) 
40   {
41     return Form("nEvents%s%s", (tr ? "Tr" : ""), (vtx ? "Vtx" : ""));
42   }
43 #if 0
44   const char* GetHitsName(UShort_t d, Char_t r) 
45   {
46     return Form("hitsFMD%d%c", d, r);
47   }
48   const char* GetStripsName(UShort_t d, Char_t r)
49   {
50     return Form("stripsFMD%d%c", d, r);
51   }
52   const char* GetPrimaryName(Char_t r, Bool_t trVtx)
53   {
54     return Form("primaries%s%s", 
55                 ((r == 'I' || r == 'i') ? "Inner" : "Outer"), 
56                 (trVtx ? "TrVtx" : "All"));
57   }
58 #endif
59 }
60
61 //====================================================================
62 AliForwardMCCorrectionsTask::AliForwardMCCorrectionsTask()
63   : AliAnalysisTaskSE(),
64     fInspector(),
65     fTrackDensity(),
66     fESDFMD(),
67     fVtxBins(0),
68     fFirstEvent(true),
69     fHEvents(0), 
70     fHEventsTr(0), 
71     fHEventsTrVtx(0),
72     fVtxAxis(),
73     fEtaAxis(),
74     fList()
75 {
76   // 
77   // Constructor 
78   // 
79   // Parameters:
80   //    name Name of task 
81   //
82 }
83
84 //____________________________________________________________________
85 AliForwardMCCorrectionsTask::AliForwardMCCorrectionsTask(const char* name)
86   : AliAnalysisTaskSE(name),
87     fInspector("eventInspector"), 
88     fTrackDensity("trackDensity"),
89     fESDFMD(),
90     fVtxBins(0),
91     fFirstEvent(true),
92     fHEvents(0), 
93     fHEventsTr(0), 
94     fHEventsTrVtx(0),
95     fVtxAxis(10,-10,10), 
96     fEtaAxis(200,-4,6),
97     fList()
98 {
99   // 
100   // Constructor 
101   // 
102   // Parameters:
103   //    name Name of task 
104   //
105   DefineOutput(1, TList::Class());
106   DefineOutput(2, TList::Class());
107   fBranchNames = 
108     "ESD:AliESDRun.,AliESDHeader.,AliMultiplicity.,"
109     "AliESDFMD.,SPDVertex.,PrimaryVertex.";
110 }
111
112 //____________________________________________________________________
113 AliForwardMCCorrectionsTask::AliForwardMCCorrectionsTask(const AliForwardMCCorrectionsTask& o)
114   : AliAnalysisTaskSE(o),
115     fInspector(o.fInspector),
116     fTrackDensity(),
117     fESDFMD(o.fESDFMD),
118     fVtxBins(0),
119     fFirstEvent(o.fFirstEvent),
120     fHEvents(o.fHEvents), 
121     fHEventsTr(o.fHEventsTr), 
122     fHEventsTrVtx(o.fHEventsTrVtx),
123     fVtxAxis(10,-10,10), 
124     fEtaAxis(200,-4,6),
125     fList(o.fList)
126 {
127   // 
128   // Copy constructor 
129   // 
130   // Parameters:
131   //    o Object to copy from 
132   //
133 }
134
135 //____________________________________________________________________
136 AliForwardMCCorrectionsTask&
137 AliForwardMCCorrectionsTask::operator=(const AliForwardMCCorrectionsTask& o)
138 {
139   // 
140   // Assignment operator 
141   // 
142   // Parameters:
143   //    o Object to assign from 
144   // 
145   // Return:
146   //    Reference to this object 
147   //
148   if (&o == this) return *this;
149   fInspector         = o.fInspector;
150   fTrackDensity      = o.fTrackDensity;
151   fESDFMD            = o.fESDFMD;
152   fVtxBins           = o.fVtxBins;
153   fFirstEvent        = o.fFirstEvent;
154   fHEvents           = o.fHEvents;
155   fHEventsTr         = o.fHEventsTr;
156   fHEventsTrVtx      = o.fHEventsTrVtx;
157   SetVertexAxis(o.fVtxAxis);
158   SetEtaAxis(o.fEtaAxis);
159
160   return *this;
161 }
162
163 //____________________________________________________________________
164 void
165 AliForwardMCCorrectionsTask::Init()
166 {
167   // 
168   // Initialize the task 
169   // 
170   //
171 }
172
173 //____________________________________________________________________
174 void
175 AliForwardMCCorrectionsTask::SetVertexAxis(Int_t nBin, Double_t min, 
176                                            Double_t max)
177 {
178   // 
179   // Set the vertex axis to use
180   // 
181   // Parameters:
182   //    nBins Number of bins
183   //    vzMin Least @f$z@f$ coordinate of interation point
184   //    vzMax Largest @f$z@f$ coordinate of interation point
185   //
186   if (max < min) { 
187     Double_t tmp = min;
188     min          = max;
189     max          = tmp;
190   }
191 /*
192   if (min < -10) 
193     AliWarning(Form("Minimum vertex %f < -10, make sure you want this",min));
194   if (max > +10) 
195     AliWarning(Form("Minimum vertex %f > +10, make sure you want this",max));
196 */
197   fVtxAxis.Set(nBin, min, max);
198 }
199 //____________________________________________________________________
200 void
201 AliForwardMCCorrectionsTask::SetVertexAxis(const TAxis& axis)
202 {
203   // 
204   // Set the vertex axis to use
205   // 
206   // Parameters:
207   //    axis Axis
208   //
209   SetVertexAxis(axis.GetNbins(),axis.GetXmin(),axis.GetXmax());
210 }
211
212 //____________________________________________________________________
213 void
214 AliForwardMCCorrectionsTask::SetEtaAxis(Int_t nBin, Double_t min, Double_t max)
215 {
216   // 
217   // Set the eta axis to use
218   // 
219   // Parameters:
220   //    nBins Number of bins
221   //    vzMin Least @f$\eta@f$ 
222   //    vzMax Largest @f$\eta@f$ 
223   //
224   if (max < min) { 
225     Double_t tmp = min;
226     min          = max;
227     max          = tmp;
228   }
229   if (min < -4) 
230     AliWarning(Form("Minimum eta %f < -4, make sure you want this",min));
231   if (max > +6) 
232     AliWarning(Form("Minimum eta %f > +6, make sure you want this",max));
233   fEtaAxis.Set(nBin, min, max);
234 }
235 //____________________________________________________________________
236 void
237 AliForwardMCCorrectionsTask::SetEtaAxis(const TAxis& axis)
238 {
239   // 
240   // Set the eta axis to use
241   // 
242   // Parameters:
243   //    axis Axis
244   //
245   SetEtaAxis(axis.GetNbins(),axis.GetXmin(),axis.GetXmax());
246 }
247
248 //____________________________________________________________________
249 void
250 AliForwardMCCorrectionsTask::DefineBins(TList* l)
251 {
252   if (!fVtxBins) fVtxBins = new TObjArray(fVtxAxis.GetNbins(), 1);
253   if (fVtxBins->GetEntries() > 0) return;
254
255   fVtxBins->SetOwner();
256   for (Int_t  i    = 1; i <= fVtxAxis.GetNbins(); i++) { 
257     Double_t  low  = fVtxAxis.GetBinLowEdge(i);
258     Double_t  high = fVtxAxis.GetBinUpEdge(i);
259     VtxBin*   bin  = new VtxBin(low, high, fEtaAxis);
260     fVtxBins->AddAt(bin, i);
261     bin->DefineOutput(l);
262   }
263 }
264
265 //____________________________________________________________________
266 void
267 AliForwardMCCorrectionsTask::UserCreateOutputObjects()
268 {
269   // 
270   // Create output objects 
271   // 
272   //
273   fList = new TList;
274   fList->SetOwner();
275   fList->SetName(Form("%sSums", GetName()));
276
277   DefineBins(fList);
278
279   fHEvents = new TH1I(GetEventName(false,false),
280                       "Number of all events", 
281                       fVtxAxis.GetNbins(), 
282                       fVtxAxis.GetXmin(), 
283                       fVtxAxis.GetXmax());
284   fHEvents->SetXTitle("v_{z} [cm]");
285   fHEvents->SetYTitle("# of events");
286   fHEvents->SetFillColor(kBlue+1);
287   fHEvents->SetFillStyle(3001);
288   fHEvents->SetDirectory(0);
289   fList->Add(fHEvents);
290
291   fHEventsTr = new TH1I(GetEventName(true, false), 
292                         "Number of triggered events",
293                         fVtxAxis.GetNbins(), 
294                         fVtxAxis.GetXmin(), 
295                         fVtxAxis.GetXmax());
296   fHEventsTr->SetXTitle("v_{z} [cm]");
297   fHEventsTr->SetYTitle("# of events");
298   fHEventsTr->SetFillColor(kRed+1);
299   fHEventsTr->SetFillStyle(3001);
300   fHEventsTr->SetDirectory(0);
301   fList->Add(fHEventsTr);
302
303   fHEventsTrVtx = new TH1I(GetEventName(true,true),
304                            "Number of events w/trigger and vertex", 
305                            fVtxAxis.GetNbins(), 
306                            fVtxAxis.GetXmin(), 
307                            fVtxAxis.GetXmax());
308   fHEventsTrVtx->SetXTitle("v_{z} [cm]");
309   fHEventsTrVtx->SetYTitle("# of events");
310   fHEventsTrVtx->SetFillColor(kBlue+1);
311   fHEventsTrVtx->SetFillStyle(3001);
312   fHEventsTrVtx->SetDirectory(0);
313   fList->Add(fHEventsTrVtx);
314
315   // Copy axis objects to output 
316   TH1* vtxAxis = new TH1D("vtxAxis", "Vertex axis", 
317                           fVtxAxis.GetNbins(), 
318                           fVtxAxis.GetXmin(), 
319                           fVtxAxis.GetXmax());
320   TH1* etaAxis = new TH1D("etaAxis", "Eta axis", 
321                           fEtaAxis.GetNbins(), 
322                           fEtaAxis.GetXmin(), 
323                           fEtaAxis.GetXmax());
324   fList->Add(vtxAxis);
325   fList->Add(etaAxis);
326
327   AliInfo(Form("Initialising sub-routines: %p, %p", 
328                &fInspector, &fTrackDensity));
329   fInspector.DefineOutput(fList);
330   fTrackDensity.DefineOutput(fList);
331
332   PostData(1, fList);
333 }
334 //____________________________________________________________________
335 void
336 AliForwardMCCorrectionsTask::UserExec(Option_t*)
337 {
338   // 
339   // Process each event 
340   // 
341   // Parameters:
342   //    option Not used
343   //  
344
345   // Get the input data - MC event
346   AliMCEvent*  mcEvent = MCEvent();
347   if (!mcEvent) { 
348     AliWarning("No MC event found");
349     return;
350   }
351
352   // Get the input data - ESD event
353   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
354   if (!esd) { 
355     AliWarning("No ESD event found for input event");
356     return;
357   }
358
359   //--- Read run information -----------------------------------------
360   if (fFirstEvent && esd->GetESDRun()) {
361     fInspector.ReadRunDetails(esd);
362     
363     AliInfo(Form("Initializing with parameters from the ESD:\n"
364                  "         AliESDEvent::GetBeamEnergy()   ->%f\n"
365                  "         AliESDEvent::GetBeamType()     ->%s\n"
366                  "         AliESDEvent::GetCurrentL3()    ->%f\n"
367                  "         AliESDEvent::GetMagneticField()->%f\n"
368                  "         AliESDEvent::GetRunNumber()    ->%d\n",
369                  esd->GetBeamEnergy(),
370                  esd->GetBeamType(),
371                  esd->GetCurrentL3(),
372                  esd->GetMagneticField(),
373                  esd->GetRunNumber()));
374
375     fInspector.Init(fVtxAxis);
376
377     Print();
378     fFirstEvent = false;
379   }
380
381   // Some variables 
382   UInt_t   triggers; // Trigger bits
383   Bool_t   lowFlux;  // Low flux flag
384   UShort_t iVz;      // Vertex bin from ESD
385   Double_t vZ;       // Z coordinate from ESD
386   Double_t cent;     // Centrality 
387   UShort_t iVzMc;    // Vertex bin from MC
388   Double_t vZMc;     // Z coordinate of IP vertex from MC
389   Double_t b;        // Impact parameter
390   Double_t cMC;      // Centrality estimate from b
391   Int_t    nPart;    // Number of participants 
392   Int_t    nBin;     // Number of binary collisions 
393   Double_t phiR;     // Reaction plane from MC
394   UShort_t nClusters;// Number of SPD clusters 
395   // Process the data 
396   UInt_t retESD = fInspector.Process(esd, triggers, lowFlux, iVz, vZ, 
397                                      cent, nClusters);
398   fInspector.ProcessMC(mcEvent, triggers, iVzMc, vZMc, 
399                        b, cMC, nPart, nBin, phiR);
400
401   Bool_t isInel   = triggers & AliAODForwardMult::kInel;
402   Bool_t hasVtx   = retESD == AliFMDMCEventInspector::kOk;
403
404   // Fill the event count histograms 
405   if (isInel)           fHEventsTr->Fill(vZMc);
406   if (isInel && hasVtx) fHEventsTrVtx->Fill(vZMc);
407   fHEvents->Fill(vZMc);
408
409   // Now find our vertex bin object 
410   VtxBin* bin = 0;
411   if (iVzMc > 0 && iVzMc <= fVtxAxis.GetNbins()) 
412     bin = static_cast<VtxBin*>(fVtxBins->At(iVzMc));
413   if (!bin) { 
414     // AliError(Form("No vertex bin object @ %d (%f)", iVzMc, vZMc));
415     return;
416   }
417
418   // Clear our ESD object 
419   fESDFMD.Clear();
420
421   // Get FMD data 
422   AliESDFMD* esdFMD = esd->GetFMDData();
423   
424   // Now process our input data and store in own ESD object 
425   fTrackDensity.Calculate(*esdFMD, *mcEvent, vZMc, fESDFMD, bin->fPrimary);
426   bin->fCounts->Fill(0.5);
427
428   // And then bin the data in our vtxbin 
429   for (UShort_t d=1; d<=3; d++) { 
430     UShort_t nr = (d == 1 ? 1 : 2);
431     for (UShort_t q=0; q<nr; q++) { 
432       Char_t      r = (q == 0 ? 'I' : 'O');
433       UShort_t    ns= (q == 0 ?  20 :  40);
434       UShort_t    nt= (q == 0 ? 512 : 256);
435       TH2D*       h = bin->fHists.Get(d,r);
436
437       for (UShort_t s=0; s<ns; s++) { 
438         for (UShort_t t=0; t<nt; t++) {
439           Float_t mult = fESDFMD.Multiplicity(d,r,s,t);
440           
441           if (mult == 0 || mult > 20) continue;
442
443           Float_t phi = fESDFMD.Phi(d,r,s,t) / 180 * TMath::Pi();
444           Float_t eta = fESDFMD.Eta(d,r,s,t);
445           h->Fill(eta,phi,mult);
446         } // for t
447       } // for s 
448     } // for q 
449   } // for d
450
451   PostData(1, fList);
452 }
453 //____________________________________________________________________
454 void
455 AliForwardMCCorrectionsTask::Terminate(Option_t*)
456 {
457   // 
458   // End of job
459   // 
460   // Parameters:
461   //    option Not used 
462   //
463
464   fList = dynamic_cast<TList*>(GetOutputData(1));
465   if (!fList) {
466     AliError("No output list defined");
467     return;
468   }
469
470   DefineBins(fList);
471
472   // Output list 
473   TList* output = new TList;
474   output->SetOwner();
475   output->SetName(Form("%sResults", GetName()));
476
477   // --- Fill correction object --------------------------------------
478   AliFMDCorrSecondaryMap* corr = new AliFMDCorrSecondaryMap;
479   corr->SetVertexAxis(fVtxAxis);
480   corr->SetEtaAxis(fEtaAxis);
481   
482   TIter     next(fVtxBins);
483   VtxBin*   bin = 0;
484   UShort_t  iVz = 1;
485   while ((bin = static_cast<VtxBin*>(next()))) 
486     bin->Finish(fList, output, iVz++, corr);
487
488   output->Add(corr);
489
490   PostData(2, output);
491 }
492
493 //____________________________________________________________________
494 void
495 AliForwardMCCorrectionsTask::Print(Option_t* option) const
496 {
497   std::cout << ClassName() << "\n"
498             << "  Vertex bins:      " << fVtxAxis.GetNbins() << '\n'
499             << "  Vertex range:     [" << fVtxAxis.GetXmin() 
500             << "," << fVtxAxis.GetXmax() << "]\n"
501             << "  Eta bins:         " << fEtaAxis.GetNbins() << '\n'
502             << "  Eta range:        [" << fEtaAxis.GetXmin() 
503             << "," << fEtaAxis.GetXmax() << "]"
504             << std::endl;
505   gROOT->IncreaseDirLevel();
506   fInspector.Print(option);
507   fTrackDensity.Print(option);
508   gROOT->DecreaseDirLevel();
509 }
510
511 //====================================================================
512 const char*
513 AliForwardMCCorrectionsTask::VtxBin::BinName(Double_t low, 
514                                              Double_t high) 
515 {
516   static TString buf;
517   buf = Form("vtx%+05.1f_%+05.1f", low, high);
518   buf.ReplaceAll("+", "p");
519   buf.ReplaceAll("-", "m");
520   buf.ReplaceAll(".", "d");
521   return buf.Data();
522 }
523
524
525 //____________________________________________________________________
526 AliForwardMCCorrectionsTask::VtxBin::VtxBin()
527   : fHists(), 
528     fPrimary(0),
529     fCounts(0)
530 {
531 }
532 //____________________________________________________________________
533 AliForwardMCCorrectionsTask::VtxBin::VtxBin(Double_t low, 
534                                             Double_t high, 
535                                             const TAxis& axis)
536   : TNamed(BinName(low, high), 
537            Form("%+5.1fcm<v_{z}<%+5.1fcm", low, high)),
538     fHists(), 
539     fPrimary(0),
540     fCounts(0)
541 {
542   fHists.Init(axis);
543
544   fPrimary = new TH2D("primary", "Primaries", 
545                       axis.GetNbins(), axis.GetXmin(), axis.GetXmax(), 
546                       40, 0, 2*TMath::Pi());
547   fPrimary->SetXTitle("#eta");
548   fPrimary->SetYTitle("#varphi [radians]");
549   fPrimary->Sumw2();
550   fPrimary->SetDirectory(0);
551
552   fCounts = new TH1D("counts", "Counts", 1, 0, 1);
553   fCounts->SetXTitle("Events");
554   fCounts->SetYTitle("# of Events");
555   fCounts->SetDirectory(0);
556 }
557
558 //____________________________________________________________________
559 AliForwardMCCorrectionsTask::VtxBin::VtxBin(const VtxBin& o)
560   : TNamed(o),
561     fHists(o.fHists),
562     fPrimary(0), 
563     fCounts(0)
564 {
565   if (o.fPrimary) {
566     fPrimary = static_cast<TH2D*>(o.fPrimary->Clone());
567     fPrimary->SetDirectory(0);
568   }
569   if (o.fCounts) {
570     fCounts = static_cast<TH1D*>(o.fCounts->Clone());
571     fCounts->SetDirectory(0);
572   }
573 }
574
575 //____________________________________________________________________
576 AliForwardMCCorrectionsTask::VtxBin&
577 AliForwardMCCorrectionsTask::VtxBin::operator=(const VtxBin& o)
578 {
579   if (&o == this) return *this;
580   TNamed::operator=(o);
581   fHists   = o.fHists;
582   fPrimary = 0;
583   fCounts  = 0;
584   if (o.fPrimary) {
585     fPrimary = static_cast<TH2D*>(o.fPrimary->Clone());
586     fPrimary->SetDirectory(0);
587   }
588   if (o.fCounts) {
589     fCounts = static_cast<TH1D*>(o.fCounts->Clone());
590     fCounts->SetDirectory(0);
591   }
592   return *this;
593 }
594
595 //____________________________________________________________________
596 void
597 AliForwardMCCorrectionsTask::VtxBin::DefineOutput(TList* l)
598 {
599   TList* d = new TList;
600   d->SetName(GetName());
601   d->SetOwner();
602   l->Add(d);
603
604   d->Add(fHists.fFMD1i);
605   d->Add(fHists.fFMD2i);
606   d->Add(fHists.fFMD2o);
607   d->Add(fHists.fFMD3i);
608   d->Add(fHists.fFMD3o);
609
610   d->Add(fPrimary);
611   d->Add(fCounts);
612 }
613
614 //____________________________________________________________________
615 TH2D*
616 AliForwardMCCorrectionsTask::VtxBin::MakeBg(const TH2D* hits, 
617                                             const TH2D* primary) const
618 {
619   TH2D* h = static_cast<TH2D*>(hits->Clone());
620   h->SetDirectory(0);
621   TString n(h->GetName());
622   n.ReplaceAll("_cache", "");
623   h->SetName(n);
624   h->Divide(primary);
625   
626   return h;
627 }
628   
629 //____________________________________________________________________
630 void
631 AliForwardMCCorrectionsTask::VtxBin::Finish(const TList* input, 
632                                             TList* output, 
633                                             UShort_t iVz,
634                                             AliFMDCorrSecondaryMap* map)
635 {
636   TList* out = new TList;
637   out->SetName(GetName());
638   out->SetOwner();
639   output->Add(out);
640
641   TList* l = static_cast<TList*>(input->FindObject(GetName()));
642   if (!l) { 
643     AliError(Form("List %s not found in %s", GetName(), input->GetName()));
644     return;
645   }
646
647   TH2D*   fmd1i = static_cast<TH2D*>(l->FindObject("FMD1I_cache"));
648   TH2D*   fmd2i = static_cast<TH2D*>(l->FindObject("FMD2I_cache"));
649   TH2D*   fmd2o = static_cast<TH2D*>(l->FindObject("FMD2O_cache"));
650   TH2D*   fmd3i = static_cast<TH2D*>(l->FindObject("FMD3I_cache"));
651   TH2D*   fmd3o = static_cast<TH2D*>(l->FindObject("FMD3O_cache"));
652   TH2D*   primO = static_cast<TH2D*>(l->FindObject("primary"));
653   if (!fmd1i || !fmd2i || !fmd2o || !fmd3i || !fmd3o || !primO) {
654     AliError(Form("Missing histogram(s): %p,%p,%p,%p,%p,%p",
655                   fmd1i, fmd2i, fmd2o, fmd3i, fmd3o, primO));
656     return;
657   }
658
659   // Half coverage in phi for inners
660   TH2D*   primI = static_cast<TH2D*>(primO->Clone());
661   primI->SetDirectory(0);
662   primI->RebinY(2); 
663
664   TH2D* bg1i = MakeBg(fmd1i, primI);
665   TH2D* bg2i = MakeBg(fmd2i, primI);
666   TH2D* bg2o = MakeBg(fmd2o, primO);
667   TH2D* bg3i = MakeBg(fmd3i, primI);
668   TH2D* bg3o = MakeBg(fmd3o, primO);
669   map->SetCorrection(1, 'I', iVz, bg1i);
670   map->SetCorrection(2, 'I', iVz, bg2i);
671   map->SetCorrection(2, 'O', iVz, bg2o);
672   map->SetCorrection(3, 'I', iVz, bg3i);
673   map->SetCorrection(3, 'O', iVz, bg3o);
674   out->Add(bg1i);
675   out->Add(bg2i);
676   out->Add(bg2o);
677   out->Add(bg3i);
678   out->Add(bg3o);
679  
680 }
681
682 //
683 // EOF
684 //