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