Actual source code: test7.c

slepc-3.22.2 2024-12-02
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[] = "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;

 42:   PetscFunctionBeginUser;
 43:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

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

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

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

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:         Create a standalone EPS with appropriate settings
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 71:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
 72:   PetscCall(EPSSetTarget(eps,1.0));
 73:   PetscCall(EPSGetST(eps,&st));
 74:   PetscCall(STSetType(st,STSINVERT));
 75:   PetscCall(STSetShift(st,1.01));
 76:   PetscCall(STGetKSP(st,&ksp));
 77:   PetscCall(KSPSetType(ksp,KSPPREONLY));
 78:   PetscCall(KSPGetPC(ksp,&pc));
 79:   PetscCall(PCSetType(pc,PCLU));
 80:   PetscCall(EPSSetFromOptions(eps));

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:         Compute singular values
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 87:   PetscCall(SVDSetOperators(svd,A,NULL));
 88:   PetscCall(SVDSetType(svd,SVDCYCLIC));
 89:   PetscCall(SVDCyclicSetEPS(svd,eps));
 90:   PetscCall(SVDCyclicSetExplicitMatrix(svd,PETSC_TRUE));
 91:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
 92:   PetscCall(SVDSetFromOptions(svd));
 93:   PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDCYCLIC,&flg));
 94:   if (flg) {
 95:     PetscCall(SVDCyclicGetExplicitMatrix(svd,&expmat));
 96:     if (expmat) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using explicit matrix with cyclic solver\n"));
 97:   }
 98:   PetscCall(SVDSolve(svd));

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:                     Display solution and clean up
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   PetscCall(SVDErrorView(svd,SVD_ERROR_ABSOLUTE,PETSC_VIEWER_STDOUT_WORLD));
104:   PetscCall(SVDDestroy(&svd));
105:   PetscCall(EPSDestroy(&eps));
106:   PetscCall(MatDestroy(&A));
107:   PetscCall(SlepcFinalize());
108:   return 0;
109: }

111: /*TEST

113:    test:
114:       suffix: 1
115:       args: -log_exclude svd

117:    testset:
118:       args: -log_exclude svd
119:       output_file: output/test7_1.out
120:       test:
121:          suffix: 2_cuda
122:          args: -mat_type aijcusparse
123:          requires: cuda
124:       test:
125:          suffix: 2_hip
126:          args: -mat_type aijhipsparse
127:          requires: hip

129: TEST*/