Actual source code: test9.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 with different matrix size.\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,B;
 36:   SVD            svd;
 37:   PetscInt       N=30,Istart,Iend,i,col[5];
 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\n",N));

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 46:         Generate the matrix of size N
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 49:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 50:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 51:   PetscCall(MatSetFromOptions(A));
 52:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 53:   for (i=Istart;i<Iend;i++) {
 54:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 55:     if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES));
 56:     else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
 57:   }
 58:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 59:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:          Create the singular value solver, set options and solve
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
 66:   PetscCall(SVDSetOperators(svd,A,NULL));
 67:   PetscCall(SVDSetTolerances(svd,1e-6,1000));
 68:   PetscCall(SVDSetFromOptions(svd));
 69:   PetscCall(SVDSolve(svd));
 70:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 73:         Generate the matrix of size 2*N
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   N *= 2;
 77:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N));

 79:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 80:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
 81:   PetscCall(MatSetFromOptions(B));
 82:   PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
 83:   for (i=Istart;i<Iend;i++) {
 84:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 85:     if (i==0) PetscCall(MatSetValues(B,1,&i,4,col+1,value+1,INSERT_VALUES));
 86:     else PetscCall(MatSetValues(B,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
 87:   }
 88:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 89:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:        Solve again, calling SVDReset() since matrix size has changed
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   PetscCall(SVDReset(svd));  /* if this is omitted, it will be called in SVDSetOperators() */
 96:   PetscCall(SVDSetOperators(svd,B,NULL));
 97:   PetscCall(SVDSolve(svd));
 98:   PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));

100:   /* Free work space */
101:   PetscCall(SVDDestroy(&svd));
102:   PetscCall(MatDestroy(&A));
103:   PetscCall(MatDestroy(&B));
104:   PetscCall(SlepcFinalize());
105:   return 0;
106: }

108: /*TEST

110:    test:
111:       suffix: 1
112:       args: -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -svd_nsv 3
113:       requires: double

115: TEST*/