Actual source code: zdtf90.c
1: #include <petsc/private/fortranimpl.h>
2: #include <petscdt.h>
3: #include <petsc/private/f90impl.h>
5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
6: #define petscquadraturegetdata_ PETSCQUADRATUREGETDATA
7: #define petscquadraturerestoredata_ PETSCQUADRATURERESTOREDATA
8: #define petscquadraturesetdata_ PETSCQUADRATURESETDATA
9: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
10: #define petscquadraturegetdata_ petscquadraturegetdata
11: #define petscquadraturerestoredata_ petscquadraturerestoredata
12: #define petscquadraturesetdata_ petscquadraturesetdata
13: #endif
15: PETSC_EXTERN void petscquadraturegetdata_(PetscQuadrature *q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, F90Array1d *ptrP, F90Array1d *ptrW, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrp) PETSC_F90_2PTR_PROTO(ptrw))
16: {
17: const PetscReal *points, *weights;
19: *ierr = PetscQuadratureGetData(*q, dim, Nc, npoints, &points, &weights);
20: if (*ierr) return;
21: *ierr = F90Array1dCreate((void *)points, MPIU_REAL, 1, (*npoints) * (*dim), ptrP PETSC_F90_2PTR_PARAM(ptrp));
22: if (*ierr) return;
23: *ierr = F90Array1dCreate((void *)weights, MPIU_REAL, 1, (*npoints) * (*Nc), ptrW PETSC_F90_2PTR_PARAM(ptrw));
24: }
26: PETSC_EXTERN void petscquadraturerestoredata_(PetscQuadrature *q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, F90Array1d *ptrP, F90Array1d *ptrW, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrp) PETSC_F90_2PTR_PROTO(ptrw))
27: {
28: *ierr = F90Array1dDestroy(ptrP, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrp));
29: if (*ierr) return;
30: *ierr = F90Array1dDestroy(ptrW, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrw));
31: }
33: PETSC_EXTERN void petscquadraturesetdata_(PetscQuadrature *q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, F90Array1d *ptrP, F90Array1d *ptrW, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrp) PETSC_F90_2PTR_PROTO(ptrw))
34: {
35: PetscReal *points, *weights;
37: *ierr = F90Array1dAccess(ptrP, MPIU_REAL, (void **)&points PETSC_F90_2PTR_PARAM(ptrp));
38: if (*ierr) return;
39: *ierr = F90Array1dAccess(ptrW, MPIU_REAL, (void **)&weights PETSC_F90_2PTR_PARAM(ptrw));
40: if (*ierr) return;
41: *ierr = PetscQuadratureSetData(*q, *dim, *Nc, *npoints, points, weights);
42: }