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