METHOD
minpack.h
Go to the documentation of this file.
1 #ifndef __MINPACK_H__
2 #define __MINPACK_H__
3 
4 #include "cminpack.h"
5 
6 /* The default floating-point type is "double" for C/C++ and "float" for CUDA,
7  but you can change this by defining one of the following symbols when
8  compiling the library, and before including cminpack.h when using it:
9  __cminpack_long_double__ for long double (requires compiler support)
10  __cminpack_double__ for double
11  __cminpack_float__ for float
12  __cminpack_half__ for half from the OpenEXR library (in this case, you must
13  compile cminpack with a C++ compiler)
14 */
15 #ifdef __cminpack_double__
16 #define __minpack_func__(func) func ## _
17 #endif
18 
19 #ifdef __cminpack_long_double__
20 #define __minpack_func__(func) ld ## func ## _
21 #endif
22 
23 #ifdef __cminpack_float__
24 #define __minpack_func__(func) s ## func ## _
25 #endif
26 
27 #ifdef __cminpack_half__
28 #define __minpack_func__(func) h ## func ## _
29 #endif
30 
31 #ifdef __cplusplus
32 extern "C" {
33 #endif /* __cplusplus */
34 
35 #define MINPACK_EXPORT CMINPACK_EXPORT
36 
37 #define __minpack_real__ __cminpack_real__
38 #define __minpack_attr__ __cminpack_attr__
39 #if defined(__CUDA_ARCH__) || defined(__CUDACC__)
40 #define __minpack_type_fcn_nn__ __minpack_attr__ void fcn_nn
41 #define __minpack_type_fcnder_nn__ __minpack_attr__ void fcnder_nn
42 #define __minpack_type_fcn_mn__ __minpack_attr__ void fcn_mn
43 #define __minpack_type_fcnder_mn__ __minpack_attr__ void fcnder_mn
44 #define __minpack_type_fcnderstr_mn__ __minpack_attr__ void fcnderstr_mn
45 #define __minpack_decl_fcn_nn__
46 #define __minpack_decl_fcnder_nn__
47 #define __minpack_decl_fcn_mn__
48 #define __minpack_decl_fcnder_mn__
49 #define __minpack_decl_fcnderstr_mn__
50 #define __minpack_param_fcn_nn__
51 #define __minpack_param_fcnder_nn__
52 #define __minpack_param_fcn_mn__
53 #define __minpack_param_fcnder_mn__
54 #define __minpack_param_fcnderstr_mn__
55 #else
56 #define __minpack_type_fcn_nn__ typedef void (*minpack_func_nn)
57 #define __minpack_type_fcnder_nn__ typedef void (*minpack_funcder_nn)
58 #define __minpack_type_fcn_mn__ typedef void (*minpack_func_mn)
59 #define __minpack_type_fcnder_mn__ typedef void (*minpack_funcder_mn)
60 #define __minpack_type_fcnderstr_mn__ typedef void (*minpack_funcderstr_mn)
61 #define __minpack_decl_fcn_nn__ minpack_func_nn fcn_nn,
62 #define __minpack_decl_fcnder_nn__ minpack_funcder_nn fcnder_nn,
63 #define __minpack_decl_fcn_mn__ minpack_func_mn fcn_mn,
64 #define __minpack_decl_fcnder_mn__ minpack_funcder_mn fcnder_mn,
65 #define __minpack_decl_fcnderstr_mn__ minpack_funcderstr_mn fcnderstr_mn,
66 #define __minpack_param_fcn_nn__ fcn_nn,
67 #define __minpack_param_fcnder_nn__ fcnder_nn,
68 #define __minpack_param_fcn_mn__ fcn_mn,
69 #define __minpack_param_fcnder_mn__ fcnder_mn,
70 #define __minpack_param_fcnderstr_mn__ fcnderstr_mn,
71 #endif
72 #undef __cminpack_type_fcn_nn__
73 #undef __cminpack_type_fcnder_nn__
74 #undef __cminpack_type_fcn_mn__
75 #undef __cminpack_type_fcnder_mn__
76 #undef __cminpack_type_fcnderstr_mn__
77 #undef __cminpack_decl_fcn_nn__
78 #undef __cminpack_decl_fcnder_nn__
79 #undef __cminpack_decl_fcn_mn__
80 #undef __cminpack_decl_fcnder_mn__
81 #undef __cminpack_decl_fcnderstr_mn__
82 #undef __cminpack_param_fcn_nn__
83 #undef __cminpack_param_fcnder_nn__
84 #undef __cminpack_param_fcn_mn__
85 #undef __cminpack_param_fcnder_mn__
86 #undef __cminpack_param_fcnderstr_mn__
87 
88 /* Declarations for minpack */
89 
90 /* Function types: */
91 /* the iflag parameter is input-only (with respect to the FORTRAN */
92 /* version), the output iflag value is the return value of the function. */
93 /* If iflag=0, the function shoulkd just print the current values (see */
94 /* the nprint parameters below). */
95 
96 /* for hybrd1 and hybrd: */
97 /* calculate the functions at x and */
98 /* return this vector in fvec. */
99 /* return a negative value to terminate hybrd1/hybrd */
100 __minpack_type_fcn_nn__(const int *n, const __minpack_real__ *x, __minpack_real__ *fvec, int *iflag );
101 
102 /* for hybrj1 and hybrj */
103 /* if iflag = 1 calculate the functions at x and */
104 /* return this vector in fvec. do not alter fjac. */
105 /* if iflag = 2 calculate the jacobian at x and */
106 /* return this matrix in fjac. do not alter fvec. */
107 /* return a negative value to terminate hybrj1/hybrj */
109  const int *ldfjac, int *iflag );
110 
111 /* for lmdif1 and lmdif */
112 /* calculate the functions at x and */
113 /* return this vector in fvec. */
114 /* if iflag = 1 the result is used to compute the residuals. */
115 /* if iflag = 2 the result is used to compute the Jacobian by finite differences. */
116 /* Jacobian computation requires exactly n function calls with iflag = 2. */
117 /* return a negative value to terminate lmdif1/lmdif */
118 __minpack_type_fcn_mn__(const int *m, const int *n, const __minpack_real__ *x, __minpack_real__ *fvec,
119  int *iflag );
120 
121 /* for lmder1 and lmder */
122 /* if iflag = 1 calculate the functions at x and */
123 /* return this vector in fvec. do not alter fjac. */
124 /* if iflag = 2 calculate the jacobian at x and */
125 /* return this matrix in fjac. do not alter fvec. */
126 /* return a negative value to terminate lmder1/lmder */
127 __minpack_type_fcnder_mn__(const int *m, const int *n, const __minpack_real__ *x, __minpack_real__ *fvec,
128  __minpack_real__ *fjac, const int *ldfjac, int *iflag );
129 
130 /* for lmstr1 and lmstr */
131 /* if iflag = 1 calculate the functions at x and */
132 /* return this vector in fvec. */
133 /* if iflag = i calculate the (i-1)-st row of the */
134 /* jacobian at x and return this vector in fjrow. */
135 /* return a negative value to terminate lmstr1/lmstr */
136 __minpack_type_fcnderstr_mn__(const int *m, const int *n, const __minpack_real__ *x, __minpack_real__ *fvec,
137  __minpack_real__ *fjrow, int *iflag );
138 
139 /* find a zero of a system of N nonlinear functions in N variables by
140  a modification of the Powell hybrid method (Jacobian calculated by
141  a forward-difference approximation) */
144  const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *tol, int *info,
145  __minpack_real__ *wa, const int *lwa );
146 
147 /* find a zero of a system of N nonlinear functions in N variables by
148  a modification of the Powell hybrid method (Jacobian calculated by
149  a forward-difference approximation, more general). */
152  const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *xtol, const int *maxfev,
153  const int *ml, const int *mu, const __minpack_real__ *epsfcn, __minpack_real__ *diag, const int *mode,
154  const __minpack_real__ *factor, const int *nprint, int *info, int *nfev,
155  __minpack_real__ *fjac, const int *ldfjac, __minpack_real__ *r, const int *lr, __minpack_real__ *qtf,
157 
158 /* find a zero of a system of N nonlinear functions in N variables by
159  a modification of the Powell hybrid method (user-supplied Jacobian) */
162  __minpack_real__ *fvec, __minpack_real__ *fjec, const int *ldfjac, const __minpack_real__ *tol,
163  int *info, __minpack_real__ *wa, const int *lwa );
164 
165 /* find a zero of a system of N nonlinear functions in N variables by
166  a modification of the Powell hybrid method (user-supplied Jacobian,
167  more general) */
170  __minpack_real__ *fvec, __minpack_real__ *fjec, const int *ldfjac, const __minpack_real__ *xtol,
171  const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor,
172  const int *nprint, int *info, int *nfev, int *njev, __minpack_real__ *r,
173  const int *lr, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2,
174  __minpack_real__ *wa3, __minpack_real__ *wa4 );
175 
176 /* minimize the sum of the squares of nonlinear functions in N
177  variables by a modification of the Levenberg-Marquardt algorithm
178  (Jacobian calculated by a forward-difference approximation) */
181  const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *tol,
182  int *info, int *iwa, __minpack_real__ *wa, const int *lwa );
183 
184 /* minimize the sum of the squares of nonlinear functions in N
185  variables by a modification of the Levenberg-Marquardt algorithm
186  (Jacobian calculated by a forward-difference approximation, more
187  general) */
190  const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *ftol,
191  const __minpack_real__ *xtol, const __minpack_real__ *gtol, const int *maxfev, const __minpack_real__ *epsfcn,
192  __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint,
193  int *info, int *nfev, __minpack_real__ *fjac, const int *ldfjac, int *ipvt,
195  __minpack_real__ *wa4 );
196 
197 /* minimize the sum of the squares of nonlinear functions in N
198  variables by a modification of the Levenberg-Marquardt algorithm
199  (user-supplied Jacobian) */
202  const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac,
203  const int *ldfjac, const __minpack_real__ *tol, int *info, int *ipvt,
204  __minpack_real__ *wa, const int *lwa );
205 
206 /* minimize the sum of the squares of nonlinear functions in N
207  variables by a modification of the Levenberg-Marquardt algorithm
208  (user-supplied Jacobian, more general) */
211  const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac,
212  const int *ldfjac, const __minpack_real__ *ftol, const __minpack_real__ *xtol, const __minpack_real__ *gtol,
213  const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor,
214  const int *nprint, int *info, int *nfev, int *njev, int *ipvt,
216  __minpack_real__ *wa4 );
217 
218 /* minimize the sum of the squares of nonlinear functions in N
219  variables by a modification of the Levenberg-Marquardt algorithm
220  (user-supplied Jacobian, minimal storage) */
222 void MINPACK_EXPORT __minpack_func__(lmstr1)( __minpack_decl_fcnderstr_mn__ const int *m, const int *n,
223  __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac,
224  const __minpack_real__ *tol, int *info, int *ipvt, __minpack_real__ *wa, const int *lwa );
225 
226 /* minimize the sum of the squares of nonlinear functions in N
227  variables by a modification of the Levenberg-Marquardt algorithm
228  (user-supplied Jacobian, minimal storage, more general) */
231  const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac,
232  const int *ldfjac, const __minpack_real__ *ftol, const __minpack_real__ *xtol, const __minpack_real__ *gtol,
233  const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor,
234  const int *nprint, int *info, int *nfev, int *njev, int *ipvt,
236  __minpack_real__ *wa4 );
237 
239 void MINPACK_EXPORT __minpack_func__(chkder)( const int *m, const int *n, const __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjec,
240  const int *ldfjac, __minpack_real__ *xp, __minpack_real__ *fvecp, const int *mode,
241  __minpack_real__ *err );
242 
245 
248 
249 /* compute a forward-difference approximation to the m by n jacobian
250  matrix associated with a specified problem of m functions in n
251  variables. */
254  const int *m, const int *n, __minpack_real__ *x, const __minpack_real__ *fvec, __minpack_real__ *fjac,
255  const int *ldfjac, int *iflag, const __minpack_real__ *epsfcn, __minpack_real__ *wa);
256 
257 /* compute a forward-difference approximation to the n by n jacobian
258  matrix associated with a specified problem of n functions in n
259  variables. if the jacobian has a banded form, then function
260  evaluations are saved by only approximating the nonzero terms. */
263  const int *n, __minpack_real__ *x, const __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac,
264  int *iflag, const int *ml, const int *mu, const __minpack_real__ *epsfcn, __minpack_real__ *wa1,
265  __minpack_real__ *wa2);
266 
267 /* internal MINPACK subroutines */
269 void __minpack_func__(dogleg)(const int *n, const __minpack_real__ *r, const int *lr,
270  const __minpack_real__ *diag, const __minpack_real__ *qtb, const __minpack_real__ *delta, __minpack_real__ *x,
273 void __minpack_func__(qrfac)(const int *m, const int *n, __minpack_real__ *a, const int *
274  lda, const int *pivot, int *ipvt, const int *lipvt, __minpack_real__ *rdiag,
275  __minpack_real__ *acnorm, __minpack_real__ *wa);
277 void __minpack_func__(qrsolv)(const int *n, __minpack_real__ *r, const int *ldr,
278  const int *ipvt, const __minpack_real__ *diag, const __minpack_real__ *qtb, __minpack_real__ *x,
279  __minpack_real__ *sdiag, __minpack_real__ *wa);
281 void __minpack_func__(qform)(const int *m, const int *n, __minpack_real__ *q, const int *
282  ldq, __minpack_real__ *wa);
284 void __minpack_func__(r1updt)(const int *m, const int *n, __minpack_real__ *s, const int *
285  ls, const __minpack_real__ *u, __minpack_real__ *v, __minpack_real__ *w, int *sing);
287 void __minpack_func__(r1mpyq)(const int *m, const int *n, __minpack_real__ *a, const int *
288  lda, const __minpack_real__ *v, const __minpack_real__ *w);
290 void __minpack_func__(lmpar)(const int *n, __minpack_real__ *r, const int *ldr,
291  const int *ipvt, const __minpack_real__ *diag, const __minpack_real__ *qtb, const __minpack_real__ *delta,
293  __minpack_real__ *wa2);
295 void __minpack_func__(rwupdt)(const int *n, __minpack_real__ *r, const int *ldr,
297  __minpack_real__ *sin);
299 void __minpack_func__(covar)(const int *n, __minpack_real__ *r, const int *ldr,
300  const int *ipvt, const __minpack_real__ *tol, __minpack_real__ *wa);
301 #ifdef __cplusplus
302 }
303 #endif /* __cplusplus */
304 
305 
306 #endif /* __MINPACK_H__ */
#define __minpack_decl_fcnder_nn__
Definition: minpack.h:62
__minpack_attr__ void __minpack_func__() dogleg(const int *n, const __minpack_real__ *r, const int *lr, const __minpack_real__ *diag, const __minpack_real__ *qtb, const __minpack_real__ *delta, __minpack_real__ *x, __minpack_real__ *wa1, __minpack_real__ *wa2)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() hybrj(__minpack_decl_fcnder_nn__ const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjec, const int *ldfjac, const __minpack_real__ *xtol, const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint, int *info, int *nfev, int *njev, __minpack_real__ *r, const int *lr, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2, __minpack_real__ *wa3, __minpack_real__ *wa4)
#define __minpack_type_fcn_nn__
Definition: minpack.h:56
__minpack_attr__ __minpack_real__ MINPACK_EXPORT __minpack_func__() dpmpar(const int *i)
__minpack_attr__ void __minpack_func__() r1updt(const int *m, const int *n, __minpack_real__ *s, const int *ls, const __minpack_real__ *u, __minpack_real__ *v, __minpack_real__ *w, int *sing)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmder1(__minpack_decl_fcnder_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, const __minpack_real__ *tol, int *info, int *ipvt, __minpack_real__ *wa, const int *lwa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmstr1(__minpack_decl_fcnderstr_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, const __minpack_real__ *tol, int *info, int *ipvt, __minpack_real__ *wa, const int *lwa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmdif(__minpack_decl_fcn_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *ftol, const __minpack_real__ *xtol, const __minpack_real__ *gtol, const int *maxfev, const __minpack_real__ *epsfcn, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint, int *info, int *nfev, __minpack_real__ *fjac, const int *ldfjac, int *ipvt, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2, __minpack_real__ *wa3, __minpack_real__ *wa4)
__minpack_attr__ void __minpack_func__() r1mpyq(const int *m, const int *n, __minpack_real__ *a, const int *lda, const __minpack_real__ *v, const __minpack_real__ *w)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() fdjac1(__minpack_decl_fcn_nn__ const int *n, __minpack_real__ *x, const __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, int *iflag, const int *ml, const int *mu, const __minpack_real__ *epsfcn, __minpack_real__ *wa1, __minpack_real__ *wa2)
__minpack_attr__ void __minpack_func__() rwupdt(const int *n, __minpack_real__ *r, const int *ldr, const __minpack_real__ *w, __minpack_real__ *b, __minpack_real__ *alpha, __minpack_real__ *cos, __minpack_real__ *sin)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmder(__minpack_decl_fcnder_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, const __minpack_real__ *ftol, const __minpack_real__ *xtol, const __minpack_real__ *gtol, const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint, int *info, int *nfev, int *njev, int *ipvt, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2, __minpack_real__ *wa3, __minpack_real__ *wa4)
__minpack_attr__ void __minpack_func__() covar(const int *n, __minpack_real__ *r, const int *ldr, const int *ipvt, const __minpack_real__ *tol, __minpack_real__ *wa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() hybrd1(__minpack_decl_fcn_nn__ const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *tol, int *info, __minpack_real__ *wa, const int *lwa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() hybrd(__minpack_decl_fcn_nn__ const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *xtol, const int *maxfev, const int *ml, const int *mu, const __minpack_real__ *epsfcn, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint, int *info, int *nfev, __minpack_real__ *fjac, const int *ldfjac, __minpack_real__ *r, const int *lr, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2, __minpack_real__ *wa3, __minpack_real__ *wa4)
#define __minpack_attr__
Definition: minpack.h:38
#define __minpack_decl_fcnderstr_mn__
Definition: minpack.h:65
#define __minpack_func__(func)
Definition: minpack.h:16
#define __minpack_type_fcnder_nn__
Definition: minpack.h:57
__minpack_attr__ __minpack_real__ MINPACK_EXPORT __minpack_func__() enorm(const int *n, const __minpack_real__ *x)
#define MINPACK_EXPORT
Definition: minpack.h:35
#define __minpack_decl_fcn_mn__
Definition: minpack.h:63
__minpack_attr__ void __minpack_func__() qform(const int *m, const int *n, __minpack_real__ *q, const int *ldq, __minpack_real__ *wa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() chkder(const int *m, const int *n, const __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjec, const int *ldfjac, __minpack_real__ *xp, __minpack_real__ *fvecp, const int *mode, __minpack_real__ *err)
__minpack_attr__ void __minpack_func__() qrsolv(const int *n, __minpack_real__ *r, const int *ldr, const int *ipvt, const __minpack_real__ *diag, const __minpack_real__ *qtb, __minpack_real__ *x, __minpack_real__ *sdiag, __minpack_real__ *wa)
#define __minpack_real__
Definition: minpack.h:37
__minpack_attr__ void __minpack_func__() lmpar(const int *n, __minpack_real__ *r, const int *ldr, const int *ipvt, const __minpack_real__ *diag, const __minpack_real__ *qtb, const __minpack_real__ *delta, __minpack_real__ *par, __minpack_real__ *x, __minpack_real__ *sdiag, __minpack_real__ *wa1, __minpack_real__ *wa2)
#define __minpack_decl_fcnder_mn__
Definition: minpack.h:64
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmstr(__minpack_decl_fcnderstr_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, const __minpack_real__ *ftol, const __minpack_real__ *xtol, const __minpack_real__ *gtol, const int *maxfev, __minpack_real__ *diag, const int *mode, const __minpack_real__ *factor, const int *nprint, int *info, int *nfev, int *njev, int *ipvt, __minpack_real__ *qtf, __minpack_real__ *wa1, __minpack_real__ *wa2, __minpack_real__ *wa3, __minpack_real__ *wa4)
#define __minpack_type_fcnder_mn__
Definition: minpack.h:59
#define __minpack_decl_fcn_nn__
Definition: minpack.h:61
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() hybrj1(__minpack_decl_fcnder_nn__ const int *n, __minpack_real__ *x, __minpack_real__ *fvec, __minpack_real__ *fjec, const int *ldfjac, const __minpack_real__ *tol, int *info, __minpack_real__ *wa, const int *lwa)
__minpack_attr__ void __minpack_func__() qrfac(const int *m, const int *n, __minpack_real__ *a, const int *lda, const int *pivot, int *ipvt, const int *lipvt, __minpack_real__ *rdiag, __minpack_real__ *acnorm, __minpack_real__ *wa)
#define __minpack_type_fcn_mn__
Definition: minpack.h:58
#define __minpack_type_fcnderstr_mn__
Definition: minpack.h:60
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() lmdif1(__minpack_decl_fcn_mn__ const int *m, const int *n, __minpack_real__ *x, __minpack_real__ *fvec, const __minpack_real__ *tol, int *info, int *iwa, __minpack_real__ *wa, const int *lwa)
__minpack_attr__ void MINPACK_EXPORT __minpack_func__() fdjac2(__minpack_decl_fcn_mn__ const int *m, const int *n, __minpack_real__ *x, const __minpack_real__ *fvec, __minpack_real__ *fjac, const int *ldfjac, int *iflag, const __minpack_real__ *epsfcn, __minpack_real__ *wa)