Actual source code: test22.c
slepc-3.17.1 2022-04-11
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 DSGSVD with compact storage.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: Mat X;
19: Vec x0;
20: SlepcSC sc;
21: PetscReal *T,*D,sigma,rnorm,aux;
22: PetscScalar *U,*V,*w,d;
23: PetscInt i,n=10,l=0,k=0,ld;
24: PetscViewer viewer;
25: PetscBool verbose,test_dsview,extrarow;
27: SlepcInitialize(&argc,&argv,(char*)0,help);
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GSVD with compact storage - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,n);
30: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
33: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
34: PetscOptionsHasName(NULL,NULL,"-test_dsview",&test_dsview);
35: PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);
37: /* Create DS object */
38: DSCreate(PETSC_COMM_WORLD,&ds);
39: DSSetType(ds,DSGSVD);
40: DSSetFromOptions(ds);
41: ld = n+2; /* test leading dimension larger than n */
42: DSAllocate(ds,ld);
43: DSSetDimensions(ds,n,l,k);
44: DSGSVDSetDimensions(ds,n,PETSC_DECIDE);
45: DSSetCompact(ds,PETSC_TRUE);
46: DSSetExtraRow(ds,extrarow);
48: /* Set up viewer */
49: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
50: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
51: DSView(ds,viewer);
52: PetscViewerPopFormat(viewer);
54: if (test_dsview) {
55: /* Fill A and B with dummy values to test DSView */
56: DSGetArrayReal(ds,DS_MAT_T,&T);
57: DSGetArrayReal(ds,DS_MAT_D,&D);
58: for (i=0;i<n;i++) { T[i] = i+1; D[i] = -i-1; }
59: for (i=0;i<n-1;i++) { T[i+ld] = -1.0; T[i+2*ld] = 1.0; }
60: DSRestoreArrayReal(ds,DS_MAT_T,&T);
61: DSRestoreArrayReal(ds,DS_MAT_D,&D);
62: DSView(ds,viewer);
63: }
65: /* Fill A and B with upper arrow-bidiagonal matrices
66: verifying that [A;B] has orthonormal columns */
67: DSGetArrayReal(ds,DS_MAT_T,&T);
68: DSGetArrayReal(ds,DS_MAT_D,&D);
69: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1)/(n+1); /* diagonal of matrix A */
70: for (i=0;i<k;i++) D[i] = PetscSqrtReal(1.0-T[i]*T[i]);
71: for (i=l;i<k;i++) {
72: T[i+ld] = PetscSqrtReal((1.0-T[k]*T[k])/(1.0+T[i]*T[i]/(D[i]*D[i])))*0.5*(1.0/k); /* upper diagonal of matrix A */
73: T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
74: }
75: aux = 1.0-T[k]*T[k];
76: for (i=l;i<k;i++) aux -= T[i+ld]*T[i+ld]+T[i+2*ld]*T[i+2*ld];
77: D[k] = PetscSqrtReal(aux);
78: for (i=k;i<n-1;i++) {
79: T[i+ld] = PetscSqrtReal((1.0-T[i+1]*T[i+1])/(1.0+T[i]*T[i]/(D[i]*D[i])))*0.5; /* upper diagonal of matrix A */
80: T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
81: D[i+1] = PetscSqrtReal(1.0-T[i+1]*T[i+1]-T[ld+i]*T[ld+i]-T[2*ld+i]*T[2*ld+i]); /* diagonal of matrix B */
82: }
83: if (extrarow) { T[n-1+ld]=-1.0; T[n-1+2*ld]=1.0; }
84: /* Fill locked eigenvalues */
85: PetscMalloc1(n,&w);
86: for (i=0;i<l;i++) w[i] = T[i]/D[i];
87: DSRestoreArrayReal(ds,DS_MAT_T,&T);
88: DSRestoreArrayReal(ds,DS_MAT_D,&D);
89: if (l==0 && k==0) DSSetState(ds,DS_STATE_INTERMEDIATE);
90: else DSSetState(ds,DS_STATE_RAW);
91: if (verbose) {
92: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
93: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
94: DSView(ds,viewer);
95: }
97: /* Solve */
98: DSGetSlepcSC(ds,&sc);
99: sc->comparison = SlepcCompareLargestReal;
100: sc->comparisonctx = NULL;
101: sc->map = NULL;
102: sc->mapobj = NULL;
103: DSSolve(ds,w,NULL);
104: DSSort(ds,w,NULL,NULL,NULL,NULL);
105: if (extrarow) DSUpdateExtraRow(ds);
106: DSSynchronize(ds,w,NULL);
107: if (verbose) {
108: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
109: DSView(ds,viewer);
110: }
112: /* Print singular values */
113: PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
114: for (i=0;i<n;i++) {
115: sigma = PetscRealPart(w[i]);
116: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma);
117: }
119: if (extrarow) {
120: /* Check that extra row is correct */
121: DSGetArrayReal(ds,DS_MAT_T,&T);
122: DSGetArray(ds,DS_MAT_U,&U);
123: DSGetArray(ds,DS_MAT_V,&V);
124: d = 0.0;
125: for (i=0;i<n;i++) d += T[i+ld]+U[n-1+i*ld];
126: if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in A's extra row of %g\n",(double)PetscAbsScalar(d));
127: d = 0.0;
128: for (i=0;i<n;i++) d += T[i+2*ld]-V[n-1+i*ld];
129: if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in B's extra row of %g\n",(double)PetscAbsScalar(d));
130: DSRestoreArrayReal(ds,DS_MAT_T,&T);
131: DSRestoreArray(ds,DS_MAT_U,&U);
132: DSRestoreArray(ds,DS_MAT_V,&V);
133: }
135: /* Singular vectors */
136: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all singular vectors */
137: DSGetMat(ds,DS_MAT_X,&X);
138: MatCreateVecs(X,NULL,&x0);
139: MatGetColumnVector(X,x0,0);
140: VecNorm(x0,NORM_2,&rnorm);
141: MatDestroy(&X);
142: VecDestroy(&x0);
143: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st X vector = %.3f\n",(double)rnorm);
145: DSGetMat(ds,DS_MAT_U,&X);
146: MatCreateVecs(X,NULL,&x0);
147: MatGetColumnVector(X,x0,0);
148: VecNorm(x0,NORM_2,&rnorm);
149: MatDestroy(&X);
150: VecDestroy(&x0);
151: if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st U vector has norm %g\n",(double)rnorm);
153: DSGetMat(ds,DS_MAT_V,&X);
154: MatCreateVecs(X,NULL,&x0);
155: MatGetColumnVector(X,x0,0);
156: VecNorm(x0,NORM_2,&rnorm);
157: MatDestroy(&X);
158: VecDestroy(&x0);
159: if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st V vector has norm %g\n",(double)rnorm);
161: PetscFree(w);
162: DSDestroy(&ds);
163: SlepcFinalize();
164: return 0;
165: }
167: /*TEST
169: testset:
170: requires: double
171: test:
172: suffix: 1
173: args: -test_dsview
174: test:
175: suffix: 2
176: args: -l 1 -k 4
177: test:
178: suffix: 2_extrarow
179: filter: sed -e "s/extrarow//"
180: args: -l 1 -k 4 -extrarow
181: output_file: output/test22_2.out
183: TEST*/