Actual source code: test26.c
slepc-3.20.1 2023-11-27
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 DSHSVD 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,*D,sigma;
20: PetscScalar *U,*w,d;
21: PetscInt i,n=10,m,l=2,k=5,p=1,ld;
22: PetscViewer viewer;
23: PetscBool verbose,extrarow;
25: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28: m = n;
29: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HSVD 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: PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
33: PetscCheck(l<=n && k<=n && l<=k && p>=0 && p<=n-l,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
34: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
35: PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));
37: /* Create DS object */
38: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
39: PetscCall(DSSetType(ds,DSHSVD));
40: PetscCall(DSSetFromOptions(ds));
41: ld = n+2; /* test leading dimension larger than n */
42: PetscCall(DSAllocate(ds,ld));
43: PetscCall(DSSetDimensions(ds,n,l,k));
44: PetscCall(DSHSVDSetDimensions(ds,m));
45: PetscCall(DSSetCompact(ds,PETSC_TRUE));
46: PetscCall(DSSetExtraRow(ds,extrarow));
48: /* Set up viewer */
49: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
50: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
51: PetscCall(DSView(ds,viewer));
52: PetscCall(PetscViewerPopFormat(viewer));
53: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
55: /* Fill upper arrow-bidiagonal matrix and signature matrix */
56: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
57: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
58: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
59: if (extrarow) T[n-1+ld] = 1.0;
60: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
61: PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
62: for (i=0;i<n;i++) D[i] = (i>=l && i<l+p)? -1.0: 1.0;
63: PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
64: if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
65: else PetscCall(DSSetState(ds,DS_STATE_RAW));
66: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
67: PetscCall(DSView(ds,viewer));
69: /* Solve */
70: PetscCall(PetscMalloc1(n,&w));
71: PetscCall(DSGetSlepcSC(ds,&sc));
72: sc->comparison = SlepcCompareLargestReal;
73: sc->comparisonctx = NULL;
74: sc->map = NULL;
75: sc->mapobj = NULL;
76: PetscCall(DSSolve(ds,w,NULL));
77: PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
78: if (extrarow) PetscCall(DSUpdateExtraRow(ds));
79: if (verbose) {
80: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
81: PetscCall(DSView(ds,viewer));
82: }
84: /* Print singular values */
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
86: for (i=0;i<n;i++) {
87: sigma = PetscRealPart(w[i]);
88: PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma));
89: }
91: if (extrarow) {
92: /* Check that extra row is correct */
93: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
94: PetscCall(DSGetArray(ds,DS_MAT_U,&U));
95: PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
96: d = 0.0;
97: for (i=l;i<n;i++) d += T[i+ld]-D[i]*U[n-1+i*ld];
98: 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)));
99: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
100: PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
101: PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
102: }
103: PetscCall(PetscFree(w));
104: PetscCall(DSDestroy(&ds));
105: PetscCall(SlepcFinalize());
106: return 0;
107: }
109: /*TEST
111: test:
112: suffix: 1
113: requires: !single
115: test:
116: args: -l 0 -k 0
117: suffix: 2
118: requires: !single
120: test:
121: args: -extrarow
122: suffix: 3
123: requires: !single
125: TEST*/