]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/AliMergeSteer.C
Put back the doxygen comments in MUONTRGda.cxx
[u/mrichter/AliRoot.git] / macros / AliMergeSteer.C
1 ////////////////////////////////////////////////////////////////////////
2 //
3 //  Class to steer several steps in merging:
4 //       - extraction of vertex from bgr event
5 //       - event simulation
6 //       - creation of sdigits
7 //       - creation of digits
8 //       - analysis
9 //                  
10 //  Author: Jiri Chudoba (CERN), 2002
11 //
12 ////////////////////////////////////////////////////////////////////////
13
14 #if !defined(__CINT__) || defined(__MAKECINT__)
15 #include <Riostream.h>
16
17 // --- ROOT system ---
18
19 #include "TTask.h"
20 #include "TString.h"
21 #include "TROOT.h"
22 #include "TFile.h"
23 #include "TTree.h"
24 #include "TStopwatch.h"
25 #include "TArrayF.h"
26 #include "TSystem.h"
27 #include "TInterpreter.h"
28
29 // --- AliRoot header files ---
30
31 #include "STEER/AliRun.h"
32 #include "STEER/AliHeader.h"
33 #include "STEER/AliGenEventHeader.h"
34 #endif
35
36 class AliTRDdigitizer;
37 class AliPHOSSDigitizer;
38 class AliTOFSDigitizer;
39
40 class AliMergeSteer: public TTask {
41
42 public:
43   AliMergeSteer(const Text_t* name="AliMergeSteer",
44                   const Text_t* title="AliMergeSteer");
45 //  AliMergeSteer(Option_t *option);
46   virtual ~AliMergeSteer();
47   Bool_t SetDetectorFlag(Option_t* option, Int_t flag);
48   Int_t GetDetectorFlag(Option_t* option);  
49 //  void Print();
50   virtual void Exec(const Option_t *option);
51   Bool_t ImportgAlice(TFile *file);
52
53   Bool_t ExtractVertex(TString fn, Int_t eventNr);
54   void PrintVertex(TArrayF &vertex);
55   void ExportVertex(TArrayF &vertex);
56
57   Bool_t Simulate();
58   Bool_t CreateSDigits();
59   Bool_t Merge();
60   Bool_t NoMerge();
61   Bool_t ITSFastPoints(const char *outputFile, const char *inputFile);
62   Bool_t RecoMerged();
63   Bool_t RecoSignalOnly();
64   Bool_t CmpMerged();
65   Bool_t CmpSignalOnly();
66   
67   Bool_t AliCopy(TFile *inputFile, TFile *outputFile);
68   Bool_t AliCopy(TString inputFileName, TString outputFileName);
69   
70   
71   
72 // file names setter/getters  
73   void SetFileNameHits(TString fileName) {fFileNameHits = fileName;}
74   void SetFileNameSDigits(TString fileName) {fFileNameSDigits = fileName;}
75   void SetFileNameBgrHits(TString fileName) {fFileNameBgrHits = fileName;}
76   void SetFileNameBgrSDigits(TString fileName) {fFileNameBgrSDigits = fileName;}
77   void SetFileNameDigitsMerged(TString fileName) {fFileNameDigitsMerged = fileName;}
78   void SetFileNameDigitsSignalOnly(TString fileName) {fFileNameDigitsSignalOnly = fileName;}
79   TString FileNameHits() {return fFileNameHits;}
80   TString FileNameSDigits() {return fFileNameSDigits;}
81   TString FileNameBgr() {return fFileNameBgrSDigits;}
82   TString FileNameDigitsMerged() {return fFileNameDigitsMerged;}
83   TString FileNameDigitsSignalOnly() {return fFileNameDigitsSignalOnly;}
84
85 // flags
86   void SetFlagSim(Bool_t flag) {fFlagSim = flag;}
87   void SetFlagSDigits(Bool_t flag) {fFlagSDigits = flag;}
88   void SetFlagMerge(Bool_t flag) {fFlagMerge = flag;}
89   void SetFlagDigitsSignalOnly(Bool_t flag) {fFlagDigitsSignalOnly = flag;}
90   void SetFlagRecoMerged(Bool_t flag) {fFlagRecoMerged = flag;}
91   void SetFlagRecoSignalOnly(Bool_t flag) {fFlagRecoSignalOnly = flag;}
92   void SetFlagCmpMerged(Bool_t flag) {fFlagCmpMerged = flag;}
93   void SetFlagCmpSignalOnly(Bool_t flag) {fFlagCmpSignalOnly = flag;}
94   
95 // event numbers
96   void SetNEvents(Int_t nEvents) {fNEvents = nEvents;}
97   void SetBgrEventNr(Int_t i) {fBgrEventNr = i;}
98   
99 private:  
100
101   TString fFileNameHits;
102   TString fFileNameSDigits;
103   TString fFileNameBgrHits;
104   TString fFileNameBgrSDigits;
105   TString fFileNameDigitsMerged;
106   TString fFileNameDigitsSignalOnly;
107   TString fFileNameCmpMerged;
108   TString fFileNameCmpSignalOnly;
109   
110   Bool_t fFlagSim;
111   Bool_t fFlagSDigits;
112   Bool_t fFlagMerge;
113   Bool_t fFlagDigitsSignalOnly;
114   Bool_t fFlagRecoMerged;
115   Bool_t fFlagRecoSignalOnly;
116   Bool_t fFlagCmpMerged;
117   Bool_t fFlagCmpSignalOnly;
118   Bool_t fCopy2S;
119   Bool_t fCopy2D;
120   
121   Int_t fNEvents;
122   Int_t fFirstEvent;
123   Int_t fBgrEventNr;
124   
125   Int_t fDEBUG;
126   
127 // flags for detectors - determine which detector will be used
128 // during merging.
129 //  0 - do not use
130 //  1 - process the whole detector
131 //  2 - process only active region (region of interest on)  
132   
133   Int_t   fFMD;
134   Int_t   fITS;
135   Int_t   fMUON;
136   Int_t   fPHOS;
137   Int_t   fPMD;
138   Int_t   fHMPID;
139   Int_t   fT0;
140   Int_t   fTOF;
141   Int_t   fTPC;
142   Int_t   fTRD;
143   Int_t   fZDC;
144   Int_t   fEMCAL;
145   
146 };
147
148
149 ////////////////////////////////////////////////////////////////////////
150 //
151 // AliMergeSteer.cxx
152 //
153 //
154 ////////////////////////////////////////////////////////////////////////
155
156 AliMergeSteer::AliMergeSteer(const Text_t *name, const Text_t* title) : TTask(name,title)
157 {
158 //
159 // default ctor
160 //
161   fFileNameHits="galice.root";
162   fFileNameSDigits = "sdigits.root";
163   fFileNameBgrSDigits = "bgr.sdigits.root";
164   fFileNameBgrHits = "bgr.hits.root";
165   fFileNameDigitsMerged = "digits.root";
166   fFileNameDigitsSignalOnly = "digits.signal.root";
167   fFileNameCmpMerged = "CmpGaRS_merged.root";
168   fFileNameCmpSignalOnly = "CmpGaRS_signal.root";
169
170   fFlagSim = kFALSE;
171   fFlagSDigits = kFALSE;
172   fFlagMerge = kFALSE;
173   fFlagDigitsSignalOnly = kFALSE;
174   fFlagRecoMerged = kFALSE;
175   fFlagRecoSignalOnly = kFALSE;
176   fFlagCmpMerged = kFALSE;
177   fFlagCmpSignalOnly = kFALSE;    
178   fCopy2S = kTRUE;
179   fCopy2D = kTRUE;
180   
181   fNEvents = 1;
182   fFirstEvent = 0;
183   fBgrEventNr = 0;
184
185   fDEBUG = 1;
186
187   fFMD = 0;
188   fITS = 0;
189   fMUON = 0;
190   fPHOS = 0;
191   fPMD = 0;
192   fHMPID = 0;
193   fT0 = 0;
194   fTOF = 0;
195   fTPC = 0;
196   fTRD = 0;
197   fZDC = 0;
198   fEMCAL = 0;
199 }
200
201 ////////////////////////////////////////////////////////////////////////
202 AliMergeSteer::~AliMergeSteer()
203 {
204 // default dtor
205
206 }
207 ////////////////////////////////////////////////////////////////////////
208 void AliMergeSteer::Exec(Option_t* option)
209 {
210   Bool_t rc = kTRUE;
211   if (gAlice) delete gAlice;
212   gAlice = 0;
213
214   if (fFlagSim) rc = Simulate();
215   if (!rc) {Error("Exec","Simulate wrong"); return;}
216
217   if (fFlagSDigits) rc = CreateSDigits();
218   if (!rc) {Error("Exec","CreateSDigits wrong"); return;}
219
220   if (fFlagSDigits && fCopy2S) rc = AliCopy(fFileNameHits,fFileNameSDigits);
221   if (!rc) {Error("Exec","AliCopy to SD wrong"); return;}
222
223   if (fFlagMerge && fCopy2D) rc = AliCopy(fFileNameHits,fFileNameDigitsMerged);
224   if (!rc) {Error("Exec","AliCopy to DigitsMerged wrong"); return;}
225
226   if (fFlagMerge) rc = Merge();
227   if (!rc) {Error("Exec","Merge wrong"); return;}
228
229   if (fFlagDigitsSignalOnly) rc = NoMerge();
230   if (!rc) {Error("Exec","NoMerge wrong"); return;}
231   
232   if (fFlagRecoMerged) rc = RecoMerged();
233   if (!rc) {Error("Exec","RecoMerged wrong"); return;}
234   
235   if (fFlagRecoSignalOnly) rc = RecoSignalOnly();
236   if (!rc) {Error("Exec","RecoSignalOnly wrong"); return;}
237
238   if (fFlagCmpMerged) rc = CmpMerged();
239   if (!rc) {Error("Exec","CmpMerged wrong"); return;}
240   
241   if (fFlagCmpSignalOnly) rc = CmpSignalOnly();
242   if (!rc) {Error("Exec","CmpSignalOnly wrong"); return;}
243   
244  
245
246
247 }
248 ////////////////////////////////////////////////////////////////////////
249 Int_t AliMergeSteer::GetDetectorFlag(Option_t* option)
250 {
251 //
252 // return current flag value for a given detector
253 //
254   if (strstr(option,"FMD")) {
255     return fFMD;
256   } else if (strstr(option,"ITS")) {
257     return fITS;
258   } else if (strstr(option,"MUON")) {
259     return fMUON;
260   } else if (strstr(option,"PHOS")) {
261     return fPHOS;
262   } else if (strstr(option,"PMD")) {
263     return fPMD;
264   } else if (strstr(option,"HMPID")) {
265     return fHMPID;
266   } else if (strstr(option,"T0")) {
267     return fT0;
268   } else if (strstr(option,"TOF")) {
269     return fTOF;
270   } else if (strstr(option,"TPC")) {
271     return fTPC;
272   } else if (strstr(option,"TRD")) {
273     return fTRD;
274   } else if (strstr(option,"ZDC")) {
275     return fZDC;
276   } else if (strstr(option,"EMCAL")) {
277     return fEMCAL;
278   } else {
279     cerr<<"Unknown detector required."<<endl;
280     return -1;
281   }
282 }
283 ////////////////////////////////////////////////////////////////////////
284 Bool_t AliMergeSteer::SetDetectorFlag(Option_t* option, Int_t flag)
285 {
286 //
287 // return current flag value for a given detector
288 //
289   if (strstr(option,"FMD")) {
290     fFMD = flag;
291   } else if (strstr(option,"ITS")) {
292     fITS = flag;
293   } else if (strstr(option,"MUON")) {
294     fMUON = flag;
295   } else if (strstr(option,"PHOS")) {
296     fPHOS = flag;
297   } else if (strstr(option,"PMD")) {
298     fPMD = flag;
299   } else if (strstr(option,"HMPID")) {
300     fHMPID = flag;
301   } else if (strstr(option,"T0")) {
302     fT0 = flag;
303   } else if (strstr(option,"TOF")) {
304     fTOF = flag;
305   } else if (strstr(option,"TPC")) {
306     fTPC = flag;
307   } else if (strstr(option,"TRD")) {
308     fTRD = flag;
309   } else if (strstr(option,"ZDC")) {
310     fZDC = flag;
311   } else if (strstr(option,"EMCAL")) {
312     fEMCAL = flag;
313   } else {
314     cerr<<"Unknown detector required."<<endl;
315     return kFALSE;
316   }
317   return kTRUE;
318 }
319 ////////////////////////////////////////////////////////////////////////
320 Bool_t AliMergeSteer::ImportgAlice(TFile *file) {
321 // read in gAlice object from the file
322   gAlice = (AliRun*)file->Get("gAlice");
323   if (!gAlice)  return kFALSE;
324   return kTRUE;
325 }
326 ////////////////////////////////////////////////////////////////////////
327 Bool_t AliMergeSteer::ExtractVertex(TString fn, Int_t eventNr)
328 {
329 // Open file with TPC geom and digits
330   TFile *file=TFile::Open(fn);
331   if (!file->IsOpen()) {cerr<<"Cannnot open "<<fn.Data()<<" !\n"; return kFALSE;}
332   if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
333     cerr<<"gAlice was not found in "<<fn.Data()<<" !\n";
334     return kFALSE;
335   }
336   
337   AliHeader *header = gAlice->GetHeader();
338   if (!header) {
339     cerr<<"header was not found in "<<fn.Data()<<" !\n";
340     return kFALSE;
341   } 
342   AliGenEventHeader* genEventHeader = header->GenEventHeader();
343   if (!genEventHeader) {
344     cerr<<"GenEventHeader was not found in "<<fn.Data()<<" !\n";
345     return kFALSE;
346   } 
347
348   TArrayF primaryVertex(3);
349   genEventHeader->PrimaryVertex(primaryVertex);
350   PrintVertex(primaryVertex);
351   ExportVertex(primaryVertex);
352 //  delete header;
353
354 // Following two lines should be there, but ....
355 //  delete genEventHeader;
356 //  delete gAlice;
357   gAlice = 0;
358   file->Close();
359     
360   return kTRUE;
361
362 }
363 ////////////////////////////////////////////////////////////////////////
364 void AliMergeSteer::PrintVertex(TArrayF &primaryVertex) 
365 {
366   cout <<"CONFIG_VERTEX: "
367        <<primaryVertex[0]<<" "
368        <<primaryVertex[1]<<" "
369        <<primaryVertex[2]<<" "<<endl;
370   return;
371
372
373 ////////////////////////////////////////////////////////////////////////
374 void AliMergeSteer::ExportVertex(TArrayF &primaryVertex) 
375 {
376   char vertexAsString[30];
377   sprintf(vertexAsString,"%f",primaryVertex[0]);
378   gSystem->Setenv("CONFIG_VERTEX_X",vertexAsString);
379   sprintf(vertexAsString,"%f",primaryVertex[1]);
380   gSystem->Setenv("CONFIG_VERTEX_Y",vertexAsString);
381   sprintf(vertexAsString,"%f",primaryVertex[2]);
382   gSystem->Setenv("CONFIG_VERTEX_Z",vertexAsString);
383   return;
384 }
385 ////////////////////////////////////////////////////////////////////////
386 Bool_t AliMergeSteer::Simulate()
387 {
388   char cmd[100];
389   TFile *f;
390
391   TDirectory *saveDir = gDirectory;
392   cout<< "The original directory is: ";
393   saveDir->pwd();
394   if (!ExtractVertex(fFileNameBgrHits,fBgrEventNr)) {
395     cerr<<" ExtractVertexAndSimulateSignal:  Error in ExtractVertex"<<endl;
396     return kFALSE;
397   }
398   saveDir->cd();
399   new AliRun("gAlice","Signal for Merging");
400   gAlice->Init("Config.C");
401 //    gSystem->Setenv("CONFIG_FILE",fileNameSigHits);
402   TStopwatch timer;
403   timer.Start();
404   gAlice->Run(fNEvents);
405   timer.Stop();
406   cout<<"Simulation of "<<fNEvents<<" of signal took: "<<endl;
407   timer.Print();
408   delete gAlice;
409   gAlice = 0;
410   f = static_cast<TFile *>(gROOT->FindObject("galice.root"));
411   if (f) f->Close();
412   f = 0;    
413   sprintf(cmd,"mv galice.root %s",fFileNameHits.Data());
414   gSystem->Exec(cmd);
415   return kTRUE;
416  
417 }
418
419 ////////////////////////////////////////////////////////////////////////
420 Bool_t AliMergeSteer::CreateSDigits()
421 {
422   char macroName[200];
423   char funcName[200]; 
424   sprintf(macroName,"AliHits2SDigits.C");
425   sprintf(funcName,
426           "AliHits2SDigits(\"%s\",\"%s\",%d,%d,%d,%d,%d,%d,%d,%d);",
427           fFileNameSDigits.Data(),fFileNameHits.Data(),
428           fNEvents,0,fITS,fTPC,fTRD,fPHOS,fTOF,0);
429   cerr<<"I'll do: "<<funcName<<endl;
430   gROOT->LoadMacro(macroName);
431   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
432   gInterpreter->ProcessLine(funcName);
433   if (fDEBUG) cerr<<"SDigits created"<<endl;
434   return kTRUE;
435 }
436
437
438 ////////////////////////////////////////////////////////////////////////
439 Bool_t AliMergeSteer::Merge()
440 {
441   char macroName[200];
442   char funcName[200]; 
443   cerr<<"Start merging"<<endl;
444   sprintf(macroName,"MergeV1.C");
445   sprintf(funcName,
446           "Merge(\"%s\",\"%s\",\"%s\",%d,%d,%d,%d,%d,%d,%d,%d);",
447           fFileNameDigitsMerged.Data(),fFileNameSDigits.Data(),
448           fFileNameBgrSDigits.Data(),
449           fNEvents,fITS,fTPC,fTRD,fPHOS,fMUON,fHMPID,0);
450   cerr<<"I'll do: "<<funcName<<endl;
451   gROOT->LoadMacro(macroName);
452   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
453   gInterpreter->ProcessLine(funcName);  
454   if (fDEBUG) cerr<<"Merging done"<<endl;
455
456 //  return kTRUE;
457 // add ITS fast points, no merging yet
458 //  return ITSFastPoints(fFileNameDigitsMerged.Data(),
459 //                     fFileNameHits.Data());
460
461 }
462 ////////////////////////////////////////////////////////////////////////
463 Bool_t AliMergeSteer::NoMerge()
464 {
465   char macroName[200];
466   char funcName[200]; 
467   cerr<<"Start NoMerging"<<endl;
468   sprintf(macroName,"AliSDigits2Digits.C");
469   sprintf(funcName,
470           "AliSDigits2Digits(\"%s\",\"%s\",%d,%d,%d,%d,%d,%d,%d,%d);",
471           fFileNameDigitsSignalOnly.Data(),fFileNameSDigits.Data(),
472           fNEvents,fITS,fTPC,fTRD,fPHOS,fMUON,fHMPID,0);
473   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
474   gROOT->LoadMacro(macroName);
475   gInterpreter->ProcessLine(funcName);  
476   if (fDEBUG) cerr<<"NoMerging done"<<endl;
477 //  return kTRUE;
478 // add ITS fast points, no merging
479 //  return ITSFastPoints(fFileNameDigitsSignalOnly.Data(),
480 //                     fFileNameHits.Data());
481
482 }
483 ////////////////////////////////////////////////////////////////////////
484 Bool_t AliMergeSteer::ITSFastPoints(const char *outputFile, const char *inputFile) {
485
486   char macroName[200];
487   char funcName[200]; 
488   sprintf(macroName,"AliITSHits2SDR.C");
489   sprintf(funcName,"AliITSH2FR2files(\"%s\",\"%s\");",
490           inputFile,  outputFile);
491   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
492   gROOT->LoadMacro(macroName);
493   gInterpreter->ProcessLine(funcName);  
494   if (fDEBUG) cerr<<"ITSFastPoints done"<<endl;
495
496   return kTRUE;
497 }
498
499 ////////////////////////////////////////////////////////////////////////
500 Bool_t AliMergeSteer::RecoMerged()
501 {
502 //
503 //
504 //
505   char macroName[200];
506   char funcName[200]; 
507   cerr<<"Start RecoMerged"<<endl;
508   sprintf(macroName,"AliBarrelRecoV3.C");
509   sprintf(funcName,"AliBarrelRecoMerged(%d);",fNEvents);
510   gROOT->LoadMacro(macroName);
511   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
512   gInterpreter->ProcessLine(funcName);
513   if (fDEBUG) cerr<<"RecoMerged done"<<endl;
514   return kTRUE;
515 }
516
517 ////////////////////////////////////////////////////////////////////////
518 Bool_t AliMergeSteer::RecoSignalOnly()
519 {
520 //
521 //
522 //
523   char macroName[200];
524   char funcName[200]; 
525   cerr<<"Start RecoSignalOnly"<<endl;
526   sprintf(macroName,"AliBarrelRecoNoITSClass.C");
527   sprintf(funcName,"AliBarrelReco(%d);",fNEvents);
528   gROOT->LoadMacro(macroName);
529   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
530   gInterpreter->ProcessLine(funcName);  
531   if (fDEBUG) cerr<<"RecoSignalOnly done"<<endl;
532   return kTRUE;
533 }
534
535 ////////////////////////////////////////////////////////////////////////
536 Bool_t AliMergeSteer::CmpMerged()
537 {
538 //
539 //
540 //
541   char macroName[200];
542   char funcName[200]; 
543   cerr<<"Start CmpMerged"<<endl;
544   sprintf(macroName,"CmpGaRS.C");
545   sprintf(funcName,
546           "CmpGaRS(%d,%d,\"%s\",\"AliTPCtracks_merged.root\",\"%s\");",
547           fNEvents, fFirstEvent, fFileNameHits.Data(),
548           fFileNameCmpMerged.Data());
549   gROOT->LoadMacro(macroName);
550   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
551   gInterpreter->ProcessLine(funcName);  
552   if (fDEBUG) cerr<<"CmpMerged done"<<endl;
553   return kTRUE;
554 }
555 ////////////////////////////////////////////////////////////////////////
556 Bool_t AliMergeSteer::CmpSignalOnly()
557 {
558 //
559 //
560 //
561   char macroName[200];
562   char funcName[200]; 
563   cerr<<"Start CmpSignalOnly"<<endl;
564   sprintf(macroName,"CmpGaRS.C");
565   sprintf(funcName,
566           "CmpGaRS(%d,%d,\"%s\",\"AliTPCtracks.root\",\"%s\");",
567           fNEvents, fFirstEvent, fFileNameHits.Data(),
568           fFileNameCmpSignalOnly.Data());
569   gROOT->LoadMacro(macroName);
570   if (fDEBUG) cerr<<"I'll do: "<<funcName<<endl;
571   gInterpreter->ProcessLine(funcName);  
572   if (fDEBUG) cerr<<"CmpSignalOnly done"<<endl;
573   return kTRUE;
574 }
575 ////////////////////////////////////////////////////////////////////////
576 Bool_t AliMergeSteer::AliCopy(TFile *inputFile, TFile *outputFile)
577 {
578 //
579 // copy gAlice object, AliceGeom and TreeE
580 //
581
582 // copy gAlice
583   if (fDEBUG) cout<<"Copy gAlice: ";
584   outputFile->cd();
585   if (gAlice) {
586     Error("AliCopy",
587           "gAlice must be deleted before AliCopy is called.");
588     return kFALSE;
589   }
590   if (!ImportgAlice(inputFile)) return kFALSE;
591   gAlice->Write();
592   if (fDEBUG) cout<<"done"<<endl;
593
594
595 // copy TreeE
596   TTree *treeE  = gAlice->TreeE();
597   if (!treeE) {
598     cerr<<"No TreeE found "<<endl;
599     return kFALSE;
600   }      
601   if (fDEBUG) cout<<"Copy TreeE: ";
602   AliHeader *header = new AliHeader();
603   treeE->SetBranchAddress("Header", &header);
604   treeE->SetBranchStatus("*",1);
605   TTree *treeENew =  treeE->CloneTree();
606   treeENew->Write();
607   if (fDEBUG) cout<<"done"<<endl;
608
609   delete gAlice;
610   gAlice = 0;
611
612   return kTRUE;
613
614 }
615
616 ////////////////////////////////////////////////////////////////////////
617 Bool_t AliMergeSteer::AliCopy(TString inputFileName, TString outputFileName)
618 {
619 //
620 // open iput and output files, 
621 // ask to copy gAlice object, AliceGeom and TreeE
622 // close input and ouput files
623 //
624   if (fDEBUG) {
625     cout<<"AliCopy: will copy gAlice from "<<inputFileName.Data()<<" to "
626         <<outputFileName.Data()<<endl;
627   }
628
629   TFile *inputFile = TFile::Open(inputFileName.Data());
630   if (!inputFile->IsOpen()) {
631     cerr<<"Can't open "<<inputFileName.Data()<<" !\n";
632     return kFALSE;
633   }
634
635   TFile *outputFile = TFile::Open(outputFileName.Data(),"UPDATE");
636   if (!outputFile->IsOpen()) {
637     cerr<<"Can't open "<<outputFileName.Data()<<" !\n";
638     return kFALSE;
639   }
640
641   AliCopy(inputFile, outputFile);
642
643   inputFile->Close();
644   outputFile->Close();
645
646   if (fDEBUG) {
647     cout<<"AliCopy copied gAlice from "<<inputFileName.Data()<<" to "
648         <<outputFileName.Data()<<endl;
649   }
650
651   return kTRUE;
652 }
653 ////////////////////////////////////////////////////////////////////////