Actual source code: test7.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 DSSVD.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   SlepcSC        sc;
 19:   PetscReal      sigma,rnorm,aux;
 20:   PetscScalar    *A,*U,*w,d;
 21:   PetscInt       i,j,k,n=15,m=10,m1,ld;
 22:   PetscViewer    viewer;
 23:   PetscBool      verbose,extrarow;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
 29:   k = PetscMin(n,m);
 30:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m));
 31:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 32:   PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));

 34:   /* Create DS object */
 35:   PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
 36:   PetscCall(DSSetType(ds,DSSVD));
 37:   PetscCall(DSSetFromOptions(ds));
 38:   ld = PetscMax(n,m)+2;  /* test leading dimension larger than n */
 39:   PetscCall(DSAllocate(ds,ld));
 40:   PetscCall(DSSetDimensions(ds,n,0,0));
 41:   PetscCall(DSSVDSetDimensions(ds,m));
 42:   PetscCall(DSSVDGetDimensions(ds,&m1));
 43:   PetscCheck(m1==m,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Inconsistent dimension value");
 44:   PetscCall(DSSetExtraRow(ds,extrarow));

 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:   /* Fill with a rectangular Toeplitz matrix */
 53:   PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 54:   for (i=0;i<k;i++) A[i+i*ld]=1.0;
 55:   for (j=1;j<3;j++) {
 56:     for (i=0;i<m-j;i++) {
 57:       if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1);
 58:     }
 59:   }
 60:   for (j=1;j<n/2;j++) {
 61:     for (i=0;i<n-j;i++) {
 62:       if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0;
 63:     }
 64:   }
 65:   if (extrarow) { A[n-2+m*ld]=1.0; A[n-1+m*ld]=1.0; }  /* really an extra column */
 66:   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 67:   PetscCall(DSSetState(ds,DS_STATE_RAW));
 68:   if (verbose) {
 69:     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
 70:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
 71:   }
 72:   PetscCall(DSView(ds,viewer));

 74:   /* Solve */
 75:   PetscCall(PetscMalloc1(k,&w));
 76:   PetscCall(DSGetSlepcSC(ds,&sc));
 77:   sc->comparison    = SlepcCompareLargestReal;
 78:   sc->comparisonctx = NULL;
 79:   sc->map           = NULL;
 80:   sc->mapobj        = NULL;
 81:   PetscCall(DSSolve(ds,w,NULL));
 82:   PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
 83:   if (extrarow) PetscCall(DSUpdateExtraRow(ds));
 84:   if (verbose) {
 85:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
 86:     PetscCall(DSView(ds,viewer));
 87:   }

 89:   /* Print singular values */
 90:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
 91:   for (i=0;i<k;i++) {
 92:     sigma = PetscRealPart(w[i]);
 93:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma));
 94:   }

 96:   if (extrarow) {
 97:     /* Check that extra column is correct */
 98:     PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 99:     PetscCall(DSGetArray(ds,DS_MAT_U,&U));
100:     d = 0.0;
101:     for (i=0;i<n;i++) d += A[i+m*ld]-U[n-2+i*ld]-U[n-1+i*ld];
102:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d)));
103:     PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
104:     PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
105:   }

107:   /* Singular vectors */
108:   PetscCall(DSVectors(ds,DS_MAT_U,NULL,NULL));  /* all singular vectors */
109:   j = 0;
110:   rnorm = 0.0;
111:   PetscCall(DSGetArray(ds,DS_MAT_U,&U));
112:   for (i=0;i<n;i++) {
113:     aux = PetscAbsScalar(U[i+j*ld]);
114:     rnorm += aux*aux;
115:   }
116:   PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
117:   rnorm = PetscSqrtReal(rnorm);
118:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st U vector = %.3f\n",(double)rnorm));

120:   PetscCall(PetscFree(w));
121:   PetscCall(DSDestroy(&ds));
122:   PetscCall(SlepcFinalize());
123:   return 0;
124: }

126: /*TEST

128:    test:
129:       args: -ds_method {{0 1}}
130:       suffix: 1
131:       filter: grep -v "solving the problem"
132:       requires: !single

134:    test:
135:       suffix: 2
136:       args: -extrarow -ds_method {{0 1}}
137:       filter: grep -v "solving the problem"
138:       requires: !single

140: TEST*/