METHOD
cminpack.h
Go to the documentation of this file.
1 /* Header file for cminpack, by Frederic Devernay.
2  The documentation for all functions can be found in the file
3  minpack-documentation.txt from the distribution, or in the source
4  code of each function. */
5 
6 #ifndef __CMINPACK_H__
7 #define __CMINPACK_H__
8 
9 /* The default floating-point type is "double" for C/C++ and "float" for CUDA,
10  but you can change this by defining one of the following symbols when
11  compiling the library, and before including cminpack.h when using it:
12  __cminpack_long_double__ for long double (requires compiler support)
13  __cminpack_double__ for double
14  __cminpack_float__ for float
15  __cminpack_half__ for half from the OpenEXR library (in this case, you must
16  compile cminpack with a C++ compiler)
17 */
18 #ifdef __cminpack_long_double__
19 #define __cminpack_real__ long double
20 #endif
21 
22 #ifdef __cminpack_double__
23 #define __cminpack_real__ double
24 #endif
25 
26 #ifdef __cminpack_float__
27 #define __cminpack_real__ float
28 #endif
29 
30 #ifdef __cminpack_half__
31 #include <OpenEXR/half.h>
32 #define __cminpack_real__ half
33 #endif
34 
35 #ifdef __cplusplus
36 extern "C" {
37 
38 #endif /* __cplusplus */
39 
40 /* Cmake will define cminpack_EXPORTS on Windows when it
41 configures to build a shared library. If you are going to use
42 another build system on windows or create the visual studio
43 projects by hand you need to define cminpack_EXPORTS when
44 building a DLL on windows.
45 */
46 #if defined (__GNUC__)
47 #define CMINPACK_DECLSPEC_EXPORT __declspec(__dllexport__)
48 #define CMINPACK_DECLSPEC_IMPORT __declspec(__dllimport__)
49 #endif
50 // #if defined (_MSC_VER) || defined (__BORLANDC__)
51 // #define CMINPACK_DECLSPEC_EXPORT __declspec(dllexport)
52 // #define CMINPACK_DECLSPEC_IMPORT __declspec(dllimport)
53 // #endif
54 // #ifdef __WATCOMC__
55 // #define CMINPACK_DECLSPEC_EXPORT __export
56 // #define CMINPACK_DECLSPEC_IMPORT __import
57 // #endif
58 // #ifdef __IBMC__
59 // #define CMINPACK_DECLSPEC_EXPORT _Export
60 // #define CMINPACK_DECLSPEC_IMPORT _Import
61 // #endif
62 
63 // #if !defined(CMINPACK_NO_DLL) && (defined(__WIN32__) || defined(WIN32) || defined (_WIN32))
64 // #if defined(cminpack_EXPORTS) || defined(CMINPACK_EXPORTS) || defined(CMINPACK_DLL_EXPORTS)
65 // #define CMINPACK_EXPORT CMINPACK_DECLSPEC_EXPORT
66 // #else
67 // #define CMINPACK_EXPORT CMINPACK_DECLSPEC_IMPORT
68 // #endif /* cminpack_EXPORTS */
69 // #else /* defined (_WIN32) */
70 #define CMINPACK_EXPORT
71 // #endif
72 
73 
74 #define __cminpack_attr__ //__host__ __device__
75 #ifndef __cminpack_real__
76 #define __cminpack_double__
77 #define __cminpack_real__ double
78 #endif
79 
80 #define __cminpack_type_fcn_nn__ typedef int (*cminpack_func_nn)
81 #define __cminpack_type_fcnder_nn__ typedef int (*cminpack_funcder_nn)
82 #define __cminpack_type_fcn_mn__ typedef int (*cminpack_func_mn)
83 #define __cminpack_type_fcnder_mn__ typedef int (*cminpack_funcder_mn)
84 #define __cminpack_type_fcnderstr_mn__ typedef int (*cminpack_funcderstr_mn)
85 #define __cminpack_decl_fcn_nn__ cminpack_func_nn fcn_nn,
86 #define __cminpack_decl_fcnder_nn__ cminpack_funcder_nn fcnder_nn,
87 #define __cminpack_decl_fcn_mn__ cminpack_func_mn fcn_mn,
88 #define __cminpack_decl_fcnder_mn__ cminpack_funcder_mn fcnder_mn,
89 #define __cminpack_decl_fcnderstr_mn__ cminpack_funcderstr_mn fcnderstr_mn,
90 #define __cminpack_param_fcn_nn__ fcn_nn,
91 #define __cminpack_param_fcnder_nn__ fcnder_nn,
92 #define __cminpack_param_fcn_mn__ fcn_mn,
93 #define __cminpack_param_fcnder_mn__ fcnder_mn,
94 #define __cminpack_param_fcnderstr_mn__ fcnderstr_mn,
95 
96 #ifdef __cminpack_double__
97 #define __cminpack_func__(func) func
98 #endif
99 
100 #ifdef __cminpack_long_double__
101 #define __cminpack_func__(func) ld ## func
102 #endif
103 
104 #ifdef __cminpack_float__
105 #define __cminpack_func__(func) s ## func
106 #endif
107 
108 #ifdef __cminpack_half__
109 #define __cminpack_func__(func) h ## func
110 #endif
111 
112 /* Declarations for minpack */
113 
114 /* Function types: */
115 /* The first argument can be used to store extra function parameters, thus */
116 /* avoiding the use of global variables. */
117 /* the iflag parameter is input-only (with respect to the FORTRAN */
118 /* version), the output iflag value is the return value of the function. */
119 /* If iflag=0, the function shoulkd just print the current values (see */
120 /* the nprint parameters below). */
121 
122 /* for hybrd1 and hybrd: */
123 /* calculate the functions at x and */
124 /* return this vector in fvec. */
125 /* return a negative value to terminate hybrd1/hybrd */
126 __cminpack_type_fcn_nn__(void *p, int n, const __cminpack_real__ *x, __cminpack_real__ *fvec, int iflag );
127 
128 /* for hybrj1 and hybrj */
129 /* if iflag = 1 calculate the functions at x and */
130 /* return this vector in fvec. do not alter fjac. */
131 /* if iflag = 2 calculate the jacobian at x and */
132 /* return this matrix in fjac. do not alter fvec. */
133 /* return a negative value to terminate hybrj1/hybrj */
135  int ldfjac, int iflag );
136 
137 
138 /* for lmdif1 and lmdif */
139 /* calculate the functions at x and */
140 /* return this vector in fvec. */
141 /* if iflag = 1 the result is used to compute the residuals. */
142 /* if iflag = 2 the result is used to compute the Jacobian by finite differences. */
143 /* Jacobian computation requires exactly n function calls with iflag = 2. */
144 /* return a negative value to terminate lmdif1/lmdif */
145 __cminpack_type_fcn_mn__(void *p, int m, int n, const __cminpack_real__ *x, __cminpack_real__ *fvec,
146  int iflag );
147 
148 /* for lmder1 and lmder */
149 /* if iflag = 1 calculate the functions at x and */
150 /* return this vector in fvec. do not alter fjac. */
151 /* if iflag = 2 calculate the jacobian at x and */
152 /* return this matrix in fjac. do not alter fvec. */
153 /* return a negative value to terminate lmder1/lmder */
154 __cminpack_type_fcnder_mn__(void *p, int m, int n, const __cminpack_real__ *x, __cminpack_real__ *fvec,
155  __cminpack_real__ *fjac, int ldfjac, int iflag );
156 
157 /* for lmstr1 and lmstr */
158 /* if iflag = 1 calculate the functions at x and */
159 /* return this vector in fvec. */
160 /* if iflag = i calculate the (i-1)-st row of the */
161 /* jacobian at x and return this vector in fjrow. */
162 /* return a negative value to terminate lmstr1/lmstr */
163 __cminpack_type_fcnderstr_mn__(void *p, int m, int n, const __cminpack_real__ *x, __cminpack_real__ *fvec,
164  __cminpack_real__ *fjrow, int iflag );
165 
166 
167 
168 
169 
170 
171 /* MINPACK functions: */
172 /* the info parameter was removed from most functions: the return */
173 /* value of the function is used instead. */
174 /* The argument 'p' can be used to store extra function parameters, thus */
175 /* avoiding the use of global variables. You can also think of it as a */
176 /* 'this' pointer a la C++. */
177 
178 /* find a zero of a system of N nonlinear functions in N variables by
179  a modification of the Powell hybrid method (Jacobian calculated by
180  a forward-difference approximation) */
183  void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ tol,
184  __cminpack_real__ *wa, int lwa );
185 
186 /* find a zero of a system of N nonlinear functions in N variables by
187  a modification of the Powell hybrid method (Jacobian calculated by
188  a forward-difference approximation, more general). */
191  void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ xtol, int maxfev,
192  int ml, int mu, __cminpack_real__ epsfcn, __cminpack_real__ *diag, int mode,
193  __cminpack_real__ factor, int nprint, int *nfev,
194  __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ *r, int lr, __cminpack_real__ *qtf,
196 
197 /* find a zero of a system of N nonlinear functions in N variables by
198  a modification of the Powell hybrid method (user-supplied Jacobian) */
201  __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ tol,
202  __cminpack_real__ *wa, int lwa );
203 
204 /* find a zero of a system of N nonlinear functions in N variables by
205  a modification of the Powell hybrid method (user-supplied Jacobian,
206  more general) */
209  __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ xtol,
210  int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor,
211  int nprint, int *nfev, int *njev, __cminpack_real__ *r,
212  int lr, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2,
214 
215 /* minimize the sum of the squares of nonlinear functions in N
216  variables by a modification of the Levenberg-Marquardt algorithm
217  (Jacobian calculated by a forward-difference approximation) */
220  void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ tol,
221  int *iwa, __cminpack_real__ *wa, int lwa );
222 
223 /* minimize the sum of the squares of nonlinear functions in N
224  variables by a modification of the Levenberg-Marquardt algorithm
225  (Jacobian calculated by a forward-difference approximation, more
226  general) */
229  void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ ftol,
230  __cminpack_real__ xtol, __cminpack_real__ gtol, int maxfev, __cminpack_real__ epsfcn,
231  __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint,
232  int *nfev, __cminpack_real__ *fjac, int ldfjac, int *ipvt,
234  __cminpack_real__ *wa4 );
235 
236 /* minimize the sum of the squares of nonlinear functions in N
237  variables by a modification of the Levenberg-Marquardt algorithm
238  (user-supplied Jacobian) */
241  void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac,
242  int ldfjac, __cminpack_real__ tol, int *ipvt,
243  __cminpack_real__ *wa, int lwa );
244 
245 /* minimize the sum of the squares of nonlinear functions in N
246  variables by a modification of the Levenberg-Marquardt algorithm
247  (user-supplied Jacobian, more general) */
250  void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac,
251  int ldfjac, __cminpack_real__ ftol, __cminpack_real__ xtol, __cminpack_real__ gtol,
252  int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor,
253  int nprint, int *nfev, int *njev, int *ipvt,
255  __cminpack_real__ *wa4 );
256 
257 /* minimize the sum of the squares of nonlinear functions in N
258  variables by a modification of the Levenberg-Marquardt algorithm
259  (user-supplied Jacobian, minimal storage) */
262  __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac,
263  __cminpack_real__ tol, int *ipvt, __cminpack_real__ *wa, int lwa );
264 
265 /* minimize the sum of the squares of nonlinear functions in N
266  variables by a modification of the Levenberg-Marquardt algorithm
267  (user-supplied Jacobian, minimal storage, more general) */
270  int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac,
271  int ldfjac, __cminpack_real__ ftol, __cminpack_real__ xtol, __cminpack_real__ gtol,
272  int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor,
273  int nprint, int *nfev, int *njev, int *ipvt,
275  __cminpack_real__ *wa4 );
276 
279  int ldfjac, __cminpack_real__ *xp, __cminpack_real__ *fvecp, int mode,
280  __cminpack_real__ *err );
281 
284 
287 
288 /* compute a forward-difference approximation to the m by n jacobian
289  matrix associated with a specified problem of m functions in n
290  variables. */
293  void *p, int m, int n, __cminpack_real__ *x, const __cminpack_real__ *fvec, __cminpack_real__ *fjac,
294  int ldfjac, __cminpack_real__ epsfcn, __cminpack_real__ *wa);
295 
296 /* compute a forward-difference approximation to the n by n jacobian
297  matrix associated with a specified problem of n functions in n
298  variables. if the jacobian has a banded form, then function
299  evaluations are saved by only approximating the nonzero terms. */
302  void *p, int n, __cminpack_real__ *x, const __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac,
303  int ml, int mu, __cminpack_real__ epsfcn, __cminpack_real__ *wa1,
304  __cminpack_real__ *wa2);
305 
306 /* compute inverse(JtJ) after a run of lmdif or lmder. The covariance matrix is obtained
307  by scaling the result by enorm(y)**2/(m-n). If JtJ is singular and k = rank(J), the
308  pseudo-inverse is computed, and the result has to be scaled by enorm(y)**2/(m-k). */
311  const int *ipvt, __cminpack_real__ tol, __cminpack_real__ *wa);
312 
313 /* covar1 estimates the variance-covariance matrix:
314  C = sigma**2 (JtJ)**+
315  where (JtJ)**+ is the inverse of JtJ or the pseudo-inverse of JtJ (in case J does not have full rank),
316  and sigma**2 = fsumsq / (m - k)
317  where fsumsq is the residual sum of squares and k is the rank of J.
318  The function returns 0 if J has full rank, else the rank of J.
319 */
321 int CMINPACK_EXPORT __cminpack_func__(covar1)(int m, int n, __cminpack_real__ fsumsq, __cminpack_real__ *r, int ldr,
322  const int *ipvt, __cminpack_real__ tol, __cminpack_real__ *wa);
323 
324 /* internal MINPACK subroutines */
326 void __cminpack_func__(dogleg)(int n, const __cminpack_real__ *r, int lr,
327  const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ delta, __cminpack_real__ *x,
330 void __cminpack_func__(qrfac)(int m, int n, __cminpack_real__ *a, int
331  lda, int pivot, int *ipvt, int lipvt, __cminpack_real__ *rdiag,
332  __cminpack_real__ *acnorm, __cminpack_real__ *wa);
334 void __cminpack_func__(qrsolv)(int n, __cminpack_real__ *r, int ldr,
335  const int *ipvt, const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ *x,
338 void __cminpack_func__(qform)(int m, int n, __cminpack_real__ *q, int
339  ldq, __cminpack_real__ *wa);
341 void __cminpack_func__(r1updt)(int m, int n, __cminpack_real__ *s, int
342  ls, const __cminpack_real__ *u, __cminpack_real__ *v, __cminpack_real__ *w, int *sing);
344 void __cminpack_func__(r1mpyq)(int m, int n, __cminpack_real__ *a, int
345  lda, const __cminpack_real__ *v, const __cminpack_real__ *w);
347 void __cminpack_func__(lmpar)(int n, __cminpack_real__ *r, int ldr,
348  const int *ipvt, const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ delta,
350  __cminpack_real__ *wa2);
352 void __cminpack_func__(rwupdt)(int n, __cminpack_real__ *r, int ldr,
354  __cminpack_real__ *sin);
355 #ifdef __cplusplus
356 }
357 #endif /* __cplusplus */
358 
359 
360 #endif /* __CMINPACK_H__ */
__cminpack_attr__ void __cminpack_func__() r1updt(int m, int n, __cminpack_real__ *s, int ls, const __cminpack_real__ *u, __cminpack_real__ *v, __cminpack_real__ *w, int *sing)
#define __cminpack_type_fcnder_mn__
Definition: cminpack.h:83
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() hybrd1(__cminpack_decl_fcn_nn__ void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ tol, __cminpack_real__ *wa, int lwa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmdif(__cminpack_decl_fcn_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ ftol, __cminpack_real__ xtol, __cminpack_real__ gtol, int maxfev, __cminpack_real__ epsfcn, __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint, int *nfev, __cminpack_real__ *fjac, int ldfjac, int *ipvt, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2, __cminpack_real__ *wa3, __cminpack_real__ *wa4)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() hybrd(__cminpack_decl_fcn_nn__ void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ xtol, int maxfev, int ml, int mu, __cminpack_real__ epsfcn, __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint, int *nfev, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ *r, int lr, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2, __cminpack_real__ *wa3, __cminpack_real__ *wa4)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmder1(__cminpack_decl_fcnder_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ tol, int *ipvt, __cminpack_real__ *wa, int lwa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() covar1(int m, int n, __cminpack_real__ fsumsq, __cminpack_real__ *r, int ldr, const int *ipvt, __cminpack_real__ tol, __cminpack_real__ *wa)
#define __cminpack_decl_fcnder_mn__
Definition: cminpack.h:88
#define __cminpack_decl_fcn_nn__
Definition: cminpack.h:85
#define __cminpack_attr__
Definition: cminpack.h:74
#define __cminpack_type_fcnderstr_mn__
Definition: cminpack.h:84
#define CMINPACK_EXPORT
Definition: cminpack.h:70
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmder(__cminpack_decl_fcnder_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ ftol, __cminpack_real__ xtol, __cminpack_real__ gtol, int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint, int *nfev, int *njev, int *ipvt, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2, __cminpack_real__ *wa3, __cminpack_real__ *wa4)
__cminpack_attr__ void CMINPACK_EXPORT __cminpack_func__() covar(int n, __cminpack_real__ *r, int ldr, const int *ipvt, __cminpack_real__ tol, __cminpack_real__ *wa)
__cminpack_attr__ void __cminpack_func__() r1mpyq(int m, int n, __cminpack_real__ *a, int lda, const __cminpack_real__ *v, const __cminpack_real__ *w)
#define __cminpack_decl_fcn_mn__
Definition: cminpack.h:87
__cminpack_attr__ void __cminpack_func__() rwupdt(int n, __cminpack_real__ *r, int ldr, const __cminpack_real__ *w, __cminpack_real__ *b, __cminpack_real__ *alpha, __cminpack_real__ *cos, __cminpack_real__ *sin)
#define __cminpack_decl_fcnder_nn__
Definition: cminpack.h:86
#define __cminpack_type_fcn_nn__
Definition: cminpack.h:80
__cminpack_attr__ void __cminpack_func__() qrsolv(int n, __cminpack_real__ *r, int ldr, const int *ipvt, const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ *x, __cminpack_real__ *sdiag, __cminpack_real__ *wa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() fdjac1(__cminpack_decl_fcn_nn__ void *p, int n, __cminpack_real__ *x, const __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, int ml, int mu, __cminpack_real__ epsfcn, __cminpack_real__ *wa1, __cminpack_real__ *wa2)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmstr1(__cminpack_decl_fcnderstr_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ tol, int *ipvt, __cminpack_real__ *wa, int lwa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() hybrj1(__cminpack_decl_fcnder_nn__ void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ tol, __cminpack_real__ *wa, int lwa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() hybrj(__cminpack_decl_fcnder_nn__ void *p, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ xtol, int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint, int *nfev, int *njev, __cminpack_real__ *r, int lr, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2, __cminpack_real__ *wa3, __cminpack_real__ *wa4)
#define __cminpack_type_fcn_mn__
Definition: cminpack.h:82
__cminpack_attr__ void CMINPACK_EXPORT __cminpack_func__() chkder(int m, int n, const __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ *xp, __cminpack_real__ *fvecp, int mode, __cminpack_real__ *err)
__cminpack_attr__ void __cminpack_func__() lmpar(int n, __cminpack_real__ *r, int ldr, const int *ipvt, const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ delta, __cminpack_real__ *par, __cminpack_real__ *x, __cminpack_real__ *sdiag, __cminpack_real__ *wa1, __cminpack_real__ *wa2)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmstr(__cminpack_decl_fcnderstr_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ ftol, __cminpack_real__ xtol, __cminpack_real__ gtol, int maxfev, __cminpack_real__ *diag, int mode, __cminpack_real__ factor, int nprint, int *nfev, int *njev, int *ipvt, __cminpack_real__ *qtf, __cminpack_real__ *wa1, __cminpack_real__ *wa2, __cminpack_real__ *wa3, __cminpack_real__ *wa4)
#define __cminpack_func__(func)
Definition: cminpack.h:97
#define __cminpack_decl_fcnderstr_mn__
Definition: cminpack.h:89
#define __cminpack_type_fcnder_nn__
Definition: cminpack.h:81
__cminpack_attr__ __cminpack_real__ CMINPACK_EXPORT __cminpack_func__() dpmpar(int i)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() lmdif1(__cminpack_decl_fcn_mn__ void *p, int m, int n, __cminpack_real__ *x, __cminpack_real__ *fvec, __cminpack_real__ tol, int *iwa, __cminpack_real__ *wa, int lwa)
#define __cminpack_real__
Definition: cminpack.h:77
__cminpack_attr__ void __cminpack_func__() qform(int m, int n, __cminpack_real__ *q, int ldq, __cminpack_real__ *wa)
__cminpack_attr__ __cminpack_real__ CMINPACK_EXPORT __cminpack_func__() enorm(int n, const __cminpack_real__ *x)
__cminpack_attr__ void __cminpack_func__() dogleg(int n, const __cminpack_real__ *r, int lr, const __cminpack_real__ *diag, const __cminpack_real__ *qtb, __cminpack_real__ delta, __cminpack_real__ *x, __cminpack_real__ *wa1, __cminpack_real__ *wa2)
__cminpack_attr__ void __cminpack_func__() qrfac(int m, int n, __cminpack_real__ *a, int lda, int pivot, int *ipvt, int lipvt, __cminpack_real__ *rdiag, __cminpack_real__ *acnorm, __cminpack_real__ *wa)
__cminpack_attr__ int CMINPACK_EXPORT __cminpack_func__() fdjac2(__cminpack_decl_fcn_mn__ void *p, int m, int n, __cminpack_real__ *x, const __cminpack_real__ *fvec, __cminpack_real__ *fjac, int ldfjac, __cminpack_real__ epsfcn, __cminpack_real__ *wa)