Actual source code: test1.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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[] = "Test the solution of a SVD without calling SVDSetFromOptions (based on ex8.c).\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n"
 14:   "  -type <svd_type> = svd type to test.\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;
 41:   char           svdtype[30] = "cross",epstype[30] = "";
 42:   PetscBool      flg;
 43:   EPS            eps;

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

 48:   PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
 49:   PetscOptionsGetString(NULL,NULL,"-type",svdtype,30,NULL);
 50:   PetscOptionsGetString(NULL,NULL,"-epstype",epstype,30,&flg);
 51:   PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%D",N);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\n\n");

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:         Generate the matrix
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   MatCreate(PETSC_COMM_WORLD,&A);
 59:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 60:   MatSetFromOptions(A);
 61:   MatSetUp(A);

 63:   MatGetOwnershipRange(A,&Istart,&Iend);
 64:   for (i=Istart;i<Iend;i++) {
 65:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 66:     if (i==0) {
 67:       MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
 68:     } else {
 69:       MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
 70:     }
 71:   }

 73:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 74:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 76:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 77:          Create the singular value solver and set the solution method
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 80:   /*
 81:      Create singular value context
 82:   */
 83:   SVDCreate(PETSC_COMM_WORLD,&svd);

 85:   /*
 86:      Set operator
 87:   */
 88:   SVDSetOperator(svd,A);

 90:   /*
 91:      Set solver parameters at runtime
 92:   */
 93:   SVDSetType(svd,svdtype);
 94:   if (flg) {
 95:     PetscObjectTypeCompare((PetscObject)svd,SVDCROSS,&flg);
 96:     if (flg) {
 97:       SVDCrossGetEPS(svd,&eps);
 98:       EPSSetType(eps,epstype);
 99:     }
100:     PetscObjectTypeCompare((PetscObject)svd,SVDCYCLIC,&flg);
101:     if (flg) {
102:       SVDCyclicGetEPS(svd,&eps);
103:       EPSSetType(eps,epstype);
104:     }
105:   }
106:   SVDSetDimensions(svd,1,PETSC_DEFAULT,PETSC_DEFAULT);
107:   SVDSetTolerances(svd,1e-6,1000);

109:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110:                       Compute the singular values
111:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

113:   /*
114:      First request the largest singular value
115:   */
116:   SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
117:   SVDSolve(svd);
118:   /*
119:      Get number of converged singular values
120:   */
121:   SVDGetConverged(svd,&nconv1);
122:   /*
123:      Get converged singular values: largest singular value is stored in sigma_1.
124:      In this example, we are not interested in the singular vectors
125:   */
126:   if (nconv1 > 0) {
127:     SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL);
128:   } else {
129:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
130:   }

132:   /*
133:      Request the smallest singular value
134:   */
135:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
136:   SVDSolve(svd);
137:   /*
138:      Get number of converged triplets
139:   */
140:   SVDGetConverged(svd,&nconv2);
141:   /*
142:      Get converged singular values: smallest singular value is stored in sigma_n.
143:      As before, we are not interested in the singular vectors
144:   */
145:   if (nconv2 > 0) {
146:     SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL);
147:   } else {
148:     PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
149:   }

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:                     Display solution and clean up
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154:   if (nconv1 > 0 && nconv2 > 0) {
155:     PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n);
156:     PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n));
157:   }

159:   /*
160:      Free work space
161:   */
162:   SVDDestroy(&svd);
163:   MatDestroy(&A);
164:   SlepcFinalize();
165:   return ierr;
166: }

168: /*TEST

170:    test:
171:       suffix: 1
172:       args: -type {{lanczos trlanczos cross cyclic lapack}}

174:    test:
175:       suffix: 1_cross_gd
176:       args: -type cross -epstype gd
177:       output_file: output/test1_1.out

179:    test:
180:       suffix: 1_cyclic_gd
181:       args: -type cyclic -epstype gd
182:       output_file: output/test1_1.out

184:    test:
185:       suffix: 1_primme
186:       args: -type primme
187:       requires: primme
188:       output_file: output/test1_1.out

190: TEST*/