Actual source code: test10.c

slepc-3.21.1 2024-04-26
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[] = "Lanczos SVD. Also illustrates the use of SVDSetBV().\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];
 35:   PetscScalar    value[] = { 1, 2 };
 36:   PetscBool      flg,oneside=PETSC_FALSE;
 37:   const char     *prefix;
 38:   BV             U,V;
 39:   Vec            u,v;

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

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

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

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:          Create standalone BV objects to illustrate use of SVDSetBV()
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   PetscCall(BVCreate(PETSC_COMM_WORLD,&U));
 71:   PetscCall(PetscObjectSetName((PetscObject)U,"U"));
 72:   PetscCall(BVSetSizesFromVec(U,u,k));
 73:   PetscCall(BVSetFromOptions(U));
 74:   PetscCall(BVCreate(PETSC_COMM_WORLD,&V));
 75:   PetscCall(PetscObjectSetName((PetscObject)V,"V"));
 76:   PetscCall(BVSetSizesFromVec(V,v,k));
 77:   PetscCall(BVSetFromOptions(V));

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:         Compute singular values
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 84:   PetscCall(SVDSetBV(svd,V,U));
 85:   PetscCall(SVDSetOptionsPrefix(svd,"check_"));
 86:   PetscCall(SVDAppendOptionsPrefix(svd,"myprefix_"));
 87:   PetscCall(SVDGetOptionsPrefix(svd,&prefix));
 88:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"SVD prefix is currently: %s\n\n",prefix));
 89:   PetscCall(PetscObjectSetName((PetscObject)svd,"SVD_solver"));

 91:   PetscCall(SVDSetOperators(svd,A,NULL));
 92:   PetscCall(SVDSetType(svd,SVDLANCZOS));
 93:   PetscCall(SVDSetFromOptions(svd));

 95:   PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDLANCZOS,&flg));
 96:   if (flg) {
 97:     PetscCall(SVDLanczosGetOneSide(svd,&oneside));
 98:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Running Lanczos %s\n\n",oneside?"(onesided)":""));
 99:   }
100:   PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg));
101:   if (flg) {
102:     PetscCall(SVDTRLanczosGetOneSide(svd,&oneside));
103:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Running thick-restart Lanczos %s\n\n",oneside?"(onesided)":""));
104:   }

106:   PetscCall(SVDSolve(svd));

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:                     Display solution and clean up
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
112:   PetscCall(BVDestroy(&U));
113:   PetscCall(BVDestroy(&V));
114:   PetscCall(VecDestroy(&u));
115:   PetscCall(VecDestroy(&v));
116:   PetscCall(SVDDestroy(&svd));
117:   PetscCall(MatDestroy(&A));
118:   PetscCall(SlepcFinalize());
119:   return 0;
120: }

122: /*TEST

124:    testset:
125:       args: -check_myprefix_svd_nsv 3
126:       requires: double
127:       test:
128:          suffix: 1
129:          args: -check_myprefix_svd_view_vectors ::ascii_info
130:       test:
131:          suffix: 2
132:          args: -check_myprefix_svd_type trlanczos -check_myprefix_svd_monitor -check_myprefix_svd_view_values ::ascii_matlab
133:          filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
134:       test:
135:          suffix: 3
136:          args: -m 22 -n 20

138: TEST*/