xref: /petsc/src/dm/dt/interface/ftn-custom/zdsf.c (revision 8434afd195968570cfdb5bc7b9cfc0a316d974ae)
1fe2efc57SMark #include <petsc/private/fortranimpl.h>
2fe2efc57SMark #include <petscds.h>
3fe2efc57SMark #include <petscviewer.h>
4fe2efc57SMark 
5fe2efc57SMark #if defined(PETSC_HAVE_FORTRAN_CAPS)
6fe2efc57SMark   #define petscdsviewfromoptions_  PETSCDSVIEWFROMOPTIONS
7e0e713ddStbridel   #define petscdsview_             PETSCDSVIEW
8e0e713ddStbridel   #define petscdssetcontext_       PETSCDSSETCONTEXT
9e0e713ddStbridel   #define petscdssetriemannsolver_ PETSCDSSETRIEMANNSOLVER
10fe2efc57SMark #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
11fe2efc57SMark   #define petscdsviewfromoptions_  petscdsviewfromoptions
12e0e713ddStbridel   #define petscdsview_             petscdsview
13e0e713ddStbridel   #define petscdssetcontext_       petscdssetcontext
14e0e713ddStbridel   #define petscdssetriemannsolver_ petscdssetriemannsolver
15fe2efc57SMark #endif
16fe2efc57SMark 
17e0e713ddStbridel static PetscFortranCallbackId riemannsolver;
18e0e713ddStbridel 
190f1503a6SJed Brown // We can't use PetscObjectUseFortranCallback() because this function returns void
200f1503a6SJed Brown static void ourriemannsolver(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx)
21e0e713ddStbridel {
220f1503a6SJed Brown   void (*func)(PetscInt *dim, PetscInt *Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], const PetscInt *numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx);
230f1503a6SJed Brown   void *_ctx;
24*8434afd1SBarry Smith   PetscCallAbort(PETSC_COMM_SELF, PetscObjectGetFortranCallback((PetscObject)ctx, PETSC_FORTRAN_CALLBACK_CLASS, riemannsolver, (PetscVoidFn **)&func, &_ctx));
255975b3b6SBarry Smith   if (func) { (*func)(&dim, &Nf, x, n, uL, uR, &numConstants, constants, flux, _ctx); }
26e0e713ddStbridel }
27e0e713ddStbridel 
2819caf8f3SSatish Balay PETSC_EXTERN void petscdsviewfromoptions_(PetscDS *ao, PetscObject obj, char *type, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
29fe2efc57SMark {
30fe2efc57SMark   char *t;
31fe2efc57SMark 
32fe2efc57SMark   FIXCHAR(type, len, t);
33b14c0cbaSBlaise Bourdin   CHKFORTRANNULLOBJECT(obj);
345975b3b6SBarry Smith   *ierr = PetscDSViewFromOptions(*ao, obj, t);
355975b3b6SBarry Smith   if (*ierr) return;
36fe2efc57SMark   FREECHAR(type, t);
37fe2efc57SMark }
38fe2efc57SMark 
39e0e713ddStbridel PETSC_EXTERN void petscdsview_(PetscDS *prob, PetscViewer *vin, PetscErrorCode *ierr)
40e0e713ddStbridel {
41e0e713ddStbridel   PetscViewer v;
42e0e713ddStbridel   PetscPatchDefaultViewers_Fortran(vin, v);
435975b3b6SBarry Smith   *ierr = PetscDSView(*prob, v);
445975b3b6SBarry Smith   if (*ierr) return;
45e0e713ddStbridel }
46e0e713ddStbridel 
47e0e713ddStbridel PETSC_EXTERN void petscdssetcontext_(PetscDS *prob, PetscInt *f, void *ctx, PetscErrorCode *ierr)
48e0e713ddStbridel {
495975b3b6SBarry Smith   *ierr = PetscDSSetContext(*prob, *f, *prob);
505975b3b6SBarry Smith   if (*ierr) return;
51e0e713ddStbridel }
52e0e713ddStbridel 
53e0e713ddStbridel PETSC_EXTERN void petscdssetriemannsolver_(PetscDS *prob, PetscInt *f, void (*rs)(PetscInt *, PetscInt *, PetscReal *, PetscReal *, PetscScalar *, PetscScalar *, PetscInt *, PetscScalar *, PetscScalar *, void *, PetscErrorCode *), PetscErrorCode *ierr)
54e0e713ddStbridel {
55*8434afd1SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*prob, PETSC_FORTRAN_CALLBACK_CLASS, &riemannsolver, (PetscVoidFn *)rs, NULL);
565975b3b6SBarry Smith   if (*ierr) return;
575975b3b6SBarry Smith   *ierr = PetscDSSetRiemannSolver(*prob, *f, ourriemannsolver);
585975b3b6SBarry Smith   if (*ierr) return;
59e0e713ddStbridel }
60