Actual source code: test20.c

slepc-3.21.1 2024-04-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test setting and getting the signature in HSVD.\n"
 12:   "Based on ex15.c. The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n"
 14:   "  -mu <mu>, where <mu> = subdiagonal value.\n\n";

 16: #include <slepcsvd.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A;               /* operator matrix */
 21:   Vec            omega,v;         /* signature */
 22:   SVD            svd;             /* singular value problem solver context */
 23:   PetscReal      mu=PETSC_SQRT_MACHINE_EPSILON;
 24:   PetscInt       n=50,i,j,Istart,Iend;

 26:   PetscFunctionBeginUser;
 27:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 30:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
 31:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLauchli hyperbolic singular value decomposition, (%" PetscInt_FMT " x %" PetscInt_FMT ") mu=%g\n\n",n+1,n,(double)mu));

 33:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34:                    Build Lauchli matrix and signature
 35:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 37:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 38:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n));
 39:   PetscCall(MatSetFromOptions(A));

 41:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 42:   for (i=Istart;i<Iend;i++) {
 43:     if (i == 0) {
 44:       for (j=0;j<n;j++) PetscCall(MatSetValue(A,0,j,1.0,INSERT_VALUES));
 45:     } else PetscCall(MatSetValue(A,i,i-1,mu,INSERT_VALUES));
 46:   }
 47:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 48:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 50:   /* signature omega = [ -1 1 1 ... 1 ] */
 51:   PetscCall(MatCreateVecs(A,NULL,&omega));
 52:   PetscCall(VecSet(omega,1.0));
 53:   PetscCall(VecGetOwnershipRange(omega,&Istart,NULL));
 54:   if (Istart==0) PetscCall(VecSetValue(omega,0,-1.0,INSERT_VALUES));
 55:   PetscCall(VecAssemblyBegin(omega));
 56:   PetscCall(VecAssemblyEnd(omega));

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:           Create the singular value solver and set various options
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 63:   PetscCall(SVDSetSignature(svd,omega));
 64:   PetscCall(SVDSetOperators(svd,A,NULL));
 65:   PetscCall(SVDSetType(svd,SVDTRLANCZOS));
 66:   PetscCall(SVDSetFromOptions(svd));

 68:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69:                 Solve the problem and display solution
 70:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 72:   PetscCall(SVDSolve(svd));
 73:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));

 75:   /* Test getting the signature */
 76:   PetscCall(MatCreateVecs(A,NULL,&v));
 77:   PetscCall(SVDGetSignature(svd,v));
 78:   PetscCall(VecView(v,NULL));

 80:   /* Free work space */
 81:   PetscCall(SVDDestroy(&svd));
 82:   PetscCall(MatDestroy(&A));
 83:   PetscCall(VecDestroy(&omega));
 84:   PetscCall(VecDestroy(&v));
 85:   PetscCall(SlepcFinalize());
 86:   return 0;
 87: }

 89: /*TEST

 91:    test:
 92:       requires: double
 93:       suffix: 1

 95: TEST*/