1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 ////////////////////////////////////////////////////////////////////////
21 // This class defines the "Standard" reconstruction for the ITS //
24 ////////////////////////////////////////////////////////////////////////
25 #include "TObjArray.h"
28 #include "AliITSresponseSDD.h"
29 #include "AliCDBManager.h"
30 #include "AliCDBStorage.h"
31 #include "AliCDBEntry.h"
33 #include "AliITSClusterFinder.h"
34 #include "AliITSClusterFinderV2.h"
35 #include "AliITSClusterFinderV2SPD.h"
36 #include "AliITSClusterFinderV2SDD.h"
37 #include "AliITSClusterFinderV2SSD.h"
38 #include "AliITSClusterFinderSPD.h"
39 #include "AliITSClusterFinderSDD.h"
40 #include "AliITSClusterFinderSSD.h"
41 #include "AliITSclusterV2.h"
42 #include "AliITSgeom.h"
43 #include "AliITSDetTypeRec.h"
44 #include "AliITSRawCluster.h"
45 #include "AliITSRawClusterSPD.h"
46 #include "AliITSRawClusterSDD.h"
47 #include "AliITSRawClusterSSD.h"
48 #include "AliITSRecPoint.h"
49 #include "AliITSCalibrationSDD.h"
50 #include "AliITSsegmentationSPD.h"
51 #include "AliITSsegmentationSDD.h"
52 #include "AliITSsegmentationSSD.h"
55 const Int_t AliITSDetTypeRec::fgkNdettypes = 3;
56 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSPD = 240;
57 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSDD = 260;
58 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSSD = 1698;
60 ClassImp(AliITSDetTypeRec)
62 //________________________________________________________________
63 AliITSDetTypeRec::AliITSDetTypeRec(): TObject(),
65 fReconstruction(),// [NDet]
66 fSegmentation(), // [NDet]
67 fCalibration(), // [NMod]
68 fPreProcess(), // [] e.g. Find Calibration values
69 fPostProcess(), // [] e.g. find primary vertex
70 fDigits(), //! [NMod][NDigits]
71 fClusterClassName(), // String with Cluster class name
72 fDigClassName(), // String with digit class name.
73 fRecPointClassName(){// String with RecPoint class name
74 // Default Constructor
80 // A properly zero-ed AliITSDetTypeRec class.
83 fReconstruction = new TObjArray(fgkNdettypes);
88 fDigits = new TObjArray(fgkNdettypes);
89 fNdtype = new Int_t[fgkNdettypes];
90 fCtype = new TObjArray(fgkNdettypes);
91 fNctype = new Int_t[fgkNdettypes];
92 fRecPoints = new TClonesArray("AliITSRecPoint",1000);
94 fClustersV2 = new TClonesArray("AliITSclusterV2",3000);
97 for(Int_t i=0;i<fgkNdettypes;i++){
105 fNMod[0] = fgkDefaultNModulesSPD;
106 fNMod[1] = fgkDefaultNModulesSDD;
107 fNMod[2] = fgkDefaultNModulesSSD;
111 //______________________________________________________________________
112 AliITSDetTypeRec::AliITSDetTypeRec(const AliITSDetTypeRec &/*rec*/):TObject(/*rec*/){
115 Error("Copy constructor","Copy constructor not allowed");
118 //______________________________________________________________________
119 AliITSDetTypeRec& AliITSDetTypeRec::operator=(const AliITSDetTypeRec& /*source*/){
120 // Assignment operator. This is a function which is not allowed to be
122 Error("operator=","Assignment operator not allowed\n");
126 //_____________________________________________________________________
127 AliITSDetTypeRec::~AliITSDetTypeRec(){
132 fReconstruction->Delete();
133 delete fReconstruction;
137 fSegmentation->Delete();
138 delete fSegmentation;
142 AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(fGeom->GetStartSPD()))->GetResponse();
143 AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(fGeom->GetStartSDD()))->GetResponse();
144 AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(fGeom->GetStartSSD()))->GetResponse();
145 if(rspd) delete rspd;
146 if(rsdd) delete rsdd;
147 if(rssd) delete rssd;
148 fCalibration->Delete();
152 if(fGeom) delete fGeom;
153 if(fPreProcess) delete fPreProcess;
154 if(fPostProcess) delete fPostProcess;
162 fRecPoints->Delete();
167 fClustersV2->Delete();
178 if(fLoader) delete fLoader;
181 //___________________________________________________________________
182 void AliITSDetTypeRec::SetReconstructionModel(Int_t dettype,AliITSClusterFinder *clf){
184 //Set reconstruction model for detector type
186 if(fReconstruction==0) fReconstruction = new TObjArray(fgkNdettypes);
187 if(fReconstruction->At(dettype)!=0) delete fReconstruction->At(dettype);
188 fReconstruction->AddAt(clf,dettype);
190 //______________________________________________________________________
191 AliITSClusterFinder* AliITSDetTypeRec::GetReconstructionModel(Int_t dettype){
193 //Get reconstruction model for detector type
194 if(fReconstruction==0) {
195 Warning("GetReconstructionModel","fReconstruction is 0!");
198 return (AliITSClusterFinder*)fReconstruction->At(dettype);
201 //______________________________________________________________________
202 void AliITSDetTypeRec::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
204 //Set segmentation model for detector type
206 if(fSegmentation==0) fSegmentation = new TObjArray(fgkNdettypes);
207 if(fSegmentation->At(dettype)!=0) delete fSegmentation->At(dettype);
208 fSegmentation->AddAt(seg,dettype);
211 //______________________________________________________________________
212 AliITSsegmentation* AliITSDetTypeRec::GetSegmentationModel(Int_t dettype){
214 //Get segmentation model for detector type
216 if(fSegmentation==0) {
217 Warning("GetSegmentationModel","fSegmentation is 0!");
220 return (AliITSsegmentation*)fSegmentation->At(dettype);
223 //_______________________________________________________________________
224 void AliITSDetTypeRec::SetCalibrationModel(Int_t iMod, AliITSCalibration *cal){
226 //Set calibration (response) for the module iMod of type dettype
227 if (fCalibration==0) {
228 fCalibration = new TObjArray(fGeom->GetIndexMax());
229 fCalibration->SetOwner(kTRUE);
230 fCalibration->Clear();
233 if (fCalibration->At(iMod) != 0)
234 delete (AliITSCalibration*) fCalibration->At(iMod);
235 fCalibration->AddAt(cal,iMod);
238 //_______________________________________________________________________
239 AliITSCalibration* AliITSDetTypeRec::GetCalibrationModel(Int_t iMod){
241 //Get calibration model for module type
243 if(fCalibration==0) {
244 Warning("GetalibrationModel","fCalibration is 0!");
248 return (AliITSCalibration*)fCalibration->At(iMod);
251 //______________________________________________________________________
252 void AliITSDetTypeRec::SetTreeAddress(){
253 // Set branch address for the Trees.
255 TTree *treeR = fLoader->TreeR();
256 TTree *treeD = fLoader->TreeD();
258 SetTreeAddressD(treeD);
259 SetTreeAddressR(treeR);
261 //______________________________________________________________________
262 void AliITSDetTypeRec::SetTreeAddressD(TTree *treeD){
263 // Set branch address for the tree of digits.
265 const char *det[4] = {"SPD","SDD","SSD","ITS"};
272 if (fDigits == 0x0) fDigits = new TObjArray(fgkNdettypes);
273 for (i=0; i<fgkNdettypes; i++) {
274 digclass = GetDigitClassName(i);
275 if(!(fDigits->At(i))) {
276 fDigits->AddAt(new TClonesArray(digclass,1000),i);
280 if (fgkNdettypes==3) sprintf(branchname,"%sDigits%s",det[3],det[i]);
281 else sprintf(branchname,"%sDigits%d",det[3],i+1);
283 branch = treeD->GetBranch(branchname);
284 if (branch) branch->SetAddress(&((*fDigits)[i]));
289 //_______________________________________________________________________
290 TBranch* AliITSDetTypeRec::MakeBranchInTree(TTree *tree, const char* name,
291 const char *classname,
292 void* address,Int_t size,
293 Int_t splitlevel, const char */*file*/)
296 // Makes branch in given tree and diverts them to a separate file
302 Error("MakeBranchInTree","Making Branch %s Tree is NULL",name);
305 TBranch *branch = tree->GetBranch(name);
310 branch = tree->Branch(name,classname,address,size,splitlevel);
313 branch = tree->Bronch(name, "TClonesArray", address, size, splitlevel);
319 //____________________________________________________________________
320 void AliITSDetTypeRec::SetDefaults(){
322 //Set defaults for segmentation and response
325 Warning("SetDefaults","fGeom is 0!");
329 AliITSsegmentation* seg;
330 if(!GetCalibration()) {AliFatal("Exit");exit(0);}
332 for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
334 seg = new AliITSsegmentationSPD(fGeom);
335 SetSegmentationModel(dettype,seg);
336 SetDigitClassName(dettype,"AliITSdigitSPD");
337 SetClusterClassName(dettype,"AliITSRawClusterSPD");
341 AliITSCalibrationSDD* res=(AliITSCalibrationSDD*) GetCalibrationModel(fGeom->GetStartSDD());
342 seg = new AliITSsegmentationSDD(fGeom,res);
343 SetSegmentationModel(dettype,seg);
344 const char *kopt = ((AliITSresponseSDD*)res->GetResponse())->ZeroSuppOption();
345 if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D"))) SetDigitClassName(dettype,"AliITSdigit");
346 else SetDigitClassName(dettype,"AliITSdigitSDD");
347 SetClusterClassName(dettype,"AliITSRawClusterSDD");
351 AliITSsegmentationSSD* seg2 = new AliITSsegmentationSSD(fGeom);
352 seg2->SetAngles(0.0075,0.0275); // strip angels rad P and N side.
353 seg2->SetAnglesLay5(0.0075,0.0275); // strip angels rad P and N side.
354 seg2->SetAnglesLay6(0.0275,0.0075); // strip angels rad P and N side.
355 SetSegmentationModel(dettype,seg2);
356 SetDigitClassName(dettype,"AliITSdigitSSD");
357 SetClusterClassName(dettype,"AliITSRawClusterSSD");
362 //______________________________________________________________________
363 Bool_t AliITSDetTypeRec::GetCalibration() {
364 // Get Default calibration if a storage is not defined.
367 AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSPD", fRunNumber);
368 AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", fRunNumber);
369 AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", fRunNumber);
370 AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", fRunNumber);
371 AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD", fRunNumber);
372 AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", fRunNumber);
374 if(!entrySPD || !entrySDD || !entrySSD || !entry2SPD || !entry2SDD || !entry2SSD){
375 AliWarning("Calibration object retrieval failed! Dummy calibration will be used.");
376 AliCDBStorage *origStorage = AliCDBManager::Instance()->GetDefaultStorage();
377 AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
379 entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSPD", fRunNumber);
380 entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", fRunNumber);
381 entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", fRunNumber);
382 entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", fRunNumber);
383 entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD", fRunNumber);
384 entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", fRunNumber);
386 AliCDBManager::Instance()->SetDefaultStorage(origStorage);
390 TObjArray *respSPD = (TObjArray *)entrySPD->GetObject();
391 entrySPD->SetObject(NULL);
392 entrySPD->SetOwner(kTRUE);
394 AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
395 entry2SPD->SetObject(NULL);
396 entry2SPD->SetOwner(kTRUE);
398 TObjArray *respSDD = (TObjArray *)entrySDD->GetObject();
399 entrySDD->SetObject(NULL);
400 entrySDD->SetOwner(kTRUE);
402 AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
403 entry2SDD->SetObject(NULL);
404 entry2SDD->SetOwner(kTRUE);
406 TObjArray *respSSD = (TObjArray *)entrySSD->GetObject();
407 entrySSD->SetObject(NULL);
408 entrySSD->SetOwner(kTRUE);
410 AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
411 entry2SSD->SetObject(NULL);
412 entry2SSD->SetOwner(kTRUE);
414 // DB entries are dleted. In this waymetadeta objects are deleted as well
423 if ((!pSPD)||(!pSDD)||(!pSSD) || (!respSPD) || (!respSDD) || (!respSSD)) {
424 AliWarning("Can not get calibration from calibration database !");
428 fNMod[0] = respSPD->GetEntries();
429 fNMod[1] = respSDD->GetEntries();
430 fNMod[2] = respSSD->GetEntries();
431 AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
432 fNMod[0], fNMod[1], fNMod[2]));
433 AliITSCalibration* res;
434 for (Int_t i=0; i<fNMod[0]; i++) {
435 res = (AliITSCalibration*) respSPD->At(i);
436 res->SetResponse((AliITSresponse*)pSPD);
437 SetCalibrationModel(i, res);
439 for (Int_t i=0; i<fNMod[1]; i++) {
440 res = (AliITSCalibration*) respSDD->At(i);
441 res->SetResponse((AliITSresponse*)pSDD);
442 Int_t iMod = i + fNMod[0];
443 SetCalibrationModel(iMod, res);
445 for (Int_t i=0; i<fNMod[2]; i++) {
446 res = (AliITSCalibration*) respSSD->At(i);
447 res->SetResponse((AliITSresponse*)pSSD);
448 Int_t iMod = i + fNMod[0] + fNMod[1];
449 SetCalibrationModel(iMod, res);
456 //________________________________________________________________
457 void AliITSDetTypeRec::SetDefaultClusterFinders(){
459 //set defaults for standard cluster finder
462 Warning("SetDefaults","fGeom is 0!");
466 AliITSClusterFinder *clf;
470 for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
473 if(!GetReconstructionModel(dettype)){
474 TClonesArray *dig0 = DigitsAddress(0);
475 TClonesArray *rec0 = ClustersAddress(0);
476 clf = new AliITSClusterFinderSPD(this,dig0,rec0);
477 SetReconstructionModel(dettype,clf);
484 if(!GetReconstructionModel(dettype)){
485 TClonesArray *dig1 = DigitsAddress(1);
486 TClonesArray *rec1 = ClustersAddress(1);
487 clf = new AliITSClusterFinderSDD(this,dig1,rec1);
488 SetReconstructionModel(dettype,clf);
494 if(!GetReconstructionModel(dettype)){
495 TClonesArray* dig2 = DigitsAddress(2);
496 clf = new AliITSClusterFinderSSD(this,dig2);
497 SetReconstructionModel(dettype,clf);
506 //________________________________________________________________
507 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata){
509 //Set defaults for cluster finder V2
512 Warning("SetDefaults","fGeom is 0!");
516 AliITSClusterFinder *clf;
519 for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
522 if(!GetReconstructionModel(dettype)){
523 clf = new AliITSClusterFinderV2SPD(this);
525 if(!rawdata) clf->SetDigits(DigitsAddress(0));
526 SetReconstructionModel(dettype,clf);
532 if(!GetReconstructionModel(dettype)){
533 clf = new AliITSClusterFinderV2SDD(this);
535 if(!rawdata) clf->SetDigits(DigitsAddress(1));
536 SetReconstructionModel(dettype,clf);
543 if(!GetReconstructionModel(dettype)){
544 clf = new AliITSClusterFinderV2SSD(this);
546 if(!rawdata) clf->SetDigits(DigitsAddress(2));
547 SetReconstructionModel(dettype,clf);
554 //______________________________________________________________________
555 void AliITSDetTypeRec::MakeBranch(Option_t* option){
557 //Creates branches for clusters and recpoints
558 Bool_t cR = (strstr(option,"R")!=0);
559 Bool_t cRF = (strstr(option,"RF")!=0);
560 Bool_t v2 = (strstr(option,"v2")!=0);
564 if(cR) MakeBranchR(0);
565 if(cRF) MakeBranchRF(0);
566 if(v2) MakeBranchR(0,"v2");
570 //_____________________________________________________________
571 void AliITSDetTypeRec::MakeTreeC(){
573 //Create a separate tree to store the clusters
575 Warning("MakeTreeC","ITS loader is null!");
578 if(fLoader->TreeC()== 0x0) fLoader->MakeTree("C");
582 //______________________________________________________________
583 void AliITSDetTypeRec::MakeBranchC(){
585 //Make branches in the tree of clusters
588 Warning("MakeBranchC","ITS loader is null!");
591 TTree* lTC = fLoader->TreeC();
593 Error("MakeTreeC","Can not get TreeC from loader");
597 Int_t buffersize = 4000;
598 Char_t branchname[30];
600 const char *det[4]={"SPD","SDD","SSD","ITS"};
602 for(Int_t i=0;i<fgkNdettypes;i++){
603 cluclass = GetClusterClassName(i);
604 if(fCtype==0x0) fCtype = new TObjArray(fgkNdettypes);
605 if(!ClustersAddress(i)){
606 fCtype->AddAt(new TClonesArray(cluclass,1000),i);
608 if(fgkNdettypes==3) sprintf(branchname,"%sClusters%s",det[3],det[i]);
609 else sprintf(branchname,"%sClusters%d",det[3],i+1);
611 if(lTC->GetBranch(branchname)){
612 Warning("MakeBranchC","Branch %s already exists in TreeC",branchname);
614 Info("MakeBranchC","Creating branch %s in TreeC",branchname);
615 lTC->Branch(branchname,&((*fCtype)[i]),buffersize);
623 //_______________________________________________________________
624 void AliITSDetTypeRec::GetTreeC(Int_t event){
626 //Get the clusters tree for this event and
627 //sets the branch address
631 Warning("GetTreeC","ITS loader is null!");
635 Char_t branchname[30];
636 const char *det[4] = {"SPD","SDD","SSD","ITS"};
637 TTree* lTC = fLoader->TreeC();
640 if(lTC) fLoader->CleanRawClusters();
645 for(Int_t i=0;i<fgkNdettypes;i++){
646 cluclass = GetClusterClassName(i);
647 if(fCtype==0x0) fCtype = new TObjArray(fgkNdettypes);
648 if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(cluclass,1000),i);
649 if(fgkNdettypes==3) sprintf(branchname,"%sClusters%s",det[3],det[i]);
650 else sprintf(branchname,"%sClusters%d",det[3],i+1);
652 branch = lTC->GetBranch(branchname);
653 if(branch) branch->SetAddress(&((*fCtype)[i]));
657 Error("GetTreeC","cannot find clusters Tree for vent %d",event);
662 //___________________________________________________________________
663 void AliITSDetTypeRec::AddCluster(Int_t id, AliITSRawCluster *c){
665 // Adds a raw cluster to the list
666 TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
669 new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
672 new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
675 new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
679 //___________________________________________________________________
680 void AliITSDetTypeRec::ResetDigits(){
681 // Reset number of digits and the digits array for the ITS detector.
684 for(Int_t i=0;i<fgkNdettypes;i++){
688 //___________________________________________________________________
689 void AliITSDetTypeRec::ResetDigits(Int_t branch){
690 // Reset number of digits and the digits array for this branch.
692 if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
693 if(fNdtype) fNdtype[branch]=0;
697 //__________________________________________________________________
698 void AliITSDetTypeRec::ResetClusters(){
700 //Resets number of clusters and the cluster array
701 for(Int_t i=0;i<fgkNdettypes;i++){
706 //__________________________________________________________________
707 void AliITSDetTypeRec::ResetClusters(Int_t i){
709 //Resets number of clusters and the cluster array for this branch
711 if (fCtype->At(i)) ((TClonesArray*)fCtype->At(i))->Clear();
712 if (fNctype) fNctype[i]=0;
714 //__________________________________________________________________
715 void AliITSDetTypeRec::MakeBranchR(const char *file, Option_t *opt){
717 //Creates tree branches for recpoints
719 // cont char *file File name where RecPoints branch is to be written
720 // to. If blank it write the SDigits to the same
721 // file in which the Hits were found.
724 Warning("MakeBranchR","ITS loader is null!");
731 // only one branch for rec points for all detector types
732 Bool_t oFast= (strstr(opt,"Fast")!=0);
733 Bool_t v2 = (strstr(opt,"v2")!=0);
735 Char_t detname[10] = "ITS";
739 sprintf(branchname,"%sRecPointsF",detname);
741 sprintf(branchname,"Clusters");
743 sprintf(branchname,"%sRecPoints",detname);
748 if(!fClustersV2)fClustersV2 = new TClonesArray("AliITSclusterV2",3000);
749 if(fLoader->TreeR()){
750 if(fClustersV2==0x0) fClustersV2 = new TClonesArray("AliITSclusterV2",3000);
751 MakeBranchInTree(fLoader->TreeR(),branchname,0,&fClustersV2,buffsz,99,file);
755 if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
756 if (fLoader->TreeR()) {
757 if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",
759 MakeBranchInTree(fLoader->TreeR(),branchname,0,&fRecPoints,buffsz,99,file);
764 //______________________________________________________________________
765 void AliITSDetTypeRec::SetTreeAddressR(TTree *treeR){
766 // Set branch address for the Reconstructed points Trees.
768 // TTree *treeR Tree containing the RecPoints.
774 Char_t namedet[10]="ITS";
777 if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
779 sprintf(branchname,"Clusters");
780 branch1 = treeR->GetBranch(branchname);
782 if(fClustersV2==0x0) fClustersV2 = new TClonesArray("AliITSclusterV2",3000);
783 branch1->SetAddress(&fClustersV2);
787 sprintf(branchname,"%sRecPoints",namedet);
788 branch = treeR->GetBranch(branchname);
790 branch->SetAddress(&fRecPoints);
792 sprintf(branchname,"%sRecPointsF",namedet);
793 branch = treeR->GetBranch(branchname);
795 branch->SetAddress(&fRecPoints);
800 //____________________________________________________________________
801 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
802 // Add a reconstructed space point to the list
804 // const AliITSRecPoint &r RecPoint class to be added to the tree
805 // of reconstructed points TreeR.
811 TClonesArray &lrecp = *fRecPoints;
812 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
814 //______________________________________________________________________
815 void AliITSDetTypeRec::AddClusterV2(const AliITSclusterV2 &r){
816 // Add a reconstructed space point to the list
818 // const AliITSClusterV2 &r class to be added to the tree
819 // of reconstructed points TreeR.
825 TClonesArray &lrecp = *fClustersV2;
826 new(lrecp[fNClustersV2++]) AliITSclusterV2(r);
829 //______________________________________________________________________
830 void AliITSDetTypeRec::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt, Bool_t v2){
831 // cluster finding and reconstruction of space points
832 // the condition below will disappear when the geom class will be
833 // initialized for all versions - for the moment it is only for v5 !
834 // 7 is the SDD beam test version
836 // Int_t evNumber Event number to be processed.
837 // Int_t lastentry Offset for module when not all of the modules
839 // Option_t *opt String indicating which ITS sub-detectors should
840 // be processed. If ="All" then all of the ITS
841 // sub detectors are processed.
844 Warning("DigitsToRecPoints","fGeom is null!");
848 Warning("DigitsToRecPoints","ITS loader is null!");
852 const char *all = strstr(opt,"All");
853 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
855 if(!v2) SetDefaultClusterFinders();
856 else SetDefaultClusterFindersV2();
859 TTree *treeC=fLoader->TreeC();
864 AliITSClusterFinder *rec = 0;
865 Int_t id,module,first=0;
866 for(module=0;module<fGeom->GetIndexMax();module++){
867 id = fGeom->GetModuleType(module);
868 if (!all && !det[id]) continue;
869 if(det[id]) first = fGeom->GetStartDet(id);
870 rec = (AliITSClusterFinder*)GetReconstructionModel(id);
871 TClonesArray *itsDigits = DigitsAddress(id);
873 Error("DigitsToRecPoints",
874 "The reconstruction class was not instanciated! event=%d",
879 TTree *lTD = fLoader->TreeD();
881 lTD->GetEvent(lastentry+module);
883 lTD->GetEvent(lastentry+(module-first));
885 Int_t ndigits = itsDigits->GetEntriesFast();
887 rec->SetDetTypeRec(this);
888 rec->SetDigits(DigitsAddress(id));
889 rec->SetClusters(ClustersAddress(id));
890 rec->FindRawClusters(module);
892 fLoader->TreeR()->Fill();
899 fLoader->WriteRecPoints("OVERWRITE");
900 fLoader->WriteRawClusters("OVERWRITE");
902 //______________________________________________________________________
903 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader){
904 // cluster finding and reconstruction of space points
905 // the condition below will disappear when the geom class will be
906 // initialized for all versions - for the moment it is only for v5 !
907 // 7 is the SDD beam test version
909 // Int_t evNumber Event number to be processed.
910 // Int_t lastentry Offset for module when not all of the modules
912 // Option_t *opt String indicating which ITS sub-detectors should
913 // be processed. If ="All" then all of the ITS
914 // sub detectors are processed.
920 Warning("DigitsToRecPoints","fGeom is null!");
924 Warning("DigitsToRecPoints","ITS loader is null!");
929 AliITSClusterFinderV2 *rec = 0;
932 if(!fLoader->TreeR()) fLoader->MakeTree("R");
933 TTree* cTree = fLoader->TreeR();
934 TClonesArray *array=new TClonesArray("AliITSclusterV2",1000);
935 cTree->Branch("Clusters",&array);
938 TClonesArray** clusters = new TClonesArray*[fGeom->GetIndexMax()];
939 for (Int_t iModule = 0; iModule < fGeom->GetIndexMax(); iModule++) {
940 clusters[iModule] = NULL;
943 rec = (AliITSClusterFinderV2*)GetReconstructionModel(id);
944 rec->SetDetTypeRec(this);
946 Error("DigitsToRecPoints",
947 "The reconstruction class was not instanciated");
950 rec->RawdataToClusters(rawReader,clusters);
953 for(Int_t iModule=0;iModule<fGeom->GetIndexMax();iModule++){
954 array = clusters[iModule];
956 Error("DigitsToRecPoints","data for module %d missing!",iModule);
957 array = new TClonesArray("AliITSclusterV2");
959 cTree->SetBranchAddress("Clusters",&array);
961 nClusters+=array->GetEntriesFast();
964 fLoader->WriteRecPoints("OVERWRITE");
967 Info("DigitsToRecPoints", "total number of found clustersV2 in ITS: %d\n",