Actual source code: test13.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 interface to external package PRIMME.\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = matrix rows.\n"
 14:   "  -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";

 16: #include <slepcsvd.h>

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

 21:               |  1  2                     |
 22:               |     1  2                  |
 23:               |        1  2               |
 24:           A = |          .  .             |
 25:               |             .  .          |
 26:               |                1  2       |
 27:               |                   1  2    |
 28:  */

 30: int main(int argc,char **argv)
 31: {
 32:   Mat             A;
 33:   SVD             svd;
 34:   PetscInt        m=20,n,Istart,Iend,i,k=6,col[2],bs;
 35:   PetscScalar     value[] = { 1, 2 };
 36:   PetscBool       flg;
 37:   SVDPRIMMEMethod meth;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 41:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 42:   if (!flg) n=m+2;
 43:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 44:   PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n);

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

 50:   MatCreate(PETSC_COMM_WORLD,&A);
 51:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 52:   MatSetFromOptions(A);
 53:   MatSetUp(A);
 54:   MatGetOwnershipRange(A,&Istart,&Iend);
 55:   for (i=Istart;i<Iend;i++) {
 56:     col[0]=i; col[1]=i+1;
 57:     if (i<n-1) MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 58:     else if (i==n-1) MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
 59:   }
 60:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:         Compute singular values
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   SVDCreate(PETSC_COMM_WORLD,&svd);
 68:   SVDSetOperators(svd,A,NULL);
 69:   SVDSetType(svd,SVDPRIMME);
 70:   SVDPRIMMESetBlockSize(svd,4);
 71:   SVDPRIMMESetMethod(svd,SVD_PRIMME_HYBRID);
 72:   SVDSetFromOptions(svd);

 74:   SVDSolve(svd);
 75:   SVDPRIMMEGetBlockSize(svd,&bs);
 76:   SVDPRIMMEGetMethod(svd,&meth);
 77:   PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,SVDPRIMMEMethods[meth]);

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:                     Display solution and clean up
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
 83:   SVDDestroy(&svd);
 84:   MatDestroy(&A);
 85:   SlepcFinalize();
 86:   return 0;
 87: }

 89: /*TEST

 91:    build:
 92:       requires: primme

 94:    test:
 95:       suffix: 1
 96:       args: -svd_nsv 4
 97:       requires: primme
 98:       filter: sed -e "s/2.99255/2.99254/"

100: TEST*/