Actual source code: test2.c
slepc-3.23.3 2025-09-08
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 with different builds with a matrix loaded from a file"
12: " (matrices available in PETSc's distribution).\n\n";
14: #include <slepcsvd.h>
16: int main(int argc,char **argv)
17: {
18: Mat A; /* operator matrix */
19: SVD svd; /* singular value problem solver context */
20: char filename[PETSC_MAX_PATH_LEN],path[PETSC_MAX_PATH_LEN];
21: const char *prefix,*scalar,*ints,*floats;
22: PetscReal tol=PETSC_SMALL;
23: PetscViewer viewer;
24: PetscBool flg;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Load the matrix for which the SVD must be computed
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32: #if defined(PETSC_USE_COMPLEX)
33: prefix = "nh";
34: scalar = "complex";
35: #else
36: prefix = "ns";
37: scalar = "real";
38: #endif
39: #if defined(PETSC_USE_64BIT_INDICES)
40: ints = "int64";
41: #else
42: ints = "int32";
43: #endif
44: #if defined(PETSC_USE_REAL_DOUBLE)
45: floats = "float64";
46: #elif defined(PETSC_USE_REAL_SINGLE)
47: floats = "float32";
48: #endif
50: PetscCall(PetscOptionsGetString(NULL,NULL,"-path",path,sizeof(path),&flg));
51: if (flg) PetscCall(PetscSNPrintf(filename,sizeof(filename),"%s/%s-%s-%s-%s",path,prefix,scalar,ints,floats));
52: else PetscCall(PetscSNPrintf(filename,sizeof(filename),"%s/share/petsc/datafiles/matrices/%s-%s-%s-%s",PETSC_DIR,prefix,scalar,ints,floats));
53: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nReading matrix from binary file...\n\n"));
54: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
55: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
56: PetscCall(MatSetFromOptions(A));
57: PetscCall(MatLoad(A,viewer));
58: PetscCall(PetscViewerDestroy(&viewer));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the SVD solver
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
64: PetscCall(SVDSetOperators(svd,A,NULL));
65: PetscCall(SVDSetTolerances(svd,tol,PETSC_CURRENT));
66: PetscCall(SVDSetFromOptions(svd));
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Compute the singular triplets and display solution
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PetscCall(SVDSolve(svd));
72: PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
73: PetscCall(SVDDestroy(&svd));
74: PetscCall(MatDestroy(&A));
75: PetscCall(SlepcFinalize());
76: return 0;
77: }
79: /*TEST
81: build:
82: requires: !__float128
84: test:
85: args: -svd_nsv 7 -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -path ${DATAFILESPATH}/matrices/petsc-small
86: requires: !single datafilespath
88: testset:
89: args: -svd_nsv 7 -svd_mpd 11 -svd_type primme -path ${DATAFILESPATH}/matrices/petsc-small
90: requires: primme !single datafilespath
91: output_file: output/test2_1.out
92: test:
93: suffix: 1_primme
94: test:
95: suffix: 1_primme_args
96: args: -svd_primme_blocksize 2 -svd_primme_method hybrid
98: TEST*/