]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSHits2Digits.C
Bug in geometry corrected
[u/mrichter/AliRoot.git] / ITS / AliITSHits2Digits.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3 #include "iostream.h"
4 #include "TDatetime.h"
5 #include "STEER/AliRun.h"
6 #include "STEER/AliRunDigitizer.h"
7 #include "ITS/AliITSDigitizer.h"
8 #include "ITS/AliITS.h"
9 #include "ITS/AliITSDetType.h"
10 #include "ITS/AliITSresponseSDD.h"
11 #include "TStopwatch.h"
12
13 #endif
14
15 TFile* AccessFile(TString inFile="galice.root", TString acctype="R");
16 void writeAR(TFile * fin, TFile *fou);
17 Int_t ChangeITSDefaults(TFile *hitfile,AliITS *ITS,TString opt="");
18 //#define DEBUG
19 Int_t AliITSHits2Digits(TString hitFile = "galice.root", 
20                         TString digFile = "galiceD.root",TString opt=""){
21     // Standeard ITS Hits to Digits, excluding creation of SDigits.
22
23     // Dynamically link some shared libs
24     if (gClassTable->GetID("AliRun") < 0) {
25         gROOT->LoadMacro("loadlibs.C");
26         loadlibs();
27     } // end if
28
29     // Connect the Root Galice file containing Geometry, Kine and Hits
30
31     TFile *hitfile = 0;   // pointer to input file.
32     TFile *digfile = 0;  // possible output file for TreeD
33     if(digFile.CompareTo(hitFile) == 0){// write output to same file as input.
34         hitfile = AccessFile(hitFile,"U");  // input file open for update.
35     }else{ // different output file then input file.
36         hitfile = AccessFile(hitFile,"R"); // input file open as read only
37         // open output file and create TreeR on it
38         digfile = gAlice->InitTreeFile("D",digFile);
39     } // end if digFile == hitFile.
40
41     AliITS *ITS = (AliITS*)gAlice->GetDetector("ITS");      
42     if (!ITS) {
43         cerr<<"AliITSHits2DigitsDefault.C : AliITS object not found on file"
44             << endl;
45         return 3;
46     }  // end if !ITS
47     if(!(ITS->GetITSgeom())){
48         cerr << " AliITSgeom not found. Can't digitize with out it." << endl;
49         return 4;
50     } // end if
51
52     ChangeITSDefaults(hitfile,ITS,opt);
53     // write the AliRun object to the output file if different from input file.
54     if(digfile) writeAR(hitfile,digfile);
55
56     TStopwatch timer;
57     Int_t evNumber1 = 0;
58     Int_t evNumber2 = gAlice->GetEventsPerRun();
59     timer.Start();
60     for(Int_t nevent = evNumber1; nevent < evNumber2; nevent++){
61         // cout<<"Producing Digits for event n."<<nevent<<endl;
62         gAlice->GetEvent(nevent);
63         if(!gAlice->TreeD() && digfile == 0){ 
64             cout << "Having to create the Digits Tree." << endl;
65             gAlice->MakeTree("D");
66         } // end if creating digits tree
67         if(digfile) gAlice->MakeTree("D",digfile);
68         ITS->MakeBranch("D");
69         ITS->SetTreeAddress();   
70         ITS->Hits2Digits();
71     } // end for nevent
72     timer.Stop();
73     timer.Print();
74     if(digfile!=0){
75         cout << digFile << " size =" << digfile->GetSize() << endl;
76     }else{
77         cout << hitFile << " size =" << hitfile->GetSize() << endl;
78     } // end if sdigfile!=0
79
80     delete gAlice; // digfile is closed by deleting gAlice if != hitfile.
81     gAlice = 0;
82     hitfile->Close(); hitfile = 0;
83 }
84 //______________________________________________________________________
85 TFile * AccessFile(TString FileName, TString acctype){
86     // Function used to open the input file and fetch the AliRun object
87
88     TFile *retfil = 0;
89     TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(FileName);
90     if(file) {
91         file->Close();
92         delete file;
93         file = 0;
94     } // end if file
95     if(acctype.Contains("U")){
96         file = new TFile(FileName,"UPDATE");
97     } // end if open for update
98     if(acctype.Contains("N") && !file){
99         file = new TFile(FileName,"RECREATE");
100     } // end if open a new file
101     if(!file) file = new TFile(FileName,"READ");   // default readonly
102     if (!file->IsOpen()) {
103         cerr << "Can't open " << FileName << " !" << endl;
104         return retfil;
105     } // end if error opeing file
106
107     // Get AliRun object from file or return if not on file
108     if (gAlice) {delete gAlice; gAlice = 0;}
109     gAlice = (AliRun*)file->Get("gAlice");
110     if (!gAlice) {
111         cerr << "AliRun object not found on file "<< FileName << "!" << endl;
112         file->Close();  // close file and return error.
113         return retfil;
114     } // end if !gAlice
115     return file;
116 }
117 //______________________________________________________________________
118 void writeAR(TFile * fin, TFile *fou) {
119     TDirectory *current = gDirectory;
120     TTree *TeOld;
121     TTree *TeNew;
122     AliHeader *alhe = new AliHeader();
123     TeOld = (TTree*)fin->Get("TE");
124     TeOld->SetBranchAddress("Header",&alhe);
125     TeOld->SetBranchStatus("*",1);
126     fou->cd();
127     TeNew = TeOld->CloneTree();
128     TeNew->Write(0,TObject::kOverwrite);
129     gAlice->Write(0,TObject::kOverwrite);
130     current->cd();
131     delete alhe;
132 #ifdef DEBUG
133     cout << "AliRun object written to file" << endl;
134 #endif
135 }
136 //______________________________________________________________________
137 Int_t ChangeITSDefaults(TFile *hitfile,AliITS *ITS,TString opt){
138
139     TDatime *ct0 = new TDatime(2002,04,26,00,00,00);
140     TDatime ct = hitfile->GetCreationDate();
141
142     if(ct0->GetDate()>ct.GetDate()){
143         // For old files, must change SDD noise.
144         AliITSresponseSDD *resp1 = (AliITSresponseSDD*)ITS->DetType(1)->
145             GetResponseModel();
146         resp1 = new AliITSresponseSDD();
147         ITS->SetResponseModel(1,resp1);
148         cout << "Changed response class for SDD:" << endl;
149         resp1->Print();
150     } // end if
151
152     if(opt.Contains("Dubna")){
153         AliITSresponseSPDdubna *resp0 = new AliITSresponseSPDdubna();
154         if(ITS->DetType(0)->GetResponseModel() !=0){
155             delete ((AliITSresponse*)ITS->DetType(0)->GetResponseModel());
156             ITS->DetType(0)->ResponseModel(0);
157         } // end if
158         ITS->DetType(0)->ResponseModel(resp0);
159         AliITSsegmentationSPD *seg0 = (AliITSsegmentationSPD*)ITS->
160             DetType(0)->GetSegmentationModel();
161         AliITSsimulationSPDdubna *sim0 = new AliITSsimulationSPDdubna(seg0,
162                                                                       resp0);
163         if(ITS->DetType(0)->GetSimulationModel() !=0){
164             delete ((AliITSsimulation*)ITS->DetType(0)->GetSimulationModel());
165             ITS->DetType(0)->SimulationModel(0);
166         } // end if
167         ITS->DetType(0)->SimulationModel(sim0);
168     } // end if Dubna
169 }