80 cout <<
"ERROR! th1morph says first input histogram doesn't exist." << endl;
84 cout <<
"ERROR! th1morph says second input histogram doesn't exist." << endl;
91 TAxis* axis1 = hist1->GetXaxis();
92 Int_t nb1 = axis1->GetNbins();
93 TAxis* axis2 = hist2->GetXaxis();
94 Int_t nb2 = axis2->GetNbins();
96 std::set<Double_t> bedgesn_tmp;
97 for(Int_t
i = 1;
i <= nb1; ++
i){
98 bedgesn_tmp.insert(axis1->GetBinLowEdge(
i));
99 bedgesn_tmp.insert(axis1->GetBinUpEdge(
i));
101 for(Int_t
i = 1;
i <= nb2; ++
i){
102 bedgesn_tmp.insert(axis2->GetBinLowEdge(
i));
103 bedgesn_tmp.insert(axis2->GetBinUpEdge(
i));
105 Int_t nbn = bedgesn_tmp.size() - 1;
106 TArrayD bedgesn(nbn+1);
108 for (std::set<Double_t>::const_iterator bedge = bedgesn_tmp.begin();
109 bedge != bedgesn_tmp.end(); ++bedge){
110 bedgesn[idx]=(*bedge);
113 Double_t xminn = bedgesn[0];
114 Double_t xmaxn = bedgesn[nbn];
125 wt1 = 1. - (parinterp-par1)/(par2-par1);
126 wt2 = 1. + (parinterp-par2)/(par2-par1);
135 if (wt1 < 0 || wt1 > 1. || wt2 < 0. || wt2 > 1. || fabs(1-(wt1+wt2))
137 cout <<
"Warning! th1fmorph: This is an extrapolation!! Weights are "
138 << wt1 <<
" and " << wt2 <<
" (sum=" << wt1+wt2 <<
")" << endl;
140 if (idebug >= 1)
cout <<
"th1morph - Weights: " << wt1 <<
" " << wt2 << endl;
142 if (idebug >= 1)
cout <<
"New hist: " << nbn <<
" " << xminn <<
" "
148 if (hist1->GetSum() <= 0 || hist2->GetSum() <=0 ) {
149 cout <<
"Warning! th1morph detects an empty input histogram. Empty interpolated histogram returned: "
150 <<endl <<
" " << chname <<
" - " << chtitle << endl;
151 TH1_t *morphedhist = (TH1_t *)gROOT->FindObject(chname);
152 if (morphedhist)
delete morphedhist;
153 morphedhist =
new TH1_t(chname,chtitle,nbn,xminn,xmaxn);
156 if (idebug >= 1)
cout <<
"Input histogram content sums: "
157 << hist1->GetSum() <<
" " << hist2->GetSum() << endl;
170 Value_t *dist1=hist1->GetArray();
171 Value_t *dist2=hist2->GetArray();
172 Double_t *sigdis1 =
new Double_t[1+nb1];
173 Double_t *sigdis2 =
new Double_t[1+nb2];
174 Double_t *sigdisn =
new Double_t[2+nb1+nb2];
175 Double_t *xdisn =
new Double_t[2+nb1+nb2];
176 Double_t *sigdisf =
new Double_t[nbn+1];
178 for(Int_t
i=0;
i<2+nb1+nb2;
i++) xdisn[
i] = 0;
179 sigdis1[0] = 0; sigdis2[0] = 0;
181 for(Int_t
i=1;
i<nb1+1;
i++) {
182 sigdis1[
i] = dist1[
i];
184 for(Int_t
i=1;
i<nb2+1;
i++) {
185 sigdis2[
i] = dist2[
i];
189 for(Int_t
i=0;
i<nb1+1;
i++) {
190 cout <<
i <<
" dist1" << dist1[
i] << endl;
192 for(Int_t
i=0;
i<nb2+1;
i++) {
193 cout <<
i <<
" dist2" << dist2[
i] << endl;
201 for(Int_t
i=0;
i<nb1+1;
i++) {
204 if (idebug >=1)
cout <<
"Total histogram 1: " << total << endl;
205 for(Int_t
i=1;
i<nb1+1;
i++) {
206 sigdis1[
i] = sigdis1[
i]/total + sigdis1[
i-1];
210 for(Int_t
i=0;
i<nb2+1;
i++) {
213 if (idebug >=1)
cout <<
"Total histogram 22: " << total << endl;
214 for(Int_t
i=1;
i<nb2+1;
i++) {
215 sigdis2[
i] = sigdis2[
i]/total + sigdis2[
i-1];
228 while(sigdis1[ix1l-1] >= sigdis1[ix1l]) {
231 while(sigdis2[ix2l-1] >= sigdis2[ix2l]) {
242 }
while(sigdis1[ix1+1] <= sigdis1[0]);
247 }
while(sigdis2[ix2+1] <= sigdis2[0]);
250 cout <<
"First and last edge of hist1: " << ix1 <<
" " << ix1l << endl;
251 cout <<
" " << sigdis1[ix1] <<
" " << sigdis1[ix1+1] << endl;
252 cout <<
"First and last edge of hist2: " << ix2 <<
" " << ix2l << endl;
253 cout <<
" " << sigdis2[ix2] <<
" " << sigdis2[ix2+1] << endl;
260 x1 = axis1->GetBinLowEdge(ix1+1);
261 x2 = axis2->GetBinLowEdge(ix2+1);
266 cout <<
"First interpolated point: " << xdisn[nx3] <<
" "
267 << sigdisn[nx3] << endl;
268 cout <<
" " << x1 <<
" <= " << x <<
" <= "
277 cout <<
"----BEFORE while with ix1=" << ix1 <<
", ix1l=" << ix1l
278 <<
", ix2=" << ix2 <<
", ix2l=" << ix2l << endl;
284 while((ix1 < ix1l) | (ix2 < ix2l)) {
285 if (idebug >= 1 )
cout <<
"----Top of while with ix1=" << ix1
286 <<
", ix1l=" << ix1l <<
", ix2=" << ix2
287 <<
", ix2l=" << ix2l << endl;
295 if ((sigdis1[ix1+1] <= sigdis2[ix2+1] || ix2 == ix2l) && ix1 < ix1l) {
297 while(sigdis1[ix1+1] <= sigdis1[ix1] && ix1 < ix1l) {
301 }
else if (ix2 < ix2l) {
303 while(sigdis2[ix2+1] <= sigdis2[ix2] && ix2 < ix2l) {
310 cout <<
"Pair for i12type=1: " << sigdis2[ix2] <<
" "
311 << sigdis1[ix1] <<
" " << sigdis2[ix2+1] << endl;
313 x1 = axis1->GetBinLowEdge(ix1+1);
315 Double_t x20 = axis2->GetBinLowEdge(ix2+1);
316 Double_t x21 = axis2->GetBinUpEdge(ix2+1);
317 Double_t y20 = sigdis2[ix2];
318 Double_t y21 = sigdis2[ix2+1];
325 x2 = x20 + (x21-x20)*(y-y20)/(y21-y20);
332 cout <<
"Pair for i12type=2: " << sigdis1[ix1] <<
" " << sigdis2[ix2]
333 <<
" " << sigdis1[ix1+1] << endl;
335 x2 = axis2->GetBinLowEdge(ix2+1);
337 Double_t x10 = axis1->GetBinLowEdge(ix1+1);
338 Double_t x11 = axis1->GetBinUpEdge(ix1+1);
339 Double_t y10 = sigdis1[ix1];
340 Double_t y11 = sigdis1[ix1+1];
347 x1 = x10 + (x11-x10)*(y-y10)/(y11-y10);
363 cout <<
" ---> y > yprev: i12type=" << i12type <<
", nx3="
364 << nx3 <<
", x= " << x <<
", y=" << y <<
", yprev=" << yprev
371 cout <<
" ix1=" << ix1 <<
", ix2= " << ix2 <<
", i12type= "
372 << i12type <<
", sigdis1[ix1]=" << sigdis1[ix1] << endl;
373 cout <<
" " <<
", nx3=" << nx3 <<
", x=" << x <<
", y= "
374 << sigdisn[nx3] << endl;
378 if (idebug >=3)
for (Int_t
i=0;
i<nx3;
i++) {
379 cout <<
" nx " <<
i <<
" " << xdisn[
i] <<
" " << sigdisn[
i] << endl;
394 if (idebug >= 1)
cout <<
"------> Any final bins to set? " << x <<
" "
395 << xdisn[nx3] << endl;
396 while(x >= xdisn[nx3]) {
397 sigdisf[ix] = sigdisn[nx3];
398 if (idebug >= 2)
cout <<
" Setting final bins" << ix <<
" " << x
399 <<
" " << sigdisf[ix] << endl;
404 if (idebug >= 1)
cout <<
" Now ixl=" << ixl <<
" ix=" << ix << endl;
416 if (idebug >= 1)
cout <<
"Start setting initial bins at x=" << x << endl;
417 while(x <= xdisn[0]) {
418 sigdisf[ix] = sigdisn[0];
419 if (idebug >= 1)
cout <<
" Setting initial bins " << ix <<
" " << x
420 <<
" " << xdisn[1] <<
" " << sigdisf[ix] << endl;
427 cout <<
"Bins left to loop over:" << ixf <<
"-" << ixl << endl;
433 for(ix=ixf;ix<ixl;ix++) {
437 }
else if (x > xdisn[nx3]) {
440 while(xdisn[ix3+1] <= x && ix3 < 2*nbn) {
443 Double_t dx2=axis2->GetBinWidth(axis2->FindBin(x));
444 if (xdisn[ix3+1]-x > 1.1*dx2) {
447 else if (xdisn[ix3+1] > xdisn[ix3]) {
448 y = sigdisn[ix3] + (sigdisn[ix3+1]-sigdisn[ix3])
449 *(x-xdisn[ix3])/(xdisn[ix3+1]-xdisn[ix3]);
452 cout <<
"Warning - th1fmorph: This probably shoudn't happen! "
454 cout <<
"Warning - th1fmorph: Zero slope solving x(y)" << endl;
459 cout << ix <<
", ix3=" << ix3 <<
", xdisn=" << xdisn[ix3] <<
", x="
460 << x <<
", next xdisn=" << xdisn[ix3+1] << endl;
461 cout <<
" cdf n=" << sigdisn[ix3] <<
", y=" << y <<
", next point="
462 << sigdisn[ix3+1] << endl;
469 TH1_t *morphedhist = (TH1_t *)gROOT->FindObject(chname);
470 if (morphedhist)
delete morphedhist;
471 morphedhist =
new TH1_t(chname,chtitle,nbn,bedgesn.GetArray());
473 for(ix=nbn-1;ix>-1;ix--) {
474 y = sigdisf[ix+1]-sigdisf[ix];
475 morphedhist->SetBinContent(ix+1,y*morphedhistnorm);
480 delete sigdis1;
delete sigdis2;
481 delete sigdisn;
delete xdisn;
delete sigdisf;