Actual source code: test11.c

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

 41:   *errest = res/norm;
 42:   return(0);
 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;

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

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

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

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

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:              Create the SVD solver and set the solution method
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   SVDCreate(PETSC_COMM_WORLD,&svd);
 84:   SVDSetOperator(svd,A);
 85:   SVDSetType(svd,SVDTRLANCZOS);
 86:   SVDSetFromOptions(svd);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                       Solve the singular value problem
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
 93:   SVDSolve(svd);
 94:   SVDGetConverged(svd,&nconv1);
 95:   if (nconv1 > 0) {
 96:     SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL);
 97:   } else {
 98:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
 99:   }

101:   /* compute smallest singular value relative to the matrix norm */
102:   SVDSetConvergenceTestFunction(svd,MyConvergedRel,&sigma_1,NULL);
103:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
104:   SVDSolve(svd);
105:   SVDGetConverged(svd,&nconv2);
106:   if (nconv2 > 0) {
107:     SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL);
108:   } else {
109:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
110:   }

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:                     Display solution and clean up
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115:   if (nconv1 > 0 && nconv2 > 0) {
116:     PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n);
117:     PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n));
118:   }

120:   SVDDestroy(&svd);
121:   MatDestroy(&A);
122:   SlepcFinalize();
123:   return ierr;
124: }

126: /*TEST

128:    test:
129:       suffix: 1

131: TEST*/