xref: /petsc/src/ksp/pc/impls/mg/ftn-custom/zmgfuncf.c (revision d0e4de75ed9992a507d929c42cd4edab7cb0b081)
1b45d2f2cSJed Brown #include <petsc-private/fortranimpl.h>
2c6db04a5SJed Brown #include <petscpc.h>
3*d0e4de75SBarry Smith #include <../src/ksp/pc/impls/mg/mgimpl.h>
4e54e4138SSatish Balay 
5e54e4138SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
6e54e4138SSatish Balay #define pcmgsetresidual_           PCMGSETRESIDUAL
7*d0e4de75SBarry Smith #define pcmgresidual_default_       PCMGRESIDUAL_DEFAULT
8e54e4138SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9e54e4138SSatish Balay #define pcmgsetresidual_           pcmgsetresidual
10*d0e4de75SBarry Smith #define pcmgresidual_default_       pcmgresidual_default
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 
21*d0e4de75SBarry Smith PETSC_EXTERN void pcmgresidual_default_(Mat *mat,Vec *b,Vec *x,Vec *r, PetscErrorCode *ierr)
221f6cc5b2SSatish Balay {
23*d0e4de75SBarry Smith   *ierr = PCMGResidual_Default(*mat,*b,*x,*r);
241f6cc5b2SSatish Balay }
25e54e4138SSatish Balay 
268cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcmgsetresidual_(PC *pc,PetscInt *l,PetscErrorCode (*residual)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*),Mat *mat, PetscErrorCode *ierr)
27e54e4138SSatish Balay {
28e54e4138SSatish Balay   MVVVV rr;
29*d0e4de75SBarry Smith   if ((PetscVoidFunction)residual == (PetscVoidFunction)pcmgresidual_default_) rr = PCMGResidual_Default;
30e54e4138SSatish Balay   else {
317850c7c0SBarry Smith     PetscObjectAllocateFortranPointers(*mat,1);
327850c7c0SBarry Smith     /*  Attach the residual computer to the Mat, this is not ideal but the only object/context passed in the residual computer */
33f68b968cSBarry Smith     ((PetscObject)*mat)->fortran_func_pointers[0] = (PetscVoidFunction)residual;
342fa5cd67SKarl Rupp 
35e54e4138SSatish Balay     rr = ourresidualfunction;
36e54e4138SSatish Balay   }
37e54e4138SSatish Balay   *ierr = PCMGSetResidual(*pc,*l,rr,*mat);
38e54e4138SSatish Balay }
39e54e4138SSatish Balay 
40