Actual source code: test8.c

slepc-3.18.2 2023-01-26
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 with compact storage.\n\n";

 13: #include <slepcds.h>

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

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   m = n;
 29:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD with compact storage - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m);
 30:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

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

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

 54:   /* Fill upper arrow-bidiagonal matrix */
 55:   DSGetArrayReal(ds,DS_MAT_T,&T);
 56:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
 57:   for (i=l;i<n-1;i++) T[i+ld] = 1.0;
 58:   if (extrarow) T[n-1+ld] = 1.0;
 59:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 60:   if (l==0 && k==0) DSSetState(ds,DS_STATE_INTERMEDIATE);
 61:   else DSSetState(ds,DS_STATE_RAW);
 62:   PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 63:   DSView(ds,viewer);

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

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

 87:   if (extrarow) {
 88:     /* Check that extra row is correct */
 89:     DSGetArrayReal(ds,DS_MAT_T,&T);
 90:     DSGetArray(ds,DS_MAT_U,&U);
 91:     d = 0.0;
 92:     for (i=l;i<n;i++) d += T[i+ld]-U[n-1+i*ld];
 93:     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));
 94:     DSRestoreArrayReal(ds,DS_MAT_T,&T);
 95:     DSRestoreArray(ds,DS_MAT_U,&U);
 96:   }
 97:   PetscFree(w);
 98:   DSDestroy(&ds);
 99:   SlepcFinalize();
100:   return 0;
101: }

103: /*TEST

105:    test:
106:       suffix: 1
107:       requires: !single

109:    test:
110:       args: -l 0 -k 0
111:       suffix: 2
112:       requires: !single

114:    test:
115:       args: -extrarow
116:       suffix: 3
117:       requires: !single

119: TEST*/