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
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
/// Get a function that returns the gradient of the provided multivariate, scalar-valued function.
///
/// The gradient is computed using forward-mode automatic differentiation.
///
/// # Arguments
///
/// * `f` - Multivariate, scalar-valued function, $f:\mathbb{R}^{n}\to\mathbb{R}$.
/// * `func_name` - Name of the function that will return the gradient of $f(\mathbf{x})$ at any
/// point $\mathbf{x}\in\mathbb{R}^{n}$.
/// * `param_type` (optional) - Type of the extra runtime parameter `p` that is passed to `f`.
/// Defaults to `[f64]` (implying that `f` accepts `p: &[f64]`).
///
/// # Warning
///
/// `f` cannot be defined as closure. It must be defined as a function.
///
/// # Note
///
/// The function produced by this macro will perform $n$ evaluations of $f(\mathbf{x})$ to evaluate
/// its gradient.
///
/// # Examples
///
/// ## Basic Example
///
/// Compute the gradient of
///
/// $$f(\mathbf{x})=x_{0}^{5}+\sin^{3}{x_{1}}$$
///
/// at $\mathbf{x}=(5,8)^{T}$, and compare the result to the true result of
///
/// $$
/// \nabla f\left((5,8)^{T}\right)=
/// \begin{bmatrix}
/// 3125 \\\\
/// 3\sin^{2}{(8)}\cos{(8)}
/// \end{bmatrix}
/// $$
///
/// #### Using standard vectors
///
/// ```
/// use linalg_traits::{Scalar, Vector};
/// use numtest::*;
///
/// use numdiff::{get_gradient, Dual, DualVector};
///
/// // Define the function, f(x).
/// fn f<S: Scalar, V: Vector<S>>(x: &V, _p: &[f64]) -> S {
/// x.vget(0).powi(5) + x.vget(1).sin().powi(3)
/// }
///
/// // Define the evaluation point.
/// let x0 = vec![5.0, 8.0];
///
/// // Parameter vector (empty for this example).
/// let p = [];
///
/// // Autogenerate the function "g" that can be used to compute the gradient of f(x) at any point
/// // x.
/// get_gradient!(f, g);
///
/// // Function defining the true gradient of f(x).
/// let g_true = |x: &Vec<f64>| vec![5.0 * x[0].powi(4), 3.0 * x[1].sin().powi(2) * x[1].cos()];
///
/// // Evaluate the gradient using "g".
/// let g_eval: Vec<f64> = g(&x0, &p);
///
/// // Verify that the gradient function obtained using get_gradient! computes the gradient
/// // correctly.
/// assert_arrays_equal_to_decimal!(g(&x0, &p), g_true(&x0), 15);
/// ```
///
/// #### Using other vector types
///
/// The function produced by `get_gradient!` can accept _any_ type for `x0`, as long as it
/// implements the `linalg_traits::Vector` trait.
///
/// ```
/// use faer::Mat;
/// use linalg_traits::{Scalar, Vector};
/// use nalgebra::{dvector, DVector, SVector};
/// use ndarray::{array, Array1};
///
/// use numdiff::{get_gradient, Dual, DualVector};
///
/// // Define the function, f(x).
/// fn f<S: Scalar, V: Vector<S>>(x: &V, _p: &[f64]) -> S {
/// x.vget(0).powi(5) + x.vget(1).sin().powi(3)
/// }
///
/// // Parameter vector (empty for this example).
/// let p = [];
///
/// // Autogenerate the function "g" that can be used to compute the gradient of f(x) at any point
/// // x.
/// get_gradient!(f, g);
///
/// // nalgebra::DVector
/// let x0: DVector<f64> = dvector![5.0, 8.0];
/// let g_eval: DVector<f64> = g(&x0, &p);
///
/// // nalgebra::SVector
/// let x0: SVector<f64, 2> = SVector::from_slice(&[5.0, 8.0]);
/// let g_eval: SVector<f64, 2> = g(&x0, &p);
///
/// // ndarray::Array1
/// let x0: Array1<f64> = array![5.0, 8.0];
/// let g_eval: Array1<f64> = g(&x0, &p);
///
/// // faer::Mat
/// let x0: Mat<f64> = Mat::from_slice(&[5.0, 8.0]);
/// let g_eval: Mat<f64> = g(&x0, &p);
/// ```
///
/// ## Example Passing Runtime Parameters
///
/// Compute the gradient of a parameterized function
///
/// $$f(\mathbf{x})=ax_{0}^{2}+bx_{1}^{2}+cx_{0}x_{1}+d$$
///
/// where $a$, $b$, $c$, and $d$ are runtime parameters. Compare the result against the true
/// gradient of
///
/// $$\nabla f=\begin{bmatrix}2ax_{0}+cx_{1}\\\\2bx_{1}+cx_{0}\end{bmatrix}$$
///
/// ```
/// use linalg_traits::{Scalar, Vector};
/// use numtest::*;
///
/// use numdiff::{get_gradient, Dual, DualVector};
///
/// // Define the function, f(x).
/// fn f<S: Scalar, V: Vector<S>>(x: &V, p: &[f64]) -> S {
/// let a = S::new(p[0]);
/// let b = S::new(p[1]);
/// let c = S::new(p[2]);
/// let d = S::new(p[3]);
/// a * x.vget(0).powi(2) + b * x.vget(1).powi(2) + c * x.vget(0) * x.vget(1) + d
/// }
///
/// // Parameter vector.
/// let a = 2.0;
/// let b = 1.5;
/// let c = 0.8;
/// let d = -3.0;
/// let p = [a, b, c, d];
///
/// // Evaluation point.
/// let x0 = vec![1.0, -2.0];
///
/// // Autogenerate the gradient function.
/// get_gradient!(f, g);
///
/// // True gradient function.
/// let g_true = |x: &Vec<f64>| vec![2.0 * a * x[0] + c * x[1], 2.0 * b * x[1] + c * x[0]];
///
/// // Compute the gradient using both the automatically generated gradient function and the true
/// // gradient function, and compare the results.
/// let g_eval: Vec<f64> = g(&x0, &p);
/// let g_eval_true: Vec<f64> = g_true(&x0);
/// assert_arrays_equal_to_decimal!(g_eval, g_eval_true, 15);
/// ```
///
/// ## Example Passing Custom Parameter Types
///
/// Use a custom parameter struct instead of `f64` values.
///
/// ```
/// use linalg_traits::{Scalar, Vector};
/// use numtest::*;
///
/// use numdiff::{get_gradient, Dual, DualVector};
///
/// struct Data {
/// a: f64,
/// b: f64,
/// c: f64,
/// d: f64,
/// }
///
/// // Define the function, f(x).
/// fn f<S: Scalar, V: Vector<S>>(x: &V, p: &Data) -> S {
/// let a = S::new(p.a);
/// let b = S::new(p.b);
/// let c = S::new(p.c);
/// let d = S::new(p.d);
/// a * x.vget(0).powi(2) + b * x.vget(1).powi(2) + c * x.vget(0) * x.vget(1) + d
/// }
///
/// // Runtime parameter struct.
/// let p = Data {
/// a: 2.0,
/// b: 1.5,
/// c: 0.8,
/// d: -3.0,
/// };
///
/// // Evaluation point.
/// let x0 = vec![1.0, -2.0];
///
/// // Autogenerate the gradient function, telling the macro to expect a runtime parameter of type
/// // &Data.
/// get_gradient!(f, g, Data);
///
/// // True gradient function.
/// let g_true = |x: &Vec<f64>| vec![
/// 2.0 * p.a * x[0] + p.c * x[1],
/// 2.0 * p.b * x[1] + p.c * x[0],
/// ];
///
/// // Compute the gradient using both the automatically generated gradient function and the true
/// // gradient function, and compare the results.
/// let g_eval: Vec<f64> = g(&x0, &p);
/// let g_eval_true: Vec<f64> = g_true(&x0);
/// assert_arrays_equal_to_decimal!(g_eval, g_eval_true, 15);
/// ```