GNU Octave  4.0.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Quad.h
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1993-2015 John W. Eaton
4 
5 This file is part of Octave.
6 
7 Octave is free software; you can redistribute it and/or modify it
8 under the terms of the GNU General Public License as published by the
9 Free Software Foundation; either version 3 of the License, or (at your
10 option) any later version.
11 
12 Octave is distributed in the hope that it will be useful, but WITHOUT
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
15 for more details.
16 
17 You should have received a copy of the GNU General Public License
18 along with Octave; see the file COPYING. If not, see
19 <http://www.gnu.org/licenses/>.
20 
21 */
22 
23 #if !defined (octave_Quad_h)
24 #define octave_Quad_h 1
25 
26 #include <cfloat>
27 
28 #include "dColVector.h"
29 #include "fColVector.h"
30 #include "lo-math.h"
31 
32 typedef double (*integrand_fcn) (double x);
33 typedef float (*float_integrand_fcn) (float x);
34 
35 // FIXME: would be nice to not have to have this global variable.
36 // Nonzero means an error occurred in the calculation of the integrand
37 // function, and the user wants us to quit.
38 extern OCTAVE_API int quad_integration_error;
39 
40 #include "Quad-opts.h"
41 
42 class
43 OCTAVE_API
44 Quad : public Quad_options
45 {
46 public:
47 
49  : Quad_options (), f (fcn), ff () { }
50 
52  : Quad_options (), f (), ff (fcn) { }
53 
54  virtual ~Quad (void) { }
55 
56  virtual double integrate (void)
57  {
58  octave_idx_type ier, neval;
59  double abserr;
60  return do_integrate (ier, neval, abserr);
61  }
62 
63  virtual float float_integrate (void)
64  {
65  octave_idx_type ier, neval;
66  float abserr;
67  return do_integrate (ier, neval, abserr);
68  }
69 
70  virtual double integrate (octave_idx_type& ier)
71  {
72  octave_idx_type neval;
73  double abserr;
74  return do_integrate (ier, neval, abserr);
75  }
76 
77  virtual float float_integrate (octave_idx_type& ier)
78  {
79  octave_idx_type neval;
80  float abserr;
81  return do_integrate (ier, neval, abserr);
82  }
83 
84  virtual double integrate (octave_idx_type& ier, octave_idx_type& neval)
85  {
86  double abserr;
87  return do_integrate (ier, neval, abserr);
88  }
89 
90  virtual float float_integrate (octave_idx_type& ier, octave_idx_type& neval)
91  {
92  float abserr;
93  return do_integrate (ier, neval, abserr);
94  }
95 
96  virtual double integrate (octave_idx_type& ier, octave_idx_type& neval,
97  double& abserr)
98  {
99  return do_integrate (ier, neval, abserr);
100  }
101 
102  virtual float float_integrate (octave_idx_type& ier, octave_idx_type& neval,
103  float& abserr)
104  {
105  return do_integrate (ier, neval, abserr);
106  }
107 
108  virtual double do_integrate (octave_idx_type& ier, octave_idx_type& neval,
109  double& abserr) = 0;
110 
111  virtual float do_integrate (octave_idx_type& ier, octave_idx_type& neval,
112  float& abserr) = 0;
113 
114 protected:
115 
118 };
119 
120 class
121 OCTAVE_API
122 DefQuad : public Quad
123 {
124 public:
125 
127  : Quad (fcn), lower_limit (0.0), upper_limit (1.0), singularities () { }
128 
129  DefQuad (integrand_fcn fcn, double ll, double ul)
130  : Quad (fcn), lower_limit (ll), upper_limit (ul), singularities () { }
131 
132  DefQuad (integrand_fcn fcn, double ll, double ul,
133  const ColumnVector& sing)
134  : Quad (fcn), lower_limit (ll), upper_limit (ul),
135  singularities (sing) { }
136 
137  DefQuad (integrand_fcn fcn, const ColumnVector& sing)
138  : Quad (fcn), lower_limit (0.0), upper_limit (1.0),
139  singularities (sing) { }
140 
141  ~DefQuad (void) { }
142 
143  double do_integrate (octave_idx_type& ier, octave_idx_type& neval,
144  double& abserr);
145 
146  float do_integrate (octave_idx_type& ier, octave_idx_type& neval,
147  float& abserr);
148 
149 private:
150 
151  double lower_limit;
152  double upper_limit;
153 
155 };
156 
157 class
158 OCTAVE_API
159 IndefQuad : public Quad
160 {
161 public:
162 
163  enum IntegralType { bound_to_inf, neg_inf_to_bound, doubly_infinite };
164 
166  : Quad (fcn), bound (0.0), type (bound_to_inf), integration_error (0) { }
167 
169  : Quad (fcn), bound (b), type (t), integration_error (0) { }
170 
171  ~IndefQuad (void) { }
172 
173  double do_integrate (octave_idx_type& ier, octave_idx_type& neval,
174  double& abserr);
175 
176  float do_integrate (octave_idx_type& ier, octave_idx_type& neval,
177  float& abserr);
178 
179 private:
180 
181  double bound;
184 };
185 
186 class
187 OCTAVE_API
188 FloatDefQuad : public Quad
189 {
190 public:
191 
193  : Quad (fcn), lower_limit (0.0), upper_limit (1.0), singularities () { }
194 
195  FloatDefQuad (float_integrand_fcn fcn, float ll, float ul)
196  : Quad (fcn), lower_limit (ll), upper_limit (ul), singularities () { }
197 
198  FloatDefQuad (float_integrand_fcn fcn, float ll, float ul,
199  const FloatColumnVector& sing)
200  : Quad (fcn), lower_limit (ll), upper_limit (ul),
201  singularities (sing) { }
202 
204  : Quad (fcn), lower_limit (0.0), upper_limit (1.0),
205  singularities (sing) { }
206 
207  ~FloatDefQuad (void) { }
208 
209  double do_integrate (octave_idx_type& ier, octave_idx_type& neval,
210  double& abserr);
211 
212  float do_integrate (octave_idx_type& ier, octave_idx_type& neval,
213  float& abserr);
214 
215 private:
216 
217  float lower_limit;
218  float upper_limit;
219 
221 };
222 
223 class
224 OCTAVE_API
225 FloatIndefQuad : public Quad
226 {
227 public:
228 
229  enum IntegralType { bound_to_inf, neg_inf_to_bound, doubly_infinite };
230 
232  : Quad (fcn), bound (0.0), type (bound_to_inf), integration_error (0) { }
233 
235  : Quad (fcn), bound (b), type (t), integration_error (0) { }
236 
237  ~FloatIndefQuad (void) { }
238 
239  double do_integrate (octave_idx_type& ier, octave_idx_type& neval,
240  double& abserr);
241 
242  float do_integrate (octave_idx_type& ier, octave_idx_type& neval,
243  float& abserr);
244 
245 private:
246 
247  float bound;
250 };
251 
252 #endif
float(* float_integrand_fcn)(float x)
Definition: Quad.h:33
DefQuad(integrand_fcn fcn, double ll, double ul, const ColumnVector &sing)
Definition: Quad.h:132
FloatIndefQuad(float_integrand_fcn fcn)
Definition: Quad.h:231
~DefQuad(void)
Definition: Quad.h:141
IndefQuad(integrand_fcn fcn, double b, IntegralType t)
Definition: Quad.h:168
virtual double integrate(void)
Definition: Quad.h:56
float_integrand_fcn ff
Definition: Quad.h:117
double(* integrand_fcn)(double x)
Definition: Quad.h:32
ColumnVector singularities
Definition: Quad.h:154
virtual float float_integrate(void)
Definition: Quad.h:63
double upper_limit
Definition: Quad.h:152
FloatDefQuad(float_integrand_fcn fcn, const FloatColumnVector &sing)
Definition: Quad.h:203
DefQuad(integrand_fcn fcn, const ColumnVector &sing)
Definition: Quad.h:137
IntegralType
Definition: Quad.h:163
~FloatIndefQuad(void)
Definition: Quad.h:237
virtual double integrate(octave_idx_type &ier, octave_idx_type &neval, double &abserr)
Definition: Quad.h:96
~FloatDefQuad(void)
Definition: Quad.h:207
FloatIndefQuad(float_integrand_fcn fcn, double b, IntegralType t)
Definition: Quad.h:234
DefQuad(integrand_fcn fcn)
Definition: Quad.h:126
FloatColumnVector singularities
Definition: Quad.h:220
int integration_error
Definition: Quad.h:183
FloatDefQuad(float_integrand_fcn fcn, float ll, float ul, const FloatColumnVector &sing)
Definition: Quad.h:198
float upper_limit
Definition: Quad.h:218
Quad(float_integrand_fcn fcn)
Definition: Quad.h:51
double bound
Definition: Quad.h:181
F77_RET_T const double const double * f
DefQuad(integrand_fcn fcn, double ll, double ul)
Definition: Quad.h:129
Definition: Quad.h:42
virtual float float_integrate(octave_idx_type &ier, octave_idx_type &neval, float &abserr)
Definition: Quad.h:102
Quad(integrand_fcn fcn)
Definition: Quad.h:48
virtual float float_integrate(octave_idx_type &ier, octave_idx_type &neval)
Definition: Quad.h:90
float lower_limit
Definition: Quad.h:217
FloatDefQuad(float_integrand_fcn fcn, float ll, float ul)
Definition: Quad.h:195
OCTAVE_API int quad_integration_error
Definition: Quad.cc:39
Definition: Quad.h:120
virtual double integrate(octave_idx_type &ier)
Definition: Quad.h:70
virtual double integrate(octave_idx_type &ier, octave_idx_type &neval)
Definition: Quad.h:84
IndefQuad(integrand_fcn fcn)
Definition: Quad.h:165
float bound
Definition: Quad.h:247
double lower_limit
Definition: Quad.h:151
IntegralType type
Definition: Quad.h:182
IntegralType type
Definition: Quad.h:248
~IndefQuad(void)
Definition: Quad.h:171
virtual float float_integrate(octave_idx_type &ier)
Definition: Quad.h:77
integrand_fcn f
Definition: Quad.h:116
FloatDefQuad(float_integrand_fcn fcn)
Definition: Quad.h:192
int integration_error
Definition: Quad.h:249
F77_RET_T const double * x
virtual ~Quad(void)
Definition: Quad.h:54