Actual source code: ex8.c

slepc-main 2024-11-15
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[] = "Estimates the 2-norm condition number of a matrix A, that is, the ratio of the largest to the smallest singular values of A. "
 12:   "The matrix is a Grcar matrix.\n\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = matrix dimension.\n\n";

 16: #include <slepcsvd.h>

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

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

 32:  */

 34: int main(int argc,char **argv)
 35: {
 36:   Mat            A;               /* Grcar matrix */
 37:   SVD            svd;             /* singular value solver context */
 38:   PetscInt       N=30,Istart,Iend,i,col[5],nconv1,nconv2;
 39:   PetscScalar    value[] = { -1, 1, 1, 1, 1 };
 40:   PetscReal      sigma_1,sigma_n;

 42:   PetscFunctionBeginUser;
 43:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

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

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:         Generate the matrix
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 53:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 54:   PetscCall(MatSetFromOptions(A));

 56:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 57:   for (i=Istart;i<Iend;i++) {
 58:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 59:     if (i==0) PetscCall(MatSetValues(A,1,&i,PetscMin(4,N-i),col+1,value+1,INSERT_VALUES));
 60:     else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
 61:   }

 63:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 64:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:              Create the singular value solver and set the solution method
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   /*
 71:      Create singular value context
 72:   */
 73:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));

 75:   /*
 76:      Set operator
 77:   */
 78:   PetscCall(SVDSetOperators(svd,A,NULL));

 80:   /*
 81:      Set solver parameters at runtime
 82:   */
 83:   PetscCall(SVDSetFromOptions(svd));
 84:   PetscCall(SVDSetDimensions(svd,1,PETSC_DETERMINE,PETSC_DETERMINE));

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                       Solve the singular value problem
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 90:   /*
 91:      First request a singular value from one end of the spectrum
 92:   */
 93:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
 94:   PetscCall(SVDSolve(svd));
 95:   /*
 96:      Get number of converged singular values
 97:   */
 98:   PetscCall(SVDGetConverged(svd,&nconv1));
 99:   /*
100:      Get converged singular values: largest singular value is stored in sigma_1.
101:      In this example, we are not interested in the singular vectors
102:   */
103:   if (nconv1 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL));
104:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n"));

106:   /*
107:      Request a singular value from the other end of the spectrum
108:   */
109:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
110:   PetscCall(SVDSolve(svd));
111:   /*
112:      Get number of converged singular triplets
113:   */
114:   PetscCall(SVDGetConverged(svd,&nconv2));
115:   /*
116:      Get converged singular values: smallest singular value is stored in sigma_n.
117:      As before, we are not interested in the singular vectors
118:   */
119:   if (nconv2 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL));
120:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n"));

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:                     Display solution and clean up
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125:   if (nconv1 > 0 && nconv2 > 0) {
126:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n));
127:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n)));
128:   }

130:   /*
131:      Free work space
132:   */
133:   PetscCall(SVDDestroy(&svd));
134:   PetscCall(MatDestroy(&A));
135:   PetscCall(SlepcFinalize());
136:   return 0;
137: }

139: /*TEST

141:    test:
142:       suffix: 1

144: TEST*/