xref: /petsc/src/ts/interface/tsregall.c (revision 2950f7e747d61362df6d6e1effb073d957ada526) !
1 #define PETSCTS_DLL
2 
3 #include "private/tsimpl.h"     /*I  "petscts.h"  I*/
4 EXTERN_C_BEGIN
5 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Euler(TS);
6 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_BEuler(TS);
7 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Pseudo(TS);
8 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS);
9 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_CN(TS);
10 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS);
11 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Alpha(TS);
12 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_GL(TS);
13 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_SSP(TS);
14 EXTERN PetscErrorCode PETSCTS_DLLEXPORT TSCreate_RK(TS);
15 EXTERN_C_END
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "TSRegisterAll"
19 /*@C
20   TSRegisterAll - Registers all of the timesteppers in the TS package.
21 
22   Not Collective
23 
24   Input parameter:
25 . path - The dynamic library path
26 
27   Level: advanced
28 
29 .keywords: TS, timestepper, register, all
30 .seealso: TSCreate(), TSRegister(), TSRegisterDestroy(), TSRegisterDynamic()
31 @*/
32 PetscErrorCode PETSCTS_DLLEXPORT TSRegisterAll(const char path[])
33 {
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   TSRegisterAllCalled = PETSC_TRUE;
38 
39   ierr = TSRegisterDynamic(TSEULER,           path, "TSCreate_Euler",    TSCreate_Euler);CHKERRQ(ierr);
40   ierr = TSRegisterDynamic(TSBEULER,          path, "TSCreate_BEuler",   TSCreate_BEuler);CHKERRQ(ierr);
41   ierr = TSRegisterDynamic(TSCN,              path, "TSCreate_CN",       TSCreate_CN);CHKERRQ(ierr);
42   ierr = TSRegisterDynamic(TSPSEUDO,          path, "TSCreate_Pseudo",   TSCreate_Pseudo);CHKERRQ(ierr);
43   ierr = TSRegisterDynamic(TSGL,              path, "TSCreate_GL",       TSCreate_GL);CHKERRQ(ierr);
44   ierr = TSRegisterDynamic(TSSSP,             path, "TSCreate_SSP",      TSCreate_SSP);CHKERRQ(ierr);
45   ierr = TSRegisterDynamic(TSTHETA,           path, "TSCreate_Theta",    TSCreate_Theta);CHKERRQ(ierr);
46   ierr = TSRegisterDynamic("alpha"/*TSALPHA*/,path, "TSCreate_Alpha",    TSCreate_Alpha);CHKERRQ(ierr);
47 #if defined(PETSC_HAVE_SUNDIALS)
48   ierr = TSRegisterDynamic(TSSUNDIALS,        path, "TSCreate_Sundials", TSCreate_Sundials);CHKERRQ(ierr);
49 #endif
50   ierr = TSRegisterDynamic(TSRK,              path, "TSCreate_RK",       TSCreate_RK);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54