Actual source code: test7.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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[] = "SVD via the cyclic matrix with a user-provided EPS.\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:   EPS                  eps;
 35:   ST                   st;
 36:   KSP                  ksp;
 37:   PC                   pc;
 38:   PetscInt             m=20,n,Istart,Iend,i,col[2];
 39:   PetscScalar          value[] = { 1, 2 };
 40:   PetscBool            flg,expmat;
 41:   PetscErrorCode       ierr;

 43:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 45:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 46:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 47:   if (!flg) n=m+2;
 48:   PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:         Generate the matrix
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:         Create a standalone EPS with appropriate settings
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   EPSCreate(PETSC_COMM_WORLD,&eps);
 75:   EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
 76:   EPSSetTarget(eps,1.0);
 77:   EPSGetST(eps,&st);
 78:   STSetType(st,STSINVERT);
 79:   STSetShift(st,1.01);
 80:   STGetKSP(st,&ksp);
 81:   KSPSetType(ksp,KSPPREONLY);
 82:   KSPGetPC(ksp,&pc);
 83:   PCSetType(pc,PCLU);
 84:   EPSSetFromOptions(eps);

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:         Compute singular values
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 90:   SVDCreate(PETSC_COMM_WORLD,&svd);
 91:   SVDSetOperator(svd,A);
 92:   SVDSetType(svd,SVDCYCLIC);
 93:   SVDCyclicSetEPS(svd,eps);
 94:   SVDCyclicSetExplicitMatrix(svd,PETSC_TRUE);
 95:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
 96:   SVDSetFromOptions(svd);
 97:   PetscObjectTypeCompare((PetscObject)svd,SVDCYCLIC,&flg);
 98:   if (flg) {
 99:     SVDCyclicGetExplicitMatrix(svd,&expmat);
100:     if (expmat) {
101:       PetscPrintf(PETSC_COMM_WORLD," Using explicit matrix with cyclic solver\n",m,n);
102:     }
103:   }
104:   SVDSolve(svd);

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:                     Display solution and clean up
108:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109:   SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
110:   SVDDestroy(&svd);
111:   EPSDestroy(&eps);
112:   MatDestroy(&A);
113:   SlepcFinalize();
114:   return ierr;
115: }

117: /*TEST

119:    test:
120:       suffix: 1
121:       args: -info_exclude svd -log_exclude svd
122:       requires: !single

124:    test:
125:       suffix: 2_cuda
126:       args: -info_exclude svd -log_exclude svd -mat_type aijcusparse
127:       requires: cuda !single
128:       output_file: output/test7_1.out

130: TEST*/