41 Int_t TSpectrum::fgIterations = 3;
42 Int_t TSpectrum::fgAverageWindow = 3;
44 #define PEAK_WINDOW 1024
50 TSpectrum::TSpectrum() :TNamed(
"Spectrum",
"Miroslav Morhac peak finder")
54 fPosition =
new Double_t[n];
55 fPositionX =
new Double_t[n];
56 fPositionY =
new Double_t[n];
70 TSpectrum::TSpectrum(Int_t maxpositions, Double_t resolution)
71 :TNamed(
"Spectrum",
"Miroslav Morhac peak finder")
73 Int_t n = maxpositions;
76 fPosition =
new Double_t[n];
77 fPositionX =
new Double_t[n];
78 fPositionY =
new Double_t[n];
81 SetResolution(resolution);
87 TSpectrum::~TSpectrum()
99 void TSpectrum::SetAverageWindow(Int_t w)
108 void TSpectrum::SetDeconIterations(Int_t n)
152 TH1 *TSpectrum::Background(
const TH1 * h,
int numberIterations,
155 if (h == 0)
return 0;
156 Int_t dimension = h->GetDimension();
158 Error(
"Search",
"Only implemented for 1-d histograms");
161 TString opt = option;
165 Int_t direction = kBackDecreasingWindow;
166 if (opt.Contains(
"backincreasingwindow")) direction = kBackIncreasingWindow;
167 Int_t filterOrder = kBackOrder2;
168 if (opt.Contains(
"backorder4")) filterOrder = kBackOrder4;
169 if (opt.Contains(
"backorder6")) filterOrder = kBackOrder6;
170 if (opt.Contains(
"backorder8")) filterOrder = kBackOrder8;
171 Bool_t smoothing = kTRUE;
172 if (opt.Contains(
"nosmoothing")) smoothing = kFALSE;
173 Int_t smoothWindow = kBackSmoothing3;
174 if (opt.Contains(
"backsmoothing5")) smoothWindow = kBackSmoothing5;
175 if (opt.Contains(
"backsmoothing7")) smoothWindow = kBackSmoothing7;
176 if (opt.Contains(
"backsmoothing9")) smoothWindow = kBackSmoothing9;
177 if (opt.Contains(
"backsmoothing11")) smoothWindow = kBackSmoothing11;
178 if (opt.Contains(
"backsmoothing13")) smoothWindow = kBackSmoothing13;
179 if (opt.Contains(
"backsmoothing15")) smoothWindow = kBackSmoothing15;
180 Bool_t compton = kFALSE;
181 if (opt.Contains(
"compton")) compton = kTRUE;
183 Int_t first = h->GetXaxis()->GetFirst();
184 Int_t last = h->GetXaxis()->GetLast();
185 Int_t size = last-first+1;
187 Double_t * source =
new Double_t[size];
188 for (i = 0; i < size; i++) source[i] = h->GetBinContent(i + first);
191 Background(source,size,numberIterations, direction, filterOrder,smoothing,
192 smoothWindow,compton);
196 Int_t nch = strlen(h->GetName());
197 char *hbname =
new char[nch+20];
198 snprintf(hbname,nch+20,
"%s_background",h->GetName());
199 TH1 *hb = (TH1*)h->Clone(hbname);
201 hb->GetListOfFunctions()->Delete();
203 for (i=0; i< size; i++) hb->SetBinContent(i+first,source[i]);
204 hb->SetEntries(size);
207 if (opt.Contains(
"same")) {
208 if (gPad)
delete gPad->GetPrimitive(hbname);
219 void TSpectrum::Print(Option_t *)
const
221 printf(
"\nNumber of positions = %d\n",fNPeaks);
222 for (Int_t i=0;i<fNPeaks;i++) {
223 printf(
" x[%d] = %g, y[%d] = %g\n",i,fPositionX[i],i,fPositionY[i]);
266 Int_t TSpectrum::Search(
const TH1 * hin, Double_t sigma, Option_t * option,
269 if (hin == 0)
return 0;
270 Int_t dimension = hin->GetDimension();
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");
279 TString opt = option;
281 Bool_t background = kTRUE;
282 if (opt.Contains(
"nobackground")) {
284 opt.ReplaceAll(
"nobackground",
"");
286 Bool_t markov = kTRUE;
287 if (opt.Contains(
"nomarkov")) {
289 opt.ReplaceAll(
"nomarkov",
"");
292 if (opt.Contains(
"nodraw")) {
294 opt.ReplaceAll(
"nodraw",
"");
296 if (dimension == 1) {
297 Int_t first = hin->GetXaxis()->GetFirst();
298 Int_t last = hin->GetXaxis()->GetLast();
299 Int_t size = last-first+1;
300 Int_t i, bin, npeaks;
301 Double_t * source =
new Double_t[size];
302 Double_t * dest =
new Double_t[size];
303 for (i = 0; i < size; i++) source[i] = hin->GetBinContent(i + first);
305 sigma = size/fMaxPeaks;
306 if (sigma < 1) sigma = 1;
307 if (sigma > 8) sigma = 8;
309 npeaks = SearchHighRes(source, dest, size, sigma, 100*threshold,
310 background, fgIterations, markov, fgAverageWindow);
312 for (i = 0; i < npeaks; i++) {
313 bin = first + Int_t(fPositionX[i] + 0.5);
314 fPositionX[i] = hin->GetBinCenter(bin);
315 fPositionY[i] = hin->GetBinContent(bin);
320 if (opt.Contains(
"goff"))
322 if (!npeaks)
return 0;
324 (TPolyMarker*)hin->GetListOfFunctions()->FindObject(
"TPolyMarker");
326 hin->GetListOfFunctions()->Remove(pm);
329 pm =
new TPolyMarker(npeaks, fPositionX, fPositionY);
330 hin->GetListOfFunctions()->Add(pm);
331 pm->SetMarkerStyle(23);
332 pm->SetMarkerColor(kRed);
333 pm->SetMarkerSize(1.3);
334 opt.ReplaceAll(
" ",
"");
335 opt.ReplaceAll(
",",
"");
337 ((TH1*)hin)->Draw(opt.Data());
351 void TSpectrum::SetResolution(Double_t resolution)
354 fResolution = resolution;
512 const char *TSpectrum::Background(Double_t *spectrum,
int ssize,
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";
526 if (smoothing == kTRUE && smoothWindow != kBackSmoothing3 && smoothWindow != kBackSmoothing5 && smoothWindow != kBackSmoothing7 && smoothWindow != kBackSmoothing9 && smoothWindow != kBackSmoothing11 && smoothWindow != kBackSmoothing13 && smoothWindow != kBackSmoothing15)
527 return "Incorrect width of smoothing window";
528 Double_t *working_space =
new Double_t[2 * ssize];
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;
534 if (direction == kBackIncreasingWindow)
536 else if(direction == kBackDecreasingWindow)
537 i = numberIterations;
538 if (filterOrder == kBackOrder2) {
540 for (j = i; j < ssize - i; j++) {
541 if (smoothing == kFALSE){
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];
586 if (direction == kBackIncreasingWindow)
588 else if(direction == kBackDecreasingWindow)
590 }
while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
593 else if (filterOrder == kBackOrder4) {
595 for (j = i; j < ssize - i; j++) {
596 if (smoothing == kFALSE){
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];
685 if (direction == kBackIncreasingWindow)
687 else if(direction == kBackDecreasingWindow)
689 }
while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
692 else if (filterOrder == kBackOrder6) {
694 for (j = i; j < ssize - i; j++) {
695 if (smoothing == kFALSE){
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];
846 if (direction == kBackIncreasingWindow)
848 else if(direction == kBackDecreasingWindow)
850 }
while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
853 else if (filterOrder == kBackOrder8) {
855 for (j = i; j < ssize - i; j++) {
856 if (smoothing == kFALSE){
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];
1087 if (direction == kBackIncreasingWindow)
1089 else if(direction == kBackDecreasingWindow)
1091 }
while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
1094 if (compton == kTRUE) {
1095 for (i = 0, b2 = 0; i < ssize; i++){
1097 a = working_space[i], b = spectrum[i];
1099 if (TMath::Abs(a - b) >= 1) {
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];
1107 if (TMath::Abs(a - b) < 1) {
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;
1195 const char* TSpectrum::SmoothMarkov(Double_t *source,
int ssize,
int averWindow)
1197 int xmin, xmax, i, l;
1198 Double_t a, b, maxch;
1199 Double_t nom, nip, nim, sp, sm, area = 0;
1201 return "Averaging Window must be positive";
1202 Double_t *working_space =
new Double_t[ssize];
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;
1233 a = TMath::Sqrt(a + nip);
1237 if((i - l + 1) < xmin)
1238 a = source[xmin] / maxch;
1241 a = source[i - l + 1] / maxch;
1246 a = TMath::Sqrt(a + nim);
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;
1459 const char *TSpectrum::Deconvolution(Double_t *source,
const Double_t *response,
1460 int ssize,
int numberIterations,
1461 int numberRepetitions, Double_t boost )
1464 return "Wrong Parameters";
1466 if (numberRepetitions <= 0)
1467 return "Wrong Parameters";
1471 Double_t *working_space =
new Double_t[4 * ssize];
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;
1665 const char *TSpectrum::DeconvolutionRL(Double_t *source,
const Double_t *response,
1666 int ssize,
int numberIterations,
1667 int numberRepetitions, Double_t boost )
1670 return "Wrong Parameters";
1672 if (numberRepetitions <= 0)
1673 return "Wrong Parameters";
1677 Double_t *working_space =
new Double_t[4 * ssize];
1678 int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
1679 Double_t lda, ldb, ldc, maximum;
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;
1889 const char *TSpectrum::Unfolding(Double_t *source,
1890 const Double_t **respMatrix,
1891 int ssizex,
int ssizey,
1892 int numberIterations,
1893 int numberRepetitions, Double_t boost)
1895 int i, j, k, lindex, lhx = 0, repet;
1896 Double_t lda, ldb, ldc, area;
1897 if (ssizex <= 0 || ssizey <= 0)
1898 return "Wrong Parameters";
1899 if (ssizex < ssizey)
1900 return "Sizex must be greater than sizey)";
1901 if (numberIterations <= 0)
1902 return "Number of iterations must be positive";
1903 Double_t *working_space =
1904 new Double_t[ssizex * ssizey + 2 * ssizey * ssizey + 4 * ssizex];
1907 for (j = 0; j < ssizey && lhx != -1; j++) {
1910 for (i = 0; i < ssizex; i++) {
1911 lda = respMatrix[j][i];
1915 working_space[j * ssizex + i] = lda;
1919 for (i = 0; i < ssizex; i++)
1920 working_space[j * ssizex + i] /= area;
1924 delete [] working_space;
1925 return (
"ZERO COLUMN IN RESPONSE MATRIX");
1929 for (i = 0; i < ssizex; i++)
1930 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1934 for (i = 0; i < ssizey; i++) {
1935 for (j = 0; j < ssizey; j++) {
1937 for (k = 0; k < ssizex; k++) {
1938 ldb = working_space[ssizex * i + k];
1939 ldc = working_space[ssizex * j + k];
1940 lda = lda + ldb * ldc;
1942 working_space[ssizex * ssizey + ssizey * i + j] = lda;
1945 for (k = 0; k < ssizex; k++) {
1946 ldb = working_space[ssizex * i + k];
1948 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
1950 lda = lda + ldb * ldc;
1952 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
1957 for (i = 0; i < ssizey; i++)
1958 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1959 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
1962 for (i = 0; i < ssizey; i++) {
1963 for (j = 0; j < ssizey; j++) {
1965 for (k = 0; k < ssizey; k++) {
1966 ldb = working_space[ssizex * ssizey + ssizey * i + k];
1967 ldc = working_space[ssizex * ssizey + ssizey * j + k];
1968 lda = lda + ldb * ldc;
1970 working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j] =
1974 for (k = 0; k < ssizey; k++) {
1975 ldb = working_space[ssizex * ssizey + ssizey * i + k];
1977 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
1979 lda = lda + ldb * ldc;
1981 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
1986 for (i = 0; i < ssizey; i++)
1987 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
1988 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
1991 for (i = 0; i < ssizey; i++)
1992 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = 1;
1995 for (repet = 0; repet < numberRepetitions; repet++) {
1997 for (i = 0; i < ssizey; i++)
1998 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = TMath::Power(working_space[ssizex * ssizey + 2 * ssizey * ssizey + i], boost);
2000 for (lindex = 0; lindex < numberIterations; lindex++) {
2001 for (i = 0; i < ssizey; i++) {
2003 for (j = 0; j < ssizey; j++) {
2005 working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j];
2006 ldc = working_space[ssizex * ssizey + 2 * ssizey * ssizey + j];
2007 lda = lda + ldb * ldc;
2010 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i];
2017 ldb = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2019 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] = lda;
2021 for (i = 0; i < ssizey; i++)
2022 working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
2023 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
2028 for (i = 0; i < ssizex; i++) {
2030 source[i] = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2035 delete[]working_space;
2126 Int_t TSpectrum::SearchHighRes(Double_t *source,Double_t *destVector,
int ssize,
2127 Double_t sigma, Double_t threshold,
2128 bool backgroundRemove,
int deconIterations,
2129 bool markov,
int averWindow)
2131 int i, j, numberIterations = (Int_t)(7 * sigma + 0.5);
2133 int k, lindex, posit, imin, imax, jmin, jmax, lh_gold, priz;
2134 Double_t lda, ldb, ldc, area, maximum, maximum_decon;
2135 int xmin, xmax, l, peak_index = 0, size_ext = ssize + 2 * numberIterations, shift = numberIterations, bw = 2, w;
2137 Double_t nom, nip, nim, sp, sm, plocha = 0;
2138 Double_t m0low=0,m1low=0,m2low=0,l0low=0,l1low=0,detlow,av,men;
2140 Error(
"SearchHighRes",
"Invalid sigma, must be greater than or equal to 1");
2144 if(threshold<=0 || threshold>=100){
2145 Error(
"SearchHighRes",
"Invalid threshold, must be positive and less than 100");
2149 j = (Int_t) (5.0 * sigma + 0.5);
2150 if (j >= PEAK_WINDOW / 2) {
2151 Error(
"SearchHighRes",
"Too large sigma");
2155 if (markov ==
true) {
2156 if (averWindow <= 0) {
2157 Error(
"SearchHighRes",
"Averaging window must be positive");
2162 if(backgroundRemove ==
true){
2163 if(ssize < 2 * numberIterations + 1){
2164 Error(
"SearchHighRes",
"Too large clipping window");
2169 k = int(2 * sigma+0.5);
2171 for(i = 0;i < k;i++){
2172 a = i,b = source[i];
2173 m0low += 1,m1low += a,m2low += a * a,l0low += b,l1low += a * b;
2175 detlow = m0low * m2low - m1low * m1low;
2177 l1low = (-l0low * m1low + l1low * m0low) / detlow;
2189 i = (Int_t)(7 * sigma + 0.5);
2191 Double_t *working_space =
new Double_t [7 * (ssize + i)];
2192 for (j=0;j<7 * (ssize + i);j++) working_space[j] = 0;
2193 for(i = 0; i < size_ext; i++){
2196 working_space[i + size_ext] = source[0] + l1low * a;
2197 if(working_space[i + size_ext] < 0)
2198 working_space[i + size_ext]=0;
2201 else if(i >= ssize + shift){
2202 a = i - (ssize - 1 + shift);
2203 working_space[i + size_ext] = source[ssize - 1];
2204 if(working_space[i + size_ext] < 0)
2205 working_space[i + size_ext]=0;
2209 working_space[i + size_ext] = source[i - shift];
2212 if(backgroundRemove ==
true){
2213 for(i = 1; i <= numberIterations; i++){
2214 for(j = i; j < size_ext - i; j++){
2215 if(markov ==
false){
2216 a = working_space[size_ext + j];
2217 b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2225 a = working_space[size_ext + j];
2228 for (w = j - bw; w <= j + bw; w++){
2229 if ( w >= 0 && w < size_ext){
2230 av += working_space[size_ext + w];
2237 for (w = j - i - bw; w <= j - i + bw; w++){
2238 if ( w >= 0 && w < size_ext){
2239 b += working_space[size_ext + w];
2246 for (w = j + i - bw; w <= j + i + bw; w++){
2247 if ( w >= 0 && w < size_ext){
2248 c += working_space[size_ext + w];
2256 working_space[j]=av;
2259 for(j = i; j < size_ext - i; j++)
2260 working_space[size_ext + j] = working_space[j];
2262 for(j = 0;j < size_ext; j++){
2265 b = source[0] + l1low * a;
2267 working_space[size_ext + j] = b - working_space[size_ext + j];
2270 else if(j >= ssize + shift){
2271 a = j - (ssize - 1 + shift);
2272 b = source[ssize - 1];
2274 working_space[size_ext + j] = b - working_space[size_ext + j];
2278 working_space[size_ext + j] = source[j - shift] - working_space[size_ext + j];
2281 for(j = 0;j < size_ext; j++){
2282 if(working_space[size_ext + j] < 0) working_space[size_ext + j] = 0;
2286 for(i = 0; i < size_ext; i++){
2287 working_space[i + 6*size_ext] = working_space[i + size_ext];
2291 for(j = 0; j < size_ext; j++)
2292 working_space[2 * size_ext + j] = working_space[size_ext + j];
2293 xmin = 0,xmax = size_ext - 1;
2294 for(i = 0, maxch = 0; i < size_ext; i++){
2295 working_space[i] = 0;
2296 if(maxch < working_space[2 * size_ext + i])
2297 maxch = working_space[2 * size_ext + i];
2298 plocha += working_space[2 * size_ext + i];
2301 delete [] working_space;
2306 working_space[xmin] = 1;
2307 for(i = xmin; i < xmax; i++){
2308 nip = working_space[2 * size_ext + i] / maxch;
2309 nim = working_space[2 * size_ext + i + 1] / maxch;
2311 for(l = 1; l <= averWindow; l++){
2313 a = working_space[2 * size_ext + xmax] / maxch;
2316 a = working_space[2 * size_ext + i + l] / maxch;
2323 a = TMath::Sqrt(a + nip);
2328 if((i - l + 1) < xmin)
2329 a = working_space[2 * size_ext + xmin] / maxch;
2332 a = working_space[2 * size_ext + i - l + 1] / maxch;
2339 a = TMath::Sqrt(a + nim);
2346 a = working_space[i + 1] = working_space[i] * a;
2349 for(i = xmin; i <= xmax; i++){
2350 working_space[i] = working_space[i] / nom;
2352 for(j = 0; j < size_ext; j++)
2353 working_space[size_ext + j] = working_space[j] * plocha;
2354 for(j = 0; j < size_ext; j++){
2355 working_space[2 * size_ext + j] = working_space[size_ext + j];
2357 if(backgroundRemove ==
true){
2358 for(i = 1; i <= numberIterations; i++){
2359 for(j = i; j < size_ext - i; j++){
2360 a = working_space[size_ext + j];
2361 b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2364 working_space[j] = a;
2366 for(j = i; j < size_ext - i; j++)
2367 working_space[size_ext + j] = working_space[j];
2369 for(j = 0; j < size_ext; j++){
2370 working_space[size_ext + j] = working_space[2 * size_ext + j] - working_space[size_ext + j];
2380 for(i = 0; i < size_ext; i++){
2381 lda = (Double_t)i - 3 * sigma;
2382 lda = lda * lda / (2 * sigma * sigma);
2383 j = (Int_t)(1000 * TMath::Exp(-lda));
2388 working_space[i] = lda;
2396 for(i = 0; i < size_ext; i++)
2397 working_space[2 * size_ext + i] = TMath::Abs(working_space[size_ext + i]);
2404 for(i = imin; i <= imax; i++){
2409 jmax = lh_gold - 1 - i;
2410 if(jmax > (lh_gold - 1))
2413 for(j = jmin;j <= jmax; j++){
2414 ldb = working_space[j];
2415 ldc = working_space[i + j];
2416 lda = lda + ldb * ldc;
2418 working_space[size_ext + i - imin] = lda;
2422 imin = -i,imax = size_ext + i - 1;
2423 for(i = imin; i <= imax; i++){
2425 for(j = 0; j <= (lh_gold - 1); j++){
2426 ldb = working_space[j];
2428 if(k >= 0 && k < size_ext){
2429 ldc = working_space[2 * size_ext + k];
2430 lda = lda + ldb * ldc;
2434 working_space[4 * size_ext + i - imin] = lda;
2437 for(i = imin; i <= imax; i++)
2438 working_space[2 * size_ext + i - imin] = working_space[4 * size_ext + i - imin];
2440 for(i = 0; i < size_ext; i++)
2441 working_space[i] = 1;
2443 for(lindex = 0; lindex < deconIterations; lindex++){
2444 for(i = 0; i < size_ext; i++){
2445 if(TMath::Abs(working_space[2 * size_ext + i]) > 0.00001 && TMath::Abs(working_space[i]) > 0.00001){
2453 if(jmax > (size_ext - 1 - i))
2456 for(j = jmin; j <= jmax; j++){
2457 ldb = working_space[j + lh_gold - 1 + size_ext];
2458 ldc = working_space[i + j];
2459 lda = lda + ldb * ldc;
2461 ldb = working_space[2 * size_ext + i];
2468 ldb = working_space[i];
2470 working_space[3 * size_ext + i] = lda;
2473 for(i = 0; i < size_ext; i++){
2474 working_space[i] = working_space[3 * size_ext + i];
2478 for(i=0;i<size_ext;i++){
2479 lda = working_space[i];
2482 working_space[size_ext + j] = lda;
2485 maximum = 0, maximum_decon = 0;
2487 for(i = 0; i < size_ext - j; i++){
2488 if(i >= shift && i < ssize + shift){
2489 working_space[i] = area * working_space[size_ext + i + j];
2490 if(maximum_decon < working_space[i])
2491 maximum_decon = working_space[i];
2492 if(maximum < working_space[6 * size_ext + i])
2493 maximum = working_space[6 * size_ext + i];
2497 working_space[i] = 0;
2505 for(i = 1; i < size_ext - 1; i++){
2506 if(working_space[i] > working_space[i - 1] && working_space[i] > working_space[i + 1]){
2507 if(i >= shift && i < ssize + shift){
2508 if(working_space[i] > lda*maximum_decon && working_space[6 * size_ext + i] > threshold * maximum / 100.0){
2509 for(j = i - 1, a = 0, b = 0; j <= i + 1; j++){
2510 a += (Double_t)(j - shift) * working_space[j];
2511 b += working_space[j];
2519 if(peak_index == 0){
2525 for(j = 0, priz = 0; j < peak_index && priz == 0; j++){
2526 if(working_space[6 * size_ext + shift + (Int_t)a] > working_space[6 * size_ext + shift + (Int_t)fPositionX[j]])
2536 for(k = peak_index; k >= j; k--){
2538 fPositionX[k] = fPositionX[k - 1];
2541 fPositionX[j - 1] = a;
2543 if(peak_index < fMaxPeaks)
2551 for(i = 0; i < ssize; i++) destVector[i] = working_space[i + shift];
2552 delete [] working_space;
2553 fNPeaks = peak_index;
2554 if(peak_index == fMaxPeaks)
2555 Warning(
"SearchHighRes",
"Peak buffer full");
2563 Int_t TSpectrum::Search1HighRes(Double_t *source,Double_t *destVector,
int ssize,
2564 Double_t sigma, Double_t threshold,
2565 bool backgroundRemove,
int deconIterations,
2566 bool markov,
int averWindow)
2570 return SearchHighRes(source,destVector,ssize,sigma,threshold,backgroundRemove,
2571 deconIterations,markov,averWindow);
2577 Int_t TSpectrum::StaticSearch(
const TH1 *hist, Double_t sigma, Option_t *option, Double_t threshold)
2581 return s.Search(hist,sigma,option,threshold);
2587 TH1 *TSpectrum::StaticBackground(
const TH1 *hist,Int_t niter, Option_t *option)
2590 return s.Background(hist,niter,option);