summaryrefslogtreecommitdiff
path: root/2.3-1/src/c/linearAlgebra/svd/dsvda.c
blob: c3bcfc29d883d2842456aee1caba2b30b6c295f2 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
/* Copyright (C) 2017 - IIT Bombay - FOSSEE

 This file must be used under the terms of the CeCILL.
 This source file is licensed as described in the file COPYING, which
 you should have received as part of this distribution.  The terms
 are also available at
 http://www.cecill.info/licences/Licence_CeCILL_V2-en.txt
 Author: Sandeep Gupta
 Organization: FOSSEE, IIT Bombay
 Email: toolbox@scilab.in

 */

/*Funtion to find singular value decomposition of given matrix */

#include "lapack.h"
#include <stdio.h>
#include <stdlib.h>
#include "string.h"
#include <math.h>
#include "svd.h"
#include "matrixTranspose.h"

int min(int a,int b);
int max(int a,int b);

extern double dgesvd_(char*,char*,int*,int*,double*,int*,double*,double*,int*,\
	double*,int*,double *,int*,int*);

#define eps 2.22044604925e-16 /* pow(2,-52) */

/*  DGESVD computes the singular value decomposition (SVD) of a real
 M-by-N matrix A, optionally computing the left and/or right singular
 vectors. The SVD is written

      A = U * SIGMA * transpose(V) */

/*Function support -

s=svd(X)
[U,S,V]=svd(X)
[U,S,V]=svd(X,0) (obsolete)
[U,S,V]=svd(X,"e")
[U,S,V,rk]=svd(X [,tol])

*/

double dsvda(double tol,double *in1,int row,int col,double in2,double nout,double *out1, \
	double *out2,double *out3){
	
	char JOBU,JOBVT;
	int i,j,k;
	int LDU=1;			/*Leading Dimension of U */
	int LDVT=1;			/*Leading Dimension of VT */
	int M = row;
	int N = col;
	double *buf;
	double *S,*U,*VT;
	double *WORK;	
	
	int rk; /*Fourth output if needed */

	/*if((nout > 1 && in2 == 1) && (M != N)){   // [U,S,VT] = svd(x,'e')
		if(M > N){		
			JOBU = 'S';
			JOBVT = 'A';
			LDVT = N;
		}
		else{
			JOBU = 'A';
			JOBVT = 'S';		
			LDVT = min(M,N);	
		}
		LDU = M; 		
		U = (double*) malloc((double) (LDU)*min(M,N)*sizeof(double));
		VT = (double*) malloc((double) (LDVT)*N*sizeof(double));
	}
	else */if(nout > 1){		/* [U,S,VT = svd(x)] */
		JOBU = 'A'; 		/*If JOBU = 'A', U contains the M-by-M orthogonal matrix U */
		JOBVT = 'A';		/*JOBVT = 'A':  all N rows of V**T are returned in the array VT;*/
		LDU = M;
		LDVT = N;
		U = (double*) malloc((double) M*M*sizeof(double));
		VT = (double*) malloc((double) N*N*sizeof(double));
	}	
	else{				/* ans = svd(x) */
		JOBU = 'N';	
		JOBVT = 'N';
	}
	int LDA = max(1,M);
	
	/* Making a copy of input matrix */
	buf = (double*) malloc((double)M*N*sizeof(double));
	memcpy(buf,in1,M*N*sizeof(double));
	
	S = (double*)malloc((double)min(col,row)*sizeof(double));

	int LWORK = 5*min(M,N);
	WORK = (double*)malloc((double)LWORK*sizeof(double));
	int INFO = 0;	/*For successful exit */

	dgesvd_(&JOBU,&JOBVT,&M,&N,buf,&LDA,S,U,&LDU,VT,&LDVT,WORK,&LWORK,&INFO);
	/*Subroutine DGESVD from Lapack lib. */

	if (nout  == 1){	/* ans = svd(x)*/
		memcpy(out1,S,min(row,col)*sizeof(double));
		//printf("%lf %lf %lf",*(S),*(S+1),*(S+2));
	}			/* [U,S,VT] = svd(x) */
	else if(in2 == 0 && nout > 1){
		memcpy(out1,U,LDU*M*sizeof(double));
		//memcpy(out3,VT,LDVT*min(row,col)*sizeof(double));
		for(j=0;j<M;j++){
			for(k=0;k<N;k++){
				if(j == k) *((out2+j*(min(M,N)))+k) = *(S+j);
				else *((out2+j*(min(M,N)))+k) = 0;					
			}
		}
							
		//dtransposea(VT,LDVT,N,out3);
		/*As there is some patch of error in SVD, these lines are added */

		for(j=1;j<=N;j++){
			for(i=j;i<=N;i++){
				*(out3+i+(j-1)*N-1) = VT[j+(i-1)*N-1];
				*(out3+j+(i-1)*N-1) = VT[i+(j-1)*N-1];
			}		
		}
		/*for(i=0;i<N;i++){
			for(j=0;j<N;j++){
				printf("%lf ",VT[i*row+j]);		
			}	
			printf("\n");
		}*/		
	}
	else{
		memcpy(out1,U,M*min(M,N)*sizeof(double));
		for(j=0;j<min(M,N);j++){
			for(k=0;k<min(M,N);k++){
				if(j == k) *((out2+j*(min(M,N)))+k) = *(S+j);
				else *((out2+j*(min(M,N)))+k) = 0;					
			}
		}
		//dtransposea(VT,LDVT,N,out3);
		/*As there is some patch of error in DGESVD, these lines are added */
		/* out3 first taken in some array then will be copied from it. */
		double *outV;
		outV = (double *)malloc(N*N*sizeof(double));
		for(j=1;j<=N;j++){
			for(i=j;i<=N;i++){
				*(outV+i+(j-1)*N-1) = VT[j+(i-1)*N-1];
				*(outV+j+(i-1)*N-1) = VT[i+(j-1)*N-1];
			}		
		}
		
		for(j=0;j<min(M,N)*N;j++){
			*(out3+j) = *(outV+j);		
		}		
	}	
	
	/* From the fortran file of scilab code - if(tol.eq.0.0d0) tol=dble(max(M,N))*eps*stk(lSV) */
	if(tol == 0){
		tol = (double)max(M,N)*eps*S[0];
	}
	if(nout == 4){	/*[U,S,VT,rk] = svd(X,tol) where tol - tolerance*/
		rk = 0;
		for(i=0;i<min(M,N);i++){
			if(S[i] > tol){
				rk = i+1;	
			}	
		}
		return rk;	
	}
	return 0;   
}

int min(int a,int b){
	if(a > b) return b;
	return a;
}

int max(int a,int b){
	if(a > b) return a;
	return b;
}