Actual source code: ex15.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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.

 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.

 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/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 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";

 27: #include <slepcsvd.h>

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            A;               /* operator matrix */
 34:   Vec            u,v;             /* left and right singular vectors */
 35:   SVD            svd;             /* singular value problem solver context */
 36:   SVDType        type;
 37:   PetscReal      error,tol,sigma,mu=PETSC_SQRT_MACHINE_EPSILON;
 38:   PetscInt       n=100,i,j,Istart,Iend,nsv,maxit,its,nconv;

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

 43:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 44:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 45:   PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%D x %D) mu=%g\n\n",n+1,n,(double)mu);

 47:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 48:                           Build the Lauchli matrix
 49:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 51:   MatCreate(PETSC_COMM_WORLD,&A);
 52:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n);
 53:   MatSetFromOptions(A);
 54:   MatSetUp(A);

 56:   MatGetOwnershipRange(A,&Istart,&Iend);
 57:   for (i=Istart;i<Iend;i++) {
 58:     if (i == 0) {
 59:       for (j=0;j<n;j++) {
 60:         MatSetValue(A,0,j,1.0,INSERT_VALUES);
 61:       }
 62:     } else {
 63:       MatSetValue(A,i,i-1,mu,INSERT_VALUES);
 64:     }
 65:   }

 67:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 69:   MatCreateVecs(A,&v,&u);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:           Create the singular value solver and set various options
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 75:   /*
 76:      Create singular value solver context
 77:   */
 78:   SVDCreate(PETSC_COMM_WORLD,&svd);

 80:   /*
 81:      Set operator
 82:   */
 83:   SVDSetOperator(svd,A);

 85:   /*
 86:      Use thick-restart Lanczos as default solver
 87:   */
 88:   SVDSetType(svd,SVDTRLANCZOS);

 90:   /*
 91:      Set solver parameters at runtime
 92:   */
 93:   SVDSetFromOptions(svd);

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                       Solve the singular value system
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   SVDSolve(svd);
100:   SVDGetIterationNumber(svd,&its);
101:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);

103:   /*
104:      Optional: Get some information from the solver and display it
105:   */
106:   SVDGetType(svd,&type);
107:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
108:   SVDGetDimensions(svd,&nsv,NULL,NULL);
109:   PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %D\n",nsv);
110:   SVDGetTolerances(svd,&tol,&maxit);
111:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:                     Display solution and clean up
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

117:   /*
118:      Get number of converged singular triplets
119:   */
120:   SVDGetConverged(svd,&nconv);
121:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %D\n\n",nconv);

123:   if (nconv>0) {
124:     /*
125:        Display singular values and relative errors
126:     */
127:     PetscPrintf(PETSC_COMM_WORLD,
128:          "          sigma           relative error\n"
129:          "  --------------------- ------------------\n");
130:     for (i=0;i<nconv;i++) {
131:       /*
132:          Get converged singular triplets: i-th singular value is stored in sigma
133:       */
134:       SVDGetSingularTriplet(svd,i,&sigma,u,v);

136:       /*
137:          Compute the error associated to each singular triplet
138:       */
139:       SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error);

141:       PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",(double)sigma);
142:       PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
143:     }
144:     PetscPrintf(PETSC_COMM_WORLD,"\n");
145:   }

147:   /*
148:      Free work space
149:   */
150:   SVDDestroy(&svd);
151:   MatDestroy(&A);
152:   VecDestroy(&u);
153:   VecDestroy(&v);
154:   SlepcFinalize();
155:   return ierr;
156: }