]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSHits2SDigits.C
First commit.
[u/mrichter/AliRoot.git] / ITS / AliITSHits2SDigits.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 AliITSHits2SDigits(TString  hitFile = "galice.root",
20                          TString sdigFile = "galice.root",TString opt=""){
21     // Standeard ITS Hits to 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 *sdigfile = 0;  // possible output file for TreeD
33     if(sdigFile.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         sdigfile = gAlice->InitTreeFile("S",sdigFile);
39     } // end if sdigFile == 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(sdigfile) writeAR(hitfile,sdigfile);
55
56     TStopwatch timer;
57     Int_t evNumber1 = 0;
58     Int_t evNumber2 = gAlice->GetEventsPerRun();
59     timer.Start();
60     for(Int_t event = evNumber1; event < evNumber2; event++){
61         gAlice->GetEvent(event);
62         if(!gAlice->TreeS() && sdigfile == 0){ 
63             cout << "Having to create the SDigits Tree." << endl;
64             gAlice->MakeTree("S");
65         } // end if !gAlice->TreeS()
66         if(sdigfile) gAlice->MakeTree("S",sdigfile);
67         //    make branch
68         ITS->MakeBranch("S");
69         ITS->SetTreeAddress();
70 #ifdef DEBUG
71         cout<<"Making ITS SDigits for event "<<event<<endl;
72 #endif
73         ITS->Hits2SDigits();
74     } // end for event
75     timer.Stop();
76     timer.Print();
77     if(sdigfile!=0){
78         cout << sdigFile << " size =" << sdigfile->GetSize() << endl;
79     }else{
80         cout << hitFile << " size =" << hitfile->GetSize() << endl;
81     } // end if sdigfile!=0
82
83     delete gAlice; // sdigfile is closed by deleting gAlice if != hitfile.
84     gAlice = 0;
85     hitfile->Close(); hitfile=0;
86 }
87 //______________________________________________________________________
88 TFile * AccessFile(TString FileName, TString acctype){
89     // Function used to open the input file and fetch the AliRun object
90
91     TFile *retfil = 0;
92     TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(FileName);
93     if(file) {
94         file->Close();
95         delete file;
96         file = 0;
97     } // end if file
98     if(acctype.Contains("U")){
99         file = new TFile(FileName,"UPDATE");
100     } // end if open for update
101     if(acctype.Contains("N") && !file){
102         file = new TFile(FileName,"RECREATE");
103     } // end if open a new file
104     if(!file) file = new TFile(FileName,"READ");   // default readonly
105     if (!file->IsOpen()) {
106         cerr << "Can't open " << FileName << " !" << endl;
107         return retfil;
108     } // end if error opeing file
109
110     // Get AliRun object from file or return if not on file
111     if (gAlice) {delete gAlice; gAlice = 0;}
112     gAlice = (AliRun*)file->Get("gAlice");
113     if (!gAlice) {
114         cerr << "AliRun object not found on file "<< FileName << "!" << endl;
115         file->Close();  // close file and return error.
116         return retfil;
117     } // end if !gAlice
118     return file;
119 }
120 //______________________________________________________________________
121 void writeAR(TFile * fin, TFile *fou) {
122     TDirectory *current = gDirectory;
123     TTree *TeOld;
124     TTree *TeNew;
125     AliHeader *alhe = new AliHeader();
126     TeOld = (TTree*)fin->Get("TE");
127     TeOld->SetBranchAddress("Header",&alhe);
128     TeOld->SetBranchStatus("*",1);
129     fou->cd();
130     TeNew = TeOld->CloneTree();
131     TeNew->Write(0,TObject::kOverwrite);
132     gAlice->Write(0,TObject::kOverwrite);
133     current->cd();
134     delete alhe;
135 #ifdef DEBUG
136     cout << "AliRun object written to file" << endl;
137 #endif
138 }
139 //______________________________________________________________________
140 Int_t ChangeITSDefaults(TFile *hitfile,AliITS *ITS,TString opt){
141
142     TDatime *ct0 = new TDatime(2002,04,26,00,00,00);
143     TDatime ct = hitfile->GetCreationDate();
144
145     if(ct0->GetDate()>ct.GetDate()){
146         // For old files, must change SDD noise.
147         AliITSresponseSDD *resp1 = (AliITSresponseSDD*)ITS->DetType(1)->
148             GetResponseModel();
149         resp1 = new AliITSresponseSDD();
150         ITS->SetResponseModel(1,resp1);
151         cout << "Changed response class for SDD:" << endl;
152         resp1->Print();
153     } // end if
154
155     if(opt.Contains("Dubna")){
156         AliITSresponseSPDdubna *resp0 = new AliITSresponseSPDdubna();
157         if(ITS->DetType(0)->GetResponseModel() !=0){
158             delete ((AliITSresponse*)ITS->DetType(0)->GetResponseModel());
159             ITS->DetType(0)->ResponseModel(0);
160         } // end if
161         ITS->DetType(0)->ResponseModel(resp0);
162         AliITSsegmentationSPD *seg0 = (AliITSsegmentationSPD*)ITS->
163             DetType(0)->GetSegmentationModel();
164         AliITSsimulationSPDdubna *sim0 = new AliITSsimulationSPDdubna(seg0,
165                                                                       resp0);
166         if(ITS->DetType(0)->GetSimulationModel() !=0){
167             delete ((AliITSsimulation*)ITS->DetType(0)->GetSimulationModel());
168             ITS->DetType(0)->SimulationModel(0);
169         } // end if
170         ITS->DetType(0)->SimulationModel(sim0);
171     } // end if Dubna
172 }