xref: /petsc/src/dm/dt/interface/ftn-custom/zdtf90.c (revision c9df829d41ec13fc316a62eed77ff350dff6e12a)
16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
26dd63270SBarry Smith #include <petscdt.h>
36dd63270SBarry Smith 
46dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
56dd63270SBarry Smith   #define petscquadraturegetdata_     PETSCQUADRATUREGETDATA
66dd63270SBarry Smith   #define petscquadraturerestoredata_ PETSCQUADRATURERESTOREDATA
76dd63270SBarry Smith   #define petscquadraturesetdata_     PETSCQUADRATURESETDATA
86dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
96dd63270SBarry Smith   #define petscquadraturegetdata_     petscquadraturegetdata
106dd63270SBarry Smith   #define petscquadraturerestoredata_ petscquadraturerestoredata
116dd63270SBarry Smith   #define petscquadraturesetdata_     petscquadraturesetdata
126dd63270SBarry Smith #endif
136dd63270SBarry Smith 
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))146dd63270SBarry Smith 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))
156dd63270SBarry Smith {
166dd63270SBarry Smith   const PetscReal *points, *weights;
176dd63270SBarry Smith 
186dd63270SBarry Smith   *ierr = PetscQuadratureGetData(*q, dim, Nc, npoints, &points, &weights);
196dd63270SBarry Smith   if (*ierr) return;
206dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)points, MPIU_REAL, 1, (*npoints) * (*dim), ptrP PETSC_F90_2PTR_PARAM(ptrp));
216dd63270SBarry Smith   if (*ierr) return;
226dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)weights, MPIU_REAL, 1, (*npoints) * (*Nc), ptrW PETSC_F90_2PTR_PARAM(ptrw));
236dd63270SBarry Smith }
246dd63270SBarry Smith 
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))256dd63270SBarry Smith 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))
266dd63270SBarry Smith {
276dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptrP, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrp));
286dd63270SBarry Smith   if (*ierr) return;
296dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptrW, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrw));
306dd63270SBarry Smith }
316dd63270SBarry Smith 
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))32*c080761bSJose E. Roman 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))
336dd63270SBarry Smith {
346dd63270SBarry Smith   PetscReal *points, *weights;
356dd63270SBarry Smith 
366dd63270SBarry Smith   *ierr = F90Array1dAccess(ptrP, MPIU_REAL, (void **)&points PETSC_F90_2PTR_PARAM(ptrp));
376dd63270SBarry Smith   if (*ierr) return;
386dd63270SBarry Smith   *ierr = F90Array1dAccess(ptrW, MPIU_REAL, (void **)&weights PETSC_F90_2PTR_PARAM(ptrw));
396dd63270SBarry Smith   if (*ierr) return;
406dd63270SBarry Smith   *ierr = PetscQuadratureSetData(*q, *dim, *Nc, *npoints, points, weights);
416dd63270SBarry Smith }
42