Actual source code: test7.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  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.\n\n";

 24:  #include slepcds.h

 28: int main( int argc, char **argv )
 29: {
 31:   DS             ds;
 32:   PetscReal      sigma;
 33:   PetscScalar    *A,*w;
 34:   PetscInt       i,j,k,n=15,m=10,ld;
 35:   PetscViewer    viewer;
 36:   PetscBool      verbose;

 38:   SlepcInitialize(&argc,&argv,(char*)0,help);
 39:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 40:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 41:   k = PetscMin(n,m);
 42:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %Dx%D.\n",n,m);
 43:   PetscOptionsHasName(PETSC_NULL,"-verbose",&verbose);

 45:   /* Create DS object */
 46:   DSCreate(PETSC_COMM_WORLD,&ds);
 47:   DSSetType(ds,DSSVD);
 48:   DSSetFromOptions(ds);
 49:   ld = n+2;  /* test leading dimension larger than n */
 50:   DSAllocate(ds,ld);
 51:   DSSetDimensions(ds,n,m,0,0);

 53:   /* Set up viewer */
 54:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 55:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 56:   DSView(ds,viewer);
 57:   PetscViewerPopFormat(viewer);
 58:   if (verbose) {
 59:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 60:   }

 62:   /* Fill with a rectangular Toeplitz matrix */
 63:   DSGetArray(ds,DS_MAT_A,&A);
 64:   for (i=0;i<k;i++) A[i+i*ld]=1.0;
 65:   for (j=1;j<3;j++) {
 66:     for (i=0;i<n-j;i++) { if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); }
 67:   }
 68:   for (j=1;j<n/2;j++) {
 69:     for (i=0;i<n-j;i++) { if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; }
 70:   }
 71:   DSRestoreArray(ds,DS_MAT_A,&A);
 72:   DSSetState(ds,DS_STATE_RAW);
 73:   if (verbose) {
 74:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 75:     DSView(ds,viewer);
 76:   }

 78:   /* Solve */
 79:   PetscMalloc(k*sizeof(PetscScalar),&w);
 80:   DSSetEigenvalueComparison(ds,SlepcCompareLargestReal,PETSC_NULL);
 81:   DSSolve(ds,w,PETSC_NULL);
 82:   DSSort(ds,w,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 83:   if (verbose) {
 84:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 85:     DSView(ds,viewer);
 86:   }
 87: 
 88:   /* Print singular values */
 89:   PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n",n);
 90:   for (i=0;i<k;i++) {
 91:     sigma = PetscRealPart(w[i]);
 92:     PetscViewerASCIIPrintf(viewer,"  %.5F\n",sigma);
 93:   }
 94:   PetscFree(w);
 95:   DSDestroy(&ds);
 96:   SlepcFinalize();
 97:   return 0;
 98: }