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
1296 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
 
1296 slepc 22
static char help[] = "Singular value decomposition of the Lauchli matrix.\n"
23
  "The command line options are:\n"
24
  "  -n <n>, where <n> = matrix dimension.\n"
25
  "  -mu <mu>, where <mu> = subdiagonal value.\n\n";
26
 
2283 jroman 27
#include <slepcsvd.h>
1296 slepc 28
 
29
#undef __FUNCT__
30
#define __FUNCT__ "main"
31
int main( int argc, char **argv )
32
{
2075 jroman 33
  Mat            A;               /* operator matrix */
34
  Vec            u,v;             /* left and right singular vectors */
35
  SVD            svd;             /* singular value problem solver context */
1507 slepc 36
  const SVDType  type;
1317 slepc 37
  PetscReal      error, tol, sigma, mu=PETSC_SQRT_MACHINE_EPSILON;
1296 slepc 38
  PetscErrorCode ierr;
1505 slepc 39
  PetscInt       n=100, i, j, Istart, Iend, nsv, maxit, its, nconv;
1296 slepc 40
 
41
  SlepcInitialize(&argc,&argv,(char*)0,help);
42
 
43
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
44
  ierr = PetscOptionsGetReal(PETSC_NULL,"-mu",&mu,PETSC_NULL);CHKERRQ(ierr);
45
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%d x %d) mu=%g\n\n",n+1,n,mu);CHKERRQ(ierr);
46
 
47
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48
                          Build the Lauchli matrix
49
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50
 
51
  ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
52
  ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n);CHKERRQ(ierr);
53
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
54
 
1505 slepc 55
  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
56
  for (i=Istart;i<Iend;i++) {
1296 slepc 57
    if (i == 0) {
58
      for (j=0;j<n;j++) {
59
        ierr = MatSetValue(A,0,j,1.0,INSERT_VALUES);CHKERRQ(ierr);
60
      }
61
    } else {
62
      ierr = MatSetValue(A,i,i-1,mu,INSERT_VALUES);CHKERRQ(ierr);
63
    }
64
  }
65
 
66
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2075 jroman 68
  ierr = MatGetVecs(A,&v,&u);CHKERRQ(ierr);
1296 slepc 69
 
70
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71
          Create the singular value solver and set various options
72
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73
 
74
  /*
75
     Create singular value solver context
76
  */
77
  ierr = SVDCreate(PETSC_COMM_WORLD,&svd);CHKERRQ(ierr);
78
 
79
  /*
80
     Set operator
81
  */
82
  ierr = SVDSetOperator(svd,A);CHKERRQ(ierr);
1391 slepc 83
 
84
  /*
85
     Use thick-restart Lanczos as default solver
86
  */
87
  ierr = SVDSetType(svd,SVDTRLANCZOS);CHKERRQ(ierr);
1296 slepc 88
 
89
  /*
90
     Set solver parameters at runtime
91
  */
92
  ierr = SVDSetFromOptions(svd);CHKERRQ(ierr);
93
 
94
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95
                      Solve the singular value system
96
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97
 
98
  ierr = SVDSolve(svd);CHKERRQ(ierr);
99
  ierr = SVDGetIterationNumber(svd, &its);CHKERRQ(ierr);
100
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
101
 
102
  /*
103
     Optional: Get some information from the solver and display it
104
  */
105
  ierr = SVDGetType(svd,&type);CHKERRQ(ierr);
106
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
1575 slepc 107
  ierr = SVDGetDimensions(svd,&nsv,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1296 slepc 108
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %d\n",nsv);CHKERRQ(ierr);
109
  ierr = SVDGetTolerances(svd,&tol,&maxit);CHKERRQ(ierr);
110
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
111
 
112
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113
                    Display solution and clean up
114
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115
 
116
  /*
117
     Get number of converged singular triplets
118
  */
119
  ierr = SVDGetConverged(svd,&nconv);CHKERRQ(ierr);
120
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %d\n\n",nconv);CHKERRQ(ierr);
121
 
122
  if (nconv>0) {
123
    /*
124
       Display singular values and relative errors
125
    */
126
    ierr = PetscPrintf(PETSC_COMM_WORLD,
127
         "          sigma           residual norm\n"
128
         "  --------------------- ------------------\n" );CHKERRQ(ierr);
129
    for( i=0; i<nconv; i++ ) {
130
      /*
131
         Get converged singular triplets: i-th singular value is stored in sigma
132
      */
2075 jroman 133
      ierr = SVDGetSingularTriplet(svd,i,&sigma,u,v);CHKERRQ(ierr);
1296 slepc 134
 
135
      /*
136
         Compute the error associated to each singular triplet
137
      */
1317 slepc 138
      ierr = SVDComputeRelativeError(svd,i,&error);CHKERRQ(ierr);
1296 slepc 139
 
140
      ierr = PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",sigma); CHKERRQ(ierr);
141
      ierr = PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);CHKERRQ(ierr);
142
    }
143
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
144
  }
145
 
146
  /*
147
     Free work space
148
  */
149
  ierr = SVDDestroy(svd);CHKERRQ(ierr);
150
  ierr = MatDestroy(A);CHKERRQ(ierr);
2075 jroman 151
  ierr = VecDestroy(u);CHKERRQ(ierr);
152
  ierr = VecDestroy(v);CHKERRQ(ierr);
1296 slepc 153
  ierr = SlepcFinalize();CHKERRQ(ierr);
154
  return 0;
155
}