xref: /petsc/src/dm/dt/interface/ftn-custom/zdsf.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 #include <petsc/private/fortranimpl.h>
2 #include <petscds.h>
3 #include <petscviewer.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6 #define petscdsviewfromoptions_   PETSCDSVIEWFROMOPTIONS
7 #define petscdsview_  PETSCDSVIEW
8 #define petscdssetcontext_  PETSCDSSETCONTEXT
9 #define petscdssetriemannsolver_ PETSCDSSETRIEMANNSOLVER
10 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
11 #define petscdsviewfromoptions_   petscdsviewfromoptions
12 #define petscdsview_  petscdsview
13 #define petscdssetcontext_ petscdssetcontext
14 #define petscdssetriemannsolver_ petscdssetriemannsolver
15 #endif
16 
17 static PetscFortranCallbackId riemannsolver;
18 
19 static PetscErrorCode ourriemannsolver(PetscInt dim,PetscInt Nf,PetscReal x[],PetscReal n[],PetscScalar uL[],PetscScalar uR[],PetscInt numConstants,PetscScalar constants[],PetscScalar flux[],void *ctx)
20 {
21   PetscObjectUseFortranCallback((PetscDS)ctx,riemannsolver,(PetscInt*,PetscInt*,PetscReal*,PetscReal*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*),(&dim,&Nf,x,n,uL,uR,&numConstants,constants,flux,_ctx,&ierr));
22 }
23 
24 PETSC_EXTERN void petscdsviewfromoptions_(PetscDS *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
25 {
26   char *t;
27 
28   FIXCHAR(type,len,t);
29   CHKFORTRANNULLOBJECT(obj);
30   *ierr = PetscDSViewFromOptions(*ao,obj,t);if (*ierr) return;
31   FREECHAR(type,t);
32 }
33 
34 PETSC_EXTERN void petscdsview_(PetscDS *prob,PetscViewer *vin,PetscErrorCode *ierr)
35 {
36   PetscViewer v;
37   PetscPatchDefaultViewers_Fortran(vin,v);
38   *ierr = PetscDSView(*prob,v);if (*ierr) return;
39 }
40 
41 PETSC_EXTERN void petscdssetcontext_(PetscDS *prob,PetscInt *f,void *ctx,PetscErrorCode *ierr)
42 {
43   *ierr = PetscDSSetContext(*prob,*f,*prob);if (*ierr) return;
44 }
45 
46 PETSC_EXTERN void petscdssetriemannsolver_(PetscDS *prob,PetscInt *f,void (*rs)(PetscInt*,PetscInt*,PetscReal*,PetscReal*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*),PetscErrorCode *ierr)
47 {
48   *ierr = PetscObjectSetFortranCallback((PetscObject)*prob,PETSC_FORTRAN_CALLBACK_CLASS,&riemannsolver,(PetscVoidFunction)rs,NULL);if (*ierr) return;
49   *ierr = PetscDSSetRiemannSolver(*prob,*f,(void*)ourriemannsolver);if (*ierr) return;
50 }
51