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