Actual source code: test4.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 an SVD problem with more columns than rows.\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 a rectangular bidiagonal matrix

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

 30: int main(int argc,char **argv)
 31: {
 32:   Mat                  A,B;
 33:   SVD                  svd;
 34:   SVDConv              conv;
 35:   SVDStop              stop;
 36:   SVDWhich             which;
 37:   SVDProblemType       ptype;
 38:   SVDConvergedReason   reason;
 39:   PetscInt             m=20,n,Istart,Iend,i,col[2],its;
 40:   PetscScalar          value[] = { 1, 2 };
 41:   PetscBool            flg,tmode;
 42:   PetscErrorCode       ierr;
 43:   PetscViewerAndFormat *vf;
 44:   const char           *ctest[] = { "absolute", "relative to the singular value", "user-defined" };
 45:   const char           *stest[] = { "basic", "user-defined" };

 47:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 49:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 50:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 51:   if (!flg) n=m+2;
 52:   PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:         Generate the matrix
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   MatCreate(PETSC_COMM_WORLD,&A);
 59:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 60:   MatSetFromOptions(A);
 61:   MatSetUp(A);

 63:   MatGetOwnershipRange(A,&Istart,&Iend);
 64:   for (i=Istart;i<Iend;i++) {
 65:     col[0]=i; col[1]=i+1;
 66:     if (i<n-1) {
 67:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 68:     } else if (i==n-1) {
 69:       MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
 70:     }
 71:   }

 73:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 74:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 76:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 77:         Compute singular values
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 80:   SVDCreate(PETSC_COMM_WORLD,&svd);
 81:   SVDSetOperators(svd,A,NULL);

 83:   /* test some interface functions */
 84:   SVDGetOperators(svd,&B,NULL);
 85:   MatView(B,PETSC_VIEWER_STDOUT_WORLD);
 86:   SVDSetConvergenceTest(svd,SVD_CONV_ABS);
 87:   SVDSetStoppingTest(svd,SVD_STOP_BASIC);
 88:   /* test monitors */
 89:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
 90:   SVDMonitorSet(svd,(PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*))SVDMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
 91:   /* SVDMonitorCancel(svd); */
 92:   SVDSetFromOptions(svd);

 94:   /* query properties and print them */
 95:   SVDGetProblemType(svd,&ptype);
 96:   PetscPrintf(PETSC_COMM_WORLD," Problem type = %d",(int)ptype);
 97:   SVDIsGeneralized(svd,&flg);
 98:   if (flg) { PetscPrintf(PETSC_COMM_WORLD," generalized"); }
 99:   SVDGetImplicitTranspose(svd,&tmode);
100:   PetscPrintf(PETSC_COMM_WORLD,"\n Transpose mode is %s\n",tmode?"implicit":"explicit");
101:   SVDGetConvergenceTest(svd,&conv);
102:   PetscPrintf(PETSC_COMM_WORLD," Convergence test is %s\n",ctest[conv]);
103:   SVDGetStoppingTest(svd,&stop);
104:   PetscPrintf(PETSC_COMM_WORLD," Stopping test is %s\n",stest[stop]);
105:   SVDGetWhichSingularTriplets(svd,&which);
106:   PetscPrintf(PETSC_COMM_WORLD," Which = %s\n",which?"smallest":"largest");

108:   /* call the solver */
109:   SVDSolve(svd);
110:   SVDGetConvergedReason(svd,&reason);
111:   SVDGetIterationNumber(svd,&its);
112:   PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
113:   /* PetscPrintf(PETSC_COMM_WORLD," its = %D\n",its); */

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:                     Display solution and clean up
117:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118:   SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
119:   SVDDestroy(&svd);
120:   MatDestroy(&A);
121:   SlepcFinalize();
122:   return ierr;
123: }

125: /*TEST

127:    testset:
128:       args: -svd_monitor_cancel
129:       filter: grep -v "Transpose mode"
130:       output_file: output/test4_1.out
131:       test:
132:          suffix: 1_lanczos
133:          args: -svd_type lanczos
134:       test:
135:          suffix: 1_randomized
136:          args: -svd_type randomized
137:       test:
138:          suffix: 1_trlanczos
139:          args: -svd_type trlanczos -svd_ncv 12
140:       test:
141:          suffix: 1_cross
142:          args: -svd_type cross
143:       test:
144:          suffix: 1_cross_exp
145:          args: -svd_type cross -svd_cross_explicitmatrix
146:       test:
147:          suffix: 1_cross_exp_imp
148:          args: -svd_type cross -svd_cross_explicitmatrix -svd_implicittranspose
149:          requires: !complex
150:       test:
151:          suffix: 1_cyclic
152:          args: -svd_type cyclic
153:       test:
154:          suffix: 1_cyclic_imp
155:          args: -svd_type cyclic -svd_implicittranspose
156:       test:
157:          suffix: 1_cyclic_exp
158:          args: -svd_type cyclic -svd_cyclic_explicitmatrix
159:       test:
160:          suffix: 1_lapack
161:          args: -svd_type lapack
162:       test:
163:          suffix: 1_scalapack
164:          args: -svd_type scalapack
165:          requires: scalapack

167:    testset:
168:       args: -svd_monitor_cancel  -mat_type aijcusparse
169:       requires: cuda !single
170:       filter: grep -v "Transpose mode" | sed -e "s/seqaijcusparse/seqaij/"
171:       output_file: output/test4_1.out
172:       test:
173:          suffix: 2_cuda_lanczos
174:          args: -svd_type lanczos
175:       test:
176:          suffix: 2_cuda_trlanczos
177:          args: -svd_type trlanczos -svd_ncv 12
178:       test:
179:          suffix: 2_cuda_cross
180:          args: -svd_type cross

182:    test:
183:       suffix: 3
184:       nsize: 2
185:       args: -svd_type trlanczos -svd_ncv 14 -svd_monitor_cancel -ds_parallel synchronized

187: TEST*/