Actual source code: test5.c
slepc-3.15.1 2021-05-28
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 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);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscPrintf(PETSC_COMM_WORLD,"\nSVD of diagonal matrix, n=%D\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++) {
39: MatSetValue(A,i,i,i+1,INSERT_VALUES);
40: }
41: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
42: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
44: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: Create the SVD solver
46: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: SVDCreate(PETSC_COMM_WORLD,&svd);
48: PetscObjectSetName((PetscObject)svd,"svd");
49: SVDSetOperators(svd,A,NULL);
50: SVDSetFromOptions(svd);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Compute the singular triplets and display solution
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: SVDSolve(svd);
56: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Check file containing the singular values
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
62: if (checkfile) {
63: SVDGetConverged(svd,&nconv);
64: PetscMalloc1(nconv,&sigma);
65: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
66: PetscViewerBinaryRead(viewer,sigma,nconv,NULL,PETSC_REAL);
67: PetscViewerDestroy(&viewer);
68: for (i=0;i<nconv;i++) {
69: SVDGetSingularTriplet(svd,i,&sval,NULL,NULL);
70: if (sval!=sigma[i]) SETERRQ(PETSC_COMM_WORLD,1,"Singular values in the file do not match");
71: }
72: PetscFree(sigma);
73: }
75: SVDDestroy(&svd);
76: MatDestroy(&A);
77: SlepcFinalize();
78: return ierr;
79: }
81: /*TEST
83: test:
84: suffix: 1
85: 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
86: 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"
87: requires: !single
89: test:
90: suffix: 2
91: args: -svd_nsv 4 -svd_view_values binary:myvalues.bin -checkfile myvalues.bin
92: requires: double
94: test:
95: suffix: 3
96: 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
97: 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"
98: requires: !single
100: test:
101: suffix: 4
102: args: -svd_monitor draw::draw_lg -svd_monitor_all draw::draw_lg -svd_view_values draw -draw_save mysingu.ppm -draw_virtual
103: requires: !single
105: TEST*/