File Coverage

levmar-2.5/Axb_core.c
Criterion Covered Total %
statement 73 78 93.5
branch 108 120 90.0
condition n/a
subroutine n/a
pod n/a
total 181 198 91.4


line stmt bran cond sub pod time code
1             /////////////////////////////////////////////////////////////////////////////////
2             //
3             // Solution of linear systems involved in the Levenberg - Marquardt
4             // minimization algorithm
5             // Copyright (C) 2004 Manolis Lourakis (lourakis at ics forth gr)
6             // Institute of Computer Science, Foundation for Research & Technology - Hellas
7             // Heraklion, Crete, Greece.
8             //
9             // This program is free software; you can redistribute it and/or modify
10             // it under the terms of the GNU General Public License as published by
11             // the Free Software Foundation; either version 2 of the License, or
12             // (at your option) any later version.
13             //
14             // This program is distributed in the hope that it will be useful,
15             // but WITHOUT ANY WARRANTY; without even the implied warranty of
16             // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17             // GNU General Public License for more details.
18             //
19             /////////////////////////////////////////////////////////////////////////////////
20              
21              
22             /* Solvers for the linear systems Ax=b. Solvers should NOT modify their A & B arguments! */
23              
24              
25             #ifndef LM_REAL // not included by Axb.c
26             #error This file should not be compiled directly!
27             #endif
28              
29              
30             #ifdef LINSOLVERS_RETAIN_MEMORY
31             #define __STATIC__ static
32             #else
33             #define __STATIC__ // empty
34             #endif /* LINSOLVERS_RETAIN_MEMORY */
35              
36             #ifdef HAVE_LAPACK
37              
38             /* prototypes of LAPACK routines */
39              
40             #define GEQRF LM_MK_LAPACK_NAME(geqrf)
41             #define ORGQR LM_MK_LAPACK_NAME(orgqr)
42             #define TRTRS LM_MK_LAPACK_NAME(trtrs)
43             #define POTF2 LM_MK_LAPACK_NAME(potf2)
44             #define POTRF LM_MK_LAPACK_NAME(potrf)
45             #define POTRS LM_MK_LAPACK_NAME(potrs)
46             #define GETRF LM_MK_LAPACK_NAME(getrf)
47             #define GETRS LM_MK_LAPACK_NAME(getrs)
48             #define GESVD LM_MK_LAPACK_NAME(gesvd)
49             #define GESDD LM_MK_LAPACK_NAME(gesdd)
50             #define SYTRF LM_MK_LAPACK_NAME(sytrf)
51             #define SYTRS LM_MK_LAPACK_NAME(sytrs)
52              
53             /* QR decomposition */
54             extern int GEQRF(int *m, int *n, LM_REAL *a, int *lda, LM_REAL *tau, LM_REAL *work, int *lwork, int *info);
55             extern int ORGQR(int *m, int *n, int *k, LM_REAL *a, int *lda, LM_REAL *tau, LM_REAL *work, int *lwork, int *info);
56              
57             /* solution of triangular systems */
58             extern int TRTRS(char *uplo, char *trans, char *diag, int *n, int *nrhs, LM_REAL *a, int *lda, LM_REAL *b, int *ldb, int *info);
59              
60             /* Cholesky decomposition and systems solution */
61             extern int POTF2(char *uplo, int *n, LM_REAL *a, int *lda, int *info);
62             extern int POTRF(char *uplo, int *n, LM_REAL *a, int *lda, int *info); /* block version of dpotf2 */
63             extern int POTRS(char *uplo, int *n, int *nrhs, LM_REAL *a, int *lda, LM_REAL *b, int *ldb, int *info);
64              
65             /* LU decomposition and systems solution */
66             extern int GETRF(int *m, int *n, LM_REAL *a, int *lda, int *ipiv, int *info);
67             extern int GETRS(char *trans, int *n, int *nrhs, LM_REAL *a, int *lda, int *ipiv, LM_REAL *b, int *ldb, int *info);
68              
69             /* Singular Value Decomposition (SVD) */
70             extern int GESVD(char *jobu, char *jobvt, int *m, int *n, LM_REAL *a, int *lda, LM_REAL *s, LM_REAL *u, int *ldu,
71             LM_REAL *vt, int *ldvt, LM_REAL *work, int *lwork, int *info);
72              
73             /* lapack 3.0 new SVD routine, faster than xgesvd().
74             * In case that your version of LAPACK does not include them, use the above two older routines
75             */
76             extern int GESDD(char *jobz, int *m, int *n, LM_REAL *a, int *lda, LM_REAL *s, LM_REAL *u, int *ldu, LM_REAL *vt, int *ldvt,
77             LM_REAL *work, int *lwork, int *iwork, int *info);
78              
79             /* LDLt/UDUt factorization and systems solution */
80             extern int SYTRF(char *uplo, int *n, LM_REAL *a, int *lda, int *ipiv, LM_REAL *work, int *lwork, int *info);
81             extern int SYTRS(char *uplo, int *n, int *nrhs, LM_REAL *a, int *lda, int *ipiv, LM_REAL *b, int *ldb, int *info);
82              
83             /* precision-specific definitions */
84             #define AX_EQ_B_QR LM_ADD_PREFIX(Ax_eq_b_QR)
85             #define AX_EQ_B_QRLS LM_ADD_PREFIX(Ax_eq_b_QRLS)
86             #define AX_EQ_B_CHOL LM_ADD_PREFIX(Ax_eq_b_Chol)
87             #define AX_EQ_B_LU LM_ADD_PREFIX(Ax_eq_b_LU)
88             #define AX_EQ_B_SVD LM_ADD_PREFIX(Ax_eq_b_SVD)
89             #define AX_EQ_B_BK LM_ADD_PREFIX(Ax_eq_b_BK)
90              
91             /*
92             * This function returns the solution of Ax = b
93             *
94             * The function is based on QR decomposition with explicit computation of Q:
95             * If A=Q R with Q orthogonal and R upper triangular, the linear system becomes
96             * Q R x = b or R x = Q^T b.
97             * The last equation can be solved directly.
98             *
99             * A is mxm, b is mx1
100             *
101             * The function returns 0 in case of error, 1 if successful
102             *
103             * This function is often called repetitively to solve problems of identical
104             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
105             * retained between calls and free'd-malloc'ed when not of the appropriate size.
106             * A call with NULL as the first argument forces this memory to be released.
107             */
108             int AX_EQ_B_QR(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
109             {
110             __STATIC__ LM_REAL *buf=NULL;
111             __STATIC__ int buf_sz=0;
112              
113             static int nb=0; /* no __STATIC__ decl. here! */
114              
115             LM_REAL *a, *tau, *r, *work;
116             int a_sz, tau_sz, r_sz, tot_sz;
117             register int i, j;
118             int info, worksz, nrhs=1;
119             register LM_REAL sum;
120              
121             if(!A)
122             #ifdef LINSOLVERS_RETAIN_MEMORY
123             {
124             if(buf) free(buf);
125             buf=NULL;
126             buf_sz=0;
127              
128             return 1;
129             }
130             #else
131             return 1; /* NOP */
132             #endif /* LINSOLVERS_RETAIN_MEMORY */
133            
134             /* calculate required memory size */
135             a_sz=m*m;
136             tau_sz=m;
137             r_sz=m*m; /* only the upper triangular part really needed */
138             if(!nb){
139             LM_REAL tmp;
140              
141             worksz=-1; // workspace query; optimal size is returned in tmp
142             GEQRF((int *)&m, (int *)&m, NULL, (int *)&m, NULL, (LM_REAL *)&tmp, (int *)&worksz, (int *)&info);
143             nb=((int)tmp)/m; // optimal worksize is m*nb
144             }
145             worksz=nb*m;
146             tot_sz=a_sz + tau_sz + r_sz + worksz;
147              
148             #ifdef LINSOLVERS_RETAIN_MEMORY
149             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
150             if(buf) free(buf); /* free previously allocated memory */
151              
152             buf_sz=tot_sz;
153             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
154             if(!buf){
155             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_QR) "() failed!\n");
156             exit(1);
157             }
158             }
159             #else
160             buf_sz=tot_sz;
161             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
162             if(!buf){
163             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_QR) "() failed!\n");
164             exit(1);
165             }
166             #endif /* LINSOLVERS_RETAIN_MEMORY */
167              
168             a=buf;
169             tau=a+a_sz;
170             r=tau+tau_sz;
171             work=r+r_sz;
172              
173             /* store A (column major!) into a */
174             for(i=0; i
175             for(j=0; j
176             a[i+j*m]=A[i*m+j];
177              
178             /* QR decomposition of A */
179             GEQRF((int *)&m, (int *)&m, a, (int *)&m, tau, work, (int *)&worksz, (int *)&info);
180             /* error treatment */
181             if(info!=0){
182             if(info<0){
183             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", GEQRF) " in ", AX_EQ_B_QR) "()\n", -info);
184             exit(1);
185             }
186             else{
187             fprintf(stderr, RCAT(RCAT("Unknown LAPACK error %d for ", GEQRF) " in ", AX_EQ_B_QR) "()\n", info);
188             #ifndef LINSOLVERS_RETAIN_MEMORY
189             free(buf);
190             #endif
191              
192             return 0;
193             }
194             }
195              
196             /* R is stored in the upper triangular part of a; copy it in r so that ORGQR() below won't destroy it */
197             for(i=0; i
198             r[i]=a[i];
199              
200             /* compute Q using the elementary reflectors computed by the above decomposition */
201             ORGQR((int *)&m, (int *)&m, (int *)&m, a, (int *)&m, tau, work, (int *)&worksz, (int *)&info);
202             if(info!=0){
203             if(info<0){
204             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", ORGQR) " in ", AX_EQ_B_QR) "()\n", -info);
205             exit(1);
206             }
207             else{
208             fprintf(stderr, RCAT("Unknown LAPACK error (%d) in ", AX_EQ_B_QR) "()\n", info);
209             #ifndef LINSOLVERS_RETAIN_MEMORY
210             free(buf);
211             #endif
212              
213             return 0;
214             }
215             }
216              
217             /* Q is now in a; compute Q^T b in x */
218             for(i=0; i
219             for(j=0, sum=0.0; j
220             sum+=a[i*m+j]*B[j];
221             x[i]=sum;
222             }
223              
224             /* solve the linear system R x = Q^t b */
225             TRTRS("U", "N", "N", (int *)&m, (int *)&nrhs, r, (int *)&m, x, (int *)&m, &info);
226             /* error treatment */
227             if(info!=0){
228             if(info<0){
229             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", TRTRS) " in ", AX_EQ_B_QR) "()\n", -info);
230             exit(1);
231             }
232             else{
233             fprintf(stderr, RCAT("LAPACK error: the %d-th diagonal element of A is zero (singular matrix) in ", AX_EQ_B_QR) "()\n", info);
234             #ifndef LINSOLVERS_RETAIN_MEMORY
235             free(buf);
236             #endif
237              
238             return 0;
239             }
240             }
241              
242             #ifndef LINSOLVERS_RETAIN_MEMORY
243             free(buf);
244             #endif
245              
246             return 1;
247             }
248              
249             /*
250             * This function returns the solution of min_x ||Ax - b||
251             *
252             * || . || is the second order (i.e. L2) norm. This is a least squares technique that
253             * is based on QR decomposition:
254             * If A=Q R with Q orthogonal and R upper triangular, the normal equations become
255             * (A^T A) x = A^T b or (R^T Q^T Q R) x = A^T b or (R^T R) x = A^T b.
256             * This amounts to solving R^T y = A^T b for y and then R x = y for x
257             * Note that Q does not need to be explicitly computed
258             *
259             * A is mxn, b is mx1
260             *
261             * The function returns 0 in case of error, 1 if successful
262             *
263             * This function is often called repetitively to solve problems of identical
264             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
265             * retained between calls and free'd-malloc'ed when not of the appropriate size.
266             * A call with NULL as the first argument forces this memory to be released.
267             */
268             int AX_EQ_B_QRLS(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m, int n)
269             {
270             __STATIC__ LM_REAL *buf=NULL;
271             __STATIC__ int buf_sz=0;
272              
273             static int nb=0; /* no __STATIC__ decl. here! */
274              
275             LM_REAL *a, *tau, *r, *work;
276             int a_sz, tau_sz, r_sz, tot_sz;
277             register int i, j;
278             int info, worksz, nrhs=1;
279             register LM_REAL sum;
280            
281             if(!A)
282             #ifdef LINSOLVERS_RETAIN_MEMORY
283             {
284             if(buf) free(buf);
285             buf=NULL;
286             buf_sz=0;
287              
288             return 1;
289             }
290             #else
291             return 1; /* NOP */
292             #endif /* LINSOLVERS_RETAIN_MEMORY */
293            
294             if(m
295             fprintf(stderr, RCAT("Normal equations require that the number of rows is greater than number of columns in ", AX_EQ_B_QRLS) "() [%d x %d]! -- try transposing\n", m, n);
296             exit(1);
297             }
298            
299             /* calculate required memory size */
300             a_sz=m*n;
301             tau_sz=n;
302             r_sz=n*n;
303             if(!nb){
304             LM_REAL tmp;
305              
306             worksz=-1; // workspace query; optimal size is returned in tmp
307             GEQRF((int *)&m, (int *)&m, NULL, (int *)&m, NULL, (LM_REAL *)&tmp, (int *)&worksz, (int *)&info);
308             nb=((int)tmp)/m; // optimal worksize is m*nb
309             }
310             worksz=nb*m;
311             tot_sz=a_sz + tau_sz + r_sz + worksz;
312              
313             #ifdef LINSOLVERS_RETAIN_MEMORY
314             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
315             if(buf) free(buf); /* free previously allocated memory */
316              
317             buf_sz=tot_sz;
318             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
319             if(!buf){
320             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_QRLS) "() failed!\n");
321             exit(1);
322             }
323             }
324             #else
325             buf_sz=tot_sz;
326             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
327             if(!buf){
328             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_QRLS) "() failed!\n");
329             exit(1);
330             }
331             #endif /* LINSOLVERS_RETAIN_MEMORY */
332              
333             a=buf;
334             tau=a+a_sz;
335             r=tau+tau_sz;
336             work=r+r_sz;
337              
338             /* store A (column major!) into a */
339             for(i=0; i
340             for(j=0; j
341             a[i+j*m]=A[i*n+j];
342              
343             /* compute A^T b in x */
344             for(i=0; i
345             for(j=0, sum=0.0; j
346             sum+=A[j*n+i]*B[j];
347             x[i]=sum;
348             }
349              
350             /* QR decomposition of A */
351             GEQRF((int *)&m, (int *)&n, a, (int *)&m, tau, work, (int *)&worksz, (int *)&info);
352             /* error treatment */
353             if(info!=0){
354             if(info<0){
355             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", GEQRF) " in ", AX_EQ_B_QRLS) "()\n", -info);
356             exit(1);
357             }
358             else{
359             fprintf(stderr, RCAT(RCAT("Unknown LAPACK error %d for ", GEQRF) " in ", AX_EQ_B_QRLS) "()\n", info);
360             #ifndef LINSOLVERS_RETAIN_MEMORY
361             free(buf);
362             #endif
363              
364             return 0;
365             }
366             }
367              
368             /* R is stored in the upper triangular part of a. Note that a is mxn while r nxn */
369             for(j=0; j
370             for(i=0; i<=j; i++)
371             r[i+j*n]=a[i+j*m];
372              
373             /* lower part is zero */
374             for(i=j+1; i
375             r[i+j*n]=0.0;
376             }
377              
378             /* solve the linear system R^T y = A^t b */
379             TRTRS("U", "T", "N", (int *)&n, (int *)&nrhs, r, (int *)&n, x, (int *)&n, &info);
380             /* error treatment */
381             if(info!=0){
382             if(info<0){
383             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", TRTRS) " in ", AX_EQ_B_QRLS) "()\n", -info);
384             exit(1);
385             }
386             else{
387             fprintf(stderr, RCAT("LAPACK error: the %d-th diagonal element of A is zero (singular matrix) in ", AX_EQ_B_QRLS) "()\n", info);
388             #ifndef LINSOLVERS_RETAIN_MEMORY
389             free(buf);
390             #endif
391              
392             return 0;
393             }
394             }
395              
396             /* solve the linear system R x = y */
397             TRTRS("U", "N", "N", (int *)&n, (int *)&nrhs, r, (int *)&n, x, (int *)&n, &info);
398             /* error treatment */
399             if(info!=0){
400             if(info<0){
401             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", TRTRS) " in ", AX_EQ_B_QRLS) "()\n", -info);
402             exit(1);
403             }
404             else{
405             fprintf(stderr, RCAT("LAPACK error: the %d-th diagonal element of A is zero (singular matrix) in ", AX_EQ_B_QRLS) "()\n", info);
406             #ifndef LINSOLVERS_RETAIN_MEMORY
407             free(buf);
408             #endif
409              
410             return 0;
411             }
412             }
413              
414             #ifndef LINSOLVERS_RETAIN_MEMORY
415             free(buf);
416             #endif
417              
418             return 1;
419             }
420              
421             /*
422             * This function returns the solution of Ax=b
423             *
424             * The function assumes that A is symmetric & postive definite and employs
425             * the Cholesky decomposition:
426             * If A=U^T U with U upper triangular, the system to be solved becomes
427             * (U^T U) x = b
428             * This amount to solving U^T y = b for y and then U x = y for x
429             *
430             * A is mxm, b is mx1
431             *
432             * The function returns 0 in case of error, 1 if successful
433             *
434             * This function is often called repetitively to solve problems of identical
435             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
436             * retained between calls and free'd-malloc'ed when not of the appropriate size.
437             * A call with NULL as the first argument forces this memory to be released.
438             */
439             int AX_EQ_B_CHOL(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
440             {
441             __STATIC__ LM_REAL *buf=NULL;
442             __STATIC__ int buf_sz=0;
443              
444             LM_REAL *a;
445             int a_sz, tot_sz;
446             register int i;
447             int info, nrhs=1;
448            
449             if(!A)
450             #ifdef LINSOLVERS_RETAIN_MEMORY
451             {
452             if(buf) free(buf);
453             buf=NULL;
454             buf_sz=0;
455              
456             return 1;
457             }
458             #else
459             return 1; /* NOP */
460             #endif /* LINSOLVERS_RETAIN_MEMORY */
461            
462             /* calculate required memory size */
463             a_sz=m*m;
464             tot_sz=a_sz;
465              
466             #ifdef LINSOLVERS_RETAIN_MEMORY
467             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
468             if(buf) free(buf); /* free previously allocated memory */
469              
470             buf_sz=tot_sz;
471             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
472             if(!buf){
473             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_CHOL) "() failed!\n");
474             exit(1);
475             }
476             }
477             #else
478             buf_sz=tot_sz;
479             buf=(LM_REAL *)malloc(buf_sz*sizeof(LM_REAL));
480             if(!buf){
481             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_CHOL) "() failed!\n");
482             exit(1);
483             }
484             #endif /* LINSOLVERS_RETAIN_MEMORY */
485              
486             a=buf;
487              
488             /* store A into a and B into x. A is assumed symmetric,
489             * hence no transposition is needed
490             */
491             for(i=0; i
492             a[i]=A[i];
493             x[i]=B[i];
494             }
495             for(i=m; i
496             a[i]=A[i];
497              
498             /* Cholesky decomposition of A */
499             //POTF2("U", (int *)&m, a, (int *)&m, (int *)&info);
500             POTRF("U", (int *)&m, a, (int *)&m, (int *)&info);
501             /* error treatment */
502             if(info!=0){
503             if(info<0){
504             fprintf(stderr, RCAT(RCAT(RCAT("LAPACK error: illegal value for argument %d of ", POTF2) "/", POTRF) " in ",
505             AX_EQ_B_CHOL) "()\n", -info);
506             exit(1);
507             }
508             else{
509             fprintf(stderr, RCAT(RCAT(RCAT("LAPACK error: the leading minor of order %d is not positive definite,\nthe factorization could not be completed for ", POTF2) "/", POTRF) " in ", AX_EQ_B_CHOL) "()\n", info);
510             #ifndef LINSOLVERS_RETAIN_MEMORY
511             free(buf);
512             #endif
513              
514             return 0;
515             }
516             }
517              
518             /* solve using the computed Cholesky in one lapack call */
519             POTRS("U", (int *)&m, (int *)&nrhs, a, (int *)&m, x, (int *)&m, &info);
520             if(info<0){
521             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", POTRS) " in ", AX_EQ_B_CHOL) "()\n", -info);
522             exit(1);
523             }
524              
525             #if 0
526             /* alternative: solve the linear system U^T y = b ... */
527             TRTRS("U", "T", "N", (int *)&m, (int *)&nrhs, a, (int *)&m, x, (int *)&m, &info);
528             /* error treatment */
529             if(info!=0){
530             if(info<0){
531             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", TRTRS) " in ", AX_EQ_B_CHOL) "()\n", -info);
532             exit(1);
533             }
534             else{
535             fprintf(stderr, RCAT("LAPACK error: the %d-th diagonal element of A is zero (singular matrix) in ", AX_EQ_B_CHOL) "()\n", info);
536             #ifndef LINSOLVERS_RETAIN_MEMORY
537             free(buf);
538             #endif
539              
540             return 0;
541             }
542             }
543              
544             /* ... solve the linear system U x = y */
545             TRTRS("U", "N", "N", (int *)&m, (int *)&nrhs, a, (int *)&m, x, (int *)&m, &info);
546             /* error treatment */
547             if(info!=0){
548             if(info<0){
549             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", TRTRS) "in ", AX_EQ_B_CHOL) "()\n", -info);
550             exit(1);
551             }
552             else{
553             fprintf(stderr, RCAT("LAPACK error: the %d-th diagonal element of A is zero (singular matrix) in ", AX_EQ_B_CHOL) "()\n", info);
554             #ifndef LINSOLVERS_RETAIN_MEMORY
555             free(buf);
556             #endif
557              
558             return 0;
559             }
560             }
561             #endif /* 0 */
562              
563             #ifndef LINSOLVERS_RETAIN_MEMORY
564             free(buf);
565             #endif
566              
567             return 1;
568             }
569              
570             /*
571             * This function returns the solution of Ax = b
572             *
573             * The function employs LU decomposition:
574             * If A=L U with L lower and U upper triangular, then the original system
575             * amounts to solving
576             * L y = b, U x = y
577             *
578             * A is mxm, b is mx1
579             *
580             * The function returns 0 in case of error, 1 if successful
581             *
582             * This function is often called repetitively to solve problems of identical
583             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
584             * retained between calls and free'd-malloc'ed when not of the appropriate size.
585             * A call with NULL as the first argument forces this memory to be released.
586             */
587             int AX_EQ_B_LU(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
588             {
589             __STATIC__ LM_REAL *buf=NULL;
590             __STATIC__ int buf_sz=0;
591              
592             int a_sz, ipiv_sz, tot_sz;
593             register int i, j;
594             int info, *ipiv, nrhs=1;
595             LM_REAL *a;
596            
597             if(!A)
598             #ifdef LINSOLVERS_RETAIN_MEMORY
599             {
600             if(buf) free(buf);
601             buf=NULL;
602             buf_sz=0;
603              
604             return 1;
605             }
606             #else
607             return 1; /* NOP */
608             #endif /* LINSOLVERS_RETAIN_MEMORY */
609            
610             /* calculate required memory size */
611             ipiv_sz=m;
612             a_sz=m*m;
613             tot_sz=a_sz*sizeof(LM_REAL) + ipiv_sz*sizeof(int); /* should be arranged in that order for proper doubles alignment */
614              
615             #ifdef LINSOLVERS_RETAIN_MEMORY
616             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
617             if(buf) free(buf); /* free previously allocated memory */
618              
619             buf_sz=tot_sz;
620             buf=(LM_REAL *)malloc(buf_sz);
621             if(!buf){
622             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_LU) "() failed!\n");
623             exit(1);
624             }
625             }
626             #else
627             buf_sz=tot_sz;
628             buf=(LM_REAL *)malloc(buf_sz);
629             if(!buf){
630             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_LU) "() failed!\n");
631             exit(1);
632             }
633             #endif /* LINSOLVERS_RETAIN_MEMORY */
634              
635             a=buf;
636             ipiv=(int *)(a+a_sz);
637              
638             /* store A (column major!) into a and B into x */
639             for(i=0; i
640             for(j=0; j
641             a[i+j*m]=A[i*m+j];
642              
643             x[i]=B[i];
644             }
645              
646             /* LU decomposition for A */
647             GETRF((int *)&m, (int *)&m, a, (int *)&m, ipiv, (int *)&info);
648             if(info!=0){
649             if(info<0){
650             fprintf(stderr, RCAT(RCAT("argument %d of ", GETRF) " illegal in ", AX_EQ_B_LU) "()\n", -info);
651             exit(1);
652             }
653             else{
654             fprintf(stderr, RCAT(RCAT("singular matrix A for ", GETRF) " in ", AX_EQ_B_LU) "()\n");
655             #ifndef LINSOLVERS_RETAIN_MEMORY
656             free(buf);
657             #endif
658              
659             return 0;
660             }
661             }
662              
663             /* solve the system with the computed LU */
664             GETRS("N", (int *)&m, (int *)&nrhs, a, (int *)&m, ipiv, x, (int *)&m, (int *)&info);
665             if(info!=0){
666             if(info<0){
667             fprintf(stderr, RCAT(RCAT("argument %d of ", GETRS) " illegal in ", AX_EQ_B_LU) "()\n", -info);
668             exit(1);
669             }
670             else{
671             fprintf(stderr, RCAT(RCAT("unknown error for ", GETRS) " in ", AX_EQ_B_LU) "()\n");
672             #ifndef LINSOLVERS_RETAIN_MEMORY
673             free(buf);
674             #endif
675              
676             return 0;
677             }
678             }
679              
680             #ifndef LINSOLVERS_RETAIN_MEMORY
681             free(buf);
682             #endif
683              
684             return 1;
685             }
686              
687             /*
688             * This function returns the solution of Ax = b
689             *
690             * The function is based on SVD decomposition:
691             * If A=U D V^T with U, V orthogonal and D diagonal, the linear system becomes
692             * (U D V^T) x = b or x=V D^{-1} U^T b
693             * Note that V D^{-1} U^T is the pseudoinverse A^+
694             *
695             * A is mxm, b is mx1.
696             *
697             * The function returns 0 in case of error, 1 if successful
698             *
699             * This function is often called repetitively to solve problems of identical
700             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
701             * retained between calls and free'd-malloc'ed when not of the appropriate size.
702             * A call with NULL as the first argument forces this memory to be released.
703             */
704             int AX_EQ_B_SVD(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
705             {
706             __STATIC__ LM_REAL *buf=NULL;
707             __STATIC__ int buf_sz=0;
708             static LM_REAL eps=LM_CNST(-1.0);
709              
710             register int i, j;
711             LM_REAL *a, *u, *s, *vt, *work;
712             int a_sz, u_sz, s_sz, vt_sz, tot_sz;
713             LM_REAL thresh, one_over_denom;
714             register LM_REAL sum;
715             int info, rank, worksz, *iwork, iworksz;
716            
717             if(!A)
718             #ifdef LINSOLVERS_RETAIN_MEMORY
719             {
720             if(buf) free(buf);
721             buf=NULL;
722             buf_sz=0;
723              
724             return 1;
725             }
726             #else
727             return 1; /* NOP */
728             #endif /* LINSOLVERS_RETAIN_MEMORY */
729            
730             /* calculate required memory size */
731             #if 1 /* use optimal size */
732             worksz=-1; // workspace query. Keep in mind that GESDD requires more memory than GESVD
733             /* note that optimal work size is returned in thresh */
734             GESVD("A", "A", (int *)&m, (int *)&m, NULL, (int *)&m, NULL, NULL, (int *)&m, NULL, (int *)&m, (LM_REAL *)&thresh, (int *)&worksz, &info);
735             //GESDD("A", (int *)&m, (int *)&m, NULL, (int *)&m, NULL, NULL, (int *)&m, NULL, (int *)&m, (LM_REAL *)&thresh, (int *)&worksz, NULL, &info);
736             worksz=(int)thresh;
737             #else /* use minimum size */
738             worksz=5*m; // min worksize for GESVD
739             //worksz=m*(7*m+4); // min worksize for GESDD
740             #endif
741             iworksz=8*m;
742             a_sz=m*m;
743             u_sz=m*m; s_sz=m; vt_sz=m*m;
744              
745             tot_sz=(a_sz + u_sz + s_sz + vt_sz + worksz)*sizeof(LM_REAL) + iworksz*sizeof(int); /* should be arranged in that order for proper doubles alignment */
746              
747             #ifdef LINSOLVERS_RETAIN_MEMORY
748             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
749             if(buf) free(buf); /* free previously allocated memory */
750              
751             buf_sz=tot_sz;
752             buf=(LM_REAL *)malloc(buf_sz);
753             if(!buf){
754             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_SVD) "() failed!\n");
755             exit(1);
756             }
757             }
758             #else
759             buf_sz=tot_sz;
760             buf=(LM_REAL *)malloc(buf_sz);
761             if(!buf){
762             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_SVD) "() failed!\n");
763             exit(1);
764             }
765             #endif /* LINSOLVERS_RETAIN_MEMORY */
766              
767             a=buf;
768             u=a+a_sz;
769             s=u+u_sz;
770             vt=s+s_sz;
771             work=vt+vt_sz;
772             iwork=(int *)(work+worksz);
773              
774             /* store A (column major!) into a */
775             for(i=0; i
776             for(j=0; j
777             a[i+j*m]=A[i*m+j];
778              
779             /* SVD decomposition of A */
780             GESVD("A", "A", (int *)&m, (int *)&m, a, (int *)&m, s, u, (int *)&m, vt, (int *)&m, work, (int *)&worksz, &info);
781             //GESDD("A", (int *)&m, (int *)&m, a, (int *)&m, s, u, (int *)&m, vt, (int *)&m, work, (int *)&worksz, iwork, &info);
782              
783             /* error treatment */
784             if(info!=0){
785             if(info<0){
786             fprintf(stderr, RCAT(RCAT(RCAT("LAPACK error: illegal value for argument %d of ", GESVD), "/" GESDD) " in ", AX_EQ_B_SVD) "()\n", -info);
787             exit(1);
788             }
789             else{
790             fprintf(stderr, RCAT("LAPACK error: dgesdd (dbdsdc)/dgesvd (dbdsqr) failed to converge in ", AX_EQ_B_SVD) "() [info=%d]\n", info);
791             #ifndef LINSOLVERS_RETAIN_MEMORY
792             free(buf);
793             #endif
794              
795             return 0;
796             }
797             }
798              
799             if(eps<0.0){
800             LM_REAL aux;
801              
802             /* compute machine epsilon */
803             for(eps=LM_CNST(1.0); aux=eps+LM_CNST(1.0), aux-LM_CNST(1.0)>0.0; eps*=LM_CNST(0.5))
804             ;
805             eps*=LM_CNST(2.0);
806             }
807              
808             /* compute the pseudoinverse in a */
809             for(i=0; i
810             for(rank=0, thresh=eps*s[0]; rankthresh; rank++){
811             one_over_denom=LM_CNST(1.0)/s[rank];
812              
813             for(j=0; j
814             for(i=0; i
815             a[i*m+j]+=vt[rank+i*m]*u[j+rank*m]*one_over_denom;
816             }
817              
818             /* compute A^+ b in x */
819             for(i=0; i
820             for(j=0, sum=0.0; j
821             sum+=a[i*m+j]*B[j];
822             x[i]=sum;
823             }
824              
825             #ifndef LINSOLVERS_RETAIN_MEMORY
826             free(buf);
827             #endif
828              
829             return 1;
830             }
831              
832             /*
833             * This function returns the solution of Ax = b for a real symmetric matrix A
834             *
835             * The function is based on UDUT factorization with the pivoting
836             * strategy of Bunch and Kaufman:
837             * A is factored as U*D*U^T where U is upper triangular and
838             * D symmetric and block diagonal (aka spectral decomposition,
839             * Banachiewicz factorization, modified Cholesky factorization)
840             *
841             * A is mxm, b is mx1.
842             *
843             * The function returns 0 in case of error, 1 if successfull
844             *
845             * This function is often called repetitively to solve problems of identical
846             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
847             * retained between calls and free'd-malloc'ed when not of the appropriate size.
848             * A call with NULL as the first argument forces this memory to be released.
849             */
850             int AX_EQ_B_BK(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
851             {
852             __STATIC__ LM_REAL *buf=NULL;
853             __STATIC__ int buf_sz=0, nb=0;
854              
855             LM_REAL *a, *work;
856             int a_sz, ipiv_sz, work_sz, tot_sz;
857             register int i, j;
858             int info, *ipiv, nrhs=1;
859            
860             if(!A)
861             #ifdef LINSOLVERS_RETAIN_MEMORY
862             {
863             if(buf) free(buf);
864             buf=NULL;
865             buf_sz=0;
866              
867             return 1;
868             }
869             #else
870             return 1; /* NOP */
871             #endif /* LINSOLVERS_RETAIN_MEMORY */
872              
873             /* calculate required memory size */
874             ipiv_sz=m;
875             a_sz=m*m;
876             if(!nb){
877             LM_REAL tmp;
878              
879             work_sz=-1; // workspace query; optimal size is returned in tmp
880             SYTRF("U", (int *)&m, NULL, (int *)&m, NULL, (LM_REAL *)&tmp, (int *)&work_sz, (int *)&info);
881             nb=((int)tmp)/m; // optimal worksize is m*nb
882             }
883             work_sz=(nb!=-1)? nb*m : 1;
884             tot_sz=(a_sz + work_sz)*sizeof(LM_REAL) + ipiv_sz*sizeof(int); /* should be arranged in that order for proper doubles alignment */
885              
886             #ifdef LINSOLVERS_RETAIN_MEMORY
887             if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
888             if(buf) free(buf); /* free previously allocated memory */
889              
890             buf_sz=tot_sz;
891             buf=(LM_REAL *)malloc(buf_sz);
892             if(!buf){
893             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_BK) "() failed!\n");
894             exit(1);
895             }
896             }
897             #else
898             buf_sz=tot_sz;
899             buf=(LM_REAL *)malloc(buf_sz);
900             if(!buf){
901             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_BK) "() failed!\n");
902             exit(1);
903             }
904             #endif /* LINSOLVERS_RETAIN_MEMORY */
905              
906             a=buf;
907             work=a+a_sz;
908             ipiv=(int *)(work+work_sz);
909              
910             /* store A into a and B into x; A is assumed to be symmetric, hence
911             * the column and row major order representations are the same
912             */
913             for(i=0; i
914             a[i]=A[i];
915             x[i]=B[i];
916             }
917             for(j=m*m; i
918             a[i]=A[i];
919              
920             /* UDUt factorization for A */
921             SYTRF("U", (int *)&m, a, (int *)&m, ipiv, work, (int *)&work_sz, (int *)&info);
922             if(info!=0){
923             if(info<0){
924             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", SYTRF) " in ", AX_EQ_B_BK) "()\n", -info);
925             exit(1);
926             }
927             else{
928             fprintf(stderr, RCAT(RCAT("LAPACK error: singular block diagonal matrix D for", SYTRF) " in ", AX_EQ_B_BK)"() [D(%d, %d) is zero]\n", info, info);
929             #ifndef LINSOLVERS_RETAIN_MEMORY
930             free(buf);
931             #endif
932              
933             return 0;
934             }
935             }
936              
937             /* solve the system with the computed factorization */
938             SYTRS("U", (int *)&m, (int *)&nrhs, a, (int *)&m, ipiv, x, (int *)&m, (int *)&info);
939             if(info<0){
940             fprintf(stderr, RCAT(RCAT("LAPACK error: illegal value for argument %d of ", SYTRS) " in ", AX_EQ_B_BK) "()\n", -info);
941             exit(1);
942             }
943              
944             #ifndef LINSOLVERS_RETAIN_MEMORY
945             free(buf);
946             #endif
947              
948             return 1;
949             }
950              
951             /* undefine all. IT MUST REMAIN IN THIS POSITION IN FILE */
952             #undef AX_EQ_B_QR
953             #undef AX_EQ_B_QRLS
954             #undef AX_EQ_B_CHOL
955             #undef AX_EQ_B_LU
956             #undef AX_EQ_B_SVD
957             #undef AX_EQ_B_BK
958              
959             #undef GEQRF
960             #undef ORGQR
961             #undef TRTRS
962             #undef POTF2
963             #undef POTRF
964             #undef POTRS
965             #undef GETRF
966             #undef GETRS
967             #undef GESVD
968             #undef GESDD
969             #undef SYTRF
970             #undef SYTRS
971              
972             #else // no LAPACK
973              
974             /* precision-specific definitions */
975             #define AX_EQ_B_LU LM_ADD_PREFIX(Ax_eq_b_LU_noLapack)
976              
977             /*
978             * This function returns the solution of Ax = b
979             *
980             * The function employs LU decomposition followed by forward/back substitution (see
981             * also the LAPACK-based LU solver above)
982             *
983             * A is mxm, b is mx1
984             *
985             * The function returns 0 in case of error, 1 if successful
986             *
987             * This function is often called repetitively to solve problems of identical
988             * dimensions. To avoid repetitive malloc's and free's, allocated memory is
989             * retained between calls and free'd-malloc'ed when not of the appropriate size.
990             * A call with NULL as the first argument forces this memory to be released.
991             */
992 102370           int AX_EQ_B_LU(LM_REAL *A, LM_REAL *B, LM_REAL *x, int m)
993             {
994             __STATIC__ void *buf=NULL;
995             __STATIC__ int buf_sz=0;
996              
997             register int i, j, k;
998 102370           int *idx, maxi=-1, idx_sz, a_sz, work_sz, tot_sz;
999             LM_REAL *a, *work, max, sum, tmp;
1000              
1001 102370 100         if(!A)
    100          
1002             #ifdef LINSOLVERS_RETAIN_MEMORY
1003             {
1004 159 50         if(buf) free(buf);
    50          
1005 159           buf=NULL;
1006 159           buf_sz=0;
1007              
1008 159           return 1;
1009             }
1010             #else
1011             return 1; /* NOP */
1012             #endif /* LINSOLVERS_RETAIN_MEMORY */
1013            
1014             /* calculate required memory size */
1015 102211           idx_sz=m;
1016 102211           a_sz=m*m;
1017 102211           work_sz=m;
1018 102211           tot_sz=(a_sz+work_sz)*sizeof(LM_REAL) + idx_sz*sizeof(int); /* should be arranged in that order for proper doubles alignment */
1019              
1020             #ifdef LINSOLVERS_RETAIN_MEMORY
1021 102211 100         if(tot_sz>buf_sz){ /* insufficient memory, allocate a "big" memory chunk at once */
    100          
1022 159 50         if(buf) free(buf); /* free previously allocated memory */
    50          
1023              
1024 159           buf_sz=tot_sz;
1025 159           buf=(void *)malloc(tot_sz);
1026 159 50         if(!buf){
    50          
1027 0           fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_LU) "() failed!\n");
1028 0           exit(1);
1029             }
1030             }
1031             #else
1032             buf_sz=tot_sz;
1033             buf=(void *)malloc(tot_sz);
1034             if(!buf){
1035             fprintf(stderr, RCAT("memory allocation in ", AX_EQ_B_LU) "() failed!\n");
1036             exit(1);
1037             }
1038             #endif /* LINSOLVERS_RETAIN_MEMORY */
1039              
1040 102211           a=buf;
1041 102211           work=a+a_sz;
1042 102211           idx=(int *)(work+work_sz);
1043              
1044             /* avoid destroying A, B by copying them to a, x resp. */
1045 427695 100         for(i=0; i
    100          
1046 325484           a[i]=A[i];
1047 325484           x[i]=B[i];
1048             }
1049 911825 100         for( ; i
    100          
1050             /****
1051             for(i=0; i
1052             for(j=0; j
1053             a[i*m+j]=A[i*m+j];
1054             x[i]=B[i];
1055             }
1056             ****/
1057              
1058             /* compute the LU decomposition of a row permutation of matrix a; the permutation itself is saved in idx[] */
1059 427695 100         for(i=0; i
    100          
1060 325484           max=0.0;
1061 1460582 100         for(j=0; j
    100          
1062 1135098 100         if((tmp=FABS(a[i*m+j]))>max)
    100          
    100          
    100          
1063 487239           max=tmp;
1064 325484 50         if(max==0.0){
    50          
1065 0           fprintf(stderr, RCAT("Singular matrix A in ", AX_EQ_B_LU) "()!\n");
1066             #ifndef LINSOLVERS_RETAIN_MEMORY
1067             free(buf);
1068             #endif
1069              
1070 0           return 0;
1071             }
1072 325484           work[i]=LM_CNST(1.0)/max;
1073             }
1074              
1075 427695 100         for(j=0; j
    100          
1076 730291 100         for(i=0; i
    100          
1077 404807           sum=a[i*m+j];
1078 646813 100         for(k=0; k
    100          
1079 242006           sum-=a[i*m+k]*a[k*m+j];
1080 404807           a[i*m+j]=sum;
1081             }
1082 325484           max=0.0;
1083 1055775 100         for(i=j; i
    100          
1084 730291           sum=a[i*m+j];
1085 1377104 100         for(k=0; k
    100          
1086 646813           sum-=a[i*m+k]*a[k*m+j];
1087 730291           a[i*m+j]=sum;
1088 730291 100         if((tmp=work[i]*FABS(sum))>=max){
    100          
    100          
    100          
1089 387470           max=tmp;
1090 387470           maxi=i;
1091             }
1092             }
1093 325484 100         if(j!=maxi){
    100          
1094 306905 100         for(k=0; k
    100          
1095 244930           tmp=a[maxi*m+k];
1096 244930           a[maxi*m+k]=a[j*m+k];
1097 244930           a[j*m+k]=tmp;
1098             }
1099 61975           work[maxi]=work[j];
1100             }
1101 325484           idx[j]=maxi;
1102 325484 50         if(a[j*m+j]==0.0)
    50          
1103 0           a[j*m+j]=LM_REAL_EPSILON;
1104 325484 100         if(j!=m-1){
    100          
1105 223273           tmp=LM_CNST(1.0)/(a[j*m+j]);
1106 628080 100         for(i=j+1; i
    100          
1107 404807           a[i*m+j]*=tmp;
1108             }
1109             }
1110              
1111             /* The decomposition has now replaced a. Solve the linear system using
1112             * forward and back substitution
1113             */
1114 427695 100         for(i=k=0; i
    100          
1115 325484           j=idx[i];
1116 325484           sum=x[j];
1117 325484           x[j]=x[i];
1118 325484 100         if(k!=0)
    100          
1119 628080 100         for(j=k-1; j
    100          
1120 404807           sum-=a[i*m+j]*x[j];
1121             else
1122 102211 50         if(sum!=0.0)
    50          
1123 102211           k=i+1;
1124 325484           x[i]=sum;
1125             }
1126              
1127 427695 100         for(i=m-1; i>=0; --i){
    100          
1128 325484           sum=x[i];
1129 730291 100         for(j=i+1; j
    100          
1130 404807           sum-=a[i*m+j]*x[j];
1131 325484           x[i]=sum/a[i*m+i];
1132             }
1133              
1134             #ifndef LINSOLVERS_RETAIN_MEMORY
1135             free(buf);
1136             #endif
1137              
1138 102211           return 1;
1139             }
1140              
1141             /* undefine all. IT MUST REMAIN IN THIS POSITION IN FILE */
1142             #undef AX_EQ_B_LU
1143              
1144             #endif /* HAVE_LAPACK */