]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
Removed.
[u/mrichter/AliRoot.git] / ITS / AliITSDetTypeSim.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17  $Id$
18 */
19
20 /////////////////////////////////////////////////////////////////////
21 // Base simulation functions for ITS                               //
22 //                                                                 //
23 //                                                                 //
24 /////////////////////////////////////////////////////////////////////          
25 #include "TBranch.h"
26 #include "TClonesArray.h"
27 #include "TObjArray.h"
28 #include "TTree.h"
29
30 #include "AliRun.h"
31 #include "AliITSdigit.h"
32 #include "AliITSdigitSPD.h"
33 #include "AliITSdigitSDD.h"
34 #include "AliITSdigitSSD.h"
35 #include "AliITSDetTypeSim.h"
36 #include "AliITSgeom.h"
37 #include "AliITSpListItem.h"
38 #include "AliITSresponseSPD.h"
39 #include "AliITSresponseSDD.h"
40 #include "AliITSresponseSSD.h"
41 #include "AliITSsegmentationSPD.h"
42 #include "AliITSsegmentationSDD.h"
43 #include "AliITSsegmentationSSD.h"
44 #include "AliITSsimulation.h"
45 #include "AliITSsimulationSPD.h"
46 #include "AliITSsimulationSDD.h"
47 #include "AliITSsimulationSSD.h"
48
49
50 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
51
52 ClassImp(AliITSDetTypeSim)
53
54 //----------------------------------------------------------------------
55 AliITSDetTypeSim::AliITSDetTypeSim():
56 TObject(),
57 fGeom(),         //
58 fSimulation(),   // [NDet]
59 fSegmentation(), // [NDet]
60 fResponse(),     // [NMod]
61 fPreProcess(),   // [] e.g. Fill fHitModule with hits
62 fPostProcess(),  // [] e.g. Wright Raw data
63 fNSDigits(0),    //! number of SDigits
64 fSDigits(),      //! [NMod][NSDigits]
65 fNDigits(0),     //! number of Digits
66 fDigits(),       //! [NMod][NDigits]
67 fHitClassName(), // String with Hit class name.
68 fSDigClassName(),// String with SDigit class name.
69 fDigClassName(){ // String with digit class name.
70     // Default Constructor
71     // Inputs:
72     //    none.
73     // Outputs:
74     //    none.
75     // Return:
76     //    A properly zero-ed AliITSDetTypeSim class.
77   fGeom = 0;
78   fSimulation = new TObjArray(fgkNdettypes);
79   fSegmentation = new TObjArray(fgkNdettypes);
80   fResponse = 0;
81   fPreProcess = 0;
82   fPostProcess = 0;
83   fNSDigits = 0;
84   fSDigits = new TClonesArray("AliITSpListItem",1000);
85   fDigits = new TObjArray(fgkNdettypes);
86   fNDigits = new Int_t[fgkNdettypes];
87   fLoader = 0;
88 }
89 //----------------------------------------------------------------------
90 AliITSDetTypeSim::~AliITSDetTypeSim(){
91     // Destructor
92     // Inputs:
93     //    none.
94     // Outputs:
95     //    none.
96     // Return:
97     //    Nothing.
98   
99     if(fGeom) delete fGeom;
100     if(fSimulation){
101       fSimulation->Delete();
102       delete fSimulation;
103       fSimulation = 0;
104     }
105     
106     if(fSegmentation){
107       fSegmentation->Delete();
108       delete fSegmentation;
109       fSegmentation = 0;
110     }
111     
112     if(fResponse){
113       fResponse->Delete();
114       delete fResponse;
115       fResponse = 0;
116     }
117     
118     if(fPreProcess){
119       fPreProcess->Delete();
120       delete fPreProcess;
121       fPreProcess = 0;
122     }
123     
124     if(fPostProcess){
125       fPostProcess->Delete();
126       delete fPostProcess;
127       fPostProcess = 0;
128     }
129     
130     if (fLoader)
131       {
132         fLoader->GetModulesFolder()->Remove(this);
133       }
134     
135         
136     if (fSDigits) {
137       fSDigits->Delete();
138       delete fSDigits;
139       fSDigits=0;
140     }
141     if (fDigits) {
142       fDigits->Delete();
143       delete fDigits;
144       fDigits=0;
145     }
146   
147 }
148 //----------------------------------------------------------------------
149 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source){
150     // Copy Constructor for object AliITSDetTypeSim not allowed
151   if(this==&source) return;
152   Error("Copy constructor",
153         "You are not allowed to make a copy of the AliITSDetTypeSim");
154   exit(1);
155
156         
157 }
158 //----------------------------------------------------------------------
159 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
160     // The = operator for object AliITSDetTypeSim
161  
162     if(&source==this) return *this;
163     Error("operator=","You are not allowed to make a copy of the AliITSDetTypeSIm");
164     exit(1);
165     return *this;
166 }
167
168
169 //______________________________________________________________________
170 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
171
172   //Set simulation model for detector type
173
174   if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
175   fSimulation->AddAt(sim,dettype);
176 }
177 //______________________________________________________________________
178 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype){
179
180   //Get simulation model for detector type
181   if(fSimulation==0)  {
182     Warning("GetSimulationModel","fSimulation is 0!");
183     return 0;     
184   }
185   return (AliITSsimulation*)(fSimulation->At(dettype));
186 }
187 //______________________________________________________________________
188 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module){
189
190   //Get simulation model by module number
191   if(fGeom==0) {
192     Warning("GetSimulationModelByModule","fGeom is 0!");
193     return 0;
194   }
195   
196   return GetSimulationModel(fGeom->GetModuleType(module));
197 }
198 //______________________________________________________________________
199 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
200    
201   //Set segmentation model for detector type
202   if(fSegmentation==0x0) fSegmentation = new TObjArray(fgkNdettypes);
203   fSegmentation->AddAt(seg,dettype);
204
205 }
206 //______________________________________________________________________
207 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype){
208
209   //Get segmentation model for detector type
210    
211    if(fSegmentation==0) {
212      Warning("GetSegmentationModel","fSegmentation is 0!");
213      return 0; 
214    } 
215    return (AliITSsegmentation*)(fSegmentation->At(dettype));
216
217 }
218 //_______________________________________________________________________
219 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module){
220   
221   //Get segmentation model by module number
222    if(fGeom==0){
223      Warning("GetSegmentationModelByModule","fGeom is 0!");
224      return 0;
225    }     
226    return GetSegmentationModel(fGeom->GetModuleType(module));
227
228 }
229 //_______________________________________________________________________
230 void AliITSDetTypeSim::SetResponseModel(Int_t dettype,AliITSresponse *resp){
231
232   
233   //Set segmentation model for module number
234   if(fResponse==0) fResponse = new TObjArray(fgkNdettypes);
235   fResponse->AddAt(resp,dettype);
236 }
237 //______________________________________________________________________
238 void AliITSDetTypeSim::ResetResponse(){
239
240   //resets response array
241   if(fResponse){
242     for(Int_t i=0;i<fResponse->GetEntries();i++){
243       if(fResponse->At(i)) delete (AliITSresponse*)fResponse->At(i);
244     }
245   }
246 }
247 //______________________________________________________________________
248 void AliITSDetTypeSim::ResetSegmentation(){
249  
250  //Resets segmentation array
251   if(fSegmentation){
252     for(Int_t i=0;i<fgkNdettypes;i++){
253       if(fSegmentation->At(i)) delete (AliITSsegmentation*)fSegmentation->At(i);
254     }
255   }
256 }
257
258 //_______________________________________________________________________
259 AliITSresponse* AliITSDetTypeSim::GetResponseModel(Int_t dettype){
260   
261   //Get segmentation model for module number
262   
263   if(fResponse==0) {
264     Warning("GetResponseModel","fResponse is 0!");
265     return 0; 
266   }  
267   return (AliITSresponse*)(fResponse->At(dettype));
268 }
269 //_______________________________________________________________________
270 void AliITSDetTypeSim::SetDefaults(){
271
272   //Set defaults for segmentation and response
273   
274
275   if(fGeom==0){
276     Warning("SetDefaults","fGeom is 0!");
277     return;
278   }
279
280   if(!fResponse) fResponse = new TObjArray(fgkNdettypes);  
281
282   AliITSsegmentation* seg;
283   AliITSresponse* res;
284   ResetResponse();
285   ResetSegmentation();
286
287   for(Int_t idet=0;idet<fgkNdettypes;idet++){
288     //SPD
289     if(idet==0){
290       if(!GetSegmentationModel(idet)){
291         seg = new AliITSsegmentationSPD(fGeom);
292         SetSegmentationModel(idet,seg);
293       }
294       if(!GetResponseModel(idet)){
295         res = new AliITSresponseSPD();
296         SetResponseModel(idet,res);
297       }
298       const char *kData0=(GetResponseModel(idet))->DataType();
299       if (strstr(kData0,"real")) {
300         SetDigitClassName(idet,"AliITSdigit");
301       } else SetDigitClassName(idet,"AliITSdigitSPD");
302     }
303     //SDD
304     if(idet==1){
305       if(!GetResponseModel(idet)){
306         SetResponseModel(idet,new AliITSresponseSDD("simulated"));
307       }
308       if(!GetSegmentationModel(idet)){
309         res = GetResponseModel(idet);
310         seg = new AliITSsegmentationSDD(fGeom,res);
311         SetSegmentationModel(idet,seg);
312       }
313       const char *kopt = GetResponseModel(idet)->ZeroSuppOption();
314       if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D"))) SetDigitClassName(idet,"AliITSdigit");
315       else SetDigitClassName(idet,"AliITSdigitSDD");
316     }
317     //SSD
318     if(idet==2){
319       if(!GetSegmentationModel(idet)){
320         seg = new AliITSsegmentationSSD(fGeom);
321         SetSegmentationModel(idet,seg);
322       }
323       if(!GetResponseModel(idet)){
324         SetResponseModel(idet,new AliITSresponseSSD("simulated"));
325       }
326       const char *kData2 = (GetResponseModel(idet))->DataType();
327       if (strstr(kData2,"real")) {
328         SetDigitClassName(idet,"AliITSdigit");
329       } else SetDigitClassName(idet,"AliITSdigitSSD");
330       
331     }
332
333   }
334
335 }
336
337 //_______________________________________________________________________
338 void AliITSDetTypeSim::SetDefaultSimulation(){
339
340   //Set default simulation for detector type
341
342  
343   if(fGeom==0){
344     Warning("SetDefaults","fGeom is 0!");
345     return;
346   }
347   
348   if(!fResponse) fResponse = new TObjArray(fgkNdettypes);
349
350   AliITSsegmentation* seg;
351   AliITSresponse* res;
352   AliITSsimulation* sim;
353
354   for(Int_t idet=0;idet<fgkNdettypes;idet++){
355    //SPD
356     if(idet==0){
357       sim = GetSimulationModel(idet);
358       if(!sim){
359         seg = (AliITSsegmentationSPD*)GetSegmentationModel(idet);
360         res = (AliITSresponseSPD*)GetResponseModel(idet);      
361         sim = new AliITSsimulationSPD(seg,res);
362         SetSimulationModel(idet,sim);
363       } else{
364         sim->SetResponseModel(GetResponseModel(idet));
365         sim->SetSegmentationModel((AliITSsegmentationSPD*)GetSegmentationModel(idet));
366         sim->Init();
367       }
368     }
369     //SDD
370     if(idet==1){
371       sim = GetSimulationModel(idet);
372       if(!sim){
373         seg = (AliITSsegmentationSDD*)GetSegmentationModel(idet);
374         res = (AliITSresponseSDD*)GetResponseModel(idet);
375         sim = new AliITSsimulationSDD(seg,res);
376         SetSimulationModel(idet,sim);
377       } else {
378         sim->SetResponseModel((AliITSresponseSDD*)GetResponseModel(idet));
379         sim->SetSegmentationModel((AliITSsegmentationSDD*)GetSegmentationModel(idet));
380         sim->Init();
381       }
382       
383     }
384     //SSD
385     if(idet==2){
386       sim = GetSimulationModel(idet);
387       if(!sim){
388         seg = (AliITSsegmentationSSD*)GetSegmentationModel(idet);
389         res = (AliITSresponseSSD*)GetResponseModel(idet);
390         sim = new AliITSsimulationSSD(seg,res);
391         SetSimulationModel(idet,sim);
392       } else{
393         sim->SetResponseModel((AliITSresponseSSD*)GetResponseModel(idet));
394         sim->SetSegmentationModel((AliITSsegmentationSSD*)GetSegmentationModel(idet));
395         sim->Init();
396       }
397
398     }
399
400   }
401 }
402
403 //___________________________________________________________________
404 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
405   // Set branch address for the ITS summable digits Trees.
406   
407   char branchname[30];
408
409   if(!treeS){
410     return;
411   }
412   if (fSDigits == 0x0){
413     fSDigits = new TClonesArray("AliITSpListItem",1000);
414   }
415   TBranch *branch;
416   sprintf(branchname,"%s",name);
417   branch = treeS->GetBranch(branchname);
418   if (branch) branch->SetAddress(&fSDigits);
419
420 }
421 //___________________________________________________________________
422 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
423   // Set branch address for the digit Trees.
424   
425   const char *det[3] = {"SPD","SDD","SSD"};
426   TBranch *branch;
427   
428   char branchname[30];
429   
430   if(!treeD){
431     return;
432   }
433   if(!fDigits){
434     fDigits = new TObjArray(fgkNdettypes); 
435   }
436   for(Int_t i=0;i<fgkNdettypes;i++){
437     Char_t* digclass = GetDigitClassName(i);
438     if(digclass==0x0){
439       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
440       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
441       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
442       digclass = GetDigitClassName(i);
443     }
444     TString classn = digclass;
445     if(!(fDigits->At(i))){
446       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
447     }else{
448       ResetDigits(i);
449     }
450     
451     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
452     else sprintf(branchname,"%sDigits%d",name,i+1);
453     if(fDigits){
454       branch = treeD->GetBranch(branchname);
455       if(branch) branch->SetAddress(&((*fDigits)[i]));
456     }
457   }
458
459 }
460 //___________________________________________________________________
461 void AliITSDetTypeSim::ResetDigits(){
462   // Reset number of digits and the digits array for the ITS detector.
463   
464
465   if(!fDigits){
466     Error("ResetDigits","fDigits is null!");
467     return;
468   }
469   for(Int_t i=0;i<fgkNdettypes;i++){
470     ResetDigits(i);
471   }
472 }
473 //___________________________________________________________________
474 void AliITSDetTypeSim::ResetDigits(Int_t branch){
475   // Reset number of digits and the digits array for this branch.
476
477   if(fDigits->At(branch)){
478     ((TClonesArray*)fDigits->At(branch))->Clear();
479   }
480   if(fNDigits) fNDigits[branch]=0;
481
482 }
483
484
485
486 //_______________________________________________________________________
487 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
488   // Standard Summable digits to Digits function.
489   if(!fGeom){
490     Warning("SDigitsToDigits","fGeom is null!!");
491     return;
492   }
493   
494   const char *all = strstr(opt,"All");
495   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
496                         strstr(opt,"SSD")};
497   if( !det[0] && !det[1] && !det[2] ) all = "All";
498   else all = 0;
499   static Bool_t setDef = kTRUE;
500   if(setDef) SetDefaultSimulation();
501   setDef = kFALSE;
502   
503   AliITSsimulation *sim =0;
504   TTree* trees = fLoader->TreeS();
505   if( !(trees && GetSDigits()) ){
506     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
507     return;
508   } 
509   sprintf(name,"%s",name);
510   TBranch* brchSDigits = trees->GetBranch(name);
511   
512   Int_t id;
513   for(Int_t module=0;module<fGeom->GetIndexMax();module++){
514      id = fGeom->GetModuleType(module);
515     if (!all && !det[id]) continue;
516     sim = (AliITSsimulation*)GetSimulationModel(id);
517     if(!sim){
518       Error("SDigit2Digits","The simulation class was not "
519             "instanciated for module %d type %s!",module,
520             fGeom->GetModuleTypeName(module));
521       exit(1);
522     }
523     sim->InitSimulationModule(module,gAlice->GetEvNumber());
524     
525     fSDigits->Clear();
526     brchSDigits->GetEvent(module);
527     sim->AddSDigitsToModule(fSDigits,0);
528     sim->FinishSDigitiseModule();
529     fLoader->TreeD()->Fill();
530     ResetDigits();
531   }
532   fLoader->TreeD()->GetEntries();
533   fLoader->TreeD()->AutoSave();
534   fLoader->TreeD()->Reset();
535 }
536
537
538
539 //_________________________________________________________
540 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){
541   
542   //Adds the module full of summable digits to the summable digits tree.
543   TClonesArray &lsdig = *fSDigits;
544   new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
545 }
546 //__________________________________________________________
547 void AliITSDetTypeSim::AddRealDigit(Int_t branch, Int_t *digits){
548   //   Add a real digit - as coming from data.
549   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
550   new(ldigits[fNDigits[branch]++]) AliITSdigit(digits); 
551 }
552 //__________________________________________________________
553 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){
554   
555   //    Add a simulated digit.
556   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
557   switch(branch){
558   case 0:
559     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
560     break;
561   case 1:
562     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
563     break;
564   case 2:
565     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
566     break;
567   } 
568   
569
570 }
571
572 //______________________________________________________________________
573 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
574                                    Int_t *tracks,Int_t *hits,Float_t *charges){
575   //   Add a simulated digit to the list.
576
577   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
578   AliITSresponseSDD *resp = 0;
579   switch(branch){
580   case 0:
581     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
582     break;
583   case 1:
584     resp = (AliITSresponseSDD*)GetResponseModel(branch);
585     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
586                                                    hits,charges,resp);
587     break;
588   case 2:
589     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
590     break;
591   } 
592 }