Actual source code: test4.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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: */
11: static char help[] = "Test an SVD problem with more columns than rows.\n\n"
12: "The command line options are:\n"
13: " -m <m>, where <m> = matrix rows.\n"
14: " -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";
16: #include <slepcsvd.h>
18: /*
19: This example computes the singular values of a rectangular bidiagonal matrix
21: | 1 2 |
22: | 1 2 |
23: | 1 2 |
24: A = | . . |
25: | . . |
26: | 1 2 |
27: | 1 2 |
28: */
30: int main(int argc,char **argv)
31: {
32: Mat A,B;
33: SVD svd;
34: SVDConv conv;
35: SVDStop stop;
36: SVDWhich which;
37: SVDConvergedReason reason;
38: PetscInt m=20,n,Istart,Iend,i,col[2],its;
39: PetscScalar value[] = { 1, 2 };
40: PetscBool flg,tmode;
41: PetscErrorCode ierr;
42: PetscViewerAndFormat *vf;
43: const char *ctest[] = { "absolute", "relative to the singular value", "user-defined" };
44: const char *stest[] = { "basic", "user-defined" };
46: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
48: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
49: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
50: if (!flg) n=m+2;
51: PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Generate the matrix
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
59: MatSetFromOptions(A);
60: MatSetUp(A);
62: MatGetOwnershipRange(A,&Istart,&Iend);
63: for (i=Istart;i<Iend;i++) {
64: col[0]=i; col[1]=i+1;
65: if (i<n-1) {
66: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
67: } else if (i==n-1) {
68: MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
69: }
70: }
72: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
73: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Compute singular values
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: SVDCreate(PETSC_COMM_WORLD,&svd);
80: SVDSetOperator(svd,A);
82: /* test some interface functions */
83: SVDGetOperator(svd,&B);
84: MatView(B,PETSC_VIEWER_STDOUT_WORLD);
85: SVDSetConvergenceTest(svd,SVD_CONV_ABS);
86: SVDSetStoppingTest(svd,SVD_STOP_BASIC);
87: /* test monitors */
88: PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
89: SVDMonitorSet(svd,(PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*))SVDMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
90: /* SVDMonitorCancel(svd); */
91: SVDSetFromOptions(svd);
93: /* query properties and print them */
94: SVDGetImplicitTranspose(svd,&tmode);
95: PetscPrintf(PETSC_COMM_WORLD," Transpose mode is %s\n",tmode?"implicit":"explicit");
96: SVDGetConvergenceTest(svd,&conv);
97: PetscPrintf(PETSC_COMM_WORLD," Convergence test is %s\n",ctest[conv]);
98: SVDGetStoppingTest(svd,&stop);
99: PetscPrintf(PETSC_COMM_WORLD," Stopping test is %s\n",stest[stop]);
100: SVDGetWhichSingularTriplets(svd,&which);
101: PetscPrintf(PETSC_COMM_WORLD," Which = %s\n",which?"largest":"smallest");
103: /* call the solver */
104: SVDSolve(svd);
105: SVDGetConvergedReason(svd,&reason);
106: SVDGetIterationNumber(svd,&its);
107: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
108: /* PetscPrintf(PETSC_COMM_WORLD," its = %D\n",its); */
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Display solution and clean up
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
114: SVDDestroy(&svd);
115: MatDestroy(&A);
116: SlepcFinalize();
117: return ierr;
118: }