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