Actual source code: zdtdsf90.c
1: #include <petsc/private/fortranimpl.h>
2: #include <petscds.h>
3: #include <petsc/private/f90impl.h>
5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
6: #define petscdsgettabulation_ PETSCDSGETTABULATION
7: #define petscdsrestoretabulation_ PETSCDSRESTORETABULATION
8: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9: #define petscdsgettabulation_ petscdsgettabulation
10: #define petscdsrestoretabulation_ petscdsrestoretabulation
11: #endif
13: PETSC_EXTERN void petscdsgettabulation_(PetscDS *prob, PetscInt *f, F90Array1d *ptrB, F90Array1d *ptrD, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrb) PETSC_F90_2PTR_PROTO(ptrd))
14: {
15: PetscFE fe;
16: PetscQuadrature q;
17: PetscInt dim, Nb, Nc, Nq;
18: PetscTabulation *T;
20: *ierr = PetscDSGetSpatialDimension(*prob, &dim);
21: if (*ierr) return;
22: *ierr = PetscDSGetDiscretization(*prob, *f, (PetscObject *)&fe);
23: if (*ierr) return;
24: *ierr = PetscFEGetDimension(fe, &Nb);
25: if (*ierr) return;
26: *ierr = PetscFEGetNumComponents(fe, &Nc);
27: if (*ierr) return;
28: *ierr = PetscFEGetQuadrature(fe, &q);
29: if (*ierr) return;
30: *ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);
31: if (*ierr) return;
32: *ierr = PetscDSGetTabulation(*prob, &T);
33: if (*ierr) return;
34: *ierr = F90Array1dCreate((void *)T[*f]->T[0], MPIU_REAL, 1, Nq * Nb * Nc, ptrB PETSC_F90_2PTR_PARAM(ptrb));
35: if (*ierr) return;
36: *ierr = F90Array1dCreate((void *)T[*f]->T[1], MPIU_REAL, 1, Nq * Nb * Nc * dim, ptrD PETSC_F90_2PTR_PARAM(ptrd));
37: }
39: PETSC_EXTERN void petscdsrestoretabulation_(PetscDS *prob, PetscInt *f, F90Array1d *ptrB, F90Array1d *ptrD, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrb) PETSC_F90_2PTR_PROTO(ptrd))
40: {
41: *ierr = F90Array1dDestroy(ptrB, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrb));
42: if (*ierr) return;
43: *ierr = F90Array1dDestroy(ptrD, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrd));
44: }