Actual source code: test8.c
slepc-3.22.1 2024-10-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: {
17: DS ds;
18: SlepcSC sc;
19: PetscReal *T,sigma;
20: PetscScalar *U,*w,d;
21: PetscInt i,n=10,m,l=2,k=5,ld;
22: PetscViewer viewer;
23: PetscBool verbose,extrarow;
25: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28: m = n;
29: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD with compact storage - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
32: PetscCheck(l<=n && k<=n && l<=k,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
33: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
34: PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));
36: /* Create DS object */
37: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
38: PetscCall(DSSetType(ds,DSSVD));
39: PetscCall(DSSetFromOptions(ds));
40: ld = n+2; /* test leading dimension larger than n */
41: PetscCall(DSAllocate(ds,ld));
42: PetscCall(DSSetDimensions(ds,n,l,k));
43: PetscCall(DSSVDSetDimensions(ds,m));
44: PetscCall(DSSetCompact(ds,PETSC_TRUE));
45: PetscCall(DSSetExtraRow(ds,extrarow));
47: /* Set up viewer */
48: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
49: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
50: PetscCall(DSView(ds,viewer));
51: PetscCall(PetscViewerPopFormat(viewer));
52: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
54: /* Fill upper arrow-bidiagonal matrix */
55: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
56: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
57: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
58: if (extrarow) T[n-1+ld] = 1.0;
59: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
60: if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
61: else PetscCall(DSSetState(ds,DS_STATE_RAW));
62: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
63: PetscCall(DSView(ds,viewer));
65: /* Solve */
66: PetscCall(PetscMalloc1(n,&w));
67: PetscCall(DSGetSlepcSC(ds,&sc));
68: sc->comparison = SlepcCompareLargestReal;
69: sc->comparisonctx = NULL;
70: sc->map = NULL;
71: sc->mapobj = NULL;
72: PetscCall(DSSolve(ds,w,NULL));
73: PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
74: if (extrarow) PetscCall(DSUpdateExtraRow(ds));
75: if (verbose) {
76: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
77: PetscCall(DSView(ds,viewer));
78: }
80: /* Print singular values */
81: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
82: for (i=0;i<n;i++) {
83: sigma = PetscRealPart(w[i]);
84: PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma));
85: }
87: if (extrarow) {
88: /* Check that extra row is correct */
89: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
90: PetscCall(DSGetArray(ds,DS_MAT_U,&U));
91: d = 0.0;
92: for (i=l;i<n;i++) d += T[i+ld]-U[n-1+i*ld];
93: if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d)));
94: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
95: PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
96: }
97: PetscCall(PetscFree(w));
98: PetscCall(DSDestroy(&ds));
99: PetscCall(SlepcFinalize());
100: return 0;
101: }
103: /*TEST
105: test:
106: args: -ds_method {{0 1}}
107: suffix: 1
108: filter: grep -v "solving the problem"
109: requires: !single
111: test:
112: args: -l 0 -k 0 -ds_method {{0 1}}
113: suffix: 2
114: filter: grep -v "solving the problem"
115: requires: !single
117: test:
118: args: -extrarow -ds_method {{0 1}}
119: suffix: 3
120: filter: grep -v "solving the problem"
121: requires: !single
123: TEST*/