44 #define PEAK_WINDOW 1024 71 :
TNamed(
"Spectrum",
"Miroslav Morhac peak finder")
155 if (h == 0)
return 0;
158 Error(
"Search",
"Only implemented for 1-d histograms");
185 Int_t size = last-first+1;
188 for (i = 0; i < size; i++) source[i] = h->
GetBinContent(i + first);
191 Background(source,size,numberIterations, direction, filterOrder,smoothing,
192 smoothWindow,compton);
197 char *
hbname =
new char[nch+20];
201 hb->GetListOfFunctions()->
Delete();
203 for (i=0; i< size; i++) hb->SetBinContent(i+first,source[i]);
204 hb->SetEntries(size);
208 if (
gPad)
delete gPad->GetPrimitive(hbname);
221 printf(
"\nNumber of positions = %d\n",
fNPeaks);
269 if (hin == 0)
return 0;
272 Error(
"Search",
"Only implemented for 1-d and 2-d histograms");
275 if (threshold <=0 || threshold >= 1) {
276 Warning(
"Search",
"threshold must 0<threshold<1, threshold=0.05 assumed");
296 if (dimension == 1) {
299 Int_t size = last-first+1;
300 Int_t i, bin, npeaks;
303 for (i = 0; i < size; i++) source[i] = hin->
GetBinContent(i + first);
306 if (sigma < 1) sigma = 1;
307 if (sigma > 8) sigma = 8;
309 npeaks =
SearchHighRes(source, dest, size, sigma, 100*threshold,
312 for (i = 0; i < npeaks; i++) {
322 if (!npeaks)
return 0;
513 int numberIterations,
514 int direction,
int filterOrder,
515 bool smoothing,
int smoothWindow,
518 int i, j, w, bw, b1, b2, priz;
519 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;
521 return "Wrong Parameters";
522 if (numberIterations < 1)
523 return "Width of Clipping Window Must Be Positive";
524 if (ssize < 2 * numberIterations + 1)
525 return "Too Large Clipping Window";
527 return "Incorrect width of smoothing window";
529 for (i = 0; i < ssize; i++){
530 working_space[i] = spectrum[i];
531 working_space[i + ssize] = spectrum[i];
533 bw=(smoothWindow-1)/2;
537 i = numberIterations;
540 for (j = i; j < ssize - i; j++) {
542 a = working_space[ssize + j];
543 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
546 working_space[j] =
a;
549 else if (smoothing ==
kTRUE){
550 a = working_space[ssize + j];
553 for (w = j - bw; w <= j + bw; w++){
554 if ( w >= 0 && w < ssize){
555 av += working_space[ssize + w];
562 for (w = j - i - bw; w <= j - i + bw; w++){
563 if ( w >= 0 && w < ssize){
564 b += working_space[ssize + w];
571 for (w = j + i - bw; w <= j + i + bw; w++){
572 if ( w >= 0 && w < ssize){
573 c += working_space[ssize + w];
584 for (j = i; j < ssize - i; j++)
585 working_space[ssize + j] = working_space[j];
595 for (j = i; j < ssize - i; j++) {
597 a = working_space[ssize + j];
598 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
601 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
602 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
603 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
604 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
609 working_space[j] =
a;
612 else if (smoothing ==
kTRUE){
613 a = working_space[ssize + j];
616 for (w = j - bw; w <= j + bw; w++){
617 if ( w >= 0 && w < ssize){
618 av += working_space[ssize + w];
625 for (w = j - i - bw; w <= j - i + bw; w++){
626 if ( w >= 0 && w < ssize){
627 b += working_space[ssize + w];
634 for (w = j + i - bw; w <= j + i + bw; w++){
635 if ( w >= 0 && w < ssize){
636 c += working_space[ssize + w];
644 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
645 if (w >= 0 && w < ssize){
646 b4 += working_space[ssize + w];
652 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
653 if (w >= 0 && w < ssize){
654 c4 += working_space[ssize + w];
660 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
661 if (w >= 0 && w < ssize){
662 d4 += working_space[ssize + w];
668 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
669 if (w >= 0 && w < ssize){
670 e4 += working_space[ssize + w];
675 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
683 for (j = i; j < ssize - i; j++)
684 working_space[ssize + j] = working_space[j];
694 for (j = i; j < ssize - i; j++) {
696 a = working_space[ssize + j];
697 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
700 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
701 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
702 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
703 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
706 d += working_space[ssize + j - (
Int_t) (3 * ai)] / 20;
707 d -= 6 * working_space[ssize + j - (
Int_t) (2 * ai)] / 20;
708 d += 15 * working_space[ssize + j - (
Int_t) ai] / 20;
709 d += 15 * working_space[ssize + j + (
Int_t) ai] / 20;
710 d -= 6 * working_space[ssize + j + (
Int_t) (2 * ai)] / 20;
711 d += working_space[ssize + j + (
Int_t) (3 * ai)] / 20;
718 working_space[j] =
a;
721 else if (smoothing ==
kTRUE){
722 a = working_space[ssize + j];
725 for (w = j - bw; w <= j + bw; w++){
726 if ( w >= 0 && w < ssize){
727 av += working_space[ssize + w];
734 for (w = j - i - bw; w <= j - i + bw; w++){
735 if ( w >= 0 && w < ssize){
736 b += working_space[ssize + w];
743 for (w = j + i - bw; w <= j + i + bw; w++){
744 if ( w >= 0 && w < ssize){
745 c += working_space[ssize + w];
753 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
754 if (w >= 0 && w < ssize){
755 b4 += working_space[ssize + w];
761 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
762 if (w >= 0 && w < ssize){
763 c4 += working_space[ssize + w];
769 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
770 if (w >= 0 && w < ssize){
771 d4 += working_space[ssize + w];
777 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
778 if (w >= 0 && w < ssize){
779 e4 += working_space[ssize + w];
784 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
787 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
788 if (w >= 0 && w < ssize){
789 b6 += working_space[ssize + w];
795 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
796 if (w >= 0 && w < ssize){
797 c6 += working_space[ssize + w];
803 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
804 if (w >= 0 && w < ssize){
805 d6 += working_space[ssize + w];
811 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
812 if (w >= 0 && w < ssize){
813 e6 += working_space[ssize + w];
819 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
820 if (w >= 0 && w < ssize){
821 f6 += working_space[ssize + w];
827 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
828 if (w >= 0 && w < ssize){
829 g6 += working_space[ssize + w];
834 b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
844 for (j = i; j < ssize - i; j++)
845 working_space[ssize + j] = working_space[j];
855 for (j = i; j < ssize - i; j++) {
857 a = working_space[ssize + j];
858 b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
861 c -= working_space[ssize + j - (
Int_t) (2 * ai)] / 6;
862 c += 4 * working_space[ssize + j - (
Int_t) ai] / 6;
863 c += 4 * working_space[ssize + j + (
Int_t) ai] / 6;
864 c -= working_space[ssize + j + (
Int_t) (2 * ai)] / 6;
867 d += working_space[ssize + j - (
Int_t) (3 * ai)] / 20;
868 d -= 6 * working_space[ssize + j - (
Int_t) (2 * ai)] / 20;
869 d += 15 * working_space[ssize + j - (
Int_t) ai] / 20;
870 d += 15 * working_space[ssize + j + (
Int_t) ai] / 20;
871 d -= 6 * working_space[ssize + j + (
Int_t) (2 * ai)] / 20;
872 d += working_space[ssize + j + (
Int_t) (3 * ai)] / 20;
875 e -= working_space[ssize + j - (
Int_t) (4 * ai)] / 70;
876 e += 8 * working_space[ssize + j - (
Int_t) (3 * ai)] / 70;
877 e -= 28 * working_space[ssize + j - (
Int_t) (2 * ai)] / 70;
878 e += 56 * working_space[ssize + j - (
Int_t) ai] / 70;
879 e += 56 * working_space[ssize + j + (
Int_t) ai] / 70;
880 e -= 28 * working_space[ssize + j + (
Int_t) (2 * ai)] / 70;
881 e += 8 * working_space[ssize + j + (
Int_t) (3 * ai)] / 70;
882 e -= working_space[ssize + j + (
Int_t) (4 * ai)] / 70;
891 working_space[j] =
a;
894 else if (smoothing ==
kTRUE){
895 a = working_space[ssize + j];
898 for (w = j - bw; w <= j + bw; w++){
899 if ( w >= 0 && w < ssize){
900 av += working_space[ssize + w];
907 for (w = j - i - bw; w <= j - i + bw; w++){
908 if ( w >= 0 && w < ssize){
909 b += working_space[ssize + w];
916 for (w = j + i - bw; w <= j + i + bw; w++){
917 if ( w >= 0 && w < ssize){
918 c += working_space[ssize + w];
926 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
927 if (w >= 0 && w < ssize){
928 b4 += working_space[ssize + w];
934 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
935 if (w >= 0 && w < ssize){
936 c4 += working_space[ssize + w];
942 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
943 if (w >= 0 && w < ssize){
944 d4 += working_space[ssize + w];
950 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
951 if (w >= 0 && w < ssize){
952 e4 += working_space[ssize + w];
957 b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
960 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
961 if (w >= 0 && w < ssize){
962 b6 += working_space[ssize + w];
968 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
969 if (w >= 0 && w < ssize){
970 c6 += working_space[ssize + w];
976 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
977 if (w >= 0 && w < ssize){
978 d6 += working_space[ssize + w];
984 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
985 if (w >= 0 && w < ssize){
986 e6 += working_space[ssize + w];
992 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
993 if (w >= 0 && w < ssize){
994 f6 += working_space[ssize + w];
1000 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
1001 if (w >= 0 && w < ssize){
1002 g6 += working_space[ssize + w];
1007 b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
1010 for (w = j - (
Int_t)(4 * ai) - bw; w <= j - (
Int_t)(4 * ai) + bw; w++){
1011 if (w >= 0 && w < ssize){
1012 b8 += working_space[ssize + w];
1018 for (w = j - (
Int_t)(3 * ai) - bw; w <= j - (
Int_t)(3 * ai) + bw; w++){
1019 if (w >= 0 && w < ssize){
1020 c8 += working_space[ssize + w];
1026 for (w = j - (
Int_t)(2 * ai) - bw; w <= j - (
Int_t)(2 * ai) + bw; w++){
1027 if (w >= 0 && w < ssize){
1028 d8 += working_space[ssize + w];
1034 for (w = j - (
Int_t)ai - bw; w <= j - (
Int_t)ai + bw; w++){
1035 if (w >= 0 && w < ssize){
1036 e8 += working_space[ssize + w];
1042 for (w = j + (
Int_t)ai - bw; w <= j + (
Int_t)ai + bw; w++){
1043 if (w >= 0 && w < ssize){
1044 f8 += working_space[ssize + w];
1050 for (w = j + (
Int_t)(2 * ai) - bw; w <= j + (
Int_t)(2 * ai) + bw; w++){
1051 if (w >= 0 && w < ssize){
1052 g8 += working_space[ssize + w];
1058 for (w = j + (
Int_t)(3 * ai) - bw; w <= j + (
Int_t)(3 * ai) + bw; w++){
1059 if (w >= 0 && w < ssize){
1060 h8 += working_space[ssize + w];
1066 for (w = j + (
Int_t)(4 * ai) - bw; w <= j + (
Int_t)(4 * ai) + bw; w++){
1067 if (w >= 0 && w < ssize){
1068 i8 += working_space[ssize + w];
1073 b8 = ( -b8 + 8 * c8 - 28 * d8 + 56 * e8 - 56 * f8 - 28 * g8 + 8 * h8 - i8)/70;
1082 working_space[j]=av;
1085 for (j = i; j < ssize - i; j++)
1086 working_space[ssize + j] = working_space[j];
1094 if (compton ==
kTRUE) {
1095 for (i = 0, b2 = 0; i < ssize; i++){
1097 a = working_space[i], b = spectrum[i];
1103 yb1 = working_space[b1];
1104 for (b2 = b1 + 1, c = 0, priz = 0; priz == 0 && b2 < ssize; b2++){
1105 a = working_space[b2], b = spectrum[b2];
1114 yb2 = working_space[b2];
1116 for (j = b1, c = 0; j <= b2; j++){
1121 c = (yb2 - yb1) / c;
1122 for (j = b1, d = 0; j <= b2 && j < ssize; j++){
1126 working_space[ssize + j] =
a;
1132 for (j = b2, c = 0; j >= b1; j--){
1137 c = (yb1 - yb2) / c;
1138 for (j = b2, d = 0;j >= b1 && j >= 0; j--){
1142 working_space[ssize + j] =
a;
1151 for (j = 0; j < ssize; j++)
1152 spectrum[j] = working_space[ssize + j];
1153 delete[]working_space;
1199 Double_t nom, nip, nim, sp, sm, area = 0;
1201 return "Averaging Window must be positive";
1203 xmin = 0,xmax = ssize - 1;
1204 for(i = 0, maxch = 0; i < ssize; i++){
1206 if(maxch < source[i])
1212 delete [] working_space;
1217 working_space[
xmin] = 1;
1218 for(i = xmin; i <
xmax; i++){
1219 nip = source[i] / maxch;
1220 nim = source[i + 1] / maxch;
1222 for(l = 1; l <= averWindow; l++){
1224 a = source[
xmax] / maxch;
1227 a = source[i +
l] / maxch;
1237 if((i - l + 1) < xmin)
1238 a = source[
xmin] / maxch;
1241 a = source[i - l + 1] / maxch;
1252 a = working_space[i + 1] = working_space[i] *
a;
1255 for(i = xmin; i <=
xmax; i++){
1256 working_space[i] = working_space[i] / nom;
1258 for(i = 0; i < ssize; i++)
1259 source[i] = working_space[i] * area;
1260 delete[]working_space;
1460 int ssize,
int numberIterations,
1464 return "Wrong Parameters";
1466 if (numberRepetitions <= 0)
1467 return "Wrong Parameters";
1472 int i, j, k, lindex, posit, lh_gold,
l, repet;
1473 Double_t lda, ldb, ldc, area, maximum;
1480 for (i = 0; i < ssize; i++) {
1484 working_space[i] = lda;
1486 if (lda > maximum) {
1491 if (lh_gold == -1) {
1492 delete [] working_space;
1493 return "ZERO RESPONSE VECTOR";
1497 for (i = 0; i < ssize; i++)
1498 working_space[2 * ssize + i] = source[i];
1501 for (i = 0; i < ssize; i++){
1503 for (j = 0; j < ssize; j++){
1504 ldb = working_space[j];
1507 ldc = working_space[k];
1508 lda = lda + ldb * ldc;
1511 working_space[ssize + i] = lda;
1513 for (k = 0; k < ssize; k++){
1516 ldb = working_space[
l];
1517 ldc = working_space[2 * ssize + k];
1518 lda = lda + ldb * ldc;
1521 working_space[3 * ssize + i]=lda;
1525 for (i = 0; i < ssize; i++){
1526 working_space[2 * ssize + i] = working_space[3 * ssize + i];
1530 for (i = 0; i < ssize; i++)
1531 working_space[i] = 1;
1534 for (repet = 0; repet < numberRepetitions; repet++) {
1536 for (i = 0; i < ssize; i++)
1537 working_space[i] =
TMath::Power(working_space[i], boost);
1539 for (lindex = 0; lindex < numberIterations; lindex++) {
1540 for (i = 0; i < ssize; i++) {
1541 if (working_space[2 * ssize + i] > 0.000001
1542 && working_space[i] > 0.000001) {
1544 for (j = 0; j < lh_gold; j++) {
1545 ldb = working_space[j + ssize];
1550 ldc = working_space[k];
1553 ldc += working_space[k];
1557 ldc = working_space[i];
1558 lda = lda + ldb * ldc;
1560 ldb = working_space[2 * ssize + i];
1566 ldb = working_space[i];
1568 working_space[3 * ssize + i] = lda;
1571 for (i = 0; i < ssize; i++)
1572 working_space[i] = working_space[3 * ssize + i];
1577 for (i = 0; i < ssize; i++) {
1578 lda = working_space[i];
1581 working_space[ssize + j] = lda;
1585 for (i = 0; i < ssize; i++)
1586 source[i] = area * working_space[ssize + i];
1587 delete[]working_space;
1666 int ssize,
int numberIterations,
1667 int numberRepetitions,
Double_t boost )
1670 return "Wrong Parameters";
1672 if (numberRepetitions <= 0)
1673 return "Wrong Parameters";
1678 int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
1685 for (i = 0; i < ssize; i++) {
1689 working_space[ssize + i] = lda;
1690 if (lda > maximum) {
1695 if (lh_gold == -1) {
1696 delete [] working_space;
1697 return "ZERO RESPONSE VECTOR";
1701 for (i = 0; i < ssize; i++)
1702 working_space[2 * ssize + i] = source[i];
1705 for (i = 0; i < ssize; i++){
1706 if (i <= ssize - lh_gold)
1707 working_space[i] = 1;
1710 working_space[i] = 0;
1714 for (repet = 0; repet < numberRepetitions; repet++) {
1716 for (i = 0; i < ssize; i++)
1717 working_space[i] =
TMath::Power(working_space[i], boost);
1719 for (lindex = 0; lindex < numberIterations; lindex++) {
1720 for (i = 0; i <= ssize - lh_gold; i++){
1722 if (working_space[i] > 0){
1723 for (j = i; j < i + lh_gold; j++){
1724 ldb = working_space[2 * ssize + j];
1728 if (kmax > lh_gold - 1)
1730 kmin = j + lh_gold - ssize;
1734 for (k = kmax; k >= kmin; k--){
1735 ldc += working_space[ssize + k] * working_space[j - k];
1743 ldb = ldb * working_space[ssize + j - i];
1747 lda = lda * working_space[i];
1749 working_space[3 * ssize + i] = lda;
1751 for (i = 0; i < ssize; i++)
1752 working_space[i] = working_space[3 * ssize + i];
1757 for (i = 0; i < ssize; i++) {
1758 lda = working_space[i];
1761 working_space[ssize + j] = lda;
1765 for (i = 0; i < ssize; i++)
1766 source[i] = working_space[ssize + i];
1767 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");
2150 j = (
Int_t) (5.0 * sigma + 0.5);
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");
2170 k = int(2 * sigma+0.5);
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;
2190 i = (
Int_t)(7 * sigma + 0.5);
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];
2294 xmin = 0,xmax = size_ext - 1;
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;
2308 for(i = xmin; i <
xmax; i++){
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;
2350 for(i = xmin; i <=
xmax; i++){
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++){
2383 lda = lda * lda / (2 * sigma *
sigma);
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)
2571 return SearchHighRes(source,destVector,ssize,sigma,threshold,backgroundRemove,
2572 deconIterations,markov,averWindow);
2582 return s.
Search(hist,sigma,option,threshold);
virtual const char * GetName() const
Returns name of object.
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.
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Int_t GetFirst() const
Return first bin on the axis i.e.
TString & ReplaceAll(const TString &s1, const TString &s2)
static TH1 * StaticBackground(const TH1 *hist, Int_t niter=20, Option_t *option="")
Static function, interface to TSpectrum::Background.
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.
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
void ToLower()
Change string to lower-case.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
virtual Int_t GetDimension() const
The TNamed class is the base class for all named ROOT classes.
Double_t * fPositionX
[fNPeaks] X position of peaks
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 void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
virtual void Delete(Option_t *option="")
Delete this object.
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
One-dimensional peak search function.
virtual ~TSpectrum()
Destructor.
Int_t GetLast() const
Return last bin on the axis i.e.
virtual TObject * Remove(TObject *obj)
Remove object from the list.
LVector boost(const LVector &v, const BoostVector &b)
Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost The only re...
TH1 * fHistogram
resulting histogram
Double_t * fPositionY
[fNPeaks] Y position of peaks
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
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 void SetAverageWindow(Int_t w=3)
Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).
static Int_t fgAverageWindow
Average window of searched peaks.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
static void SetDeconIterations(Int_t n=3)
Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::Search...
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
A PolyMarker is defined by an array on N points in a 2-D space.
virtual void Print(Option_t *option="") const
Print the array of positions.
virtual TH1 * Background(const TH1 *hist, Int_t niter=20, Option_t *option="")
One-dimensional background estimation function.
static Int_t fgIterations
Maximum number of decon iterations (default=3)
Int_t fMaxPeaks
Maximum number of peaks to be found.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
void SetResolution(Double_t resolution=1)
NOT USED resolution: determines resolution of the neighbouring peaks default value is 1 correspond to...
Advanced Spectra Processing.
Double_t * fPosition
[fNPeaks] array of current peak positions
virtual void Add(TObject *obj)
#define dest(otri, vertexptr)
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.
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Int_t fNPeaks
number of peaks found
Double_t Sqrt(Double_t x)
TList * GetListOfFunctions() const
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.
Double_t fResolution
NOT USED resolution of the neighboring peaks
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
const char * Data() const