Actual source code: test11.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[] = "Tests a user-defined convergence test (based on ex8.c).\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n\n";

 15: #include <slepcsvd.h>

 17: /*
 18:    This example computes the singular values of an nxn Grcar matrix,
 19:    which is a nonsymmetric Toeplitz matrix:

 21:               |  1  1  1  1               |
 22:               | -1  1  1  1  1            |
 23:               |    -1  1  1  1  1         |
 24:               |       .  .  .  .  .       |
 25:           A = |          .  .  .  .  .    |
 26:               |            -1  1  1  1  1 |
 27:               |               -1  1  1  1 |
 28:               |                  -1  1  1 |
 29:               |                     -1  1 |

 31:  */

 33: /*
 34:   MyConvergedRel - Convergence test relative to the norm of A (given in ctx).
 35: */
 36: PetscErrorCode MyConvergedRel(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
 37: {
 38:   PetscReal norm = *(PetscReal*)ctx;

 40:   PetscFunctionBegin;
 41:   *errest = res/norm;
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: int main(int argc,char **argv)
 46: {
 47:   Mat            A;               /* Grcar matrix */
 48:   SVD            svd;             /* singular value solver context */
 49:   PetscInt       N=30,Istart,Iend,i,col[5],nconv1,nconv2;
 50:   PetscScalar    value[] = { -1, 1, 1, 1, 1 };
 51:   PetscReal      sigma_1,sigma_n;

 53:   PetscFunctionBeginUser;
 54:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 56:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
 57:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N));

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:         Generate the matrix
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 64:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 65:   PetscCall(MatSetFromOptions(A));
 66:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 67:   for (i=Istart;i<Iend;i++) {
 68:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 69:     if (i==0) PetscCall(MatSetValues(A,1,&i,PetscMin(4,N-i),col+1,value+1,INSERT_VALUES));
 70:     else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
 71:   }
 72:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 73:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:              Create the SVD solver and set the solution method
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 80:   PetscCall(SVDSetOperators(svd,A,NULL));
 81:   PetscCall(SVDSetType(svd,SVDTRLANCZOS));
 82:   PetscCall(SVDSetFromOptions(svd));

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                       Solve the singular value problem
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
 89:   PetscCall(SVDSolve(svd));
 90:   PetscCall(SVDGetConverged(svd,&nconv1));
 91:   if (nconv1 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL));
 92:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n"));

 94:   /* compute smallest singular value relative to the matrix norm */
 95:   PetscCall(SVDSetConvergenceTestFunction(svd,MyConvergedRel,&sigma_1,NULL));
 96:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
 97:   PetscCall(SVDSolve(svd));
 98:   PetscCall(SVDGetConverged(svd,&nconv2));
 99:   if (nconv2 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL));
100:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n"));

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:                     Display solution and clean up
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105:   if (nconv1 > 0 && nconv2 > 0) {
106:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n));
107:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n)));
108:   }

110:   PetscCall(SVDDestroy(&svd));
111:   PetscCall(MatDestroy(&A));
112:   PetscCall(SlepcFinalize());
113:   return 0;
114: }

116: /*TEST

118:    test:
119:       suffix: 1

121: TEST*/