Actual source code: test5.c
slepc-3.18.3 2023-03-24
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 SVD view and monitor functionality.\n\n";
13: #include <slepcsvd.h>
15: int main(int argc,char **argv)
16: {
17: Mat A;
18: SVD svd;
19: PetscReal *sigma,sval;
20: PetscInt n=6,Istart,Iend,i,nconv;
21: PetscBool checkfile;
22: char filename[PETSC_MAX_PATH_LEN];
23: PetscViewer viewer;
26: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscPrintf(PETSC_COMM_WORLD,"\nSVD of diagonal matrix, n=%" PetscInt_FMT "\n\n",n);
30: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: Generate the matrix
32: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33: MatCreate(PETSC_COMM_WORLD,&A);
34: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
35: MatSetFromOptions(A);
36: MatSetUp(A);
37: MatGetOwnershipRange(A,&Istart,&Iend);
38: for (i=Istart;i<Iend;i++) MatSetValue(A,i,i,i+1,INSERT_VALUES);
39: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
40: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Create the SVD solver
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45: SVDCreate(PETSC_COMM_WORLD,&svd);
46: PetscObjectSetName((PetscObject)svd,"svd");
47: SVDSetOperators(svd,A,NULL);
48: SVDSetFromOptions(svd);
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Compute the singular triplets and display solution
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: SVDSolve(svd);
54: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Check file containing the singular values
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
60: if (checkfile) {
61: SVDGetConverged(svd,&nconv);
62: PetscMalloc1(nconv,&sigma);
63: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
64: PetscViewerBinaryRead(viewer,sigma,nconv,NULL,PETSC_REAL);
65: PetscViewerDestroy(&viewer);
66: for (i=0;i<nconv;i++) {
67: SVDGetSingularTriplet(svd,i,&sval,NULL,NULL);
69: }
70: PetscFree(sigma);
71: }
73: SVDDestroy(&svd);
74: MatDestroy(&A);
75: SlepcFinalize();
76: return 0;
77: }
79: /*TEST
81: test:
82: suffix: 1
83: args: -svd_error_relative ::ascii_info_detail -svd_view_values -svd_monitor_conv -svd_error_absolute ::ascii_matlab -svd_monitor_all -svd_converged_reason -svd_view
84: filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/1.999999/2.000000/" | sed -e "s/2.000001/2.000000/" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"
85: requires: !single
87: test:
88: suffix: 2
89: args: -svd_nsv 4 -svd_view_values binary:myvalues.bin -checkfile myvalues.bin
90: requires: double
92: test:
93: suffix: 3
94: args: -svd_type trlanczos -svd_error_relative ::ascii_info_detail -svd_view_values -svd_monitor_conv -svd_error_absolute ::ascii_matlab -svd_monitor_all -svd_converged_reason -svd_view
95: filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/1.999999/2.000000/" | sed -e "s/2.000001/2.000000/" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"
96: requires: !single
98: test:
99: suffix: 4
100: args: -svd_monitor draw::draw_lg -svd_monitor_all draw::draw_lg -svd_view_values draw -draw_save mysingu.ppm -draw_virtual
101: requires: x !single
103: TEST*/