Actual source code: test8.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[] = "Tests multiple calls to SVDSolve changing ncv.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n\n";

 15: #include <slepcsvd.h>

 17: /*
 18:    This example computes the singular values of an nxn Grcar matrix,
 19:    which is a nonsymmetric Toeplitz matrix:

 21:               |  1  1  1  1               |
 22:               | -1  1  1  1  1            |
 23:               |    -1  1  1  1  1         |
 24:               |       .  .  .  .  .       |
 25:           A = |          .  .  .  .  .    |
 26:               |            -1  1  1  1  1 |
 27:               |               -1  1  1  1 |
 28:               |                  -1  1  1 |
 29:               |                     -1  1 |

 31:  */

 33: int main(int argc,char **argv)
 34: {
 35:   Mat            A;
 36:   SVD            svd;
 37:   PetscInt       N=30,Istart,Iend,i,col[5],nsv,ncv;
 38:   PetscScalar    value[] = { -1, 1, 1, 1, 1 };

 40:   PetscFunctionBeginUser;
 41:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 42:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
 43:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT,N));
 44:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));

 46:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 47:         Generate the matrix
 48:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 50:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 51:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 52:   PetscCall(MatSetFromOptions(A));

 54:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 55:   for (i=Istart;i<Iend;i++) {
 56:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 57:     if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES));
 58:     else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
 59:   }

 61:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 62:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 64:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 65:          Create the singular value solver and set the solution method
 66:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 68:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 69:   PetscCall(SVDSetOperators(svd,A,NULL));
 70:   PetscCall(SVDSetTolerances(svd,1e-6,1000));
 71:   PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
 72:   PetscCall(SVDSetFromOptions(svd));

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:                       Compute the singular values
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   /* First solve */
 79:   PetscCall(SVDSolve(svd));
 80:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," - - - First solve, default subspace dimension - - -\n"));
 81:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));

 83:   /* Second solve */
 84:   PetscCall(SVDGetDimensions(svd,&nsv,&ncv,NULL));
 85:   PetscCall(SVDSetDimensions(svd,nsv,ncv+2,PETSC_DEFAULT));
 86:   PetscCall(SVDSolve(svd));
 87:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," - - - Second solve, subspace of increased size - - -\n"));
 88:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));

 90:   /* Free work space */
 91:   PetscCall(SVDDestroy(&svd));
 92:   PetscCall(MatDestroy(&A));
 93:   PetscCall(SlepcFinalize());
 94:   return 0;
 95: }

 97: /*TEST

 99:    test:
100:       suffix: 1
101:       args: -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -svd_nsv 3 -svd_ncv 12
102:       requires: !single

104: TEST*/