Actual source code: test9.c

slepc-3.17.0 2022-03-31
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:   SlepcInitialize(&argc,&argv,(char*)0,help);
 41:   PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
 42:   PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N);

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

 48:   MatCreate(PETSC_COMM_WORLD,&A);
 49:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 50:   MatSetFromOptions(A);
 51:   MatSetUp(A);
 52:   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) MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
 56:     else MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
 57:   }
 58:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

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

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

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

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

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

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

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

109: /*TEST

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

116: TEST*/