Actual source code: svdlapack.c
slepc-3.6.1 2015-09-03
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-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepc/private/svdimpl.h>
28: PetscErrorCode SVDSetUp_LAPACK(SVD svd)
29: {
31: PetscInt M,N;
34: SVDMatGetSize(svd,&M,&N);
35: svd->ncv = N;
36: if (svd->mpd) { PetscInfo(svd,"Warning: parameter mpd ignored\n"); }
37: svd->max_it = 1;
38: svd->leftbasis = PETSC_TRUE;
39: SVDAllocateSolution(svd,0);
40: DSSetType(svd->ds,DSSVD);
41: DSAllocate(svd->ds,PetscMax(M,N));
42: return(0);
43: }
47: PetscErrorCode SVDSolve_LAPACK(SVD svd)
48: {
50: PetscInt M,N,n,i,j,k,ld;
51: Mat mat;
52: Vec u,v;
53: PetscScalar *pU,*pVT,*pmat,*pu,*pv,*A,*w;
56: DSGetLeadingDimension(svd->ds,&ld);
57: MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
58: MatGetSize(mat,&M,&N);
59: DSSetDimensions(svd->ds,M,N,0,0);
60: MatDenseGetArray(mat,&pmat);
61: DSGetArray(svd->ds,DS_MAT_A,&A);
62: for (i=0;i<M;i++)
63: for (j=0;j<N;j++)
64: A[i+j*ld] = pmat[i+j*M];
65: DSRestoreArray(svd->ds,DS_MAT_A,&A);
66: MatDenseRestoreArray(mat,&pmat);
67: DSSetState(svd->ds,DS_STATE_RAW);
69: n = PetscMin(M,N);
70: PetscMalloc1(n,&w);
71: DSSolve(svd->ds,w,NULL);
72: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
74: /* copy singular vectors */
75: DSGetArray(svd->ds,DS_MAT_U,&pU);
76: DSGetArray(svd->ds,DS_MAT_VT,&pVT);
77: for (i=0;i<n;i++) {
78: if (svd->which == SVD_SMALLEST) k = n - i - 1;
79: else k = i;
80: svd->sigma[k] = PetscRealPart(w[i]);
81: BVGetColumn(svd->U,k,&u);
82: BVGetColumn(svd->V,k,&v);
83: VecGetArray(u,&pu);
84: VecGetArray(v,&pv);
85: if (M>=N) {
86: for (j=0;j<M;j++) pu[j] = pU[i*ld+j];
87: for (j=0;j<N;j++) pv[j] = PetscConj(pVT[j*ld+i]);
88: } else {
89: for (j=0;j<N;j++) pu[j] = PetscConj(pVT[j*ld+i]);
90: for (j=0;j<M;j++) pv[j] = pU[i*ld+j];
91: }
92: VecRestoreArray(u,&pu);
93: VecRestoreArray(v,&pv);
94: BVRestoreColumn(svd->U,k,&u);
95: BVRestoreColumn(svd->V,k,&v);
96: }
97: DSRestoreArray(svd->ds,DS_MAT_U,&pU);
98: DSRestoreArray(svd->ds,DS_MAT_VT,&pVT);
100: svd->nconv = n;
101: svd->reason = SVD_CONVERGED_TOL;
103: MatDestroy(&mat);
104: PetscFree(w);
105: return(0);
106: }
110: PETSC_EXTERN PetscErrorCode SVDCreate_LAPACK(SVD svd)
111: {
113: svd->ops->setup = SVDSetUp_LAPACK;
114: svd->ops->solve = SVDSolve_LAPACK;
115: return(0);
116: }