1ce0a2cd1SBarry Smith #include "private/fortranimpl.h" 2f4e70085SSatish Balay #include "petscmat.h" 3f4e70085SSatish Balay 4f4e70085SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS) 5f4e70085SSatish Balay #define matshellsetoperation_ MATSHELLSETOPERATION 6f4e70085SSatish Balay #define matcreateshell_ MATCREATESHELL 7c6866cfdSSatish Balay #define matshellgetcontext_ MATSHELLGETCONTEXT 8f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 9f4e70085SSatish Balay #define matcreateshell_ matcreateshell 10f4e70085SSatish Balay #define matshellsetoperation_ matshellsetoperation 11c6866cfdSSatish Balay #define matshellgetcontext_ matshellgetcontext 12f4e70085SSatish Balay #endif 13f4e70085SSatish Balay 14f4e70085SSatish Balay EXTERN_C_BEGIN 15f4e70085SSatish Balay 16f4e70085SSatish Balay /* 17f4e70085SSatish Balay The MatShell Matrix Vector product requires a C routine. 18f4e70085SSatish Balay This C routine then calls the corresponding Fortran routine that was 19f4e70085SSatish Balay set by the user. 20f4e70085SSatish Balay */ 21f4e70085SSatish Balay void PETSC_STDCALL matcreateshell_(MPI_Comm *comm,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N,void **ctx,Mat *mat,PetscErrorCode *ierr) 22f4e70085SSatish Balay { 23a542b6e8SBarry Smith *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint *)&*comm),*m,*n,*M,*N,*ctx,mat); 245db8bc65SBarry Smith 25f4e70085SSatish Balay } 26f4e70085SSatish Balay 27f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat,Vec x,Vec y) 28f4e70085SSatish Balay { 29f4e70085SSatish Balay PetscErrorCode ierr = 0; 30f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[0]))(&mat,&x,&y,&ierr); 31f4e70085SSatish Balay return ierr; 32f4e70085SSatish Balay } 33f4e70085SSatish Balay 34f4e70085SSatish Balay static PetscErrorCode ourmulttranspose(Mat mat,Vec x,Vec y) 35f4e70085SSatish Balay { 36f4e70085SSatish Balay PetscErrorCode ierr = 0; 37f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[2]))(&mat,&x,&y,&ierr); 38f4e70085SSatish Balay return ierr; 39f4e70085SSatish Balay } 40f4e70085SSatish Balay 41f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat,Vec x,Vec y,Vec z) 42f4e70085SSatish Balay { 43f4e70085SSatish Balay PetscErrorCode ierr = 0; 44f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[1]))(&mat,&x,&y,&z,&ierr); 45f4e70085SSatish Balay return ierr; 46f4e70085SSatish Balay } 47f4e70085SSatish Balay 48f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat,Vec x,Vec y,Vec z) 49f4e70085SSatish Balay { 50f4e70085SSatish Balay PetscErrorCode ierr = 0; 51f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[3]))(&mat,&x,&y,&z,&ierr); 52f4e70085SSatish Balay return ierr; 53f4e70085SSatish Balay } 54f4e70085SSatish Balay 5522612f2fSMatthew Knepley static PetscErrorCode ourgetdiagonal(Mat mat,Vec x) 5622612f2fSMatthew Knepley { 5722612f2fSMatthew Knepley PetscErrorCode ierr = 0; 5822612f2fSMatthew Knepley (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[4]))(&mat,&x,&ierr); 5922612f2fSMatthew Knepley return ierr; 6022612f2fSMatthew Knepley } 6122612f2fSMatthew Knepley 62160922c2SBarry Smith static PetscErrorCode ourdiagonalscale(Mat mat,Vec l,Vec r) 63160922c2SBarry Smith { 64160922c2SBarry Smith PetscErrorCode ierr = 0; 65160922c2SBarry Smith if (!l) { 6635b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,(Vec*)PETSC_NULL_OBJECT_Fortran,&r,&ierr); 67160922c2SBarry Smith } else if (!r) { 6835b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,(Vec*)PETSC_NULL_OBJECT_Fortran,&ierr); 69160922c2SBarry Smith } else { 70160922c2SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,&r,&ierr); 71160922c2SBarry Smith } 72160922c2SBarry Smith return ierr; 73160922c2SBarry Smith } 74160922c2SBarry Smith 757911a512SBarry Smith static PetscErrorCode ourgetvecs(Mat mat,Vec *l,Vec *r) 767911a512SBarry Smith { 777911a512SBarry Smith PetscErrorCode ierr = 0; 78501d9185SBarry Smith PetscInt none = -1; 797911a512SBarry Smith if (!l) { 80501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,(Vec*)&none,r,&ierr); 817911a512SBarry Smith } else if (!r) { 82501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,(Vec*)&none,&ierr); 837911a512SBarry Smith } else { 847911a512SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,r,&ierr); 857911a512SBarry Smith } 867911a512SBarry Smith return ierr; 877911a512SBarry Smith } 887911a512SBarry Smith 89f5a4496aSBarry Smith static PetscErrorCode ourdiagonalset(Mat mat,Vec x,InsertMode ins) 90f5a4496aSBarry Smith { 91f5a4496aSBarry Smith PetscErrorCode ierr = 0; 92f5a4496aSBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,InsertMode*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[6]))(&mat,&x,&ins,&ierr); 93f5a4496aSBarry Smith return ierr; 94f5a4496aSBarry Smith } 95f5a4496aSBarry Smith 96f4e70085SSatish Balay void PETSC_STDCALL matshellsetoperation_(Mat *mat,MatOperation *op,PetscErrorCode (PETSC_STDCALL *f)(Mat*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode *ierr) 97f4e70085SSatish Balay { 98e32f2f54SBarry Smith MPI_Comm comm; 99e32f2f54SBarry Smith 100e32f2f54SBarry Smith *ierr = PetscObjectGetComm((PetscObject)*mat,&comm);if (*ierr) return; 1017911a512SBarry Smith PetscObjectAllocateFortranPointers(*mat,8); 102f4e70085SSatish Balay if (*op == MATOP_MULT) { 103f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmult); 104f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[0] = (PetscVoidFunction)f; 105f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE) { 106f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttranspose); 107f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[2] = (PetscVoidFunction)f; 108f4e70085SSatish Balay } else if (*op == MATOP_MULT_ADD) { 109f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmultadd); 110f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[1] = (PetscVoidFunction)f; 111f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE_ADD) { 112f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttransposeadd); 113f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[3] = (PetscVoidFunction)f; 11422612f2fSMatthew Knepley } else if (*op == MATOP_GET_DIAGONAL) { 11522612f2fSMatthew Knepley *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetdiagonal); 11622612f2fSMatthew Knepley ((PetscObject)*mat)->fortran_func_pointers[4] = (PetscVoidFunction)f; 117160922c2SBarry Smith } else if (*op == MATOP_DIAGONAL_SCALE) { 118160922c2SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalscale); 119160922c2SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[5] = (PetscVoidFunction)f; 12035153367SBarry Smith } else if (*op == MATOP_DIAGONAL_SET) { 121f5a4496aSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalset); 122f5a4496aSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[6] = (PetscVoidFunction)f; 1237911a512SBarry Smith } else if (*op == MATOP_GET_VECS) { 1247911a512SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetvecs); 1257911a512SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[7] = (PetscVoidFunction)f; 126f4e70085SSatish Balay } else { 127*d736bfebSBarry Smith PetscError(comm,__LINE__,"MatShellSetOperation_Fortran",__FILE__,__SDIR__,PETSC_ERR_ARG_WRONG,PETSC_ERROR_INITIAL, 128f4e70085SSatish Balay "Cannot set that matrix operation"); 129f4e70085SSatish Balay *ierr = 1; 130f4e70085SSatish Balay } 131f4e70085SSatish Balay } 132f4e70085SSatish Balay 133c6866cfdSSatish Balay void PETSC_STDCALL matshellgetcontext_(Mat *mat,void **ctx,PetscErrorCode *ierr) 134c6866cfdSSatish Balay { 135c6866cfdSSatish Balay *ierr = MatShellGetContext(*mat,ctx); 136c6866cfdSSatish Balay } 137c6866cfdSSatish Balay 138c6866cfdSSatish Balay 139f4e70085SSatish Balay EXTERN_C_END 140