Actual source code: test8.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test DSSVD with compact storage.\n\n";

 24: #include <slepcds.h>

 28: int main(int argc,char **argv)
 29: {
 31:   DS             ds;
 32:   SlepcSC        sc;
 33:   PetscReal      *T,sigma;
 34:   PetscScalar    *w;
 35:   PetscInt       i,n=10,m,l=2,k=5,ld;
 36:   PetscViewer    viewer;
 37:   PetscBool      verbose;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 41:   m = n;
 42:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD with compact storage - dimension %Dx%D.\n",n,m);
 43:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 44:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 45:   if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
 46:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

 48:   /* Create DS object */
 49:   DSCreate(PETSC_COMM_WORLD,&ds);
 50:   DSSetType(ds,DSSVD);
 51:   DSSetFromOptions(ds);
 52:   ld = n+2;  /* test leading dimension larger than n */
 53:   DSAllocate(ds,ld);
 54:   DSSetDimensions(ds,n,m,l,k);
 55:   DSSetCompact(ds,PETSC_TRUE);

 57:   /* Set up viewer */
 58:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 59:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 60:   DSView(ds,viewer);
 61:   PetscViewerPopFormat(viewer);
 62:   if (verbose) {
 63:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 64:   }

 66:   /* Fill upper arrow-tridiagonal matrix */
 67:   DSGetArrayReal(ds,DS_MAT_T,&T);
 68:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
 69:   for (i=l;i<n-1;i++) T[i+ld] = 1.0;
 70:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 71:   if (l==0 && k==0) {
 72:     DSSetState(ds,DS_STATE_INTERMEDIATE);
 73:   } else {
 74:     DSSetState(ds,DS_STATE_RAW);
 75:   }
 76:   if (verbose) {
 77:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 78:     DSView(ds,viewer);
 79:   }

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

 95:   /* Print singular values */
 96:   PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n");
 97:   for (i=0;i<n;i++) {
 98:     sigma = PetscRealPart(w[i]);
 99:     PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma);
100:   }
101:   PetscFree(w);
102:   DSDestroy(&ds);
103:   SlepcFinalize();
104:   return ierr;
105: }