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