xref: /petsc/src/ksp/pc/impls/mg/ftn-custom/zmgfuncf.c (revision f68b968ce39302dfa79eb1a6cfa1998ce074e829)
1e54e4138SSatish Balay #include "zpetsc.h"
2e54e4138SSatish Balay #include "petscpc.h"
3e54e4138SSatish Balay #include "petscmg.h"
4e54e4138SSatish Balay 
5e54e4138SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
6e54e4138SSatish Balay #define pcmgsetresidual_           PCMGSETRESIDUAL
7e54e4138SSatish Balay #define pcmgdefaultresidual_       PCMGDEFAULTRESIDUAL
8e54e4138SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9e54e4138SSatish Balay #define pcmgsetresidual_           pcmgsetresidual
10e54e4138SSatish Balay #define pcmgdefaultresidual_       pcmgdefaultresidual
11e54e4138SSatish Balay #endif
12e54e4138SSatish Balay 
13e54e4138SSatish Balay typedef PetscErrorCode (*MVVVV)(Mat,Vec,Vec,Vec);
14e54e4138SSatish Balay static PetscErrorCode ourresidualfunction(Mat mat,Vec b,Vec x,Vec R)
15e54e4138SSatish Balay {
16e54e4138SSatish Balay   PetscErrorCode ierr = 0;
17e54e4138SSatish Balay   (*(void (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[0]))(&mat,&b,&x,&R,&ierr);
18e54e4138SSatish Balay   return 0;
19e54e4138SSatish Balay }
20e54e4138SSatish Balay 
21e54e4138SSatish Balay EXTERN_C_BEGIN
22e54e4138SSatish Balay extern void PETSC_STDCALL pcmgdefaultresidual_(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*);
23e54e4138SSatish Balay 
24e54e4138SSatish Balay void PETSC_STDCALL pcmgsetresidual_(PC *pc,PetscInt *l,PetscErrorCode (*residual)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*),Mat *mat, PetscErrorCode *ierr)
25e54e4138SSatish Balay {
26e54e4138SSatish Balay   MVVVV rr;
27*f68b968cSBarry Smith   if ((PetscVoidFunction)residual == (PetscVoidFunction)pcmgdefaultresidual_) rr = PCMGDefaultResidual;
28e54e4138SSatish Balay   else {
29e54e4138SSatish Balay     if (!((PetscObject)*mat)->fortran_func_pointers) {
30e54e4138SSatish Balay       *ierr = PetscMalloc(1*sizeof(void*),&((PetscObject)*mat)->fortran_func_pointers);
31e54e4138SSatish Balay     }
32*f68b968cSBarry Smith     ((PetscObject)*mat)->fortran_func_pointers[0] = (PetscVoidFunction)residual;
33e54e4138SSatish Balay     rr = ourresidualfunction;
34e54e4138SSatish Balay   }
35e54e4138SSatish Balay   *ierr = PCMGSetResidual(*pc,*l,rr,*mat);
36e54e4138SSatish Balay }
37e54e4138SSatish Balay 
38e54e4138SSatish Balay EXTERN_C_END
39