Actual source code: test8.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 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;

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

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

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

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

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

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

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

102: /*TEST

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

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

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

118: TEST*/