Actual source code: ex48.c

slepc-3.18.3 2023-03-24
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[] = "Solves a GSVD problem with matrices loaded from a file.\n"
 12:   "The command line options are:\n"
 13:   "  -f1 <filename>, PETSc binary file containing matrix A.\n"
 14:   "  -f2 <filename>, PETSc binary file containing matrix B (optional). Instead of"
 15:   "     a file it is possible to specify one of 'identity', 'bidiagonal' or 'tridiagonal'"
 16:   "  -p <p>, in case B is not taken from a file.\n\n";

 18: #include <slepcsvd.h>

 20: int main(int argc,char **argv)
 21: {
 22:   Mat            A,B;             /* matrices */
 23:   SVD            svd;             /* singular value problem solver context */
 24:   PetscInt       i,m,n,p,Istart,Iend,col[3];
 25:   PetscScalar    vals[3];
 26:   char           filename[PETSC_MAX_PATH_LEN];
 27:   PetscViewer    viewer;
 28:   PetscBool      flg,terse;

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

 33:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34:         Load matrices that define the generalized singular value problem
 35:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 37:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value problem stored in file.\n\n");
 38:   PetscOptionsGetString(NULL,NULL,"-f1",filename,sizeof(filename),&flg);

 41: #if defined(PETSC_USE_COMPLEX)
 42:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
 43: #else
 44:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
 45: #endif
 46:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 47:   MatCreate(PETSC_COMM_WORLD,&A);
 48:   MatSetFromOptions(A);
 49:   MatLoad(A,viewer);
 50:   PetscViewerDestroy(&viewer);

 52:   MatGetSize(A,&m,&n);

 54:   PetscOptionsGetString(NULL,NULL,"-f2",filename,sizeof(filename),&flg);
 56:   PetscStrcmp(filename,"identity",&flg);
 57:   if (flg) {
 58:     p = n;
 59:     PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
 60:     PetscPrintf(PETSC_COMM_WORLD," Using B=I with p=%" PetscInt_FMT "\n\n",p);
 61:     MatCreate(PETSC_COMM_WORLD,&B);
 62:     MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 63:     MatSetFromOptions(B);
 64:     MatSetUp(B);
 65:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 66:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 67:     MatShift(B,1.0);
 68:   } else {
 69:     PetscStrcmp(filename,"bidiagonal",&flg);
 70:     if (flg) {
 71:       p = n+1;
 72:       PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
 73:       vals[0]=-1; vals[1]=1;
 74:       PetscPrintf(PETSC_COMM_WORLD," Using B=bidiag(1,-1) with p=%" PetscInt_FMT "\n\n",p);
 75:       MatCreate(PETSC_COMM_WORLD,&B);
 76:       MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 77:       MatSetFromOptions(B);
 78:       MatSetUp(B);
 79:       MatGetOwnershipRange(B,&Istart,&Iend);
 80:       for (i=Istart;i<Iend;i++) {
 81:         col[0]=i-1; col[1]=i;
 82:         if (i==0) MatSetValue(B,i,col[1],vals[1],INSERT_VALUES);
 83:         else if (i<n) MatSetValues(B,1,&i,2,col,vals,INSERT_VALUES);
 84:         else if (i==n) MatSetValue(B,i,col[0],vals[0],INSERT_VALUES);
 85:       }
 86:       MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 87:       MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 88:     } else {
 89:       PetscStrcmp(filename,"tridiagonal",&flg);
 90:       if (flg) {
 91:         p = n-2;
 92:         PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
 93:         vals[0]=-1; vals[1]=2; vals[2]=-1;
 94:         PetscPrintf(PETSC_COMM_WORLD," Using B=tridiag(-1,2,-1) with p=%" PetscInt_FMT "\n\n",p);
 95:         MatCreate(PETSC_COMM_WORLD,&B);
 96:         MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 97:         MatSetFromOptions(B);
 98:         MatSetUp(B);
 99:         MatGetOwnershipRange(B,&Istart,&Iend);
100:         for (i=Istart;i<Iend;i++) {
101:           col[0]=i; col[1]=i+1; col[2]=i+2;
102:           MatSetValues(B,1,&i,3,col,vals,INSERT_VALUES);
103:         }
104:         MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
105:         MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
106:       } else {  /* load file */
107:         PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
108:         MatCreate(PETSC_COMM_WORLD,&B);
109:         MatSetFromOptions(B);
110:         MatLoad(B,viewer);
111:         PetscViewerDestroy(&viewer);
112:       }
113:     }
114:   }

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:            Create the singular value solver and set various options
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   /*
121:      Create singular value solver context
122:   */
123:   SVDCreate(PETSC_COMM_WORLD,&svd);

125:   /*
126:      Set operators of GSVD problem
127:   */
128:   SVDSetOperators(svd,A,B);
129:   SVDSetProblemType(svd,SVD_GENERALIZED);

131:   /*
132:      Set solver parameters at runtime
133:   */
134:   SVDSetFromOptions(svd);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                       Solve the problem and print solution
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   SVDSolve(svd);

142:   /* show detailed info unless -terse option is given by user */
143:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
144:   if (terse) SVDErrorView(svd,SVD_ERROR_NORM,NULL);
145:   else {
146:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
147:     SVDConvergedReasonView(svd,PETSC_VIEWER_STDOUT_WORLD);
148:     SVDErrorView(svd,SVD_ERROR_NORM,PETSC_VIEWER_STDOUT_WORLD);
149:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
150:   }
151:   SVDDestroy(&svd);
152:   MatDestroy(&A);
153:   MatDestroy(&B);
154:   SlepcFinalize();
155:   return 0;
156: }
157: /*TEST

159:    testset:
160:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
161:       args: -f1 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62a.petsc -f2 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62b.petsc -svd_nsv 3 -terse
162:       output_file: output/ex48_1.out
163:       test:
164:          suffix: 1
165:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}} -svd_trlanczos_scale 1e5 -svd_trlanczos_ksp_rtol 1e-13
166:       test:
167:          suffix: 1_spqr
168:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix -svd_trlanczos_pc_type qr -svd_trlanczos_scale 1e5 -svd_trlanczos_oneside {{0 1}}
169:          requires: suitesparse
170:       test:
171:          suffix: 1_autoscale
172:          args: -svd_type trlanczos -svd_trlanczos_gbidiag {{lower upper}} -svd_trlanczos_scale -5 -svd_trlanczos_ksp_rtol 1e-14 -svd_trlanczos_oneside {{0 1}}
173:       test:
174:          suffix: 1_cross
175:          args: -svd_type cross -svd_cross_explicitmatrix
176:       test:
177:          suffix: 1_cyclic
178:          args: -svd_type cyclic -svd_cyclic_explicitmatrix

180:    testset:
181:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
182:       args: -f1 ${DATAFILESPATH}/matrices/complex/qc324.petsc -f2 bidiagonal -svd_nsv 3 -terse
183:       output_file: output/ex48_2.out
184:       filter: sed -e "s/30749/30748/"
185:       timeoutfactor: 2
186:       test:
187:          suffix: 2
188:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}} -svd_trlanczos_ksp_rtol 1e-10 -svd_trlanczos_scale 100
189:          requires: !defined(PETSCTEST_VALGRIND)
190:       test:
191:          suffix: 2_spqr
192:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix -svd_trlanczos_pc_type qr -svd_trlanczos_ksp_rtol 1e-10
193:          requires: suitesparse
194:       test:
195:          suffix: 2_cross
196:          args: -svd_type cross -svd_cross_explicitmatrix
197:       test:
198:          suffix: 2_cyclic
199:          args: -svd_type cyclic -svd_cyclic_explicitmatrix

201:    test:
202:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES) !defined(PETSCTEST_VALGRIND)
203:       args: -f1 ${DATAFILESPATH}/matrices/complex/qc324.petsc -f2 bidiagonal -p 320 -svd_nsv 3 -svd_type trlanczos -svd_trlanczos_ksp_rtol 1e-14 -svd_trlanczos_scale 100 -terse
204:       timeoutfactor: 2
205:       suffix: 3

207:    testset:
208:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
209:       args: -f1 ${SLEPC_DIR}/share/slepc/datafiles/matrices/rdb200.petsc -f2 identity -svd_nsv 3 -svd_ncv 24 -svd_smallest -terse
210:       output_file: output/ex48_4.out
211:       test:
212:          suffix: 4
213:          args: -svd_type trlanczos
214:       test:
215:          suffix: 4_spqr
216:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix -svd_trlanczos_pc_type qr
217:          requires: suitesparse
218:       test:
219:          suffix: 4_cross
220:          args: -svd_type cross -svd_cross_explicitmatrix
221:       test:
222:          suffix: 4_cross_implicit
223:          args: -svd_type cross -svd_cross_eps_type lobpcg -svd_cross_st_ksp_type cg -svd_cross_st_pc_type jacobi -svd_max_it 1000
224:       test:
225:          suffix: 4_cyclic
226:          args: -svd_type cyclic -svd_cyclic_explicitmatrix
227:       test:
228:          suffix: 4_hpddm
229:          nsize: 4
230:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix -svd_trlanczos_pc_type hpddm
231:          args: -prefix_push svd_trlanczos_pc_hpddm_ -levels_1_st_share_sub_ksp -levels_1_eps_nev 10 -levels_1_eps_threshold 0.005 -levels_1_pc_asm_type basic -define_subdomains -levels_1_pc_asm_sub_mat_type sbaij -coarse_pc_type cholesky -levels_1_sub_pc_type cholesky -prefix_pop
232:          requires: hpddm

234: TEST*/