Actual source code: test1.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[] = "Test the solution of a SVD without calling SVDSetFromOptions (based on ex8.c).\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = matrix dimension.\n"
 25:   "  -type <svd_type> = svd type to test.\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;
 54:   char           svdtype[30] = "cross", epstype[30] = "";
 55:   PetscBool      flg;
 56:   EPS            eps;

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

 61:   PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
 62:   PetscOptionsGetString(PETSC_NULL,"-type",svdtype,30,PETSC_NULL);
 63:   PetscOptionsGetString(PETSC_NULL,"-epstype",epstype,30,&flg);
 64:   PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%D",N);
 65:   PetscPrintf(PETSC_COMM_WORLD,"\nSVD type: %s",svdtype);
 66:   if (flg) {
 67:     PetscPrintf(PETSC_COMM_WORLD,"\nEPS type: %s",epstype);
 68:   }
 69:   PetscPrintf(PETSC_COMM_WORLD,"\n\n");
 70: 

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 73:         Generate the matrix 
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   MatCreate(PETSC_COMM_WORLD,&A);
 77:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 78:   MatSetFromOptions(A);

 80:   MatGetOwnershipRange(A,&Istart,&Iend);
 81:   for (i=Istart;i<Iend;i++) {
 82:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 83:     if (i==0) {
 84:       MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
 85:     } else {
 86:       MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
 87:     }
 88:   }

 90:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 91:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 94:              Create the singular value solver and set the solution method
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 97:   /* 
 98:      Create singular value context
 99:   */
100:   SVDCreate(PETSC_COMM_WORLD,&svd);

102:   /* 
103:      Set operator
104:   */
105:   SVDSetOperator(svd,A);

107:   /*
108:      Set solver parameters at runtime
109:   */
110:   SVDSetType(svd,svdtype);
111:   if (flg) {
112:     PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg);
113:     if (flg) {
114:       SVDCrossGetEPS(svd,&eps);
115:       EPSSetType(eps,epstype);
116:     }
117:     PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg);
118:     if (flg) {
119:       SVDCyclicGetEPS(svd,&eps);
120:       EPSSetType(eps,epstype);
121:     }
122:   }
123:   SVDSetDimensions(svd,1,PETSC_IGNORE,PETSC_IGNORE);

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
126:                       Solve the eigensystem
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

129:   /*
130:      First request an eigenvalue from one end of the spectrum
131:   */
132:   SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
133:   SVDSolve(svd);
134:   /* 
135:      Get number of converged singular values
136:   */
137:   SVDGetConverged(svd,&nconv1);
138:   /* 
139:      Get converged singular values: largest singular value is stored in sigma_1.
140:      In this example, we are not interested in the singular vectors
141:   */
142:   if (nconv1 > 0) {
143:     SVDGetSingularTriplet(svd,0,&sigma_1,PETSC_NULL,PETSC_NULL);
144:   } else {
145:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
146:   }

148:   /*
149:      Request an eigenvalue from the other end of the spectrum
150:   */
151:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
152:   SVDSolve(svd);
153:   /* 
154:      Get number of converged eigenpairs
155:   */
156:   SVDGetConverged(svd,&nconv2);
157:   /* 
158:      Get converged singular values: smallest singular value is stored in sigma_n. 
159:      As before, we are not interested in the singular vectors
160:   */
161:   if (nconv2 > 0) {
162:     SVDGetSingularTriplet(svd,0,&sigma_n,PETSC_NULL,PETSC_NULL);
163:   } else {
164:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
165:   }

167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
168:                     Display solution and clean up
169:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170:   if (nconv1 > 0 && nconv2 > 0) {
171:     PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%6F, sigma_n=%6F\n",sigma_1,sigma_n);
172:     PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%6F\n\n",sigma_1/sigma_n);
173:   }
174: 
175:   /* 
176:      Free work space
177:   */
178:   SVDDestroy(&svd);
179:   MatDestroy(&A);
180:   SlepcFinalize();
181:   return 0;
182: }