]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/macros/testAnyPadIterators.C
348284a112d177161730cf1e9d407e5a2c9a6fd0
[u/mrichter/AliRoot.git] / MUON / mapping / macros / testAnyPadIterators.C
1 // $Id$
2 // $MpId: testAnyPadIterators.C,v 1.17 2006/03/15 13:07:07 ivana Exp $
3 //
4 // Test macro for reading  sector, and iterate over it
5
6 #if !defined(__CINT__) || defined(__MAKECINT__)
7
8 #include "AliMpStation12Type.h"
9 #include "AliMpPlaneType.h"
10 #include "AliMpDataProcessor.h"
11 #include "AliMpDataMap.h"
12 #include "AliMpDataStreams.h"
13 #include "AliMpSector.h"
14 #include "AliMpSectorReader.h"
15 #include "AliMpSectorSegmentation.h" 
16 #include "AliMpMotifType.h"
17 #include "AliMpMotifMap.h"
18 #include "AliMpVMotif.h"
19 #include "AliMpVPadIterator.h"
20 #include "AliMpSectorPadIterator.h"
21 #include "AliMpNeighboursPadIterator.h"
22 #include "AliMpMotifPositionPadIterator.h"
23 #include "AliMpConstants.h"
24
25 #include <Riostream.h>
26 #include <TCanvas.h>
27 #include <TMarker.h>
28 #include <TH2.h>
29
30 #endif
31
32 class AliMpVPadIterator;
33 void MarkPads(AliMpVPadIterator& it,Int_t xmax,Int_t ymax,Bool_t print=kTRUE)
34 {
35 // This function works with pad iterator, no matter which kind of iterator
36 // it is. So it can be used for drawing all pad of the sector (AliMpSectorPadIterator)
37 // or all pad around a given pad (AliMpNeighboursPadIterator), as with pads
38 // of a given motif type (AliMpMotifTypePadIterator)
39
40   Int_t num=0;
41
42   for (it.First(); ! it.IsDone(); it.Next()){
43     AliMpIntPair indices = it.CurrentItem().GetIndices();
44     if (print) cout<<"Iterator number "<< num << " at "<< indices <<endl;
45     num++;
46     TMarker* marker = new TMarker( (Double_t)indices.GetFirst() /xmax,
47                                    (Double_t)indices.GetSecond()/ymax,
48                                    2);
49     marker->Draw();
50   }
51 }
52
53 void testAnyPadIterators(AliMq::Station12Type station, AliMp::PlaneType plane,
54                          Int_t i=50, Int_t j=50)
55 {
56   AliMpDataProcessor mp;
57   AliMpDataMap* dataMap = mp.CreateDataMap("data");
58   AliMpDataStreams dataStreams(dataMap);
59
60   AliMpSectorReader r(dataStreams, station, plane);
61   AliMpSector* sector = r.BuildSector();
62   AliMpSectorSegmentation segmentation(sector);
63     
64   TCanvas *canv = new TCanvas("canv");
65   canv->Divide(2,2);
66   
67   canv->cd(1);
68   AliMpSectorPadIterator its(sector);
69   MarkPads(its, 150, 250, kFALSE);
70
71   canv->cd(2);
72   AliMpVMotif* motif = sector->FindMotif(TVector2(30,3));
73   if ( motif ) {
74     AliMpMotifType* motifType = motif->GetMotifType();
75     AliMpMotifTypePadIterator itm(motifType);
76     MarkPads(itm,15,15);
77     cout << "______________ MotifType "  << motifType->GetID() 
78          << "__________________________" << endl;
79   } 
80   else 
81     cout << "No motif found at given position..." << endl;
82   
83   canv->cd(3);
84   AliMpPad pad = segmentation.PadByIndices(AliMpIntPair(i,j));
85   AliMpNeighboursPadIterator itn(&segmentation,pad);
86   MarkPads(itn,i+8,j+8);
87   cout<<"________________ Neighbours __________________________"<<endl;
88   
89   canv->cd(4);
90   Int_t motifPosId = 20 | AliMpConstants::ManuMask(plane); 
91   if ( plane == AliMp::kNonBendingPlane ) motifPosId = 19;
92   AliMpMotifPosition* motifPos = sector->GetMotifMap()->FindMotifPosition(motifPosId);
93   if ( motifPos ){
94     AliMpMotifPositionPadIterator itmp(motifPos);
95     MarkPads(itmp,15,15);
96     cout<<"_________________ MotifPosition _________________________"<<endl;
97   }
98 }
99
100 void testAnyPadIterators()
101 {
102   AliMq::Station12Type  station[2] = { AliMq::kStation1, AliMq::kStation2 }; 
103   AliMp::PlaneType      plane[2]   = { AliMp::kBendingPlane, AliMp::kNonBendingPlane };
104   
105   for ( Int_t is = 0; is < 2; is++ ) {
106     for ( Int_t ip = 0; ip < 2; ip++ ) {
107     
108       cout << "Running testAnyPadIterators for " 
109            << AliMq::Station12TypeName(station[is]) << "  "
110            << AliMp::PlaneTypeName(plane[ip])  << " ... " << endl;
111        
112       testAnyPadIterators(station[is], plane[ip]);
113     
114       cout << "... end running " << endl << endl;
115     }  
116   }   
117 }  
118