This file is indexed.

/usr/include/root/Math/GSLIntegrator.h is in libroot-math-mathmore-dev 5.34.30-0ubuntu8.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
// @(#)root/mathmore:$Id$
// Authors: L. Moneta, A. Zsenei   08/2005

 /**********************************************************************
  *                                                                    *
  * Copyright (c) 2004 ROOT Foundation,  CERN/PH-SFT                   *
  *                                                                    *
  * This library is free software; you can redistribute it and/or      *
  * modify it under the terms of the GNU General Public License        *
  * as published by the Free Software Foundation; either version 2     *
  * of the License, or (at your option) any later version.             *
  *                                                                    *
  * This library is distributed in the hope that it will be useful,    *
  * but WITHOUT ANY WARRANTY; without even the implied warranty of     *
  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU   *
  * General Public License for more details.                           *
  *                                                                    *
  * You should have received a copy of the GNU General Public License  *
  * along with this library (see file COPYING); if not, write          *
  * to the Free Software Foundation, Inc., 59 Temple Place, Suite      *
  * 330, Boston, MA 02111-1307 USA, or contact the author.             *
  *                                                                    *
  **********************************************************************/

// Header file for class GSLIntegrator
//
// Created by: Lorenzo Moneta  at Thu Nov 11 14:22:32 2004
//
// Last update: Thu Nov 11 14:22:32 2004
//
#ifndef ROOT_Math_GSLIntegrator
#define ROOT_Math_GSLIntegrator


#ifndef ROOT_Math_VirtualIntegrator
#include "Math/VirtualIntegrator.h"
#endif

#ifndef ROOT_Math_IntegrationTypes
#include "Math/IntegrationTypes.h"
#endif

#ifndef ROOT_Math_IFunctionfwd
#include "Math/IFunctionfwd.h"
#endif




#ifndef ROOT_Math_GSLFunctionAdapter
#include "Math/GSLFunctionAdapter.h"
#endif

#include <vector>


/**

@defgroup Integration Numerical Integration
@ingroup NumAlgo 
*/


namespace ROOT {
namespace Math {


   
   class GSLIntegrationWorkspace;
   class GSLFunctionWrapper;
   
   //_________________________________________________________________________
   /**
      
   Class for performing numerical integration of a function in one dimension.
   It uses the numerical integration algorithms of GSL, which reimplements the
   algorithms used in the QUADPACK, a numerical integration package written in Fortran.
    
   Various types of adaptive and non-adaptive integration are supported. These include
   integration over infinite and semi-infinite ranges and singular integrals.
    
   The integration type is selected using the Integration::type enumeration
   in the class constructor.
   The default type is adaptive integration with singularity
   (ADAPTIVESINGULAR or QAGS in the QUADPACK convention) applying a Gauss-Kronrod 21-point integration rule.
   In the case of ADAPTIVE type, the integration rule can also be specified via the
   Integration::GKRule. The default rule is 31 points.
    
   In the case of integration over infinite and semi-infinite ranges, the type used is always
   ADAPTIVESINGULAR applying a transformation from the original interval into (0,1).
    
   The ADAPTIVESINGULAR type is the most sophicticated type. When performances are
   important, it is then recommened to use the NONADAPTIVE type in case of smooth functions or
   ADAPTIVE with a lower Gauss-Kronrod rule.
    
   For detailed description on GSL integration algorithms see the
   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Numerical-Integration.html">GSL Manual</A>.
    
    
   @ingroup Integration
   */
   
   
   class GSLIntegrator : public VirtualIntegratorOneDim  {
      
   public:
      
      
      
      // constructors
      
      
      /** Default constructor of GSL Integrator for Adaptive Singular integration
      
      @param absTol desired absolute Error
      @param relTol desired relative Error
      @param size maximum number of sub-intervals
      */
      
      GSLIntegrator(double absTol = 1.E-9, double relTol = 1E-6, size_t size = 1000);
      
      
      
      
      /** constructor of GSL Integrator. In the case of Adaptive integration the Gauss-Krond rule of 31 points is used
         
         @param type type of integration. The possible types are defined in the Integration::Type enumeration
         @param absTol desired absolute Error
         @param relTol desired relative Error
         @param size maximum number of sub-intervals
         */
      
      
      GSLIntegrator(const Integration::Type type, double absTol = 1.E-9, double relTol = 1E-6, size_t size = 1000);
      
      
      /**
         generic constructor for GSL Integrator
       
       @param type type of integration. The possible types are defined in the Integration::Type enumeration
       @param rule Gauss-Kronrod rule. It is used only for ADAPTIVE::Integration types. The possible rules are defined in the Integration::GKRule enumeration
       @param absTol desired absolute Error
       @param relTol desired relative Error
       @param size maximum number of sub-intervals
       
       */
      
      GSLIntegrator(const Integration::Type type, const Integration::GKRule rule, double absTol = 1.E-9, double relTol = 1E-6, size_t size = 1000);
      

      /** constructor of GSL Integrator. In the case of Adaptive integration the Gauss-Krond rule of 31 points is used
          This is used by the plug-in manager (need a char * instead of enumerations)
         
         @param type type of integration. The possible types are defined in the Integration::Type enumeration
         @param rule Gauss-Kronrod rule (from 1 to 6)
         @param absTol desired absolute Error
         @param relTol desired relative Error
         @param size maximum number of sub-intervals
         */            
      GSLIntegrator(const char *  type, int rule, double absTol, double relTol, size_t size );
      
      virtual ~GSLIntegrator();
      //~GSLIntegrator();
      
      // disable copy ctrs
   private:
         
      GSLIntegrator(const GSLIntegrator &);
      GSLIntegrator & operator=(const GSLIntegrator &);
      
   public:
         
         
         // template methods for generic functors
         
         /**
         method to set the a generic integration function
          
          @param f integration function. The function type must implement the assigment operator, <em>  double  operator() (  double  x ) </em>
          
          */
         
         
      void SetFunction(const IGenFunction &f);

      /**
         Set function from a GSL pointer function type 
       */
      void SetFunction( GSLFuncPointer f, void * p = 0); 
      
      // methods using IGenFunction
      
      /**
         evaluate the Integral of a function f over the defined interval (a,b)
       @param f integration function. The function type must implement the mathlib::IGenFunction interface
       @param a lower value of the integration interval
       @param b upper value of the integration interval
       */
      
      double Integral(const IGenFunction & f, double a, double b);
      
      
      /**
         evaluate the Integral of a function f over the infinite interval (-inf,+inf)
       @param f integration function. The function type must implement the mathlib::IGenFunction interface
       */

      double Integral(const IGenFunction & f);
   
      /**
       evaluate the Cauchy principal value of the integral of  a previously defined function f over
        the defined interval (a,b) with a singularity at c 
        @param a lower interval value
        @param b lower interval value
        @param c singular value of f
        */   
      double IntegralCauchy(double a, double b, double c);

      /**
       evaluate the Cauchy principal value of the integral of  a function f over the defined interval (a,b)
        with a singularity at c 
        @param f integration function. The function type must implement the mathlib::IGenFunction interface
        @param a lower interval value
        @param b lower interval value
        @param c singular value of f
      */
      double IntegralCauchy(const IGenFunction & f, double a, double b, double c);
      
      /**
         evaluate the Integral of a function f over the semi-infinite interval (a,+inf)
       @param f integration function. The function type must implement the mathlib::IGenFunction interface
       @param a lower value of the integration interval
       
       */
      double IntegralUp(const IGenFunction & f, double a );
      
      /**
         evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b)
       @param f integration function. The function type must implement the mathlib::IGenFunction interface
       @param b upper value of the integration interval
       */
      double IntegralLow(const IGenFunction & f, double b );
      
      /**
         evaluate the Integral of a function f with known singular points over the defined Integral (a,b)
       @param f integration function. The function type must implement the mathlib::IGenFunction interface
       @param pts vector containing both the function singular points and the lower/upper edges of the interval. The vector must have as first element the lower edge of the integration Integral ( \a a) and last element the upper value.
       
       */
      double Integral(const IGenFunction & f, const std::vector<double> & pts );
      
      // evaluate using cached function
      
      /**
         evaluate the Integral over the defined interval (a,b) using the function previously set with GSLIntegrator::SetFunction method
       @param a lower value of the integration interval
       @param b upper value of the integration interval
       */
      
      double Integral(double a, double b);

      
      /**
         evaluate the Integral over the infinite interval (-inf,+inf) using the function previously set with GSLIntegrator::SetFunction method.
       */
      double Integral( );
      
      /**
         evaluate the Integral of a function f over the semi-infinite interval (a,+inf) using the function previously set with GSLIntegrator::SetFunction method.
       @param a lower value of the integration interval
       */
      double IntegralUp(double a );
      
      /**
         evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b) using the function previously set with GSLIntegrator::SetFunction method.
       @param b upper value of the integration interval
       */
      double IntegralLow( double b );
      
      /**
         evaluate the Integral over the defined interval (a,b) using the function previously set with GSLIntegrator::SetFunction method. The function has known singular points.
       @param pts vector containing both the function singular points and the lower/upper edges of the interval. The vector must have as first element the lower edge of the integration Integral ( \a a) and last element the upper value.
       
       */
      double Integral( const std::vector<double> & pts);
      
      // evaluate using free function pointer (same GSL signature)
      
      /**
         signature for function pointers used by GSL
       */
      //typedef double ( * GSLFuncPointer ) ( double, void * );
      
      /**
         evaluate the Integral of  of a function f over the defined interval (a,b) passing a free function pointer
       The integration function must be a free function and have a signature consistent with GSL functions:
       
       <em>double my_function ( double x, void * p ) { ...... } </em>
       
       This method is the most efficient since no internal adapter to GSL function is created.
       @param f pointer to the integration function
       @param p pointer to the Parameters of the function
       @param a lower value of the integration interval
       @param b upper value of the integration interval
       
       */
      double Integral(GSLFuncPointer f, void * p, double a, double b);
      
      /**
         evaluate the Integral  of a function f over the infinite interval (-inf,+inf) passing a free function pointer
       */
      double Integral(GSLFuncPointer f, void * p);
      
      /**
         evaluate the Integral of a function f over the semi-infinite interval (a,+inf) passing a free function pointer
       */
      double IntegralUp(GSLFuncPointer f, void * p, double a);
      
      /**
         evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b) passing a free function pointer
       */
      double IntegralLow(GSLFuncPointer f, void * p, double b);
      
      /**
         evaluate the Integral of a function f with knows singular points over the over a defined interval passing a free function pointer
       */
      double Integral(GSLFuncPointer f, void * p, const std::vector<double> & pts);
      
      /**
         return  the Result of the last Integral calculation
       */
      double Result() const;
      
      /**
         return the estimate of the absolute Error of the last Integral calculation
       */
      double Error() const;
      
      /**
         return the Error Status of the last Integral calculation
       */
      int Status() const;
            
      /** 
          return number of function evaluations in calculating the integral 
      */
      int NEval() const { return fNEval; }
      
      // setter for control Parameters  (getters are not needed so far )
      
      /**
         set the desired relative Error
       */
      void SetRelTolerance(double relTolerance);
      
      
      /**
         set the desired absolute Error
       */
      void SetAbsTolerance(double absTolerance);
      
      /**
         set the integration rule (Gauss-Kronrod rule).
       The possible rules are defined in the Integration::GKRule enumeration.
       The integration rule can be modified only for ADAPTIVE type integrations
       */
      void SetIntegrationRule(Integration::GKRule );
      
      /// set the options 
      virtual void SetOptions(const ROOT::Math::IntegratorOneDimOptions & opt);

      ///  get the option used for the integration 
      virtual ROOT::Math::IntegratorOneDimOptions Options() const;
      
      /// get type name
      IntegrationOneDim::Type GetType() const { return fType; }
      
      /** 
          return the name
      */
      const char * GetTypeName() const; 


   protected:
         
      // internal method to check validity of GSL function pointer
      bool CheckFunction(); 

      
   private:
         
      Integration::Type fType;
      Integration::GKRule fRule;
      double fAbsTol;
      double fRelTol;
      size_t fSize;
      size_t fMaxIntervals;
      
      // cache Error, Result and Status of integration
      
      double fResult;
      double fError;
      int fStatus;
      int fNEval; 
      
      // GSLIntegrationAlgorithm * fAlgorithm;
      
      GSLFunctionWrapper  *     fFunction;
      GSLIntegrationWorkspace * fWorkspace;
     
   };
   
   



} // namespace Math
} // namespace ROOT


#endif /* ROOT_Math_GSLIntegrator */