Actual source code: test7.c
slepc-3.7.0 2016-05-16
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test DSSVD.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: SlepcSC sc;
33: PetscReal sigma;
34: PetscScalar *A,*w;
35: PetscInt i,j,k,n=15,m=10,ld;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
42: k = PetscMin(n,m);
43: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %Dx%D.\n",n,m);
44: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
46: /* Create DS object */
47: DSCreate(PETSC_COMM_WORLD,&ds);
48: DSSetType(ds,DSSVD);
49: DSSetFromOptions(ds);
50: ld = n+2; /* test leading dimension larger than n */
51: DSAllocate(ds,ld);
52: DSSetDimensions(ds,n,m,0,0);
54: /* Set up viewer */
55: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
56: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
57: DSView(ds,viewer);
58: PetscViewerPopFormat(viewer);
59: if (verbose) {
60: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
61: }
63: /* Fill with a rectangular Toeplitz matrix */
64: DSGetArray(ds,DS_MAT_A,&A);
65: for (i=0;i<k;i++) A[i+i*ld]=1.0;
66: for (j=1;j<3;j++) {
67: for (i=0;i<n-j;i++) { if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); }
68: }
69: for (j=1;j<n/2;j++) {
70: for (i=0;i<n-j;i++) { if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; }
71: }
72: DSRestoreArray(ds,DS_MAT_A,&A);
73: DSSetState(ds,DS_STATE_RAW);
74: if (verbose) {
75: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
76: DSView(ds,viewer);
77: }
79: /* Solve */
80: PetscMalloc1(k,&w);
81: DSGetSlepcSC(ds,&sc);
82: sc->comparison = SlepcCompareLargestReal;
83: sc->comparisonctx = NULL;
84: sc->map = NULL;
85: sc->mapobj = NULL;
86: DSSolve(ds,w,NULL);
87: DSSort(ds,w,NULL,NULL,NULL,NULL);
88: if (verbose) {
89: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
90: DSView(ds,viewer);
91: }
93: /* Print singular values */
94: PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
95: for (i=0;i<k;i++) {
96: sigma = PetscRealPart(w[i]);
97: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma);
98: }
99: PetscFree(w);
100: DSDestroy(&ds);
101: SlepcFinalize();
102: return ierr;
103: }