Actual source code: test4.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 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" };

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

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

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

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

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

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

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

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

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

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

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

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

121: /*TEST

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

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

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

183: TEST*/