Actual source code: test26.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 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,reorthog,flg;
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 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));
36: PetscCall(PetscOptionsHasName(NULL,NULL,"-reorthog",&reorthog));
38: /* Create DS object */
39: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
40: PetscCall(DSSetType(ds,DSHSVD));
41: PetscCall(DSSetFromOptions(ds));
42: ld = n+2; /* test leading dimension larger than n */
43: PetscCall(DSAllocate(ds,ld));
44: PetscCall(DSSetDimensions(ds,n,l,k));
45: PetscCall(DSHSVDSetDimensions(ds,m));
46: PetscCall(DSSetCompact(ds,PETSC_TRUE));
47: PetscCall(DSSetExtraRow(ds,extrarow));
48: PetscCall(DSHSVDSetReorthogonalize(ds,reorthog));
49: PetscCall(DSHSVDGetReorthogonalize(ds,&flg));
50: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"reorthogonalizing\n"));
52: /* Set up viewer */
53: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
54: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
55: PetscCall(DSView(ds,viewer));
56: PetscCall(PetscViewerPopFormat(viewer));
57: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
59: /* Fill upper arrow-bidiagonal matrix and signature matrix */
60: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
61: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
62: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
63: if (extrarow) T[n-1+ld] = 1.0;
64: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
65: PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
66: for (i=0;i<n;i++) D[i] = (i>=l && i<l+p)? -1.0: 1.0;
67: PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
68: if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
69: else PetscCall(DSSetState(ds,DS_STATE_RAW));
70: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
71: PetscCall(DSView(ds,viewer));
73: /* Solve */
74: PetscCall(PetscMalloc1(n,&w));
75: PetscCall(DSGetSlepcSC(ds,&sc));
76: sc->comparison = SlepcCompareLargestReal;
77: sc->comparisonctx = NULL;
78: sc->map = NULL;
79: sc->mapobj = NULL;
80: PetscCall(DSSolve(ds,w,NULL));
81: PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
82: if (extrarow) PetscCall(DSUpdateExtraRow(ds));
83: if (verbose) {
84: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
85: PetscCall(DSView(ds,viewer));
86: }
88: /* Print singular values */
89: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
90: for (i=0;i<n;i++) {
91: sigma = PetscRealPart(w[i]);
92: PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma));
93: }
95: if (extrarow) {
96: /* Check that extra row is correct */
97: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
98: PetscCall(DSGetArray(ds,DS_MAT_U,&U));
99: PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
100: d = 0.0;
101: for (i=l;i<n;i++) d += T[i+ld]-D[i]*U[n-1+i*ld];
102: 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)));
103: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
104: PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
105: PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
106: }
107: PetscCall(PetscFree(w));
108: PetscCall(DSDestroy(&ds));
109: PetscCall(SlepcFinalize());
110: return 0;
111: }
113: /*TEST
115: test:
116: suffix: 1
117: requires: !single
119: test:
120: args: -l 0 -k 0
121: suffix: 2
122: requires: !single
124: test:
125: args: -extrarow -reorthog {{0 1}}
126: suffix: 3
127: requires: !single
128: filter: grep -v reorthogonalizing
130: TEST*/