Actual source code: test8.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 with compact storage.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: SlepcSC sc;
20: PetscReal *T,sigma;
21: PetscScalar *w;
22: PetscInt i,n=10,m,l=2,k=5,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: m = n;
29: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD with compact storage - dimension %Dx%D.\n",n,m);
30: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
32: if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
33: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
35: /* Create DS object */
36: DSCreate(PETSC_COMM_WORLD,&ds);
37: DSSetType(ds,DSSVD);
38: DSSetFromOptions(ds);
39: ld = n+2; /* test leading dimension larger than n */
40: DSAllocate(ds,ld);
41: DSSetDimensions(ds,n,m,l,k);
42: DSSetCompact(ds,PETSC_TRUE);
44: /* Set up viewer */
45: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
46: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
47: DSView(ds,viewer);
48: PetscViewerPopFormat(viewer);
49: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
51: /* Fill upper arrow-tridiagonal matrix */
52: DSGetArrayReal(ds,DS_MAT_T,&T);
53: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
54: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
55: DSRestoreArrayReal(ds,DS_MAT_T,&T);
56: if (l==0 && k==0) {
57: DSSetState(ds,DS_STATE_INTERMEDIATE);
58: } else {
59: DSSetState(ds,DS_STATE_RAW);
60: }
61: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
62: DSView(ds,viewer);
64: /* Solve */
65: PetscMalloc1(n,&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<n;i++) {
81: sigma = PetscRealPart(w[i]);
82: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma);
83: }
84: PetscFree(w);
85: DSDestroy(&ds);
86: SlepcFinalize();
87: return ierr;
88: }
90: /*TEST
92: test:
93: suffix: 1
94: requires: !single
96: TEST*/