source: trunk/MagicSoft/Mars/mbase/MSpline3.cc@ 9448

Last change on this file since 9448 was 9424, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 7.7 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 1/2009 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: Software Development, 2000-2009
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MSpline3
28//
29// This is a extension of TSpline3. In addition to TSpline3 it allows access
30// to Xmin, Xman and Np. The construction is a bit simplified because no
31// title hase to be given (it can be given later by SetTitle anyway)
32// and is provides constructors which allow to scale the x-values by
33// pre-defined multiplier (e.g. frequency) to create the spline.
34//
35//////////////////////////////////////////////////////////////////////////////
36#include "MSpline3.h"
37
38#include <TF1.h>
39
40#include "MArrayD.h"
41
42ClassImp(MSpline3);
43
44using namespace std;
45
46// --------------------------------------------------------------------------
47//
48// Constructor.
49//
50MSpline3::MSpline3(const TF1 &f, const char *opt, Double_t valbeg, Double_t valend)
51 : TSpline3("MSpline3", f.GetXmin(), f.GetXmax(), &f, f.GetNpx(), opt, valbeg, valend)
52{
53}
54
55MSpline3::MSpline3(const TF1 &f, Double_t freq, const char *opt,Double_t valbeg, Double_t valend)
56 : TSpline3("MSpline3", f.GetXmin()*freq, f.GetXmax()*freq, ConvertFunc(f, freq).GetArray(), f.GetNpx(), opt, valbeg, valend)
57{
58}
59
60// --------------------------------------------------------------------------
61//
62// This is a helper to convert the x-values by multiplying with freq
63// before initializing the spline
64//
65TGraph *MSpline3::ConvertSpline(const TSpline &s, Float_t freq) const
66{
67 const UInt_t npx = s.GetNpx();
68
69 // WARNING: This is a stupid workaround because the TSpline3-
70 // constructor takes a pointer as input! It is not thread-safe!
71 static TGraph g;
72 g.Set(npx);
73
74 for (UInt_t i=0; i<npx; i++)
75 {
76 Double_t x, y;
77 s.GetKnot(i, x, y);
78 g.SetPoint(i, x*freq, y);
79 }
80
81 return &g;
82}
83
84// --------------------------------------------------------------------------
85//
86// This is a helper to convert the x-values by multiplying with freq
87// before initializing the spline
88//
89TGraph *MSpline3::ConvertGraph(const TGraph &s, Float_t freq) const
90{
91 const UInt_t npx = s.GetN();
92
93 // WARNING: This is a stupid workaround because the TSpline3-
94 // constructor takes a pointer as input! It is not thread-safe!
95 static TGraph g;
96 g.Set(npx);
97
98 for (UInt_t i=0; i<npx; i++)
99 {
100 Double_t x, y;
101 s.GetPoint(i, x, y);
102 g.SetPoint(i, x*freq, y);
103 }
104
105 return &g;
106}
107
108// --------------------------------------------------------------------------
109//
110// This is a helper to convert the x-values by multiplying with freq
111// before initializing the spline. The conversion from the function to
112// a discrete binning is done similar to the constructor of TSpline
113//
114MArrayD &MSpline3::ConvertFunc(const TF1 &f, Float_t freq) const
115{
116 const UInt_t npx = f.GetNpx();
117
118 // WARNING: This is a stupid workaround because the TSpline3-
119 // constructor takes a pointer as input! It is not thread-safe!
120 static MArrayD g;
121 g.Set(npx);
122
123 const Double_t step = (f.GetXmax()-f.GetXmin())/(npx-1);
124
125 for (UInt_t i=0; i<npx; ++i)
126 {
127 const Double_t x = f.GetXmin() + i*step;
128 g[i] = f.Eval(x);
129 }
130
131 return g;
132}
133
134// --------------------------------------------------------------------------
135//
136// Return the integral in the splines bin i up to x.
137//
138// The TSpline3 in the Interval [fX[i], fX[i+1]] is defined as:
139//
140// dx = x-fX[i]
141// y = fY + dx*fB + dx*dx*fC + dx*dx*dx*fD
142//
143// This yields the integral:
144//
145// int(y) = dx*fY + 1/2*dx*dx*fB + 1/3*dx*dx*dx*fC + 1/4*dx*dx*dx*dx*fD
146// = dx*(fY + dx*(1/2*fB + dx*(1/3*fC + dx*(1/4*fD))))
147//
148// Which gives for the integral range [fX[i], fX[i]+w]:
149// int(fX[i]+w)-int(fX[i]) = w*(fY + w*(1/2*fB + w*(1/3*fC + w*(1/4*fD))))
150//
151// and for the integral range [fX[i]+w, fX[i+1]]:
152// int(fX[i+1])-int(fX[i]+w) = `
153// W*(fY + W*(1/2*fB + W*(1/3*fC + W*(1/4*fD)))) -
154// w*(fY + w*(1/2*fB + w*(1/3*fC + w*(1/4*fD))))
155// with
156// W := fX[i+1]-fX[i]
157//
158Double_t MSpline3::Integral(Int_t i, Double_t x) const
159{
160 Double_t x0, y, b, c, d;
161 const_cast<MSpline3*>(this)->GetCoeff(i, x0, y, b, c, d);
162
163 const Double_t w = x-x0;
164
165 return w*(y + w*(b/2 + w*(c/3 + w*d/4)));
166}
167
168// --------------------------------------------------------------------------
169//
170// Return the integral of the spline's bin i.
171//
172Double_t MSpline3::Integral(Int_t i) const
173{
174 Double_t x, y;
175
176 GetKnot(i+1, x, y);
177
178 return Integral(i, x);
179}
180
181// --------------------------------------------------------------------------
182//
183// Return the integral from a to b
184//
185Double_t MSpline3::Integral(Double_t a, Double_t b) const
186{
187 const Int_t n = FindX(a);
188 const Int_t m = FindX(b);
189
190 Double_t sum = -Integral(n, a);
191
192 for (int i=n; i<=m-1; i++)
193 sum += Integral(i);
194
195 sum += Integral(m, b);
196
197 return sum;
198}
199
200// --------------------------------------------------------------------------
201//
202// Return the integral between Xmin and Xmax
203//
204Double_t MSpline3::Integral() const
205{
206 Double_t sum = 0;
207
208 for (int i=0; i<GetNp()-1; i++)
209 sum += Integral(i);
210
211 return sum;
212}
213
214// --------------------------------------------------------------------------
215//
216// Return the integral between Xmin and Xmax of int( f(x)*sin(x) )
217//
218// The x-axis is assumed to be in degrees
219//
220Double_t MSpline3::IntegralSolidAngle() const
221{
222 const Int_t n = GetNp();
223
224 MArrayD x(n);
225 MArrayD y(n);
226
227 for (int i=0; i<n; i++)
228 {
229 GetKnot(i, x[i], y[i]);
230
231 x[i] *= TMath::DegToRad();
232 y[i] *= TMath::Sin(x[i]);
233 }
234
235 return TMath::TwoPi()*MSpline3(x.GetArray(), y.GetArray(), n).Integral();
236}
237
238
239// FIXME: As soon as TSpline3 allows access to fPoly we can implement
240// a much faster evaluation of the spline, especially in
241// special conditions like in MAnalogSignal::AddPulse
242// This will be the case for root > 5.22/00
243
244/*
245Double_t MSpline3::EvalFast(Double_t x) const
246{
247 // Eval this spline at x
248 const Int_t klow=FindFast(x);
249 return fPoly[klow].Eval(x);
250}
251
252Int_t MSpline3::FindFast(Double_t x) const
253{
254 //
255 // If out of boundaries, extrapolate
256 // It may be badly wrong
257
258 // if (x<=fXmin)
259 // return 0;
260 //
261 // if (x>=fXmax)
262 // return fNp-1;
263
264 //
265 // Equidistant knots, use histogramming
266 if (fKstep)
267 return TMath::Min(Int_t((x-fXmin)/fDelta),fNp-1);
268
269 //
270 // Non equidistant knots, binary search
271 Int_t klow = 0;
272 Int_t khig = fNp-1;
273
274 Int_t khalf;
275 while (khig-klow>1)
276 if(x>fPoly[khalf=(klow+khig)/2].X())
277 klow=khalf;
278 else
279 khig=khalf;
280
281 // This could be removed, sanity check
282 //if(!(fPoly[klow].X()<=x && x<=fPoly[klow+1].X()))
283 // Error("Eval",
284 // "Binary search failed x(%d) = %f < %f < x(%d) = %f\n",
285 // klow,fPoly[klow].X(),x,fPoly[klow+1].X());
286
287 return klow;
288}
289*/
Note: See TracBrowser for help on using the repository browser.