]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONContourHandler.cxx
Attempt to fix Coverity defect 14695
[u/mrichter/AliRoot.git] / MUON / AliMUONContourHandler.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
18 ///
19 /// \class AliMUONContourHandler
20 /// 
21 /// Class used to create contours of the muon tracker parts :
22 /// manu, bus patches, detection elements, chambers.
23 ///
24 /// \author Laurent Aphecetche, Subatech
25 ///
26
27 #include "AliMUONContourHandler.h"
28
29 #include "AliCodeTimer.h"
30 #include "AliLog.h"
31 #include "AliMUONContour.h"
32 #include "AliMUONContourMaker.h"
33 #include "AliMUONGeometryDetElement.h"
34 #include "AliMUONGeometryTransformer.h"
35 #include "AliMUONManuContourMaker.h"
36 #include "AliMUONPolygon.h"
37 #include "AliMUONSegment.h"
38 #include "AliMpArea.h"
39 #include "AliMpCDB.h"
40 #include "AliMpDDLStore.h"
41 #include "AliMpDEManager.h"
42 #include "AliMpExMap.h"
43 #include <float.h>
44 #include "Riostream.h"
45 #include "TArrow.h"
46 #include "TCanvas.h"
47 #include "TFile.h"
48 #include "TGeoMatrix.h"
49 #include "TLine.h"
50 #include "TObjArray.h"
51 #include "TPolyLine.h"
52 #include "TSystem.h"
53
54 ///\cond CLASSIMP
55 ClassImp(AliMUONContourHandler)
56 ///\endcond
57
58 //_____________________________________________________________________________
59 AliMUONContourHandler::AliMUONContourHandler(Bool_t explodedView) 
60 : TObject(),
61 fTransformations(0x0),
62 fAllContourMap(0x0),
63 fAllContourArray(0x0)
64 {
65   /// ctor. If explodedView=kTRUE, we generate views that will look 
66   /// fine in 2D (i.e. won't show overlaps of DE as in reality)
67   /// Use kFALSE if you want realistic geometry, though.
68   ///
69   /// IMPORTANT : as we many MUON classes, this one cannot work
70   /// before you've loaded the mapping in memory (see e.g. AliMpCDB::LoadDDLStore)
71   ///
72   
73   if ( explodedView ) 
74   {
75     fTransformations = GenerateTransformations(explodedView);
76   }
77   
78   AliMUONManuContourMaker manuMaker(fTransformations);
79   
80   TObjArray* manus = manuMaker.GenerateManuContours(kTRUE);
81   
82   if (manus)
83   {
84     manus->SetOwner(kFALSE);
85     GenerateAllContours(*manus);
86   }
87   delete manus;
88 }
89
90 //_____________________________________________________________________________
91 AliMUONContourHandler::~AliMUONContourHandler()
92 {
93   /// Dtor
94   delete fTransformations;
95   delete fAllContourMap;
96   delete fAllContourArray;
97 }
98
99 //______________________________________________________________________________
100 Bool_t 
101 AliMUONContourHandler::Adopt(AliMUONContour* contour)
102 {
103   /// Adopt the given contour
104   if ( GetContour(contour->GetName()) ) 
105   {
106     // contour already exists
107     return kFALSE;
108   }
109   fAllContourMap->Add(new TObjString(contour->GetName()),contour);
110   return kTRUE;
111 }
112
113 //______________________________________________________________________________
114 TObjArray*
115 AliMUONContourHandler::CreateContourList(const TObjArray& manuContours)
116 {    
117   /// Create an array of maps of contour names
118   ///
119   /// Assyming that key is something like station#/chamber#/de#/buspatch#/manu#
120   /// the idea here is to put one TMap for each level in mapArray :
121   ///
122   /// mapArray[0].key = station0
123   /// mapArray[0].value = map of strings { station0/chamber0, station0/chamber1 }
124   ///
125   /// Then each entry in mapArray will be converted into a contour by
126   /// merging its children (e.g. station0 contour will be made from the merging
127   /// of station0/chamber0 and station0/chamber1 in the example above).
128   ///
129   
130   AliCodeTimerAuto("",0);
131   
132   Int_t start(0);
133   
134   if ( !fTransformations )
135   {
136     start = 2; // skip chamber and station contours, as we seem to be in 3D
137   }
138     
139   TIter next(&manuContours);
140   AliMUONContour* contour;
141   TObjArray* mapArray = new TObjArray;
142   
143   while ( ( contour = static_cast<AliMUONContour*>(next()) ) )
144   {
145     // Key is something like station#/chamber#/de#/buspatch#/manu#
146     
147     TString key(contour->GetName());
148     TObjArray* s = key.Tokenize("/");
149     for ( Int_t i = start; i < s->GetLast(); ++i ) 
150     {
151       TMap* m = static_cast<TMap*>(mapArray->At(i));
152       if (!m)
153       {
154         m = new TMap;
155         if ( i > mapArray->GetSize() ) mapArray->Expand(i);
156         mapArray->AddAt(m,i);
157       }
158       TString parent;
159       for ( Int_t k = 0; k <= i; ++k )
160       {
161         TObjString* str = static_cast<TObjString*>(s->At(k));
162         parent += str->String();
163         if ( k < i ) parent += "/";
164       }
165       TString child(parent);
166       child += "/";
167       child += static_cast<TObjString*>(s->At(i+1))->String();
168       
169       TObjArray* ma = static_cast<TObjArray*>(m->GetValue(parent.Data()));
170       if (!ma)
171       {
172         ma = new TObjArray;
173         m->Add(new TObjString(parent.Data()),ma);
174       }
175       TPair* p = static_cast<TPair*>(ma->FindObject(child.Data()));
176       if ( !p ) 
177       {
178         ma->Add(new TObjString(child.Data()));
179       }
180     }
181     delete s;
182   }
183   
184   return mapArray;
185 }
186
187 //______________________________________________________________________________
188 void
189 AliMUONContourHandler::GenerateAllContours(const TObjArray& manuContours)
190 {
191   /// From a map of manu contours, generate the compound contours (bp, de, etc...)
192   /// by merging them.
193   /// Note that manuContours should NOT be the owner of its contours,
194   /// as they are adopted by the array returned by this method.
195   
196   AliCodeTimerAuto("",0);
197   
198   // Get the list of contours to create
199   TObjArray* mapArray = CreateContourList(manuContours);
200   
201   fAllContourMap = new TMap(20000,1);
202   fAllContourMap->SetOwnerKeyValue(kTRUE,kTRUE);
203   
204   fAllContourArray = new TObjArray;
205   fAllContourArray->SetOwner(kFALSE); // the map is the real owner of the contours
206   
207   TIter nextContour(&manuContours);  
208   AliMUONContour* contour(0x0);
209   
210   // start by adding the manu contours we begin with
211   while ( ( contour = static_cast<AliMUONContour*>(nextContour()) ) )
212   {
213     fAllContourMap->Add(new TObjString(contour->GetName()),contour);
214   }
215   
216   AliMUONContourMaker maker;
217   
218   for ( Int_t i = mapArray->GetLast(); i >= 1; --i ) 
219     // end at 1 to avoid merging different cathodes together, which
220     // would not work...
221   {
222     TMap* a = static_cast<TMap*>(mapArray->At(i));
223     TIter next3(a);
224     TObjString* str;
225     while ( ( str = static_cast<TObjString*>(next3()) ) )
226     {
227       TObjArray* m = static_cast<TObjArray*>(a->GetValue(str->String().Data()));
228       TIter next4(m);
229       TObjString* k;
230       TObjArray subcontours;
231       subcontours.SetOwner(kFALSE);
232       while ( ( k = static_cast<TObjString*>(next4()) ) )
233       {
234         contour = static_cast<AliMUONContour*>(fAllContourMap->GetValue(k->String().Data()));
235         if ( contour ) 
236         {
237           subcontours.Add(contour);
238         }
239         else
240         {
241           AliError(Form("Did not find contour %s",k->String().Data()))
242           continue;
243         }
244       }
245       
246       contour = maker.MergeContour(subcontours,str->String().Data());
247       
248       bool error(kFALSE);
249       
250       if (!contour)
251       {
252         error=kTRUE;
253         AliError(Form("ERROR : could not merge into %s",str->String().Data()));
254       }
255       else
256       {
257         if ( contour->Area().IsValid() == kFALSE ) 
258         {
259           error=kTRUE;
260           AliError(Form("ERROR : area of contour %s is invalid",str->String().Data()));
261         }
262       }
263       
264       if (!error)
265       {
266         fAllContourMap->Add(new TObjString(str->String().Data()),contour);
267         fAllContourArray->Add(contour);
268       }
269     }
270   }
271 }
272
273 //_____________________________________________________________________________
274 AliMpExMap* 
275 AliMUONContourHandler::GenerateTransformations(Bool_t exploded)
276 {
277   /// Generate geometric transformations to be used to compute the contours
278   /// If exploded=kFALSE then we generate real transformations, otherwise
279   /// we generate tweaked ones that look fine on screen.
280   
281   AliCodeTimerAuto("",0);
282   
283   AliMUONGeometryTransformer transformer;
284   Bool_t ok = transformer.LoadGeometryData("transform.dat");
285   //  transformer.LoadGeometryData("geometry.root"); //FIXME: add a protection if geometry.root file does not exist
286   if (!ok)
287   {
288     cout << "ERROR : cannot get geometry !" << endl;
289     return 0x0;
290   }
291   AliMpExMap* transformations = new AliMpExMap;
292   AliMpDEIterator deIt;
293   deIt.First();
294   while ( !deIt.IsDone() )
295   {
296     Int_t detElemId = deIt.CurrentDEId();
297     const AliMUONGeometryDetElement* de = transformer.GetDetElement(detElemId);
298     
299     TGeoHMatrix* matrix = static_cast<TGeoHMatrix*>(de->GetGlobalTransformation()->Clone());
300
301     if (exploded)
302     {
303       Double_t* translation = matrix->GetTranslation();
304       Double_t xscale = 1.0;
305       Double_t yscale = 1.5;
306       Double_t shift = 5.0; // cm
307           
308       if ( AliMpDEManager::GetStationType(detElemId) == AliMp::kStation345 ) 
309       {
310         translation[0] *= xscale;
311         translation[1] *= yscale; 
312       }
313       else
314       {
315         Double_t xshift[] = { shift, -shift, -shift, shift };
316         Double_t yshift[] = { shift, shift, -shift, -shift };
317         Int_t ishift = detElemId % 100;
318         
319         translation[0] += xshift[ishift];
320         translation[1] += yshift[ishift];
321       }
322       matrix->SetTranslation(translation);
323     }
324     transformations->Add(detElemId,matrix);
325     deIt.Next();
326   }
327   return transformations;
328 }
329
330 //_____________________________________________________________________________
331 AliMUONContour* 
332 AliMUONContourHandler::GetContour(const char* contourname) const
333 {
334   /// Get a given contour
335   return static_cast<AliMUONContour*>(fAllContourMap->GetValue(contourname));
336 }
337
338 //_____________________________________________________________________________
339 void
340 AliMUONContourHandler::Print(Option_t* opt) const
341 {
342   /// printout
343   
344   if ( ! fAllContourMap )  return;
345   
346   cout << Form("Contour map : collisions = %5.3f size = %d capacity = %d", 
347                fAllContourMap->AverageCollisions(),
348                fAllContourMap->GetSize(),
349                fAllContourMap->Capacity()) << endl;
350
351   TString sopt(opt);
352   sopt.ToUpper();
353   
354   if ( sopt.Contains("ALL") || sopt.Contains("FULL") )
355   {
356     fAllContourMap->Print();
357   }
358 }