xref: /petsc/src/dm/dt/interface/ftn-custom/zdtf90.c (revision 6dd63270497ad23dcf16ae500a87ff2b2a0b7474)
1*6dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2*6dd63270SBarry Smith #include <petscdt.h>
3*6dd63270SBarry Smith 
4*6dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
5*6dd63270SBarry Smith   #define petscquadraturegetdata_     PETSCQUADRATUREGETDATA
6*6dd63270SBarry Smith   #define petscquadraturerestoredata_ PETSCQUADRATURERESTOREDATA
7*6dd63270SBarry Smith   #define petscquadraturesetdata_     PETSCQUADRATURESETDATA
8*6dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9*6dd63270SBarry Smith   #define petscquadraturegetdata_     petscquadraturegetdata
10*6dd63270SBarry Smith   #define petscquadraturerestoredata_ petscquadraturerestoredata
11*6dd63270SBarry Smith   #define petscquadraturesetdata_     petscquadraturesetdata
12*6dd63270SBarry Smith #endif
13*6dd63270SBarry Smith 
14*6dd63270SBarry 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))
15*6dd63270SBarry Smith {
16*6dd63270SBarry Smith   const PetscReal *points, *weights;
17*6dd63270SBarry Smith 
18*6dd63270SBarry Smith   *ierr = PetscQuadratureGetData(*q, dim, Nc, npoints, &points, &weights);
19*6dd63270SBarry Smith   if (*ierr) return;
20*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)points, MPIU_REAL, 1, (*npoints) * (*dim), ptrP PETSC_F90_2PTR_PARAM(ptrp));
21*6dd63270SBarry Smith   if (*ierr) return;
22*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)weights, MPIU_REAL, 1, (*npoints) * (*Nc), ptrW PETSC_F90_2PTR_PARAM(ptrw));
23*6dd63270SBarry Smith }
24*6dd63270SBarry Smith 
25*6dd63270SBarry 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))
26*6dd63270SBarry Smith {
27*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptrP, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrp));
28*6dd63270SBarry Smith   if (*ierr) return;
29*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptrW, MPIU_REAL PETSC_F90_2PTR_PARAM(ptrw));
30*6dd63270SBarry Smith }
31*6dd63270SBarry Smith 
32*6dd63270SBarry Smith PETSC_EXTERN void petscquadraturesetdata_NOTTODAY(PetscQuadrature *q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, F90Array1d *ptrP, F90Array1d *ptrW, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrp) PETSC_F90_2PTR_PROTO(ptrw))
33*6dd63270SBarry Smith {
34*6dd63270SBarry Smith   PetscReal *points, *weights;
35*6dd63270SBarry Smith 
36*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptrP, MPIU_REAL, (void **)&points PETSC_F90_2PTR_PARAM(ptrp));
37*6dd63270SBarry Smith   if (*ierr) return;
38*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptrW, MPIU_REAL, (void **)&weights PETSC_F90_2PTR_PARAM(ptrw));
39*6dd63270SBarry Smith   if (*ierr) return;
40*6dd63270SBarry Smith   *ierr = PetscQuadratureSetData(*q, *dim, *Nc, *npoints, points, weights);
41*6dd63270SBarry Smith }
42