]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSAlignMille2.cxx
Millepede implementation operating with AliMillePede2 (R. Shahoyan) added to the...
[u/mrichter/AliRoot.git] / ITS / AliITSAlignMille2.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 /// \class AliITSAlignMille
20 /// Alignment class fro the ALICE ITS detector
21 ///
22 /// ITS specific alignment class which interface to AliMillepede.   
23 /// For each track ProcessTrack calculates the local and global derivatives
24 /// at each hit and fill the corresponding local equations. Provide methods for
25 /// fixing or constraining detection elements for best results. 
26 ///
27 /// \author M. Lunardon (thanks to J. Castillo)
28 //-----------------------------------------------------------------------------
29
30 #include <TF1.h>
31 #include <TFile.h>
32 #include <TClonesArray.h>
33 #include <TGraph.h>
34 #include <TGeoMatrix.h>
35 #include <TMath.h>
36 #include <TGraphErrors.h>
37 #include <TVirtualFitter.h>
38
39 #include "AliITSAlignMille2.h"
40 #include "AliITSgeomTGeo.h"
41 #include "AliGeomManager.h"
42 #include "AliMillePede2.h"
43 #include "AliTrackPointArray.h"
44 #include "AliAlignObjParams.h"
45 #include "AliLog.h"
46 #include "TSystem.h"  // come si fa?
47 #include "AliTrackFitterRieman.h"
48
49 /// \cond CLASSIMP
50 ClassImp(AliITSAlignMille2)
51 /// \endcond
52
53 AliITSAlignMille2* AliITSAlignMille2::fgInstance = 0;  
54 Int_t              AliITSAlignMille2::fgInstanceID = 0;
55
56 AliITSAlignMille2::AliITSAlignMille2(const Char_t *configFilename, Bool_t initmille) 
57 : TObject(),
58   fMillepede(0),
59   fStartFac(16.), 
60   fResCutInitial(100.), 
61   fResCut(100.),
62   fNGlobal(0),
63   fNLocal(4),
64   fNStdDev(3),
65   fIsMilleInit(kFALSE),
66   fSensorsIn(kFALSE),
67   fParSigTranslations(0.0100),
68   fParSigRotations(0.1),
69 //
70   fCurrentModule(0),
71   fTrack(0),
72   fCluster(),
73   fGlobalDerivatives(0), 
74 //
75   fMinNPtsPerTrack(3),
76   fInitTrackParamsMeth(1),
77   fTotBadLocEqPoints(0),
78   fRieman(0),
79   //
80   fUseGlobalDelta(kFALSE),
81   fRequirePoints(kFALSE),
82   fTempExcludedModule(-1),
83   //
84   fGeometryFileName("geometry.root"),
85   fPreAlignmentFileName(""),
86   fGeoManager(0),
87   fIsConfigured(kFALSE),
88   fPreAlignQF(0),
89 //
90   fPrealignment(0),
91   fMilleModule(2),
92   fSuperModule(2),
93   fNModules(0),
94   fNSuperModules(0),
95   fUsePreAlignment(kFALSE),
96   fUseSortedTracks(kTRUE),
97   fBOn(kFALSE),
98   fBField(0.0),
99   fBug(0)
100 {
101   /// main constructor that takes input from configuration file
102   //  
103   fMillepede = new AliMillePede2();
104   for (int i=3;i--;) fSigmaFactor[i] = 1.0;
105   //
106   // new RS
107   for (Int_t i=0; i<6; i++) {
108     fNReqLayUp[i]=0;
109     fNReqLayDown[i]=0;
110     fNReqLay[i]=0;
111   }
112   for (Int_t i=0; i<3; i++) {
113     fNReqDetUp[i]=0;
114     fNReqDetDown[i]=0;
115     fNReqDet[i]=0;
116   }
117   //
118   Int_t lc=LoadConfig(configFilename);
119   if (lc) {
120     AliInfo(Form("Error %d loading configuration from %s",lc,configFilename));
121   }
122   else {    
123     fIsConfigured=kTRUE;
124     if (initmille) {
125       AliInfo(Form("Initializing Millepede with %d gpar, %d lpar and %d stddev ...",fNGlobal, fNLocal, fNStdDev));
126       Init(fNGlobal, fNLocal, fNStdDev);      
127       ResetLocalEquation();    
128       AliInfo("Parameters initialized to zero");
129     }
130     else {
131       AliInfo("Millepede has not been initialized ... ");
132     }
133   }
134   //
135   fgInstance = this;
136   fgInstanceID++;
137   //
138 }
139
140 AliITSAlignMille2::~AliITSAlignMille2() {
141   /// Destructor
142   if (fMillepede)         delete fMillepede;
143   if (fGlobalDerivatives) delete[] fGlobalDerivatives;
144   if (fRieman)            delete fRieman;
145   if (fPrealignment)      delete fPrealignment;
146   fMilleModule.Delete();
147   fSuperModule.Delete();
148   if (--fgInstanceID==0) fgInstance = 0;
149 }
150
151 ///////////////////////////////////////////////////////////////////////
152
153 Int_t AliITSAlignMille2::LoadConfig(const Char_t *cfile) {
154   /// return 0 if success
155   ///        1 if error in module index or voluid
156   
157   FILE *pfc=fopen(cfile,"r");
158   if (!pfc) return -1;
159   
160   Char_t st[200],st2[200];
161   Char_t tmp[100];
162   Int_t idx,itx,ity,itz,ith,ips,iph;
163   Float_t f1,f2,f3;
164   UShort_t voluid;
165   Int_t nmod=0;
166   //
167   while (fgets(st,200,pfc)) {
168     //
169     for (int i=0; i<int(strlen(st)); i++) if (st[i]=='#') st[i]=0; // skip comments
170     //
171     if (strstr(st,"GEOMETRY_FILE")) {
172       sscanf(st,"%s %s",tmp,st2);
173       if (gSystem->AccessPathName(st2)) { AliInfo("*** WARNING! *** geometry file not found! "); return -1;}  
174       fGeometryFileName=st2;
175       InitGeometry();
176     }
177     //
178     if (strstr(st,"PREALIGNMENT_FILE")) {
179       sscanf(st,"%s %s",tmp,st2);
180       if (gSystem->AccessPathName(st2)) { AliInfo("*** WARNING! *** prealignment file not found! "); return -1;}  
181       fPreAlignmentFileName=st2;
182       itx=ApplyToGeometry();
183       if (itx) { AliInfo(Form("*** WARNING! *** error %d reading prealignment file! ",itx)); return -6;}
184     }
185     //
186     if (strstr(st,"SUPERMODULE_FILE")) {
187       sscanf(st,"%s %s",tmp,st2);
188       if (gSystem->AccessPathName(st2)) { AliInfo("*** WARNING! *** supermodule file not found! "); return -1;}  
189       if (LoadSuperModuleFile(st2)) return -1;
190     }
191     //
192     if (strstr(st,"SET_B_FIELD")) {
193       sscanf(st,"%s %f",tmp,&f1);
194       if (f1>0) {
195         fBField = f1;
196         fBOn = kTRUE;
197         fNLocal = 5; // helices
198         fRieman = new AliTrackFitterRieman();
199       }  
200       else {
201         fBField = 0.0;
202         fBOn = kFALSE;
203         fNLocal = 4;
204       }
205     }
206     //
207     if (strstr(st,"SET_MINPNT_TRA")) {
208       sscanf(st,"%s %d",tmp,&idx);
209       fMinNPtsPerTrack=idx;
210     }
211     //
212     if (strstr(st,"SET_PARSIG_TRA")) {
213       sscanf(st,"%s %f",tmp,&f1);
214       fParSigTranslations=f1;
215     }
216     //
217     if (strstr(st,"SET_PARSIG_ROT")) {
218       sscanf(st,"%s %f",tmp,&f1);
219       fParSigRotations=f1;
220     }
221     //
222     if (strstr(st,"SET_NSTDDEV")) {
223       sscanf(st,"%s %d",tmp,&idx);
224       fNStdDev=idx;
225     }
226     //
227     if (strstr(st,"SET_RESCUT_INIT")) {
228       sscanf(st,"%s %f",tmp,&f1);
229       fResCutInitial=f1;
230     }
231     //
232     if (strstr(st,"SET_RESCUT_OTHER")) {
233       sscanf(st,"%s %f",tmp,&f1);
234       fResCut=f1;
235     }
236     //
237     if (strstr(st,"SET_LOCALSIGMAFACTOR")) {
238       f1=f2=f3=0;
239       sscanf(st,"%s %f %f %f",tmp,&f1,&f2,&f3);
240       if (f1>0) fSigmaFactor[0] = f1;
241       if (f2>0) fSigmaFactor[1] = f2; else fSigmaFactor[1]=fSigmaFactor[0];
242       if (f3>0) fSigmaFactor[2] = f3; else fSigmaFactor[2]=fSigmaFactor[1];
243     }
244     //
245     if (strstr(st,"SET_STARTFAC")) {
246       sscanf(st,"%s %f",tmp,&f1);
247       fStartFac=f1;
248     }
249     //
250     // >> RS
251     if (strstr(st,"REQUIRE_POINT")) {
252       // syntax:   REQUIRE_POINT where ndet updw nreqpts
253       //    where = LAYER or DETECTOR
254       //    ndet = detector number: 1-6 for LAYER and 1-3 for DETECTOR (SPD=1, SDD=2, SSD=3)
255       //    updw = 1 for Y>0, -1 for Y<0, 0 if not specified
256       //    nreqpts = minimum number of points of that type
257       sscanf(st,"%s %s %d %d %d",tmp,st2,&itx,&ity,&itz);
258       itx--;
259       if (strstr(st2,"LAYER")) {
260         if (itx<0 || itx>5) return -7;
261         if (ity>0) fNReqLayUp[itx]=itz;
262         else if (ity<0) fNReqLayDown[itx]=itz;
263         else fNReqLay[itx]=itz;
264         fRequirePoints=kTRUE;
265       }
266       else if (strstr(st2,"DETECTOR")) { // DETECTOR
267         if (itx<0 || itx>2) return -7;
268         if (ity>0) fNReqDetUp[itx]=itz;
269         else if (ity<0) fNReqDetDown[itx]=itz;
270         else fNReqDet[itx]=itz; 
271         fRequirePoints=kTRUE;
272       }
273     }
274     // << RS
275     
276     //
277     if (strstr(st,"MODULE_INDEX") || strstr(st,"MODULE_VOLUID")) { 
278       f1=f2=f3=0;
279       sscanf(st,"%s %d %d %d %d %d %d %d %f %f %f",tmp,&idx,&itx,&ity,&itz,&iph,&ith,&ips,&f1,&f2,&f3);
280       //
281       if (idx<=kMaxITSSensID) voluid=GetModuleVolumeID(idx);
282       else voluid = UShort_t(idx);
283       //
284       if (voluid>=kMinITSSupeModuleID) { // custom supermodule
285         int ism=-1;
286         for (int j=0; j<fNSuperModules; j++) if (voluid==GetSuperModule(j)->GetVolumeID()) ism=j;
287         if (ism<0) return -1; // bad volid
288         fMilleModule.AddAtAndExpand(new AliITSAlignMille2Module(*GetSuperModule(ism)),nmod);
289         // >> RS
290 //      if (f1>0) {
291 //        for (int kk=0; kk<GetMilleModule(nmod)->GetNSensitiveVolumes(); kk++) {
292 //          idx=AliITSAlignMille2Module::GetIndexFromVolumeID(GetMilleModule(nmod)->GetSensitiveVolumeVolumeID()[kk]);
293 //          if (idx>=0) fSensVolSigmaXfactor[idx]=f1;
294 //        }
295 //      }
296 //      if (f2>0) {
297 //        for (int kk=0; kk<GetMilleModule(nmod)->GetNSensitiveVolumes(); kk++) {
298 //          idx=AliITSAlignMille2Module::GetIndexFromVolumeID(GetMilleModule(nmod)->GetSensitiveVolumeVolumeID()[kk]);
299 //          if (idx>=0) fSensVolSigmaZfactor[idx]=f2;
300 //        }
301 //      }
302         // << RS
303       }
304       else if (idx<=kMaxITSSensVID) {
305         fMilleModule.AddAtAndExpand(new AliITSAlignMille2Module(voluid),nmod);
306         AliITSAlignMille2Module* md = (AliITSAlignMille2Module*) fMilleModule[nmod];
307         fSensorsIn = kTRUE;
308         md->SetSensorsProvided();
309       }
310       else return -1;  // bad volid
311       //
312       AliITSAlignMille2Module* mod = GetMilleModule(nmod);
313       mod->SetFreeDOF(kDOFTX,itx!=0);
314       mod->SetFreeDOF(kDOFTY,ity!=0);
315       mod->SetFreeDOF(kDOFTZ,itz!=0);
316       mod->SetFreeDOF(kDOFPH,iph!=0);
317       mod->SetFreeDOF(kDOFTH,ith!=0);
318       mod->SetFreeDOF(kDOFPS,ips!=0);
319       //
320       mod->SetUniqueID(nmod);
321       if (f1>0) mod->SetSigmaXFactor(f1);
322       if (f2>0) mod->SetSigmaYFactor(f2); else mod->SetSigmaYFactor(mod->GetSigmaXFactor());
323       if (f3>0) mod->SetSigmaZFactor(f3); else mod->SetSigmaZFactor(mod->GetSigmaYFactor());
324       nmod++;
325     }
326     //
327   } // end while
328   //
329   fNModules = nmod;
330   fNGlobal = fNModules*kNParCh;
331   //
332   fclose(pfc);
333   //
334   // set parent/child relationship for modules to align
335   printf("Setting parent/child relationships\n");
336   //
337   for (int ipar=0;ipar<nmod;ipar++) {
338     AliITSAlignMille2Module* parent = GetMilleModule(ipar);
339     if (parent->GetIndex()<=kMaxITSSensID) continue; // sensor cannot be a parent
340     //
341     for (int icld=0;icld<nmod;icld++) {
342       if (icld==ipar) continue;
343       AliITSAlignMille2Module* child = GetMilleModule(icld);
344       if (!child->BelongsTo(parent)) continue;
345       //
346       AliITSAlignMille2Module* parOld = child->GetParent();
347       if (parOld && parOld->GetNSensitiveVolumes()<parent->GetNSensitiveVolumes()) continue; // parOld is closer
348       child->SetParent(parent);
349     }
350     //
351   }
352   //
353   // reorder the modules in such a way that parents come first
354   for (int icld=0;icld<nmod;icld++) {
355     AliITSAlignMille2Module* child  = GetMilleModule(icld);
356     AliITSAlignMille2Module* parent; 
357     while ( (parent=child->GetParent()) &&  (parent->GetUniqueID()<child->GetUniqueID()) ) {
358       // swap
359       fMilleModule[icld] = parent;
360       fMilleModule[parent->GetUniqueID()] = child;
361       child->SetUniqueID(parent->GetUniqueID());
362       parent->SetUniqueID(icld);
363       child = parent;
364     }
365     //
366   } 
367   //
368   // go over the child->parent chain and mark modules with explicitly provided sensors
369   for (int icld=nmod;icld--;) {
370     AliITSAlignMille2Module* child = GetMilleModule(icld);
371     AliITSAlignMille2Module* parent = child->GetParent();
372     if (!parent) continue;
373     parent->SetSensorsProvided( child->AreSensorsProvided() );
374     if (!parent->AreSensorsProvided()) continue;
375     // suppress unused sensors
376     for (int isn=0;isn<parent->GetNSensitiveVolumes();isn++) {
377       int snVID = parent->GetSensVolVolumeID(isn);
378       // check if this sensor is explicitly requested
379       for (int imd=nmod;imd--;) if (GetMilleModule(imd)->GetVolumeID() == snVID) {snVID = -1; break;}
380       //
381       if (snVID==-1) continue; // found this sensor, do nothing
382       //
383       // otherwise, remove this sensor from the module list
384       AliInfo(Form("Removing sensor %d from %s",snVID,parent->GetName()));
385       parent->DelSensitiveVolume(isn--);
386     }
387   }
388   //
389   fGlobalDerivatives = new Double_t[fNGlobal];
390   memset(fGlobalDerivatives,0,fNGlobal*sizeof(Double_t));
391   //
392   return 0;
393 }
394
395
396 void AliITSAlignMille2::SetRequiredPoint(Char_t* where, Int_t ndet, Int_t updw, Int_t nreqpts) 
397 {
398   // set minimum number of points in specific detector or layer
399   // where = LAYER or DETECTOR
400   // ndet = detector number: 1-6 for LAYER and 1-3 for DETECTOR (SPD=1, SDD=2, SSD=3)
401   // updw = 1 for Y>0, -1 for Y<0, 0 if not specified
402   // nreqpts = minimum number of points of that type
403   ndet--;
404   if (strstr(where,"LAYER")) {
405     if (ndet<0 || ndet>5) return;
406     if (updw>0) fNReqLayUp[ndet]=nreqpts;
407     else if (updw<0) fNReqLayDown[ndet]=nreqpts;
408     else fNReqLay[ndet]=nreqpts;
409     fRequirePoints=kTRUE;
410   }
411   else if (strstr(where,"DETECTOR")) {
412     if (ndet<0 || ndet>2) return;
413     if (updw>0) fNReqDetUp[ndet]=nreqpts;
414     else if (updw<0) fNReqDetDown[ndet]=nreqpts;
415     else fNReqDet[ndet]=nreqpts;        
416     fRequirePoints=kTRUE;
417   }
418 }
419
420
421 Int_t AliITSAlignMille2::GetModuleIndex(const Char_t *symname) {
422   /// index from symname
423   if (!symname) return -1;
424   for (Int_t i=0;i<=kMaxITSSensID; i++) {
425     if (!strcmp(symname,AliITSgeomTGeo::GetSymName(i))) return i;
426   }
427   return -1;
428 }
429
430 Int_t AliITSAlignMille2::GetModuleIndex(UShort_t voluid) {
431   /// index from volume ID
432   AliGeomManager::ELayerID lay = AliGeomManager::VolUIDToLayer(voluid);
433   if (lay<1|| lay>6) return -1;
434   Int_t idx=Int_t(voluid)-2048*lay;
435   if (idx>=AliGeomManager::LayerSize(lay)) return -1;
436   for (Int_t ilay=1; ilay<lay; ilay++) 
437     idx += AliGeomManager::LayerSize(ilay);
438   return idx;
439 }
440
441 UShort_t AliITSAlignMille2::GetModuleVolumeID(const Char_t *symname) {
442   /// volume ID from symname
443   /// works for sensitive volumes only
444   if (!symname) return 0;
445
446   for (UShort_t voluid=2000; voluid<13300; voluid++) {
447     Int_t modId;
448     AliGeomManager::ELayerID layerId = AliGeomManager::VolUIDToLayer(voluid,modId);
449     if (layerId>0 && layerId<7 && modId>=0 && modId<AliGeomManager::LayerSize(layerId)) {
450       if (!strcmp(symname,AliGeomManager::SymName(layerId,modId))) return voluid;
451     }
452   }
453
454   return 0;
455 }
456
457 UShort_t AliITSAlignMille2::GetModuleVolumeID(Int_t index) {
458   /// volume ID from index
459   if (index<0) return 0;
460   if (index<2198)
461     return GetModuleVolumeID(AliITSgeomTGeo::GetSymName(index));
462   else {
463     for (int i=0; i<fNSuperModules; i++) {
464       if (GetSuperModule(i)->GetIndex()==index) return GetSuperModule(i)->GetVolumeID();
465     }
466   }
467   return 0;
468 }
469
470 void AliITSAlignMille2::InitGeometry() {
471   /// initialize geometry
472   AliGeomManager::LoadGeometry(fGeometryFileName.Data());
473   fGeoManager = AliGeomManager::GetGeometry();
474   if (!fGeoManager) {
475     AliInfo("Couldn't initialize geometry");
476     return;
477   }
478 }
479
480 void AliITSAlignMille2::Init(Int_t nGlobal,  /* number of global paramers */
481                            Int_t nLocal,   /* number of local parameters */
482                            Int_t nStdDev   /* std dev cut */ )
483 {
484   /// Initialization of AliMillepede. Fix parameters, define constraints ...
485   fMillepede->InitMille(nGlobal,nLocal,nStdDev,fResCut,fResCutInitial);
486   fIsMilleInit = kTRUE;
487   
488   /// Fix non free parameters
489   for (Int_t i=0; i<fNModules; i++) {
490     for (Int_t j=0; j<kNParCh; j++) {
491       if (!GetMilleModule(i)->IsFreeDOF(j)) FixParameter(i*kNParCh+j,0.0);
492       else {
493         // pepopepo: da verificare il settaggio delle sigma, ma forse va bene...
494         Double_t parsig=0;
495         if (j<3) parsig = fParSigTranslations; // translations (0.0100 cm)
496         else     parsig = fParSigRotations; // rotations (1/10 deg)
497         FixParameter(i*kNParCh+j,parsig);
498       }
499     }    
500   }
501   //
502   // Set iterations
503   if (fStartFac>1) fMillepede->SetIterations(fStartFac);          
504 }
505
506
507 void AliITSAlignMille2::AddConstraint(Double_t *par, Double_t value) {
508   /// Constrain equation defined by par to value
509   if (!fIsMilleInit) {
510     AliInfo("Millepede has not been initialized!");
511     return;
512   }
513   fMillepede->SetGlobalConstraint(par, value);
514   AliInfo("Adding constraint");
515 }
516
517 void AliITSAlignMille2::InitGlobalParameters(Double_t *par) {
518   /// Initialize global parameters with par array
519   if (!fIsMilleInit) {
520     AliInfo("Millepede has not been initialized!");
521     return;
522   }
523   fMillepede->SetGlobalParameters(par);
524   AliInfo("Init Global Parameters");
525 }
526  
527 void AliITSAlignMille2::FixParameter(Int_t iPar, Double_t value) {
528   /// Parameter iPar is encourage to vary in [-value;value]. 
529   /// If value == 0, parameter is fixed
530   if (!fIsMilleInit) {
531     AliInfo("Millepede has not been initialized!");
532     return;
533   }
534   fMillepede->SetParSigma(iPar, value);
535   if (value==0) AliInfo(Form("Parameter %i Fixed", iPar));
536 }
537
538 void AliITSAlignMille2::ResetLocalEquation()
539 {
540   /// Reset the derivative vectors
541   for(int i=fNLocal;i--;)  fLocalDerivatives[i] = 0.0;
542   for(int i=fNGlobal;i--;) fGlobalDerivatives[i] = 0.0;
543 }
544
545 Int_t AliITSAlignMille2::ApplyToGeometry() 
546 {
547   // apply starting realignment to ideal geometry
548   //
549   if (!fGeoManager) return -1; 
550   TFile *pref = new TFile(fPreAlignmentFileName.Data());
551   if (!pref->IsOpen()) return -2;
552   fPrealignment = (TClonesArray*)pref->Get("ITSAlignObjs");
553   if (!fPrealignment) return -3;  
554   Int_t nprea = fPrealignment->GetEntriesFast();
555   AliInfo(Form("Array of input misalignments with %d entries",nprea));
556   //
557   for (int ix=0; ix<nprea; ix++) {
558     AliAlignObjParams *preo=(AliAlignObjParams*) fPrealignment->At(ix);
559     Int_t index=AliITSAlignMille2Module::GetIndexFromVolumeID(preo->GetVolUID());
560     if (index>=0) {
561       if (index>=fPreAlignQF.GetSize()) fPreAlignQF.Set(index+10);
562       fPreAlignQF[index] = (int) preo->GetUniqueID()+1;
563     }
564     //TString nms = preo->GetSymName();
565     //if (!nms.Contains("Ladder")) continue; //RRR
566     if (!preo->ApplyToGeometry()) return -4;
567   }
568   //
569   pref->Close();
570   delete pref;
571   //
572   fUsePreAlignment = kTRUE;
573   return 0;
574 }
575
576 Int_t AliITSAlignMille2::GetPreAlignmentQualityFactor(Int_t index) const
577 {
578   if (!fUsePreAlignment || index<0 || index>=fPreAlignQF.GetSize()) return -1;
579   return fPreAlignQF[index]-1;
580 }
581
582 AliTrackPointArray *AliITSAlignMille2::PrepareTrack(const AliTrackPointArray *atp) {
583   /// create a new AliTrackPointArray keeping only defined modules
584   /// move points according to a given prealignment, if any
585   /// sort alitrackpoints w.r.t. global Y direction, if selected
586
587   AliTrackPointArray *atps=NULL;
588   Int_t idx[20];
589   Int_t npts=atp->GetNPoints();
590
591   /// checks if AliTrackPoints belong to defined modules
592   Int_t ngoodpts=0;
593   Int_t intidx[20];
594   
595   for (int j=0; j<npts; j++) {
596     intidx[j] = IsContained(atp->GetVolumeID()[j]);
597     if (intidx[j]>=0) ngoodpts++;
598   }
599   AliDebug(3,Form("Number of points in defined modules: %d",ngoodpts));
600
601   // reject track if not enough points are left
602   if (ngoodpts<fMinNPtsPerTrack) {
603     AliInfo("Track with not enough points!");
604     return NULL;
605   }
606   // >> RS
607   AliTrackPoint p;
608   // check points in specific places
609   if (fRequirePoints) {
610     Int_t nlayup[6],nlaydown[6],nlay[6];
611     Int_t ndetup[3],ndetdown[3],ndet[3];
612     for (Int_t j=0; j<6; j++) {nlayup[j]=0; nlaydown[j]=0; nlay[j]=0;}
613     for (Int_t j=0; j<3; j++) {ndetup[j]=0; ndetdown[j]=0; ndet[j]=0;}
614     
615     for (int i=0; i<npts; i++) {
616       // skip not defined points
617       if (intidx[i]<0) continue;
618       Float_t xx=atp->GetX()[i];
619       Float_t yy=atp->GetY()[i];
620       Float_t r=TMath::Sqrt(xx*xx + yy*yy);
621       int lay=-1;
622       if (r<5) lay=0;
623       else if (r>5 && r<10) lay=1;
624       else if (r>10 && r<18) lay=2;
625       else if (r>18 && r<30) lay=3;
626       else if (r>30 && r<40) lay=4;
627       else if (r>40) lay=5;
628       if (lay<0) continue;
629       int det=lay/2;
630       //printf("Point %d - x=%f  y=%f  R=%f  lay=%d  det=%d\n",i,xx,yy,r,lay,det);
631
632       if (yy>=0.0) { // UP point
633         nlayup[lay]++;
634         nlay[lay]++;
635         ndetup[det]++;
636         ndet[det]++;
637       }
638       else {
639         nlaydown[lay]++;
640         nlay[lay]++;
641         ndetdown[det]++;
642         ndet[det]++;
643       }
644     }
645     
646     // checks minimum values
647     Bool_t isok=kTRUE;
648     for (Int_t j=0; j<6; j++) {
649       if (nlayup[j]<fNReqLayUp[j]) isok=kFALSE; 
650       if (nlaydown[j]<fNReqLayDown[j]) isok=kFALSE; 
651       if (nlay[j]<fNReqLay[j]) isok=kFALSE; 
652     }
653     for (Int_t j=0; j<3; j++) {
654       if (ndetup[j]<fNReqDetUp[j]) isok=kFALSE; 
655       if (ndetdown[j]<fNReqDetDown[j]) isok=kFALSE; 
656       if (ndet[j]<fNReqDet[j]) isok=kFALSE; 
657     }
658     if (!isok) {
659       AliDebug(2,Form("Track does not meet all location point requirements!"));
660       return NULL;
661     }
662   }
663   
664   // << RS
665   // build a new track with (sorted) (prealigned) good points
666   atps=new AliTrackPointArray(ngoodpts);
667
668   for (int i=0; i<npts; i++) idx[i]=i;
669   // sort track if required
670   if (fUseSortedTracks) TMath::Sort(npts,atp->GetY(),idx); // sort descending...
671
672   Int_t npto=0;
673   for (int i=0; i<npts; i++) {
674     // skip not defined points
675     if (intidx[idx[i]]<0) continue;
676     atp->GetPoint(p,idx[i]);
677
678     // prealign point if required
679     // get IDEAL matrix
680     TGeoHMatrix *svOrigMatrix = GetMilleModule(intidx[idx[i]])->GetSensitiveVolumeOrigGlobalMatrix(p.GetVolumeID());
681     // get back real local coordinates: use OriginalGlobalMatrix because AliTrackPoints were written
682     // with idel geometry  
683     Double_t pg[3],pl[3];
684     pg[0]=p.GetX();
685     pg[1]=p.GetY();
686     pg[2]=p.GetZ();
687     //    printf("Global coordinates of measured point : X=%f  Y=%f  Z=%f \n",pg[0],pg[1],pg[2]);
688     AliDebug(3,Form("Global coordinates of measured point : X=%f  Y=%f  Z=%f \n",pg[0],pg[1],pg[2]));
689     svOrigMatrix->MasterToLocal(pg,pl);
690
691     AliDebug(3,Form("Local coordinates of measured point : X=%f  Y=%f  Z=%f \n",pl[0],pl[1],pl[2]));
692
693     // update covariance matrix
694     TGeoHMatrix hcov;
695     Double_t hcovel[9];
696     hcovel[0]=double(p.GetCov()[0]);
697     hcovel[1]=double(p.GetCov()[1]);
698     hcovel[2]=double(p.GetCov()[2]);
699     hcovel[3]=double(p.GetCov()[1]);
700     hcovel[4]=double(p.GetCov()[3]);
701     hcovel[5]=double(p.GetCov()[4]);
702     hcovel[6]=double(p.GetCov()[2]);
703     hcovel[7]=double(p.GetCov()[4]);
704     hcovel[8]=double(p.GetCov()[5]);
705     hcov.SetRotation(hcovel);
706     // now rotate in local system
707     hcov.Multiply(svOrigMatrix);
708     hcov.MultiplyLeft(&svOrigMatrix->Inverse());
709     // now hcov is LOCAL COVARIANCE MATRIX
710
711     // >> RS
712     if (fBug==1) {
713       // correzione bug LAYER 5  SSD temporanea..
714       int ssdidx=AliITSAlignMille2Module::GetIndexFromVolumeID(p.GetVolumeID());
715       if (ssdidx>=500 && ssdidx<1248) {
716         int ladder=(ssdidx-500)%22;
717         if (ladder==18) p.SetVolumeID(AliITSAlignMille2Module::GetVolumeIDFromIndex(ssdidx+1));
718         if (ladder==19) p.SetVolumeID(AliITSAlignMille2Module::GetVolumeIDFromIndex(ssdidx-1));
719       }
720     }
721     
722     // << RS
723
724     /// get (evenctually prealigned) matrix of sens. vol.
725     TGeoHMatrix *svMatrix = GetMilleModule(intidx[idx[i]])->GetSensitiveVolumeMatrix(p.GetVolumeID());
726     // modify global coordinates according with pre-aligment
727     svMatrix->LocalToMaster(pl,pg);
728     // now rotate in local system
729     hcov.Multiply(&svMatrix->Inverse());
730     hcov.MultiplyLeft(svMatrix);
731     // hcov is back in GLOBAL RF
732     Float_t pcov[6];
733     pcov[0]=hcov.GetRotationMatrix()[0];
734     pcov[1]=hcov.GetRotationMatrix()[1];
735     pcov[2]=hcov.GetRotationMatrix()[2];
736     pcov[3]=hcov.GetRotationMatrix()[4];
737     pcov[4]=hcov.GetRotationMatrix()[5];
738     pcov[5]=hcov.GetRotationMatrix()[8];
739
740     p.SetXYZ(pg[0],pg[1],pg[2],pcov);
741     //    printf("New Gl coordinates of measured point : X=%f  Y=%f  Z=%f \n",pg[0],pg[1],pg[2]);
742     AliDebug(3,Form("New global coordinates of measured point : X=%f  Y=%f  Z=%f \n",pg[0],pg[1],pg[2]));
743     atps->AddPoint(npto,&p);
744     AliDebug(2,Form("Adding point[%d] = ( %f , %f , %f )     volid = %d",npto,atps->GetX()[npto],atps->GetY()[npto],atps->GetZ()[npto],atps->GetVolumeID()[npto] ));
745
746     npto++;
747   }
748
749   return atps;
750 }
751
752
753
754 AliTrackPointArray *AliITSAlignMille2::SortTrack(const AliTrackPointArray *atp) {
755   /// sort alitrackpoints w.r.t. global Y direction
756   AliTrackPointArray *atps=NULL;
757   Int_t idx[20];
758   Int_t npts=atp->GetNPoints();
759   AliTrackPoint p;
760   atps=new AliTrackPointArray(npts);
761
762   TMath::Sort(npts,atp->GetY(),idx);
763
764   for (int i=0; i<npts; i++) {
765     atp->GetPoint(p,idx[i]);
766     atps->AddPoint(i,&p);
767     AliDebug(2,Form("Point[%d] = ( %f , %f , %f )     volid = %d",i,atps->GetX()[i],atps->GetY()[i],atps->GetZ()[i],atps->GetVolumeID()[i] ));
768   }
769   return atps;
770 }
771
772
773 Int_t AliITSAlignMille2::InitModuleParams() {
774   /// initialize geometry parameters for a given detector
775   /// for current cluster (fCluster)
776   /// fGlobalInitParam[] is set as:
777   ///    [tx,ty,tz,psi,theta,phi]
778   ///    (old was [tx,ty,tz,theta,psi,phi] ROOT's angles...)
779   /// *** At the moment: using Raffalele's angles definition ***
780   ///
781   /// return 0 if success
782
783   if (!fGeoManager) {
784     AliInfo("ITS geometry not initialized!");
785     return -1;
786   }
787
788   // now 'voluid' is the volumeID of a SENSITIVE VOLUME (coming from a cluster)
789
790   // set the internal index (index in module list)
791   UShort_t voluid=fCluster.GetVolumeID();
792   //
793   // IT IS VERY IMPORTANT: start from the end of the list, where the childs are located !!!
794   Int_t k=fNModules-1;
795   fCurrentModule = 0;
796   while (k>=0 && ! (fCurrentModule=GetMilleModule(k))->IsIn(voluid)) {
797      // VERY IMPORTANT: if the sensors were explicitly provided, don't look in the supermodules  
798     if (fSensorsIn && fCurrentModule->GetVolumeID() > kMaxITSSensVID) {k=-1; break;} 
799     k--; 
800   }
801   if (k<0) return -3;
802   //  fCurrentModule = GetMilleModule(k);
803   //
804   fModuleInitParam[0] = 0.0;
805   fModuleInitParam[1] = 0.0;
806   fModuleInitParam[2] = 0.0;
807   fModuleInitParam[3] = 0.0; // psi   (X)
808   fModuleInitParam[4] = 0.0; // theta (Y)
809   fModuleInitParam[5] = 0.0; // phi   (Z)
810   
811   /// get (evenctually prealigned) matrix of sens. vol.
812   TGeoHMatrix *svMatrix = fCurrentModule->GetSensitiveVolumeMatrix(voluid);
813   
814   fMeasGlo[0] = fCluster.GetX();
815   fMeasGlo[1] = fCluster.GetY();
816   fMeasGlo[2] = fCluster.GetZ(); 
817   svMatrix->MasterToLocal(fMeasGlo,fMeasLoc);  
818   AliDebug(2,Form("Local coordinates of measured point : X=%f  Y=%f  Z=%f \n",fMeasLoc[0] ,fMeasLoc[1] ,fMeasLoc[2] ));
819   
820   // set stdev from cluster
821   TGeoHMatrix hcov;
822   Double_t hcovel[9];
823   hcovel[0]=double(fCluster.GetCov()[0]);
824   hcovel[1]=double(fCluster.GetCov()[1]);
825   hcovel[2]=double(fCluster.GetCov()[2]);
826   hcovel[3]=double(fCluster.GetCov()[1]);
827   hcovel[4]=double(fCluster.GetCov()[3]);
828   hcovel[5]=double(fCluster.GetCov()[4]);
829   hcovel[6]=double(fCluster.GetCov()[2]);
830   hcovel[7]=double(fCluster.GetCov()[4]);
831   hcovel[8]=double(fCluster.GetCov()[5]);
832   hcov.SetRotation(hcovel);
833   // now rotate in local system
834   hcov.Multiply(svMatrix);
835   hcov.MultiplyLeft(&svMatrix->Inverse());
836
837   // set local sigmas
838   fSigmaLoc[0] = TMath::Sqrt(TMath::Abs(hcov.GetRotationMatrix()[0]));
839   fSigmaLoc[1] = TMath::Sqrt(TMath::Abs(hcov.GetRotationMatrix()[4])); // RS
840   fSigmaLoc[2] = TMath::Sqrt(TMath::Abs(hcov.GetRotationMatrix()[8]));
841
842   // set minimum value for SigmaLoc to 10 micron 
843   if (fSigmaLoc[0]<0.0010) fSigmaLoc[0]=0.0010;
844   if (fSigmaLoc[2]<0.0010) fSigmaLoc[2]=0.0010;
845
846   // multiply local sigmas by global and module specific factor 
847   for (int i=3;i--;) fSigmaLoc[i] *= fSigmaFactor[i]*fCurrentModule->GetSigmaFactor(i);
848   //
849   AliDebug(2,Form("Setting StDev from CovMat : fSigmaLocX=%g  fSigmaLocY=%g fSigmaLocZ=%g \n",fSigmaLoc[0] ,fSigmaLoc[1] ,fSigmaLoc[2] ));
850    
851   return 0;
852 }
853
854 void AliITSAlignMille2::Print(Option_t*) const 
855 {
856   ///
857   printf("*** AliMillepede for ITS ***\n");
858   printf("    Number of defined super modules: %d\n",fNModules);
859   printf("    Obtained parameters refer to %s Deltas\n",fUseGlobalDelta ? "GLOBAL":"LOCAL");
860   //
861   if (fGeoManager)
862     printf("    geometry loaded from %s\n",fGeometryFileName.Data());
863   else
864     printf("    geometry not loaded\n");
865   //  
866   if (fUsePreAlignment) 
867     printf("    using prealignment from %s \n",fPreAlignmentFileName.Data());
868   else
869     printf("    prealignment not used\n");    
870   //
871   //
872   if (fBOn) 
873     printf("    B Field set to %f T - using Riemann's helices\n",fBField);
874   else
875     printf("    B Field OFF - using straight lines \n");
876   //
877   if (fRequirePoints) printf("    Required points in tracks:\n");
878   for (Int_t i=0; i<6; i++) {
879     if (fNReqLayUp[i]>0) printf("        Layer %d : %d points with Y>0\n",i+1,fNReqLayUp[i]);
880     if (fNReqLayDown[i]>0) printf("        Layer %d : %d points with Y<0\n",i+1,fNReqLayDown[i]);
881     if (fNReqLay[i]>0) printf("        Layer %d : %d points \n",i+1,fNReqLay[i]);
882   }
883   for (Int_t i=0; i<3; i++) {
884     if (fNReqDetUp[i]>0) printf("        Detector %d : %d points with Y>0\n",i+1,fNReqDetUp[i]);
885     if (fNReqDetDown[i]>0) printf("        Detector %d : %d points with Y<0\n",i+1,fNReqDetDown[i]);
886     if (fNReqDet[i]>0) printf("        Detector %d : %d points \n",i+1,fNReqDet[i]);
887   }
888   //  
889   printf("\n    Millepede configuration parameters:\n");
890   printf("        init parsig for translations  : %.4f\n",fParSigTranslations);
891   printf("        init parsig for rotations     : %.4f\n",fParSigRotations);
892   printf("        init value for chi2 cut       : %.4f\n",fStartFac);
893   printf("        first iteration cut value     : %.4f\n",fResCutInitial);
894   printf("        other iterations cut value    : %.4f\n",fResCut);
895   printf("        number of stddev for chi2 cut : %d\n",fNStdDev);
896   printf("        mult. fact. for local sigmas  : %.4f %.4f %.4f\n",fSigmaFactor[0],fSigmaFactor[1],fSigmaFactor[2]);
897
898   printf("List of defined modules:\n");
899   printf("  intidx\tindex\tvoluid\tname\n");
900   for (int i=0; i<fNModules; i++) {
901     AliITSAlignMille2Module* md = GetMilleModule(i); 
902     printf("  %d\t%d\t%d\t%s\n",i,md->GetIndex(),md->GetVolumeID(),md->GetName());
903   }
904 }
905
906 AliITSAlignMille2Module  *AliITSAlignMille2::GetMilleModuleByVID(UShort_t voluid) const
907 {
908   // return pointer to a define supermodule
909   // return NULL if error
910   Int_t i=IsDefined(voluid);
911   if (i<0) return NULL;
912   return GetMilleModule(i);
913 }
914
915 Bool_t AliITSAlignMille2::InitRiemanFit() 
916 {
917   // Initialize Riemann Fitter for current track
918   // return kFALSE if error
919
920   if (!fBOn) return kFALSE;
921
922   Int_t npts=0;
923   AliTrackPoint ap;
924   npts = fTrack->GetNPoints();
925   AliDebug(3,Form("Fitting track with %d points",npts));
926
927   fRieman->Reset();
928   fRieman->SetTrackPointArray(fTrack);
929
930   TArrayI ai(npts);
931   for (Int_t ipt=0; ipt<npts; ipt++) ai[ipt]=fTrack->GetVolumeID()[ipt];
932   
933   // fit track with 5 params in his own tracking-rotated reference system
934   // xc = -p[1]/p[0];
935   // yc = 1/p[0];
936   // R  = sqrt( x0*x0 + y0*y0 - y0*p[2]);
937   if (!fRieman->Fit(&ai,NULL,(AliGeomManager::ELayerID)1,(AliGeomManager::ELayerID)6)) {
938     return kFALSE;
939   }
940
941   for (int i=0; i<5; i++)
942     fLocalInitParam[i] = fRieman->GetParam()[i];
943   
944   return kTRUE;
945 }
946
947 Bool_t fullErr2D = kTRUE;
948
949 void trackFit2D(Int_t &, Double_t *, double &chi2, double *par, int)
950 {
951   const double kTiny = 1.e-14;
952   chi2 = 0;
953   static AliTrackPoint pnt;
954   //
955   enum {kAX,kAZ,kBX,kBZ};
956   enum {kXX=0,kXY=1,kXZ=2,kYX=kXY,kYY=3,kYZ=4,kZX=kXZ,kZY=kYZ,kZZ=5};
957   //
958   AliITSAlignMille2* alig = AliITSAlignMille2::GetInstance();
959   AliTrackPointArray* track = alig->GetCurrentTrack();
960   //
961   int npts = track->GetNPoints();
962   for (int ip=0;ip<npts;ip++) {
963     track->GetPoint(pnt,ip);
964     const float *cov = pnt.GetCov();
965     double y  = pnt.GetY();
966     double dx = pnt.GetX() - (par[kAX]+y*par[kBX]);
967     double dz = pnt.GetZ() - (par[kAZ]+y*par[kBZ]);
968     double xxe = cov[kXX];
969     double zze = cov[kZZ];
970     double xze = cov[kXZ];
971     //
972     if (fullErr2D) {
973       xxe += par[kBX]*par[kBX]*cov[kYY]-2.*par[kBX]*cov[kXY];
974       zze += par[kBZ]*par[kBZ]*cov[kYY]-2.*par[kBZ]*cov[kZY];
975       xze += par[kBX]*par[kBZ]*cov[kYY]-cov[kYZ]*par[kBZ]-cov[kXY]*par[kBX];
976     }
977     //
978     double det = xxe*zze - xze*xze;
979     if (det<kTiny) {
980       printf("Negative diag. error (det=%+e) |sxx:%+e szz:%+e sxz:%+e| bx:%+e bz:%+e|\n"
981              "Discarding correlation term\n",det,xxe,zze,xze,par[kBX],par[kBZ]);
982       xxe = cov[kXX];
983       zze = cov[kZZ];
984       xze = cov[kXZ];
985       fullErr2D = kFALSE;
986     }
987     double xxeI = zze/det;
988     double zzeI = xxe/det;
989     double xzeI =-xze/det;
990     //
991     chi2 += dx*dx*xxeI + dz*dz*zzeI + 2.*dx*dz*xzeI;
992     // 
993     //    printf("%d | %+e %+e %+e %+e %+e -> %+e\n",ip,dx,dz,xxeI,zzeI,xzeI,  chi2);
994   }
995   //
996 }
997
998 void AliITSAlignMille2::InitTrackParams(int meth) 
999 {
1000   /// initialize local parameters with different methods
1001   /// for current track (fTrack)
1002   Int_t npts=0;
1003   AliTrackPoint ap;
1004   double sX=0,sXY=0,sZ=0,sZY=0,sY=0,sYY=0,det=0;
1005   // simple linear interpolation
1006   // get local starting parameters (to be substituted by ESD track parms)
1007   // local parms (fLocalInitParam[]) are:
1008   //      [0] = global x coord. of straight line intersection at y=0 plane
1009   //      [1] = global z coord. of straight line intersection at y=0 plane
1010   //      [2] = px/py  
1011   //      [3] = pz/py
1012   // test #1: linear fit in x(y) and z(y)
1013   npts = fTrack->GetNPoints();
1014   AliDebug(3,Form("*** initializing track with %d points ***",npts));
1015   for (int i=npts;i--;) {
1016     sY  += fTrack->GetY()[i];
1017     sYY += fTrack->GetY()[i]*fTrack->GetY()[i];
1018     sX  += fTrack->GetX()[i];
1019     sXY += fTrack->GetX()[i]*fTrack->GetY()[i];
1020     sZ  += fTrack->GetZ()[i];
1021     sZY += fTrack->GetZ()[i]*fTrack->GetY()[i];
1022   }
1023   det = sYY*npts-sY*sY;
1024   if (det==0) det = 1E-20;
1025   fLocalInitParam[0] = (sX*sYY-sY*sXY)/det;
1026   fLocalInitParam[2] = (sXY*npts-sY*sX)/det;
1027   //
1028   fLocalInitParam[1] = (sZ*sYY-sY*sZY)/det;
1029   fLocalInitParam[3] = (sZY*npts-sY*sZ)/det;
1030   //AliDebug(2,Form("X = p0gx + ugx*Y : p0gx = %f +- %f    ugx = %f +- %f\n",fLocalInitParam[0],f1->GetParError(0),fLocalInitParam[2],f1->GetParError(1)));
1031   AliDebug(2,Form("X = p0gx + ugx*Y : p0gx = %f    ugx = %f\n",fLocalInitParam[0],fLocalInitParam[2]));
1032   //
1033   if (meth==1) return;
1034   //
1035   // perform full fit accounting for cov.matrix
1036   static TVirtualFitter *minuit = 0;
1037   static Double_t step[5]   = {1E-3,1E-3,1E-4,1E-4,1E-5};
1038   static Double_t arglist[10];
1039   //
1040   if (!minuit) {
1041     minuit = TVirtualFitter::Fitter(0,4);
1042     minuit->SetFCN(trackFit2D);
1043     arglist[0] = 1;
1044     minuit->ExecuteCommand("SET ERR",arglist, 1);
1045     //
1046     arglist[0] = -1;
1047     minuit->ExecuteCommand("SET PRINT",arglist,1);
1048     //
1049   }
1050   //
1051   minuit->SetParameter(0, "ax",   fLocalInitParam[0], step[0], 0,0);
1052   minuit->SetParameter(1, "az",   fLocalInitParam[1], step[1], 0,0);
1053   minuit->SetParameter(2, "bx",   fLocalInitParam[2], step[2], 0,0);
1054   minuit->SetParameter(3, "bz",   fLocalInitParam[3], step[3], 0,0);
1055   //
1056   arglist[0] = 1000; // number of function calls 
1057   arglist[1] = 0.001; // tolerance 
1058   fullErr2D = kTRUE;
1059   minuit->ExecuteCommand("MIGRAD",arglist,2);
1060   //
1061   for (int i=0;i<4;i++) fLocalInitParam[i] = minuit->GetParameter(i);
1062   for (int i=0;i<4;i++) for (int j=0;j<4;j++) fLocalInitParEr[i][j] = minuit->GetCovarianceMatrixElement(i,j);
1063   //
1064 }
1065
1066
1067 Int_t AliITSAlignMille2::IsDefined(UShort_t voluid) const
1068 {
1069   // checks if supermodule 'voluid' is defined and return the internal index
1070   // return -1 if error
1071   Int_t k = fNModules-1;
1072   while ( k>=0 && !(voluid==GetMilleModule(k)->GetVolumeID()) ) k--;  
1073   if (k<0) return -1; 
1074   return k;
1075 }
1076
1077 Int_t AliITSAlignMille2::IsContained(UShort_t voluid) const
1078 {
1079   // checks if the sensitive module 'voluid' is contained inside a supermodule and return the internal index of the last identified supermodule
1080   // return -1 if error
1081   if (AliITSAlignMille2Module::GetIndexFromVolumeID(voluid)<0) return -1;
1082   Int_t k=fNModules-1;
1083   while (k>=0 && !(GetMilleModule(k)->IsIn(voluid)) ) k--;  
1084   if (k<0) return -1; 
1085   return k;
1086 }
1087
1088 Bool_t AliITSAlignMille2::CheckVolumeID(UShort_t voluid) const 
1089 {
1090   /// check if a sensitive volume is contained inside one of the defined supermodules
1091   Int_t k=fNModules-1;
1092   while (k>=0 && !(GetMilleModule(k)->IsIn(voluid)) ) k--;  
1093   if (k>=0) return kTRUE;
1094   return kFALSE;
1095 }
1096
1097 Int_t AliITSAlignMille2::CheckCurrentTrack() {
1098   /// checks if AliTrackPoints belongs to defined modules
1099   /// return number of good poins
1100   /// return 0 if not enough points
1101
1102   Int_t npts = fTrack->GetNPoints();
1103   Int_t ngoodpts=0;
1104   // debug points
1105   for (int j=0; j<npts; j++) {
1106     UShort_t voluid = fTrack->GetVolumeID()[j];    
1107     if (CheckVolumeID(voluid)) {
1108       ngoodpts++;
1109     }
1110   }
1111
1112   if (ngoodpts<fMinNPtsPerTrack) return 0;
1113
1114   return ngoodpts;
1115 }
1116
1117 Int_t AliITSAlignMille2::ProcessTrack(const AliTrackPointArray *track) {
1118   /// Process track; Loop over hits and set local equations
1119   /// here 'track' is a AliTrackPointArray
1120   /// return 0 if success;
1121   
1122   if (!fIsMilleInit) {
1123     AliInfo("Millepede has not been initialized!");
1124     return -1;
1125   }
1126
1127   Int_t npts = track->GetNPoints();
1128   AliDebug(2,Form("*** Input track with %d points ***",npts));
1129
1130   // preprocessing of the input track: keep only points in defined volumes,
1131   // move points if prealignment is set, sort by Yglo if required
1132   
1133   fTrack=PrepareTrack(track);
1134   if (!fTrack) return -1;
1135
1136   npts = fTrack->GetNPoints();
1137   AliDebug(2,Form("*** Processing prepared track with %d points ***",npts));
1138   
1139   if (!fBOn) { // straight lines  
1140     // set local starting parameters (to be substituted by ESD track parms)
1141     // local parms (fLocalInitParam[]) are:
1142     //      [0] = global x coord. of straight line intersection at y=0 plane
1143     //      [1] = global z coord. of straight line intersection at y=0 plane
1144     //      [2] = px/py  
1145     //      [3] = pz/py
1146     InitTrackParams(fInitTrackParamsMeth);  
1147   } 
1148   else {
1149     // local parms (fLocalInitParam[]) are the Riemann Fitter params
1150     if (!InitRiemanFit()) {
1151       AliInfo("Riemann fit failed! skipping this track...");
1152       delete fTrack;
1153       fTrack=NULL;
1154       return -5;
1155     }
1156   }
1157   
1158   Int_t nloceq=0;
1159   static Mille2Data md[100];
1160   
1161   for (Int_t ipt=0; ipt<npts; ipt++) {
1162     fTrack->GetPoint(fCluster,ipt);
1163     AliDebug(2,Form("\n--- processing point %d --- \n",ipt));    
1164
1165     // set geometry parameters for the the current module
1166     if (InitModuleParams()) continue;
1167     AliDebug(2,Form("    VolID=%d  Index=%d  InternalIdx=%d  symname=%s\n", 
1168                     track->GetVolumeID()[ipt], fCurrentModule->GetIndex(),
1169                     fCurrentModule->GetUniqueID(), AliGeomManager::SymName(track->GetVolumeID()[ipt]) ));
1170     AliDebug(2,Form("    Preprocessed Point = ( %f , %f , %f ) \n",fCluster.GetX(),fCluster.GetY(),fCluster.GetZ()));
1171     
1172     if (!AddLocalEquation(md[nloceq])) nloceq++;    
1173     else fTotBadLocEqPoints++;
1174   } // end loop over points
1175   //
1176   delete fTrack;
1177   fTrack=NULL;
1178   // not enough good points!
1179   if (nloceq<fMinNPtsPerTrack) return -1;
1180   //
1181   // finally send local equations to millepede
1182   SetLocalEquations(md,nloceq);
1183   fMillepede->SaveRecordData(); // RRR
1184   //
1185   return 0;
1186 }
1187
1188 Int_t AliITSAlignMille2::CalcIntersectionPoint(Double_t *lpar, Double_t *gpar) {
1189   /// calculate track intersection point in local coordinates
1190   /// according with a given set of parameters (local(4) and global(6))
1191   /// and fill fPintLoc/Glo
1192   ///    local are:   pgx0, pgz0, ugx, ugz   OR   riemann fitters pars
1193   ///    global are:  tx,ty,tz,psi,theta,phi (Raff's delta angles in deg.)
1194   /// return 0 if success
1195   
1196   AliDebug(3,Form("lpar = %g %g %g %g %g\ngpar= %g %g %g %g %g %g\n",lpar[0],lpar[1],lpar[2],lpar[3],lpar[4],gpar[0],gpar[1],gpar[2],gpar[3],gpar[4],gpar[5]));
1197   AliDebug(3,Form("deltalpar = %g %g %g %g %g\n",lpar[0]-fLocalInitParam[0],lpar[1]-fLocalInitParam[1],lpar[2]-fLocalInitParam[2],lpar[3]-fLocalInitParam[3],lpar[4]-fLocalInitParam[4]));
1198
1199   
1200   // prepare the TGeoHMatrix
1201   TGeoHMatrix *fTempHMat = fCurrentModule->GetSensitiveVolumeModifiedMatrix(fCluster.GetVolumeID(),gpar,
1202                                                                             !fUseGlobalDelta);
1203   if (!fTempHMat) return -1;
1204   
1205   Double_t v0g[3]; // vector with straight line direction in global coord.
1206   Double_t p0g[3]; // point of the straight line (glo)
1207   
1208   if (fBOn) { // B FIELD!
1209     AliTrackPoint prf; 
1210     for (int ip=0; ip<5; ip++)
1211       fRieman->SetParam(ip,lpar[ip]);
1212
1213     if (!fRieman->GetPCA(fCluster,prf))  {
1214       AliInfo(Form("error in GetPCA for point %d",fCluster.GetVolumeID()));
1215       return -3;
1216     }
1217     // now determine straight line passing tangent to fit curve at prf
1218     // ugx = dX/dY_glo = DeltaX/DeltaY_glo
1219     // mo' P1=(X,Y,Z)_glo_prf
1220     //       => (x,y,Z)_trk_prf ruotando di alpha...
1221     Double_t alpha=fRieman->GetAlpha();
1222     Double_t x1g = prf.GetX();
1223     Double_t y1g = prf.GetY();
1224     Double_t z1g = prf.GetZ();
1225     Double_t x1t =  x1g*TMath::Cos(alpha) + y1g*TMath::Sin(alpha);
1226     Double_t y1t = -x1g*TMath::Sin(alpha) + y1g*TMath::Cos(alpha);
1227     Double_t z1t =  z1g;    
1228
1229     Double_t x2t = x1t+1.0;
1230     Double_t y2t = y1t+fRieman->GetDYat(x1t);
1231     Double_t z2t = z1t+fRieman->GetDZat(x1t);
1232     Double_t x2g =  x2t*TMath::Cos(alpha) - y2t*TMath::Sin(alpha);
1233     Double_t y2g =  x2t*TMath::Sin(alpha) + y2t*TMath::Cos(alpha);
1234     Double_t z2g =  z2t;  
1235
1236     AliDebug(3,Form("Riemann frame:  fAlpha = %f  =  %f  ",alpha,alpha*180./TMath::Pi()));
1237     AliDebug(3,Form("   prf_glo=( %f , %f , %f )  prf_rf=( %f , %f , %f )\n", x1g,y1g,z1g, x1t,y1t,z1t));
1238     AliDebug(3,Form("   mov_glo=( %f , %f , %f )      rf=( %f , %f , %f )\n",x2g,y2g,z2g, x2t,y2t,z2t));
1239         
1240     if (TMath::Abs(y2g-y1g)<1e-15) {
1241       AliInfo("DeltaY=0! Cannot proceed...");
1242       return -1;
1243     }
1244     // ugx,1,ugz
1245     v0g[0] = (x2g-x1g)/(y2g-y1g);
1246     v0g[1]=1.0;
1247     v0g[2] = (z2g-z1g)/(y2g-y1g);
1248     
1249     // point: just keep prf
1250     p0g[0]=x1g;
1251     p0g[1]=y1g;
1252     p0g[2]=z1g;
1253   }  
1254   else { // staight line
1255     // vector of initial straight line direction in glob. coord
1256     v0g[0]=lpar[2];
1257     v0g[1]=1.0;
1258     v0g[2]=lpar[3];
1259     
1260     // intercept in yg=0 plane in glob coord
1261     p0g[0]=lpar[0];
1262     p0g[1]=0.0;
1263     p0g[2]=lpar[1];
1264   }
1265   AliDebug(3,Form("Line vector: ( %f , %f , %f )  point:( %f , %f , %f )\n",v0g[0],v0g[1],v0g[2],p0g[0],p0g[1],p0g[2]));
1266   
1267   // same in local coord.
1268   Double_t p0l[3],v0l[3];
1269   fTempHMat->MasterToLocalVect(v0g,v0l);
1270   fTempHMat->MasterToLocal(p0g,p0l);
1271   
1272   if (TMath::Abs(v0l[1])<1e-15) {
1273     AliInfo("Track Y direction in local frame is zero! Cannot proceed...");
1274     return -1;
1275   }
1276   
1277   // local intersection point
1278   fPintLoc[0] = p0l[0] - (v0l[0]/v0l[1])*p0l[1];
1279   fPintLoc[1] = 0;
1280   fPintLoc[2] = p0l[2] - (v0l[2]/v0l[1])*p0l[1];
1281   
1282   // global intersection point
1283   fTempHMat->LocalToMaster(fPintLoc,fPintGlo);
1284   AliDebug(3,Form("Intesect. point: L( %f , %f , %f )  G( %f , %f , %f )\n",fPintLoc[0],fPintLoc[1],fPintLoc[2],fPintGlo[0],fPintGlo[1],fPintGlo[2]));
1285   
1286   return 0;
1287 }
1288
1289 Int_t AliITSAlignMille2::CalcDerivatives(Int_t paridx, Bool_t islpar) {
1290   /// calculate numerically (ROOT's style) the derivatives for
1291   /// local X intersection  and local Z intersection
1292   /// parlist: local  (islpar=kTRUE)  pgx0, pgz0, ugx0, ugz0  OR riemann's params
1293   ///          global (islpar=kFALSE) tx, ty, tz, psi, theta, phi (Raf's angles in deg)
1294   /// return 0 if success
1295   
1296   // copy initial parameters
1297   Double_t lpar[ITSMILLE2_NLOCAL];
1298   Double_t gpar[ITSMILLE2_NPARCH];
1299   for (Int_t i=0; i<ITSMILLE2_NLOCAL; i++) lpar[i]=fLocalInitParam[i];
1300   for (Int_t i=0; i<ITSMILLE2_NPARCH; i++) gpar[i]=fModuleInitParam[i];
1301
1302   // trial with fixed dpar...
1303   Double_t dpar=0.0;
1304
1305   if (islpar) { // track parameters
1306     //dpar=fLocalInitParam[paridx]*0.001;
1307     // set minimum dpar
1308     if (!fBOn) {
1309       if (paridx<2) dpar=1.0e-4; // translations
1310       else dpar=1.0e-6; // direction
1311     }
1312     else { // B Field
1313       // pepo: proviamo con 1/1000, poi evenctually 1/100...
1314       Double_t dfrac=0.01;
1315       switch(paridx) {
1316       case 0:
1317         // RMS cosmics: 1e-4
1318         dpar = TMath::Max(1.0e-6,TMath::Abs(fLocalInitParam[paridx]*dfrac)); 
1319         break;
1320       case 1: 
1321         // RMS cosmics: 0.2
1322         dpar = TMath::Max(0.002,TMath::Abs(fLocalInitParam[paridx]*dfrac)); 
1323         break;
1324       case 2: 
1325         // RMS cosmics: 9
1326         dpar = TMath::Max(0.09,TMath::Abs(fLocalInitParam[paridx]*dfrac)); 
1327         break;
1328       case 3: 
1329         // RMS cosmics: 7
1330         dpar = TMath::Max(0.07,TMath::Abs(fLocalInitParam[paridx]*dfrac)); 
1331         break;
1332       case 4: 
1333         // RMS cosmics: 0.3
1334         dpar = TMath::Max(0.003,TMath::Abs(fLocalInitParam[paridx]*dfrac)); 
1335         break;
1336       }
1337     }
1338   }
1339   else { // alignment global parameters
1340     //dpar=fModuleInitParam[paridx]*0.001;
1341     if (paridx<3) dpar=1.0e-4; // translations
1342     else dpar=1.0e-2; // angles    
1343   }
1344
1345   AliDebug(3,Form("+++ using dpar=%g",dpar));
1346   
1347   // calculate derivative ROOT's like:
1348   //  using f(x+h),f(x-h),f(x+h/2),f(x-h2)...
1349   Double_t pintl1[3]; // f(x-h)
1350   Double_t pintl2[3]; // f(x-h/2)
1351   Double_t pintl3[3]; // f(x+h/2)
1352   Double_t pintl4[3]; // f(x+h)
1353     
1354   // first values
1355   if (islpar) lpar[paridx] -= dpar;
1356   else gpar[paridx] -= dpar;
1357   if (CalcIntersectionPoint(lpar, gpar)) return -2;
1358   for (Int_t i=0; i<3; i++) pintl1[i]=fPintLoc[i];
1359
1360   // second values
1361   if (islpar) lpar[paridx] += dpar/2;
1362   else gpar[paridx] += dpar/2;
1363   if (CalcIntersectionPoint(lpar, gpar)) return -2;
1364   for (Int_t i=0; i<3; i++) pintl2[i]=fPintLoc[i];
1365
1366   // third values
1367   if (islpar) lpar[paridx] += dpar;
1368   else gpar[paridx] += dpar;
1369   if (CalcIntersectionPoint(lpar, gpar)) return -2;
1370   for (Int_t i=0; i<3; i++) pintl3[i]=fPintLoc[i];
1371
1372   // fourth values
1373   if (islpar) lpar[paridx] += dpar/2;
1374   else gpar[paridx] += dpar/2;
1375   if (CalcIntersectionPoint(lpar, gpar)) return -2;
1376   for (Int_t i=0; i<3; i++) pintl4[i]=fPintLoc[i];
1377
1378   Double_t h2 = 1./(2.*dpar);
1379   Double_t d0 = pintl4[0]-pintl1[0];
1380   Double_t d2 = 2.*(pintl3[0]-pintl2[0]);
1381   fDerivativeLoc[0] = h2*(4*d2 - d0)/3.;
1382   if (TMath::Abs(fDerivativeLoc[0]) < 1.0e-9) fDerivativeLoc[0] = 0.0;
1383
1384   d0 = pintl4[2]-pintl1[2];
1385   d2 = 2.*(pintl3[2]-pintl2[2]);
1386   fDerivativeLoc[2] = h2*(4*d2 - d0)/3.;
1387   if (TMath::Abs(fDerivativeLoc[2]) < 1.0e-9) fDerivativeLoc[2]=0.0;
1388
1389   AliDebug(3,Form("\n+++ derivatives +++ \n"));
1390   AliDebug(3,Form("+++ dXLoc/dpar = %g +++\n",fDerivativeLoc[0]));
1391   AliDebug(3,Form("+++ dZLoc/dpar = %g +++\n\n",fDerivativeLoc[0]));
1392   
1393   return 0;
1394 }
1395
1396
1397 Int_t AliITSAlignMille2::AddLocalEquation(Mille2Data &m) {
1398   /// Define local equation for current cluster in X and Z coor.
1399   /// and store them to memory
1400   /// return 0 if success
1401   int nLev = 0; 
1402   // store first interaction point
1403   if (CalcIntersectionPoint(fLocalInitParam, fModuleInitParam)) return -4;  
1404   for (Int_t i=0; i<3; i++) fPintLoc0[i]=fPintLoc[i];
1405   AliDebug(2,Form("Intesect. point: L( %f , %f , %f )",fPintLoc[0],fPintLoc[1],fPintLoc[2]));
1406   
1407   // calculate local derivatives numerically
1408   Bool_t zeroX = kTRUE;
1409   Bool_t zeroZ = kTRUE;
1410   //
1411   for (Int_t i=0; i<fNLocal; i++) {
1412     if (CalcDerivatives(i,kTRUE)) return -1;
1413     m.derlocX[i] = fDerivativeLoc[0];
1414     m.derlocZ[i] = fDerivativeLoc[2];
1415     if (zeroX) zeroX = fDerivativeLoc[0]==0;
1416     if (zeroZ) zeroZ = fDerivativeLoc[2]==0;
1417   }
1418   //  for (Int_t i=0; i<fNLocal; i++) AliDebug(2,Form("Local parameter %d - dXdpar = %g  - dZdpar = %g\n",i,dXdL[i],dZdL[i]));
1419   if (zeroX) {AliInfo("Aborting: zero local X derivatives!"); return -2;}
1420   if (zeroZ) {AliInfo("Aborting: zero local Z derivatives!"); return -2;}
1421   //
1422   //
1423   AliITSAlignMille2Module* endModule = fCurrentModule;
1424   //
1425   do {
1426     if (nLev==0 || !fUseGlobalDelta) zeroX = zeroZ = kTRUE;
1427     int shiftL = nLev*ITSMILLE2_NPARCH;
1428     for (Int_t i=0; i<ITSMILLE2_NPARCH; i++) {
1429       if (nLev==0 || !fUseGlobalDelta) {
1430         if (CalcDerivatives(i,kFALSE)) return -1;
1431         m.dergloX[shiftL + i] = fDerivativeLoc[0];
1432         m.dergloZ[shiftL + i] = fDerivativeLoc[2];
1433         if (zeroX) zeroX = fDerivativeLoc[0]==0;
1434         if (zeroZ) zeroZ = fDerivativeLoc[2]==0;      
1435       } 
1436       else { // for the global delta the deravites of different levels are the same
1437         m.dergloX[shiftL + i] = m.dergloX[shiftL + i - ITSMILLE2_NPARCH];
1438         m.dergloZ[shiftL + i] = m.dergloZ[shiftL + i - ITSMILLE2_NPARCH];
1439       }
1440     }
1441     //  for (Int_t i=0; i<ITSMILLE2_NPARCH; i++) AliDebug(2,Form("Global parameter %d - dXdpar = %g  - dZdpar = %g\n",i,dXdG[i],dZdG[i]));
1442     if (zeroX) {AliInfo("Aborting: zero global X derivatives!");return -2;}
1443     if (zeroZ) {AliInfo("Aborting: zero global Z derivatives!");return -2;}
1444     // set equation for Xloc coordinate
1445     m.moduleIDX[nLev] = fCurrentModule->GetUniqueID();
1446     nLev++;
1447     //
1448   } while( (fCurrentModule=fCurrentModule->GetParent()) );
1449   //
1450   // ok, can copy to m
1451   AliDebug(2,Form("Adding local equation X with fMeas=%.6f  and fSigma=%.6f",(fMeasLoc[0]-fPintLoc0[0]), fSigmaLoc[0]));
1452   m.measX = fMeasLoc[0]-fPintLoc0[0];
1453   m.sigmaX = fSigmaLoc[0];
1454   //
1455   AliDebug(2,Form("Adding local equation Z with fMeas=%.6f  and fSigma=%.6f",(fMeasLoc[2]-fPintLoc0[2]), fSigmaLoc[2]));
1456   m.measZ = fMeasLoc[2]-fPintLoc0[2];
1457   m.sigmaZ = fSigmaLoc[2];
1458   //
1459   m.levFilled = nLev;
1460   fCurrentModule = endModule;
1461   //
1462   return 0;
1463 }
1464
1465 void AliITSAlignMille2::SetLocalEquations(const Mille2Data *marr, Int_t neq) {
1466   /// Set local equations with data stored in m
1467   /// return 0 if success
1468   //
1469   for (Int_t j=0; j<neq; j++) {
1470     //
1471     const Mille2Data &m = marr[j];
1472     //
1473     // set equation for Xloc coordinate
1474     AliDebug(2,Form("setting local equation X with fMeas=%.6f  and fSigma=%.6f",m.measX, m.sigmaX));
1475     for (int i=fNLocal; i--;) SetLocalDerivative( i, m.derlocX[i] );
1476     for (int il=m.levFilled;il--;) {
1477       GetMilleModule(m.moduleIDX[il])->IncNProcessedPoints();
1478       int hlev = m.moduleIDX[il]*ITSMILLE2_NPARCH;       // id of the supermodule 
1479       int llev = il*ITSMILLE2_NPARCH;
1480       for (int i=ITSMILLE2_NPARCH; i--;) SetGlobalDerivative( hlev+i, m.dergloX[llev+i] );
1481     }
1482     //
1483     fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, m.measX, m.sigmaX);  
1484     //
1485     // set equation for Zloc coordinate
1486     AliDebug(2,Form("setting local equation Z with fMeas=%.6f  and fSigma=%.6f",m.measZ, m.sigmaZ));
1487     for (int i=fNLocal; i--;) SetLocalDerivative( i, m.derlocZ[i] );
1488     for (int il=m.levFilled;il--;) {
1489       int hlev = m.moduleIDX[il]*ITSMILLE2_NPARCH;       // id of the supermodule 
1490       int llev = il*ITSMILLE2_NPARCH;
1491       for (int i=ITSMILLE2_NPARCH; i--;) SetGlobalDerivative(hlev+i, m.dergloZ[llev+i] );
1492     }
1493     fMillepede->SetLocalEquation(fGlobalDerivatives, fLocalDerivatives, m.measZ, m.sigmaZ);  
1494   }
1495 }
1496
1497 Int_t AliITSAlignMille2::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls) {
1498   /// Call global fit; Global parameters are stored in parameters
1499   if (!fIsMilleInit) {
1500     AliInfo("Millepede has not been initialized!");
1501     return 0;
1502   }
1503   int res = fMillepede->GlobalFit(parameters,errors,pulls);
1504   AliInfo(Form("%s fitting global parameters!",res ? "Done":"Failed"));
1505   return res;
1506 }
1507
1508 Double_t AliITSAlignMille2::GetParError(Int_t iPar) {
1509   /// Get error of parameter iPar
1510   if (!fIsMilleInit) {
1511     AliInfo("Millepede has not been initialized!");
1512     return 0;
1513   }
1514   Double_t lErr = fMillepede->GetParError(iPar);
1515   return lErr;
1516 }
1517
1518 void AliITSAlignMille2::PrintGlobalParameters() {
1519   /// Print global parameters
1520   if (!fIsMilleInit) {
1521     AliInfo("Millepede has not been initialized!");
1522     return;
1523   }
1524   fMillepede->PrintGlobalParameters();
1525 }
1526
1527 // //_________________________________________________________________________
1528 Int_t AliITSAlignMille2::LoadSuperModuleFile(const Char_t *sfile)
1529
1530   // load definitions of supermodules from a root file
1531   // return 0 if success
1532
1533   TFile *smf=TFile::Open(sfile);
1534   if (!smf->IsOpen()) {
1535     AliInfo(Form("Cannot open supermodule file %s",sfile));
1536     return -1;
1537   }
1538
1539   TClonesArray *sma=(TClonesArray*)smf->Get("ITSMilleSuperModules");
1540   if (!sma) {
1541     AliInfo(Form("Cannot find ITSMilleSuperModules array in file"));
1542     return -2;  
1543   }  
1544   Int_t nsma=sma->GetEntriesFast();
1545   AliInfo(Form("Array of SuperModules with %d entries\n",nsma));
1546   //
1547   Char_t st[250];
1548   char symname[150];
1549   UShort_t volid;
1550   TGeoHMatrix m;
1551   //
1552   for (Int_t i=0; i<nsma; i++) {
1553     AliAlignObjParams *a = (AliAlignObjParams*)sma->UncheckedAt(i);
1554     volid=a->GetVolUID();
1555     strcpy(st,a->GetSymName());
1556     a->GetMatrix(m);
1557     //
1558     sscanf(st,"%s",symname);
1559     // decode module list
1560     char *stp=strstr(st,"ModuleList:");
1561     if (!stp) return -3;
1562     stp += 11;
1563     int idx[2200];
1564     char spp[200]; int jp=0;
1565     char cl[20];
1566     strcpy(st,stp);
1567     int l=strlen(st);
1568     int j=0;
1569     int n=0;
1570     //
1571     while (j<=l) {
1572       if (st[j]==9 || st[j]==32 || st[j]==10 || st[j]==0) {
1573         spp[jp]=0;
1574         jp=0;
1575         if (strlen(spp)) {
1576           int k=strcspn(spp,"-");
1577           if (k<int(strlen(spp))) { // c'e' il -
1578             strcpy(cl,&(spp[k+1]));
1579             spp[k]=0;
1580             int ifrom=atoi(spp); int ito=atoi(cl);
1581             for (int b=ifrom; b<=ito; b++) {
1582               idx[n]=b;
1583               n++;
1584             }
1585           }
1586           else { // numerillo singolo
1587             idx[n]=atoi(spp);
1588             n++;
1589           }
1590         }
1591       }
1592       else {
1593         spp[jp]=st[j];
1594         jp++;
1595       }
1596       j++;
1597     }
1598     UShort_t volidsv[2198];
1599     for (j=0;j<n;j++) {
1600       volidsv[j]=AliITSAlignMille2Module::GetVolumeIDFromIndex(idx[j]);
1601       if (!volidsv[j]) {
1602         AliInfo(Form("Index %d not valid (range 0->%d)",idx[j],kMaxITSSensID));
1603         return -5;
1604       }
1605     }
1606     Int_t smindex=int(2198+volid-14336); // virtual index
1607     fSuperModule.AddAtAndExpand(new AliITSAlignMille2Module(smindex,volid,symname,&m,n,volidsv),fNSuperModules);
1608     //
1609     fNSuperModules++;
1610   }
1611
1612   smf->Close();
1613   //
1614   return 0;
1615 }
1616
1617 //_________________________________________________________________________
1618 void AliITSAlignMille2::ConstrainModuleSubUnits(Int_t idm, Double_t val, UInt_t pattern)
1619 {
1620   // require that sum of modifications for the childs of this module is = val, i.e.
1621   // the internal corrections moves the module as a whole by fixed value (0 by default).
1622   // pattern is the bit pattern for the parameters to constrain
1623   //
1624   static TObjArray childs;
1625   childs.Clear();
1626   //
1627   // build list of childs for this module
1628   int nChilds = 0;
1629   AliITSAlignMille2Module* parent = GetMilleModule(idm);
1630   if (!parent) return;
1631   //
1632   for (int i=fNModules;i--;) {
1633     AliITSAlignMille2Module* child = GetMilleModule(i);
1634     if (child->GetParent() == parent) childs.AddAtAndExpand(child,nChilds++);
1635   }
1636   if (nChilds<1) return;
1637   //
1638   int npc = 0;
1639   for (int ip=0;ip<kNParCh;ip++) {
1640     if ( !((pattern>>ip)&0x1) /*|| !parent->IsFreeDOF(ip)*/) continue;
1641     ResetLocalEquation();
1642     for (int ich=nChilds;ich--;) fGlobalDerivatives[childs[ich]->GetUniqueID()*kNParCh+ip] = 1.0;
1643     AddConstraint(fGlobalDerivatives,val);
1644     npc++;
1645   }
1646   //
1647   AliInfo(Form("Constrained %d params for %d submodules of module #%d: %s",npc,nChilds,idm,parent->GetName()));
1648   //
1649 }
1650
1651 //_________________________________________________________________________
1652 void AliITSAlignMille2::PostConstrainModuleSubUnitsMedian(Int_t idm, Double_t val, UInt_t pattern)
1653 {
1654   // require that median of the modifications for the childs of this module is = val, i.e.
1655   // the internal corrections moves the module as a whole by fixed value (0 by default) 
1656   // module the outliers.
1657   // pattern is the bit pattern for the parameters to constrain
1658   // The difference between the mean and the median will be transfered to the parent
1659   //
1660   static TObjArray childs;
1661   childs.Clear();
1662   //
1663   // build list of childs for this module
1664   int nChilds = 0;
1665   AliITSAlignMille2Module* parent = GetMilleModule(idm);
1666   if (!parent) return;
1667   //
1668   for (int i=fNModules;i--;) {
1669     AliITSAlignMille2Module* child = GetMilleModule(i);
1670     if (child->GetParent() == parent) childs.AddAtAndExpand(child,nChilds++);
1671   }
1672   if (nChilds<1) return;
1673   //
1674   int npc = 0;
1675   double *deltas = fMillepede->GetDeltaPars();
1676   double *tmpArr = new double[nChilds]; 
1677   //
1678   for (int ip=0;ip<kNParCh;ip++) {
1679     if ( !((pattern>>ip)&0x1) /*|| !parent->IsFreeDOF(ip)*/) continue;
1680     // compute the median of the deltas
1681     for (int ich=nChilds;ich--;) tmpArr[ich] = deltas[childs[ich]->GetUniqueID()*kNParCh+ip];
1682     for (int ic0=0;ic0<nChilds;ic0++) // order the deltas 
1683       for (int ic1=ic0+1;ic1<nChilds;ic1++) 
1684         if (tmpArr[ic0]>tmpArr[ic1]) {double tv=tmpArr[ic0]; tmpArr[ic0]=tmpArr[ic1]; tmpArr[ic1]=tv;}
1685     //
1686     int kmed = nChilds/2;
1687     double median = (tmpArr[kmed]+tmpArr[nChilds-kmed-1])/2.;
1688     //
1689     for (int ich=nChilds;ich--;) deltas[childs[ich]->GetUniqueID()*kNParCh+ip] -= median - val;
1690     deltas[parent->GetUniqueID()*kNParCh+ip] += median - val;
1691     npc++;
1692   }
1693   delete[] tmpArr;  
1694   //
1695   AliInfo(Form("Applied median constraint to %d params for %d submodules of module #%d: %s",npc,nChilds,idm,parent->GetName()));
1696   //
1697 }
1698
1699 //_________________________________________________________________________
1700 void AliITSAlignMille2::ConstrainOrphans(Double_t val, UInt_t pattern)
1701 {
1702   // require that median of the modifications for the supermodules which have no parents is = val, i.e.
1703   // the corrections moves the whole setup by fixed value (0 by default) modulo the outliers.
1704   // pattern is the bit pattern for the parameters to constrain
1705   //
1706   static TObjArray modSet;
1707   modSet.Clear();
1708   //
1709   // build list of childs for this module
1710   int nModules = 0;
1711   int *nFree = new int[kNParCh];
1712   for (int i=0;i<kNParCh;i++) nFree[i] = 0;
1713   //
1714   for (int i=fNModules;i--;) {
1715     AliITSAlignMille2Module* module = GetMilleModule(i);
1716     if (module->GetParent()) continue;   // skip this    
1717     for (int ip=0;ip<kNParCh;ip++) if ( ((pattern>>ip)&0x1) && module->IsFreeDOF(ip)) nFree[ip]++;
1718     modSet.AddAtAndExpand(module,nModules++);
1719   }
1720   if (nModules<1) return;
1721   //
1722   int npc = 0;
1723   for (int ip=0;ip<kNParCh;ip++) {
1724     if (nFree[ip]<1) continue; // nothing to do
1725     ResetLocalEquation();
1726     for (int ich=nModules;ich--;) {
1727       AliITSAlignMille2Module* module = (AliITSAlignMille2Module*) modSet[ich];
1728       if ( !((pattern>>ip)&0x1) || !module->IsFreeDOF(ip)) continue;
1729       fGlobalDerivatives[module->GetUniqueID()*kNParCh+ip] = 1.0;
1730     }
1731     AddConstraint(fGlobalDerivatives,val);
1732     npc++;
1733   }
1734   //
1735   delete[] nFree;
1736   AliInfo(Form("Constrained %d params for %d orphan modules",npc,nModules));
1737   //
1738 }
1739 //_________________________________________________________________________
1740 void AliITSAlignMille2::PostConstrainOrphansMedian(Double_t val, UInt_t pattern)
1741 {
1742   // require that sum of modifications for the supermodules which have no parents is = val, i.e.
1743   // the corrections moves the whole setup by fixed value (0 by default).
1744   // pattern is the bit pattern for the parameters to constrain
1745   //
1746   static TObjArray modSet;
1747   modSet.Clear();
1748   //
1749   // build list of childs for this module
1750   int nModules = 0;
1751   int *nFree = new int[kNParCh];
1752   for (int i=0;i<kNParCh;i++) nFree[i] = 0;
1753   //
1754   for (int i=fNModules;i--;) {
1755     AliITSAlignMille2Module* module = GetMilleModule(i);
1756     if (module->GetParent()) continue;   // skip this    
1757     for (int ip=0;ip<kNParCh;ip++) if ( ((pattern>>ip)&0x1) && module->IsFreeDOF(ip)) nFree[ip]++;
1758     modSet.AddAtAndExpand(module,nModules++);
1759   }
1760   if (nModules<1) return;
1761   //
1762   int npc = 0;
1763   double *deltas = fMillepede->GetDeltaPars();
1764   double *tmpArr = new double[nModules]; 
1765   //
1766   for (int ip=0;ip<kNParCh;ip++) {
1767     if (nFree[ip]<1) continue; // nothing to do
1768     // compute the median of the deltas
1769     for (int ich=nModules;ich--;) tmpArr[ich] = deltas[modSet[ich]->GetUniqueID()*kNParCh+ip];
1770     for (int ic0=0;ic0<nModules;ic0++) // order the deltas 
1771       for (int ic1=ic0+1;ic1<nModules;ic1++) 
1772         if (tmpArr[ic0]>tmpArr[ic1])  {double tv=tmpArr[ic0]; tmpArr[ic0]=tmpArr[ic1]; tmpArr[ic1]=tv;};
1773     //
1774     int kmed = nModules/2;
1775     double median = (tmpArr[kmed]+tmpArr[nModules-kmed-1])/2.;
1776     //    
1777     for (int ich=nModules;ich--;) {
1778       AliITSAlignMille2Module* module = (AliITSAlignMille2Module*) modSet[ich];
1779       if ( !((pattern>>ip)&0x1) || !module->IsFreeDOF(ip)) continue;
1780       deltas[module->GetUniqueID()*kNParCh+ip] -= median - val;
1781     }
1782     npc++;
1783   }
1784   //
1785   delete[] nFree;
1786   delete[] tmpArr;
1787   AliInfo(Form("Applied median constraint to %d params for %d orphan modules",npc,nModules));
1788   //
1789 }
1790
1791 //_________________________________________________________________________
1792 void AliITSAlignMille2::ConstrainLinComb(const Int_t *vidLst, const Float_t *wghLst, Int_t nmd, Double_t val, UInt_t pattern)
1793 {
1794   // require that the linear combinations of the nmd modules (refered by their volume ID) from the 
1795   // modList with the coefficients wghLst adds up to val.
1796   // pattern is the bit pattern for the parameters to constrain.
1797   //
1798   static TObjArray modSet;
1799   modSet.Clear();
1800   //
1801   // build list of childs for this module
1802   int nModules = 0;
1803   int *nFree = new int[kNParCh];
1804   for (int i=0;i<kNParCh;i++) nFree[i] = 0;
1805   //
1806   for (int imd=nmd;imd--;) {
1807     UShort_t vid = (UShort_t)vidLst[imd];
1808     for (int i=fNModules;i--;) {
1809       AliITSAlignMille2Module* module = GetMilleModule(i);
1810       if (module->GetVolumeID() == vid) {
1811         modSet.AddAtAndExpand(module,nModules++);
1812         for (int ip=0;ip<kNParCh;ip++) if ( ((pattern>>ip)&0x1) && module->IsFreeDOF(ip)) nFree[ip]++;
1813         break;
1814       }
1815     }
1816   }
1817   //
1818   if (nModules != nmd) {
1819     AliInfo(Form("Error: constraint for %d modules requested but %d are found",nmd,nModules));
1820     delete[] nFree;
1821     return;
1822   }
1823   int npc = 0;
1824   for (int ip=0;ip<kNParCh;ip++) {
1825     if (nFree[ip]<1) continue; // nothing to do
1826     ResetLocalEquation();
1827     for (int ich=nModules;ich--;) {
1828       AliITSAlignMille2Module* module = (AliITSAlignMille2Module*) modSet[ich];
1829       if ( !((pattern>>ip)&0x1) || !module->IsFreeDOF(ip)) continue;
1830       fGlobalDerivatives[module->GetUniqueID()*kNParCh+ip] = wghLst[ich];
1831     }
1832     AddConstraint(fGlobalDerivatives,val);
1833     npc++;
1834   }
1835   //
1836   delete[] nFree;
1837   AliInfo(Form("Constrained %d params for linerar combination of %d  modules",npc,nModules));
1838   //
1839 }