Actual source code: test3.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 SVD with user-provided initial vectors.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = row dimension.\n"
 14:   "  -m <m>, where <m> = column dimension.\n\n";

 16: #include <slepcsvd.h>

 18: /*
 19:    This example computes the singular values of a rectangular nxm Grcar matrix:

 21:               |  1  1  1  1               |
 22:               | -1  1  1  1  1            |
 23:               |    -1  1  1  1  1         |
 24:           A = |       .  .  .  .  .       |
 25:               |          .  .  .  .  .    |
 26:               |            -1  1  1  1  1 |
 27:               |               -1  1  1  1 |

 29:  */

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            A;               /* Grcar matrix */
 34:   SVD            svd;             /* singular value solver context */
 35:   Vec            v0,w0;           /* initial vectors */
 36:   PetscInt       N=35,M=30,Istart,Iend,i,col[5];
 37:   PetscScalar    value[] = { -1, 1, 1, 1, 1 };

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 41:   PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL);
 43:   PetscPrintf(PETSC_COMM_WORLD,"\nSVD of a rectangular Grcar matrix, %Dx%D\n\n",N,M);

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 46:         Generate the matrix
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 49:   MatCreate(PETSC_COMM_WORLD,&A);
 50:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M);
 51:   MatSetFromOptions(A);
 52:   MatSetUp(A);

 54:   MatGetOwnershipRange(A,&Istart,&Iend);
 55:   for (i=Istart;i<Iend;i++) {
 56:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 57:     if (i==0) {
 58:       MatSetValues(A,1,&i,PetscMin(4,M-i+1),col+1,value+1,INSERT_VALUES);
 59:     } else {
 60:       MatSetValues(A,1,&i,PetscMin(5,M-i+1),col,value,INSERT_VALUES);
 61:     }
 62:   }
 63:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 64:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:              Create the SVD context and solve the problem
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   SVDCreate(PETSC_COMM_WORLD,&svd);
 71:   SVDSetOperator(svd,A);
 72:   SVDSetFromOptions(svd);

 74:   /*
 75:      Set the initial vectors. This is optional, if not done the initial
 76:      vectors are set to random values
 77:   */
 78:   MatCreateVecs(A,&v0,&w0);
 79:   VecSet(v0,1.0);
 80:   VecSet(w0,1.0);
 81:   SVDSetInitialSpaces(svd,1,&v0,1,&w0);

 83:   /*
 84:      Compute solution
 85:   */
 86:   SVDSolve(svd);
 87:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);

 89:   /*
 90:      Free work space
 91:   */
 92:   VecDestroy(&v0);
 93:   VecDestroy(&w0);
 94:   SVDDestroy(&svd);
 95:   MatDestroy(&A);
 96:   SlepcFinalize();
 97:   return ierr;
 98: }