]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/AliMpPlaneSegmentation.cxx
Added copy constructor and assignement operator (I. Hrivnacova)
[u/mrichter/AliRoot.git] / MUON / mapping / AliMpPlaneSegmentation.cxx
1 // $Id$
2 // Category: plane
3 //
4 // Class AliMpPlaneSegmentation
5 // ----------------------------
6 // Class describing the segmentation of the plane.
7 //
8 // Transformation of pad characteristics according to sectors:
9 //
10 //   I.  ( posId,  Guassi ), ( i, j), ( x, y)         II. |  I.
11 //  II.  ( posId', Guassi'), (-i, j), (-x, y)       _____ | ____
12 // III.  (-posId,  Guassi),  (-i,-j), (-x,-y)             |
13 //  IV.  (-posId', Guassi'), ( i,-j), ( x,-y)        III. |  IV.
14 //   
15 // Where (posId', Guassi') is the location of the pad
16 // in the clipped sector.
17 //
18 // Authors: David Guez, Ivana Hrivnacova; IPN Orsay
19
20 #include <Riostream.h>
21 #include <TMath.h>
22 #include <TError.h>
23
24 #include "AliMpPlaneSegmentation.h"
25 #include "AliMpPlaneAreaPadIterator.h"
26 #include "AliMpPlane.h"
27 #include "AliMpSectorPosition.h"
28 #include "AliMpSectorSegmentation.h"
29
30 ClassImp(AliMpPlaneSegmentation)
31
32 //_____________________________________________________________________________
33 AliMpPlaneSegmentation::AliMpPlaneSegmentation(const AliMpPlane* plane) 
34   : AliMpVSegmentation(),
35     fkPlane(plane),
36     fFrontSectorSegmentation(0),
37     fBackSectorSegmentation(0)
38 {
39 //
40   fFrontSectorSegmentation = new AliMpSectorSegmentation(plane->GetFrontSector());
41   fBackSectorSegmentation = new AliMpSectorSegmentation(plane->GetBackSector());
42   
43   for (Int_t i=0; i<fkPlane->GetNofSectorPositions(); i++) {
44
45 #ifdef WITH_STL
46     fTransformers.push_back(
47       new AliMpTransformer(fkPlane->GetSectorPosition(i)->GetOffset(),
48                            fkPlane->GetSectorPosition(i)->GetScale()));
49 #endif
50
51 #ifdef WITH_ROOT
52     fTransformers.Add(
53       new AliMpTransformer(fkPlane->GetSectorPosition(i)->GetOffset(),
54                            fkPlane->GetSectorPosition(i)->GetScale()));
55 #endif
56  }                     
57 }
58
59 ///_____________________________________________________________________________
60 AliMpPlaneSegmentation::AliMpPlaneSegmentation() 
61   : AliMpVSegmentation(),
62     fkPlane(0),
63     fFrontSectorSegmentation(0),
64     fBackSectorSegmentation(0)
65 {
66 //
67 }
68
69 //_____________________________________________________________________________
70 AliMpPlaneSegmentation::AliMpPlaneSegmentation(
71                                   const AliMpPlaneSegmentation& right) 
72   : AliMpVSegmentation(right) {
73 // 
74   Fatal("AliMpPlaneSegmentation", "Copy constructor not provided.");
75 }
76
77 //_____________________________________________________________________________
78 AliMpPlaneSegmentation::~AliMpPlaneSegmentation() {
79 // 
80   delete fFrontSectorSegmentation;
81   delete fBackSectorSegmentation;
82
83   for (Int_t i=0; i<GetNofTransformers(); i++) 
84     delete GetTransformer(i);
85 }
86
87 //
88 // operators
89 //
90
91 //_____________________________________________________________________________
92 AliMpPlaneSegmentation& 
93 AliMpPlaneSegmentation::operator=(const AliMpPlaneSegmentation& right)
94 {
95   // check assignement to self
96   if (this == &right) return *this;
97
98   Fatal("operator =", "Assignement operator not provided.");
99     
100   return *this;  
101 }    
102
103 //
104 // private methods
105 //
106
107 //_____________________________________________________________________________
108 const AliMpTransformer* 
109 AliMpPlaneSegmentation::GetTransformer(const AliMpIntPair& scale) const
110 {
111 // Returns the sector position specified by scale.
112 // ---
113
114   for (Int_t i=0; i<GetNofTransformers(); i++) 
115     if (GetTransformer(i)->GetScale() == scale) return GetTransformer(i);
116
117   Fatal("GetTransformer", "Wrong scale");
118   return 0; 
119 }
120
121 //_____________________________________________________________________________
122 AliMpIntPair AliMpPlaneSegmentation::GetScale(const AliMpIntPair& pair) const
123 {
124 // Returns pair of the signs of the values of the given pair.
125 // ---
126
127   AliMpIntPair scale(1, 1);
128   
129   if (pair.GetFirst() < 0)  scale.SetFirst(-1);
130   if (pair.GetSecond() < 0) scale.SetSecond(-1);
131   
132   return scale;
133 }
134
135 //_____________________________________________________________________________
136 AliMpIntPair AliMpPlaneSegmentation::GetScale(const TVector2& vector) const
137 {
138 // Returns pair of the signs of the values of the given vector.
139 // ---
140
141   AliMpIntPair scale(1, 1);
142   
143   if (vector.X() < 0) scale.SetFirst(-1);
144   if (vector.Y() < 0) scale.SetSecond(-1);
145   
146   return scale;
147 }
148
149 //_____________________________________________________________________________
150 AliMpIntPair 
151 AliMpPlaneSegmentation::GetLocationScale(const AliMpIntPair& location) const
152 {
153 // Returns the scale transformation of the specified location. 
154 // ---
155
156   // Find the sector
157   Bool_t inFront;
158   if (fFrontSectorSegmentation
159         ->HasMotifPosition(TMath::Abs(location.GetFirst())))
160     inFront = true;
161   else if (fBackSectorSegmentation
162              ->HasMotifPosition(TMath::Abs(location.GetFirst())))  
163     inFront = false;
164   else {
165     Fatal("GetLocationScale", "Motif position not found.");
166     return AliMpIntPair();
167   }  
168     
169   if      (inFront  && location.GetFirst() > 0) return  AliMpIntPair(1, 1); 
170   else if (inFront  && location.GetFirst() < 0) return  AliMpIntPair(-1, -1);
171   else if (!inFront && location.GetFirst() > 0) return  AliMpIntPair(-1, 1);
172   else if (!inFront && location.GetFirst() < 0) return  AliMpIntPair( 1,-1);  
173
174   // cannot get there
175   Fatal("GetLocationScale", "Condition failed.");
176   return AliMpIntPair();
177 }
178
179
180 //_____________________________________________________________________________
181 AliMpSectorSegmentation* 
182 AliMpPlaneSegmentation::GetSectorSegmentation(const AliMpIntPair& scale) const
183 {    
184 // Returns front sector or back sector segmentation
185 // according to quadrant specified by scale.
186 // ---
187
188   if (scale.GetFirst()*scale.GetSecond() > 0) {
189     // quadrant I or III
190     return fFrontSectorSegmentation;
191   }  
192   else  {
193     // quadrant II or IV
194     return fBackSectorSegmentation;
195   }  
196 }
197
198 //_____________________________________________________________________________
199 AliMpSectorSegmentation* 
200 AliMpPlaneSegmentation::GetSectorSegmentation(Int_t motifPositionId) const
201 {    
202 // Returns front sector or back sector segmentation
203 // according to specified motifPositionId
204 // ---
205
206   if (fFrontSectorSegmentation->HasMotifPosition(motifPositionId))
207     return fFrontSectorSegmentation;
208   else if (fBackSectorSegmentation->HasMotifPosition(motifPositionId)) 
209     return fBackSectorSegmentation;
210   else {
211     Fatal("GetSectorSegmentation", "Motif position not found.");
212     return 0;
213   }
214 }
215
216 //
217 // public methods
218 //
219
220 //_____________________________________________________________________________
221 AliMpVPadIterator* 
222 AliMpPlaneSegmentation::CreateIterator(const AliMpArea& area) const
223 {
224 // Creates the are iterator. 
225 // (The inherited method cannot be used)
226 // ---
227
228   return new AliMpPlaneAreaPadIterator(this, area);  
229 }   
230   
231 //______________________________________________________________________________
232 AliMpPad AliMpPlaneSegmentation::PadByLocation(const AliMpIntPair& location, 
233                                        Bool_t warning) const
234 {
235 // Find the pad which corresponds to the given location
236 // ---
237
238   // Get segmentation
239   AliMpSectorSegmentation* segmentation 
240     = GetSectorSegmentation(TMath::Abs(location.GetFirst()));
241
242   // Get pad in the segmentation
243   AliMpPad pad
244      = segmentation
245        ->PadByLocation(
246            AliMpIntPair(TMath::Abs(location.GetFirst()),location.GetSecond()), 
247                         warning);
248
249   // Get transformation
250   AliMpIntPair scale  = GetLocationScale(location);
251   const AliMpTransformer* kTransformer = GetTransformer(scale);
252   
253   // Transform pad characteristics
254   return kTransformer->Transform(pad);        
255 }
256
257 //______________________________________________________________________________
258 AliMpPad AliMpPlaneSegmentation::PadByIndices (const AliMpIntPair& indices,
259                                                Bool_t warning ) const
260 {
261 // Find the pad which corresponds to the given indices  
262 //
263
264   AliMpIntPair scale = GetScale(indices);
265   const AliMpTransformer* kTransformer = GetTransformer(scale);
266
267   AliMpIntPair scaledIndices = kTransformer->Scale(indices);
268   AliMpPad pad 
269     = GetSectorSegmentation(scale)->PadByIndices(scaledIndices, warning);
270     
271   return kTransformer->Transform(pad);        
272 }
273
274 //_____________________________________________________________________________
275 AliMpPad AliMpPlaneSegmentation::PadByPosition(const TVector2& position,
276                                                Bool_t warning) const
277 {
278 // Find the pad which corresponds to the given position
279 // ---
280
281   AliMpIntPair scale = GetScale(position);
282   const AliMpTransformer* kTransformer = GetTransformer(scale);
283
284   TVector2 scaledPosition = kTransformer->ITransform(position);  
285   AliMpPad pad 
286     = GetSectorSegmentation(scale)->PadByPosition(scaledPosition, warning);
287   
288   return kTransformer->Transform(pad);        
289 }
290
291 //_____________________________________________________________________________
292 Bool_t AliMpPlaneSegmentation::HasPad(const AliMpIntPair& indices) const
293 {
294 // Does the pad located by <indices> exists ?
295 // ---
296
297   AliMpIntPair scale = GetScale(indices);
298   const AliMpTransformer* kTransformer = GetTransformer(scale);
299
300   AliMpIntPair scaledIndices = kTransformer->Scale(indices);
301
302   return GetSectorSegmentation(scale)->HasPad(scaledIndices);
303 }
304
305 //_____________________________________________________________________________
306 Int_t AliMpPlaneSegmentation::Zone(const AliMpPad& pad, Bool_t warning) const
307 {
308 // Returns the zone index of the zone containing the specified pad.
309 // This zone index is different from the zone ID,
310 // as it is unique for each pad dimensions.
311 // It is composed in this way:
312 //   sectorID*100 + zoneID*10 + specific index 
313 // Where sectorID = 0,1 for front/back sector.
314 // Specific index is present only for zones containing special motifs.
315 // ---
316
317   if (!pad.IsValid()) {
318     if (warning) Warning("Zone(AliMpPad)", "Invalid pad");
319     return 0;
320   }  
321
322   AliMpIntPair scale = GetScale(pad.GetIndices());  
323   const AliMpTransformer* kTransformer = GetTransformer(scale);
324
325   AliMpPad scaledPad = kTransformer->ITransform(pad);
326   
327   AliMpSectorSegmentation* segmentation = GetSectorSegmentation(scale);
328   Int_t zoneID = segmentation->Zone(scaledPad, warning);
329   
330   // Distinguish zones from front/back sector
331   // For back sector - add 10
332   if (segmentation == fBackSectorSegmentation) zoneID += 100;
333
334   return zoneID;
335 }  
336
337 //_____________________________________________________________________________
338 TVector2 
339 AliMpPlaneSegmentation::PadDimensions(Int_t zone, Bool_t warning) const
340 {
341 // Returns the pad dimensions for the zone with the specified zone index.
342 // ---
343
344   if (zone < 100)
345     return fFrontSectorSegmentation->PadDimensions(zone, warning);
346   else  
347     return fBackSectorSegmentation->PadDimensions(zone - 100, warning);
348 }  
349
350 //_____________________________________________________________________________
351 Bool_t AliMpPlaneSegmentation::CircleTest(const AliMpIntPair& indices) const
352 {
353 // Verifies that all methods for retrieving pads are consistents between them.
354 // Returns true if the pad with specified indices was found and verified,
355 // false otherwise.
356 // ---
357
358   if (!HasPad(indices)) return false;
359
360   // Verify the indice->location->position->indice way
361   AliMpIntPair location = PadByIndices(indices).GetLocation();
362   TVector2 position = PadByLocation(location).Position();
363   AliMpIntPair retIndices = PadByPosition(position).GetIndices();
364     
365   if (retIndices != indices) {
366     cout << "Pad " << indices << " lead to inconsistency" << endl;
367     cout << "in indice->location->position->indice way..." << endl;
368     cout << "starting from " << indices << "-->" << location << "-->" 
369          << '(' << position.X() << ',' << position.Y() << ')'
370          << " and retIndices: " << retIndices << endl;
371   }
372     
373     
374   // Verify the indice->position->location->indice way    
375   position = PadByIndices(indices).Position();
376   location = PadByPosition(position).GetLocation();
377   retIndices = PadByLocation(location).GetIndices();
378
379   if (retIndices != indices) {
380     cout << "Pad " << indices << " lead to inconsistency" << endl;
381     cout << "in indice->position->location->indice way..." <<endl;
382     cout << "starting from " << indices 
383          << " and retIndices: " << retIndices << endl;
384   }
385   
386   return true;
387
388 }
389
390 //_____________________________________________________________________________
391 Int_t AliMpPlaneSegmentation::GetNofTransformers() const
392 {
393 // Returns number of transformers.
394 // ---
395
396 #ifdef WITH_STL
397   return fTransformers.size();
398 #endif
399
400 #ifdef WITH_ROOT
401   return fTransformers.GetEntriesFast();
402 #endif
403 }  
404
405
406 //_____________________________________________________________________________
407 AliMpTransformer* AliMpPlaneSegmentation::GetTransformer(Int_t i) const
408 {
409 // Returns i-th transformer.
410 // ---
411  
412 #ifdef WITH_STL
413   return  fTransformers[i];
414 #endif
415
416 #ifdef WITH_ROOT
417   return  (AliMpTransformer*)fTransformers[i];
418 #endif
419 }     
420
421