Actual source code: zrkf.c
1: #include <petsc/private/ftnimpl.h>
2: #include <petsc/private/tsimpl.h>
3: #include <../src/ts/impls/explicit/rk/rk.h>
5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
6: #define tsrkgettableau_ TSRKGETTABLEAU
7: #define tsrkrestoretableau_ TSRKRESTORETABLEAU
8: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9: #define tsrkgettableau_ tsrkgettableau
10: #define tsrkrestoretableau_ tsrkrestoretableau
11: #endif
13: PETSC_EXTERN void tsrkgettableau_(TS *ts, PetscInt *s, F90Array1d *A, F90Array1d *b, F90Array1d *c, F90Array1d *bembed, PetscInt *p, F90Array1d *binterp, PetscBool *FSAL, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(p1) PETSC_F90_2PTR_PROTO(p2) PETSC_F90_2PTR_PROTO(p3) PETSC_F90_2PTR_PROTO(p4) PETSC_F90_2PTR_PROTO(p5))
14: {
15: TS_RK *rk = (TS_RK *)(*ts)->data;
16: RKTableau tab = rk->tableau;
18: CHKFORTRANNULLINTEGER(s);
19: CHKFORTRANNULLINTEGER(p);
20: CHKFORTRANNULLBOOL(FSAL);
21: if (s) *s = tab->s;
22: if (!FORTRANNULLSCALARPOINTER(A)) {
23: *ierr = F90Array1dCreate(tab->A, MPIU_REAL, 1, tab->s * tab->s, A PETSC_F90_2PTR_PARAM(p1));
24: if (*ierr) return;
25: }
26: if (!FORTRANNULLSCALARPOINTER(b)) {
27: *ierr = F90Array1dCreate(tab->b, MPIU_REAL, 1, tab->s, b PETSC_F90_2PTR_PARAM(p2));
28: if (*ierr) return;
29: }
30: if (!FORTRANNULLSCALARPOINTER(c)) {
31: *ierr = F90Array1dCreate(tab->c, MPIU_REAL, 1, tab->s, c PETSC_F90_2PTR_PARAM(p3));
32: if (*ierr) return;
33: }
34: if (!FORTRANNULLSCALARPOINTER(bembed) && tab->bembed) {
35: *ierr = F90Array1dCreate(tab->bembed, MPIU_REAL, 1, tab->s, bembed PETSC_F90_2PTR_PARAM(p4));
36: if (*ierr) return;
37: }
38: if (p) *p = tab->p;
39: if (!FORTRANNULLSCALARPOINTER(binterp)) {
40: *ierr = F90Array1dCreate(tab->binterp, MPIU_REAL, 1, tab->s * tab->p, binterp PETSC_F90_2PTR_PARAM(p5));
41: if (*ierr) return;
42: }
43: if (FSAL) *FSAL = tab->FSAL;
44: *ierr = PETSC_SUCCESS;
45: }
47: PETSC_EXTERN void tsrkrestoretableau_(TS *ts, PetscInt *s, F90Array1d *A, F90Array1d *b, F90Array1d *c, F90Array1d *bembed, PetscInt *p, F90Array1d *binterp, PetscBool *FSAL, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(p1) PETSC_F90_2PTR_PROTO(p2) PETSC_F90_2PTR_PROTO(p3) PETSC_F90_2PTR_PROTO(p4) PETSC_F90_2PTR_PROTO(p5))
48: {
49: if (!FORTRANNULLSCALARPOINTER(A)) {
50: *ierr = F90Array1dDestroy(A, MPIU_REAL PETSC_F90_2PTR_PARAM(p1));
51: if (*ierr) return;
52: }
53: if (!FORTRANNULLSCALARPOINTER(b)) {
54: *ierr = F90Array1dDestroy(b, MPIU_REAL PETSC_F90_2PTR_PARAM(p2));
55: if (*ierr) return;
56: }
57: if (!FORTRANNULLSCALARPOINTER(c)) {
58: *ierr = F90Array1dDestroy(c, MPIU_REAL PETSC_F90_2PTR_PARAM(p3));
59: if (*ierr) return;
60: }
61: if (!FORTRANNULLSCALARPOINTER(bembed)) {
62: *ierr = F90Array1dDestroy(bembed, MPIU_REAL PETSC_F90_2PTR_PARAM(p4));
63: if (*ierr) return;
64: }
65: if (!FORTRANNULLSCALARPOINTER(binterp)) {
66: *ierr = F90Array1dDestroy(binterp, MPIU_REAL PETSC_F90_2PTR_PARAM(p5));
67: if (*ierr) return;
68: }
69: *ierr = PETSC_SUCCESS;
70: }