Actual source code: test7.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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 DSSVD.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: SlepcSC sc;
20: PetscReal sigma,rnorm,aux;
21: PetscScalar *A,*U,*w;
22: PetscInt i,j,k,n=15,m=10,ld;
23: PetscViewer viewer;
24: PetscBool verbose;
26: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
29: k = PetscMin(n,m);
30: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %Dx%D.\n",n,m);
31: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: /* Create DS object */
34: DSCreate(PETSC_COMM_WORLD,&ds);
35: DSSetType(ds,DSSVD);
36: DSSetFromOptions(ds);
37: ld = n+2; /* test leading dimension larger than n */
38: DSAllocate(ds,ld);
39: DSSetDimensions(ds,n,m,0,0);
41: /* Set up viewer */
42: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
43: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
44: DSView(ds,viewer);
45: PetscViewerPopFormat(viewer);
47: /* Fill with a rectangular Toeplitz matrix */
48: DSGetArray(ds,DS_MAT_A,&A);
49: for (i=0;i<k;i++) A[i+i*ld]=1.0;
50: for (j=1;j<3;j++) {
51: for (i=0;i<n-j;i++) { if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); }
52: }
53: for (j=1;j<n/2;j++) {
54: for (i=0;i<n-j;i++) { if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; }
55: }
56: DSRestoreArray(ds,DS_MAT_A,&A);
57: DSSetState(ds,DS_STATE_RAW);
58: if (verbose) {
59: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
60: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
61: }
62: DSView(ds,viewer);
64: /* Solve */
65: PetscMalloc1(k,&w);
66: DSGetSlepcSC(ds,&sc);
67: sc->comparison = SlepcCompareLargestReal;
68: sc->comparisonctx = NULL;
69: sc->map = NULL;
70: sc->mapobj = NULL;
71: DSSolve(ds,w,NULL);
72: DSSort(ds,w,NULL,NULL,NULL,NULL);
73: if (verbose) {
74: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
75: DSView(ds,viewer);
76: }
78: /* Print singular values */
79: PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
80: for (i=0;i<k;i++) {
81: sigma = PetscRealPart(w[i]);
82: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma);
83: }
85: /* Singular vectors */
86: DSVectors(ds,DS_MAT_U,NULL,NULL); /* all singular vectors */
87: j = 0;
88: rnorm = 0.0;
89: DSGetArray(ds,DS_MAT_U,&U);
90: for (i=0;i<n;i++) {
91: aux = PetscAbsScalar(U[i+j*ld]);
92: rnorm += aux*aux;
93: }
94: DSRestoreArray(ds,DS_MAT_U,&U);
95: rnorm = PetscSqrtReal(rnorm);
96: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st U vector = %.3f\n",(double)rnorm);
98: PetscFree(w);
99: DSDestroy(&ds);
100: SlepcFinalize();
101: return ierr;
102: }
104: /*TEST
106: test:
107: suffix: 1
108: requires: !single
110: TEST*/