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.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 12 : int main(int argc,char **argv)
16 : {
17 12 : DS ds;
18 12 : SlepcSC sc;
19 12 : Mat X;
20 12 : Vec x0;
21 12 : PetscReal sigma,rnorm,cond;
22 12 : PetscScalar *A,*B,*w;
23 12 : PetscInt i,j,k,n=15,m=10,p=10,m1,p1,ld;
24 12 : PetscViewer viewer;
25 12 : PetscBool verbose;
26 :
27 12 : PetscFunctionBeginUser;
28 12 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
31 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
32 12 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GSVD - dimension (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT ".\n",n,p,m));
33 12 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
34 :
35 : /* Create DS object */
36 12 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
37 12 : PetscCall(DSSetType(ds,DSGSVD));
38 12 : PetscCall(DSSetFromOptions(ds));
39 12 : ld = PetscMax(PetscMax(p,m),n)+2; /* test leading dimension larger than n */
40 12 : PetscCall(DSAllocate(ds,ld));
41 12 : PetscCall(DSSetDimensions(ds,n,0,0));
42 12 : PetscCall(DSGSVDSetDimensions(ds,m,p));
43 12 : PetscCall(DSGSVDGetDimensions(ds,&m1,&p1));
44 12 : PetscCheck(m1==m && p1==p,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Inconsistent dimension values");
45 :
46 : /* Set up viewer */
47 12 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
48 12 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
49 12 : PetscCall(DSView(ds,viewer));
50 12 : PetscCall(PetscViewerPopFormat(viewer));
51 :
52 12 : k = PetscMin(n,m);
53 : /* Fill A with a rectangular Toeplitz matrix */
54 12 : PetscCall(DSGetArray(ds,DS_MAT_A,&A));
55 132 : for (i=0;i<k;i++) A[i+i*ld]=1.0;
56 36 : for (j=1;j<3;j++) {
57 348 : for (i=0;i<n-j;i++) { if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); }
58 : }
59 84 : for (j=1;j<n/2;j++) {
60 900 : for (i=0;i<n-j;i++) { if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; }
61 : }
62 12 : PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
63 :
64 12 : k = PetscMin(p,m);
65 : /* Fill B with a shifted bidiagonal matrix */
66 12 : PetscCall(DSGetArray(ds,DS_MAT_B,&B));
67 132 : for (i=m-k;i<m;i++) {
68 120 : B[i-m+k+i*ld]=2.0-1.0/(PetscScalar)(i+1);
69 120 : if (i) B[i-1-m+k+i*ld]=1.0;
70 : }
71 12 : PetscCall(DSRestoreArray(ds,DS_MAT_B,&B));
72 :
73 12 : PetscCall(DSSetState(ds,DS_STATE_RAW));
74 12 : if (verbose) {
75 0 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
76 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
77 : }
78 12 : PetscCall(DSView(ds,viewer));
79 :
80 : /* Condition number */
81 12 : PetscCall(DSCond(ds,&cond));
82 12 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Condition number = %.3f\n",(double)cond));
83 :
84 : /* Solve */
85 12 : PetscCall(PetscMalloc1(m,&w));
86 12 : PetscCall(DSGetSlepcSC(ds,&sc));
87 12 : sc->comparison = SlepcCompareLargestReal;
88 12 : sc->comparisonctx = NULL;
89 12 : sc->map = NULL;
90 12 : sc->mapobj = NULL;
91 12 : PetscCall(DSSolve(ds,w,NULL));
92 12 : PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
93 12 : PetscCall(DSSynchronize(ds,w,NULL));
94 12 : if (verbose) {
95 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
96 0 : PetscCall(DSView(ds,viewer));
97 : }
98 : /* Print singular values */
99 12 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
100 12 : PetscCall(DSGetDimensions(ds,NULL,NULL,NULL,&k));
101 132 : for (i=0;i<k;i++) {
102 120 : sigma = PetscRealPart(w[i]);
103 120 : PetscCall(PetscViewerASCIIPrintf(viewer," %g\n",(double)sigma));
104 : }
105 :
106 : /* Singular vectors */
107 12 : PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all singular vectors */
108 12 : PetscCall(DSGetMat(ds,DS_MAT_X,&X));
109 12 : PetscCall(MatCreateVecs(X,NULL,&x0));
110 12 : PetscCall(MatGetColumnVector(X,x0,0));
111 12 : PetscCall(VecNorm(x0,NORM_2,&rnorm));
112 12 : PetscCall(DSRestoreMat(ds,DS_MAT_X,&X));
113 12 : PetscCall(VecDestroy(&x0));
114 12 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st X vector = %.3f\n",(double)rnorm));
115 :
116 12 : PetscCall(DSGetMat(ds,DS_MAT_U,&X));
117 12 : PetscCall(MatCreateVecs(X,NULL,&x0));
118 12 : PetscCall(MatGetColumnVector(X,x0,0));
119 12 : PetscCall(VecNorm(x0,NORM_2,&rnorm));
120 12 : PetscCall(DSRestoreMat(ds,DS_MAT_U,&X));
121 12 : PetscCall(VecDestroy(&x0));
122 12 : 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));
123 :
124 12 : PetscCall(DSGetMat(ds,DS_MAT_V,&X));
125 12 : PetscCall(MatCreateVecs(X,NULL,&x0));
126 12 : PetscCall(MatGetColumnVector(X,x0,0));
127 12 : PetscCall(VecNorm(x0,NORM_2,&rnorm));
128 12 : PetscCall(DSRestoreMat(ds,DS_MAT_V,&X));
129 12 : PetscCall(VecDestroy(&x0));
130 12 : 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));
131 :
132 12 : PetscCall(PetscFree(w));
133 12 : PetscCall(DSDestroy(&ds));
134 12 : PetscCall(SlepcFinalize());
135 : return 0;
136 : }
137 :
138 : /*TEST
139 :
140 : testset:
141 : output_file: output/test21_1.out
142 : requires: !single
143 : nsize: {{1 2 3}}
144 : filter: grep -v "parallel operation mode" | grep -v " MPI process"
145 : test:
146 : suffix: 1
147 : args: -ds_parallel redundant
148 : test:
149 : suffix: 2
150 : args: -ds_parallel synchronized
151 :
152 : TEST*/
|