1 //_____________________________________________________________________________
2 Int_t AliTRDdisplayDigits3D(Int_t event = 0, Int_t thresh = 2
3 , Bool_t sdigits = kFALSE)
9 // <event> : Event number
10 // <thresh> : Threshold to suppress the noise
11 // <sdigits> : If kTRUE it will display summable digits, normal digits otherwise.
12 // The signal event is displayed in yellow.
15 Char_t *inputFile = "galice.root";
17 const Int_t kNdict = 3;
26 TString evfoldname = AliConfig::GetDefaultEventFolderName();
27 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
29 fRunLoader = AliRunLoader::Open(inputFile
30 ,AliConfig::GetDefaultEventFolderName()
34 Printf("Can not open session for file %s.",inputFile);
38 if (!fRunLoader->GetAliRun()) {
39 fRunLoader->LoadgAlice();
41 gAlice = fRunLoader->GetAliRun();
43 printf("Could not find AliRun object.\n");
47 fRunLoader->GetEvent(event);
49 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
51 printf("Can not get TRD loader from Run Loader");
55 // Get the pointer to the detector object
56 trd = (AliTRDv1*) gAlice->GetDetector("TRD");
58 // Get the pointer to the geometry object
60 geo = trd->GetGeometry();
63 printf("Cannot find the geometry\n");
67 AliCDBManager *cdbManager = AliCDBManager::Instance();
68 cdbManager->SetDefaultStorage("local://$ALICE_ROOT");
69 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
70 calibration->SetRun(0);
72 TCanvas *c1 = new TCanvas("digits","TRD digits display",0,0,700,730);
73 TView *v = new TView(1);
74 v->SetRange(-430,-560,-430,430,560,1710);
80 Int_t markerColorSignal = 2;
81 Int_t markerColorBgnd = 7;
82 Int_t markerColorMerged = 5;
83 Int_t mask = 10000000;
85 // Create the digits manager
86 AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
87 digitsManager->SetSDigits(sdigits);
89 // Read the digits from the file
91 digitsManager->ReadDigits(loader->TreeS());
94 if (!loader->TreeD()) {
98 digitsManager->ReadDigits(loader->TreeD());
101 Int_t totalsignal = 0;
103 Int_t totalmerged = 0;
104 Int_t timeMax = calibration->GetNumberOfTimeBins();
106 // Loop through all detectors
107 for (Int_t idet = 0; idet < geo->Ndet(); idet++) {
109 printf("<AliTRDdisplayDigits3D> Loading detector %d\n",idet);
110 AliTRDdataArrayI *digits = digitsManager->GetDigits(idet);
112 AliTRDdataArrayI *tracks[kNdict];
113 for (Int_t idict = 0; idict < kNdict; idict++) {
114 tracks[idict] = digitsManager->GetDictionary(idet,idict);
115 tracks[idict]->Expand();
118 Int_t isec = geo->GetSector(idet);
119 Int_t icha = geo->GetChamber(idet);
120 Int_t ipla = geo->GetPlane(idet);
121 AliTRDpadPlane *padPlane = new AliTRDpadPlane(ipla,icha);
122 Int_t rowMax = padPlane->GetNrows();
123 Int_t colMax = padPlane->GetNcols();
125 Int_t ndigits = digits->GetOverThreshold(thresh);
129 TPolyMarker3D *pmSignal = new TPolyMarker3D(ndigits);
130 TPolyMarker3D *pmBgnd = new TPolyMarker3D(ndigits);
131 TPolyMarker3D *pmMerged = new TPolyMarker3D(ndigits);
137 for (Int_t time = 0; time < timeMax; time++) {
138 for (Int_t col = 0; col < colMax; col++) {
139 for (Int_t row = 0; row < rowMax; row++) {
143 Int_t amp = digits->GetDataUnchecked(row,col,time);
144 for (idict = 0; idict < kNdict; idict++) {
145 Int_t trk = tracks[idict]->GetDataUnchecked(row,col,time) - 1;
146 if ((idict == 0) && (trk >= mask)) {
149 if ((type == 1) && (trk >= mask)) {
162 geo->Local2Global(idet,loc,glb);
168 pmSignal->SetPoint(isignal,x,y,z);
172 else if (type == 2) {
173 pmBgnd->SetPoint(ibgnd,x,y,z);
177 else if (type == 3) {
178 pmMerged->SetPoint(imerged,x,y,z);
189 digits->Compress(1,0);
190 for (idict = 0; idict < kNdict; idict++) {
191 tracks[idict]->Compress(1,0);
194 pmMerged->SetMarkerSize(1);
195 pmMerged->SetMarkerColor(markerColorMerged);
196 pmMerged->SetMarkerStyle(1);
199 pmBgnd->SetMarkerSize(1);
200 pmBgnd->SetMarkerColor(markerColorBgnd);
201 pmBgnd->SetMarkerStyle(1);
204 pmSignal->SetMarkerSize(1);
205 pmSignal->SetMarkerColor(markerColorSignal);
206 pmSignal->SetMarkerStyle(1);
215 TGeometry *geoAlice = gAlice->GetGeometry();
216 TNode *main = (TNode *) ((geoAlice->GetListOfNodes())->First());
217 TIter next(main->GetListOfNodes());
219 while ((module = (TNode *) next())) {
221 sprintf(ch,"%s\n",module->GetTitle());
222 if ((ch[0] == 'T') && ((ch[1] == 'R') || (ch[1] == 'P'))) {
223 module->SetVisibility( 3);
226 module->SetVisibility(-1);
230 geoAlice->Draw("same");
235 printf("<AliTRDdisplayDigits3D> Number of digits:\n");
236 printf(" signal = %d, bgnd = %d, merged = %d\n"
237 ,totalsignal,totalbgnd,totalmerged);