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-2009, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include private/svdimpl.h
25: #include slepcblaslapack.h
29: PetscErrorCode SVDSetUp_LAPACK(SVD svd)
30: {
31: PetscErrorCode ierr;
32: PetscInt N,i,nloc;
33: PetscScalar *pU;
36: SVDMatGetSize(svd,PETSC_NULL,&N);
37: svd->ncv = N;
38: if (svd->mpd) PetscInfo(svd,"Warning: parameter mpd ignored\n");
39: svd->max_it = 1;
40: if (svd->ncv!=svd->n) {
41: if (svd->U) {
42: VecGetArray(svd->U[0],&pU);
43: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
44: PetscFree(pU);
45: PetscFree(svd->U);
46: }
47: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
48: SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
49: PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
50: for (i=0;i<svd->ncv;i++) {
51: VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
52: }
53: }
54: return(0);
55: }
59: PetscErrorCode SVDSolve_LAPACK(SVD svd)
60: {
61: PetscErrorCode ierr;
62: PetscInt M,N,n,i,j,k;
63: Mat mat;
64: PetscScalar *pU,*pVT,*pmat,*pu,*pv;
65: PetscReal *sigma;
66:
68: MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
69: MatGetArray(mat,&pmat);
70: MatGetSize(mat,&M,&N);
71: if (M>=N) {
72: n = N;
73: pU = PETSC_NULL;
74: PetscMalloc(sizeof(PetscScalar)*N*N,&pVT);
75: } else {
76: n = M;
77: PetscMalloc(sizeof(PetscScalar)*M*M,&pU);
78: pVT = PETSC_NULL;
79: }
80: PetscMalloc(sizeof(PetscReal)*n,&sigma);
81:
82: SVDDense(M,N,pmat,sigma,pU,pVT);
84: /* copy singular vectors */
85: for (i=0;i<n;i++) {
86: if (svd->which == SVD_SMALLEST) k = n - i - 1;
87: else k = i;
88: svd->sigma[k] = sigma[i];
89: VecGetArray(svd->U[k],&pu);
90: VecGetArray(svd->V[k],&pv);
91: if (M>=N) {
92: for (j=0;j<M;j++) pu[j] = pmat[i*M+j];
93: for (j=0;j<N;j++) pv[j] = pVT[j*N+i];
94: } else {
95: for (j=0;j<N;j++) pu[j] = pmat[j*M+i];
96: for (j=0;j<M;j++) pv[j] = pU[i*M+j];
97: }
98: VecRestoreArray(svd->U[k],&pu);
99: VecRestoreArray(svd->V[k],&pv);
100: }
102: svd->nconv = n;
103: svd->reason = SVD_CONVERGED_TOL;
105: MatRestoreArray(mat,&pmat);
106: MatDestroy(mat);
107: PetscFree(sigma);
108: if (M>=N) {
109: PetscFree(pVT);
110: } else {
111: PetscFree(pU);
112: }
113: return(0);
114: }
119: PetscErrorCode SVDCreate_LAPACK(SVD svd)
120: {
122: svd->ops->setup = SVDSetUp_LAPACK;
123: svd->ops->solve = SVDSolve_LAPACK;
124: svd->ops->destroy = SVDDestroy_Default;
125: if (svd->transmode == PETSC_DECIDE)
126: svd->transmode = SVD_TRANSPOSE_IMPLICIT; /* don't build the transpose */
127: return(0);
128: }