Line data Source code
1 : /* 2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 3 : SLEPc - Scalable Library for Eigenvalue Problem Computations 4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain 5 : 6 : This file is part of SLEPc. 7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE). 8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 9 : */ 10 : 11 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/ 12 : 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 31 : 32 : /*@C 33 : SVDRegisterAll - Registers all the singular value solvers in the SVD package. 34 : 35 : Not Collective 36 : 37 : Level: advanced 38 : 39 : .seealso: SVDRegister() 40 : @*/ 41 471 : PetscErrorCode SVDRegisterAll(void) 42 : { 43 471 : PetscFunctionBegin; 44 471 : if (SVDRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 45 240 : SVDRegisterAllCalled = PETSC_TRUE; 46 240 : PetscCall(SVDRegister(SVDCROSS,SVDCreate_Cross)); 47 240 : PetscCall(SVDRegister(SVDCYCLIC,SVDCreate_Cyclic)); 48 240 : PetscCall(SVDRegister(SVDLAPACK,SVDCreate_LAPACK)); 49 240 : PetscCall(SVDRegister(SVDLANCZOS,SVDCreate_Lanczos)); 50 240 : PetscCall(SVDRegister(SVDTRLANCZOS,SVDCreate_TRLanczos)); 51 240 : 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 240 : PetscCall(SVDRegister(SVDPRIMME,SVDCreate_PRIMME)); 63 : #endif 64 240 : PetscFunctionReturn(PETSC_SUCCESS); 65 : } 66 : 67 : /*@C 68 : SVDMonitorRegisterAll - Registers all the monitors in the SVD package. 69 : 70 : Not Collective 71 : 72 : Level: advanced 73 : 74 : .seealso: SVDMonitorRegister() 75 : @*/ 76 240 : PetscErrorCode SVDMonitorRegisterAll(void) 77 : { 78 240 : PetscFunctionBegin; 79 240 : if (SVDMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 80 240 : SVDMonitorRegisterAllCalled = PETSC_TRUE; 81 : 82 240 : PetscCall(SVDMonitorRegister("first_approximation",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorFirst,NULL,NULL)); 83 240 : PetscCall(SVDMonitorRegister("first_approximation",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorFirstDrawLG,SVDMonitorFirstDrawLGCreate,NULL)); 84 240 : PetscCall(SVDMonitorRegister("all_approximations",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorAll,NULL,NULL)); 85 240 : PetscCall(SVDMonitorRegister("all_approximations",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorAllDrawLG,SVDMonitorAllDrawLGCreate,NULL)); 86 240 : PetscCall(SVDMonitorRegister("convergence_history",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorConverged,SVDMonitorConvergedCreate,SVDMonitorConvergedDestroy)); 87 240 : PetscCall(SVDMonitorRegister("convergence_history",PETSCVIEWERDRAW,PETSC_VIEWER_DRAW_LG,SVDMonitorConvergedDrawLG,SVDMonitorConvergedDrawLGCreate,SVDMonitorConvergedDestroy)); 88 240 : PetscCall(SVDMonitorRegister("conditioning",PETSCVIEWERASCII,PETSC_VIEWER_DEFAULT,SVDMonitorConditioning,NULL,NULL)); 89 240 : PetscFunctionReturn(PETSC_SUCCESS); 90 : }