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