Actual source code: test22.c

slepc-3.16.0 2021-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {
 18:   DS             ds;
 19:   Mat            X;
 20:   Vec            x0;
 21:   SlepcSC        sc;
 22:   PetscReal      *T,*D,sigma,rnorm,aux;
 23:   PetscScalar    *U,*V,*w,d;
 24:   PetscInt       i,n=10,l=0,k=0,ld;
 25:   PetscViewer    viewer;
 26:   PetscBool      verbose,test_dsview,extrarow;

 28:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 29:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 30:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GSVD with compact storage - dimension %Dx%D.\n",n,n);
 31:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 33:   if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
 34:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 35:   PetscOptionsHasName(NULL,NULL,"-test_dsview",&test_dsview);
 36:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

 38:   /* Create DS object */
 39:   DSCreate(PETSC_COMM_WORLD,&ds);
 40:   DSSetType(ds,DSGSVD);
 41:   DSSetFromOptions(ds);
 42:   ld = n+2;  /* test leading dimension larger than n */
 43:   DSAllocate(ds,ld);
 44:   DSSetDimensions(ds,n,l,k);
 45:   DSGSVDSetDimensions(ds,n,PETSC_DECIDE);
 46:   DSSetCompact(ds,PETSC_TRUE);
 47:   DSSetExtraRow(ds,extrarow);

 49:   /* Set up viewer */
 50:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 51:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 52:   DSView(ds,viewer);
 53:   PetscViewerPopFormat(viewer);

 55:   if (test_dsview) {
 56:     /* Fill A and B with dummy values to test DSView */
 57:     DSGetArrayReal(ds,DS_MAT_T,&T);
 58:     DSGetArrayReal(ds,DS_MAT_D,&D);
 59:     for (i=0;i<n;i++) { T[i] = i+1; D[i] = -i-1; }
 60:     for (i=0;i<n-1;i++) { T[i+ld] = -1.0; T[i+2*ld] = 1.0; }
 61:     DSRestoreArrayReal(ds,DS_MAT_T,&T);
 62:     DSRestoreArrayReal(ds,DS_MAT_D,&D);
 63:     DSView(ds,viewer);
 64:   }

 66:   /* Fill A and B with upper arrow-bidiagonal matrices
 67:      verifying that [A;B] has orthonormal columns */
 68:   DSGetArrayReal(ds,DS_MAT_T,&T);
 69:   DSGetArrayReal(ds,DS_MAT_D,&D);
 70:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1)/(n+1); /* diagonal of matrix A */
 71:   for (i=0;i<k;i++) D[i] = PetscSqrtReal(1.0-T[i]*T[i]);
 72:   for (i=l;i<k;i++) {
 73:     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 */
 74:     T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
 75:   }
 76:   aux = 1.0-T[k]*T[k];
 77:   for (i=l;i<k;i++) aux -= T[i+ld]*T[i+ld]+T[i+2*ld]*T[i+2*ld];
 78:   D[k] = PetscSqrtReal(aux);
 79:   for (i=k;i<n-1;i++) {
 80:     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 */
 81:     T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
 82:     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 */
 83:   }
 84:   if (extrarow) { T[n-1+ld]=-1.0; T[n-1+2*ld]=1.0; }
 85:   /* Fill locked eigenvalues */
 86:   PetscMalloc1(n,&w);
 87:   for (i=0;i<l;i++) w[i] = T[i]/D[i];
 88:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 89:   DSRestoreArrayReal(ds,DS_MAT_D,&D);
 90:   if (l==0 && k==0) {
 91:     DSSetState(ds,DS_STATE_INTERMEDIATE);
 92:   } else {
 93:     DSSetState(ds,DS_STATE_RAW);
 94:   }
 95:   if (verbose) {
 96:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 97:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 98:     DSView(ds,viewer);
 99:   }

101:   /* Solve */
102:   DSGetSlepcSC(ds,&sc);
103:   sc->comparison    = SlepcCompareLargestReal;
104:   sc->comparisonctx = NULL;
105:   sc->map           = NULL;
106:   sc->mapobj        = NULL;
107:   DSSolve(ds,w,NULL);
108:   DSSort(ds,w,NULL,NULL,NULL,NULL);
109:   if (extrarow) { DSUpdateExtraRow(ds); }
110:   DSSynchronize(ds,w,NULL);
111:   if (verbose) {
112:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
113:     DSView(ds,viewer);
114:   }

116:   /* Print singular values */
117:   PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
118:   for (i=0;i<n;i++) {
119:     sigma = PetscRealPart(w[i]);
120:     PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma);
121:   }

123:   if (extrarow) {
124:     /* Check that extra row is correct */
125:     DSGetArrayReal(ds,DS_MAT_T,&T);
126:     DSGetArray(ds,DS_MAT_U,&U);
127:     DSGetArray(ds,DS_MAT_V,&V);
128:     d = 0.0;
129:     for (i=0;i<n;i++) d += T[i+ld]+U[n-1+i*ld];
130:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) {
131:       PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in A's extra row of %g\n",(double)PetscAbsScalar(d));
132:     }
133:     d = 0.0;
134:     for (i=0;i<n;i++) d += T[i+2*ld]-V[n-1+i*ld];
135:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) {
136:       PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in B's extra row of %g\n",(double)PetscAbsScalar(d));
137:     }
138:     DSRestoreArrayReal(ds,DS_MAT_T,&T);
139:     DSRestoreArray(ds,DS_MAT_U,&U);
140:     DSRestoreArray(ds,DS_MAT_V,&V);
141:   }

143:   /* Singular vectors */
144:   DSVectors(ds,DS_MAT_X,NULL,NULL);  /* all singular vectors */
145:   DSGetMat(ds,DS_MAT_X,&X);
146:   MatCreateVecs(X,NULL,&x0);
147:   MatGetColumnVector(X,x0,0);
148:   VecNorm(x0,NORM_2,&rnorm);
149:   MatDestroy(&X);
150:   VecDestroy(&x0);
151:   PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st X vector = %.3f\n",(double)rnorm);

153:   DSGetMat(ds,DS_MAT_U,&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) {
160:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st U vector has norm %g\n",(double)rnorm);
161:   }

163:   DSGetMat(ds,DS_MAT_V,&X);
164:   MatCreateVecs(X,NULL,&x0);
165:   MatGetColumnVector(X,x0,0);
166:   VecNorm(x0,NORM_2,&rnorm);
167:   MatDestroy(&X);
168:   VecDestroy(&x0);
169:   if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) {
170:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st V vector has norm %g\n",(double)rnorm);
171:   }

173:   PetscFree(w);
174:   DSDestroy(&ds);
175:   SlepcFinalize();
176:   return ierr;
177: }

179: /*TEST

181:    testset:
182:       requires: double
183:       test:
184:          suffix: 1
185:          args: -test_dsview
186:       test:
187:          suffix: 2
188:          args: -l 1 -k 4
189:       test:
190:          suffix: 2_extrarow
191:          filter: sed -e "s/extrarow//"
192:          args: -l 1 -k 4 -extrarow
193:          output_file: output/test22_2.out

195: TEST*/