Actual source code: ex8.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: 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. "
 23:   "The matrix is a Grcar matrix.\n\n"
 24:   "The command line options are:\n"
 25:   "  -n <n>, where <n> = matrix dimension.\n\n";

 27: #include <slepcsvd.h>

 29: /*
 30:    This example computes the singular values of an nxn Grcar matrix,
 31:    which is a nonsymmetric Toeplitz matrix:

 33:               |  1  1  1  1               |
 34:               | -1  1  1  1  1            |
 35:               |    -1  1  1  1  1         |
 36:               |       .  .  .  .  .       |
 37:           A = |          .  .  .  .  .    |
 38:               |            -1  1  1  1  1 |
 39:               |               -1  1  1  1 |
 40:               |                  -1  1  1 |
 41:               |                     -1  1 |

 43:  */

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

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

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

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 62:         Generate the matrix 
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   MatCreate(PETSC_COMM_WORLD,&A);
 66:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 67:   MatSetFromOptions(A);

 69:   MatGetOwnershipRange(A,&Istart,&Iend);
 70:   for (i=Istart;i<Iend;i++) {
 71:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 72:     if (i==0) {
 73:       MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
 74:     } else {
 75:       MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
 76:     }
 77:   }

 79:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 83:              Create the singular value solver and set the solution method
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   /* 
 87:      Create singular value context
 88:   */
 89:   SVDCreate(PETSC_COMM_WORLD,&svd);

 91:   /* 
 92:      Set operator
 93:   */
 94:   SVDSetOperator(svd,A);

 96:   /*
 97:      Set solver parameters at runtime
 98:   */
 99:   SVDSetFromOptions(svd);
100:   SVDSetDimensions(svd,1,PETSC_IGNORE,PETSC_IGNORE);

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
103:                       Solve the singular value problem
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

106:   /*
107:      First request a singular value from one end of the spectrum
108:   */
109:   SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
110:   SVDSolve(svd);
111:   /* 
112:      Get number of converged singular values
113:   */
114:   SVDGetConverged(svd,&nconv1);
115:   /* 
116:      Get converged singular values: largest singular value is stored in sigma_1.
117:      In this example, we are not interested in the singular vectors
118:   */
119:   if (nconv1 > 0) {
120:     SVDGetSingularTriplet(svd,0,&sigma_1,PETSC_NULL,PETSC_NULL);
121:   } else {
122:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
123:   }

125:   /*
126:      Request a singular value from the other end of the spectrum
127:   */
128:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
129:   SVDSolve(svd);
130:   /* 
131:      Get number of converged eigenpairs
132:   */
133:   SVDGetConverged(svd,&nconv2);
134:   /* 
135:      Get converged singular values: smallest singular value is stored in sigma_n. 
136:      As before, we are not interested in the singular vectors
137:   */
138:   if (nconv2 > 0) {
139:     SVDGetSingularTriplet(svd,0,&sigma_n,PETSC_NULL,PETSC_NULL);
140:   } else {
141:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
142:   }

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
145:                     Display solution and clean up
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147:   if (nconv1 > 0 && nconv2 > 0) {
148:     PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%6F, sigma_n=%6F\n",sigma_1,sigma_n);
149:     PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%6F\n\n",sigma_1/sigma_n);
150:   }
151: 
152:   /* 
153:      Free work space
154:   */
155:   SVDDestroy(&svd);
156:   MatDestroy(&A);
157:   SlepcFinalize();
158:   return 0;
159: }