Actual source code: ex8.c

  2: static char help[] = "Wstimates the 2-norm condition number of a matrix A, that is, the ratio of the largest  to the smallest singular values of A. "
  3:   "The matrix is a Grcar matrix.\n\n"
  4:   "The command line options are:\n\n"
  5:   "  -n <n>, where <n> = matrix dimension.\n\n";

 7:  #include slepceps.h

  9: /*
 10:    This example computes the singular values of A by computing the eigenvalues
 11:    of A^T*A, where A^T denotes the transpose of A. 

 13:    An nxn Grcar matrix is a nonsymmetric Toeplitz matrix:

 15:               |  1  1  1  1               |
 16:               | -1  1  1  1  1            |
 17:               |    -1  1  1  1  1         |
 18:               |       .  .  .  .  .       |
 19:           A = |          .  .  .  .  .    |
 20:               |            -1  1  1  1  1 |
 21:               |               -1  1  1  1 |
 22:               |                  -1  1  1 |
 23:               |                     -1  1 |

 25:  */


 28: /* 
 29:    Matrix multiply routine
 30: */
 33: int MatSVD_Mult(Mat H,Vec x,Vec y)
 34: {
 35:   Mat      A;
 36:   Vec      w;
 37:   int      n, m, N, M, ierr;
 38:   MPI_Comm comm;

 40:   MatShellGetContext(H,(void**)&A);
 41:   PetscObjectGetComm((PetscObject)A,&comm);
 42:   MatGetLocalSize(A,&n,&m);
 43:   MatGetSize(A,&N,&M);
 44:   VecCreate(comm,&w);
 45:   VecSetSizes(w,n,N);
 46:   VecSetFromOptions(w);
 47:   MatMult(A,x,w);
 48:   MatMultTranspose(A,w,y);
 49:   VecDestroy(w);

 51:   return 0;
 52: }

 56: int main( int argc, char **argv )
 57: {
 58:   Mat         A;               /* Grcar matrix */
 59:   Mat         H;               /* eigenvalue problem matrix, H=A^T*A */
 60:   EPS         eps;             /* eigenproblem solver context */
 61:   int         N=30, n, ierr, i, nconv1, nconv2, col[5], Istart, Iend;
 62:   PetscScalar kl, ks, sigma_1, sigma_n, value[] = { -1, 1, 1, 1, 1 };

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

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

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 70:         Generate the matrix 
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
 74:   MatSetFromOptions(A);

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

 87:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 88:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 90:   /* 
 91:      Now create a symmetric shell matrix H=A^T*A
 92:   */
 93:   MatGetLocalSize(A,PETSC_NULL,&n);
 94:   MatCreateShell(PETSC_COMM_WORLD,n,n,PETSC_DETERMINE,PETSC_DETERMINE,(void*)A,&H);
 95:   MatShellSetOperation(H,MATOP_MULT,(void(*)())MatSVD_Mult);
 96:   MatShellSetOperation(H,MATOP_MULT_TRANSPOSE,(void(*)())MatSVD_Mult);

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 99:              Create the eigensolver and set the solution method
100:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

102:   /* 
103:      Create eigensolver context
104:   */
105:   EPSCreate(PETSC_COMM_WORLD,&eps);

107:   /* 
108:      Set operators. In this case, it is a standard symmetric eigenvalue problem
109:   */
110:   EPSSetOperators(eps,H,PETSC_NULL);

112:   /*
113:      Set solver parameters at runtime
114:   */
115: #ifdef SLEPC_HAVE_ARPACK
116:   EPSSetType(eps, EPSARPACK);
117: #else
118:   EPSSetType(eps, EPSLAPACK);
119: #endif
120:   EPSSetFromOptions(eps);
121:   EPSSetDimensions(eps,1,PETSC_DEFAULT);
122:   EPSSetTolerances(eps,PETSC_DEFAULT,1000);

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

128:   /*
129:      First request an eigenvalue from one end of the spectrum
130:   */
131:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
132:   EPSSolve(eps);
133:   /* 
134:      Get number of converged eigenpairs
135:   */
136:   EPSGetConverged(eps,&nconv1);
137:     /* 
138:        Get converged eigenpairs: largest eigenvalue is stored in kl. In this
139:        example, we are not interested in the eigenvectors
140:     */
141:   if (nconv1 > 0) {
142:     EPSGetEigenpair(eps,0,&kl,PETSC_NULL,PETSC_NULL,PETSC_NULL);
143:   }

145:   /*
146:      Request an eigenvalue from the other end of the spectrum
147:   */
148:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
149:   EPSSolve(eps);
150:   /* 
151:      Get number of converged eigenpairs
152:   */
153:   EPSGetConverged(eps,&nconv2);
154:   /* 
155:      Get converged eigenpairs: smallest eigenvalue is stored in ks. In this
156:      example, we are not interested in the eigenvectors
157:   */
158:   if (nconv2 > 0) {
159:     EPSGetEigenpair(eps,0,&ks,PETSC_NULL,PETSC_NULL,PETSC_NULL);
160:   }

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
163:                     Display solution and clean up
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   if (nconv1 > 0 && nconv2 > 0) {
166:     sigma_1 = PetscSqrtScalar(kl);
167:     sigma_n = PetscSqrtScalar(ks);

169:     PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%6f, sigma_n=%6f\n",sigma_1,sigma_n);
170:     PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%6f\n\n",sigma_1/sigma_n);
171:   } else {
172:     PetscPrintf(PETSC_COMM_WORLD," Process did not converge!\n\n");
173:   }
174: 
175:   /* 
176:      Free work space
177:   */
178:   EPSDestroy(eps);
179:   MatDestroy(A);
180:   MatDestroy(H);
181:   SlepcFinalize();
182:   return 0;
183: }