]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMillePedeRecord.cxx
Fix for raw tag file creation.
[u/mrichter/AliRoot.git] / STEER / AliMillePedeRecord.cxx
1 #include "AliMillePedeRecord.h"
2 #include <TMath.h>
3 #include "AliLog.h"
4
5 /**********************************************************************************************/
6 /* AliMillePedeRecords: class to store the data of single track processing                    */
7 /* Format: for each measured point the data is stored consequtively                           */
8 /* INDEX                                                      VALUE                           */
9 /* -1                                                         residual                        */
10 /* Local_param_id                                             dResidual/dLocal_param          */
11 /* ...                                                        ...                             */
12 /* -2                                                         weight of the measurement       */
13 /* Global_param_od                                            dResidual/dGlobal_param         */
14 /* ...                                                        ...                             */
15 /*                                                                                            */
16 /* The records for all processed tracks are stored in the temporary tree in orgder to be      */
17 /* reused for multiple iterations of MillePede                                                */
18 /*                                                                                            */ 
19 /* Author: ruben.shahoyan@cern.ch                                                             */
20 /*                                                                                            */ 
21 /**********************************************************************************************/
22
23 ClassImp(AliMillePedeRecord)
24
25 //_____________________________________________________________________________________________
26 AliMillePedeRecord::AliMillePedeRecord() : 
27 fSize(0),fNGroups(0),fRunID(0),fGroupID(0),fIndex(0),fValue(0),fWeight(1) {SetUniqueID(0);}
28
29 //_____________________________________________________________________________________________
30 AliMillePedeRecord::AliMillePedeRecord(const AliMillePedeRecord& src) : 
31   TObject(src),fSize(src.fSize),fNGroups(src.fNGroups),fRunID(src.fRunID),fGroupID(0),fIndex(0),fValue(0),fWeight(src.fWeight)
32 {
33   fIndex = new Int_t[GetDtBufferSize()];
34   memcpy(fIndex,src.fIndex,fSize*sizeof(Int_t));
35   fValue = new Double_t[GetDtBufferSize()];
36   memcpy(fValue,src.fValue,fSize*sizeof(Double_t));
37   fGroupID = new UShort_t[GetGrBufferSize()];
38   memcpy(fGroupID,src.fGroupID,GetGrBufferSize()*sizeof(UShort_t));
39 }
40
41 //_____________________________________________________________________________________________
42 AliMillePedeRecord& AliMillePedeRecord::operator=(const AliMillePedeRecord& rhs)
43 {
44   if (this!=&rhs) {
45     Reset();
46     for (int i=0;i<rhs.GetSize();i++) {
47       Double_t val;
48       Int_t    ind;
49       rhs.GetIndexValue(i,ind,val);
50       AddIndexValue(ind,val);
51     }
52     fWeight = rhs.fWeight;
53     fRunID = rhs.fRunID;
54     for (int i=0;i<rhs.GetNGroups();i++) MarkGroup(rhs.GetGroupID(i));
55   }
56   return *this;
57 }
58
59 //_____________________________________________________________________________________________
60 AliMillePedeRecord::~AliMillePedeRecord() {delete[] fIndex; delete[] fValue; delete[] fGroupID;}
61
62 //_____________________________________________________________________________________________
63 void AliMillePedeRecord::Reset()
64 {
65   fSize = 0;
66   for (int i=fNGroups;i--;) fGroupID[i] = 0;
67   fNGroups = 0;
68   fRunID = 0;
69   fWeight = 1.;
70 }
71
72 //_____________________________________________________________________________________________
73 void AliMillePedeRecord::Print(const Option_t *) const
74 {
75   if (!fSize) {AliInfo("No data"); return;}
76   int cnt=0,point=0;
77   //  
78   if (fNGroups) printf("Groups: ");
79   for (int i=0;i<fNGroups;i++) printf("%4d |",GetGroupID(i)); 
80   printf("Run: %9d Weight: %+.2e\n",fRunID,fWeight);
81   while(cnt<fSize) {
82     //
83     Double_t resid = fValue[cnt++];
84     Double_t *derLoc = GetValue()+cnt;
85     int    *indLoc = GetIndex()+cnt;
86     int     nLoc = 0;
87     while(!IsWeight(cnt)) {nLoc++;cnt++;}
88     Double_t weight = GetValue(cnt++);
89     Double_t *derGlo = GetValue()+cnt;
90     int    *indGlo = GetIndex()+cnt;
91     int     nGlo = 0;
92     while(!IsResidual(cnt) && cnt<fSize) {nGlo++; cnt++;} 
93     //
94     printf("\n*** Point#%2d | Residual = %+.4e | Weight = %+.4e\n",point++,resid,weight);
95     printf("Locals : "); 
96     for (int i=0;i<nLoc;i++) printf("[%5d] %+.4e|",indLoc[i],derLoc[i]); printf("\n");
97     printf("Globals: "); 
98     for (int i=0;i<nGlo;i++) printf("[%5d] %+.4e|",indGlo[i],derGlo[i]); printf("\n");
99     //
100   }
101   //
102 }
103
104 //_____________________________________________________________________________________________
105 void AliMillePedeRecord::ExpandDtBuffer(Int_t bfsize)
106 {
107   // add extra space for derivatives data
108   bfsize = TMath::Max(bfsize,GetDtBufferSize());
109   Int_t *tmpI = new Int_t[bfsize];
110   memcpy(tmpI,fIndex, fSize*sizeof(Int_t));
111   delete [] fIndex;
112   fIndex = tmpI;
113   //
114   Double_t *tmpD = new Double_t[bfsize];
115   memcpy(tmpD,fValue, fSize*sizeof(Double_t));
116   delete [] fValue;
117   fValue = tmpD;
118   //
119   SetDtBufferSize(bfsize);
120 }
121
122 //_____________________________________________________________________________________________
123 void AliMillePedeRecord::ExpandGrBuffer(Int_t bfsize)
124 {
125   // add extra space for groupID data 
126   bfsize = TMath::Max(bfsize,GetGrBufferSize());
127   UShort_t *tmpI = new UShort_t[bfsize];
128   memcpy(tmpI,fGroupID, fNGroups*sizeof(UShort_t));
129   delete [] fGroupID;
130   fGroupID = tmpI;
131   for (int i=fNGroups;i<bfsize;i++) fGroupID[i] = 0;
132   //
133   SetGrBufferSize(bfsize);
134 }
135
136 //_____________________________________________________________________________________________
137 void AliMillePedeRecord::MarkGroup(Int_t id)
138 {
139   // mark the presence of the detector group
140   id++; // groupID is stored as realID+1
141   if (fNGroups>0 && fGroupID[fNGroups-1]==id) return; // already there
142   if (fNGroups>=GetGrBufferSize()) ExpandGrBuffer(2*(fNGroups+1));
143   fGroupID[fNGroups++] = id;  
144 }
145