xref: /petsc/src/mat/matfd/ftn-custom/zfdmatrixf.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1edaa7c33SBarry Smith #include <petsc/private/f90impl.h>
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
3fcfc5002SJed Brown 
4fcfc5002SJed Brown /* Declare these pointer types instead of void* for clarity, but do not include petscts.h so that this code does have an actual reverse dependency. */
5fcfc5002SJed Brown typedef struct _p_TS *TS;
6fcfc5002SJed Brown typedef struct _p_SNES *SNES;
7f4e70085SSatish Balay 
8f4e70085SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
9f4e70085SSatish Balay #define matfdcoloringsetfunctionts_              MATFDCOLORINGSETFUNCTIONTS
1089d42083SSatish Balay #define matfdcoloringsetfunction_                MATFDCOLORINGSETFUNCTION
111d2e4005SSatish Balay #define matfdcoloringview_                       MATFDCOLORINGVIEW
1261ab5df0SBarry Smith #define matfdcoloingsettype_                     MATFDCOLORINGSETTYPE
13edaa7c33SBarry Smith #define matfdcoloringgetperturbedcolumnsf90_     MATFDCOLORINGGETPERTURBEDCOLUMNSF90
14edaa7c33SBarry Smith #define matfdcoloringrestoreperturbedcolumnsf90_ MATFDCOLORINGRESTOREPERTURBEDCOLUMNSF90
15f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
16f4e70085SSatish Balay #define matfdcoloringsetfunctionts_              matfdcoloringsetfunctionts
17372a5eeaSSatish Balay #define matfdcoloringsetfunction_                matfdcoloringsetfunction
181d2e4005SSatish Balay #define matfdcoloringview_                       matfdcoloringview
1961ab5df0SBarry Smith #define matfdcoloingsettype_                     matfdcoloringsettype
20edaa7c33SBarry Smith #define matfdcoloringgetperturbedcolumnsf90_     matfdcoloringgetperturbedcolumnsf90
21edaa7c33SBarry Smith #define matfdcoloringrestoreperturbedcolumnsf90_ matfdcoloringrestoreperturbedcolumnsf90
22f4e70085SSatish Balay #endif
23f4e70085SSatish Balay 
2419caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringgetperturbedcolumnsf90_(MatFDColoring *x,F90Array1d *ptr,int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
25edaa7c33SBarry Smith {
26edaa7c33SBarry Smith   const PetscInt *fa;
27edaa7c33SBarry Smith   PetscInt       len;
28edaa7c33SBarry Smith 
29edaa7c33SBarry Smith   *__ierr = MatFDColoringGetPerturbedColumns(*x,&len,&fa);   if (*__ierr) return;
30b7b8f77aSBarry Smith   *__ierr = F90Array1dCreate((void*)fa,MPIU_INT,1,len,ptr PETSC_F90_2PTR_PARAM(ptrd));
31edaa7c33SBarry Smith }
3219caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringrestoreperturbedcolumnsf90_(MatFDColoring *x,F90Array1d *ptr,int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
33edaa7c33SBarry Smith {
34b7b8f77aSBarry Smith   *__ierr = F90Array1dDestroy(ptr,MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
35edaa7c33SBarry Smith }
36f4e70085SSatish Balay 
37f4e70085SSatish Balay /* These are not extern C because they are passed into non-extern C user level functions */
387850c7c0SBarry Smith static PetscErrorCode ourmatfdcoloringfunctionts(TS ts, PetscReal t, Vec x, Vec y, MatFDColoring fd)
39f4e70085SSatish Balay {
40*3ba16761SJacob Faibussowitsch   PetscErrorCode ierr = PETSC_SUCCESS;
4119caf8f3SSatish Balay   (*(void (*)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *))(fd->ftn_func_pointer))(&ts, &t, &x, &y, fd->ftn_func_cntx, &ierr);
42f4e70085SSatish Balay   return ierr;
43f4e70085SSatish Balay }
44f4e70085SSatish Balay 
457850c7c0SBarry Smith static PetscErrorCode ourmatfdcoloringfunctionsnes(SNES snes, Vec x, Vec y, MatFDColoring fd)
46f4e70085SSatish Balay {
47*3ba16761SJacob Faibussowitsch   PetscErrorCode ierr = PETSC_SUCCESS;
4819caf8f3SSatish Balay   (*(void (*)(SNES *, Vec *, Vec *, void *, PetscErrorCode *))(fd->ftn_func_pointer))(&snes, &x, &y, fd->ftn_func_cntx, &ierr);
49f4e70085SSatish Balay   return ierr;
50f4e70085SSatish Balay }
51f4e70085SSatish Balay 
52f4e70085SSatish Balay /*
537850c7c0SBarry Smith         MatFDColoringSetFunction sticks the Fortran function and its context into the MatFDColoring structure and passes the MatFDColoring object
547850c7c0SBarry Smith     in as the function context. ourmafdcoloringfunctionsnes() and ourmatfdcoloringfunctionts()  then access the function and its context from the
557850c7c0SBarry Smith     MatFDColoring that is passed in. This is the same way that fortran_func_pointers is used in PETSc objects.
56f4e70085SSatish Balay 
57f4e70085SSatish Balay    NOTE: FORTRAN USER CANNOT PUT IN A NEW J OR B currently.
58f4e70085SSatish Balay */
59f4e70085SSatish Balay 
6019caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringsetfunctionts_(MatFDColoring *fd,void (*f)(TS*,double*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
61f4e70085SSatish Balay {
6278bd23c2SBarry Smith   (*fd)->ftn_func_pointer =  (void (*)(void)) f;
637850c7c0SBarry Smith   (*fd)->ftn_func_cntx    = ctx;
648865f1eaSKarl Rupp 
657850c7c0SBarry Smith   *ierr = MatFDColoringSetFunction(*fd,(PetscErrorCodeFunction)ourmatfdcoloringfunctionts,*fd);
66f4e70085SSatish Balay }
67f4e70085SSatish Balay 
6819caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringsetfunction_(MatFDColoring *fd,void (*f)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
69f4e70085SSatish Balay {
7078bd23c2SBarry Smith   (*fd)->ftn_func_pointer = (void (*)(void)) f;
717850c7c0SBarry Smith   (*fd)->ftn_func_cntx    = ctx;
728865f1eaSKarl Rupp 
737850c7c0SBarry Smith   *ierr = MatFDColoringSetFunction(*fd,(PetscErrorCodeFunction)ourmatfdcoloringfunctionsnes,*fd);
74f4e70085SSatish Balay }
75f4e70085SSatish Balay 
7619caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringview_(MatFDColoring *c,PetscViewer *vin,PetscErrorCode *ierr)
771d2e4005SSatish Balay {
781d2e4005SSatish Balay   PetscViewer v;
791d2e4005SSatish Balay 
801d2e4005SSatish Balay   PetscPatchDefaultViewers_Fortran(vin,v);
811d2e4005SSatish Balay   *ierr = MatFDColoringView(*c,v);
821d2e4005SSatish Balay }
831d2e4005SSatish Balay 
8419caf8f3SSatish Balay PETSC_EXTERN void matfdcoloringsettype_(MatFDColoring *matfdcoloring,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
8561ab5df0SBarry Smith {
8661ab5df0SBarry Smith   char *t;
8761ab5df0SBarry Smith 
8861ab5df0SBarry Smith   FIXCHAR(type,len,t);
89d49bb8f9SBarry Smith   *ierr = MatFDColoringSetType(*matfdcoloring,t);if (*ierr) return;
9061ab5df0SBarry Smith   FREECHAR(type,t);
9161ab5df0SBarry Smith }
92