]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMillePedeRecord.cxx
in addition to 41626: liReconstruction will store in the AliESDRun the average intens...
[u/mrichter/AliRoot.git] / STEER / AliMillePedeRecord.cxx
CommitLineData
de34b538 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
23ClassImp(AliMillePedeRecord)
24
25//_____________________________________________________________________________________________
26AliMillePedeRecord::AliMillePedeRecord() :
bce37084 27fSize(0),fNGroups(0),fGroupID(0),fIndex(0),fValue(0),fWeight(1) {SetUniqueID(0);}
de34b538 28
29//_____________________________________________________________________________________________
30AliMillePedeRecord::AliMillePedeRecord(const AliMillePedeRecord& src) :
a8c5b94c 31 TObject(src),fSize(src.fSize),fNGroups(src.fNGroups),fGroupID(0),fIndex(0),fValue(0),fWeight(src.fWeight)
de34b538 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));
a8c5b94c 37 fGroupID = new UShort_t[GetGrBufferSize()];
38 memcpy(fGroupID,src.fGroupID,GetGrBufferSize()*sizeof(UShort_t));
de34b538 39}
40
41//_____________________________________________________________________________________________
42AliMillePedeRecord& 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 }
a8c5b94c 52 fWeight = rhs.fWeight;
de34b538 53 for (int i=0;i<rhs.GetNGroups();i++) MarkGroup(rhs.GetGroupID(i));
54 }
55 return *this;
56}
57
58//_____________________________________________________________________________________________
59AliMillePedeRecord::~AliMillePedeRecord() {delete[] fIndex; delete[] fValue; delete[] fGroupID;}
60
61//_____________________________________________________________________________________________
62void AliMillePedeRecord::Reset()
63{
64 fSize = 0;
a8c5b94c 65 for (int i=fNGroups;i--;) fGroupID[i] = 0;
de34b538 66 fNGroups = 0;
a8c5b94c 67 fWeight = 1.;
de34b538 68}
69
70//_____________________________________________________________________________________________
71void AliMillePedeRecord::Print(const Option_t *) const
72{
73 if (!fSize) {AliInfo("No data"); return;}
74 int cnt=0,point=0;
75 //
76 if (fNGroups) printf("Groups: ");
a8c5b94c 77 for (int i=0;i<fNGroups;i++) printf("%4d |",GetGroupID(i));
78 printf(" Weight: %+.2e\n",fWeight);
de34b538 79 while(cnt<fSize) {
80 //
81 Double_t resid = fValue[cnt++];
82 Double_t *derLoc = GetValue()+cnt;
83 int *indLoc = GetIndex()+cnt;
84 int nLoc = 0;
85 while(!IsWeight(cnt)) {nLoc++;cnt++;}
86 Double_t weight = GetValue(cnt++);
87 Double_t *derGlo = GetValue()+cnt;
88 int *indGlo = GetIndex()+cnt;
89 int nGlo = 0;
90 while(!IsResidual(cnt) && cnt<fSize) {nGlo++; cnt++;}
91 //
92 printf("\n*** Point#%2d | Residual = %+.4e | Weight = %+.4e\n",point++,resid,weight);
93 printf("Locals : ");
94 for (int i=0;i<nLoc;i++) printf("[%5d] %+.4e|",indLoc[i],derLoc[i]); printf("\n");
95 printf("Globals: ");
96 for (int i=0;i<nGlo;i++) printf("[%5d] %+.4e|",indGlo[i],derGlo[i]); printf("\n");
97 //
98 }
99 //
100}
101
102//_____________________________________________________________________________________________
103void AliMillePedeRecord::ExpandDtBuffer(Int_t bfsize)
104{
105 // add extra space for derivatives data
106 bfsize = TMath::Max(bfsize,GetDtBufferSize());
107 Int_t *tmpI = new Int_t[bfsize];
108 memcpy(tmpI,fIndex, fSize*sizeof(Int_t));
109 delete fIndex;
110 fIndex = tmpI;
111 //
112 Double_t *tmpD = new Double_t[bfsize];
113 memcpy(tmpD,fValue, fSize*sizeof(Double_t));
114 delete fValue;
115 fValue = tmpD;
116 //
117 SetDtBufferSize(bfsize);
118}
119
120//_____________________________________________________________________________________________
121void AliMillePedeRecord::ExpandGrBuffer(Int_t bfsize)
122{
123 // add extra space for groupID data
124 bfsize = TMath::Max(bfsize,GetGrBufferSize());
a8c5b94c 125 UShort_t *tmpI = new UShort_t[bfsize];
126 memcpy(tmpI,fGroupID, fNGroups*sizeof(UShort_t));
127 delete[] fGroupID;
de34b538 128 fGroupID = tmpI;
a8c5b94c 129 for (int i=fNGroups;i<bfsize;i++) fGroupID[i] = 0;
de34b538 130 //
131 SetGrBufferSize(bfsize);
132}
133
134//_____________________________________________________________________________________________
135void AliMillePedeRecord::MarkGroup(Int_t id)
136{
137 // mark the presence of the detector group
a8c5b94c 138 id++; // groupID is stored as realID+1
de34b538 139 if (fNGroups>0 && fGroupID[fNGroups-1]==id) return; // already there
140 if (fNGroups>=GetGrBufferSize()) ExpandGrBuffer(2*(fNGroups+1));
141 fGroupID[fNGroups++] = id;
142}
143