Actual source code: ex15.c

slepc-3.20.1 2023-11-27
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[] = "Singular value decomposition of the Lauchli matrix.\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n"
 14:   "  -mu <mu>, where <mu> = subdiagonal value.\n\n";

 16: #include <slepcsvd.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A;               /* operator matrix */
 21:   Vec            u,v;             /* left and right singular vectors */
 22:   SVD            svd;             /* singular value problem solver context */
 23:   SVDType        type;
 24:   PetscReal      error,tol,sigma,mu=PETSC_SQRT_MACHINE_EPSILON;
 25:   PetscInt       n=100,i,j,Istart,Iend,nsv,maxit,its,nconv;

 27:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 31:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
 32:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%" PetscInt_FMT " x %" PetscInt_FMT ") mu=%g\n\n",n+1,n,(double)mu));

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35:                           Build the Lauchli matrix
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 39:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n));
 40:   PetscCall(MatSetFromOptions(A));
 41:   PetscCall(MatSetUp(A));

 43:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 44:   for (i=Istart;i<Iend;i++) {
 45:     if (i == 0) {
 46:       for (j=0;j<n;j++) PetscCall(MatSetValue(A,0,j,1.0,INSERT_VALUES));
 47:     } else PetscCall(MatSetValue(A,i,i-1,mu,INSERT_VALUES));
 48:   }

 50:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 52:   PetscCall(MatCreateVecs(A,&v,&u));

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:           Create the singular value solver and set various options
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   /*
 59:      Create singular value solver context
 60:   */
 61:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));

 63:   /*
 64:      Set operators and problem type
 65:   */
 66:   PetscCall(SVDSetOperators(svd,A,NULL));
 67:   PetscCall(SVDSetProblemType(svd,SVD_STANDARD));

 69:   /*
 70:      Use thick-restart Lanczos as default solver
 71:   */
 72:   PetscCall(SVDSetType(svd,SVDTRLANCZOS));

 74:   /*
 75:      Set solver parameters at runtime
 76:   */
 77:   PetscCall(SVDSetFromOptions(svd));

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:                       Solve the singular value system
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   PetscCall(SVDSolve(svd));
 84:   PetscCall(SVDGetIterationNumber(svd,&its));
 85:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));

 87:   /*
 88:      Optional: Get some information from the solver and display it
 89:   */
 90:   PetscCall(SVDGetType(svd,&type));
 91:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
 92:   PetscCall(SVDGetDimensions(svd,&nsv,NULL,NULL));
 93:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %" PetscInt_FMT "\n",nsv));
 94:   PetscCall(SVDGetTolerances(svd,&tol,&maxit));
 95:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:                     Display solution and clean up
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   /*
102:      Get number of converged singular triplets
103:   */
104:   PetscCall(SVDGetConverged(svd,&nconv));
105:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %" PetscInt_FMT "\n\n",nconv));

107:   if (nconv>0) {
108:     /*
109:        Display singular values and relative errors
110:     */
111:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,
112:          "          sigma           relative error\n"
113:          "  --------------------- ------------------\n"));
114:     for (i=0;i<nconv;i++) {
115:       /*
116:          Get converged singular triplets: i-th singular value is stored in sigma
117:       */
118:       PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));

120:       /*
121:          Compute the error associated to each singular triplet
122:       */
123:       PetscCall(SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error));

125:       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",(double)sigma));
126:       PetscCall(PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error));
127:     }
128:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
129:   }

131:   /*
132:      Free work space
133:   */
134:   PetscCall(SVDDestroy(&svd));
135:   PetscCall(MatDestroy(&A));
136:   PetscCall(VecDestroy(&u));
137:   PetscCall(VecDestroy(&v));
138:   PetscCall(SlepcFinalize());
139:   return 0;
140: }

142: /*TEST

144:    testset:
145:       filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
146:       requires: double
147:       test:
148:          suffix: 1
149:       test:
150:          suffix: 1_scalapack
151:          nsize: {{1 2}}
152:          args: -svd_type scalapack
153:          requires: scalapack

155: TEST*/