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