Actual source code: test4.c

slepc-3.18.0 2022-10-01
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 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:   PetscViewerAndFormat *vf;
 43:   const char           *ctest[] = { "absolute", "relative to the singular value", "user-defined" };
 44:   const char           *stest[] = { "basic", "user-defined" };

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

 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=%" PetscInt_FMT " n=%" PetscInt_FMT "\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) MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 67:     else if (i==n-1) MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
 68:   }

 70:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 73:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74:         Compute singular values
 75:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 77:   SVDCreate(PETSC_COMM_WORLD,&svd);
 78:   SVDSetOperators(svd,A,NULL);

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

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

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

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:                     Display solution and clean up
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115:   SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
116:   SVDDestroy(&svd);
117:   MatDestroy(&A);
118:   SlepcFinalize();
119:   return 0;
120: }

122: /*TEST

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

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

179:    test:
180:       suffix: 3
181:       nsize: 2
182:       args: -svd_type trlanczos -svd_ncv 14 -svd_monitor_cancel -ds_parallel synchronized

184: TEST*/