Actual source code: test24.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 and rectangular matrix A.\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,m,l=0,k=0,ld;
24: PetscViewer viewer;
25: PetscBool verbose,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+1,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,"-extrarow",&extrarow);
35: m = n+1;
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,m,l,k);
44: DSGSVDSetDimensions(ds,n,n);
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: /* Fill A and B with lower/upper arrow-bidiagonal matrices
55: verifying that [A;B] has orthonormal columns */
56: DSGetArrayReal(ds,DS_MAT_T,&T);
57: DSGetArrayReal(ds,DS_MAT_D,&D);
58: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1)/(n+1); /* diagonal of matrix A */
59: for (i=0;i<k;i++) D[i] = PetscSqrtReal(1.0-T[i]*T[i]);
60: for (i=l;i<k;i++) {
61: 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 */
62: T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
63: }
64: aux = 1.0-T[k]*T[k];
65: for (i=l;i<k;i++) aux -= T[i+ld]*T[i+ld]+T[i+2*ld]*T[i+2*ld];
66: T[k+ld] = PetscSqrtReal((1.0-aux)*.1);
67: aux -= T[k+ld]*T[k+ld];
68: D[k] = PetscSqrtReal(aux);
69: for (i=k+1;i<n;i++) {
70: T[i-1+2*ld] = -T[i-1+ld]*T[i]/D[i-1]; /* upper diagonal of matrix B */
71: aux = 1.0-T[i]*T[i]-T[2*ld+i-1]*T[2*ld+i-1];
72: T[i+ld] = PetscSqrtReal(aux)*.1; /* upper diagonal of matrix A */
73: D[i] = PetscSqrtReal(aux-T[i+ld]*T[i+ld]);
74: }
75: if (extrarow) { T[n]=-1.0; T[n-1+2*ld]=1.0; }
76: /* Fill locked eigenvalues */
77: PetscMalloc1(n,&w);
78: for (i=0;i<l;i++) w[i] = T[i]/D[i];
79: DSRestoreArrayReal(ds,DS_MAT_T,&T);
80: DSRestoreArrayReal(ds,DS_MAT_D,&D);
81: if (l==0 && k==0) DSSetState(ds,DS_STATE_INTERMEDIATE);
82: else DSSetState(ds,DS_STATE_RAW);
83: if (verbose) {
84: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
85: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
86: DSView(ds,viewer);
87: }
89: /* Solve */
90: DSGetSlepcSC(ds,&sc);
91: sc->comparison = SlepcCompareLargestReal;
92: sc->comparisonctx = NULL;
93: sc->map = NULL;
94: sc->mapobj = NULL;
95: DSSolve(ds,w,NULL);
96: DSSort(ds,w,NULL,NULL,NULL,NULL);
97: if (extrarow) DSUpdateExtraRow(ds);
98: DSSynchronize(ds,w,NULL);
99: if (verbose) {
100: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
101: DSView(ds,viewer);
102: }
104: /* Print singular values */
105: PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
106: for (i=0;i<n;i++) {
107: sigma = PetscRealPart(w[i]);
108: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma);
109: }
111: if (extrarow) {
112: /* Check that extra row is correct */
113: DSGetArrayReal(ds,DS_MAT_T,&T);
114: DSGetArray(ds,DS_MAT_U,&U);
115: DSGetArray(ds,DS_MAT_V,&V);
116: d = 0.0;
117: for (i=0;i<n;i++) d += T[i+ld]+U[n+i*ld];
118: 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));
119: d = 0.0;
120: for (i=0;i<n;i++) d += T[i+2*ld]-V[n-1+i*ld];
121: 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));
122: DSRestoreArrayReal(ds,DS_MAT_T,&T);
123: DSRestoreArray(ds,DS_MAT_U,&U);
124: DSRestoreArray(ds,DS_MAT_V,&V);
125: }
127: /* Singular vectors */
128: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all singular vectors */
129: DSGetMat(ds,DS_MAT_X,&X);
130: MatCreateVecs(X,NULL,&x0);
131: MatGetColumnVector(X,x0,0);
132: VecNorm(x0,NORM_2,&rnorm);
133: MatDestroy(&X);
134: VecDestroy(&x0);
135: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st X vector = %.3f\n",(double)rnorm);
137: DSGetMat(ds,DS_MAT_U,&X);
138: MatCreateVecs(X,NULL,&x0);
139: MatGetColumnVector(X,x0,0);
140: VecNorm(x0,NORM_2,&rnorm);
141: MatDestroy(&X);
142: VecDestroy(&x0);
143: if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st U vector has norm %g\n",(double)rnorm);
145: DSGetMat(ds,DS_MAT_V,&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 V vector has norm %g\n",(double)rnorm);
153: PetscFree(w);
154: DSDestroy(&ds);
155: SlepcFinalize();
156: return 0;
157: }
159: /*TEST
161: testset:
162: requires: double
163: output_file: output/test24_1.out
164: test:
165: suffix: 1
166: args: -l 1 -k 4
167: test:
168: suffix: 1_extrarow
169: filter: sed -e "s/extrarow//"
170: args: -l 1 -k 4 -extrarow
172: TEST*/