Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

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