Actual source code: svdlapack.c
1: /*
2: This file implements a wrapper to the LAPACK SVD subroutines.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/svd/svdimpl.h
14: #include slepcblaslapack.h
18: PetscErrorCode SVDSetup_LAPACK(SVD svd)
19: {
20: PetscErrorCode ierr;
21: PetscInt N;
22: int i;
25: SVDMatGetSize(svd,PETSC_NULL,&N);
26: svd->ncv = N;
27: svd->max_it = 1;
28: if (svd->ncv!=svd->n) {
29: if (svd->U) {
30: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
31: PetscFree(svd->U);
32: }
33: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
34: for (i=0;i<svd->ncv;i++) { SVDMatGetVecs(svd,PETSC_NULL,svd->U+i); }
35: }
36: return(0);
37: }
41: PetscErrorCode SVDSolve_LAPACK(SVD svd)
42: {
43: PetscErrorCode ierr;
44: PetscInt M,N,n;
45: Mat mat;
46: PetscScalar *pU,*pVT,*pmat,*pu,*pv;
47: PetscReal *sigma;
48: int i,j,k;
49:
51: MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
52: MatGetArray(mat,&pmat);
53: MatGetSize(mat,&M,&N);
54: if (M>=N) {
55: n = N;
56: pU = PETSC_NULL;
57: PetscMalloc(sizeof(PetscScalar)*N*N,&pVT);
58: } else {
59: n = M;
60: PetscMalloc(sizeof(PetscScalar)*M*M,&pU);
61: pVT = PETSC_NULL;
62: }
63: PetscMalloc(sizeof(PetscReal)*n,&sigma);
64:
65: SVDDense(M,N,pmat,sigma,pU,pVT);
67: /* copy singular vectors */
68: for (i=0;i<n;i++) {
69: if (svd->which == SVD_SMALLEST) k = n - i - 1;
70: else k = i;
71: svd->sigma[k] = sigma[i];
72: VecGetArray(svd->U[k],&pu);
73: VecGetArray(svd->V[k],&pv);
74: if (M>=N) {
75: for (j=0;j<M;j++) pu[j] = pmat[i*M+j];
76: for (j=0;j<N;j++) pv[j] = pVT[j*N+i];
77: } else {
78: for (j=0;j<N;j++) pu[j] = pmat[j*M+i];
79: for (j=0;j<M;j++) pv[j] = pU[i*M+j];
80: }
81: VecRestoreArray(svd->U[k],&pu);
82: VecRestoreArray(svd->V[k],&pv);
83: }
85: svd->nconv = n;
86: svd->reason = SVD_CONVERGED_TOL;
88: MatRestoreArray(mat,&pmat);
89: MatDestroy(mat);
90: PetscFree(sigma);
91: if (M>=N) {
92: PetscFree(pVT);
93: } else {
94: PetscFree(pU);
95: }
96: return(0);
97: }
102: PetscErrorCode SVDCreate_LAPACK(SVD svd)
103: {
105: svd->ops->setup = SVDSetup_LAPACK;
106: svd->ops->solve = SVDSolve_LAPACK;
107: svd->ops->destroy = SVDDestroy_Default;
108: if (svd->transmode == PETSC_DECIDE)
109: svd->transmode = SVD_TRANSPOSE_IMPLICIT; /* don't build the transpose */
110: return(0);
111: }