Actual source code: ex15.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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;

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

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

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

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

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

 54:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 55:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 56:   MatCreateVecs(A,&v,&u);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:           Create the singular value solver and set various options
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   /*
 63:      Create singular value solver context
 64:   */
 65:   SVDCreate(PETSC_COMM_WORLD,&svd);

 67:   /*
 68:      Set operators and problem type
 69:   */
 70:   SVDSetOperators(svd,A,NULL);
 71:   SVDSetProblemType(svd,SVD_STANDARD);

 73:   /*
 74:      Use thick-restart Lanczos as default solver
 75:   */
 76:   SVDSetType(svd,SVDTRLANCZOS);

 78:   /*
 79:      Set solver parameters at runtime
 80:   */
 81:   SVDSetFromOptions(svd);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                       Solve the singular value system
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 87:   SVDSolve(svd);
 88:   SVDGetIterationNumber(svd,&its);
 89:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);

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

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                     Display solution and clean up
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   /*
106:      Get number of converged singular triplets
107:   */
108:   SVDGetConverged(svd,&nconv);
109:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %D\n\n",nconv);

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

124:       /*
125:          Compute the error associated to each singular triplet
126:       */
127:       SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error);

129:       PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",(double)sigma);
130:       PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
131:     }
132:     PetscPrintf(PETSC_COMM_WORLD,"\n");
133:   }

135:   /*
136:      Free work space
137:   */
138:   SVDDestroy(&svd);
139:   MatDestroy(&A);
140:   VecDestroy(&u);
141:   VecDestroy(&v);
142:   SlepcFinalize();
143:   return ierr;
144: }

146: /*TEST

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

159: TEST*/