]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDAlignFaker.cxx
Added new 'status' histogram to the DQM data maker and checker. This histogram
[u/mrichter/AliRoot.git] / FMD / AliFMDAlignFaker.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 /* $Id$ */
16 /** 
17  * @file    AliFMDAlignFaker.cxx
18  * @author  Christian Holm Christensen <cholm@nbi.dk>
19  * @date    Sun Mar 26 17:57:55 2006
20  * @brief   Implementation of AliFMDAlignFaker 
21  */
22 //____________________________________________________________________
23 //
24 //  Class 
25 //  to 
26 //  make 
27 //  fake 
28 //  alignment
29 //  parameters 
30 //
31 //____________________________________________________________________
32 //                                                                          
33 // Forward Multiplicity Detector based on Silicon wafers. 
34 //
35 // This task creates fake alignment. Which alignment, depends on
36 // the bit mask passed to the constructor, or added by `AddAlign'.
37 //
38 // The default is to write all alignment parameters to a local
39 // storage `local://cdb' which is a directory in the current
40 // directory. 
41 //                                                       
42 #include "AliFMDDebug.h"                // ALIFMDDEBUG_H ALILOG_H
43 #include "AliFMDAlignFaker.h"      // ALIFMDALIGNFAKER_H
44 #include <AliCDBManager.h>         // ALICDBMANAGER_H
45 #include <AliCDBStorage.h>         // ALICDBSTORAGE_H
46 #include <AliCDBEntry.h>           // ALICDBMANAGER_H
47 // #include <AliAlignObj.h>
48 #include <AliAlignObjParams.h>
49 // #include <Riostream.h>
50 #include <TSystem.h>
51 // #include <TMath.h>
52 #include <TRandom.h>
53 #include <TClonesArray.h>
54 #include <TString.h>
55 #include <TFile.h>
56 #include <TGeoManager.h>
57 #include <TGeoNode.h>
58 // #include <TGeoVolume.h>
59 #include <TROOT.h>
60 #include <TClass.h>
61
62 //====================================================================
63 ClassImp(AliFMDAlignFaker)
64 #if 0
65   ; // This is here to keep Emacs for indenting the next line
66 #endif
67
68 //____________________________________________________________________
69 AliFMDAlignFaker::AliFMDAlignFaker(Int_t mask, const char* geo, 
70                                    const char* loc) 
71   : TTask(geo, loc),
72     fMask(mask),
73     fSensorTransMin(0,0,0),
74     fSensorTransMax(0,0,0),
75     fSensorRotMin(0,0,0),
76     fSensorRotMax(0,0,0),
77     fHalfTransMin(0,0,0),
78     fHalfTransMax(0,0,0),
79     fHalfRotMin(0,0,0),
80     fHalfRotMax(0,0,0),
81     fRunMin(0),
82     fRunMax(AliCDBRunRange::Infinity()), 
83     fArray(0),
84     fComment("")
85 {
86   // Default constructor 
87   if (!loc) { 
88     AliCDBStorage* storage = AliCDBManager::Instance()->GetDefaultStorage();
89     if (!storage) AliFatal("Default Storage not set");
90     const TString& uri = storage->GetURI();
91     fTitle = uri;
92   }
93   SetSensorDisplacement();
94   SetSensorRotation();
95   SetHalfDisplacement();
96   SetHalfRotation();
97   SetComment();
98 }
99
100 //__________________________________________________________________
101 void
102 AliFMDAlignFaker::SetSensorDisplacement(Double_t x1, Double_t y1, Double_t z1,
103                                         Double_t x2, Double_t y2, Double_t z2)
104 {
105   // Set sensor displacement (unit is centimeters)
106   fSensorTransMin.SetXYZ(x1, y1, z1);
107   fSensorTransMax.SetXYZ(x2, y2, z2);
108 }
109
110 //__________________________________________________________________
111 void
112 AliFMDAlignFaker::SetSensorRotation(Double_t x1, Double_t y1, Double_t z1,
113                                     Double_t x2, Double_t y2, Double_t z2)
114 {
115   // Set sensor rotations (unit is degrees)
116   fSensorRotMin.SetXYZ(x1, y1, z1);
117   fSensorRotMax.SetXYZ(x2, y2, z2);
118 }
119
120 //__________________________________________________________________
121 void
122 AliFMDAlignFaker::SetHalfDisplacement(Double_t x1, Double_t y1, Double_t z1,
123                                       Double_t x2, Double_t y2, Double_t z2)
124 {
125   // Set half ring/cone displacement (unit is centimeters)
126   fHalfTransMin.SetXYZ(x1, y1, z1);
127   fHalfTransMax.SetXYZ(x2, y2, z2);
128 }
129
130 //__________________________________________________________________
131 void
132 AliFMDAlignFaker::SetHalfRotation(Double_t x1, Double_t y1, Double_t z1,
133                                   Double_t x2, Double_t y2, Double_t z2)
134 {
135   // Set half ring/cone rotations (unit is degrees)
136   fHalfRotMin.SetXYZ(x1, y1, z1);
137   fHalfRotMax.SetXYZ(x2, y2, z2);
138 }
139
140 //__________________________________________________________________
141 #define IS_NODE_HALF(name) \
142   (name[0] == 'F' && name[2] == 'M' && (name[3] == 'T' || name[3] == 'B'))
143 #define IS_NODE_SENSOR(name) \
144   (name[0] == 'F' && name[2] == 'S' && name[3] == 'E')
145
146 //__________________________________________________________________
147 Bool_t
148 AliFMDAlignFaker::GetGeometry(Bool_t toCdb, const TString& storage)
149 {
150   // 
151   // Get the geometry
152   // 
153   // Parameters:
154   //    toCdb   Whether to store in CDB
155   //    storage Storage element to use 
156   // 
157   // Return:
158   //    true on success 
159   //
160   if (!toCdb) { 
161      //load geom from default CDB storage
162     AliGeomManager::LoadGeometry(); 
163     return kTRUE;
164   }
165   if(!storage.BeginsWith("local://") && 
166      !storage.BeginsWith("alien://")) {
167     AliErrorClass(Form("STORAGE=\"%s\" invalid. Exiting\n", storage.Data()));
168     return kFALSE;
169   }
170
171   AliCDBManager* cdb   = AliCDBManager::Instance();
172   AliCDBStorage* store = cdb->GetStorage(storage.Data());
173   if(!store){
174     AliErrorClass(Form("Unable to open storage %s\n", storage.Data()));
175     return kFALSE;
176   }
177
178   AliCDBPath   path("GRP","Geometry","Data");
179   AliCDBEntry* entry = store->Get(path.GetPath(),cdb->GetRun());
180   if(!entry) {
181     AliErrorClass("Could not get the specified CDB entry!");
182     return kFALSE;
183   }
184   
185
186   entry->SetOwner(0);
187   TGeoManager* geom = static_cast<TGeoManager*>(entry->GetObject());
188   AliGeomManager::SetGeometry(geom);
189   return kTRUE;
190 }
191   
192 //__________________________________________________________________
193 void
194 AliFMDAlignFaker::Exec(Option_t*)
195 {
196   // Make the objects. 
197
198   // Get geometry 
199   if (!gGeoManager) {
200     if (!TGeoManager::Import(GetName())) {
201       AliFatal(Form("Failed to import geometry from %s", GetName()));
202       return;
203     }
204   }
205   // Get top volume 
206   TGeoVolume* topVolume = gGeoManager->GetTopVolume();
207   if (!topVolume) {
208     AliFatal("No top-level volume defined");
209     return;
210   }
211   // Make container of transforms 
212   if (!fArray) fArray = new TClonesArray("AliAlignObjParams");
213   fArray->Clear();
214   
215   // Make an iterator
216   TGeoIterator next(topVolume);
217   next.SetTopName(Form("/%s_1", topVolume->GetName()));
218   TGeoNode* node = 0;
219   
220   Char_t currentDet  = '\0';
221   Char_t currentHalf = '\0';
222   // Loop over all entries in geometry to find our nodes. 
223   while ((node = static_cast<TGeoNode*>(next()))) {
224     const char* name =  node->GetName();
225     if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
226         !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors))) 
227       continue;
228
229     TString path, alignName;
230     next.GetPath(path);
231     Int_t id = node->GetVolume()->GetNumber();
232     if (IS_NODE_HALF(name))   { 
233       currentDet    = name[1];
234       currentHalf   = name[3];
235       alignName     = Form("FMD/FMD%c_%c", currentDet, currentHalf);
236     }
237     if (IS_NODE_SENSOR(name)) {
238       Char_t ring      = name[1];
239       Int_t  lvl       = next.GetLevel();
240       TGeoNode* parent = next.GetNode(lvl-1);
241       Int_t     copy   = parent->GetNumber();
242       alignName        = Form("FMD/FMD%c_%c/FMD%c_%02d", 
243                           currentDet, currentHalf, ring, copy);
244     }
245     if (alignName.IsNull()) continue;
246     if (!gGeoManager->GetAlignableEntry(alignName.Data())) {
247       AliWarning(Form("No alignable entry for %s, using path %s",
248                       alignName.Data(), path.Data()));
249       alignName = path;
250     }
251     AliFMDDebug(1, ("Making alignment for %s -> %s (%d)", 
252                      alignName.Data(), path.Data(), id));
253     if (IS_NODE_HALF(name))   MakeAlignHalf(alignName, id);
254     if (IS_NODE_SENSOR(name)) MakeAlignSensor(alignName, id);
255 #if 0    
256     if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
257         !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors))) 
258       continue;
259
260     // Get the path 
261     TString path(Form("/%s", gGeoManager->GetNode(0)->GetName()));
262     Int_t nLevel = next.GetLevel();
263     for (Int_t lvl = 0; lvl <= nLevel; lvl++) {
264       TGeoNode* p = next.GetNode(lvl);
265       if (!p) {
266         if (lvl != 0)
267           AliWarning(Form("No node at level %d in path %s",lvl,path.Data()));
268         continue;
269       }
270       if (!path.IsNull()) path.Append("/");
271       path.Append(p->GetName());
272     }
273     Int_t id = node->GetVolume()->GetNumber();
274     if (IS_NODE_HALF(name))   MakeAlignHalf(path, id);
275     if (IS_NODE_SENSOR(name)) MakeAlignSensor(path, id);
276 #endif 
277   }
278
279   TString t(GetTitle());
280   if (t.IsNull() || t.Contains("local://") || t.Contains("alien://")) 
281     WriteToCDB();
282   else 
283     WriteToFile();
284 }
285   
286 //__________________________________________________________________
287 Bool_t
288 AliFMDAlignFaker::MakeAlign(const TString& path, Int_t id, 
289                             Double_t transX, Double_t transY, Double_t transZ,
290                             Double_t rotX,   Double_t rotY, Double_t rotZ)
291 {
292   // make alignment for a path 
293   // Params: 
294   //   path      Path to node 
295   //   id        Volume number 
296   //   transX    Translation in X
297   //   transZ    Translation in Y
298   //   transZ    Translation in Z
299   //   rotX      Rotation about X-axis 
300   //   rotY      Rotation about Y-axis
301   //   rotZ      Rotation about Z-axis 
302   AliFMDDebug(3, ("Make alignment for %s (volume %d): (%f,%f,%f) (%f,%f,%f)", 
303                    path.Data(), id, transX, transY, transZ, rotX, rotY, rotZ));
304   Int_t nAlign = fArray->GetEntries();
305   id = 0;
306   AliAlignObjParams* obj = 
307     new ((*fArray)[nAlign]) AliAlignObjParams(path.Data(),
308                                               id,0,0,0,0,0,0,kTRUE);
309   if (!obj) {
310     AliError(Form("Failed to create alignment object for %s", path.Data()));
311     return kFALSE;
312   }
313   if (!obj->SetLocalPars(transX, transY, transZ, rotX, rotY, rotZ)) {
314     AliError(Form("Failed to set local transforms on %s", path.Data()));
315     return kTRUE;
316   }
317   return kTRUE;
318 }
319
320 //__________________________________________________________________
321 Bool_t
322 AliFMDAlignFaker::MakeAlignHalf(const TString& path, Int_t id)
323 {
324   // Make alignment of a half ring/cone 
325   AliFMDDebug(15, ("Make alignment for half-ring/cone %s", path.Data()));
326   Double_t transX = gRandom->Uniform(fHalfTransMin.X(), fHalfTransMax.X());
327   Double_t transY = gRandom->Uniform(fHalfTransMin.Y(), fHalfTransMax.Y());
328   Double_t transZ = gRandom->Uniform(fHalfTransMin.Z(), fHalfTransMax.Z());
329   Double_t rotX   = gRandom->Uniform(fHalfRotMin.X(),   fHalfRotMax.X());
330   Double_t rotY   = gRandom->Uniform(fHalfRotMin.Y(),   fHalfRotMax.Y());
331   Double_t rotZ   = gRandom->Uniform(fHalfRotMin.Z(),   fHalfRotMax.Z());
332   return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
333 }
334
335   
336 //__________________________________________________________________
337 Bool_t
338 AliFMDAlignFaker::MakeAlignSensor(const TString& path, Int_t id)
339 {
340   // Make alignment of a sensor 
341   AliFMDDebug(15, ("Make alignment for sensor %s", path.Data()));
342   Double_t transX = gRandom->Uniform(fSensorTransMin.X(), fSensorTransMax.X());
343   Double_t transY = gRandom->Uniform(fSensorTransMin.Y(), fSensorTransMax.Y());
344   Double_t transZ = gRandom->Uniform(fSensorTransMin.Z(), fSensorTransMax.Z());
345   Double_t rotX   = gRandom->Uniform(fSensorRotMin.X(),   fSensorRotMax.X());
346   Double_t rotY   = gRandom->Uniform(fSensorRotMin.Y(),   fSensorRotMax.Y());
347   Double_t rotZ   = gRandom->Uniform(fSensorRotMin.Z(),   fSensorRotMax.Z());
348   return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
349 }
350
351 //__________________________________________________________________
352 void
353 AliFMDAlignFaker::WriteToCDB()
354 {
355   // Make the objects. 
356   const char* loc = GetTitle();
357   AliCDBManager*     cdb  = AliCDBManager::Instance();
358   AliCDBStorage*  storage = cdb->GetStorage(!loc ? "" : loc);
359   AliCDBMetaData*    meta = new AliCDBMetaData; 
360   meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data()); 
361   meta->SetAliRootVersion(gROOT->GetVersion()); 
362   meta->SetBeamPeriod(1); 
363   meta->SetComment(fComment.Data());
364
365   AliCDBId id("FMD/Align/Data", fRunMin, fRunMax);
366   storage->Put(fArray, id, meta);
367 }
368
369 //__________________________________________________________________
370 void
371 AliFMDAlignFaker::WriteToFile()
372 {
373   // Write to a local file 
374   TFile* file = TFile::Open(GetTitle(), "RECREATE");
375   if (!file) {
376     AliFatal(Form("Failed to open file '%s' for output", GetTitle()));
377     return;
378   }
379   file->cd();
380   fArray->Write("FMDAlignment",TObject::kSingleKey);
381   file->Write();
382   file->Close();
383 }
384
385   
386   
387 //____________________________________________________________________
388 //
389 // EOF
390 //