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