Actual source code: test16.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 equal matrix size (GSVD).\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = number of rows of A.\n"
 14:   "  -n <n>, where <n> = number of columns of A.\n"
 15:   "  -p <p>, where <p> = number of rows of B.\n\n";

 17: #include <slepcsvd.h>

 19: /*
 20:    This example solves two GSVD problems for the bidiagonal matrices

 22:               |  1  2                     |       |  1                        |
 23:               |     1  2                  |       |  2  1                     |
 24:               |        1  2               |       |     2  1                  |
 25:          A1 = |          .  .             |  A2 = |       .  .                |
 26:               |             .  .          |       |          .  .             |
 27:               |                1  2       |       |             2  1          |
 28:               |                   1  2    |       |                2  1       |

 30:    with B = tril(ones(p,n))
 31:  */

 33: int main(int argc,char **argv)
 34: {
 35:   Mat            A1,A2,B;
 36:   SVD            svd;
 37:   PetscInt       m=15,n=20,p=21,Istart,Iend,i,j,d,col[2];
 38:   PetscScalar    valsa[] = { 1, 2 }, valsb[] = { 2, 1 };

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);
 41:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
 44:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n);

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

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

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

 77:   MatCreate(PETSC_COMM_WORLD,&B);
 78:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 79:   MatSetFromOptions(B);
 80:   MatSetUp(B);
 81:   MatGetOwnershipRange(B,&Istart,&Iend);
 82:   d = PetscMax(0,n-p);
 83:   for (i=Istart;i<Iend;i++) {
 84:     for (j=0;j<=PetscMin(i,n-1);j++) MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
 85:   }
 86:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 87:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:          Create the singular value solver, set options and solve
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 93:   SVDCreate(PETSC_COMM_WORLD,&svd);
 94:   SVDSetOperators(svd,A1,B);
 95:   SVDSetFromOptions(svd);
 96:   SVDSolve(svd);
 97:   SVDErrorView(svd,SVD_ERROR_NORM,NULL);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:                        Solve second problem
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

103:   SVDSetOperators(svd,A2,B);
104:   SVDSolve(svd);
105:   SVDErrorView(svd,SVD_ERROR_NORM,NULL);

107:   /* Free work space */
108:   SVDDestroy(&svd);
109:   MatDestroy(&A1);
110:   MatDestroy(&A2);
111:   MatDestroy(&B);
112:   SlepcFinalize();
113:   return 0;
114: }

116: /*TEST

118:    testset:
119:       args: -svd_nsv 3
120:       requires: !single
121:       output_file: output/test16_1.out
122:       test:
123:          suffix: 1_lapack
124:          args: -svd_type lapack
125:       test:
126:          suffix: 1_cross
127:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
128:       test:
129:          suffix: 1_cyclic
130:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
131:       test:
132:          suffix: 1_trlanczos
133:          args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}}

135: TEST*/