Actual source code: test18.c

slepc-3.17.2 2022-08-09
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[] = "Test GSVD with user-provided initial vectors.\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 a GSVD problem for the bidiagonal matrices

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

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            A,B;
 34:   SVD            svd;
 35:   Vec            v0,w0;           /* initial vectors */
 36:   VecType        vtype;
 37:   PetscInt       m=22,n=20,p=22,Istart,Iend,i,col[2];
 38:   PetscScalar    valsa[] = { 1, 2 }, valsb[] = { -1, 2 };

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

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

 77:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 78:          Create the singular value solver, set options and solve
 79:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 81:   SVDCreate(PETSC_COMM_WORLD,&svd);
 82:   SVDSetOperators(svd,A,B);
 83:   SVDSetFromOptions(svd);

 85:   /*
 86:      Set the initial vectors. This is optional, if not done the initial
 87:      vectors are set to random values
 88:   */
 89:   MatCreateVecs(A,&v0,NULL);        /* right initial vector, length n */
 90:   VecCreate(PETSC_COMM_WORLD,&w0);  /* left initial vector, length m+p */
 91:   VecSetSizes(w0,PETSC_DECIDE,m+p);
 92:   VecGetType(v0,&vtype);
 93:   VecSetType(w0,vtype);
 94:   VecSet(v0,1.0);
 95:   VecSet(w0,1.0);
 96:   SVDSetInitialSpaces(svd,1,&v0,1,&w0);

 98:   /*
 99:      Compute solution
100:   */
101:   SVDSolve(svd);
102:   SVDErrorView(svd,SVD_ERROR_NORM,NULL);

104:   /* Free work space */
105:   VecDestroy(&v0);
106:   VecDestroy(&w0);
107:   SVDDestroy(&svd);
108:   MatDestroy(&A);
109:   MatDestroy(&B);
110:   SlepcFinalize();
111:   return 0;
112: }

114: /*TEST

116:    testset:
117:       args: -svd_nsv 3
118:       requires: !single
119:       output_file: output/test18_1.out
120:       test:
121:          suffix: 1
122:          args: -svd_type {{lapack cross cyclic}}
123:       test:
124:          suffix: 1_trlanczos
125:          args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}}
126:          requires: !__float128

128: TEST*/