21 extern void Yields(Int_t &, Double_t *, Double_t &f, Double_t *x, Int_t iflag);
316 TSPlot::TSPlot(Int_t nx, Int_t ny, Int_t ne, Int_t ns, TTree *tree) :
327 fXvar.ResizeTo(fNevents, fNx);
328 fYvar.ResizeTo(fNevents, fNy);
329 fYpdf.ResizeTo(fNevents, fNSpecies*fNy);
330 fSWeights.ResizeTo(fNevents, fNSpecies*(fNy+1));
332 fNumbersOfEvents = 0;
340 if (fNumbersOfEvents)
341 delete [] fNumbersOfEvents;
342 if (!fXvarHists.IsEmpty())
344 if (!fYvarHists.IsEmpty())
346 if (!fYpdfHists.IsEmpty())
353 void TSPlot::Browse(TBrowser *b)
355 if (!fSWeightsHists.IsEmpty()) {
356 TIter next(&fSWeightsHists);
358 while ((h = (TH1D*)next()))
359 b->Add(h,h->GetName());
362 if (!fYpdfHists.IsEmpty()) {
363 TIter next(&fYpdfHists);
365 while ((h = (TH1D*)next()))
366 b->Add(h,h->GetName());
368 if (!fYvarHists.IsEmpty()) {
369 TIter next(&fYvarHists);
371 while ((h = (TH1D*)next()))
372 b->Add(h,h->GetName());
374 if (!fXvarHists.IsEmpty()) {
375 TIter next(&fXvarHists);
377 while ((h = (TH1D*)next()))
378 b->Add(h,h->GetName());
380 b->Add(&fSWeights,
"sWeights");
387 void TSPlot::SetInitialNumbersOfSpecies(Int_t *numbers)
389 if (!fNumbersOfEvents)
390 fNumbersOfEvents =
new Double_t[fNSpecies];
391 for (Int_t i=0; i<fNSpecies; i++)
392 fNumbersOfEvents[i]=numbers[i];
403 void TSPlot::MakeSPlot(Option_t *option)
406 if (!fNumbersOfEvents){
407 Error(
"MakeSPlot",
"Initial numbers of events in species have not been set");
410 Int_t i, j, ispecies;
412 TString opt = option;
414 opt.ReplaceAll(
"VV",
"W");
418 if (TVirtualFitter::GetFitter()){
419 Int_t strdiff=strcmp(TVirtualFitter::GetFitter()->IsA()->GetName(), s);
421 delete TVirtualFitter::GetFitter();
425 TVirtualFitter *minuit = TVirtualFitter::Fitter(0, 2);
426 fPdfTot.ResizeTo(fNevents, fNSpecies);
430 for (Int_t iplot=-1; iplot<fNy; iplot++){
431 for (i=0; i<fNevents; i++){
432 for (ispecies=0; ispecies<fNSpecies; ispecies++){
433 fPdfTot(i, ispecies)=1;
434 for (j=0; j<fNy; j++){
436 fPdfTot(i, ispecies)*=fYpdf(i, ispecies*fNy+j);
441 minuit->SetFCN(Yields);
442 Double_t arglist[10];
444 if (opt.Contains(
"Q")||opt.Contains(
"V")){
447 if (opt.Contains(
"W"))
449 minuit->ExecuteCommand(
"SET PRINT", arglist, 1);
451 minuit->SetObjectFit(&fPdfTot);
452 for (ispecies=0; ispecies<fNSpecies; ispecies++)
453 minuit->SetParameter(ispecies,
"", fNumbersOfEvents[ispecies], 1, 0, 0);
455 minuit->ExecuteCommand(
"MIGRAD", arglist, 0);
456 for (ispecies=0; ispecies<fNSpecies; ispecies++){
457 fNumbersOfEvents[ispecies]=minuit->GetParameter(ispecies);
458 if (!opt.Contains(
"Q"))
459 printf(
"estimated #of events in species %d = %f\n", ispecies, fNumbersOfEvents[ispecies]);
461 if (!opt.Contains(
"Q"))
463 Double_t *covmat = minuit->GetCovarianceMatrix();
464 SPlots(covmat, iplot);
466 if (opt.Contains(
"W")){
467 Double_t *sumweight =
new Double_t[fNSpecies];
468 for (i=0; i<fNSpecies; i++){
470 for (j=0; j<fNevents; j++)
471 sumweight[i]+=fSWeights(j, (iplot+1)*fNSpecies + i);
472 printf(
"checking sum of weights[%d]=%f\n", i, sumweight[i]);
483 void TSPlot::SPlots(Double_t *covmat, Int_t i_excl)
485 Int_t i, ispecies, k;
486 Double_t numerator, denominator;
487 for (i=0; i<fNevents; i++){
489 for (ispecies=0; ispecies<fNSpecies; ispecies++)
490 denominator+=fNumbersOfEvents[ispecies]*fPdfTot(i, ispecies);
491 for (ispecies=0; ispecies<fNSpecies; ispecies++){
493 for (k=0; k<fNSpecies; k++)
494 numerator+=covmat[ispecies*fNSpecies+k]*fPdfTot(i, k);
495 fSWeights(i, (i_excl+1)*fNSpecies + ispecies)=numerator/denominator;
504 void TSPlot::GetSWeights(TMatrixD &weights)
506 if (weights.GetNcols()!=fNSpecies*(fNy+1) || weights.GetNrows()!=fNevents)
507 weights.ResizeTo(fNevents, fNSpecies*(fNy+1));
515 void TSPlot::GetSWeights(Double_t *weights)
517 for (Int_t i=0; i<fNevents; i++){
518 for (Int_t j=0; j<fNSpecies; j++){
519 weights[i*fNSpecies+j]=fSWeights(i, j);
527 void TSPlot::FillXvarHists(Int_t nbins)
531 if (!fXvarHists.IsEmpty()){
532 if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
540 for (i=0; i<fNx; i++){
541 snprintf(name,
sizeof(name),
"x%d", i);
542 TH1D *h =
new TH1D(name, name, nbins, fMinmax(0, i), fMinmax(1, i));
543 for (j=0; j<fNevents; j++)
544 h->Fill(fXvar(j, i));
555 TObjArray* TSPlot::GetXvarHists()
558 if (fXvarHists.IsEmpty())
559 FillXvarHists(nbins);
560 else if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
561 FillXvarHists(nbins);
570 TH1D *TSPlot::GetXvarHist(Int_t ixvar)
573 if (fXvarHists.IsEmpty())
574 FillXvarHists(nbins);
575 else if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
576 FillXvarHists(nbins);
578 return (TH1D*)fXvarHists.UncheckedAt(ixvar);
584 void TSPlot::FillYvarHists(Int_t nbins)
588 if (!fYvarHists.IsEmpty()){
589 if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
597 for (i=0; i<fNy; i++){
598 snprintf(name,
sizeof(name),
"y%d", i);
599 TH1D *h=
new TH1D(name, name, nbins, fMinmax(0, fNx+i), fMinmax(1, fNx+i));
600 for (j=0; j<fNevents; j++)
601 h->Fill(fYvar(j, i));
610 TObjArray* TSPlot::GetYvarHists()
613 if (fYvarHists.IsEmpty())
614 FillYvarHists(nbins);
615 else if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
616 FillYvarHists(nbins);
624 TH1D *TSPlot::GetYvarHist(Int_t iyvar)
627 if (fYvarHists.IsEmpty())
628 FillYvarHists(nbins);
629 else if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
630 FillYvarHists(nbins);
631 return (TH1D*)fYvarHists.UncheckedAt(iyvar);
637 void TSPlot::FillYpdfHists(Int_t nbins)
639 Int_t i, j, ispecies;
641 if (!fYpdfHists.IsEmpty()){
642 if (((TH1D*)fYpdfHists.First())->GetNbinsX()!=nbins)
649 for (ispecies=0; ispecies<fNSpecies; ispecies++){
650 for (i=0; i<fNy; i++){
651 snprintf(name,
sizeof(name),
"pdf_species%d_y%d", ispecies, i);
653 TH1D *h =
new TH1D(name, name, nbins, fMinmax(0, fNx+fNy+ispecies*fNy+i), fMinmax(1, fNx+fNy+ispecies*fNy+i));
654 for (j=0; j<fNevents; j++)
655 h->Fill(fYpdf(j, ispecies*fNy+i));
666 TObjArray* TSPlot::GetYpdfHists()
669 if (fYpdfHists.IsEmpty())
670 FillYpdfHists(nbins);
680 TH1D *TSPlot::GetYpdfHist(Int_t iyvar, Int_t ispecies)
683 if (fYpdfHists.IsEmpty())
684 FillYpdfHists(nbins);
686 return (TH1D*)fYpdfHists.UncheckedAt(fNy*ispecies+iyvar);
697 void TSPlot::FillSWeightsHists(Int_t nbins)
699 if (fSWeights.GetNoElements()==0){
700 Error(
"GetSWeightsHists",
"SWeights were not computed");
704 if (!fSWeightsHists.IsEmpty()){
705 if (((TH1D*)fSWeightsHists.First())->GetNbinsX()!=nbins)
706 fSWeightsHists.Delete();
714 for (Int_t ivar=0; ivar<fNx; ivar++){
715 for (Int_t ispecies=0; ispecies<fNSpecies; ispecies++){
716 snprintf(name,30,
"x%d_species%d", ivar, ispecies);
717 TH1D *h =
new TH1D(name, name, nbins, fMinmax(0, ivar), fMinmax(1, ivar));
719 for (Int_t ievent=0; ievent<fNevents; ievent++)
720 h->Fill(fXvar(ievent, ivar), fSWeights(ievent, ispecies));
721 fSWeightsHists.AddLast(h);
726 for (Int_t iexcl=0; iexcl<fNy; iexcl++){
727 for(Int_t ispecies=0; ispecies<fNSpecies; ispecies++){
728 snprintf(name,30,
"y%d_species%d", iexcl, ispecies);
729 TH1D *h =
new TH1D(name, name, nbins, fMinmax(0, fNx+iexcl), fMinmax(1, fNx+iexcl));
731 for (Int_t ievent=0; ievent<fNevents; ievent++)
732 h->Fill(fYvar(ievent, iexcl), fSWeights(ievent, fNSpecies*(iexcl+1)+ispecies));
733 fSWeightsHists.AddLast(h);
745 TObjArray *TSPlot::GetSWeightsHists()
748 if (fSWeightsHists.IsEmpty())
749 FillSWeightsHists(nbins);
751 return &fSWeightsHists;
766 void TSPlot::RefillHist(Int_t type, Int_t nvar, Int_t nbins, Double_t min, Double_t max, Int_t nspecies)
768 if (type<1 || type>5){
769 Error(
"RefillHist",
"type must lie between 1 and 5");
776 hremove = (TH1D*)fXvarHists.RemoveAt(nvar);
778 snprintf(name,20,
"x%d",nvar);
779 TH1D *h =
new TH1D(name, name, nbins, min, max);
780 for (j=0; j<fNevents;j++)
781 h->Fill(fXvar(j, nvar));
782 fXvarHists.AddAt(h, nvar);
785 hremove = (TH1D*)fYvarHists.RemoveAt(nvar);
787 snprintf(name,20,
"y%d", nvar);
788 TH1D *h =
new TH1D(name, name, nbins, min, max);
789 for (j=0; j<fNevents;j++)
790 h->Fill(fYvar(j, nvar));
791 fXvarHists.AddAt(h, nvar);
794 hremove = (TH1D*)fYpdfHists.RemoveAt(nspecies*fNy+nvar);
796 snprintf(name,20,
"pdf_species%d_y%d", nspecies, nvar);
797 TH1D *h=
new TH1D(name, name, nbins, min, max);
798 for (j=0; j<fNevents; j++)
799 h->Fill(fYpdf(j, nspecies*fNy+nvar));
800 fYpdfHists.AddAt(h, nspecies*fNy+nvar);
803 hremove = (TH1D*)fSWeightsHists.RemoveAt(fNSpecies*nvar+nspecies);
805 snprintf(name,20,
"x%d_species%d", nvar, nspecies);
806 TH1D *h =
new TH1D(name, name, nbins, min, max);
808 for (Int_t ievent=0; ievent<fNevents; ievent++)
809 h->Fill(fXvar(ievent, nvar), fSWeights(ievent, nspecies));
810 fSWeightsHists.AddAt(h, fNSpecies*nvar+nspecies);
813 hremove = (TH1D*)fSWeightsHists.RemoveAt(fNx*fNSpecies + fNSpecies*nvar+nspecies);
815 snprintf(name,20,
"y%d_species%d", nvar, nspecies);
816 TH1D *h =
new TH1D(name, name, nbins, min, max);
818 for (Int_t ievent=0; ievent<fNevents; ievent++)
819 h->Fill(fYvar(ievent, nvar), fSWeights(ievent, nspecies));
820 fSWeightsHists.AddAt(h, fNx*fNSpecies + fNSpecies*nvar+nspecies);
831 TH1D *TSPlot::GetSWeightsHist(Int_t ixvar, Int_t ispecies,Int_t iyexcl)
835 if (fSWeightsHists.IsEmpty())
836 FillSWeightsHists(nbins);
839 return (TH1D*)fSWeightsHists.UncheckedAt(fNx*fNSpecies + fNSpecies*iyexcl+ispecies);
841 return (TH1D*)fSWeightsHists.UncheckedAt(fNSpecies*ixvar + ispecies);
849 void TSPlot::SetTree(TTree *tree)
867 void TSPlot::SetTreeSelection(
const char* varexp,
const char *selection, Long64_t firstentry)
870 std::vector<TString> cnames;
871 TList *formulaList =
new TList();
872 TSelectorDraw *selector = (TSelectorDraw*)(((TTreePlayer*)fTree->GetPlayer())->GetSelector());
874 Long64_t entry, entryNumber;
877 TObjArray *leaves = fTree->GetListOfLeaves();
879 fTreename=
new TString(fTree->GetName());
881 fVarexp =
new TString(varexp);
883 fSelection=
new TString(selection);
885 nch = varexp ? strlen(varexp) : 0;
889 TTreeFormula *select = 0;
890 if (selection && strlen(selection)) {
891 select =
new TTreeFormula(
"Selection",selection,fTree);
893 if (!select->GetNdim()) {
delete select;
return; }
894 formulaList->Add(select);
899 ncols = fNx + fNy + fNy*fNSpecies;
900 for (i=0;i<ncols;i++) {
901 cnames.push_back( leaves->At(i)->GetName() );
905 ncols = selector->SplitNames(varexp,cnames);
907 var =
new TTreeFormula* [ncols];
908 Double_t *xvars =
new Double_t[ncols];
910 fMinmax.ResizeTo(2, ncols);
911 for (i=0; i<ncols; i++){
917 for (i=0;i<ncols;i++) {
918 var[i] =
new TTreeFormula(
"Var1",cnames[i].Data(),fTree);
919 formulaList->Add(var[i]);
923 TTreeFormulaManager *manager=0;
924 if (formulaList->LastIndex()>=0) {
925 manager =
new TTreeFormulaManager;
926 for(i=0;i<=formulaList->LastIndex();i++) {
927 manager->Add((TTreeFormula*)formulaList->At(i));
934 Long64_t selectedrows=0;
935 for (entry=firstentry;entry<firstentry+fNevents;entry++) {
936 entryNumber = fTree->GetEntryNumber(entry);
937 if (entryNumber < 0)
break;
938 Long64_t localEntry = fTree->LoadTree(entryNumber);
939 if (localEntry < 0)
break;
940 if (tnumber != fTree->GetTreeNumber()) {
941 tnumber = fTree->GetTreeNumber();
942 if (manager) manager->UpdateFormulaLeaves();
945 if (manager && manager->GetMultiplicity()) {
946 ndata = manager->GetNdata();
949 for(Int_t inst=0;inst<ndata;inst++) {
950 Bool_t loaded = kFALSE;
952 if (select->EvalInstance(inst) == 0) {
957 if (inst==0) loaded = kTRUE;
961 for (i=0;i<ncols;i++) {
962 var[i]->EvalInstance(0);
967 for (i=0;i<ncols;i++) {
968 xvars[i] = var[i]->EvalInstance(inst);
976 for (i=0; i<fNx; i++){
977 fXvar(selectedrows, i) = xvars[i];
978 if (fXvar(selectedrows, i) < fMinmax(0, i))
979 fMinmax(0, i)=fXvar(selectedrows, i);
980 if (fXvar(selectedrows, i) > fMinmax(1, i))
981 fMinmax(1, i)=fXvar(selectedrows, i);
983 for (i=0; i<fNy; i++){
984 fYvar(selectedrows, i) = xvars[i+fNx];
986 if (fYvar(selectedrows, i) < fMinmax(0, i+fNx))
987 fMinmax(0, i+fNx) = fYvar(selectedrows, i);
988 if (fYvar(selectedrows, i) > fMinmax(1, i+fNx))
989 fMinmax(1, i+fNx) = fYvar(selectedrows, i);
990 for (Int_t j=0; j<fNSpecies; j++){
991 fYpdf(selectedrows, j*fNy + i)=xvars[j*fNy + i+fNx+fNy];
992 if (fYpdf(selectedrows, j*fNy+i) < fMinmax(0, j*fNy+i+fNx+fNy))
993 fMinmax(0, j*fNy+i+fNx+fNy) = fYpdf(selectedrows, j*fNy+i);
994 if (fYpdf(selectedrows, j*fNy+i) > fMinmax(1, j*fNy+i+fNx+fNy))
995 fMinmax(1, j*fNy+i+fNx+fNy) = fYpdf(selectedrows, j*fNy+i);
1001 fNevents=selectedrows;
1018 void Yields(Int_t &, Double_t *, Double_t &f, Double_t *x, Int_t )
1023 TVirtualFitter *fitter = TVirtualFitter::GetFitter();
1024 TMatrixD *pdftot = (TMatrixD*)fitter->GetObjectFit();
1025 Int_t nev = pdftot->GetNrows();
1026 Int_t nes = pdftot->GetNcols();
1028 for (i=0; i<nev; i++){
1030 for (ispecies=0; ispecies<nes; ispecies++)
1031 lik+=x[ispecies]*(*pdftot)(i, ispecies);
1037 for (i=0; i<nes; i++)