Subversion Repositories slepc-dev

Rev

Go to most recent revision | 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
2575 eromero 4
   Copyright (c) 2002-2011, Universitat 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"
2331 jroman 31
int main(int argc,char **argv)
1296 slepc 32
{
2316 jroman 33
  Mat            A;               /* operator matrix */
2075 jroman 34
  Vec            u,v;             /* left and right singular vectors */
2316 jroman 35
  SVD            svd;             /* singular value problem solver context */
1507 slepc 36
  const SVDType  type;
2331 jroman 37
  PetscReal      error,tol,sigma,mu=PETSC_SQRT_MACHINE_EPSILON;
38
  PetscInt       n=100,i,j,Istart,Iend,nsv,maxit,its,nconv;
2317 jroman 39
  PetscErrorCode ierr;
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);
2515 jroman 45
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%D x %D) mu=%G\n\n",n+1,n,mu);CHKERRQ(ierr);
1296 slepc 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);
2717 jroman 54
  ierr = MatSetUp(A);CHKERRQ(ierr);
1296 slepc 55
 
1505 slepc 56
  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
57
  for (i=Istart;i<Iend;i++) {
1296 slepc 58
    if (i == 0) {
59
      for (j=0;j<n;j++) {
60
        ierr = MatSetValue(A,0,j,1.0,INSERT_VALUES);CHKERRQ(ierr);
61
      }
62
    } else {
63
      ierr = MatSetValue(A,i,i-1,mu,INSERT_VALUES);CHKERRQ(ierr);
64
    }
65
  }
66
 
67
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2075 jroman 69
  ierr = MatGetVecs(A,&v,&u);CHKERRQ(ierr);
1296 slepc 70
 
71
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72
          Create the singular value solver and set various options
73
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74
 
75
  /*
76
     Create singular value solver context
77
  */
78
  ierr = SVDCreate(PETSC_COMM_WORLD,&svd);CHKERRQ(ierr);
79
 
80
  /*
81
     Set operator
82
  */
83
  ierr = SVDSetOperator(svd,A);CHKERRQ(ierr);
1391 slepc 84
 
85
  /*
86
     Use thick-restart Lanczos as default solver
87
  */
88
  ierr = SVDSetType(svd,SVDTRLANCZOS);CHKERRQ(ierr);
1296 slepc 89
 
90
  /*
91
     Set solver parameters at runtime
92
  */
93
  ierr = SVDSetFromOptions(svd);CHKERRQ(ierr);
94
 
95
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96
                      Solve the singular value system
97
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98
 
99
  ierr = SVDSolve(svd);CHKERRQ(ierr);
2331 jroman 100
  ierr = SVDGetIterationNumber(svd,&its);CHKERRQ(ierr);
2515 jroman 101
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);CHKERRQ(ierr);
1296 slepc 102
 
103
  /*
104
     Optional: Get some information from the solver and display it
105
  */
106
  ierr = SVDGetType(svd,&type);CHKERRQ(ierr);
107
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
1575 slepc 108
  ierr = SVDGetDimensions(svd,&nsv,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2515 jroman 109
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %D\n",nsv);CHKERRQ(ierr);
1296 slepc 110
  ierr = SVDGetTolerances(svd,&tol,&maxit);CHKERRQ(ierr);
2515 jroman 111
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);CHKERRQ(ierr);
1296 slepc 112
 
113
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114
                    Display solution and clean up
115
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116
 
117
  /*
118
     Get number of converged singular triplets
119
  */
120
  ierr = SVDGetConverged(svd,&nconv);CHKERRQ(ierr);
2515 jroman 121
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %D\n\n",nconv);CHKERRQ(ierr);
1296 slepc 122
 
123
  if (nconv>0) {
124
    /*
125
       Display singular values and relative errors
126
    */
127
    ierr = PetscPrintf(PETSC_COMM_WORLD,
2515 jroman 128
         "          sigma           relative error\n"
2331 jroman 129
         "  --------------------- ------------------\n");CHKERRQ(ierr);
130
    for (i=0;i<nconv;i++) {
1296 slepc 131
      /*
132
         Get converged singular triplets: i-th singular value is stored in sigma
133
      */
2075 jroman 134
      ierr = SVDGetSingularTriplet(svd,i,&sigma,u,v);CHKERRQ(ierr);
1296 slepc 135
 
136
      /*
137
         Compute the error associated to each singular triplet
138
      */
1317 slepc 139
      ierr = SVDComputeRelativeError(svd,i,&error);CHKERRQ(ierr);
1296 slepc 140
 
2515 jroman 141
      ierr = PetscPrintf(PETSC_COMM_WORLD,"       % 6F      ",sigma);CHKERRQ(ierr);
142
      ierr = PetscPrintf(PETSC_COMM_WORLD," % 12G\n",error);CHKERRQ(ierr);
1296 slepc 143
    }
2331 jroman 144
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr);
1296 slepc 145
  }
146
 
147
  /*
148
     Free work space
149
  */
2312 jroman 150
  ierr = SVDDestroy(&svd);CHKERRQ(ierr);
2305 jroman 151
  ierr = MatDestroy(&A);CHKERRQ(ierr);
152
  ierr = VecDestroy(&u);CHKERRQ(ierr);
153
  ierr = VecDestroy(&v);CHKERRQ(ierr);
1296 slepc 154
  ierr = SlepcFinalize();CHKERRQ(ierr);
155
  return 0;
156
}