Actual source code: svdregis.c
slepc-3.22.2 2024-12-02
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: #include <slepc/private/svdimpl.h>
13: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD);
14: SLEPC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD);
15: SLEPC_EXTERN PetscErrorCode SVDCreate_LAPACK(SVD);
16: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD);
17: SLEPC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD);
18: SLEPC_EXTERN PetscErrorCode SVDCreate_Randomized(SVD);
19: #if defined(SLEPC_HAVE_SCALAPACK)
20: SLEPC_EXTERN PetscErrorCode SVDCreate_ScaLAPACK(SVD);
21: #endif
22: #if defined(SLEPC_HAVE_KSVD)
23: SLEPC_EXTERN PetscErrorCode SVDCreate_KSVD(SVD);
24: #endif
25: #if defined(SLEPC_HAVE_ELEMENTAL)
26: SLEPC_EXTERN PetscErrorCode SVDCreate_Elemental(SVD);
27: #endif
28: #if defined(SLEPC_HAVE_PRIMME)
29: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD);
30: #endif
32: /*@C
33: SVDRegisterAll - Registers all the singular value solvers in the SVD package.
35: Not Collective
37: Level: advanced
39: .seealso: SVDRegister()
40: @*/
41: PetscErrorCode SVDRegisterAll(void)
42: {
43: PetscFunctionBegin;
44: if (SVDRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
45: SVDRegisterAllCalled = PETSC_TRUE;
46: PetscCall(SVDRegister(SVDCROSS,SVDCreate_Cross));
47: PetscCall(SVDRegister(SVDCYCLIC,SVDCreate_Cyclic));
48: PetscCall(SVDRegister(SVDLAPACK,SVDCreate_LAPACK));
49: PetscCall(SVDRegister(SVDLANCZOS,SVDCreate_Lanczos));
50: PetscCall(SVDRegister(SVDTRLANCZOS,SVDCreate_TRLanczos));
51: PetscCall(SVDRegister(SVDRANDOMIZED,SVDCreate_Randomized));
52: #if defined(SLEPC_HAVE_SCALAPACK)
53: PetscCall(SVDRegister(SVDSCALAPACK,SVDCreate_ScaLAPACK));
54: #endif
55: #if defined(SLEPC_HAVE_KSVD)
56: PetscCall(SVDRegister(SVDKSVD,SVDCreate_KSVD));
57: #endif
58: #if defined(SLEPC_HAVE_ELEMENTAL)
59: PetscCall(SVDRegister(SVDELEMENTAL,SVDCreate_Elemental));
60: #endif
61: #if defined(SLEPC_HAVE_PRIMME)
62: PetscCall(SVDRegister(SVDPRIMME,SVDCreate_PRIMME));
63: #endif
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: /*@C
68: SVDMonitorRegisterAll - Registers all the monitors in the SVD package.
70: Not Collective
72: Level: advanced
74: .seealso: SVDMonitorRegister()
75: @*/
76: PetscErrorCode SVDMonitorRegisterAll(void)
77: {
78: PetscFunctionBegin;
79: if (SVDMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
80: SVDMonitorRegisterAllCalled = PETSC_TRUE;
82: PetscCall(SVDMonitorRegister("first_approximation",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorFirst,NULL,NULL));
83: PetscCall(SVDMonitorRegister("first_approximation",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorFirstDrawLG,SVDMonitorFirstDrawLGCreate,NULL));
84: PetscCall(SVDMonitorRegister("all_approximations",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorAll,NULL,NULL));
85: PetscCall(SVDMonitorRegister("all_approximations",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorAllDrawLG,SVDMonitorAllDrawLGCreate,NULL));
86: PetscCall(SVDMonitorRegister("convergence_history",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorConverged,SVDMonitorConvergedCreate,SVDMonitorConvergedDestroy));
87: PetscCall(SVDMonitorRegister("convergence_history",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorConvergedDrawLG,SVDMonitorConvergedDrawLGCreate,SVDMonitorConvergedDestroy));
88: PetscCall(SVDMonitorRegister("conditioning",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorConditioning,NULL,NULL));
89: PetscFunctionReturn(PETSC_SUCCESS);
90: }