00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 #include "TTreeIndex.h"
00019 #include "TTree.h"
00020 #include "TMath.h"
00021
00022 ClassImp(TTreeIndex)
00023
00024
00025 TTreeIndex::TTreeIndex(): TVirtualIndex()
00026 {
00027
00028
00029 fTree = 0;
00030 fN = 0;
00031 fIndexValues = 0;
00032 fIndex = 0;
00033 fMajorFormula = 0;
00034 fMinorFormula = 0;
00035 fMajorFormulaParent = 0;
00036 fMinorFormulaParent = 0;
00037 }
00038
00039
00040 TTreeIndex::TTreeIndex(const TTree *T, const char *majorname, const char *minorname)
00041 : TVirtualIndex()
00042 {
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 fTree = (TTree*)T;
00101 fN = 0;
00102 fIndexValues = 0;
00103 fIndex = 0;
00104 fMajorFormula = 0;
00105 fMinorFormula = 0;
00106 fMajorFormulaParent = 0;
00107 fMinorFormulaParent = 0;
00108 fMajorName = majorname;
00109 fMinorName = minorname;
00110 if (!T) return;
00111 fN = T->GetEntries();
00112 if (fN <= 0) {
00113 MakeZombie();
00114 Error("TreeIndex","Cannot build a TreeIndex with a Tree having no entries");
00115 return;
00116 }
00117
00118 GetMajorFormula();
00119 GetMinorFormula();
00120 if (!fMajorFormula || !fMinorFormula) {
00121 MakeZombie();
00122 Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
00123 return;
00124 }
00125 if ((fMajorFormula->GetNdim() != 1) || (fMinorFormula->GetNdim() != 1)) {
00126 MakeZombie();
00127 Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
00128 return;
00129 }
00130
00131
00132
00133
00134
00135
00136
00137 Long64_t *w = new Long64_t[fN];
00138 Long64_t i;
00139 Long64_t oldEntry = fTree->GetReadEntry();
00140 Int_t current = -1;
00141 for (i=0;i<fN;i++) {
00142 Long64_t centry = fTree->LoadTree(i);
00143 if (centry < 0) break;
00144 if (fTree->GetTreeNumber() != current) {
00145 current = fTree->GetTreeNumber();
00146 fMajorFormula->UpdateFormulaLeaves();
00147 fMinorFormula->UpdateFormulaLeaves();
00148 }
00149 Double_t majord = fMajorFormula->EvalInstance();
00150 Double_t minord = fMinorFormula->EvalInstance();
00151 Long64_t majorv = (Long64_t)majord;
00152 Long64_t minorv = (Long64_t)minord;
00153 w[i] = majorv<<31;
00154 w[i] += minorv;
00155 }
00156 fIndex = new Long64_t[fN];
00157 TMath::Sort(fN,w,fIndex,0);
00158 fIndexValues = new Long64_t[fN];
00159 for (i=0;i<fN;i++) {
00160 fIndexValues[i] = w[fIndex[i]];
00161 }
00162
00163 delete [] w;
00164 fTree->LoadTree(oldEntry);
00165 }
00166
00167
00168 TTreeIndex::~TTreeIndex()
00169 {
00170
00171
00172 if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(0);
00173 delete [] fIndexValues; fIndexValues = 0;
00174 delete [] fIndex; fIndex = 0;
00175 delete fMajorFormula; fMajorFormula = 0;
00176 delete fMinorFormula; fMinorFormula = 0;
00177 delete fMajorFormulaParent; fMajorFormulaParent = 0;
00178 delete fMinorFormulaParent; fMinorFormulaParent = 0;
00179 }
00180
00181
00182 void TTreeIndex::Append(const TVirtualIndex *add, Bool_t delaySort )
00183 {
00184
00185
00186
00187
00188
00189 if (add && add->GetN()) {
00190
00191
00192 const TTreeIndex *ti_add = dynamic_cast<const TTreeIndex*>(add);
00193 if (ti_add == 0) {
00194 Error("Append","Can only Append a TTreeIndex to a TTreeIndex but got a %s",
00195 add->IsA()->GetName());
00196 }
00197
00198 Long64_t oldn = fN;
00199 fN += add->GetN();
00200
00201 Long64_t *oldIndex = fIndex;
00202 Long64_t *oldValues = fIndexValues;
00203
00204 fIndex = new Long64_t[fN];
00205 fIndexValues = new Long64_t[fN];
00206
00207
00208
00209 Long_t size = sizeof(Long64_t) * oldn;
00210 Long_t add_size = sizeof(Long64_t) * add->GetN();
00211
00212 memcpy(fIndex,oldIndex, size);
00213 memcpy(fIndexValues,oldValues, size);
00214
00215 Long64_t *addIndex = ti_add->GetIndex();
00216 Long64_t *addValues = ti_add->GetIndexValues();
00217
00218 memcpy(fIndex + oldn, addIndex, add_size);
00219 memcpy(fIndexValues + oldn, addValues, add_size);
00220 for(Int_t i = 0; i < add->GetN(); i++) {
00221 fIndex[oldn + i] += oldn;
00222 }
00223
00224 delete [] oldIndex;
00225 delete [] oldValues;
00226 }
00227
00228
00229 if (!delaySort) {
00230 Long64_t *w = fIndexValues;
00231 Long64_t *ind = fIndex;
00232 Long64_t *conv = new Long64_t[fN];
00233
00234 TMath::Sort(fN,w,conv,0);
00235
00236 fIndex = new Long64_t[fN];
00237 fIndexValues = new Long64_t[fN];
00238
00239 for (Int_t i=0;i<fN;i++) {
00240 fIndex[i] = ind[conv[i]];
00241 fIndexValues[i] = w[conv[i]];
00242 }
00243 delete [] w;
00244 delete [] ind;
00245 delete [] conv;
00246 }
00247 }
00248
00249
00250 Int_t TTreeIndex::GetEntryNumberFriend(const TTree *T)
00251 {
00252
00253
00254
00255
00256
00257
00258 if (!T) return -3;
00259 GetMajorFormulaParent(T);
00260 GetMinorFormulaParent(T);
00261 if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
00262 if (!fMajorFormulaParent->GetNdim() || !fMinorFormulaParent->GetNdim()) {
00263
00264
00265
00266 Int_t pentry = T->GetReadEntry();
00267 if (pentry >= (Int_t)fTree->GetEntries()) return -2;
00268
00269
00270 return pentry;
00271 }
00272
00273
00274
00275 Double_t majord = fMajorFormulaParent->EvalInstance();
00276 Double_t minord = fMinorFormulaParent->EvalInstance();
00277 Long64_t majorv = (Long64_t)majord;
00278 Long64_t minorv = (Long64_t)minord;
00279
00280
00281
00282 return fTree->GetEntryNumberWithIndex(majorv,minorv);
00283 }
00284
00285
00286 Long64_t TTreeIndex::GetEntryNumberWithBestIndex(Int_t major, Int_t minor) const
00287 {
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303 if (fN == 0) return -1;
00304 Long64_t value = Long64_t(major)<<31;
00305 value += minor;
00306 Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
00307 if (i < 0) return -1;
00308 return fIndex[i];
00309 }
00310
00311
00312
00313 Long64_t TTreeIndex::GetEntryNumberWithIndex(Int_t major, Int_t minor) const
00314 {
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326 if (fN == 0) return -1;
00327 Long64_t value = Long64_t(major)<<31;
00328 value += minor;
00329 Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
00330 if (i < 0) return -1;
00331 if (fIndexValues[i] != value) return -1;
00332 return fIndex[i];
00333 }
00334
00335
00336 TTreeFormula *TTreeIndex::GetMajorFormula()
00337 {
00338
00339
00340 if (!fMajorFormula) {
00341 fMajorFormula = new TTreeFormula("Major",fMajorName.Data(),fTree);
00342 fMajorFormula->SetQuickLoad(kTRUE);
00343 }
00344 return fMajorFormula;
00345 }
00346
00347
00348 TTreeFormula *TTreeIndex::GetMinorFormula()
00349 {
00350
00351
00352 if (!fMinorFormula) {
00353 fMinorFormula = new TTreeFormula("Minor",fMinorName.Data(),fTree);
00354 fMinorFormula->SetQuickLoad(kTRUE);
00355 }
00356 return fMinorFormula;
00357 }
00358
00359
00360 TTreeFormula *TTreeIndex::GetMajorFormulaParent(const TTree *T)
00361 {
00362
00363
00364 if (!fMajorFormulaParent) {
00365 fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),(TTree*)T);
00366 fMajorFormulaParent->SetQuickLoad(kTRUE);
00367 }
00368 if (fMajorFormulaParent->GetTree() != T) {
00369 fMajorFormulaParent->SetTree((TTree*)T);
00370 fMajorFormulaParent->UpdateFormulaLeaves();
00371 }
00372 return fMajorFormulaParent;
00373 }
00374
00375
00376 TTreeFormula *TTreeIndex::GetMinorFormulaParent(const TTree *T)
00377 {
00378
00379
00380 if (!fMinorFormulaParent) {
00381 fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),(TTree*)T);
00382 fMinorFormulaParent->SetQuickLoad(kTRUE);
00383 }
00384 if (fMinorFormulaParent->GetTree() != T) {
00385 fMinorFormulaParent->SetTree((TTree*)T);
00386 fMinorFormulaParent->UpdateFormulaLeaves();
00387 }
00388
00389 return fMinorFormulaParent;
00390 }
00391
00392
00393 void TTreeIndex::Print(Option_t * option) const
00394 {
00395
00396
00397
00398
00399
00400 TString opt = option;
00401 Bool_t printEntry = kFALSE;
00402 Long64_t n = fN;
00403 if (opt.Contains("10")) n = 10;
00404 if (opt.Contains("100")) n = 100;
00405 if (opt.Contains("1000")) n = 1000;
00406 if (opt.Contains("all")) {
00407 printEntry = kTRUE;
00408 }
00409
00410 if (printEntry) {
00411
00412 Printf("\n*****************************************************************");
00413 Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
00414 Printf("*****************************************************************");
00415 Printf("%8s : %16s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data(),"entry number");
00416 Printf("*****************************************************************");
00417 for (Long64_t i=0;i<n;i++) {
00418 Long64_t minor = fIndexValues[i] & 0xffff;
00419 Long64_t major = fIndexValues[i]>>31;
00420 Printf("%8lld : %8lld : %8lld : %8lld",i,major,minor,fIndex[i]);
00421 }
00422
00423 } else {
00424
00425 Printf("\n**********************************************");
00426 Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
00427 Printf("**********************************************");
00428 Printf("%8s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data());
00429 Printf("**********************************************");
00430 for (Long64_t i=0;i<n;i++) {
00431 Long64_t minor = fIndexValues[i] & 0xffff;
00432 Long64_t major = fIndexValues[i]>>31;
00433 Printf("%8lld : %8lld : %8lld",i,major,minor);
00434 }
00435 }
00436 }
00437
00438
00439 void TTreeIndex::Streamer(TBuffer &R__b)
00440 {
00441
00442
00443
00444
00445 UInt_t R__s, R__c;
00446 if (R__b.IsReading()) {
00447 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
00448 TVirtualIndex::Streamer(R__b);
00449 fMajorName.Streamer(R__b);
00450 fMinorName.Streamer(R__b);
00451 R__b >> fN;
00452 fIndexValues = new Long64_t[fN];
00453 R__b.ReadFastArray(fIndexValues,fN);
00454 fIndex = new Long64_t[fN];
00455 R__b.ReadFastArray(fIndex,fN);
00456 R__b.CheckByteCount(R__s, R__c, TTreeIndex::IsA());
00457 } else {
00458 R__c = R__b.WriteVersion(TTreeIndex::IsA(), kTRUE);
00459 TVirtualIndex::Streamer(R__b);
00460 fMajorName.Streamer(R__b);
00461 fMinorName.Streamer(R__b);
00462 R__b << fN;
00463 R__b.WriteFastArray(fIndexValues, fN);
00464 R__b.WriteFastArray(fIndex, fN);
00465 R__b.SetByteCount(R__c, kTRUE);
00466 }
00467 }
00468
00469
00470 void TTreeIndex::UpdateFormulaLeaves(const TTree *parent)
00471 {
00472
00473
00474 if (fMajorFormula) { fMajorFormula->UpdateFormulaLeaves();}
00475 if (fMinorFormula) { fMinorFormula->UpdateFormulaLeaves();}
00476 if (fMajorFormulaParent) {
00477 if (parent) fMajorFormulaParent->SetTree((TTree*)parent);
00478 fMajorFormulaParent->UpdateFormulaLeaves();
00479 }
00480 if (fMinorFormulaParent) {
00481 if (parent) fMinorFormulaParent->SetTree((TTree*)parent);
00482 fMinorFormulaParent->UpdateFormulaLeaves();
00483 }
00484 }
00485
00486 void TTreeIndex::SetTree(const TTree *T)
00487 {
00488
00489
00490
00491
00492
00493 fTree = (TTree*)T;
00494 }
00495