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 "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>
51 #include <TClonesArray.h>
54 #include <TGeoManager.h>
56 // #include <TGeoVolume.h>
59 //====================================================================
60 ClassImp(AliFMDAlignFaker)
62 ; // This is here to keep Emacs for indenting the next line
65 //____________________________________________________________________
66 AliFMDAlignFaker::AliFMDAlignFaker(Int_t mask, const char* geo,
74 // Default constructor
75 SetSensorDisplacement();
77 SetHalfDisplacement();
81 //__________________________________________________________________
83 AliFMDAlignFaker::SetSensorDisplacement(Double_t x1, Double_t y1, Double_t z1,
84 Double_t x2, Double_t y2, Double_t z2)
86 // Set sensor displacement (unit is centimeters)
87 fSensorTransMin.SetXYZ(x1, y1, z1);
88 fSensorTransMax.SetXYZ(x2, y2, z2);
91 //__________________________________________________________________
93 AliFMDAlignFaker::SetSensorRotation(Double_t x1, Double_t y1, Double_t z1,
94 Double_t x2, Double_t y2, Double_t z2)
96 // Set sensor rotations (unit is degrees)
97 fSensorRotMin.SetXYZ(x1, y1, z1);
98 fSensorRotMax.SetXYZ(x2, y2, z2);
101 //__________________________________________________________________
103 AliFMDAlignFaker::SetHalfDisplacement(Double_t x1, Double_t y1, Double_t z1,
104 Double_t x2, Double_t y2, Double_t z2)
106 // Set half ring/cone displacement (unit is centimeters)
107 fHalfTransMin.SetXYZ(x1, y1, z1);
108 fHalfTransMax.SetXYZ(x2, y2, z2);
111 //__________________________________________________________________
113 AliFMDAlignFaker::SetHalfRotation(Double_t x1, Double_t y1, Double_t z1,
114 Double_t x2, Double_t y2, Double_t z2)
116 // Set half ring/cone rotations (unit is degrees)
117 fHalfRotMin.SetXYZ(x1, y1, z1);
118 fHalfRotMax.SetXYZ(x2, y2, z2);
121 //__________________________________________________________________
122 #define IS_NODE_HALF(name) \
123 (name[0] == 'F' && name[2] == 'M' && (name[3] == 'T' || name[3] == 'B'))
124 #define IS_NODE_SENSOR(name) \
125 (name[0] == 'F' && name[2] == 'S' && name[3] == 'E')
127 //__________________________________________________________________
129 AliFMDAlignFaker::Exec(Option_t*)
135 if (!TGeoManager::Import(GetName())) {
136 AliFatal(Form("Failed to import geometry from %s", GetName()));
141 TGeoVolume* topVolume = gGeoManager->GetTopVolume();
143 AliFatal("No top-level volume defined");
146 // Make container of transforms
147 if (!fArray) fArray = new TClonesArray("AliAlignObjAngles");
151 TGeoIterator next(topVolume);
154 // Loop over all entries in geometry to find our nodes.
155 while ((node = static_cast<TGeoNode*>(next()))) {
156 const char* name = node->GetName();
157 if (!(IS_NODE_HALF(name) && TESTBIT(fMask, kHalves)) &&
158 !(IS_NODE_SENSOR(name) && TESTBIT(fMask, kSensors)))
162 TString path(Form("/%s", gGeoManager->GetNode(0)->GetName()));
163 Int_t nLevel = next.GetLevel();
164 for (Int_t lvl = 0; lvl <= nLevel; lvl++) {
165 TGeoNode* p = next.GetNode(lvl);
168 AliWarning(Form("No node at level %d in path %s",lvl,path.Data()));
171 if (!path.IsNull()) path.Append("/");
172 path.Append(p->GetName());
174 Int_t id = node->GetVolume()->GetNumber();
175 if (IS_NODE_HALF(name)) MakeAlignHalf(path, id);
176 if (IS_NODE_SENSOR(name)) MakeAlignSensor(path, id);
179 TString t(GetTitle());
180 if (t.IsNull() || t.Contains("local://") || t.Contains("alien://"))
186 //__________________________________________________________________
188 AliFMDAlignFaker::MakeAlign(const TString& path, Int_t id,
189 Double_t transX, Double_t transY, Double_t transZ,
190 Double_t rotX, Double_t rotY, Double_t rotZ)
192 // make alignment for a path
196 // transX Translation in X
197 // transZ Translation in Y
198 // transZ Translation in Z
199 // rotX Rotation about X-axis
200 // rotY Rotation about Y-axis
201 // rotZ Rotation about Z-axis
202 AliDebug(1, Form("Make alignment for %s (volume %d): (%f,%f,%f) (%f,%f,%f)",
203 path.Data(), id, transX, transY, transZ, rotX, rotY, rotZ));
204 Int_t nAlign = fArray->GetEntries();
206 AliAlignObjAngles* obj =
207 new ((*fArray)[nAlign]) AliAlignObjAngles(path.Data(), id,0,0,0,0,0,0);
209 AliError(Form("Failed to create alignment object for %s", path.Data()));
212 if (!obj->SetLocalPars(transX, transY, transZ, rotX, rotY, rotZ)) {
213 AliError(Form("Failed to set local transforms on %s", path.Data()));
219 //__________________________________________________________________
221 AliFMDAlignFaker::MakeAlignHalf(const TString& path, Int_t id)
223 // Make alignment of a half ring/cone
224 AliDebug(15, Form("Make alignment for half-ring/cone %s", path.Data()));
225 Double_t transX = gRandom->Uniform(fHalfTransMin.X(), fHalfTransMax.X());
226 Double_t transY = gRandom->Uniform(fHalfTransMin.Y(), fHalfTransMax.Y());
227 Double_t transZ = gRandom->Uniform(fHalfTransMin.Z(), fHalfTransMax.Z());
228 Double_t rotX = gRandom->Uniform(fHalfRotMin.X(), fHalfRotMax.X());
229 Double_t rotY = gRandom->Uniform(fHalfRotMin.Y(), fHalfRotMax.Y());
230 Double_t rotZ = gRandom->Uniform(fHalfRotMin.Z(), fHalfRotMax.Z());
231 return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
235 //__________________________________________________________________
237 AliFMDAlignFaker::MakeAlignSensor(const TString& path, Int_t id)
239 // Make alignment of a sensor
240 AliDebug(15, Form("Make alignment for sensor %s", path.Data()));
241 Double_t transX = gRandom->Uniform(fSensorTransMin.X(), fSensorTransMax.X());
242 Double_t transY = gRandom->Uniform(fSensorTransMin.Y(), fSensorTransMax.Y());
243 Double_t transZ = gRandom->Uniform(fSensorTransMin.Z(), fSensorTransMax.Z());
244 Double_t rotX = gRandom->Uniform(fSensorRotMin.X(), fSensorRotMax.X());
245 Double_t rotY = gRandom->Uniform(fSensorRotMin.Y(), fSensorRotMax.Y());
246 Double_t rotZ = gRandom->Uniform(fSensorRotMin.Z(), fSensorRotMax.Z());
247 return MakeAlign(path, id, transX, transY, transZ, rotX, rotY, rotZ);
250 //__________________________________________________________________
252 AliFMDAlignFaker::WriteToCDB()
255 AliCDBManager* cdb = AliCDBManager::Instance();
256 if (GetTitle() && GetTitle()[0] != '\0')
257 cdb->SetDefaultStorage(GetTitle());
259 AliCDBMetaData* meta = new AliCDBMetaData;
260 meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data());
261 meta->SetAliRootVersion(gROOT->GetVersion());
262 meta->SetBeamPeriod(1);
263 meta->SetComment("Dummy data for testing");
265 AliCDBId id("FMD/Align/Data", fRunMin, fRunMax);
266 cdb->Put(fArray, id, meta);
270 //__________________________________________________________________
272 AliFMDAlignFaker::WriteToFile()
274 // Write to a local file
275 TFile* file = TFile::Open(GetTitle(), "RECREATE");
277 AliFatal(Form("Failed to open file '%s' for output", GetTitle()));
281 fArray->Write("FMDAlignment");
288 //____________________________________________________________________