Actual source code: test14.c

slepc-3.17.1 2022-04-11
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 equal matrix size.\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = matrix rows.\n"
 14:   "  -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";

 16: #include <slepcsvd.h>

 18: /*
 19:    This example computes the singular values of two rectangular bidiagonal matrices

 21:               |  1  2                     |       |  1                        |
 22:               |     1  2                  |       |  2  1                     |
 23:               |        1  2               |       |     2  1                  |
 24:           A = |          .  .             |   B = |       .  .                |
 25:               |             .  .          |       |          .  .             |
 26:               |                1  2       |       |             2  1          |
 27:               |                   1  2    |       |                2  1       |
 28:  */

 30: int main(int argc,char **argv)
 31: {
 32:   Mat            A,B;
 33:   SVD            svd;
 34:   PetscInt       m=20,n,Istart,Iend,i,col[2];
 35:   PetscScalar    valsa[] = { 1, 2 }, valsb[] = { 2, 1 };
 36:   PetscBool      flg;

 38:   SlepcInitialize(&argc,&argv,(char*)0,help);
 39:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 40:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 41:   if (!flg) n=m+2;
 42:   PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n);

 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:                      Generate the matrices
 46:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 48:   MatCreate(PETSC_COMM_WORLD,&A);
 49:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 50:   MatSetFromOptions(A);
 51:   MatSetUp(A);
 52:   MatGetOwnershipRange(A,&Istart,&Iend);
 53:   for (i=Istart;i<Iend;i++) {
 54:     col[0]=i; col[1]=i+1;
 55:     if (i<n-1) MatSetValues(A,1,&i,2,col,valsa,INSERT_VALUES);
 56:     else if (i==n-1) MatSetValue(A,i,col[0],valsa[0],INSERT_VALUES);
 57:   }
 58:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 61:   MatCreate(PETSC_COMM_WORLD,&B);
 62:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
 63:   MatSetFromOptions(B);
 64:   MatSetUp(B);
 65:   MatGetOwnershipRange(B,&Istart,&Iend);
 66:   for (i=Istart;i<Iend;i++) {
 67:     col[0]=i-1; col[1]=i;
 68:     if (i==0) MatSetValue(B,i,col[1],valsb[1],INSERT_VALUES);
 69:     else if (i<n) MatSetValues(B,1,&i,2,col,valsb,INSERT_VALUES);
 70:     else if (i==n) MatSetValue(B,i,col[0],valsb[0],INSERT_VALUES);
 71:   }
 72:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:          Create the singular value solver, set options and solve
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   SVDCreate(PETSC_COMM_WORLD,&svd);
 80:   SVDSetOperators(svd,A,NULL);
 81:   SVDSetTolerances(svd,PETSC_DEFAULT,1000);
 82:   SVDSetFromOptions(svd);
 83:   SVDSolve(svd);
 84:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                        Solve with second matrix
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 90:   SVDSetOperators(svd,B,NULL);
 91:   SVDSolve(svd);
 92:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);

 94:   /* Free work space */
 95:   SVDDestroy(&svd);
 96:   MatDestroy(&A);
 97:   MatDestroy(&B);
 98:   SlepcFinalize();
 99:   return 0;
100: }

102: /*TEST

104:    testset:
105:       args: -svd_nsv 3
106:       requires: !single
107:       output_file: output/test14_1.out
108:       test:
109:          suffix: 1
110:          args: -svd_type {{lanczos trlanczos lapack}}
111:       test:
112:          suffix: 1_cross
113:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
114:       test:
115:          suffix: 1_cyclic
116:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}

118:    testset:
119:       args: -n 18 -svd_nsv 3
120:       requires: !single
121:       output_file: output/test14_2.out
122:       test:
123:          suffix: 2
124:          args: -svd_type {{lanczos trlanczos lapack}}
125:       test:
126:          suffix: 2_cross
127:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
128:       test:
129:          suffix: 2_cyclic
130:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}

132: TEST*/