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