Actual source code: test1.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, 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);
 79:   MatSetUp(A);

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

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

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

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

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

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

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
128:                       Solve the eigensystem
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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