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