// @(#)root/spectrum:$Id$ // Author: Miroslav Morhac 27/05/99 #include "TSpectrum.h" #include "TPolyMarker.h" #include "TVirtualPad.h" #include "TList.h" #include "TH1.h" #include "TMath.h" //______________________________________________________________________________ /* Begin_Html
Author:
Miroslav Morhac
Institute of Physics
Slovak Academy of Sciences
Dubravska cesta 9, 842 28 BRATISLAVA
SLOVAKIA
email:fyzimiro@savba.sk, fax:+421 7 54772479
The original code in C has been repackaged as a C++ class by R.Brun.
The algorithms in this class have been published in the following references:
This function calculates the background spectrum in the input histogram h. The background is returned as a histogram.
Function parameters:
This function searches for peaks in source spectrum in hin The number of found peaks and their positions are written into the members fNpeaks and fPositionX. The search is performed in the current histogram range.
Function parameters:
By default the "Markov" chain algorithm is used. Specify the option "noMarkov" to disable this algorithm Note that by default the source spectrum is replaced by a new spectrum
By default a polymarker object is created and added to the list of functions of the histogram. The histogram is drawn with the specified option and the polymarker object drawn on top of the histogram. The polymarker coordinates correspond to the npeaks peaks found in the histogram.
A pointer to the polymarker object can be retrieved later via:
TList *functions = hin->GetListOfFunctions(); TPolyMarker *pm = (TPolyMarker*)functions->FindObject("TPolyMarker");Specify the option "goff" to disable the storage and drawing of the polymarker.
To disable the final drawing of the histogram with the search results (in case
you want to draw it yourself) specify "nodraw" in the options parameter.
End_Html */
if (hin == 0) return 0;
Int_t dimension = hin->GetDimension();
if (dimension > 2) {
Error("Search", "Only implemented for 1-d and 2-d histograms");
return 0;
}
if (threshold <=0 || threshold >= 1) {
Warning("Search","threshold must 0
Figure 1 Example of the estimation of background for number of iterations=6.
Original spectrum is shown in black color, estimated background in red color.
Script:
In Figure 1. one can notice that at the edges of the peaks the estimated
background goes under the peaks. An alternative approach is to decrease the
clipping window from a given value numberIterations to the value of one, which
is presented in this example.
Figure 2 Example of the estimation of background for numberIterations=6 using
decreasing clipping window algorithm. Original spectrum is shown in black
color, estimated background in red color.
Script:
The question is how to choose the width of the clipping window, i.e.,
numberIterations parameter. The influence of this parameter on the estimated
background is illustrated in Figure 3.
Figure 3 Example of the influence of clipping window width on the estimated
background for numberIterations=4 (red line), 6 (blue line) 8 (green line) using
decreasing clipping window algorithm.
in general one should set this parameter so that the value
2*numberIterations+1 was greater than the widths of preserved objects (peaks).
Script:
another example for very complex spectrum is given in Figure 4.
Figure 4 Example of the influence of clipping window width on the estimated
background for numberIterations=10 (red line), 20 (blue line), 30 (green line)
and 40 (magenta line) using decreasing clipping window algorithm.
Script:
Second order difference filter removes linear (quasi-linear) background and
preserves symmetrical peaks. However if the shape of the background is more
complex one can employ higher-order clipping filters (see example in Figure 5)
Figure 5 Example of the influence of clipping filter difference order on the
estimated background for fNnumberIterations=40, 2-nd order red line, 4-th order
blue line, 6-th order green line and 8-th order magenta line, and using
decreasing clipping window algorithm.
Script:
The estimate of the background can be influenced by noise present in the
spectrum. We proposed the algorithm of the background estimate with
simultaneous smoothing. In the original algorithm without smoothing, the
estimated background snatches the lower spikes in the noise. Consequently,
the areas of peaks are biased by this error.
Figure 7 Principle of background estimation algorithm with simultaneous
smoothing.
Figure 8 Illustration of non-smoothing (red line) and smoothing algorithm of
background estimation (blue line).
Script:
Sometimes it is necessary to include also the Compton edges into the estimate of
the background. In Figure 8 we present the example of the synthetic spectrum
with Compton edges. The background was estimated using the 8-th order filter
with the estimation of the Compton edges using decreasing
clipping window algorithm (numberIterations=10) with smoothing
(smoothingWindow=5).
Figure 8 Example of the estimate of the background with Compton edges (red
line) for numberIterations=10, 8-th order difference filter, using decreasing
clipping window algorithm and smoothing (smoothingWindow=5).
Script:
This function calculates smoothed spectrum from source spectrum based on
Markov chain method. The result is placed in the array pointed by source
pointer. On successful completion it returns 0. On error it returns pointer
to the string describing error.
Function parameters:
being defined
from the normalization condition
.
n is the length of the smoothed spectrum and
Reference:
Example 14 - script Smoothing.c
Fig. 23 Original noisy spectrum
Fig. 24 Smoothed spectrum m=3
Fig. 25 Smoothed spectrum
Fig.26 Smoothed spectrum m=10
Script:
This function calculates deconvolution from source spectrum according to
response spectrum using Gold deconvolution algorithm. The result is placed
in the vector pointed by source pointer. On successful completion it
returns 0. On error it returns pointer to the string describing error. If
desired after every numberIterations one can apply boosting operation
(exponential function with exponent given by boost coefficient) and repeat
it numberRepetitions times.
Function parameters:
where h(i) is the impulse response function, x, y are input and output
vectors, respectively, N is the length of x and h vectors. In matrix form
we have:
Let us assume that we know the response and the output vector (spectrum) of
the above given system. The deconvolution represents solution of the
overdetermined system of linear equations, i.e., the calculation of the
vector x. From numerical stability point of view the operation of
deconvolution is extremely critical (ill-posed problem) as well as time
consuming operation. The Gold deconvolution algorithm proves to work very
well, other methods (Fourier, VanCittert etc) oscillate. It is suitable to
process positive definite data (e.g. histograms).
Gold deconvolution algorithm:
Where L is given number of iterations (numberIterations parameter).
Boosted deconvolution:
References:
Example 8 - script Deconvolution.c :
response function (usually peak) should be shifted left to the first
non-zero channel (bin) (see Figure 9)
Figure 9 Response spectrum.
Figure 10 Principle how the response matrix is composed inside of the
Deconvolution function.
Figure 11 Example of Gold deconvolution. The original source spectrum is
drawn with black color, the spectrum after the deconvolution (10000
iterations) with red color.
Script:
Examples of Gold deconvolution method:
First let us study the influence of the number of iterations on the
deconvolved spectrum (Figure 12).
Figure 12 Study of Gold deconvolution algorithm.The original source spectrum
is drawn with black color, spectrum after 100 iterations with red color,
spectrum after 1000 iterations with blue color, spectrum after 10000
iterations with green color and spectrum after 100000 iterations with
magenta color.
For relatively narrow peaks in the above given example the Gold
deconvolution method is able to decompose overlapping peaks practically to
delta - functions. In the next example we have chosen a synthetic data
(spectrum, 256 channels) consisting of 5 very closely positioned, relatively
wide peaks (sigma =5), with added noise (Figure 13). Thin lines represent
pure Gaussians (see Table 1); thick line is a resulting spectrum with
additive noise (10% of the amplitude of small peaks).
Figure 13 Testing example of synthetic spectrum composed of 5 Gaussians with
added noise.
Table 1 Positions, heights and areas of peaks in the spectrum shown in
Figure 13.
In ideal case, we should obtain the result given in Figure 14. The areas of
the Gaussian components of the spectrum are concentrated completely to
delta-functions. When solving the overdetermined system of linear equations
with data from Figure 13 in the sense of minimum least squares criterion
without any regularization we obtain the result with large oscillations
(Figure 15). From mathematical point of view, it is the optimal solution in
the unconstrained space of independent variables. From physical point of
view we are interested only in a meaningful solution. Therefore, we have to
employ regularization techniques (e.g. Gold deconvolution) and/or to
confine the space of allowed solutions to subspace of positive solutions.
Figure 14 The same spectrum like in Figure 13, outlined bars show the
contents of present components (peaks).
Figure 15 Least squares solution of the system of linear equations without
regularization.
Example 9 - script Deconvolution_wide.c
When we employ Gold deconvolution algorithm we obtain the result given in
Fig. 16. One can observe that the resulting spectrum is smooth. On the
other hand the method is not able to decompose completely the peaks in the
spectrum.
Figure 16 Example of Gold deconvolution for closely positioned wide peaks.
The original source spectrum is drawn with black color, the spectrum after
the deconvolution (10000 iterations) with red color.
Script:
Example 10 - script Deconvolution_wide_boost.c :
Further let us employ boosting operation into deconvolution (Fig. 17).
Figure 17 The original source spectrum is drawn with black color, the
spectrum after the deconvolution with red color. Number of iterations = 200,
number of repetitions = 50 and boosting coefficient = 1.2.
Table 2 Results of the estimation of peaks in spectrum shown in Figure 17.
One can observe that peaks are decomposed practically to delta functions.
Number of peaks is correct, positions of big peaks as well as their areas
are relatively well estimated. However there is a considerable error in
the estimation of the position of small right hand peak.
Script:
This function calculates deconvolution from source spectrum according to
response spectrum using Richardson-Lucy deconvolution algorithm. The result
is placed in the vector pointed by source pointer. On successful completion
it returns 0. On error it returns pointer to the string describing error.
If desired after every numberIterations one can apply boosting operation
(exponential function with exponent given by boost coefficient) and repeat
it numberRepetitions times (see Gold deconvolution).
Function parameters:
Richardson-Lucy deconvolution algorithm:
For discrete systems it has the form:
for positive input data and response matrix this iterative method forces
the deconvoluted spectra to be non-negative. The Richardson-Lucy
iteration converges to the maximum likelihood solution for Poisson statistics
in the data.
References:
Examples of Richardson-Lucy deconvolution method:
Example 11 - script DeconvolutionRL_wide.c :
When we employ Richardson-Lucy deconvolution algorithm to our data from
Fig. 13 we obtain the result given in Fig. 18. One can observe improvements
as compared to the result achieved by Gold deconvolution. Neverthless it is
unable to decompose the multiplet.
Figure 18 Example of Richardson-Lucy deconvolution for closely positioned
wide peaks. The original source spectrum is drawn with black color, the
spectrum after the deconvolution (10000 iterations) with red color.
Script:
Example 12 - script DeconvolutionRL_wide_boost.c :
Further let us employ boosting operation into deconvolution (Fig. 19).
Figure 19 The original source spectrum is drawn with black color, the
spectrum after the deconvolution with red color. Number of iterations = 200,
number of repetitions = 50 and boosting coefficient = 1.2.
Table 3 Results of the estimation of peaks in spectrum shown in Figure 19.
One can observe improvements in the estimation of peak positions as compared
to the results achieved by Gold deconvolution.
Script:
This function unfolds source spectrum according to response matrix columns.
The result is placed in the vector pointed by source pointer.
The coefficients of the resulting vector represent contents of the columns
(weights) in the input vector. On successful completion it returns 0. On
error it returns pointer to the string describing error. If desired after
every numberIterations one can apply boosting operation (exponential
function with exponent given by boost coefficient) and repeat it
numberRepetitions times. For details we refer to [1].
Function parameters:
Unfolding:
The goal is the decomposition of spectrum to a given set of component
spectra.
The mathematical formulation of the discrete linear system is:
References:
Example of unfolding:
Example 13 - script Unfolding.c:
Fig. 20 Response matrix composed of neutron spectra of pure
chemical elements.
Fig. 21 Source neutron spectrum to be decomposed
Fig. 22 Spectrum after decomposition, contains 10 coefficients, which
correspond to contents of chemical components (dominant 8-th and 10-th
components, i.e. O, Si)
Script:
This function searches for peaks in source spectrum. It is based on
deconvolution method. First the background is removed (if desired), then
Markov smoothed spectrum is calculated (if desired), then the response
function is generated according to given sigma and deconvolution is
carried out. The order of peaks is arranged according to their heights in
the spectrum after background elimination. The highest peak is the first in
the list. On success it returns number of found peaks.
Function parameters:
Peaks searching:
The goal of this function is to identify automatically the peaks in spectrum
with the presence of the continuous background and statistical
fluctuations - noise.
The common problems connected with correct peak identification are:
Fig. 27 An example of one-dimensional synthetic spectrum with found peaks
denoted by markers.
References:
Examples of peak searching method:
The SearchHighRes function provides users with the possibility to vary the
input parameters and with the access to the output deconvolved data in the
destination spectrum. Based on the output data one can tune the parameters.
Example 15 - script SearchHR1.c:
Fig. 28 One-dimensional spectrum with found peaks denoted by markers, 3
iterations steps in the deconvolution.
Fig. 29 One-dimensional spectrum with found peaks denoted by markers, 8
iterations steps in the deconvolution.
Script:
Example 16 - script SearchHR3.c:
Table 4 Positions and sigma of peaks in the following examples.
Fig. 30 Influence of number of iterations (3-red, 10-blue, 100- green,
1000-magenta), sigma=8, smoothing width=3.
Fig. 31 Influence of sigma (3-red, 8-blue, 20- green, 43-magenta),
num. iter.=10, sm. width=3.
Fig. 32 Influence smoothing width (0-red, 3-blue, 7- green, 20-magenta), num.
iter.=10, sigma=8.
Script:
where p = 1, 2, ..., numberIterations. In fact it represents second order
difference filter (-1,2,-1).
One can also change the
direction of the change of the clipping window, the order of the clipping
filter, to include smoothing, to set width of smoothing window and to include
the estimation of Compton edges. On successful completion it returns 0. On
error it returns pointer to the string describing error.
Parameters:
References:
Example 1 script Background_incr.c:
// Example to illustrate the background estimator (class TSpectrum).
// To execute this example, do
// root > .x Background_incr.C
#include
Example 2 script Background_decr.c:
// Example to illustrate the background estimator (class TSpectrum).
// To execute this example, do
// root > .x Background_decr.C
#include
Example 3 script Background_width.c:
// Example to illustrate the influence of the clipping window width on the
// estimated background. To execute this example, do:
// root > .x Background_width.C
#include
Example 4 script Background_width2.c:
// Example to illustrate the influence of the clipping window width on the
// estimated background. To execute this example, do:
// root > .x Background_width2.C
#include
Example 5 script Background_order.c:
// Example to illustrate the influence of the clipping filter difference order
// on the estimated background. To execute this example, do
// root > .x Background_order.C
#include
Example 6 script Background_smooth.c:
// Example to illustrate the background estimator (class TSpectrum) including
// Compton edges. To execute this example, do:
// root > .x Background_smooth.C
#include
Example 8 script Background_compton.c:
// Example to illustrate the background estimator (class TSpectrum) including
// Compton edges. To execute this example, do:
// root > .x Background_compton.C
#include
End_Html */
int i, j, w, bw, b1, b2, priz;
float 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;
if (ssize <= 0)
return "Wrong Parameters";
if (numberIterations < 1)
return "Width of Clipping Window Must Be Positive";
if (ssize < 2 * numberIterations + 1)
return "Too Large Clipping Window";
if (smoothing == kTRUE && smoothWindow != kBackSmoothing3 && smoothWindow != kBackSmoothing5 && smoothWindow != kBackSmoothing7 && smoothWindow != kBackSmoothing9 && smoothWindow != kBackSmoothing11 && smoothWindow != kBackSmoothing13 && smoothWindow != kBackSmoothing15)
return "Incorrect width of smoothing window";
float *working_space = new float[2 * ssize];
for (i = 0; i < ssize; i++){
working_space[i] = spectrum[i];
working_space[i + ssize] = spectrum[i];
}
bw=(smoothWindow-1)/2;
if (direction == kBackIncreasingWindow)
i = 1;
else if(direction == kBackDecreasingWindow)
i = numberIterations;
if (filterOrder == kBackOrder2) {
do{
for (j = i; j < ssize - i; j++) {
if (smoothing == kFALSE){
a = working_space[ssize + j];
b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
if (b < a)
a = b;
working_space[j] = a;
}
else if (smoothing == kTRUE){
a = working_space[ssize + j];
av = 0;
men = 0;
for (w = j - bw; w <= j + bw; w++){
if ( w >= 0 && w < ssize){
av += working_space[ssize + w];
men +=1;
}
}
av = av / men;
b = 0;
men = 0;
for (w = j - i - bw; w <= j - i + bw; w++){
if ( w >= 0 && w < ssize){
b += working_space[ssize + w];
men +=1;
}
}
b = b / men;
c = 0;
men = 0;
for (w = j + i - bw; w <= j + i + bw; w++){
if ( w >= 0 && w < ssize){
c += working_space[ssize + w];
men +=1;
}
}
c = c / men;
b = (b + c) / 2;
if (b < a)
av = b;
working_space[j]=av;
}
}
for (j = i; j < ssize - i; j++)
working_space[ssize + j] = working_space[j];
if (direction == kBackIncreasingWindow)
i+=1;
else if(direction == kBackDecreasingWindow)
i-=1;
}while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
}
else if (filterOrder == kBackOrder4) {
do{
for (j = i; j < ssize - i; j++) {
if (smoothing == kFALSE){
a = working_space[ssize + j];
b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
c = 0;
ai = i / 2;
c -= working_space[ssize + j - (int) (2 * ai)] / 6;
c += 4 * working_space[ssize + j - (int) ai] / 6;
c += 4 * working_space[ssize + j + (int) ai] / 6;
c -= working_space[ssize + j + (int) (2 * ai)] / 6;
if (b < c)
b = c;
if (b < a)
a = b;
working_space[j] = a;
}
else if (smoothing == kTRUE){
a = working_space[ssize + j];
av = 0;
men = 0;
for (w = j - bw; w <= j + bw; w++){
if ( w >= 0 && w < ssize){
av += working_space[ssize + w];
men +=1;
}
}
av = av / men;
b = 0;
men = 0;
for (w = j - i - bw; w <= j - i + bw; w++){
if ( w >= 0 && w < ssize){
b += working_space[ssize + w];
men +=1;
}
}
b = b / men;
c = 0;
men = 0;
for (w = j + i - bw; w <= j + i + bw; w++){
if ( w >= 0 && w < ssize){
c += working_space[ssize + w];
men +=1;
}
}
c = c / men;
b = (b + c) / 2;
ai = i / 2;
b4 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b4 += working_space[ssize + w];
men +=1;
}
}
b4 = b4 / men;
c4 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
c4 += working_space[ssize + w];
men +=1;
}
}
c4 = c4 / men;
d4 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
d4 += working_space[ssize + w];
men +=1;
}
}
d4 = d4 / men;
e4 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
e4 += working_space[ssize + w];
men +=1;
}
}
e4 = e4 / men;
b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
if (b < b4)
b = b4;
if (b < a)
av = b;
working_space[j]=av;
}
}
for (j = i; j < ssize - i; j++)
working_space[ssize + j] = working_space[j];
if (direction == kBackIncreasingWindow)
i+=1;
else if(direction == kBackDecreasingWindow)
i-=1;
}while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
}
else if (filterOrder == kBackOrder6) {
do{
for (j = i; j < ssize - i; j++) {
if (smoothing == kFALSE){
a = working_space[ssize + j];
b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
c = 0;
ai = i / 2;
c -= working_space[ssize + j - (int) (2 * ai)] / 6;
c += 4 * working_space[ssize + j - (int) ai] / 6;
c += 4 * working_space[ssize + j + (int) ai] / 6;
c -= working_space[ssize + j + (int) (2 * ai)] / 6;
d = 0;
ai = i / 3;
d += working_space[ssize + j - (int) (3 * ai)] / 20;
d -= 6 * working_space[ssize + j - (int) (2 * ai)] / 20;
d += 15 * working_space[ssize + j - (int) ai] / 20;
d += 15 * working_space[ssize + j + (int) ai] / 20;
d -= 6 * working_space[ssize + j + (int) (2 * ai)] / 20;
d += working_space[ssize + j + (int) (3 * ai)] / 20;
if (b < d)
b = d;
if (b < c)
b = c;
if (b < a)
a = b;
working_space[j] = a;
}
else if (smoothing == kTRUE){
a = working_space[ssize + j];
av = 0;
men = 0;
for (w = j - bw; w <= j + bw; w++){
if ( w >= 0 && w < ssize){
av += working_space[ssize + w];
men +=1;
}
}
av = av / men;
b = 0;
men = 0;
for (w = j - i - bw; w <= j - i + bw; w++){
if ( w >= 0 && w < ssize){
b += working_space[ssize + w];
men +=1;
}
}
b = b / men;
c = 0;
men = 0;
for (w = j + i - bw; w <= j + i + bw; w++){
if ( w >= 0 && w < ssize){
c += working_space[ssize + w];
men +=1;
}
}
c = c / men;
b = (b + c) / 2;
ai = i / 2;
b4 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b4 += working_space[ssize + w];
men +=1;
}
}
b4 = b4 / men;
c4 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
c4 += working_space[ssize + w];
men +=1;
}
}
c4 = c4 / men;
d4 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
d4 += working_space[ssize + w];
men +=1;
}
}
d4 = d4 / men;
e4 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
e4 += working_space[ssize + w];
men +=1;
}
}
e4 = e4 / men;
b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
ai = i / 3;
b6 = 0, men = 0;
for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b6 += working_space[ssize + w];
men +=1;
}
}
b6 = b6 / men;
c6 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
c6 += working_space[ssize + w];
men +=1;
}
}
c6 = c6 / men;
d6 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
d6 += working_space[ssize + w];
men +=1;
}
}
d6 = d6 / men;
e6 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
e6 += working_space[ssize + w];
men +=1;
}
}
e6 = e6 / men;
f6 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
f6 += working_space[ssize + w];
men +=1;
}
}
f6 = f6 / men;
g6 = 0, men = 0;
for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
g6 += working_space[ssize + w];
men +=1;
}
}
g6 = g6 / men;
b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
if (b < b6)
b = b6;
if (b < b4)
b = b4;
if (b < a)
av = b;
working_space[j]=av;
}
}
for (j = i; j < ssize - i; j++)
working_space[ssize + j] = working_space[j];
if (direction == kBackIncreasingWindow)
i+=1;
else if(direction == kBackDecreasingWindow)
i-=1;
}while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
}
else if (filterOrder == kBackOrder8) {
do{
for (j = i; j < ssize - i; j++) {
if (smoothing == kFALSE){
a = working_space[ssize + j];
b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
c = 0;
ai = i / 2;
c -= working_space[ssize + j - (int) (2 * ai)] / 6;
c += 4 * working_space[ssize + j - (int) ai] / 6;
c += 4 * working_space[ssize + j + (int) ai] / 6;
c -= working_space[ssize + j + (int) (2 * ai)] / 6;
d = 0;
ai = i / 3;
d += working_space[ssize + j - (int) (3 * ai)] / 20;
d -= 6 * working_space[ssize + j - (int) (2 * ai)] / 20;
d += 15 * working_space[ssize + j - (int) ai] / 20;
d += 15 * working_space[ssize + j + (int) ai] / 20;
d -= 6 * working_space[ssize + j + (int) (2 * ai)] / 20;
d += working_space[ssize + j + (int) (3 * ai)] / 20;
e = 0;
ai = i / 4;
e -= working_space[ssize + j - (int) (4 * ai)] / 70;
e += 8 * working_space[ssize + j - (int) (3 * ai)] / 70;
e -= 28 * working_space[ssize + j - (int) (2 * ai)] / 70;
e += 56 * working_space[ssize + j - (int) ai] / 70;
e += 56 * working_space[ssize + j + (int) ai] / 70;
e -= 28 * working_space[ssize + j + (int) (2 * ai)] / 70;
e += 8 * working_space[ssize + j + (int) (3 * ai)] / 70;
e -= working_space[ssize + j + (int) (4 * ai)] / 70;
if (b < e)
b = e;
if (b < d)
b = d;
if (b < c)
b = c;
if (b < a)
a = b;
working_space[j] = a;
}
else if (smoothing == kTRUE){
a = working_space[ssize + j];
av = 0;
men = 0;
for (w = j - bw; w <= j + bw; w++){
if ( w >= 0 && w < ssize){
av += working_space[ssize + w];
men +=1;
}
}
av = av / men;
b = 0;
men = 0;
for (w = j - i - bw; w <= j - i + bw; w++){
if ( w >= 0 && w < ssize){
b += working_space[ssize + w];
men +=1;
}
}
b = b / men;
c = 0;
men = 0;
for (w = j + i - bw; w <= j + i + bw; w++){
if ( w >= 0 && w < ssize){
c += working_space[ssize + w];
men +=1;
}
}
c = c / men;
b = (b + c) / 2;
ai = i / 2;
b4 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b4 += working_space[ssize + w];
men +=1;
}
}
b4 = b4 / men;
c4 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
c4 += working_space[ssize + w];
men +=1;
}
}
c4 = c4 / men;
d4 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
d4 += working_space[ssize + w];
men +=1;
}
}
d4 = d4 / men;
e4 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
e4 += working_space[ssize + w];
men +=1;
}
}
e4 = e4 / men;
b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
ai = i / 3;
b6 = 0, men = 0;
for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b6 += working_space[ssize + w];
men +=1;
}
}
b6 = b6 / men;
c6 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
c6 += working_space[ssize + w];
men +=1;
}
}
c6 = c6 / men;
d6 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
d6 += working_space[ssize + w];
men +=1;
}
}
d6 = d6 / men;
e6 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
e6 += working_space[ssize + w];
men +=1;
}
}
e6 = e6 / men;
f6 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
f6 += working_space[ssize + w];
men +=1;
}
}
f6 = f6 / men;
g6 = 0, men = 0;
for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
g6 += working_space[ssize + w];
men +=1;
}
}
g6 = g6 / men;
b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
ai = i / 4;
b8 = 0, men = 0;
for (w = j - (int)(4 * ai) - bw; w <= j - (int)(4 * ai) + bw; w++){
if (w >= 0 && w < ssize){
b8 += working_space[ssize + w];
men +=1;
}
}
b8 = b8 / men;
c8 = 0, men = 0;
for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
c8 += working_space[ssize + w];
men +=1;
}
}
c8 = c8 / men;
d8 = 0, men = 0;
for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
d8 += working_space[ssize + w];
men +=1;
}
}
d8 = d8 / men;
e8 = 0, men = 0;
for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
if (w >= 0 && w < ssize){
e8 += working_space[ssize + w];
men +=1;
}
}
e8 = e8 / men;
f8 = 0, men = 0;
for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
if (w >= 0 && w < ssize){
f8 += working_space[ssize + w];
men +=1;
}
}
f8 = f8 / men;
g8 = 0, men = 0;
for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
if (w >= 0 && w < ssize){
g8 += working_space[ssize + w];
men +=1;
}
}
g8 = g8 / men;
h8 = 0, men = 0;
for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
if (w >= 0 && w < ssize){
h8 += working_space[ssize + w];
men +=1;
}
}
h8 = h8 / men;
i8 = 0, men = 0;
for (w = j + (int)(4 * ai) - bw; w <= j + (int)(4 * ai) + bw; w++){
if (w >= 0 && w < ssize){
i8 += working_space[ssize + w];
men +=1;
}
}
i8 = i8 / men;
b8 = ( -b8 + 8 * c8 - 28 * d8 + 56 * e8 - 56 * f8 - 28 * g8 + 8 * h8 - i8)/70;
if (b < b8)
b = b8;
if (b < b6)
b = b6;
if (b < b4)
b = b4;
if (b < a)
av = b;
working_space[j]=av;
}
}
for (j = i; j < ssize - i; j++)
working_space[ssize + j] = working_space[j];
if (direction == kBackIncreasingWindow)
i += 1;
else if(direction == kBackDecreasingWindow)
i -= 1;
}while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
}
if (compton == kTRUE) {
for (i = 0, b2 = 0; i < ssize; i++){
b1 = b2;
a = working_space[i], b = spectrum[i];
j = i;
if (TMath::Abs(a - b) >= 1) {
b1 = i - 1;
if (b1 < 0)
b1 = 0;
yb1 = working_space[b1];
for (b2 = b1 + 1, c = 0, priz = 0; priz == 0 && b2 < ssize; b2++){
a = working_space[b2], b = spectrum[b2];
c = c + b - yb1;
if (TMath::Abs(a - b) < 1) {
priz = 1;
yb2 = b;
}
}
if (b2 == ssize)
b2 -= 1;
yb2 = working_space[b2];
if (yb1 <= yb2){
for (j = b1, c = 0; j <= b2; j++){
b = spectrum[j];
c = c + b - yb1;
}
if (c > 1){
c = (yb2 - yb1) / c;
for (j = b1, d = 0; j <= b2 && j < ssize; j++){
b = spectrum[j];
d = d + b - yb1;
a = c * d + yb1;
working_space[ssize + j] = a;
}
}
}
else{
for (j = b2, c = 0; j >= b1; j--){
b = spectrum[j];
c = c + b - yb2;
}
if (c > 1){
c = (yb1 - yb2) / c;
for (j = b2, d = 0;j >= b1 && j >= 0; j--){
b = spectrum[j];
d = d + b - yb2;
a = c * d + yb2;
working_space[ssize + j] = a;
}
}
}
i=b2;
}
}
}
for (j = 0; j < ssize; j++)
spectrum[j] = working_space[ssize + j];
delete[]working_space;
return 0;
}
//______________________________________________________________________________
const char* TSpectrum::SmoothMarkov(float *source, int ssize, int averWindow)
{
/* Begin_Html
One-dimensional markov spectrum smoothing function
The goal of this function is the suppression of the statistical fluctuations.
The algorithm is based on discrete Markov chain, which has very simple
invariant distribution:
// Example to illustrate smoothing using Markov algorithm (class TSpectrum).
// To execute this example, do
// root > .x Smoothing.C
void Smoothing() {
Int_t i;
Double_t nbins = 1024;
Double_t xmin = 0;
Double_t xmax = (Double_t)nbins;
Float_t * source = new float[nbins];
TH1F *h = new TH1F("h","Smoothed spectrum for m=3",nbins,xmin,xmax);
TFile *f = new TFile("spectra\\TSpectrum.root");
h=(TH1F*) f->Get("smooth1;1");
for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
TCanvas *Smooth1 = gROOT->GetListOfCanvases()->FindObject("Smooth1");
if (!Smooth1) Smooth1 = new TCanvas("Smooth1","Smooth1",10,10,1000,700);
TSpectrum *s = new TSpectrum();
s->SmoothMarkov(source,1024,3); //3, 7, 10
for (i = 0; i < nbins; i++) h->SetBinContent(i + 1,source[i]);
h->SetAxisRange(330,880);
h->Draw("L");
}
End_Html */
int xmin, xmax, i, l;
float a, b, maxch;
float nom, nip, nim, sp, sm, area = 0;
if(averWindow <= 0)
return "Averaging Window must be positive";
float *working_space = new float[ssize];
xmin = 0,xmax = ssize - 1;
for(i = 0, maxch = 0; i < ssize; i++){
working_space[i]=0;
if(maxch < source[i])
maxch = source[i];
area += source[i];
}
if(maxch == 0) {
delete [] working_space;
return 0 ;
}
nom = 1;
working_space[xmin] = 1;
for(i = xmin; i < xmax; i++){
nip = source[i] / maxch;
nim = source[i + 1] / maxch;
sp = 0,sm = 0;
for(l = 1; l <= averWindow; l++){
if((i + l) > xmax)
a = source[xmax] / maxch;
else
a = source[i + l] / maxch;
b = a - nip;
if(a + nip <= 0)
a = 1;
else
a = TMath::Sqrt(a + nip);
b = b / a;
b = TMath::Exp(b);
sp = sp + b;
if((i - l + 1) < xmin)
a = source[xmin] / maxch;
else
a = source[i - l + 1] / maxch;
b = a - nim;
if(a + nim <= 0)
a = 1;
else
a = TMath::Sqrt(a + nim);
b = b / a;
b = TMath::Exp(b);
sm = sm + b;
}
a = sp / sm;
a = working_space[i + 1] = working_space[i] * a;
nom = nom + a;
}
for(i = xmin; i <= xmax; i++){
working_space[i] = working_space[i] / nom;
}
for(i = 0; i < ssize; i++)
source[i] = working_space[i] * area;
delete[]working_space;
return 0;
}
//______________________________________________________________________________
const char *TSpectrum::Deconvolution(float *source, const float *response,
int ssize, int numberIterations,
int numberRepetitions, double boost )
{
/* Begin_Html
One-dimensional deconvolution function
The goal of this function is the improvement of the resolution in spectra,
decomposition of multiplets. The mathematical formulation of
the convolution system is:
// Example to illustrate deconvolution function (class TSpectrum).
// To execute this example, do
// root > .x Deconvolution.C
#include
Peak # Position Height Area
1 50 500 10159
2 70 3000 60957
3 80 1000 20319
4 100 5000 101596
5 110 500 10159
// Example to illustrate deconvolution function (class TSpectrum).
// To execute this example, do
// root > .x Deconvolution_wide.C
#include
Peak # Original/Estimated (max) position Original/Estimated area
1 50/49 10159/10419
2 70/70 60957/58933
3 80/79 20319/19935
4 100/100 101596/105413
5 110/117 10159/6676
// Example to illustrate deconvolution function (class TSpectrum).
// To execute this example, do
// root > .x Deconvolution_wide_boost.C
#include
End_Html */
if (ssize <= 0)
return "Wrong Parameters";
if (numberRepetitions <= 0)
return "Wrong Parameters";
// working_space-pointer to the working vector
// (its size must be 4*ssize of source spectrum)
double *working_space = new double[4 * ssize];
int i, j, k, lindex, posit, lh_gold, l, repet;
double lda, ldb, ldc, area, maximum;
area = 0;
lh_gold = -1;
posit = 0;
maximum = 0;
//read response vector
for (i = 0; i < ssize; i++) {
lda = response[i];
if (lda != 0)
lh_gold = i + 1;
working_space[i] = lda;
area += lda;
if (lda > maximum) {
maximum = lda;
posit = i;
}
}
if (lh_gold == -1) {
delete [] working_space;
return "ZERO RESPONSE VECTOR";
}
//read source vector
for (i = 0; i < ssize; i++)
working_space[2 * ssize + i] = source[i];
// create matrix at*a and vector at*y
for (i = 0; i < ssize; i++){
lda = 0;
for (j = 0; j < ssize; j++){
ldb = working_space[j];
k = i + j;
if (k < ssize){
ldc = working_space[k];
lda = lda + ldb * ldc;
}
}
working_space[ssize + i] = lda;
lda = 0;
for (k = 0; k < ssize; k++){
l = k - i;
if (l >= 0){
ldb = working_space[l];
ldc = working_space[2 * ssize + k];
lda = lda + ldb * ldc;
}
}
working_space[3 * ssize + i]=lda;
}
// move vector at*y
for (i = 0; i < ssize; i++){
working_space[2 * ssize + i] = working_space[3 * ssize + i];
}
//initialization of resulting vector
for (i = 0; i < ssize; i++)
working_space[i] = 1;
//**START OF ITERATIONS**
for (repet = 0; repet < numberRepetitions; repet++) {
if (repet != 0) {
for (i = 0; i < ssize; i++)
working_space[i] = TMath::Power(working_space[i], boost);
}
for (lindex = 0; lindex < numberIterations; lindex++) {
for (i = 0; i < ssize; i++) {
if (working_space[2 * ssize + i] > 0.000001
&& working_space[i] > 0.000001) {
lda = 0;
for (j = 0; j < lh_gold; j++) {
ldb = working_space[j + ssize];
if (j != 0){
k = i + j;
ldc = 0;
if (k < ssize)
ldc = working_space[k];
k = i - j;
if (k >= 0)
ldc += working_space[k];
}
else
ldc = working_space[i];
lda = lda + ldb * ldc;
}
ldb = working_space[2 * ssize + i];
if (lda != 0)
lda = ldb / lda;
else
lda = 0;
ldb = working_space[i];
lda = lda * ldb;
working_space[3 * ssize + i] = lda;
}
}
for (i = 0; i < ssize; i++)
working_space[i] = working_space[3 * ssize + i];
}
}
//shift resulting spectrum
for (i = 0; i < ssize; i++) {
lda = working_space[i];
j = i + posit;
j = j % ssize;
working_space[ssize + j] = lda;
}
//write back resulting spectrum
for (i = 0; i < ssize; i++)
source[i] = area * working_space[ssize + i];
delete[]working_space;
return 0;
}
//______________________________________________________________________________
const char *TSpectrum::DeconvolutionRL(float *source, const float *response,
int ssize, int numberIterations,
int numberRepetitions, double boost )
{
/* Begin_Html
One-dimensional deconvolution function.
// Example to illustrate deconvolution function (class TSpectrum).
// To execute this example, do
// root > .x DeconvolutionRL_wide.C
#include
Peak # Original/Estimated (max) position Original/Estimated area 1 50/51 10159/11426 2 70/71 60957/65003 3 80/81 20319/12813 4 100/100 101596/101851 5 110/111 10159/8920
// Example to illustrate deconvolution function (class TSpectrum).
// To execute this example, do
// root > .x DeconvolutionRL_wide_boost.C
#include
End_Html */
if (ssize <= 0)
return "Wrong Parameters";
if (numberRepetitions <= 0)
return "Wrong Parameters";
// working_space-pointer to the working vector
// (its size must be 4*ssize of source spectrum)
double *working_space = new double[4 * ssize];
int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
double lda, ldb, ldc, maximum;
lh_gold = -1;
posit = 0;
maximum = 0;
//read response vector
for (i = 0; i < ssize; i++) {
lda = response[i];
if (lda != 0)
lh_gold = i + 1;
working_space[ssize + i] = lda;
if (lda > maximum) {
maximum = lda;
posit = i;
}
}
if (lh_gold == -1) {
delete [] working_space;
return "ZERO RESPONSE VECTOR";
}
//read source vector
for (i = 0; i < ssize; i++)
working_space[2 * ssize + i] = source[i];
//initialization of resulting vector
for (i = 0; i < ssize; i++){
if (i <= ssize - lh_gold)
working_space[i] = 1;
else
working_space[i] = 0;
}
//**START OF ITERATIONS**
for (repet = 0; repet < numberRepetitions; repet++) {
if (repet != 0) {
for (i = 0; i < ssize; i++)
working_space[i] = TMath::Power(working_space[i], boost);
}
for (lindex = 0; lindex < numberIterations; lindex++) {
for (i = 0; i <= ssize - lh_gold; i++){
lda = 0;
if (working_space[i] > 0){//x[i]
for (j = i; j < i + lh_gold; j++){
ldb = working_space[2 * ssize + j];//y[j]
if (j < ssize){
if (ldb > 0){//y[j]
kmax = j;
if (kmax > lh_gold - 1)
kmax = lh_gold - 1;
kmin = j + lh_gold - ssize;
if (kmin < 0)
kmin = 0;
ldc = 0;
for (k = kmax; k >= kmin; k--){
ldc += working_space[ssize + k] * working_space[j - k];//h[k]*x[j-k]
}
if (ldc > 0)
ldb = ldb / ldc;
else
ldb = 0;
}
ldb = ldb * working_space[ssize + j - i];//y[j]*h[j-i]/suma(h[j][k]x[k])
}
lda += ldb;
}
lda = lda * working_space[i];
}
working_space[3 * ssize + i] = lda;
}
for (i = 0; i < ssize; i++)
working_space[i] = working_space[3 * ssize + i];
}
}
//shift resulting spectrum
for (i = 0; i < ssize; i++) {
lda = working_space[i];
j = i + posit;
j = j % ssize;
working_space[ssize + j] = lda;
}
//write back resulting spectrum
for (i = 0; i < ssize; i++)
source[i] = working_space[ssize + i];
delete[]working_space;
return 0;
}
//______________________________________________________________________________
const char *TSpectrum::Unfolding(float *source,
const float **respMatrix,
int ssizex, int ssizey,
int numberIterations,
int numberRepetitions, double boost)
{
/* Begin_Html
One-dimensional unfolding function
// Example to illustrate unfolding function (class TSpectrum).
// To execute this example, do
// root > .x Unfolding.C
#include
End_Html */
int i, j, k, lindex, lhx = 0, repet;
double lda, ldb, ldc, area;
if (ssizex <= 0 || ssizey <= 0)
return "Wrong Parameters";
if (ssizex < ssizey)
return "Sizex must be greater than sizey)";
if (numberIterations <= 0)
return "Number of iterations must be positive";
double *working_space =
new double[ssizex * ssizey + 2 * ssizey * ssizey + 4 * ssizex];
/*read response matrix*/
for (j = 0; j < ssizey && lhx != -1; j++) {
area = 0;
lhx = -1;
for (i = 0; i < ssizex; i++) {
lda = respMatrix[j][i];
if (lda != 0) {
lhx = i + 1;
}
working_space[j * ssizex + i] = lda;
area = area + lda;
}
if (lhx != -1) {
for (i = 0; i < ssizex; i++)
working_space[j * ssizex + i] /= area;
}
}
if (lhx == -1) {
delete [] working_space;
return ("ZERO COLUMN IN RESPONSE MATRIX");
}
/*read source vector*/
for (i = 0; i < ssizex; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
source[i];
/*create matrix at*a + at*y */
for (i = 0; i < ssizey; i++) {
for (j = 0; j < ssizey; j++) {
lda = 0;
for (k = 0; k < ssizex; k++) {
ldb = working_space[ssizex * i + k];
ldc = working_space[ssizex * j + k];
lda = lda + ldb * ldc;
}
working_space[ssizex * ssizey + ssizey * i + j] = lda;
}
lda = 0;
for (k = 0; k < ssizex; k++) {
ldb = working_space[ssizex * i + k];
ldc =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
k];
lda = lda + ldb * ldc;
}
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
lda;
}
/*move vector at*y*/
for (i = 0; i < ssizey; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
/*create matrix at*a*at*a + vector at*a*at*y */
for (i = 0; i < ssizey; i++) {
for (j = 0; j < ssizey; j++) {
lda = 0;
for (k = 0; k < ssizey; k++) {
ldb = working_space[ssizex * ssizey + ssizey * i + k];
ldc = working_space[ssizex * ssizey + ssizey * j + k];
lda = lda + ldb * ldc;
}
working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j] =
lda;
}
lda = 0;
for (k = 0; k < ssizey; k++) {
ldb = working_space[ssizex * ssizey + ssizey * i + k];
ldc =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
k];
lda = lda + ldb * ldc;
}
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
lda;
}
/*move at*a*at*y*/
for (i = 0; i < ssizey; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
/*initialization in resulting vector */
for (i = 0; i < ssizey; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = 1;
/***START OF ITERATIONS***/
for (repet = 0; repet < numberRepetitions; repet++) {
if (repet != 0) {
for (i = 0; i < ssizey; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = TMath::Power(working_space[ssizex * ssizey + 2 * ssizey * ssizey + i], boost);
}
for (lindex = 0; lindex < numberIterations; lindex++) {
for (i = 0; i < ssizey; i++) {
lda = 0;
for (j = 0; j < ssizey; j++) {
ldb =
working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j];
ldc = working_space[ssizex * ssizey + 2 * ssizey * ssizey + j];
lda = lda + ldb * ldc;
}
ldb =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i];
if (lda != 0) {
lda = ldb / lda;
}
else
lda = 0;
ldb = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
lda = lda * ldb;
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] = lda;
}
for (i = 0; i < ssizey; i++)
working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
}
}
/*write back resulting spectrum*/
for (i = 0; i < ssizex; i++) {
if (i < ssizey)
source[i] = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
else
source[i] = 0;
}
delete[]working_space;
return 0;
}
//______________________________________________________________________________
Int_t TSpectrum::SearchHighRes(float *source,float *destVector, int ssize,
float sigma, double threshold,
bool backgroundRemove,int deconIterations,
bool markov, int averWindow)
{
/* Begin_Html
One-dimensional high-resolution peak search function
// Example to illustrate high resolution peak searching function (class TSpectrum).
// To execute this example, do
// root > .x SearchHR1.C
#include
Peak # Position Sigma 1 118 26 2 162 41 3 310 4 4 330 8 5 482 22 6 491 26 7 740 21 8 852 15 9 954 12 10 989 13
// Example to illustrate the influence of number of iterations in deconvolution in high resolution peak searching function (class TSpectrum).
// To execute this example, do
// root > .x SearchHR3.C
#include
End_Html */
int i, j, numberIterations = (int)(7 * sigma + 0.5);
double a, b, c;
int k, lindex, posit, imin, imax, jmin, jmax, lh_gold, priz;
double lda, ldb, ldc, area, maximum, maximum_decon;
int xmin, xmax, l, peak_index = 0, size_ext = ssize + 2 * numberIterations, shift = numberIterations, bw = 2, w;
double maxch;
double nom, nip, nim, sp, sm, plocha = 0;
double m0low=0,m1low=0,m2low=0,l0low=0,l1low=0,detlow,av,men;
if (sigma < 1) {
Error("SearchHighRes", "Invalid sigma, must be greater than or equal to 1");
return 0;
}
if(threshold<=0 || threshold>=100){
Error("SearchHighRes", "Invalid threshold, must be positive and less than 100");
return 0;
}
j = (int) (5.0 * sigma + 0.5);
if (j >= PEAK_WINDOW / 2) {
Error("SearchHighRes", "Too large sigma");
return 0;
}
if (markov == true) {
if (averWindow <= 0) {
Error("SearchHighRes", "Averanging window must be positive");
return 0;
}
}
if(backgroundRemove == true){
if(ssize < 2 * numberIterations + 1){
Error("SearchHighRes", "Too large clipping window");
return 0;
}
}
k = int(2 * sigma+0.5);
if(k >= 2){
for(i = 0;i < k;i++){
a = i,b = source[i];
m0low += 1,m1low += a,m2low += a * a,l0low += b,l1low += a * b;
}
detlow = m0low * m2low - m1low * m1low;
if(detlow != 0)
l1low = (-l0low * m1low + l1low * m0low) / detlow;
else
l1low = 0;
if(l1low > 0)
l1low=0;
}
else{
l1low = 0;
}
i = (int)(7 * sigma + 0.5);
i = 2 * i;
double *working_space = new double [7 * (ssize + i)];
for (j=0;j<7 * (ssize + i);j++) working_space[j] = 0;
for(i = 0; i < size_ext; i++){
if(i < shift){
a = i - shift;
working_space[i + size_ext] = source[0] + l1low * a;
if(working_space[i + size_ext] < 0)
working_space[i + size_ext]=0;
}
else if(i >= ssize + shift){
a = i - (ssize - 1 + shift);
working_space[i + size_ext] = source[ssize - 1];
if(working_space[i + size_ext] < 0)
working_space[i + size_ext]=0;
}
else
working_space[i + size_ext] = source[i - shift];
}
if(backgroundRemove == true){
for(i = 1; i <= numberIterations; i++){
for(j = i; j < size_ext - i; j++){
if(markov == false){
a = working_space[size_ext + j];
b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
if(b < a)
a = b;
working_space[j]=a;
}
else{
a = working_space[size_ext + j];
av = 0;
men = 0;
for (w = j - bw; w <= j + bw; w++){
if ( w >= 0 && w < size_ext){
av += working_space[size_ext + w];
men +=1;
}
}
av = av / men;
b = 0;
men = 0;
for (w = j - i - bw; w <= j - i + bw; w++){
if ( w >= 0 && w < size_ext){
b += working_space[size_ext + w];
men +=1;
}
}
b = b / men;
c = 0;
men = 0;
for (w = j + i - bw; w <= j + i + bw; w++){
if ( w >= 0 && w < size_ext){
c += working_space[size_ext + w];
men +=1;
}
}
c = c / men;
b = (b + c) / 2;
if (b < a)
av = b;
working_space[j]=av;
}
}
for(j = i; j < size_ext - i; j++)
working_space[size_ext + j] = working_space[j];
}
for(j = 0;j < size_ext; j++){
if(j < shift){
a = j - shift;
b = source[0] + l1low * a;
if (b < 0) b = 0;
working_space[size_ext + j] = b - working_space[size_ext + j];
}
else if(j >= ssize + shift){
a = j - (ssize - 1 + shift);
b = source[ssize - 1];
if (b < 0) b = 0;
working_space[size_ext + j] = b - working_space[size_ext + j];
}
else{
working_space[size_ext + j] = source[j - shift] - working_space[size_ext + j];
}
}
for(j = 0;j < size_ext; j++){
if(working_space[size_ext + j] < 0) working_space[size_ext + j] = 0;
}
}
for(i = 0; i < size_ext; i++){
working_space[i + 6*size_ext] = working_space[i + size_ext];
}
if(markov == true){
for(j = 0; j < size_ext; j++)
working_space[2 * size_ext + j] = working_space[size_ext + j];
xmin = 0,xmax = size_ext - 1;
for(i = 0, maxch = 0; i < size_ext; i++){
working_space[i] = 0;
if(maxch < working_space[2 * size_ext + i])
maxch = working_space[2 * size_ext + i];
plocha += working_space[2 * size_ext + i];
}
if(maxch == 0) {
delete [] working_space;
return 0;
}
nom = 1;
working_space[xmin] = 1;
for(i = xmin; i < xmax; i++){
nip = working_space[2 * size_ext + i] / maxch;
nim = working_space[2 * size_ext + i + 1] / maxch;
sp = 0,sm = 0;
for(l = 1; l <= averWindow; l++){
if((i + l) > xmax)
a = working_space[2 * size_ext + xmax] / maxch;
else
a = working_space[2 * size_ext + i + l] / maxch;
b = a - nip;
if(a + nip <= 0)
a=1;
else
a = TMath::Sqrt(a + nip);
b = b / a;
b = TMath::Exp(b);
sp = sp + b;
if((i - l + 1) < xmin)
a = working_space[2 * size_ext + xmin] / maxch;
else
a = working_space[2 * size_ext + i - l + 1] / maxch;
b = a - nim;
if(a + nim <= 0)
a = 1;
else
a = TMath::Sqrt(a + nim);
b = b / a;
b = TMath::Exp(b);
sm = sm + b;
}
a = sp / sm;
a = working_space[i + 1] = working_space[i] * a;
nom = nom + a;
}
for(i = xmin; i <= xmax; i++){
working_space[i] = working_space[i] / nom;
}
for(j = 0; j < size_ext; j++)
working_space[size_ext + j] = working_space[j] * plocha;
for(j = 0; j < size_ext; j++){
working_space[2 * size_ext + j] = working_space[size_ext + j];
}
if(backgroundRemove == true){
for(i = 1; i <= numberIterations; i++){
for(j = i; j < size_ext - i; j++){
a = working_space[size_ext + j];
b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
if(b < a)
a = b;
working_space[j] = a;
}
for(j = i; j < size_ext - i; j++)
working_space[size_ext + j] = working_space[j];
}
for(j = 0; j < size_ext; j++){
working_space[size_ext + j] = working_space[2 * size_ext + j] - working_space[size_ext + j];
}
}
}
//deconvolution starts
area = 0;
lh_gold = -1;
posit = 0;
maximum = 0;
//generate response vector
for(i = 0; i < size_ext; i++){
lda = (double)i - 3 * sigma;
lda = lda * lda / (2 * sigma * sigma);
j = (int)(1000 * TMath::Exp(-lda));
lda = j;
if(lda != 0)
lh_gold = i + 1;
working_space[i] = lda;
area = area + lda;
if(lda > maximum){
maximum = lda;
posit = i;
}
}
//read source vector
for(i = 0; i < size_ext; i++)
working_space[2 * size_ext + i] = TMath::Abs(working_space[size_ext + i]);
//create matrix at*a(vector b)
i = lh_gold - 1;
if(i > size_ext)
i = size_ext;
imin = -i,imax = i;
for(i = imin; i <= imax; i++){
lda = 0;
jmin = 0;
if(i < 0)
jmin = -i;
jmax = lh_gold - 1 - i;
if(jmax > (lh_gold - 1))
jmax = lh_gold - 1;
for(j = jmin;j <= jmax; j++){
ldb = working_space[j];
ldc = working_space[i + j];
lda = lda + ldb * ldc;
}
working_space[size_ext + i - imin] = lda;
}
//create vector p
i = lh_gold - 1;
imin = -i,imax = size_ext + i - 1;
for(i = imin; i <= imax; i++){
lda = 0;
for(j = 0; j <= (lh_gold - 1); j++){
ldb = working_space[j];
k = i + j;
if(k >= 0 && k < size_ext){
ldc = working_space[2 * size_ext + k];
lda = lda + ldb * ldc;
}
}
working_space[4 * size_ext + i - imin] = lda;
}
//move vector p
for(i = imin; i <= imax; i++)
working_space[2 * size_ext + i - imin] = working_space[4 * size_ext + i - imin];
//initialization of resulting vector
for(i = 0; i < size_ext; i++)
working_space[i] = 1;
//START OF ITERATIONS
for(lindex = 0; lindex < deconIterations; lindex++){
for(i = 0; i < size_ext; i++){
if(TMath::Abs(working_space[2 * size_ext + i]) > 0.00001 && TMath::Abs(working_space[i]) > 0.00001){
lda=0;
jmin = lh_gold - 1;
if(jmin > i)
jmin = i;
jmin = -jmin;
jmax = lh_gold - 1;
if(jmax > (size_ext - 1 - i))
jmax=size_ext-1-i;
for(j = jmin; j <= jmax; j++){
ldb = working_space[j + lh_gold - 1 + size_ext];
ldc = working_space[i + j];
lda = lda + ldb * ldc;
}
ldb = working_space[2 * size_ext + i];
if(lda != 0)
lda = ldb / lda;
else
lda = 0;
ldb = working_space[i];
lda = lda * ldb;
working_space[3 * size_ext + i] = lda;
}
}
for(i = 0; i < size_ext; i++){
working_space[i] = working_space[3 * size_ext + i];
}
}
//shift resulting spectrum
for(i=0;i