Actual source code: svdlapack.c
slepc-3.7.3 2016-09-29
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-2016, 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: if (svd->stop!=SVD_STOP_BASIC) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined stopping test not supported in this solver");
38: svd->max_it = 1;
39: svd->leftbasis = PETSC_TRUE;
40: SVDAllocateSolution(svd,0);
41: DSSetType(svd->ds,DSSVD);
42: DSAllocate(svd->ds,PetscMax(M,N));
43: return(0);
44: }
48: PetscErrorCode SVDSolve_LAPACK(SVD svd)
49: {
51: PetscInt M,N,n,i,j,k,ld;
52: Mat mat;
53: Vec u,v;
54: PetscScalar *pU,*pVT,*pmat,*pu,*pv,*A,*w;
57: DSGetLeadingDimension(svd->ds,&ld);
58: MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
59: MatGetSize(mat,&M,&N);
60: DSSetDimensions(svd->ds,M,N,0,0);
61: MatDenseGetArray(mat,&pmat);
62: DSGetArray(svd->ds,DS_MAT_A,&A);
63: for (i=0;i<M;i++)
64: for (j=0;j<N;j++)
65: A[i+j*ld] = pmat[i+j*M];
66: DSRestoreArray(svd->ds,DS_MAT_A,&A);
67: MatDenseRestoreArray(mat,&pmat);
68: DSSetState(svd->ds,DS_STATE_RAW);
70: n = PetscMin(M,N);
71: PetscMalloc1(n,&w);
72: DSSolve(svd->ds,w,NULL);
73: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
75: /* copy singular vectors */
76: DSGetArray(svd->ds,DS_MAT_U,&pU);
77: DSGetArray(svd->ds,DS_MAT_VT,&pVT);
78: for (i=0;i<n;i++) {
79: if (svd->which == SVD_SMALLEST) k = n - i - 1;
80: else k = i;
81: svd->sigma[k] = PetscRealPart(w[i]);
82: BVGetColumn(svd->U,k,&u);
83: BVGetColumn(svd->V,k,&v);
84: VecGetArray(u,&pu);
85: VecGetArray(v,&pv);
86: if (M>=N) {
87: for (j=0;j<M;j++) pu[j] = pU[i*ld+j];
88: for (j=0;j<N;j++) pv[j] = PetscConj(pVT[j*ld+i]);
89: } else {
90: for (j=0;j<N;j++) pu[j] = PetscConj(pVT[j*ld+i]);
91: for (j=0;j<M;j++) pv[j] = pU[i*ld+j];
92: }
93: VecRestoreArray(u,&pu);
94: VecRestoreArray(v,&pv);
95: BVRestoreColumn(svd->U,k,&u);
96: BVRestoreColumn(svd->V,k,&v);
97: }
98: DSRestoreArray(svd->ds,DS_MAT_U,&pU);
99: DSRestoreArray(svd->ds,DS_MAT_VT,&pVT);
101: svd->nconv = n;
102: svd->reason = SVD_CONVERGED_TOL;
104: MatDestroy(&mat);
105: PetscFree(w);
106: return(0);
107: }
111: PETSC_EXTERN PetscErrorCode SVDCreate_LAPACK(SVD svd)
112: {
114: svd->ops->setup = SVDSetUp_LAPACK;
115: svd->ops->solve = SVDSolve_LAPACK;
116: return(0);
117: }