Actual source code: test11.c

slepc-3.17.1 2022-04-11
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:   *errest = res/norm;
 41:   PetscFunctionReturn(0);
 42: }

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

 52:   SlepcInitialize(&argc,&argv,(char*)0,help);

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

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:         Generate the matrix
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

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

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

115: /*TEST

117:    test:
118:       suffix: 1

120: TEST*/