Actual source code: ex14.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Solves a singular value problem with the matrix loaded from a file.\n"
12: "This example works for both real and complex numbers.\n\n"
13: "The command line options are:\n"
14: " -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n";
16: #include slepcsvd.h
20: int main( int argc, char **argv )
21: {
22: Mat A; /* operator matrix */
23: SVD svd; /* singular value problem solver context */
24: SVDType type;
25: PetscReal error, tol, sigma;
27: int nsv, maxit, i, its, nconv;
28: char filename[256];
29: PetscViewer viewer;
30: PetscTruth flg;
33: SlepcInitialize(&argc,&argv,(char*)0,help);
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Load the operator matrix that defines the singular value problem
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: PetscPrintf(PETSC_COMM_WORLD,"\nSingular value problem stored in file.\n\n");
40: PetscOptionsGetString(PETSC_NULL,"-file",filename,256,&flg);
41: if (!flg) {
42: SETERRQ(1,"Must indicate a file name with the -file option.");
43: }
45: #if defined(PETSC_USE_COMPLEX)
46: PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
47: #else
48: PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
49: #endif
50: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
51: MatLoad(viewer,MATAIJ,&A);
52: PetscViewerDestroy(viewer);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Create the singular value solver and set various options
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: /*
59: Create singular value solver context
60: */
61: SVDCreate(PETSC_COMM_WORLD,&svd);
63: /*
64: Set operator
65: */
66: SVDSetOperator(svd,A);
68: /*
69: Set solver parameters at runtime
70: */
71: SVDSetFromOptions(svd);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Solve the singular value system
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: SVDSolve(svd);
78: SVDGetIterationNumber(svd, &its);
79: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
81: /*
82: Optional: Get some information from the solver and display it
83: */
84: SVDGetType(svd,&type);
85: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
86: SVDGetDimensions(svd,&nsv,PETSC_NULL);
87: PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %d\n",nsv);
88: SVDGetTolerances(svd,&tol,&maxit);
89: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Display solution and clean up
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: /*
96: Get number of converged singular triplets
97: */
98: SVDGetConverged(svd,&nconv);
99: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %d\n\n",nconv);
101: if (nconv>0) {
102: /*
103: Display singular values and relative errors
104: */
105: PetscPrintf(PETSC_COMM_WORLD,
106: " sigma residual norm\n"
107: " --------------------- ------------------\n" );
108: for( i=0; i<nconv; i++ ) {
109: /*
110: Get converged singular triplets: i-th singular value is stored in sigma
111: */
112: SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
114: /*
115: Compute the error associated to each singular triplet
116: */
117: SVDComputeRelativeError(svd,i,&error);
118: PetscPrintf(PETSC_COMM_WORLD," % 6f ",sigma);
119: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
120: }
121: PetscPrintf(PETSC_COMM_WORLD,"\n" );
122: }
123:
124: /*
125: Free work space
126: */
127: SVDDestroy(svd);
128: MatDestroy(A);
129: SlepcFinalize();
130: return 0;
131: }