1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
21 //____________________________________________________________________
30 //____________________________________________________________________
32 // Forward Multiplicity Detector based on Silicon wafers.
34 // This task creates fake alignment. Which alignment, depends on
35 // the bit mask passed to the constructor, or added by `AddAlign'.
37 // The default is to write all alignment parameters to a local
38 // storage `local://cdb' which is a directory in the current
41 #include "AliFMDDebug.h" // ALIFMDDEBUG_H ALILOG_H
42 #include "AliFMDAlignFaker.h" // ALIFMDALIGNFAKER_H
43 #include <AliCDBManager.h> // ALICDBMANAGER_H
44 #include <AliCDBStorage.h> // ALICDBSTORAGE_H
45 #include <AliCDBEntry.h> // ALICDBMANAGER_H
46 // #include <AliAlignObj.h>
47 #include <AliAlignObjParams.h>
48 // #include <Riostream.h>
52 #include <TClonesArray.h>
55 #include <TGeoManager.h>
57 // #include <TGeoVolume.h>
61 //====================================================================
62 ClassImp(AliFMDAlignFaker)
64 ; // This is here to keep Emacs for indenting the next line
67 //____________________________________________________________________
68 AliFMDAlignFaker::AliFMDAlignFaker(Int_t mask, const char* geo,
72 fSensorTransMin(0,0,0),
73 fSensorTransMax(0,0,0),
81 fRunMax(AliCDBRunRange::Infinity()),
85 // Default constructor
87 AliCDBStorage* storage = AliCDBManager::Instance()->GetDefaultStorage();
88 if (!storage) AliFatal("Default Storage not set");
89 const TString& uri = storage->GetURI();
92 SetSensorDisplacement();
94 SetHalfDisplacement();
99 //__________________________________________________________________
101 AliFMDAlignFaker::SetSensorDisplacement(Double_t x1, Double_t y1, Double_t z1,
102 Double_t x2, Double_t y2, Double_t z2)
104 // Set sensor displacement (unit is centimeters)
105 fSensorTransMin.SetXYZ(x1, y1, z1);
106 fSensorTransMax.SetXYZ(x2, y2, z2);
109 //__________________________________________________________________
111 AliFMDAlignFaker::SetSensorRotation(Double_t x1, Double_t y1, Double_t z1,
112 Double_t x2, Double_t y2, Double_t z2)
114 // Set sensor rotations (unit is degrees)
115 fSensorRotMin.SetXYZ(x1, y1, z1);
116 fSensorRotMax.SetXYZ(x2, y2, z2);
119 //__________________________________________________________________
121 AliFMDAlignFaker::SetHalfDisplacement(Double_t x1, Double_t y1, Double_t z1,
122 Double_t x2, Double_t y2, Double_t z2)
124 // Set half ring/cone displacement (unit is centimeters)
125 fHalfTransMin.SetXYZ(x1, y1, z1);
126 fHalfTransMax.SetXYZ(x2, y2, z2);
129 //__________________________________________________________________
131 AliFMDAlignFaker::SetHalfRotation(Double_t x1, Double_t y1, Double_t z1,
132 Double_t x2, Double_t y2, Double_t z2)
134 // Set half ring/cone rotations (unit is degrees)
135 fHalfRotMin.SetXYZ(x1, y1, z1);
136 fHalfRotMax.SetXYZ(x2, y2, z2);
139 //__________________________________________________________________
140 #define IS_NODE_HALF(name) \
141 (name[0] == 'F' && name[2] == 'M' && (name[3] == 'T' || name[3] == 'B'))
142 #define IS_NODE_SENSOR(name) \
143 (name[0] == 'F' && name[2] == 'S' && name[3] == 'E')
145 //__________________________________________________________________
147 AliFMDAlignFaker::GetGeometry(Bool_t toCdb, const TString& storage)
150 //load geom from default CDB storage
151 AliGeomManager::LoadGeometry();
154 if(!storage.BeginsWith("local://") &&
155 !storage.BeginsWith("alien://")) {
156 AliErrorClass(Form("STORAGE=\"%s\" invalid. Exiting\n", storage.Data()));
160 AliCDBManager* cdb = AliCDBManager::Instance();
161 AliCDBStorage* store = cdb->GetStorage(storage.Data());
163 AliErrorClass(Form("Unable to open storage %s\n", storage.Data()));
167 AliCDBPath path("GRP","Geometry","Data");
168 AliCDBEntry* entry = store->Get(path.GetPath(),cdb->GetRun());
170 AliErrorClass("Could not get the specified CDB entry!");
176 TGeoManager* geom = static_cast<TGeoManager*>(entry->GetObject());
177 AliGeomManager::SetGeometry(geom);
181 //__________________________________________________________________
183 AliFMDAlignFaker::Exec(Option_t*)
189 if (!TGeoManager::Import(GetName())) {
190 AliFatal(Form("Failed to import geometry from %s", GetName()));
195 TGeoVolume* topVolume = gGeoManager->GetTopVolume();
197 AliFatal("No top-level volume defined");
200 // Make container of transforms
201 if (!fArray) fArray = new TClonesArray("AliAlignObjParams");
205 TGeoIterator next(topVolume);
206 next.SetTopName(Form("/%s_1", topVolume->GetName()));
209 Char_t currentDet = '\0';
210 Char_t currentHalf = '\0';
211 // Loop over all entries in geometry to find our nodes.
212 while ((node = static_cast<TGeoNode*>(next()))) {
213 const char* name = node->GetName();
214 if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
215 !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors)))
218 TString path, alignName;
220 Int_t id = node->GetVolume()->GetNumber();
221 if (IS_NODE_HALF(name)) {
222 currentDet = name[1];
223 currentHalf = name[3];
224 alignName = Form("FMD/FMD%c_%c", currentDet, currentHalf);
226 if (IS_NODE_SENSOR(name)) {
227 Char_t ring = name[1];
228 Int_t lvl = next.GetLevel();
229 TGeoNode* parent = next.GetNode(lvl-1);
230 Int_t copy = parent->GetNumber();
231 alignName = Form("FMD/FMD%c_%c/FMD%c_%02d",
232 currentDet, currentHalf, ring, copy);
234 if (alignName.IsNull()) continue;
235 if (!gGeoManager->GetAlignableEntry(alignName.Data())) {
236 AliWarning(Form("No alignable entry for %s, using path %s",
237 alignName.Data(), path.Data()));
240 AliFMDDebug(1, ("Making alignment for %s -> %s (%d)",
241 alignName.Data(), path.Data(), id));
242 if (IS_NODE_HALF(name)) MakeAlignHalf(alignName, id);
243 if (IS_NODE_SENSOR(name)) MakeAlignSensor(alignName, id);
245 if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
246 !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors)))
250 TString path(Form("/%s", gGeoManager->GetNode(0)->GetName()));
251 Int_t nLevel = next.GetLevel();
252 for (Int_t lvl = 0; lvl <= nLevel; lvl++) {
253 TGeoNode* p = next.GetNode(lvl);
256 AliWarning(Form("No node at level %d in path %s",lvl,path.Data()));
259 if (!path.IsNull()) path.Append("/");
260 path.Append(p->GetName());
262 Int_t id = node->GetVolume()->GetNumber();
263 if (IS_NODE_HALF(name)) MakeAlignHalf(path, id);
264 if (IS_NODE_SENSOR(name)) MakeAlignSensor(path, id);
268 TString t(GetTitle());
269 if (t.IsNull() || t.Contains("local://") || t.Contains("alien://"))
275 //__________________________________________________________________
277 AliFMDAlignFaker::MakeAlign(const TString& path, Int_t id,
278 Double_t transX, Double_t transY, Double_t transZ,
279 Double_t rotX, Double_t rotY, Double_t rotZ)
281 // make alignment for a path
285 // transX Translation in X
286 // transZ Translation in Y
287 // transZ Translation in Z
288 // rotX Rotation about X-axis
289 // rotY Rotation about Y-axis
290 // rotZ Rotation about Z-axis
291 AliFMDDebug(3, ("Make alignment for %s (volume %d): (%f,%f,%f) (%f,%f,%f)",
292 path.Data(), id, transX, transY, transZ, rotX, rotY, rotZ));
293 Int_t nAlign = fArray->GetEntries();
295 AliAlignObjParams* obj =
296 new ((*fArray)[nAlign]) AliAlignObjParams(path.Data(),
297 id,0,0,0,0,0,0,kTRUE);
299 AliError(Form("Failed to create alignment object for %s", path.Data()));
302 if (!obj->SetLocalPars(transX, transY, transZ, rotX, rotY, rotZ)) {
303 AliError(Form("Failed to set local transforms on %s", path.Data()));
309 //__________________________________________________________________
311 AliFMDAlignFaker::MakeAlignHalf(const TString& path, Int_t id)
313 // Make alignment of a half ring/cone
314 AliFMDDebug(15, ("Make alignment for half-ring/cone %s", path.Data()));
315 Double_t transX = gRandom->Uniform(fHalfTransMin.X(), fHalfTransMax.X());
316 Double_t transY = gRandom->Uniform(fHalfTransMin.Y(), fHalfTransMax.Y());
317 Double_t transZ = gRandom->Uniform(fHalfTransMin.Z(), fHalfTransMax.Z());
318 Double_t rotX = gRandom->Uniform(fHalfRotMin.X(), fHalfRotMax.X());
319 Double_t rotY = gRandom->Uniform(fHalfRotMin.Y(), fHalfRotMax.Y());
320 Double_t rotZ = gRandom->Uniform(fHalfRotMin.Z(), fHalfRotMax.Z());
321 return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
325 //__________________________________________________________________
327 AliFMDAlignFaker::MakeAlignSensor(const TString& path, Int_t id)
329 // Make alignment of a sensor
330 AliFMDDebug(15, ("Make alignment for sensor %s", path.Data()));
331 Double_t transX = gRandom->Uniform(fSensorTransMin.X(), fSensorTransMax.X());
332 Double_t transY = gRandom->Uniform(fSensorTransMin.Y(), fSensorTransMax.Y());
333 Double_t transZ = gRandom->Uniform(fSensorTransMin.Z(), fSensorTransMax.Z());
334 Double_t rotX = gRandom->Uniform(fSensorRotMin.X(), fSensorRotMax.X());
335 Double_t rotY = gRandom->Uniform(fSensorRotMin.Y(), fSensorRotMax.Y());
336 Double_t rotZ = gRandom->Uniform(fSensorRotMin.Z(), fSensorRotMax.Z());
337 return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
340 //__________________________________________________________________
342 AliFMDAlignFaker::WriteToCDB()
345 const char* loc = GetTitle();
346 AliCDBManager* cdb = AliCDBManager::Instance();
347 AliCDBStorage* storage = cdb->GetStorage(!loc ? "" : loc);
348 AliCDBMetaData* meta = new AliCDBMetaData;
349 meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data());
350 meta->SetAliRootVersion(gROOT->GetVersion());
351 meta->SetBeamPeriod(1);
352 meta->SetComment(fComment.Data());
354 AliCDBId id("FMD/Align/Data", fRunMin, fRunMax);
355 storage->Put(fArray, id, meta);
358 //__________________________________________________________________
360 AliFMDAlignFaker::WriteToFile()
362 // Write to a local file
363 TFile* file = TFile::Open(GetTitle(), "RECREATE");
365 AliFatal(Form("Failed to open file '%s' for output", GetTitle()));
369 fArray->Write("FMDAlignment",TObject::kSingleKey);
376 //____________________________________________________________________