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-2011, Universitat 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> /*I "slepcsvd.h" I*/
25: #include <slepcblaslapack.h>
29: PetscErrorCode SVDSetUp_LAPACK(SVD svd)
30: {
32: PetscInt N;
35: SVDMatGetSize(svd,PETSC_NULL,&N);
36: svd->ncv = N;
37: if (svd->mpd) PetscInfo(svd,"Warning: parameter mpd ignored\n");
38: svd->max_it = 1;
39: if (svd->ncv!=svd->n) {
40: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
41: }
42: return(0);
43: }
47: PetscErrorCode SVDSolve_LAPACK(SVD svd)
48: {
50: PetscInt M,N,n,i,j,k;
51: Mat mat;
52: PetscScalar *pU,*pVT,*pmat,*pu,*pv;
53: PetscReal *sigma;
54:
56: MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
57: MatGetArray(mat,&pmat);
58: MatGetSize(mat,&M,&N);
59: if (M>=N) {
60: n = N;
61: pU = PETSC_NULL;
62: PetscMalloc(sizeof(PetscScalar)*N*N,&pVT);
63: } else {
64: n = M;
65: PetscMalloc(sizeof(PetscScalar)*M*M,&pU);
66: pVT = PETSC_NULL;
67: }
68: PetscMalloc(sizeof(PetscReal)*n,&sigma);
69:
70: SVDDense(M,N,pmat,sigma,pU,pVT);
72: /* copy singular vectors */
73: for (i=0;i<n;i++) {
74: if (svd->which == SVD_SMALLEST) k = n - i - 1;
75: else k = i;
76: svd->sigma[k] = sigma[i];
77: VecGetArray(svd->U[k],&pu);
78: VecGetArray(svd->V[k],&pv);
79: if (M>=N) {
80: for (j=0;j<M;j++) pu[j] = pmat[i*M+j];
81: for (j=0;j<N;j++) pv[j] = pVT[j*N+i];
82: } else {
83: for (j=0;j<N;j++) pu[j] = pmat[j*M+i];
84: for (j=0;j<M;j++) pv[j] = pU[i*M+j];
85: }
86: VecRestoreArray(svd->U[k],&pu);
87: VecRestoreArray(svd->V[k],&pv);
88: }
90: svd->nconv = n;
91: svd->reason = SVD_CONVERGED_TOL;
93: MatRestoreArray(mat,&pmat);
94: MatDestroy(&mat);
95: PetscFree(sigma);
96: if (M>=N) {
97: PetscFree(pVT);
98: } else {
99: PetscFree(pU);
100: }
101: return(0);
102: }
106: PetscErrorCode SVDReset_LAPACK(SVD svd)
107: {
111: VecDestroyVecs(svd->n,&svd->U);
112: return(0);
113: }
117: PetscErrorCode SVDDestroy_LAPACK(SVD svd)
118: {
122: PetscFree(svd->data);
123: return(0);
124: }
126: EXTERN_C_BEGIN
129: PetscErrorCode SVDCreate_LAPACK(SVD svd)
130: {
132: svd->ops->setup = SVDSetUp_LAPACK;
133: svd->ops->solve = SVDSolve_LAPACK;
134: svd->ops->destroy = SVDDestroy_LAPACK;
135: svd->ops->reset = SVDReset_LAPACK;
136: if (svd->transmode == PETSC_DECIDE)
137: svd->transmode = SVD_TRANSPOSE_IMPLICIT; /* don't build the transpose */
138: return(0);
139: }
140: EXTERN_C_END