Actual source code: svdlapack.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This file implements a wrapper to the LAPACK SVD subroutines
12: */
14: #include <slepc/private/svdimpl.h>
15: #include <slepcblaslapack.h>
17: PetscErrorCode SVDSetUp_LAPACK(SVD svd)
18: {
20: PetscInt M,N,P=0;
23: MatGetSize(svd->A,&M,&N);
24: if (!svd->isgeneralized) svd->ncv = N;
25: else {
26: MatGetSize(svd->OPb,&P,NULL);
27: svd->ncv = PetscMin(M,PetscMin(N,P));
28: }
29: if (svd->mpd!=PETSC_DEFAULT) { PetscInfo(svd,"Warning: parameter mpd ignored\n"); }
30: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
31: if (svd->max_it==PETSC_DEFAULT) svd->max_it = 1;
32: svd->leftbasis = PETSC_TRUE;
33: SVDAllocateSolution(svd,0);
34: DSSetType(svd->ds,svd->isgeneralized?DSGSVD:DSSVD);
35: DSAllocate(svd->ds,PetscMax(N,PetscMax(M,P)));
36: return(0);
37: }
39: PetscErrorCode SVDSolve_LAPACK(SVD svd)
40: {
41: PetscErrorCode ierr;
42: PetscInt M,N,n,i,j,k,ld,lowu,lowv,highu,highv;
43: Mat Ar,mat;
44: Vec u,v;
45: PetscScalar *pU,*pV,*pu,*pv,*A,*w;
46: const PetscScalar *pmat;
49: DSGetLeadingDimension(svd->ds,&ld);
50: MatCreateRedundantMatrix(svd->OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
51: MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
52: MatDestroy(&Ar);
53: MatGetSize(mat,&M,&N);
54: DSSetDimensions(svd->ds,M,0,0);
55: DSSVDSetDimensions(svd->ds,N);
56: MatDenseGetArrayRead(mat,&pmat);
57: DSGetArray(svd->ds,DS_MAT_A,&A);
58: for (i=0;i<M;i++)
59: for (j=0;j<N;j++)
60: A[i+j*ld] = pmat[i+j*M];
61: DSRestoreArray(svd->ds,DS_MAT_A,&A);
62: MatDenseRestoreArrayRead(mat,&pmat);
63: DSSetState(svd->ds,DS_STATE_RAW);
65: n = PetscMin(M,N);
66: PetscMalloc1(n,&w);
67: DSSolve(svd->ds,w,NULL);
68: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
69: DSSynchronize(svd->ds,w,NULL);
71: /* copy singular vectors */
72: DSGetArray(svd->ds,DS_MAT_U,&pU);
73: DSGetArray(svd->ds,DS_MAT_V,&pV);
74: for (i=0;i<n;i++) {
75: if (svd->which == SVD_SMALLEST) k = n - i - 1;
76: else k = i;
77: svd->sigma[k] = PetscRealPart(w[i]);
78: BVGetColumn(svd->U,k,&u);
79: BVGetColumn(svd->V,k,&v);
80: VecGetOwnershipRange(u,&lowu,&highu);
81: VecGetOwnershipRange(v,&lowv,&highv);
82: VecGetArray(u,&pu);
83: VecGetArray(v,&pv);
84: if (M>=N) {
85: for (j=lowu;j<highu;j++) pu[j-lowu] = pU[i*ld+j];
86: for (j=lowv;j<highv;j++) pv[j-lowv] = pV[i*ld+j];
87: } else {
88: for (j=lowu;j<highu;j++) pu[j-lowu] = pV[i*ld+j];
89: for (j=lowv;j<highv;j++) pv[j-lowv] = pU[i*ld+j];
90: }
91: VecRestoreArray(u,&pu);
92: VecRestoreArray(v,&pv);
93: BVRestoreColumn(svd->U,k,&u);
94: BVRestoreColumn(svd->V,k,&v);
95: }
96: DSRestoreArray(svd->ds,DS_MAT_U,&pU);
97: DSRestoreArray(svd->ds,DS_MAT_V,&pV);
99: svd->nconv = n;
100: svd->its = 1;
101: svd->reason = SVD_CONVERGED_TOL;
103: MatDestroy(&mat);
104: PetscFree(w);
105: return(0);
106: }
108: PetscErrorCode SVDSolve_LAPACK_GSVD(SVD svd)
109: {
110: PetscErrorCode ierr;
111: PetscInt nsv,m,n,p,i,j,mlocal,plocal,ld,lowx,lowu,lowv,highx;
112: Mat Ar,A,Br,B;
113: Vec uv,x;
114: PetscScalar *Ads,*Bds,*U,*V,*X,*px,*puv,*w;
115: const PetscScalar *pA,*pB;
118: DSGetLeadingDimension(svd->ds,&ld);
119: MatCreateRedundantMatrix(svd->OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
120: MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&A);
121: MatDestroy(&Ar);
122: MatCreateRedundantMatrix(svd->OPb,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Br);
123: MatConvert(Br,MATSEQDENSE,MAT_INITIAL_MATRIX,&B);
124: MatDestroy(&Br);
125: MatGetSize(A,&m,&n);
126: MatGetLocalSize(svd->OP,&mlocal,NULL);
127: MatGetLocalSize(svd->OPb,&plocal,NULL);
128: MatGetSize(B,&p,NULL);
129: DSSetDimensions(svd->ds,m,0,0);
130: DSGSVDSetDimensions(svd->ds,n,p);
131: MatDenseGetArrayRead(A,&pA);
132: MatDenseGetArrayRead(B,&pB);
133: DSGetArray(svd->ds,DS_MAT_A,&Ads);
134: DSGetArray(svd->ds,DS_MAT_B,&Bds);
135: for (j=0;j<n;j++) {
136: for (i=0;i<m;i++) Ads[i+j*ld] = pA[i+j*m];
137: for (i=0;i<p;i++) Bds[i+j*ld] = pB[i+j*p];
138: }
139: DSRestoreArray(svd->ds,DS_MAT_B,&Bds);
140: DSRestoreArray(svd->ds,DS_MAT_A,&Ads);
141: MatDenseRestoreArrayRead(B,&pB);
142: MatDenseRestoreArrayRead(A,&pA);
143: DSSetState(svd->ds,DS_STATE_RAW);
145: nsv = PetscMin(n,PetscMin(p,m));
146: PetscMalloc1(nsv,&w);
147: DSSolve(svd->ds,w,NULL);
148: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
149: DSSynchronize(svd->ds,w,NULL);
150: DSGetDimensions(svd->ds,NULL,NULL,NULL,&nsv);
152: /* copy singular vectors */
153: MatGetOwnershipRange(svd->OP,&lowu,NULL);
154: MatGetOwnershipRange(svd->OPb,&lowv,NULL);
155: DSGetArray(svd->ds,DS_MAT_X,&X);
156: DSGetArray(svd->ds,DS_MAT_U,&U);
157: DSGetArray(svd->ds,DS_MAT_V,&V);
158: for (j=0;j<nsv;j++) {
159: svd->sigma[j] = PetscRealPart(w[j]);
160: BVGetColumn(svd->V,j,&x);
161: VecGetOwnershipRange(x,&lowx,&highx);
162: VecGetArrayWrite(x,&px);
163: for (i=lowx;i<highx;i++) px[i-lowx] = X[i+j*ld];
164: VecRestoreArrayWrite(x,&px);
165: BVRestoreColumn(svd->V,j,&x);
166: BVGetColumn(svd->U,j,&uv);
167: VecGetArrayWrite(uv,&puv);
168: for (i=0;i<mlocal;i++) puv[i] = U[i+lowu+j*ld];
169: for (i=0;i<plocal;i++) puv[i+mlocal] = V[i+lowv+j*ld];
170: VecRestoreArrayWrite(uv,&puv);
171: BVRestoreColumn(svd->U,j,&uv);
172: }
173: DSRestoreArray(svd->ds,DS_MAT_X,&X);
174: DSRestoreArray(svd->ds,DS_MAT_U,&U);
175: DSRestoreArray(svd->ds,DS_MAT_V,&V);
177: svd->nconv = nsv;
178: svd->its = 1;
179: svd->reason = SVD_CONVERGED_TOL;
181: MatDestroy(&A);
182: MatDestroy(&B);
183: PetscFree(w);
184: return(0);
185: }
187: SLEPC_EXTERN PetscErrorCode SVDCreate_LAPACK(SVD svd)
188: {
190: svd->ops->setup = SVDSetUp_LAPACK;
191: svd->ops->solve = SVDSolve_LAPACK;
192: svd->ops->solveg = SVDSolve_LAPACK_GSVD;
193: return(0);
194: }