]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONContourHandler.cxx
Fixing previous change, empty argument macros not supported everywhere
[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   fTransformations = GenerateTransformations(explodedView);
74   AliMUONManuContourMaker manuMaker(fTransformations);
75   TObjArray* manus = manuMaker.GenerateManuContours(kTRUE);
76   if (manus)
77   {
78     manus->SetOwner(kFALSE);
79     GenerateAllContours(*manus);
80   }
81   delete manus;
82 }
83
84 //_____________________________________________________________________________
85 AliMUONContourHandler::~AliMUONContourHandler()
86 {
87   /// Dtor
88   delete fAllContourMap;
89   delete fAllContourArray;
90 }
91
92 //______________________________________________________________________________
93 Bool_t 
94 AliMUONContourHandler::Adopt(AliMUONContour* contour)
95 {
96   /// Adopt the given contour
97   if ( GetContour(contour->GetName()) ) 
98   {
99     // contour already exists
100     return kFALSE;
101   }
102   fAllContourMap->Add(new TObjString(contour->GetName()),contour);
103   return kTRUE;
104 }
105
106 //______________________________________________________________________________
107 TObjArray*
108 AliMUONContourHandler::CreateContourList(const TObjArray& manuContours)
109 {    
110   /// Create an array of maps of contour names
111   ///
112   /// Assyming that key is something like station#/chamber#/de#/buspatch#/manu#
113   /// the idea here is to put one TMap for each level in mapArray :
114   ///
115   /// mapArray[0].key = station0
116   /// mapArray[0].value = map of strings { station0/chamber0, station0/chamber1 }
117   ///
118   /// Then each entry in mapArray will be converted into a contour by
119   /// merging its children (e.g. station0 contour will be made from the merging
120   /// of station0/chamber0 and station0/chamber1 in the example above).
121   ///
122   
123   AliCodeTimerAuto("",0);
124   
125   TIter next(&manuContours);
126   AliMUONContour* contour;
127   TObjArray* mapArray = new TObjArray;
128   
129   while ( ( contour = static_cast<AliMUONContour*>(next()) ) )
130   {
131     // Key is something like station#/chamber#/de#/buspatch#/manu#
132     
133     TString key(contour->GetName());
134     TObjArray* s = key.Tokenize("/");
135     for ( Int_t i = 0; i < s->GetLast(); ++i ) 
136     {
137       TMap* m = static_cast<TMap*>(mapArray->At(i));
138       if (!m)
139       {
140         m = new TMap;
141         if ( i > mapArray->GetSize() ) mapArray->Expand(i);
142         mapArray->AddAt(m,i);
143       }
144       TString parent;
145       for ( Int_t k = 0; k <= i; ++k )
146       {
147         TObjString* str = static_cast<TObjString*>(s->At(k));
148         parent += str->String();
149         if ( k < i ) parent += "/";
150       }
151       TString child(parent);
152       child += "/";
153       child += static_cast<TObjString*>(s->At(i+1))->String();
154       
155       TObjArray* ma = static_cast<TObjArray*>(m->GetValue(parent.Data()));
156       if (!ma)
157       {
158         ma = new TObjArray;
159         m->Add(new TObjString(parent.Data()),ma);
160       }
161       TPair* p = static_cast<TPair*>(ma->FindObject(child.Data()));
162       if ( !p ) 
163       {
164         ma->Add(new TObjString(child.Data()));
165       }
166     }
167     delete s;
168   }
169   
170   return mapArray;
171 }
172
173 //______________________________________________________________________________
174 void
175 AliMUONContourHandler::GenerateAllContours(const TObjArray& manuContours)
176 {
177   /// From a map of manu contours, generate the compound contours (bp, de, etc...)
178   /// by merging them.
179   /// Note that manuContours should NOT be the owner of its contours,
180   /// as they are adopted by the array returned by this method.
181   
182   AliCodeTimerAuto("",0);
183   
184   // Get the list of contours to create
185   TObjArray* mapArray = CreateContourList(manuContours);
186   
187   // Now loop over the mapArray to actually create the contours
188   TIter next2(mapArray,kIterBackward);
189   
190   fAllContourMap = new TMap(20000,1);
191   fAllContourMap->SetOwnerKeyValue(kTRUE,kTRUE);
192   
193   fAllContourArray = new TObjArray;
194   fAllContourArray->SetOwner(kFALSE); // the map is the real owner of the contours
195   
196   TIter nextContour(&manuContours);  
197   AliMUONContour* contour(0x0);
198   
199   // start by adding the manu contours we begin with
200   while ( ( contour = static_cast<AliMUONContour*>(nextContour()) ) )
201   {
202     fAllContourMap->Add(new TObjString(contour->GetName()),contour);
203   }
204   
205   AliMUONContourMaker maker;
206   
207   for ( Int_t i = mapArray->GetLast(); i >= 1; --i ) 
208     // end at 1 to avoid merging different cathodes together, which
209     // would not work...
210   {
211     TMap* a = static_cast<TMap*>(mapArray->At(i));
212     TIter next3(a);
213     TObjString* str;
214     while ( ( str = static_cast<TObjString*>(next3()) ) )
215     {
216       TObjArray* m = static_cast<TObjArray*>(a->GetValue(str->String().Data()));
217       TIter next4(m);
218       TObjString* k;
219       TObjArray subcontours;
220       subcontours.SetOwner(kFALSE);
221       while ( ( k = static_cast<TObjString*>(next4()) ) )
222       {
223         contour = static_cast<AliMUONContour*>(fAllContourMap->GetValue(k->String().Data()));
224         if ( contour ) 
225         {
226           subcontours.Add(contour);
227         }
228         else
229         {
230           AliError(Form("Did not find contour %s",k->String().Data()))
231           continue;
232         }
233       }
234       
235       contour = maker.MergeContour(subcontours,str->String().Data());
236       
237       bool error(kFALSE);
238       
239       if (!contour)
240       {
241         error=kTRUE;
242         AliError(Form("ERROR : could not merge into %s",str->String().Data()));
243       }
244       else
245       {
246         if ( contour->Area().IsValid() == kFALSE ) 
247         {
248           error=kTRUE;
249           AliError(Form("ERROR : area of contour %s is invalid",str->String().Data()));
250         }
251       }
252       
253       if (!error)
254       {
255         fAllContourMap->Add(new TObjString(str->String().Data()),contour);
256         fAllContourArray->Add(contour);
257       }
258     }
259   }
260 }
261
262 //_____________________________________________________________________________
263 AliMpExMap* 
264 AliMUONContourHandler::GenerateTransformations(Bool_t exploded)
265 {
266   /// Generate geometric transformations to be used to compute the contours
267   /// If exploded=kFALSE then we generate real transformations, otherwise
268   /// we generate tweaked ones that look fine on screen.
269   
270   AliCodeTimerAuto("",0);
271   
272   AliMUONGeometryTransformer transformer;
273   Bool_t ok = transformer.LoadGeometryData("transform.dat");
274   //  transformer.LoadGeometryData("geometry.root"); //FIXME: add a protection if geometry.root file does not exist
275   if (!ok)
276   {
277     cout << "ERROR : cannot get geometry !" << endl;
278     return 0x0;
279   }
280   AliMpExMap* transformations = new AliMpExMap;
281   AliMpDEIterator deIt;
282   deIt.First();
283   while ( !deIt.IsDone() )
284   {
285     Int_t detElemId = deIt.CurrentDEId();
286     const AliMUONGeometryDetElement* de = transformer.GetDetElement(detElemId);
287     
288     TGeoHMatrix* matrix = static_cast<TGeoHMatrix*>(de->GetGlobalTransformation()->Clone());
289
290     if (exploded)
291     {
292       Double_t* translation = matrix->GetTranslation();
293       Double_t xscale = 1.0;
294       Double_t yscale = 1.5;
295       Double_t shift = 5.0; // cm
296           
297       if ( AliMpDEManager::GetStationType(detElemId) == AliMp::kStation345 ) 
298       {
299         translation[0] *= xscale;
300         translation[1] *= yscale; 
301       }
302       else
303       {
304         Double_t xshift[] = { shift, -shift, -shift, shift };
305         Double_t yshift[] = { shift, shift, -shift, -shift };
306         Int_t ishift = detElemId % 100;
307         
308         translation[0] += xshift[ishift];
309         translation[1] += yshift[ishift];
310       }
311       matrix->SetTranslation(translation);
312     }
313     transformations->Add(detElemId,matrix);
314     deIt.Next();
315   }
316   return transformations;
317 }
318
319 //_____________________________________________________________________________
320 AliMUONContour* 
321 AliMUONContourHandler::GetContour(const char* contourname) const
322 {
323   /// Get a given contour
324   return static_cast<AliMUONContour*>(fAllContourMap->GetValue(contourname));
325 }
326
327 //_____________________________________________________________________________
328 void
329 AliMUONContourHandler::Print(Option_t* /*opt*/) const
330 {
331   /// printout
332   
333   if ( fAllContourMap ) 
334   {
335     cout << Form("Contour map : collisions = %5.3f size = %d capacity = %d", 
336                  fAllContourMap->AverageCollisions(),
337                  fAllContourMap->GetSize(),
338                  fAllContourMap->Capacity()) << endl;
339   }
340 }