45#define PEAK_WINDOW 1024
72 :
TNamed(
"Spectrum",
"Miroslav Morhac peak finder")
156 if (
h == 0)
return 0;
157 Int_t dimension =
h->GetDimension();
159 Error(
"Search",
"Only implemented for 1-d histograms");
185 Int_t last =
h->GetXaxis()->GetLast();
189 for (i = 0; i < size; i++) source[i] =
h->GetBinContent(i +
first);
192 Background(source,size,numberIterations, direction, filterOrder,smoothing,
193 smoothWindow,compton);
197 Int_t nch = strlen(
h->GetName());
198 char *
hbname =
new char[nch+20];
222 printf(
"\nNumber of positions = %d\n",
fNPeaks);
270 if (hin == 0)
return 0;
273 Error(
"Search",
"Only implemented for 1-d and 2-d histograms");
276 if (threshold <=0 || threshold >= 1) {
277 Warning(
"Search",
"threshold must 0<threshold<1, threshold=0.05 assumed");
297 if (dimension == 1) {
301 Int_t i, bin, npeaks;
313 for (i = 0; i < npeaks; i++) {
323 if (!npeaks)
return 0;
514 int numberIterations,
515 int direction,
int filterOrder,
516 bool smoothing,
int smoothWindow,
519 int i, j, w, bw, b1, b2, priz;
520 Double_t a,
b,
c,
d,
e, yb1, yb2, ai, av, men, b4, c4, d4, e4, b6, c6, d6, e6, f6, g6, b8, c8, d8, e8, f8, g8, h8, i8;
522 return "Wrong Parameters";
523 if (numberIterations < 1)
524 return "Width of Clipping Window Must Be Positive";
525 if (ssize < 2 * numberIterations + 1)
526 return "Too Large Clipping Window";
528 return "Incorrect width of smoothing window";
530 for (i = 0; i < ssize; i++){
531 working_space[i] = spectrum[i];
532 working_space[i + ssize] = spectrum[i];
534 bw=(smoothWindow-1)/2;
538 i = numberIterations;
541 for (j = i; j < ssize - i; j++) {
543 a = working_space[ssize + j];
544 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
547 working_space[j] =
a;
550 else if (smoothing ==
kTRUE){
551 a = working_space[ssize + j];
554 for (w = j - bw; w <= j + bw; w++){
555 if ( w >= 0 && w < ssize){
556 av += working_space[ssize + w];
563 for (w = j - i - bw; w <= j - i + bw; w++){
564 if ( w >= 0 && w < ssize){
565 b += working_space[ssize + w];
572 for (w = j + i - bw; w <= j + i + bw; w++){
573 if ( w >= 0 && w < ssize){
574 c += working_space[ssize + w];
585 for (j = i; j < ssize - i; j++)
586 working_space[ssize + j] = working_space[j];
596 for (j = i; j < ssize - i; j++) {
598 a = working_space[ssize + j];
599 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
602 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
603 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
604 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
605 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
610 working_space[j] =
a;
613 else if (smoothing ==
kTRUE){
614 a = working_space[ssize + j];
617 for (w = j - bw; w <= j + bw; w++){
618 if ( w >= 0 && w < ssize){
619 av += working_space[ssize + w];
626 for (w = j - i - bw; w <= j - i + bw; w++){
627 if ( w >= 0 && w < ssize){
628 b += working_space[ssize + w];
635 for (w = j + i - bw; w <= j + i + bw; w++){
636 if ( w >= 0 && w < ssize){
637 c += working_space[ssize + w];
645 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
646 if (w >= 0 && w < ssize){
647 b4 += working_space[ssize + w];
653 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
654 if (w >= 0 && w < ssize){
655 c4 += working_space[ssize + w];
661 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
662 if (w >= 0 && w < ssize){
663 d4 += working_space[ssize + w];
669 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
670 if (w >= 0 && w < ssize){
671 e4 += working_space[ssize + w];
676 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
684 for (j = i; j < ssize - i; j++)
685 working_space[ssize + j] = working_space[j];
695 for (j = i; j < ssize - i; j++) {
697 a = working_space[ssize + j];
698 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
701 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
702 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
703 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
704 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
707 d += working_space[ssize + j - (
Int_t) (3 * ai)] / 20;
708 d -= 6 * working_space[ssize + j - (
Int_t) (2 * ai)] / 20;
709 d += 15 * working_space[ssize + j - (
Int_t) ai] / 20;
710 d += 15 * working_space[ssize + j + (
Int_t) ai] / 20;
711 d -= 6 * working_space[ssize + j + (
Int_t) (2 * ai)] / 20;
712 d += working_space[ssize + j + (
Int_t) (3 * ai)] / 20;
719 working_space[j] =
a;
722 else if (smoothing ==
kTRUE){
723 a = working_space[ssize + j];
726 for (w = j - bw; w <= j + bw; w++){
727 if ( w >= 0 && w < ssize){
728 av += working_space[ssize + w];
735 for (w = j - i - bw; w <= j - i + bw; w++){
736 if ( w >= 0 && w < ssize){
737 b += working_space[ssize + w];
744 for (w = j + i - bw; w <= j + i + bw; w++){
745 if ( w >= 0 && w < ssize){
746 c += working_space[ssize + w];
754 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
755 if (w >= 0 && w < ssize){
756 b4 += working_space[ssize + w];
762 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
763 if (w >= 0 && w < ssize){
764 c4 += working_space[ssize + w];
770 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
771 if (w >= 0 && w < ssize){
772 d4 += working_space[ssize + w];
778 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
779 if (w >= 0 && w < ssize){
780 e4 += working_space[ssize + w];
785 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
788 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
789 if (w >= 0 && w < ssize){
790 b6 += working_space[ssize + w];
796 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
797 if (w >= 0 && w < ssize){
798 c6 += working_space[ssize + w];
804 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
805 if (w >= 0 && w < ssize){
806 d6 += working_space[ssize + w];
812 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
813 if (w >= 0 && w < ssize){
814 e6 += working_space[ssize + w];
820 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
821 if (w >= 0 && w < ssize){
822 f6 += working_space[ssize + w];
828 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
829 if (w >= 0 && w < ssize){
830 g6 += working_space[ssize + w];
835 b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
845 for (j = i; j < ssize - i; j++)
846 working_space[ssize + j] = working_space[j];
856 for (j = i; j < ssize - i; j++) {
858 a = working_space[ssize + j];
859 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
862 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
863 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
864 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
865 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
868 d += working_space[ssize + j - (
Int_t) (3 * ai)] / 20;
869 d -= 6 * working_space[ssize + j - (
Int_t) (2 * ai)] / 20;
870 d += 15 * working_space[ssize + j - (
Int_t) ai] / 20;
871 d += 15 * working_space[ssize + j + (
Int_t) ai] / 20;
872 d -= 6 * working_space[ssize + j + (
Int_t) (2 * ai)] / 20;
873 d += working_space[ssize + j + (
Int_t) (3 * ai)] / 20;
876 e -= working_space[ssize + j - (
Int_t) (4 * ai)] / 70;
877 e += 8 * working_space[ssize + j - (
Int_t) (3 * ai)] / 70;
878 e -= 28 * working_space[ssize + j - (
Int_t) (2 * ai)] / 70;
879 e += 56 * working_space[ssize + j - (
Int_t) ai] / 70;
880 e += 56 * working_space[ssize + j + (
Int_t) ai] / 70;
881 e -= 28 * working_space[ssize + j + (
Int_t) (2 * ai)] / 70;
882 e += 8 * working_space[ssize + j + (
Int_t) (3 * ai)] / 70;
883 e -= working_space[ssize + j + (
Int_t) (4 * ai)] / 70;
892 working_space[j] =
a;
895 else if (smoothing ==
kTRUE){
896 a = working_space[ssize + j];
899 for (w = j - bw; w <= j + bw; w++){
900 if ( w >= 0 && w < ssize){
901 av += working_space[ssize + w];
908 for (w = j - i - bw; w <= j - i + bw; w++){
909 if ( w >= 0 && w < ssize){
910 b += working_space[ssize + w];
917 for (w = j + i - bw; w <= j + i + bw; w++){
918 if ( w >= 0 && w < ssize){
919 c += working_space[ssize + w];
927 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
928 if (w >= 0 && w < ssize){
929 b4 += working_space[ssize + w];
935 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
936 if (w >= 0 && w < ssize){
937 c4 += working_space[ssize + w];
943 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
944 if (w >= 0 && w < ssize){
945 d4 += working_space[ssize + w];
951 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
952 if (w >= 0 && w < ssize){
953 e4 += working_space[ssize + w];
958 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
961 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
962 if (w >= 0 && w < ssize){
963 b6 += working_space[ssize + w];
969 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
970 if (w >= 0 && w < ssize){
971 c6 += working_space[ssize + w];
977 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
978 if (w >= 0 && w < ssize){
979 d6 += working_space[ssize + w];
985 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
986 if (w >= 0 && w < ssize){
987 e6 += working_space[ssize + w];
993 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
994 if (w >= 0 && w < ssize){
995 f6 += working_space[ssize + w];
1001 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
1002 if (w >= 0 && w < ssize){
1003 g6 += working_space[ssize + w];
1008 b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
1011 for (w = j - (
Int_t)(4 * ai) - bw; w <= j - (
Int_t)(4 * ai) + bw; w++){
1012 if (w >= 0 && w < ssize){
1013 b8 += working_space[ssize + w];
1019 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
1020 if (w >= 0 && w < ssize){
1021 c8 += working_space[ssize + w];
1027 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
1028 if (w >= 0 && w < ssize){
1029 d8 += working_space[ssize + w];
1035 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
1036 if (w >= 0 && w < ssize){
1037 e8 += working_space[ssize + w];
1043 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
1044 if (w >= 0 && w < ssize){
1045 f8 += working_space[ssize + w];
1051 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
1052 if (w >= 0 && w < ssize){
1053 g8 += working_space[ssize + w];
1059 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
1060 if (w >= 0 && w < ssize){
1061 h8 += working_space[ssize + w];
1067 for (w = j + (
Int_t)(4 * ai) - bw; w <= j + (
Int_t)(4 * ai) + bw; w++){
1068 if (w >= 0 && w < ssize){
1069 i8 += working_space[ssize + w];
1074 b8 = ( -b8 + 8 * c8 - 28 * d8 + 56 * e8 - 56 * f8 - 28 * g8 + 8 * h8 - i8)/70;
1083 working_space[j]=av;
1086 for (j = i; j < ssize - i; j++)
1087 working_space[ssize + j] = working_space[j];
1095 if (compton ==
kTRUE) {
1096 for (i = 0, b2 = 0; i < ssize; i++){
1098 a = working_space[i],
b = spectrum[i];
1104 yb1 = working_space[b1];
1105 for (b2 = b1 + 1,
c = 0, priz = 0; priz == 0 && b2 < ssize; b2++){
1106 a = working_space[b2],
b = spectrum[b2];
1115 yb2 = working_space[b2];
1117 for (j = b1,
c = 0; j <= b2; j++){
1122 c = (yb2 - yb1) /
c;
1123 for (j = b1,
d = 0; j <= b2 && j < ssize; j++){
1127 working_space[ssize + j] =
a;
1133 for (j = b2,
c = 0; j >= b1; j--){
1138 c = (yb1 - yb2) /
c;
1139 for (j = b2,
d = 0;j >= b1 && j >= 0; j--){
1143 working_space[ssize + j] =
a;
1152 for (j = 0; j < ssize; j++)
1153 spectrum[j] = working_space[ssize + j];
1154 delete[]working_space;
1200 Double_t nom, nip, nim, sp, sm, area = 0;
1202 return "Averaging Window must be positive";
1205 for(i = 0, maxch = 0; i < ssize; i++){
1207 if(maxch < source[i])
1213 delete [] working_space;
1218 working_space[
xmin] = 1;
1220 nip = source[i] / maxch;
1221 nim = source[i + 1] / maxch;
1223 for(
l = 1;
l <= averWindow;
l++){
1225 a = source[
xmax] / maxch;
1228 a = source[i +
l] / maxch;
1238 if((i -
l + 1) <
xmin)
1239 a = source[
xmin] / maxch;
1242 a = source[i -
l + 1] / maxch;
1253 a = working_space[i + 1] = working_space[i] *
a;
1257 working_space[i] = working_space[i] / nom;
1259 for(i = 0; i < ssize; i++)
1260 source[i] = working_space[i] * area;
1261 delete[]working_space;
1461 int ssize,
int numberIterations,
1462 int numberRepetitions,
Double_t boost )
1465 return "Wrong Parameters";
1467 if (numberRepetitions <= 0)
1468 return "Wrong Parameters";
1473 int i, j, k, lindex, posit, lh_gold,
l, repet;
1474 Double_t lda, ldb, ldc, area, maximum;
1481 for (i = 0; i < ssize; i++) {
1485 working_space[i] = lda;
1487 if (lda > maximum) {
1492 if (lh_gold == -1) {
1493 delete [] working_space;
1494 return "ZERO RESPONSE VECTOR";
1498 for (i = 0; i < ssize; i++)
1499 working_space[2 * ssize + i] = source[i];
1502 for (i = 0; i < ssize; i++){
1504 for (j = 0; j < ssize; j++){
1505 ldb = working_space[j];
1508 ldc = working_space[k];
1509 lda = lda + ldb * ldc;
1512 working_space[ssize + i] = lda;
1514 for (k = 0; k < ssize; k++){
1517 ldb = working_space[
l];
1518 ldc = working_space[2 * ssize + k];
1519 lda = lda + ldb * ldc;
1522 working_space[3 * ssize + i]=lda;
1526 for (i = 0; i < ssize; i++){
1527 working_space[2 * ssize + i] = working_space[3 * ssize + i];
1531 for (i = 0; i < ssize; i++)
1532 working_space[i] = 1;
1535 for (repet = 0; repet < numberRepetitions; repet++) {
1537 for (i = 0; i < ssize; i++)
1538 working_space[i] =
TMath::Power(working_space[i], boost);
1540 for (lindex = 0; lindex < numberIterations; lindex++) {
1541 for (i = 0; i < ssize; i++) {
1542 if (working_space[2 * ssize + i] > 0.000001
1543 && working_space[i] > 0.000001) {
1545 for (j = 0; j < lh_gold; j++) {
1546 ldb = working_space[j + ssize];
1551 ldc = working_space[k];
1554 ldc += working_space[k];
1558 ldc = working_space[i];
1559 lda = lda + ldb * ldc;
1561 ldb = working_space[2 * ssize + i];
1567 ldb = working_space[i];
1569 working_space[3 * ssize + i] = lda;
1572 for (i = 0; i < ssize; i++)
1573 working_space[i] = working_space[3 * ssize + i];
1578 for (i = 0; i < ssize; i++) {
1579 lda = working_space[i];
1582 working_space[ssize + j] = lda;
1586 for (i = 0; i < ssize; i++)
1587 source[i] = area * working_space[ssize + i];
1588 delete[]working_space;
1667 int ssize,
int numberIterations,
1668 int numberRepetitions,
Double_t boost )
1671 return "Wrong Parameters";
1673 if (numberRepetitions <= 0)
1674 return "Wrong Parameters";
1679 int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
1686 for (i = 0; i < ssize; i++) {
1690 working_space[ssize + i] = lda;
1691 if (lda > maximum) {
1696 if (lh_gold == -1) {
1697 delete [] working_space;
1698 return "ZERO RESPONSE VECTOR";
1702 for (i = 0; i < ssize; i++)
1703 working_space[2 * ssize + i] = source[i];
1706 for (i = 0; i < ssize; i++){
1707 if (i <= ssize - lh_gold)
1708 working_space[i] = 1;
1711 working_space[i] = 0;
1715 for (repet = 0; repet < numberRepetitions; repet++) {
1717 for (i = 0; i < ssize; i++)
1718 working_space[i] =
TMath::Power(working_space[i], boost);
1720 for (lindex = 0; lindex < numberIterations; lindex++) {
1721 for (i = 0; i <= ssize - lh_gold; i++){
1723 if (working_space[i] > 0){
1724 for (j = i; j < i + lh_gold; j++){
1725 ldb = working_space[2 * ssize + j];
1729 if (kmax > lh_gold - 1)
1731 kmin = j + lh_gold - ssize;
1735 for (k = kmax; k >= kmin; k--){
1736 ldc += working_space[ssize + k] * working_space[j - k];
1744 ldb = ldb * working_space[ssize + j - i];
1748 lda = lda * working_space[i];
1750 working_space[3 * ssize + i] = lda;
1752 for (i = 0; i < ssize; i++)
1753 working_space[i] = working_space[3 * ssize + i];
1758 for (i = 0; i < ssize; i++) {
1759 lda = working_space[i];
1762 working_space[ssize + j] = lda;
1766 for (i = 0; i < ssize; i++)
1767 source[i] = working_space[ssize + i];
1768 delete[]working_space;
1892 int ssizex,
int ssizey,
1893 int numberIterations,
1894 int numberRepetitions,
Double_t boost)
1896 int i, j, k, lindex, lhx = 0, repet;
1898 if (ssizex <= 0 || ssizey <= 0)
1899 return "Wrong Parameters";
1900 if (ssizex < ssizey)
1901 return "Sizex must be greater than sizey)";
1902 if (numberIterations <= 0)
1903 return "Number of iterations must be positive";
1905 new Double_t[ssizex * ssizey + 2 * ssizey * ssizey + 4 * ssizex];
1908 for (j = 0; j < ssizey && lhx != -1; j++) {
1911 for (i = 0; i < ssizex; i++) {
1912 lda = respMatrix[j][i];
1916 working_space[j * ssizex + i] = lda;
1920 for (i = 0; i < ssizex; i++)
1921 working_space[j * ssizex + i] /= area;
1925 delete [] working_space;
1926 return (
"ZERO COLUMN IN RESPONSE MATRIX");
1930 for (i = 0; i < ssizex; i++)
1931 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1935 for (i = 0; i < ssizey; i++) {
1936 for (j = 0; j < ssizey; j++) {
1938 for (k = 0; k < ssizex; k++) {
1939 ldb = working_space[ssizex * i + k];
1940 ldc = working_space[ssizex * j + k];
1941 lda = lda + ldb * ldc;
1943 working_space[ssizex * ssizey + ssizey * i + j] = lda;
1946 for (k = 0; k < ssizex; k++) {
1947 ldb = working_space[ssizex * i + k];
1949 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
1951 lda = lda + ldb * ldc;
1953 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
1958 for (i = 0; i < ssizey; i++)
1959 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1960 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
1963 for (i = 0; i < ssizey; i++) {
1964 for (j = 0; j < ssizey; j++) {
1966 for (k = 0; k < ssizey; k++) {
1967 ldb = working_space[ssizex * ssizey + ssizey * i + k];
1968 ldc = working_space[ssizex * ssizey + ssizey * j + k];
1969 lda = lda + ldb * ldc;
1971 working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j] =
1975 for (k = 0; k < ssizey; k++) {
1976 ldb = working_space[ssizex * ssizey + ssizey * i + k];
1978 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
1980 lda = lda + ldb * ldc;
1982 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
1987 for (i = 0; i < ssizey; i++)
1988 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1989 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
1992 for (i = 0; i < ssizey; i++)
1993 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = 1;
1996 for (repet = 0; repet < numberRepetitions; repet++) {
1998 for (i = 0; i < ssizey; i++)
1999 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
TMath::Power(working_space[ssizex * ssizey + 2 * ssizey * ssizey + i], boost);
2001 for (lindex = 0; lindex < numberIterations; lindex++) {
2002 for (i = 0; i < ssizey; i++) {
2004 for (j = 0; j < ssizey; j++) {
2006 working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j];
2007 ldc = working_space[ssizex * ssizey + 2 * ssizey * ssizey + j];
2008 lda = lda + ldb * ldc;
2011 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i];
2018 ldb = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2020 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] = lda;
2022 for (i = 0; i < ssizey; i++)
2023 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
2024 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
2029 for (i = 0; i < ssizex; i++) {
2031 source[i] = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2036 delete[]working_space;
2129 bool backgroundRemove,
int deconIterations,
2130 bool markov,
int averWindow)
2132 int i, j, numberIterations = (
Int_t)(7 *
sigma + 0.5);
2134 int k, lindex, posit, imin, imax, jmin, jmax, lh_gold, priz;
2135 Double_t lda, ldb, ldc, area, maximum, maximum_decon;
2136 int xmin,
xmax,
l, peak_index = 0, size_ext = ssize + 2 * numberIterations, shift = numberIterations, bw = 2, w;
2138 Double_t nom, nip, nim, sp, sm, plocha = 0;
2139 Double_t m0low=0,m1low=0,m2low=0,l0low=0,l1low=0,detlow,av,men;
2141 Error(
"SearchHighRes",
"Invalid sigma, must be greater than or equal to 1");
2145 if(threshold<=0 || threshold>=100){
2146 Error(
"SearchHighRes",
"Invalid threshold, must be positive and less than 100");
2152 Error(
"SearchHighRes",
"Too large sigma");
2156 if (markov ==
true) {
2157 if (averWindow <= 0) {
2158 Error(
"SearchHighRes",
"Averaging window must be positive");
2163 if(backgroundRemove ==
true){
2164 if(ssize < 2 * numberIterations + 1){
2165 Error(
"SearchHighRes",
"Too large clipping window");
2172 for(i = 0;i < k;i++){
2173 a = i,
b = source[i];
2174 m0low += 1,m1low +=
a,m2low +=
a *
a,l0low +=
b,l1low +=
a *
b;
2176 detlow = m0low * m2low - m1low * m1low;
2178 l1low = (-l0low * m1low + l1low * m0low) / detlow;
2193 for (j=0;j<7 * (ssize + i);j++) working_space[j] = 0;
2194 for(i = 0; i < size_ext; i++){
2197 working_space[i + size_ext] = source[0] + l1low *
a;
2198 if(working_space[i + size_ext] < 0)
2199 working_space[i + size_ext]=0;
2202 else if(i >= ssize + shift){
2203 a = i - (ssize - 1 + shift);
2204 working_space[i + size_ext] = source[ssize - 1];
2205 if(working_space[i + size_ext] < 0)
2206 working_space[i + size_ext]=0;
2210 working_space[i + size_ext] = source[i - shift];
2213 if(backgroundRemove ==
true){
2214 for(i = 1; i <= numberIterations; i++){
2215 for(j = i; j < size_ext - i; j++){
2216 if(markov ==
false){
2217 a = working_space[size_ext + j];
2218 b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2226 a = working_space[size_ext + j];
2229 for (w = j - bw; w <= j + bw; w++){
2230 if ( w >= 0 && w < size_ext){
2231 av += working_space[size_ext + w];
2238 for (w = j - i - bw; w <= j - i + bw; w++){
2239 if ( w >= 0 && w < size_ext){
2240 b += working_space[size_ext + w];
2247 for (w = j + i - bw; w <= j + i + bw; w++){
2248 if ( w >= 0 && w < size_ext){
2249 c += working_space[size_ext + w];
2257 working_space[j]=av;
2260 for(j = i; j < size_ext - i; j++)
2261 working_space[size_ext + j] = working_space[j];
2263 for(j = 0;j < size_ext; j++){
2266 b = source[0] + l1low *
a;
2268 working_space[size_ext + j] =
b - working_space[size_ext + j];
2271 else if(j >= ssize + shift){
2272 a = j - (ssize - 1 + shift);
2273 b = source[ssize - 1];
2275 working_space[size_ext + j] =
b - working_space[size_ext + j];
2279 working_space[size_ext + j] = source[j - shift] - working_space[size_ext + j];
2282 for(j = 0;j < size_ext; j++){
2283 if(working_space[size_ext + j] < 0) working_space[size_ext + j] = 0;
2287 for(i = 0; i < size_ext; i++){
2288 working_space[i + 6*size_ext] = working_space[i + size_ext];
2292 for(j = 0; j < size_ext; j++)
2293 working_space[2 * size_ext + j] = working_space[size_ext + j];
2295 for(i = 0, maxch = 0; i < size_ext; i++){
2296 working_space[i] = 0;
2297 if(maxch < working_space[2 * size_ext + i])
2298 maxch = working_space[2 * size_ext + i];
2299 plocha += working_space[2 * size_ext + i];
2302 delete [] working_space;
2307 working_space[
xmin] = 1;
2309 nip = working_space[2 * size_ext + i] / maxch;
2310 nim = working_space[2 * size_ext + i + 1] / maxch;
2312 for(
l = 1;
l <= averWindow;
l++){
2314 a = working_space[2 * size_ext +
xmax] / maxch;
2317 a = working_space[2 * size_ext + i +
l] / maxch;
2329 if((i -
l + 1) <
xmin)
2330 a = working_space[2 * size_ext +
xmin] / maxch;
2333 a = working_space[2 * size_ext + i -
l + 1] / maxch;
2347 a = working_space[i + 1] = working_space[i] *
a;
2351 working_space[i] = working_space[i] / nom;
2353 for(j = 0; j < size_ext; j++)
2354 working_space[size_ext + j] = working_space[j] * plocha;
2355 for(j = 0; j < size_ext; j++){
2356 working_space[2 * size_ext + j] = working_space[size_ext + j];
2358 if(backgroundRemove ==
true){
2359 for(i = 1; i <= numberIterations; i++){
2360 for(j = i; j < size_ext - i; j++){
2361 a = working_space[size_ext + j];
2362 b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2365 working_space[j] =
a;
2367 for(j = i; j < size_ext - i; j++)
2368 working_space[size_ext + j] = working_space[j];
2370 for(j = 0; j < size_ext; j++){
2371 working_space[size_ext + j] = working_space[2 * size_ext + j] - working_space[size_ext + j];
2381 for(i = 0; i < size_ext; i++){
2389 working_space[i] = lda;
2397 for(i = 0; i < size_ext; i++)
2398 working_space[2 * size_ext + i] =
TMath::Abs(working_space[size_ext + i]);
2405 for(i = imin; i <= imax; i++){
2410 jmax = lh_gold - 1 - i;
2411 if(jmax > (lh_gold - 1))
2414 for(j = jmin;j <= jmax; j++){
2415 ldb = working_space[j];
2416 ldc = working_space[i + j];
2417 lda = lda + ldb * ldc;
2419 working_space[size_ext + i - imin] = lda;
2423 imin = -i,imax = size_ext + i - 1;
2424 for(i = imin; i <= imax; i++){
2426 for(j = 0; j <= (lh_gold - 1); j++){
2427 ldb = working_space[j];
2429 if(k >= 0 && k < size_ext){
2430 ldc = working_space[2 * size_ext + k];
2431 lda = lda + ldb * ldc;
2435 working_space[4 * size_ext + i - imin] = lda;
2438 for(i = imin; i <= imax; i++)
2439 working_space[2 * size_ext + i - imin] = working_space[4 * size_ext + i - imin];
2441 for(i = 0; i < size_ext; i++)
2442 working_space[i] = 1;
2444 for(lindex = 0; lindex < deconIterations; lindex++){
2445 for(i = 0; i < size_ext; i++){
2446 if(
TMath::Abs(working_space[2 * size_ext + i]) > 0.00001 &&
TMath::Abs(working_space[i]) > 0.00001){
2454 if(jmax > (size_ext - 1 - i))
2457 for(j = jmin; j <= jmax; j++){
2458 ldb = working_space[j + lh_gold - 1 + size_ext];
2459 ldc = working_space[i + j];
2460 lda = lda + ldb * ldc;
2462 ldb = working_space[2 * size_ext + i];
2469 ldb = working_space[i];
2471 working_space[3 * size_ext + i] = lda;
2474 for(i = 0; i < size_ext; i++){
2475 working_space[i] = working_space[3 * size_ext + i];
2479 for(i=0;i<size_ext;i++){
2480 lda = working_space[i];
2483 working_space[size_ext + j] = lda;
2486 maximum = 0, maximum_decon = 0;
2488 for(i = 0; i < size_ext - j; i++){
2489 if(i >= shift && i < ssize + shift){
2490 working_space[i] = area * working_space[size_ext + i + j];
2491 if(maximum_decon < working_space[i])
2492 maximum_decon = working_space[i];
2493 if(maximum < working_space[6 * size_ext + i])
2494 maximum = working_space[6 * size_ext + i];
2498 working_space[i] = 0;
2506 for(i = 1; i < size_ext - 1; i++){
2507 if(working_space[i] > working_space[i - 1] && working_space[i] > working_space[i + 1]){
2508 if(i >= shift && i < ssize + shift){
2509 if(working_space[i] > lda*maximum_decon && working_space[6 * size_ext + i] > threshold * maximum / 100.0){
2510 for(j = i - 1,
a = 0,
b = 0; j <= i + 1; j++){
2511 a += (
Double_t)(j - shift) * working_space[j];
2512 b += working_space[j];
2520 if(peak_index == 0){
2526 for(j = 0, priz = 0; j < peak_index && priz == 0; j++){
2527 if(working_space[6 * size_ext + shift + (
Int_t)
a] > working_space[6 * size_ext + shift + (
Int_t)
fPositionX[j]])
2537 for(k = peak_index; k >= j; k--){
2552 for(i = 0; i < ssize; i++) destVector[i] = working_space[i + shift];
2553 delete [] working_space;
2556 Warning(
"SearchHighRes",
"Peak buffer full");
2566 bool backgroundRemove,
int deconIterations,
2567 bool markov,
int averWindow)
2572 deconIterations,markov,averWindow);
virtual void SetLineColor(Color_t lcolor)
Set the line color.
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Int_t GetLast() const
Return last bin on the axis i.e.
Int_t GetFirst() const
Return first bin on the axis i.e.
TH1 is the base class of all histogram classes in ROOT.
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
virtual Int_t GetDimension() const
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
TList * GetListOfFunctions() const
virtual void Draw(Option_t *option="")
Draw this histogram with options.
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
virtual void SetEntries(Double_t n)
virtual void Add(TObject *obj)
virtual TObject * Remove(TObject *obj)
Remove object from the list.
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
The TNamed class is the base class for all named ROOT classes.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
A PolyMarker is defined by an array on N points in a 2-D space.
Advanced Spectra Processing.
Double_t fResolution
NOT USED resolution of the neighboring peaks
static TH1 * StaticBackground(const TH1 *hist, Int_t niter=20, Option_t *option="")
Static function, interface to TSpectrum::Background.
Int_t fMaxPeaks
Maximum number of peaks to be found.
Double_t * fPositionX
[fNPeaks] X position of peaks
const char * SmoothMarkov(Double_t *source, Int_t ssize, Int_t averWindow)
One-dimensional markov spectrum smoothing function.
const char * Unfolding(Double_t *source, const Double_t **respMatrix, Int_t ssizex, Int_t ssizey, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional unfolding function.
void SetResolution(Double_t resolution=1)
NOT USED resolution: determines resolution of the neighbouring peaks default value is 1 correspond to...
virtual ~TSpectrum()
Destructor.
virtual TH1 * Background(const TH1 *hist, Int_t niter=20, Option_t *option="")
One-dimensional background estimation function.
Int_t SearchHighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
One-dimensional high-resolution peak search function.
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
One-dimensional peak search function.
Int_t fNPeaks
number of peaks found
Int_t Search1HighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
Old name of SearcHighRes introduced for back compatibility.
const char * Deconvolution(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
static Int_t fgAverageWindow
Average window of searched peaks.
static void SetDeconIterations(Int_t n=3)
Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::Search...
Double_t * fPosition
[fNPeaks] array of current peak positions
static void SetAverageWindow(Int_t w=3)
Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).
const char * DeconvolutionRL(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
static Int_t fgIterations
Maximum number of decon iterations (default=3)
TH1 * fHistogram
resulting histogram
virtual void Print(Option_t *option="") const
Print the array of positions.
static Int_t StaticSearch(const TH1 *hist, Double_t sigma=2, Option_t *option="goff", Double_t threshold=0.05)
Static function, interface to TSpectrum::Search.
Double_t * fPositionY
[fNPeaks] Y position of peaks
void ToLower()
Change string to lower-case.
const char * Data() const
TString & ReplaceAll(const TString &s1, const TString &s2)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Double_t Sqrt(Double_t x)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
#define dest(otri, vertexptr)