]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ITS/AliITSDetTypeSim.cxx
Modified treatment of SDD gain in simulation: single anode gain normalized to average...
[u/mrichter/AliRoot.git] / ITS / AliITSDetTypeSim.cxx
... / ...
CommitLineData
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/////////////////////////////////////////////////////////////////////
22// Base simulation functions for ITS //
23// //
24// //
25/////////////////////////////////////////////////////////////////////
26#include "TBranch.h"
27#include "TClonesArray.h"
28#include "TObjArray.h"
29#include "TTree.h"
30
31#include "AliRun.h"
32
33#include "AliCDBManager.h"
34#include "AliCDBId.h"
35#include "AliCDBStorage.h"
36#include "AliCDBEntry.h"
37#include "AliCDBMetaData.h"
38#include "AliITSdigit.h"
39#include "AliITSdigitSPD.h"
40#include "AliITSdigitSDD.h"
41#include "AliITSdigitSSD.h"
42#include "AliITSgeom.h"
43#include "AliITSDetTypeSim.h"
44#include "AliITSpListItem.h"
45#include "AliITSCalibration.h"
46#include "AliITSCalibrationSDD.h"
47#include "AliITSMapSDD.h"
48#include "AliITSCorrMapSDD.h"
49#include "AliITSDriftSpeedArraySDD.h"
50#include "AliITSDriftSpeedSDD.h"
51#include "AliITSCalibrationSSD.h"
52#include "AliITSNoiseSSDv2.h"
53#include "AliITSGainSSDv2.h"
54#include "AliITSBadChannelsSSDv2.h"
55#include "AliITSNoiseSSD.h"
56#include "AliITSGainSSD.h"
57#include "AliITSBadChannelsSSD.h"
58#include "AliITSCalibrationSSD.h"
59#include "AliITSsegmentationSPD.h"
60#include "AliITSsegmentationSDD.h"
61#include "AliITSsegmentationSSD.h"
62#include "AliITSsimulation.h"
63#include "AliITSsimulationSPD.h"
64#include "AliITSsimulationSDD.h"
65#include "AliITSsimulationSSD.h"
66#include "AliITSDDLModuleMapSDD.h"
67#include "AliITSTriggerConditions.h"
68#include "AliBaseLoader.h"
69
70const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
71const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD = 240;
72const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD = 260;
73const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
74
75ClassImp(AliITSDetTypeSim)
76
77//----------------------------------------------------------------------
78AliITSDetTypeSim::AliITSDetTypeSim():
79TObject(),
80fSimulation(), // [NDet]
81fSegmentation(), // [NDet]
82fCalibration(), // [NMod]
83fSSDCalibration(0),
84fSPDNoisy(0),
85fNSDigits(0), //! number of SDigits
86fSDigits("AliITSpListItem",1000),
87fNDigits(0), //! number of Digits
88fRunNumber(0), //! Run number (to access DB)
89fDigits(), //! [NMod][NDigits]
90fSimuPar(0),
91fDDLMapSDD(0),
92fAveGainSDD(0),
93fkDigClassName(), // String with digit class name.
94fLoader(0), // local pointer to loader
95fFirstcall(kTRUE),
96fFOGenerator(),
97fTriggerConditions(NULL)
98{
99 // Default Constructor
100 // Inputs:
101 // none.
102 // Outputs:
103 // none.
104 // Return:
105 // A properly zero-ed AliITSDetTypeSim class.
106
107 fSimulation = new TObjArray(fgkNdettypes);
108 fSegmentation = new TObjArray(fgkNdettypes);
109 fSegmentation->SetOwner(kTRUE);
110 fDigits = new TObjArray(fgkNdettypes);
111 fNDigits = new Int_t[fgkNdettypes];
112 fDDLMapSDD=new AliITSDDLModuleMapSDD();
113 fSimuPar= new AliITSSimuParam();
114 fSSDCalibration=new AliITSCalibrationSSD();
115 fNMod[0] = fgkDefaultNModulesSPD;
116 fNMod[1] = fgkDefaultNModulesSDD;
117 fNMod[2] = fgkDefaultNModulesSSD;
118 SetRunNumber();
119}
120//----------------------------------------------------------------------
121AliITSDetTypeSim::~AliITSDetTypeSim(){
122 // Destructor
123 // Inputs:
124 // none.
125 // Outputs:
126 // none.
127 // Return:
128 // Nothing.
129
130 if(fSimulation){
131 fSimulation->Delete();
132 delete fSimulation;
133 }
134 fSimulation = 0;
135 if(fSegmentation){
136 fSegmentation->Delete();
137 delete fSegmentation;
138 }
139 fSegmentation = 0;
140 if(fCalibration && fRunNumber<0){
141
142 fCalibration->Delete();
143 delete fCalibration;
144 }
145 fCalibration = 0;
146 if(fSSDCalibration) delete fSSDCalibration;
147 if(fSPDNoisy){
148 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
149 fSPDNoisy->Delete();
150 delete fSPDNoisy;
151 fSPDNoisy = 0;
152 }
153 }
154 if(fSimuPar) delete fSimuPar;
155 if(fDDLMapSDD) delete fDDLMapSDD;
156 if(fNDigits) delete [] fNDigits;
157 fNDigits = 0;
158 if (fLoader)fLoader->GetModulesFolder()->Remove(this);
159 fLoader = 0; // Not deleting it.
160 fSDigits.Delete();
161 if (fDigits) {
162 fDigits->Delete();
163 delete fDigits;
164 }
165 fDigits=0;
166}
167//----------------------------------------------------------------------
168AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
169fSimulation(source.fSimulation), // [NDet]
170fSegmentation(source.fSegmentation), // [NDet]
171fCalibration(source.fCalibration), // [NMod]
172fSSDCalibration(source.fSSDCalibration),
173fSPDNoisy(source.fSPDNoisy),
174fNSDigits(source.fNSDigits), //! number of SDigits
175fSDigits(*((TClonesArray*)source.fSDigits.Clone())),
176fNDigits(source.fNDigits), //! number of Digits
177fRunNumber(source.fRunNumber), //! Run number (to access DB)
178fDigits(source.fDigits), //! [NMod][NDigits]
179fSimuPar(source.fSimuPar),
180fDDLMapSDD(source.fDDLMapSDD),
181fAveGainSDD(source.fAveGainSDD),
182fkDigClassName(), // String with digit class name.
183fLoader(source.fLoader), // local pointer to loader
184fFirstcall(source.fFirstcall),
185fFOGenerator(source.fFOGenerator),
186fTriggerConditions(source.fTriggerConditions)
187{
188 // Copy Constructor for object AliITSDetTypeSim not allowed
189 for(Int_t i=0;i<fgkNdettypes;i++){
190 fkDigClassName[i] = source.fkDigClassName[i];
191 }
192}
193//----------------------------------------------------------------------
194AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
195 // The = operator for object AliITSDetTypeSim
196
197 this->~AliITSDetTypeSim();
198 new(this) AliITSDetTypeSim(source);
199 return *this;
200}
201
202//______________________________________________________________________
203void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
204 // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
205 //
206 // Inputs:
207 // AliITSgoem *geom The AliITSgeom object to be used.
208 // Output:
209 // none.
210 // Return:
211 // none.
212 if(!fLoader){
213 Error("SetITSgeom","No pointer to loader - nothing done");
214 return;
215 }
216 else {
217 fLoader->SetITSgeom(geom); // protections in AliITSLoader::SetITSgeom
218 }
219
220}
221//______________________________________________________________________
222void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
223 // Sets the local copy of the AliITSLoader, and passes on the
224 // AliITSgeom object as needed.
225 // Inputs
226 // AliITSLoader *loader pointer to AliITSLoader for local use
227 // Outputs:
228 // none.
229 // Return:
230 // none.
231
232 if(fLoader==loader) return; // Same do nothing
233 if(fLoader){ // alread have an existing loader
234 Error("SetLoader",
235 "Already have an exisiting loader ptr=%p Nothing done",
236 fLoader);
237 } // end if
238 fLoader = loader;
239}
240//______________________________________________________________________
241void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
242
243 //Set simulation model for detector type
244
245 if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
246 fSimulation->AddAt(sim,dettype);
247}
248//______________________________________________________________________
249AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype) const {
250
251 //Get simulation model for detector type
252 if(fSimulation==0) {
253 Warning("GetSimulationModel","fSimulation is 0!");
254 return 0;
255 }
256 return (AliITSsimulation*)(fSimulation->At(dettype));
257}
258//______________________________________________________________________
259AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module) const {
260
261 //Get simulation model by module number
262 if(GetITSgeom()==0) {
263 Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
264 return 0;
265 }
266
267 return GetSimulationModel(GetITSgeom()->GetModuleType(module));
268}
269//_______________________________________________________________________
270void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
271 // Set default segmentation model objects
272 AliITSsegmentation *seg;
273
274 if(fSegmentation==0x0){
275 fSegmentation = new TObjArray(fgkNdettypes);
276 fSegmentation->SetOwner(kTRUE);
277 }
278 if(GetSegmentationModel(idet))
279 delete (AliITSsegmentation*)fSegmentation->At(idet);
280 if(idet==0){
281 seg = new AliITSsegmentationSPD();
282 }else if(idet==1){
283 seg = new AliITSsegmentationSDD();
284 AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
285 if(cal->IsAMAt20MHz()){
286 seg->SetPadSize(seg->Dpz(0),20.);
287 seg->SetNPads(seg->Npz()/2,128);
288 }
289 }else {
290 seg = new AliITSsegmentationSSD();
291 }
292 SetSegmentationModel(idet,seg);
293}
294//______________________________________________________________________
295void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
296 AliITSsegmentation *seg){
297
298 //Set segmentation model for detector type
299 if(fSegmentation==0x0){
300 fSegmentation = new TObjArray(fgkNdettypes);
301 fSegmentation->SetOwner(kTRUE);
302 }
303 fSegmentation->AddAt(seg,dettype);
304}
305//______________________________________________________________________
306AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype) const{
307 //Get segmentation model for detector type
308
309 if(fSegmentation==0) {
310 Warning("GetSegmentationModel","fSegmentation is 0!");
311 return 0;
312 }
313 return (AliITSsegmentation*)(fSegmentation->At(dettype));
314}
315//_______________________________________________________________________
316AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module) const{
317 //Get segmentation model by module number
318 if(GetITSgeom()==0){
319 Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
320 return 0;
321 }
322 return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
323}
324//_______________________________________________________________________
325void AliITSDetTypeSim::CreateCalibrationArray() {
326 //Create the container of calibration functions with correct size
327 if (fCalibration) {
328 Warning("CreateCalibration","pointer to calibration object exists\n");
329 fCalibration->Delete();
330 delete fCalibration;
331 }
332
333 Int_t nModTot = GetITSgeom()->GetIndexMax();
334 fCalibration = new TObjArray(nModTot);
335 fCalibration->SetOwner(kTRUE);
336 fCalibration->Clear();
337}
338//_______________________________________________________________________
339void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
340 //Set response model for modules
341
342 if (fCalibration==0) CreateCalibrationArray();
343
344 if (fCalibration->At(iMod)!=0)
345 delete (AliITSCalibration*) fCalibration->At(iMod);
346 fCalibration->AddAt(resp, iMod);
347}
348//_______________________________________________________________________
349void AliITSDetTypeSim::SetSPDNoisyModel(Int_t iMod, AliITSCalibration *cal){
350 //Set noisy pixel info for the SPD module iMod
351 if (fSPDNoisy==0) {
352 fSPDNoisy = new TObjArray(fgkDefaultNModulesSPD);
353 fSPDNoisy->SetOwner(kTRUE);
354 fSPDNoisy->Clear();
355 }
356
357 if (fSPDNoisy->At(iMod) != 0)
358 delete (AliITSCalibration*) fSPDNoisy->At(iMod);
359 fSPDNoisy->AddAt(cal,iMod);
360}
361//______________________________________________________________________
362void AliITSDetTypeSim::ResetCalibrationArray(){
363 //resets response array
364 if(fCalibration && fRunNumber<0){ // if fRunNumber<0 fCalibration is owner
365 /*
366 AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
367 GetITSgeom()->GetStartSPD()))->GetResponse();
368 AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
369 GetITSgeom()->GetStartSSD()))->GetResponse();
370 if(rspd) delete rspd;
371 if(rssd) delete rssd;
372 */
373 fCalibration->Clear();
374 }else if (fCalibration && fRunNumber>=0){
375 fCalibration->Clear();
376 }
377}
378//______________________________________________________________________
379void AliITSDetTypeSim::ResetSegmentation(){
380 //Resets segmentation array
381 if(fSegmentation) fSegmentation->Clear();
382}
383//_______________________________________________________________________
384AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod) const {
385 //Get response model for module number iMod
386
387 if(fCalibration==0) {
388 AliError("fCalibration is 0!");
389 return 0;
390 }
391 if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
392 return (AliITSCalibration*)fCalibration->At(iMod);
393 }else{
394 Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
395 fSSDCalibration->SetModule(i);
396 return (AliITSCalibration*)fSSDCalibration;
397 }
398
399}
400//_______________________________________________________________________
401AliITSCalibration* AliITSDetTypeSim::GetSPDNoisyModel(Int_t iMod) const {
402 //Get SPD noisy calib for module iMod
403 if(fSPDNoisy==0) {
404 AliWarning("fSPDNoisy is 0!");
405 return 0;
406 }
407 return (AliITSCalibration*)fSPDNoisy->At(iMod);
408}
409//_______________________________________________________________________
410void AliITSDetTypeSim::SetDefaults(){
411 //Set defaults for segmentation and response
412
413 if(GetITSgeom()==0){
414 Warning("SetDefaults","GetITSgeom() is 0!");
415 return;
416 } // end if
417 if (fCalibration==0) {
418 CreateCalibrationArray();
419 } // end if
420
421 ResetSegmentation();
422 if(!GetCalibration()){AliFatal("Exit"); exit(0);}
423
424 SetDigitClassName(0,"AliITSdigitSPD");
425 SetDigitClassName(1,"AliITSdigitSDD");
426 SetDigitClassName(2,"AliITSdigitSSD");
427
428 for(Int_t idet=0;idet<fgkNdettypes;idet++){
429 if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
430 }
431}
432//______________________________________________________________________
433Bool_t AliITSDetTypeSim::GetCalibration() {
434 // Get Default calibration if a storage is not defined.
435
436 if(!fFirstcall){
437 AliITSCalibration* cal = GetCalibrationModel(0);
438 if(cal)return kTRUE;
439 }
440 else {
441 fFirstcall = kFALSE;
442 }
443
444 SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
445 Int_t run=GetRunNumber();
446
447 Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
448 Bool_t isCacheActive = kTRUE;
449 if(GetRunNumber()<0){
450 isCacheActive=kFALSE;
451 fCalibration->SetOwner(kTRUE);
452 }
453 else{
454 fCalibration->SetOwner(kFALSE);
455 }
456
457 AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
458
459 AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead", run);
460 AliCDBEntry *noisySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy", run);
461 AliCDBEntry *foEffSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDFOEfficiency", run);
462 AliCDBEntry *foNoiSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDFONoise", run);
463 AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
464 AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD",run);
465 AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD",run);
466 //AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
467 AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
468 // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
469 AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
470 AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
471 AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
472
473if(!deadSPD || !noisySPD || !foEffSPD || !foNoiSPD
474 || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD
475 || !drSpSDD || !ddlMapSDD || !mapTSDD){
476 AliFatal("Calibration object retrieval failed! ");
477 return kFALSE;
478 }
479
480
481 TObjArray *calDeadSPD = (TObjArray*) deadSPD->GetObject();
482 if (!isCacheActive) deadSPD->SetObject(NULL);
483 deadSPD->SetOwner(kTRUE);
484
485 TObjArray *calNoisySPD = (TObjArray*) noisySPD->GetObject();
486 if (!isCacheActive) noisySPD->SetObject(NULL);
487 noisySPD->SetOwner(kTRUE);
488
489 AliITSFOEfficiencySPD *calFoEffSPD = (AliITSFOEfficiencySPD*) foEffSPD->GetObject();
490 if (!isCacheActive) foEffSPD->SetObject(NULL);
491 foEffSPD->SetOwner(kTRUE);
492
493 AliITSFONoiseSPD *calFoNoiSPD = (AliITSFONoiseSPD*) foNoiSPD->GetObject();
494 if (!isCacheActive) foNoiSPD->SetObject(NULL);
495 foNoiSPD->SetOwner(kTRUE);
496
497 TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
498 if(!isCacheActive)entrySDD->SetObject(NULL);
499 entrySDD->SetOwner(kTRUE);
500
501 TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
502 if(!isCacheActive)drSpSDD->SetObject(NULL);
503 drSpSDD->SetOwner(kTRUE);
504
505 AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
506 if(!isCacheActive)ddlMapSDD->SetObject(NULL);
507 ddlMapSDD->SetOwner(kTRUE);
508
509// TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
510// if(!isCacheActive)mapASDD->SetObject(NULL);
511// mapASDD->SetOwner(kTRUE);
512
513 TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
514 if(!isCacheActive)mapTSDD->SetObject(NULL);
515 mapTSDD->SetOwner(kTRUE);
516
517 /*
518 TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
519 if(!isCacheActive)entrySSD->SetObject(NULL);
520 entrySSD->SetOwner(kTRUE);
521 */
522
523 TObject *emptyssd = 0; TString ssdobjectname;
524 AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();
525 emptyssd = (TObject *)entryNoiseSSD->GetObject();
526 ssdobjectname = emptyssd->GetName();
527 if(ssdobjectname=="TObjArray") {
528 TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
529 ReadOldSSDNoise(noiseSSDOld, noiseSSD);
530 }
531 else if(ssdobjectname=="AliITSNoiseSSDv2")
532 noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
533 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
534 entryNoiseSSD->SetOwner(kTRUE);
535
536 AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
537 emptyssd = (TObject *)entryGainSSD->GetObject();
538 ssdobjectname = emptyssd->GetName();
539 if(ssdobjectname=="Gain") {
540 TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
541 ReadOldSSDGain(gainSSDOld, gainSSD);
542 }
543 else if(ssdobjectname=="AliITSGainSSDv2")
544 gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
545 if(!isCacheActive)entryGainSSD->SetObject(NULL);
546 entryGainSSD->SetOwner(kTRUE);
547
548 AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
549 emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
550 ssdobjectname = emptyssd->GetName();
551 if(ssdobjectname=="TObjArray") {
552 TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
553 ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
554 }
555 else if(ssdobjectname=="AliITSBadChannelsSSDv2")
556 badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
557 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
558 entryBadChannelsSSD->SetOwner(kTRUE);
559
560 /*AliITSNoiseSSDv2 *noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
561 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
562 entryNoiseSSD->SetOwner(kTRUE);
563
564 AliITSGainSSDv2 *gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
565 if(!isCacheActive)entryGainSSD->SetObject(NULL);
566 entryGainSSD->SetOwner(kTRUE);
567
568 AliITSBadChannelsSSDv2 *badchannelsSSD =
569 (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
570 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
571 entryBadChannelsSSD->SetOwner(kTRUE);*/
572
573 // DB entries are deleted. In this way metadeta objects are deleted as well
574 if(!isCacheActive){
575 delete deadSPD;
576 delete noisySPD;
577 delete foEffSPD;
578 delete foNoiSPD;
579 delete entrySDD;
580 delete entryNoiseSSD;
581 delete entryGainSSD;
582 delete entryBadChannelsSSD;
583// delete mapASDD;
584 delete mapTSDD;
585 delete drSpSDD;
586 delete ddlMapSDD;
587 }
588
589 AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
590
591 if ((!calDeadSPD) || (!calNoisySPD) || (!calFoEffSPD) || (!calFoNoiSPD)
592 || (!calSDD) || (!drSp) || (!ddlsdd)
593 || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
594 AliWarning("Can not get calibration from calibration database !");
595 return kFALSE;
596 }
597
598
599 fNMod[0] = calDeadSPD->GetEntries();
600 fNMod[1] = calSDD->GetEntries();
601 // fNMod[2] = noiseSSD->GetEntries();
602 AliDebug(1,Form("%i SPD, %i SDD and %i SSD in calibration database",
603 fNMod[0], fNMod[1], fNMod[2]));
604 AliITSCalibration* cal;
605 for (Int_t i=0; i<fNMod[0]; i++) {
606 cal = (AliITSCalibration*) calDeadSPD->At(i);
607 SetCalibrationModel(i, cal);
608 cal = (AliITSCalibration*) calNoisySPD->At(i);
609 SetSPDNoisyModel(i, cal);
610 }
611
612 fFOGenerator.SetEfficiency(calFoEffSPD); // this cal object is used only by the generator
613 fFOGenerator.SetNoise(calFoNoiSPD); // this cal object is used only by the generator
614
615 fDDLMapSDD->SetDDLMap(ddlsdd);
616 Float_t avegain=0.;
617 Float_t nGdAnodes=0;
618 Bool_t oldMapFormat=kFALSE;
619 TObject* objmap=(TObject*)mapT->At(0);
620 TString cname(objmap->ClassName());
621 if(cname.CompareTo("AliITSMapSDD")==0){
622 oldMapFormat=kTRUE;
623 AliInfo("SDD Maps converted to new format");
624 }
625
626 for (Int_t i=0; i<fgkDefaultNModulesSDD; i++) {
627 Int_t iddl,icarlos;
628 Int_t iMod = i + fgkDefaultNModulesSPD;
629 fDDLMapSDD->FindInDDLMap(iMod,iddl,icarlos);
630 if(iddl<0){
631 AliITSCalibrationSDD* calsdddead=new AliITSCalibrationSDD();
632 calsdddead->SetBad();
633 AliITSDriftSpeedSDD* driftspdef = new AliITSDriftSpeedSDD();
634 AliITSDriftSpeedArraySDD* arrdrsp=new AliITSDriftSpeedArraySDD(1);
635 arrdrsp->AddDriftSpeed(driftspdef);
636 calsdddead->SetDriftSpeed(0,arrdrsp);
637 calsdddead->SetDriftSpeed(1,arrdrsp);
638 SetCalibrationModel(iMod, calsdddead);
639 AliWarning(Form("SDD module %d not present in DDL map: set it as dead",iMod));
640 }else{
641 cal = (AliITSCalibration*) calSDD->At(i);
642 for(Int_t iAnode=0;iAnode< ((AliITSCalibrationSDD*)cal)->NOfAnodes(); iAnode++){
643 if(((AliITSCalibrationSDD*)cal)->IsBadChannel(iAnode)) continue;
644 avegain+= ((AliITSCalibrationSDD*)cal)->GetChannelGain(iAnode);
645 nGdAnodes++;
646 }
647 Int_t i0=2*i;
648 Int_t i1=1+2*i;
649 AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
650 AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
651
652 AliITSCorrMapSDD* mt0 = 0;
653 AliITSCorrMapSDD* mt1 = 0;
654 if(oldMapFormat){
655 AliITSMapSDD* oldmap0=(AliITSMapSDD*)mapT->At(i0);
656 AliITSMapSDD* oldmap1=(AliITSMapSDD*)mapT->At(i1);
657 mt0=oldmap0->ConvertToNewFormat();
658 mt1=oldmap1->ConvertToNewFormat();
659 }else{
660 mt0=(AliITSCorrMapSDD*)mapT->At(i0);
661 mt1=(AliITSCorrMapSDD*)mapT->At(i1);
662 }
663 cal->SetDriftSpeed(0,arr0);
664 cal->SetDriftSpeed(1,arr1);
665 cal->SetMapT(0,mt0);
666 cal->SetMapT(1,mt1);
667 SetCalibrationModel(iMod, cal);
668 }
669 }
670 if(nGdAnodes) fAveGainSDD=avegain/nGdAnodes;
671 AliDebug(3,Form("SDD average gain=%f\n",fAveGainSDD));
672 fSSDCalibration->SetNoise(noiseSSD);
673 fSSDCalibration->SetGain(gainSSD);
674 fSSDCalibration->SetBadChannels(badChannelsSSD);
675 //fSSDCalibration->FillBadChipMap();
676
677
678
679 return kTRUE;
680}
681//_______________________________________________________________________
682void AliITSDetTypeSim::SetDefaultSimulation(){
683 //Set default simulation for detector type
684
685 if(GetITSgeom()==0){
686 Warning("SetDefaultSimulation","GetITSgeom() is 0!");
687 return;
688 }
689 if(fCalibration==0){
690 Warning("SetDefaultSimulation","fCalibration is 0!");
691 return;
692 }
693 if(fSegmentation==0){
694 Warning("SetDefaultSimulation","fSegmentation is 0!");
695 for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
696 }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
697 Warning("SetDefaultSimulation",
698 "Segmentation not defined for det %d - Default taken\n!",i);
699 SetDefaultSegmentation(i);
700 }
701 AliITSsimulation* sim;
702
703 for(Int_t idet=0;idet<fgkNdettypes;idet++){
704 //SPD
705 if(idet==0){
706 sim = GetSimulationModel(idet);
707 if(!sim){
708 sim = new AliITSsimulationSPD(this);
709 SetSimulationModel(idet,sim);
710 }
711 }
712 //SDD
713 if(idet==1){
714 sim = GetSimulationModel(idet);
715 if(!sim){
716 sim = new AliITSsimulationSDD(this);
717 SetSimulationModel(idet,sim);
718 }
719 }
720 //SSD
721 if(idet==2){
722 sim = GetSimulationModel(idet);
723 if(!sim){
724 sim = new AliITSsimulationSSD(this);
725 SetSimulationModel(idet,sim);
726 }
727 }
728
729 }
730}
731//___________________________________________________________________
732void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, const Char_t* name){
733 // Set branch address for the ITS summable digits Trees.
734 char branchname[30];
735
736 if(!treeS){
737 return;
738 }
739 TBranch *branch;
740 sprintf(branchname,"%s",name);
741 branch = treeS->GetBranch(branchname);
742 TClonesArray *sdigi = &fSDigits;
743 if (branch) branch->SetAddress(&sdigi);
744
745}
746//___________________________________________________________________
747void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, const Char_t* name){
748 // Set branch address for the digit Trees.
749
750 const char *det[3] = {"SPD","SDD","SSD"};
751 TBranch *branch;
752
753 char branchname[30];
754
755 if(!treeD){
756 return;
757 }
758 if(!fDigits){
759 fDigits = new TObjArray(fgkNdettypes);
760 }
761 for(Int_t i=0;i<fgkNdettypes;i++){
762 const Char_t* digclass = GetDigitClassName(i);
763 if(digclass==0x0){
764 if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
765 if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
766 if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
767 digclass = GetDigitClassName(i);
768 }
769 TString classn = digclass;
770 if(!(fDigits->At(i))){
771 fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
772 }else{
773 ResetDigits(i);
774 }
775
776 if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
777 else sprintf(branchname,"%sDigits%d",name,i+1);
778 if(fDigits){
779 branch = treeD->GetBranch(branchname);
780 if(branch) branch->SetAddress(&((*fDigits)[i]));
781 }
782 }
783
784}
785//___________________________________________________________________
786void AliITSDetTypeSim::ResetDigits(){
787 // Reset number of digits and the digits array for the ITS detector.
788
789 if(!fDigits){
790 Error("ResetDigits","fDigits is null!");
791 return;
792 }
793 for(Int_t i=0;i<fgkNdettypes;i++){
794 ResetDigits(i);
795 }
796}
797//___________________________________________________________________
798void AliITSDetTypeSim::ResetDigits(Int_t branch){
799 // Reset number of digits and the digits array for this branch.
800
801 if(fDigits->At(branch)){
802 ((TClonesArray*)fDigits->At(branch))->Clear();
803 }
804 if(fNDigits) fNDigits[branch]=0;
805
806}
807
808
809
810//_______________________________________________________________________
811void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
812 // Standard Summable digits to Digits function.
813 if(!GetITSgeom()){
814 Warning("SDigitsToDigits","GetITSgeom() is null!!");
815 return;
816 }
817
818 const char *all = strstr(opt,"All");
819 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
820 strstr(opt,"SSD")};
821 if( !det[0] && !det[1] && !det[2] ) all = "All";
822 else all = 0;
823 static Bool_t setDef = kTRUE;
824 if(setDef) SetDefaultSimulation();
825 setDef = kFALSE;
826
827 AliITSsimulation *sim =0;
828 TTree* trees = fLoader->TreeS();
829 if( !(trees && GetSDigits()) ){
830 Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
831 return;
832 }
833 sprintf(name,"%s",name);
834 TBranch* brchSDigits = trees->GetBranch(name);
835
836 Int_t id;
837 for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
838 id = GetITSgeom()->GetModuleType(module);
839 if (!all && !det[id]) continue;
840 sim = (AliITSsimulation*)GetSimulationModel(id);
841 if(!sim){
842 Error("SDigit2Digits","The simulation class was not "
843 "instanciated for module %d type %s!",module,
844 GetITSgeom()->GetModuleTypeName(module));
845 exit(1);
846 }
847 sim->InitSimulationModule(module,gAlice->GetEvNumber());
848
849 fSDigits.Clear();
850 brchSDigits->GetEvent(module);
851 sim->AddSDigitsToModule(&fSDigits,0);
852 sim->FinishSDigitiseModule();
853 fLoader->TreeD()->Fill();
854 ResetDigits();
855 }
856 fLoader->TreeD()->GetEntries();
857 fLoader->TreeD()->AutoSave();
858 fLoader->TreeD()->Reset();
859}
860//_________________________________________________________
861void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){
862 //Adds the module full of summable digits to the summable digits tree.
863
864 new(fSDigits[fNSDigits++]) AliITSpListItem(sdig);
865}
866//__________________________________________________________
867void AliITSDetTypeSim::AddSimDigit(Int_t branch, const AliITSdigit* d){
868 // Add a simulated digit.
869
870 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
871 switch(branch){
872 case 0:
873 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
874 break;
875 case 1:
876 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
877 break;
878 case 2:
879 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
880 break;
881 }
882}
883//______________________________________________________________________
884void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
885 Int_t *tracks,Int_t *hits,Float_t *charges,
886 Int_t sigexpanded){
887 // Add a simulated digit to the list.
888
889 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
890 switch(branch){
891 case 0:
892 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
893 break;
894 case 1:
895 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
896 hits,charges,sigexpanded);
897 break;
898 case 2:
899 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
900 break;
901 }
902}
903//______________________________________________________________________
904void AliITSDetTypeSim::ReadOldSSDNoise(const TObjArray *array,
905 AliITSNoiseSSDv2 *noiseSSD) {
906 //Reads the old SSD calibration object and converts it to the new format
907 const Int_t fgkSSDSTRIPSPERMODULE = 1536;
908 const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
909
910 Int_t gNMod = array->GetEntries();
911 cout<<"Converting old calibration object for noise..."<<endl;
912
913 //NOISE
914 Double_t noise = 0.0;
915 for (Int_t iModule = 0; iModule < gNMod; iModule++) {
916 AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
917 for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
918 noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
919 if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
920 noiseSSD->AddNoiseP(iModule,iStrip,noise);
921 if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
922 noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
923 }//loop over strips
924 }//loop over modules
925}
926
927//______________________________________________________________________
928void AliITSDetTypeSim::ReadOldSSDBadChannels(const TObjArray *array,
929 AliITSBadChannelsSSDv2 *badChannelsSSD) {
930 //Reads the old SSD calibration object and converts it to the new format
931 Int_t nMod = array->GetEntries();
932 cout<<"Converting old calibration object for bad channels..."<<endl;
933 for (Int_t iModule = 0; iModule < nMod; iModule++) {
934 //for (Int_t iModule = 0; iModule < 1; iModule++) {
935 AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
936 TArrayI arrayPSide = bad->GetBadPChannelsList();
937 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
938 badChannelsSSD->AddBadChannelP(iModule,
939 iPCounter,
940 (Char_t)arrayPSide.At(iPCounter));
941
942 TArrayI arrayNSide = bad->GetBadNChannelsList();
943 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
944 badChannelsSSD->AddBadChannelN(iModule,
945 iNCounter,
946 (Char_t)arrayNSide.At(iNCounter));
947
948 }//loop over modules
949}
950
951//______________________________________________________________________
952void AliITSDetTypeSim::ReadOldSSDGain(const TObjArray *array,
953 AliITSGainSSDv2 *gainSSD) {
954 //Reads the old SSD calibration object and converts it to the new format
955
956 Int_t nMod = array->GetEntries();
957 cout<<"Converting old calibration object for gain..."<<endl;
958
959 //GAIN
960 for (Int_t iModule = 0; iModule < nMod; iModule++) {
961 AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
962 TArrayF arrayPSide = gainModule->GetGainP();
963 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
964 gainSSD->AddGainP(iModule,
965 iPCounter,
966 arrayPSide.At(iPCounter));
967 TArrayF arrayNSide = gainModule->GetGainN();
968 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
969 gainSSD->AddGainN(iModule,
970 iNCounter,
971 arrayNSide.At(iNCounter));
972 }//loop over modules
973}
974//______________________________________________________________________
975void AliITSDetTypeSim::ProcessSPDDigitForFastOr(UInt_t module, UInt_t colM, UInt_t rowM) {
976 // Processes wether a single fired pixel will give rise to a fast-or signal
977 fFOGenerator.ProcessPixelHitM(module,colM,rowM);
978}
979//_______________________________________________________________________
980AliITSTriggerConditions* AliITSDetTypeSim::GetTriggerConditions() {
981 // Get Pixel Trigger Conditions (separate method since it is used only when simulating trigger)
982 if (fTriggerConditions==NULL) { // read from db
983 fRunNumber = ((Int_t)AliCDBManager::Instance()->GetRun());
984 Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
985 Bool_t isCacheActive;
986 if (fRunNumber<0) isCacheActive=kFALSE;
987 else isCacheActive=kTRUE;
988 AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
989 AliCDBEntry *pitCond = AliCDBManager::Instance()->Get("TRIGGER/SPD/PITConditions", fRunNumber);
990 if (!pitCond) {
991 AliError("Trigger conditions retrieval failed! ");
992 return NULL;
993 }
994 fTriggerConditions = (AliITSTriggerConditions*) pitCond->GetObject();
995 if (!isCacheActive) pitCond->SetObject(NULL);
996 pitCond->SetOwner(kTRUE);
997 if (!isCacheActive) {
998 delete pitCond;
999 }
1000 AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
1001 if (fTriggerConditions==NULL) {
1002 AliWarning("fTriggerConditions is NULL!");
1003 }
1004 }
1005 return fTriggerConditions;
1006}
1007//_______________________________________________________________________
1008void AliITSDetTypeSim::WriteFOSignals() {
1009 // write fo signals to event
1010
1011 if (!fLoader) {
1012 AliError("ITS loader is NULL.");
1013 return;
1014 }
1015
1016 if(!fLoader->TreeD()){
1017 AliError("No TreeD available");
1018 return;
1019 }
1020 TTree *tree = fLoader->TreeD();
1021 AliITSFOSignalsSPD *foSignals = new AliITSFOSignalsSPD(*GetFOSignals());
1022 tree->GetUserInfo()->Add(foSignals);
1023}
1024