]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/AliMpSectorSegmentation.cxx
Fixing part of the Coding violation
[u/mrichter/AliRoot.git] / MUON / mapping / AliMpSectorSegmentation.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 // $Id$
17 // $MpId: AliMpSectorSegmentation.cxx,v 1.15 2006/05/24 13:58:46 ivana Exp $
18 // Category: sector
19
20 //-----------------------------------------------------------------------------
21 // Class AliMpSectorSegmentation
22 // -----------------------------
23 // Class describing the segmentation of the sector.        
24 // Provides methods related to pads:
25 // conversion between pad indices, pad location, pad position;
26 // finding pad neighbour.
27 //
28 // Authors: David Guez, Ivana Hrivnacova; IPN Orsay
29 //-----------------------------------------------------------------------------
30
31 #include "AliMpSectorSegmentation.h"
32 #include "AliMpSector.h"
33 #include "AliMpZone.h"
34 #include "AliMpSubZone.h"
35 #include "AliMpRow.h"
36 #include "AliMpVRowSegment.h"
37 #include "AliMpMotifMap.h"
38 #include "AliMpVMotif.h"
39 #include "AliMpMotifPosition.h"
40 #include "AliMpConnection.h"
41 #include "AliMpNeighboursPadIterator.h"
42 #include "AliMpSectorAreaHPadIterator.h"
43 #include "AliMpSectorAreaVPadIterator.h"
44 #include "AliMpSectorPadIterator.h"
45 #include "AliMpIntPair.h"
46 #include "AliMpArea.h"
47 #include "AliMpConstants.h"
48
49 #include "AliLog.h"
50
51 #include <Riostream.h>
52 #include <TMath.h>
53
54 /// \cond CLASSIMP
55 ClassImp(AliMpSectorSegmentation)
56 /// \endcond
57
58 #ifdef WITH_ROOT
59 const Double_t AliMpSectorSegmentation::fgkS1 = 100000.;
60 const Double_t AliMpSectorSegmentation::fgkS2 = 1000.;
61 #endif
62
63 //______________________________________________________________________________
64 AliMpSectorSegmentation::AliMpSectorSegmentation(
65                             const AliMpSector* sector, Bool_t own) 
66   : AliMpVSegmentation(),
67     fkSector(sector),
68     fIsOwner(own),
69     fPadBuffer(0),
70     fPadDimensionsMap(),
71     fMaxIndexInX(0),
72     fMaxIndexInY(0)
73 {
74 /// Standard constructor
75
76   AliDebugStream(1) << "this = " << this << endl;
77
78   fPadBuffer = new AliMpPad(AliMpPad::Invalid());
79   
80   FillPadDimensionsMap();
81 }
82
83 //______________________________________________________________________________
84 AliMpSectorSegmentation::AliMpSectorSegmentation() 
85   : AliMpVSegmentation(),
86     fkSector(0),
87     fIsOwner(false),
88     fPadBuffer(0),
89     fPadDimensionsMap(),      
90     fMaxIndexInX(0),
91     fMaxIndexInY(0)
92 {
93 /// Default constructor
94
95   AliDebugStream(1) << "this = " << this << endl;
96 }
97
98 //______________________________________________________________________________
99 AliMpSectorSegmentation::~AliMpSectorSegmentation() 
100 {
101 /// Destructor 
102
103   AliDebugStream(1) << "this = " << this << endl;
104
105   if ( fIsOwner ) delete fkSector;
106
107   delete fPadBuffer;
108   
109 }
110
111 //
112 // private methods
113 //
114
115 //_____________________________________________________________________________
116 void 
117 AliMpSectorSegmentation::GetAllElectronicCardIDs(TArrayI& ecn) const
118 {
119   /// Fill the array ecn with all manuIds
120
121   GetSector()->GetAllMotifPositionsIDs(ecn);
122 }
123
124 #ifdef WITH_ROOT
125 //______________________________________________________________________________
126 Long_t AliMpSectorSegmentation::GetIndex(const TVector2& vector2) const
127 {
128 /// Convert the two vector to long.
129
130   return Long_t(TMath::Floor((vector2.X()*fgkS1 + vector2.Y())*fgkS2));
131 }  
132
133 //______________________________________________________________________________
134 TVector2  AliMpSectorSegmentation::GetVector(Long_t index) const
135 {
136 /// Convert the long index to twovector.
137
138   return TVector2( TMath::Floor(index/fgkS1)/fgkS2,
139                    (index - TMath::Floor(index/fgkS1)*fgkS1)/fgkS2 );
140 }  
141 #endif
142
143 //______________________________________________________________________________
144 void AliMpSectorSegmentation::FillPadDimensionsMap()
145 {
146 /// Fill the maps between zone ids and pad dimensions.
147
148   for (Int_t i=0; i<fkSector->GetNofZones(); i++) {
149     AliMpZone* zone   = fkSector->GetZone(i+1);
150     Int_t  zoneID = zone->GetID();
151     
152     if (!AliMpConstants::IsEqual(zone->GetPadDimensions(), TVector2())) {
153
154       // regular zone
155 #ifdef WITH_STL
156       fPadDimensionsMap[zoneID*10] = zone->GetPadDimensions();
157 #endif
158 #ifdef WITH_ROOT
159      AliDebugStream(3)
160        << "Filling fPadDimensions[" << zoneID*10 << "] = ("
161        << zone->GetPadDimensions().X() << ", "
162        << zone->GetPadDimensions().Y() << ")" << endl;
163
164      fPadDimensionsMap.Add((Long_t)(zoneID*10), 
165                             GetIndex(zone->GetPadDimensions()));
166 #endif
167     }
168     else {
169       // special zone
170       Int_t subIndex = 0;
171       for (Int_t j=0; j<zone->GetNofSubZones(); j++) {
172         AliMpSubZone* subZone = zone->GetSubZone(j);
173         AliMpVMotif*  motif = subZone->GetMotif();
174         
175         for (Int_t k=0; k<motif->GetNofPadDimensions(); k++) {
176           Int_t index = zoneID*10 +  subIndex++;
177 #ifdef WITH_STL
178           fPadDimensionsMap[index] = motif->GetPadDimensions(k);
179 #endif
180 #ifdef WITH_ROOT
181           AliDebugStream(3)
182             << "Filling fPadDimensions[" << index << "] = ("
183             << motif->GetPadDimensions(k).X() << ", "
184             << motif->GetPadDimensions(k).Y() << ") motif "
185             << motif->GetID().Data() << "-" << k << endl;
186
187           fPadDimensionsMap.Add((Long_t)(index), 
188                             GetIndex(motif->GetPadDimensions(k)));
189 #endif
190         }
191       }   
192     }     
193   }      
194 }
195
196 //______________________________________________________________________________
197 AliMpMotifPosition* 
198 AliMpSectorSegmentation::FindMotifPosition(const AliMpIntPair& indices) const
199 {
200 /// Find the motif position which contains the given pad indices
201 /// return 0 if not found
202
203   switch (fkSector->GetDirection()) {
204     case AliMp::kX : {
205     // Case where all the pads have the same size along X direction
206
207       for (Int_t irow=0; irow<fkSector->GetNofRows(); ++irow) {
208         AliMpRow* row = fkSector->GetRow(irow);
209         if (row->GetLowIndicesLimit().GetFirst()<=indices.GetFirst() &&
210             row->GetHighIndicesLimit().GetFirst()>=indices.GetFirst()) {
211             
212            for (Int_t iseg=0;iseg<row->GetNofRowSegments();++iseg){
213             AliMpVRowSegment* seg = row->GetRowSegment(iseg);
214             if (seg->GetLowIndicesLimit().GetFirst()<=indices.GetFirst() &&
215                 seg->GetHighIndicesLimit().GetFirst()>=indices.GetFirst()) {
216
217               AliMpMotifPosition* motifPos;
218               for (Int_t imot=0;imot<seg->GetNofMotifs();++imot) {
219                 motifPos 
220                   = fkSector->GetMotifMap()
221                     ->FindMotifPosition(seg->GetMotifPositionId(imot));
222                 if (motifPos && motifPos->HasPad(indices)) return motifPos;
223               }
224             }
225           }
226         }
227       }
228       return 0;
229     }
230     break;
231     ////////////////////////////////////////////////////////////////////////////////
232     case AliMp::kY : {
233       // Case where all the pads have the same size along Y direction   
234       // look for the row which contains the indices
235       AliMpRow* row=0;
236       Int_t irow;
237       for (irow=0; irow<fkSector->GetNofRows(); ++irow) {
238         row = fkSector->GetRow(irow);
239         AliMpVRowSegment* lastSeg = row->GetRowSegment(row->GetNofRowSegments()-1);
240         if (lastSeg->GetLowIndicesLimit().GetSecond()<=indices.GetSecond() &&
241             lastSeg->GetHighIndicesLimit().GetSecond()>=indices.GetSecond()) break;
242         // NOTE : We use the last row segment in order to ensure that
243         // we are not on a special motif
244       }
245       if (irow==fkSector->GetNofRows()) return 0;
246       // look for the row segment, in the found row, which contains the indices
247       AliMpVRowSegment* seg=0;
248       Int_t iseg;
249       for (iseg=0;iseg<row->GetNofRowSegments();++iseg){
250         seg = row->GetRowSegment(iseg);
251         if (seg->HasIndices(indices)) break;
252       }
253       if (iseg==row->GetNofRowSegments()) return 0;
254   
255       // look for the motif position which contains the indices
256       AliMpMotifPosition* motifPos=0;
257       Int_t imot=0;
258       for (imot=0;imot<seg->GetNofMotifs();++imot) {
259         motifPos 
260           = fkSector->GetMotifMap()
261             ->FindMotifPosition(seg->GetMotifPositionId(imot));
262         if (motifPos && motifPos->HasPad(indices)) break;
263       }      
264       if (imot==seg->GetNofMotifs()) return 0;
265    
266       return motifPos;      
267     }
268     default: return 0;
269   }
270 }
271
272 //______________________________________________________________________________
273 AliMpPad 
274 AliMpSectorSegmentation::PadByXDirection(const TVector2& startPosition, 
275                                          Double_t maxX) const
276 {
277 /// Find the first valid pad from starting position in the
278 /// direction of pad lines up to distance dx.
279
280   // Define step limits
281   Double_t  stepX = fkSector->GetMinPadDimensions().X();
282  
283   // Search in X direction
284   AliMpPad pad;
285   TVector2 position(startPosition);    
286   do {
287     pad = PadByPosition(position, false);
288     position += TVector2(stepX, 0.);
289   }   
290   while ( !pad.IsValid() && 
291           position.X() - fkSector->GetMaxPadDimensions().X() < maxX ); 
292   
293   // Invalidate pad if it is outside limits
294   if ( (pad.Position().X() - pad.Dimensions().X()) > maxX ) 
295     pad = AliMpPad::Invalid();
296
297   return pad;
298 }
299
300 //______________________________________________________________________________
301 AliMpPad 
302 AliMpSectorSegmentation::PadByYDirection(const TVector2& startPosition, 
303                                          Double_t maxY) const
304 {
305 /// Find the first valid pad from starting position in the
306 /// direction of pad columns up to distance dx.
307   
308   // Define step limits
309   Double_t stepY = fkSector->GetMinPadDimensions().Y();
310  
311   // Search in Y direction
312   AliMpPad pad;
313   TVector2 position(startPosition);    
314   do {
315     pad = PadByPosition(position, false);
316     position += TVector2(0., stepY);
317   }   
318   while ( !pad.IsValid() && 
319           position.Y() - fkSector->GetMaxPadDimensions().Y()< maxY ); 
320   
321   // Invalidate pad if it is outside limits
322   if ((pad.Position().Y() - pad.Dimensions().Y()) > maxY) 
323     pad = AliMpPad::Invalid();
324
325   return pad;
326 }
327
328 //
329 // public methods
330 //
331
332 //______________________________________________________________________________
333 AliMpVPadIterator* 
334 AliMpSectorSegmentation::CreateIterator() const
335 {
336 /// Create the sector iterator
337
338   return new AliMpSectorPadIterator(fkSector);
339 }
340
341 //______________________________________________________________________________
342 AliMpVPadIterator* 
343 AliMpSectorSegmentation::CreateIterator(const AliMpArea& area) const
344 {
345 /// Create the area iterator. 
346
347   switch (fkSector->GetDirection()) {
348   
349     case AliMp::kX: return new AliMpSectorAreaVPadIterator(this, area);
350              ;;
351     case AliMp::kY: return new AliMpSectorAreaHPadIterator(this, area);
352              ;;
353   }
354   
355   Fatal("CreateIterator", "Incomplete switch on Sector direction");
356   return 0;  
357 }   
358   
359 //______________________________________________________________________________
360 Int_t 
361 AliMpSectorSegmentation::GetNeighbours(const AliMpPad& pad, TObjArray& neighbours,
362                                        Bool_t includeSelf,
363                                        Bool_t includeVoid) const
364 {
365   /// Uses default implementation
366   return AliMpVSegmentation::GetNeighbours(pad,neighbours,includeSelf,includeVoid);
367 }
368
369 //______________________________________________________________________________
370 AliMpVPadIterator* 
371 AliMpSectorSegmentation::CreateIterator(const AliMpPad& centerPad,
372                                         Bool_t includeCenter) const
373 {
374   /// Create the neighbours pad iterator.
375
376   return new AliMpNeighboursPadIterator(this, centerPad, includeCenter);
377 }   
378   
379 //______________________________________________________________________________
380 TVector2
381 AliMpSectorSegmentation::Dimensions() const
382 {
383   return GetSector()->Dimensions();
384 }
385
386 //______________________________________________________________________________
387 AliMp::PlaneType
388 AliMpSectorSegmentation::PlaneType() const
389 {
390   return GetSector()->GetPlaneType();
391 }
392
393 //______________________________________________________________________________
394 AliMpPad 
395 AliMpSectorSegmentation::PadByLocation(const AliMpIntPair& location, 
396                                        Bool_t warning) const
397 {
398 /// Find the pad which corresponds to the given location
399   
400   if ((*fPadBuffer).GetLocation()==location) return (*fPadBuffer);
401   
402   AliMpMotifPosition* motifPos = 
403     fkSector->GetMotifMap()->FindMotifPosition(location.GetFirst());
404   if (!motifPos){
405     if (warning) Warning("PadByLocation","The pad motif position ID doesn't exists");
406     return AliMpPad::Invalid();
407   }
408   
409   AliMpVMotif* motif = motifPos->GetMotif();
410   AliMpIntPair localIndices = 
411     motif->GetMotifType()->FindLocalIndicesByGassiNum(location.GetSecond());
412   if (! localIndices.IsValid()) {
413     if (warning) Warning("PadByLocation","The pad number doesn't exists");
414     return AliMpPad::Invalid();
415   }
416   TVector2 delta = motif->PadPositionLocal(localIndices);
417   return (*fPadBuffer) = AliMpPad(location,
418               motifPos->GlobalIndices(localIndices),
419               motifPos->Position()+delta,
420               motif->GetPadDimensions(localIndices));
421
422 }
423 //______________________________________________________________________________
424 AliMpPad 
425 AliMpSectorSegmentation::PadByIndices(const AliMpIntPair& indices,
426                                       Bool_t warning ) const
427 {
428 /// Find the pad which corresponds to the given indices  
429
430   if ((*fPadBuffer).GetIndices()==indices) return (*fPadBuffer);    
431    
432   AliMpMotifPosition* motifPos = FindMotifPosition(indices);
433   if (!motifPos) {    
434     if (warning) 
435       Warning("PadByIndices","Pad indices not contained in any motif!");
436     return AliMpPad::Invalid();
437   }
438   
439   // retrieve the local indices in the found motif
440   AliMpVMotif* motif = motifPos->GetMotif();
441   AliMpIntPair localIndices = indices - motifPos->GetLowIndicesLimit();
442   
443   AliMpConnection* connection=
444     motif->GetMotifType()->FindConnectionByLocalIndices(localIndices);
445     
446   if (!connection){
447     if (warning) Warning("PadByIndices","No connection with the given indices!");
448     return AliMpPad::Invalid();
449   }
450
451   TVector2 localPos = motif->PadPositionLocal(localIndices);
452
453   return (*fPadBuffer) 
454     = AliMpPad(AliMpIntPair(motifPos->GetID(),connection->GetGassiNum()),
455                indices,
456                motifPos->Position()+localPos,
457                motif->GetPadDimensions(localIndices)); 
458
459 }
460 //______________________________________________________________________________
461 AliMpPad 
462 AliMpSectorSegmentation::PadByPosition(const TVector2& position,
463                                        Bool_t warning) const
464 {
465 /// Find the pad which corresponds to the given position
466
467   if ((*fPadBuffer).Position().X()==position.X() && 
468       (*fPadBuffer).Position().Y()==position.Y()) return (*fPadBuffer);  
469
470   Int_t motifPosID = fkSector->FindMotifPositionId(position);
471   AliMpMotifPosition* motifPos 
472     = fkSector->GetMotifMap()
473         ->FindMotifPosition(motifPosID);
474     
475   if (!motifPos){
476     if (warning) Warning("PadByPosition","Position outside limits");
477     return AliMpPad::Invalid();
478   }
479
480   AliMpVMotif* motif =  motifPos->GetMotif();  
481   AliMpIntPair localIndices 
482     = motif->PadIndicesLocal(position-motifPos->Position());
483     
484   AliMpConnection* connect = 
485     motif->GetMotifType()->FindConnectionByLocalIndices(localIndices);
486
487    if (!connect){
488     if (warning) Warning("PadByPosition","Position outside motif limits");
489     return AliMpPad::Invalid();
490   }
491   
492   return (*fPadBuffer)
493     = AliMpPad(AliMpIntPair(motifPosID,connect->GetGassiNum()),
494                motifPos->GlobalIndices(localIndices),
495                motifPos->Position()+motif->PadPositionLocal(localIndices),
496                motif->GetPadDimensions(localIndices));
497
498 }
499
500 //______________________________________________________________________________
501 AliMpPad 
502 AliMpSectorSegmentation::PadByDirection(const TVector2& startPosition, 
503                                         Double_t distance) const
504 {
505 /// Find the first valid pad from starting position in the
506 /// direction of pad lines/columns up to the specified distance.
507 /// Pad lines are the lines of pads in the sector with constant pad y size,
508 /// pad columns are the columns of pads in the sector with constant pad x size. 
509
510   switch (fkSector->GetDirection()) {
511   
512     case AliMp::kX: return PadByYDirection(startPosition, distance);
513              ;;
514     case AliMp::kY: return PadByXDirection(startPosition, distance);
515              ;;
516   }
517   
518   Fatal("PadByDirection", "Incomplete switch on Sector direction");
519   return AliMpPad::Invalid();  
520 }
521
522 //______________________________________________________________________________
523 Int_t  AliMpSectorSegmentation::MaxPadIndexX() const
524 {
525 /// Return maximum pad index in x
526
527   return fkSector->GetMaxPadIndices().GetFirst();
528 }
529
530 //______________________________________________________________________________
531 Int_t  AliMpSectorSegmentation::MaxPadIndexY() const
532 {
533 /// Return maximum pad index in y
534
535   return fkSector->GetMaxPadIndices().GetSecond();
536 }
537
538 //______________________________________________________________________________
539 Int_t  AliMpSectorSegmentation::NofPads() const
540 {
541 /// Return number of pads defined in the sector
542
543   return fkSector->GetNofPads();
544 }
545
546 //______________________________________________________________________________
547 Bool_t AliMpSectorSegmentation::HasPad(const AliMpIntPair& indices) const
548 {
549 /// Does the pad specified by \a indices exist ?
550
551   return PadByIndices(indices,kFALSE) != AliMpPad::Invalid();
552 }
553
554 //______________________________________________________________________________
555 Bool_t AliMpSectorSegmentation::HasMotifPosition(Int_t motifPositionID) const
556 {
557 /// Does the motif position specified by motifPositionID exist ?
558
559   return (fkSector->GetMotifMap()->FindMotifPosition(motifPositionID) != 0);
560 }
561
562 //______________________________________________________________________________
563 TVector2  AliMpSectorSegmentation::GetMinPadDimensions() const
564 {
565 /// Returne the dimensions of the smallest pad.
566
567   return fkSector->GetMinPadDimensions();
568 }  
569
570 //______________________________________________________________________________
571 Int_t AliMpSectorSegmentation::Zone(const AliMpPad& pad, Bool_t warning) const
572 {
573 /// Return the zone index of the zone containing the specified pad.
574 /// This zone index is different from the zone ID,
575 /// as it is unique for each pad dimensions.
576 /// It is composed in this way:
577 ///   zoneID*10 + specific index 
578 /// Specific index is present only for zones containing special motifs.
579
580   if (!pad.IsValid()) {
581     if (warning) Warning("Zone(AliMpPad)", "Invalid pad");
582     return 0;
583   }  
584
585 #ifdef WITH_STL
586   PadDimensionsMapCIterator it;
587   for (it = fPadDimensionsMap.begin(); it != fPadDimensionsMap.end(); ++it) {
588     if (AliMpConstants::IsEqual(it->second, pad.Dimensions()))
589       return it->first;
590   }
591 #endif
592
593 #ifdef WITH_ROOT
594   PadDimensionsMapCIterator it(&fPadDimensionsMap);
595   Long_t key, value;
596   while ( it.Next(key, value) ) {
597     TVector2 dimensions =  GetVector(value);
598     if (AliMpConstants::IsEqual(dimensions, pad.Dimensions()))
599       return (Int_t)key;
600   } 
601   
602   AliError(Form("fPadDimensionsMap size is %d",fPadDimensionsMap.GetSize()));
603   
604 #endif
605
606   // Should never happen
607   AliErrorStream() 
608     << "Zone(AliMpPad pad) not found, where pad is: " << pad << endl;
609   return 0;
610 }  
611
612 //______________________________________________________________________________
613 TVector2 
614 AliMpSectorSegmentation::PadDimensions(Int_t zone, Bool_t warning) const
615 {
616 /// Return the pad dimensions for the zone with the specified zone index.
617
618 #ifdef WITH_STL
619   PadDimensionsMapCIterator it = fPadDimensionsMap.find(zone);
620   if (it != fPadDimensionsMap.end()) return it->second;
621 #endif
622
623 #ifdef WITH_ROOT
624   Long_t value = fPadDimensionsMap.GetValue(zone);
625   if (value) return GetVector(value);
626 #endif
627
628   if (warning) Warning("PadDimensions(zone)", "not found");
629   return TVector2();
630 }  
631
632 //______________________________________________________________________________
633 Bool_t AliMpSectorSegmentation::CircleTest(const AliMpIntPair& indices) const
634 {
635 /// Verify that all methods for retrieving pads are consistents between them.
636 /// Return true if the pad with specified indices was found and verified,
637 /// false otherwise.
638
639   if (!HasPad(indices)) return false;
640
641   // Verify the indice->location->position->indice way
642   AliMpIntPair location = PadByIndices(indices).GetLocation();
643   TVector2 position = PadByLocation(location).Position();
644   AliMpIntPair retIndices = PadByPosition(position).GetIndices();
645     
646   if (retIndices != indices) {
647     cout << "Pad " << indices << " lead to inconsistency" << endl;
648     cout << "in indice->location->position->indice way..." << endl;
649     cout << "starting from " << indices << "-->" << location << "-->" 
650          << '(' << position.X() << ',' << position.Y() << ')'
651          << " and retIndices: " << retIndices << endl;
652   }
653     
654     
655   // Verify the indice->position->location->indice way    
656   position = PadByIndices(indices).Position();
657   location = PadByPosition(position).GetLocation();
658   retIndices = PadByLocation(location).GetIndices();
659
660   if (retIndices != indices) {
661     cout << "Pad " << indices << " lead to inconsistency" << endl;
662     cout << "in indice->position->location->indice way..." <<endl;
663     cout << "starting from " << indices 
664          << " and retIndices: " << retIndices << endl;
665   }
666   
667   return true;
668 }
669
670 //______________________________________________________________________________
671 void
672 AliMpSectorSegmentation::Print(Option_t* opt) const
673 {
674 /// Print the sector
675
676   fkSector->Print(opt);
677 }
678
679 //______________________________________________________________________________
680 void AliMpSectorSegmentation::PrintZones() const
681 {
682 /// Print all zones and pads dimensions from the map.
683
684   cout << "Zones: " << endl;
685
686 #ifdef WITH_STL
687   PadDimensionsMapCIterator it;
688   for (it = fPadDimensionsMap.begin(); it != fPadDimensionsMap.end(); ++it) {
689     cout << "    zone: " <<   setw(4) << it->first;
690     cout << "    pad dimensions: ( " 
691          << it->second.X() << ", " << it->second.Y() << ")" << endl; 
692   }
693 #endif
694
695 #ifdef WITH_ROOT
696   PadDimensionsMapCIterator it(&fPadDimensionsMap);
697   Long_t key, value;
698   while ( it.Next(key, value) ) {
699     //cout << "Iterating over: " << key << ", " << value << endl;
700     TVector2 dimensions =  GetVector(value);
701
702     cout << "    zone: " <<   setw(4) << key;
703     cout << "    pad dimensions: ( " 
704          << dimensions.X() << ", " << dimensions.Y() << ")" << endl; 
705   }
706 #endif
707 }
708