Actual source code: test3.c

slepc-3.7.0 2016-05-16
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 SVD with user-provided initial vectors.\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = row dimension.\n"
 25:   "  -m <m>, where <m> = column dimension.\n\n";

 27: #include <slepcsvd.h>

 29: /*
 30:    This example computes the singular values of a rectangular nxm Grcar matrix:

 32:               |  1  1  1  1               |
 33:               | -1  1  1  1  1            |
 34:               |    -1  1  1  1  1         |
 35:           A = |       .  .  .  .  .       |
 36:               |          .  .  .  .  .    |
 37:               |            -1  1  1  1  1 |
 38:               |               -1  1  1  1 |

 40:  */

 44: int main(int argc,char **argv)
 45: {
 46:   Mat            A;               /* Grcar matrix */
 47:   SVD            svd;             /* singular value solver context */
 48:   Vec            v0,w0;           /* initial vectors */
 49:   PetscInt       N=35,M=30,Istart,Iend,i,col[5];
 50:   PetscScalar    value[] = { -1, 1, 1, 1, 1 };

 53:   SlepcInitialize(&argc,&argv,(char*)0,help);
 54:   PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
 55:   PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL);
 56:   PetscPrintf(PETSC_COMM_WORLD,"\nSVD of a rectangular Grcar matrix, %Dx%D\n\n",N,M);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:         Generate the matrix
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   MatCreate(PETSC_COMM_WORLD,&A);
 63:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M);
 64:   MatSetFromOptions(A);
 65:   MatSetUp(A);

 67:   MatGetOwnershipRange(A,&Istart,&Iend);
 68:   for (i=Istart;i<Iend;i++) {
 69:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 70:     if (i==0) {
 71:       MatSetValues(A,1,&i,PetscMin(4,M-i+1),col+1,value+1,INSERT_VALUES);
 72:     } else {
 73:       MatSetValues(A,1,&i,PetscMin(5,M-i+1),col,value,INSERT_VALUES);
 74:     }
 75:   }
 76:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:              Create the SVD context and solve the problem
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   SVDCreate(PETSC_COMM_WORLD,&svd);
 84:   SVDSetOperator(svd,A);
 85:   SVDSetFromOptions(svd);

 87:   /*
 88:      Set the initial vectors. This is optional, if not done the initial
 89:      vectors are set to random values
 90:   */
 91:   MatCreateVecs(A,&v0,&w0);
 92:   VecSet(v0,1.0);
 93:   VecSet(w0,1.0);
 94:   SVDSetInitialSpace(svd,1,&v0);
 95:   SVDSetInitialSpaceLeft(svd,1,&w0);

 97:   /*
 98:      Compute solution
 99:   */
100:   SVDSolve(svd);
101:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);

103:   /*
104:      Free work space
105:   */
106:   VecDestroy(&v0);
107:   VecDestroy(&w0);
108:   SVDDestroy(&svd);
109:   MatDestroy(&A);
110:   SlepcFinalize();
111:   return ierr;
112: }