Actual source code: test7.c

slepc-3.17.2 2022-08-09
Report Typos and Errors
  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:   SlepcInitialize(&argc,&argv,(char*)0,help);
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 27:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 28:   k = PetscMin(n,m);
 29:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m);
 30:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 31:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

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

 45:   /* Set up viewer */
 46:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 47:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 48:   DSView(ds,viewer);
 49:   PetscViewerPopFormat(viewer);

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

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

 84:   /* Print singular values */
 85:   PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
 86:   for (i=0;i<k;i++) {
 87:     sigma = PetscRealPart(w[i]);
 88:     PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma);
 89:   }

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

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

115:   PetscFree(w);
116:   DSDestroy(&ds);
117:   SlepcFinalize();
118:   return 0;
119: }

121: /*TEST

123:    test:
124:       suffix: 1
125:       requires: !single

127:    test:
128:       suffix: 2
129:       args: -extrarow
130:       requires: !single

132: TEST*/