]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/mapping/macros/testAllIndices.C
In mapping/macros:
[u/mrichter/AliRoot.git] / MUON / mapping / macros / testAllIndices.C
CommitLineData
66e0997c 1// $Id$
a387ee7c 2// $MpId: testAllIndices.C,v 1.7 2005/08/24 08:53:27 ivana Exp $
66e0997c 3//
4// Test macro for testing which pad is seen as "existing" by AliMpSector.
5
f05d3eb1 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 "AliMpVPainter.h"
17#include "AliMpRow.h"
18#include "AliMpVRowSegment.h"
19#include "AliMpMotifMap.h"
20#include "AliMpMotifPosition.h"
21#include "AliMpMotifType.h"
22
23#include <Riostream.h>
24#include <TCanvas.h>
25#include <TH2.h>
26
27#endif
28
29void testAllIndices(AliMq::Station12Type station, AliMp::PlaneType plane)
66e0997c 30{
f05d3eb1 31 AliMpDataProcessor mp;
32 AliMpDataMap* dataMap = mp.CreateDataMap("data");
33 AliMpDataStreams dataStreams(dataMap);
34
35 AliMpSectorReader r(dataStreams, station, plane);
66e0997c 36
37 AliMpSector *sector=r.BuildSector();
72d8f376 38 AliMpSectorSegmentation segmentation(sector);
66e0997c 39 AliMpVPainter* painter = AliMpVPainter::CreatePainter(sector);
40
41 TCanvas* c1 = new TCanvas("view",
42 "MSectorPainter::Draw() output (view per pad)");
43 painter->Draw("ZSSMP");
44 c1->Update();
45
72d8f376 46 Int_t maxPadIndexX = segmentation.MaxPadIndexX();
47 Int_t maxPadIndexY = segmentation.MaxPadIndexY();
48
49 // Define histogram limits
50 Int_t nx = (maxPadIndexX/10 + 1)*10;
51 Int_t ny = (maxPadIndexY/10 + 1)*10;
52 TH2C* histo = new TH2C("histo","Existing pads",
53 nx, -0.5, nx-0.5, ny, -0.5, ny-0.5);
54
f05d3eb1 55 Int_t nx2 = 95/2;
56 Int_t ny2 = 95/2;
57 if (station == AliMq::kStation2) {
58 nx2 = 120/2;
59 ny2 = 120/2;
72d8f376 60 }
61 TH2F* histo2 = new TH2F("histo2","Existing positions",
62 nx2, 0, nx2*2, ny2, 0, ny2*2);
63
64 // Define canvas
66e0997c 65 TCanvas* c2 = new TCanvas("c2","Only existing pads are filled");
66 TCanvas* c3 = new TCanvas("c3","Positions");
66e0997c 67
f05d3eb1 68 for ( Int_t irow=0; irow<sector->GetNofRows(); irow++ ) {
66e0997c 69 AliMpRow* row = sector->GetRow(irow);
70
f05d3eb1 71 for ( Int_t iseg=0; iseg<row->GetNofRowSegments(); iseg++ ) {
66e0997c 72 AliMpVRowSegment* seg = row->GetRowSegment(iseg);
73
f05d3eb1 74 for ( Int_t imot=0; imot<seg->GetNofMotifs(); imot++) {
66e0997c 75 AliMpMotifPosition* motifPos
76 = sector->GetMotifMap()->FindMotifPosition(seg->GetMotifPositionId(imot));
77
f05d3eb1 78 for ( Int_t gassNum=0; gassNum<64; gassNum++ ) {
66e0997c 79 if (motifPos->GetMotif()->GetMotifType()->FindConnectionByGassiNum(gassNum)){
80
81 AliMpPad pad = segmentation.PadByLocation(AliMpIntPair(motifPos->GetID(),gassNum));
82 if (pad != AliMpPad::Invalid()) {
f05d3eb1 83 histo->Fill(pad.GetIndices().GetFirst(),
66e0997c 84 pad.GetIndices().GetSecond());
85 histo2->Fill(pad.Position().X(),
86 pad.Position().Y());
87 }
88 }
89 }
90 }
91 }
92 }
93 c2->cd();
94 histo->Draw("col");
95 c3->cd();
96 histo2->Draw("box");
97}
f05d3eb1 98
99void testAllIndices()
100{
101 AliMq::Station12Type station[2] = { AliMq::kStation1, AliMq::kStation2 };
102 AliMp::PlaneType plane[2] = { AliMp::kBendingPlane, AliMp::kNonBendingPlane };
103
104 for ( Int_t is = 0; is < 2; is++ ) {
105 for ( Int_t ip = 0; ip < 2; ip++ ) {
106
107 cout << "Running testAllIndices for "
108 << AliMq::Station12TypeName(station[is]) << " "
109 << AliMp::PlaneTypeName(plane[ip]) << " ... " << endl;
110
111 testAllIndices(station[is], plane[ip]);
112
113 cout << "... end running " << endl << endl;
114 }
115 }
116}
117