Actual source code: test27.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 DSHSVD with dense storage.\n\n";

 13: #include <slepcds.h>

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

 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 HSVD - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
 32:   PetscCheck(p>=0 && p<=m,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value %" PetscInt_FMT ", must be 0<=p<=%" PetscInt_FMT,p,m);
 33:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));

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

 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:   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 66:   /* Fill signature matrix */
 67:   PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
 68:   for (i=0;i<n;i++) D[i] = (i<p)? -1.0: 1.0;
 69:   PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
 70:   PetscCall(DSSetState(ds,DS_STATE_RAW));
 71:   if (verbose) {
 72:     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
 73:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
 74:   }
 75:   PetscCall(DSView(ds,viewer));

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

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

 98:   /* Singular vectors */
 99:   PetscCall(DSVectors(ds,DS_MAT_U,NULL,NULL));  /* all singular vectors */
100:   j = 0;
101:   rnorm = 0.0;
102:   PetscCall(DSGetArray(ds,DS_MAT_U,&U));
103:   for (i=0;i<n;i++) {
104:     aux = PetscAbsScalar(U[i+j*ld]);
105:     rnorm += aux*aux;
106:   }
107:   PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
108:   rnorm = PetscSqrtReal(rnorm);
109:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st U vector = %.3f\n",(double)rnorm));

111:   PetscCall(PetscFree(w));
112:   PetscCall(DSDestroy(&ds));
113:   PetscCall(SlepcFinalize());
114:   return 0;
115: }

117: /*TEST

119:    test:
120:       suffix: 1
121:       requires: !single

123: TEST*/