]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDAlignFaker.cxx
Bug fixed (Christian)
[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 /** @file    AliFMDAlignFaker.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Sun Mar 26 17:57:55 2006
19     @brief   Implementation of AliFMDAlignFaker 
20 */
21 //____________________________________________________________________
22 //
23 //  Class 
24 //  to 
25 //  make 
26 //  fake 
27 //  alignment
28 //  parameters 
29 //
30 //____________________________________________________________________
31 //                                                                          
32 // Forward Multiplicity Detector based on Silicon wafers. 
33 //
34 // This task creates fake alignment. Which alignment, depends on
35 // the bit mask passed to the constructor, or added by `AddAlign'.
36 //
37 // The default is to write all alignment parameters to a local
38 // storage `local://cdb' which is a directory in the current
39 // directory. 
40 //                                                       
41 #include "AliLog.h"                // ALILOG_H
42 #include "AliFMDAlignFaker.h"      // ALIFMDALIGNFAKER_H
43 #include <AliCDBManager.h>         // ALICDBMANAGER_H
44 #include <AliCDBEntry.h>           // ALICDBMANAGER_H
45 // #include <AliAlignObj.h>
46 #include <AliAlignObjAngles.h>
47 // #include <Riostream.h>
48 // #include <TSystem.h>
49 // #include <TMath.h>
50 #include <TRandom.h>
51 #include <TClonesArray.h>
52 #include <TString.h>
53 #include <TFile.h>
54 #include <TGeoManager.h>
55 #include <TGeoNode.h>
56 // #include <TGeoVolume.h>
57 // #include <TROOT.h>
58
59 //====================================================================
60 ClassImp(AliFMDAlignFaker)
61 #if 0
62   ; // This is here to keep Emacs for indenting the next line
63 #endif
64
65 //____________________________________________________________________
66 AliFMDAlignFaker::AliFMDAlignFaker(Int_t mask, const char* geo, 
67                                    const char* loc) 
68   : TTask(geo, loc),
69     fMask(mask),
70     fSensorTransMin(0,0,0),
71     fSensorTransMax(0,0,0),
72     fSensorRotMin(0,0,0),
73     fSensorRotMax(0,0,0),
74     fHalfTransMin(0,0,0),
75     fHalfTransMax(0,0,0),
76     fHalfRotMin(0,0,0),
77     fHalfRotMax(0,0,0),
78     fRunMin(0),
79     fRunMax(10), 
80     fArray(0),
81     fComment("")
82 {
83   // Default constructor 
84   SetSensorDisplacement();
85   SetSensorRotation();
86   SetHalfDisplacement();
87   SetHalfRotation();
88   SetComment();
89 }
90
91 //__________________________________________________________________
92 void
93 AliFMDAlignFaker::SetSensorDisplacement(Double_t x1, Double_t y1, Double_t z1,
94                                         Double_t x2, Double_t y2, Double_t z2)
95 {
96   // Set sensor displacement (unit is centimeters)
97   fSensorTransMin.SetXYZ(x1, y1, z1);
98   fSensorTransMax.SetXYZ(x2, y2, z2);
99 }
100
101 //__________________________________________________________________
102 void
103 AliFMDAlignFaker::SetSensorRotation(Double_t x1, Double_t y1, Double_t z1,
104                                     Double_t x2, Double_t y2, Double_t z2)
105 {
106   // Set sensor rotations (unit is degrees)
107   fSensorRotMin.SetXYZ(x1, y1, z1);
108   fSensorRotMax.SetXYZ(x2, y2, z2);
109 }
110
111 //__________________________________________________________________
112 void
113 AliFMDAlignFaker::SetHalfDisplacement(Double_t x1, Double_t y1, Double_t z1,
114                                       Double_t x2, Double_t y2, Double_t z2)
115 {
116   // Set half ring/cone displacement (unit is centimeters)
117   fHalfTransMin.SetXYZ(x1, y1, z1);
118   fHalfTransMax.SetXYZ(x2, y2, z2);
119 }
120
121 //__________________________________________________________________
122 void
123 AliFMDAlignFaker::SetHalfRotation(Double_t x1, Double_t y1, Double_t z1,
124                                   Double_t x2, Double_t y2, Double_t z2)
125 {
126   // Set half ring/cone rotations (unit is degrees)
127   fHalfRotMin.SetXYZ(x1, y1, z1);
128   fHalfRotMax.SetXYZ(x2, y2, z2);
129 }
130
131 //__________________________________________________________________
132 #define IS_NODE_HALF(name) \
133   (name[0] == 'F' && name[2] == 'M' && (name[3] == 'T' || name[3] == 'B'))
134 #define IS_NODE_SENSOR(name) \
135   (name[0] == 'F' && name[2] == 'S' && name[3] == 'E')
136
137 //__________________________________________________________________
138 void
139 AliFMDAlignFaker::Exec(Option_t*)
140 {
141   // Make the objects. 
142
143   // Get geometry 
144   if (!gGeoManager) {
145     if (!TGeoManager::Import(GetName())) {
146       AliFatal(Form("Failed to import geometry from %s", GetName()));
147       return;
148     }
149   }
150   // Get top volume 
151   TGeoVolume* topVolume = gGeoManager->GetTopVolume();
152   if (!topVolume) {
153     AliFatal("No top-level volume defined");
154     return;
155   }
156   // Make container of transforms 
157   if (!fArray) fArray = new TClonesArray("AliAlignObjAngles");
158   fArray->Clear();
159   
160   // Make an iterator
161   TGeoIterator next(topVolume);
162   TGeoNode* node = 0;
163
164   // Loop over all entries in geometry to find our nodes. 
165   while ((node = static_cast<TGeoNode*>(next()))) {
166     const char* name =  node->GetName();
167     if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
168         !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors))) 
169       continue;
170     
171     // Get the path 
172     TString path(Form("/%s", gGeoManager->GetNode(0)->GetName()));
173     Int_t nLevel = next.GetLevel();
174     for (Int_t lvl = 0; lvl <= nLevel; lvl++) {
175       TGeoNode* p = next.GetNode(lvl);
176       if (!p) {
177         if (lvl != 0)
178           AliWarning(Form("No node at level %d in path %s",lvl,path.Data()));
179         continue;
180       }
181       if (!path.IsNull()) path.Append("/");
182       path.Append(p->GetName());
183     }
184     Int_t id = node->GetVolume()->GetNumber();
185     if (IS_NODE_HALF(name))   MakeAlignHalf(path, id);
186     if (IS_NODE_SENSOR(name)) MakeAlignSensor(path, id);
187   }
188
189   TString t(GetTitle());
190   if (t.IsNull() || t.Contains("local://") || t.Contains("alien://")) 
191     WriteToCDB();
192   else 
193     WriteToFile();
194 }
195   
196 //__________________________________________________________________
197 Bool_t
198 AliFMDAlignFaker::MakeAlign(const TString& path, Int_t id, 
199                             Double_t transX, Double_t transY, Double_t transZ,
200                             Double_t rotX,   Double_t rotY, Double_t rotZ)
201 {
202   // make alignment for a path 
203   // Params: 
204   //   path      Path to node 
205   //   id        Volume number 
206   //   transX    Translation in X
207   //   transZ    Translation in Y
208   //   transZ    Translation in Z
209   //   rotX      Rotation about X-axis 
210   //   rotY      Rotation about Y-axis
211   //   rotZ      Rotation about Z-axis 
212   AliDebug(1, Form("Make alignment for %s (volume %d): (%f,%f,%f) (%f,%f,%f)", 
213                    path.Data(), id, transX, transY, transZ, rotX, rotY, rotZ));
214   Int_t nAlign = fArray->GetEntries();
215   id = 0;
216   AliAlignObjAngles* obj = 
217     new ((*fArray)[nAlign]) AliAlignObjAngles(path.Data(), id,0,0,0,0,0,0,kTRUE);
218   if (!obj) {
219     AliError(Form("Failed to create alignment object for %s", path.Data()));
220     return kFALSE;
221   }
222   if (!obj->SetLocalPars(transX, transY, transZ, rotX, rotY, rotZ)) {
223     AliError(Form("Failed to set local transforms on %s", path.Data()));
224     return kTRUE;
225   }
226   return kTRUE;
227 }
228
229 //__________________________________________________________________
230 Bool_t
231 AliFMDAlignFaker::MakeAlignHalf(const TString& path, Int_t id)
232 {
233   // Make alignment of a half ring/cone 
234   AliDebug(15, Form("Make alignment for half-ring/cone %s", path.Data()));
235   Double_t transX = gRandom->Uniform(fHalfTransMin.X(), fHalfTransMax.X());
236   Double_t transY = gRandom->Uniform(fHalfTransMin.Y(), fHalfTransMax.Y());
237   Double_t transZ = gRandom->Uniform(fHalfTransMin.Z(), fHalfTransMax.Z());
238   Double_t rotX   = gRandom->Uniform(fHalfRotMin.X(),   fHalfRotMax.X());
239   Double_t rotY   = gRandom->Uniform(fHalfRotMin.Y(),   fHalfRotMax.Y());
240   Double_t rotZ   = gRandom->Uniform(fHalfRotMin.Z(),   fHalfRotMax.Z());
241   return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
242 }
243
244   
245 //__________________________________________________________________
246 Bool_t
247 AliFMDAlignFaker::MakeAlignSensor(const TString& path, Int_t id)
248 {
249   // Make alignment of a sensor 
250   AliDebug(15, Form("Make alignment for sensor %s", path.Data()));
251   Double_t transX = gRandom->Uniform(fSensorTransMin.X(), fSensorTransMax.X());
252   Double_t transY = gRandom->Uniform(fSensorTransMin.Y(), fSensorTransMax.Y());
253   Double_t transZ = gRandom->Uniform(fSensorTransMin.Z(), fSensorTransMax.Z());
254   Double_t rotX   = gRandom->Uniform(fSensorRotMin.X(),   fSensorRotMax.X());
255   Double_t rotY   = gRandom->Uniform(fSensorRotMin.Y(),   fSensorRotMax.Y());
256   Double_t rotZ   = gRandom->Uniform(fSensorRotMin.Z(),   fSensorRotMax.Z());
257   return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
258 }
259
260 //__________________________________________________________________
261 void
262 AliFMDAlignFaker::WriteToCDB()
263 {
264   // Make the objects. 
265   AliCDBManager*     cdb  = AliCDBManager::Instance();
266   AliCDBMetaData*    meta = new AliCDBMetaData; 
267   meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data()); 
268   meta->SetAliRootVersion(gROOT->GetVersion()); 
269   meta->SetBeamPeriod(1); 
270   meta->SetComment(fComment.Data());
271
272   AliCDBId id("FMD/Align/Data", fRunMin, fRunMax);
273   cdb->Put(fArray, id, meta);
274   cdb->Destroy();
275 }
276
277 //__________________________________________________________________
278 void
279 AliFMDAlignFaker::WriteToFile()
280 {
281   // Write to a local file 
282   TFile* file = TFile::Open(GetTitle(), "RECREATE");
283   if (!file) {
284     AliFatal(Form("Failed to open file '%s' for output", GetTitle()));
285     return;
286   }
287   file->cd();
288   fArray->Write("FMDAlignment");
289   file->Write();
290   file->Close();
291 }
292
293   
294   
295 //____________________________________________________________________
296 //
297 // EOF
298 //