1*b45d2f2cSJed Brown #include <petsc-private/fortranimpl.h> 2c6db04a5SJed Brown #include <petscmat.h> 3f4e70085SSatish Balay 4f4e70085SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS) 5f4e70085SSatish Balay #define matshellsetoperation_ MATSHELLSETOPERATION 6f4e70085SSatish Balay #define matcreateshell_ MATCREATESHELL 7f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 8f4e70085SSatish Balay #define matcreateshell_ matcreateshell 9f4e70085SSatish Balay #define matshellsetoperation_ matshellsetoperation 10f4e70085SSatish Balay #endif 11f4e70085SSatish Balay 12f4e70085SSatish Balay EXTERN_C_BEGIN 13f4e70085SSatish Balay 14f4e70085SSatish Balay /* 15f4e70085SSatish Balay The MatShell Matrix Vector product requires a C routine. 16f4e70085SSatish Balay This C routine then calls the corresponding Fortran routine that was 17f4e70085SSatish Balay set by the user. 18f4e70085SSatish Balay */ 192e843561SJed Brown void PETSC_STDCALL matcreateshell_(MPI_Comm *comm,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N,void *ctx,Mat *mat,PetscErrorCode *ierr) 20f4e70085SSatish Balay { 212e843561SJed Brown *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint *)&*comm),*m,*n,*M,*N,ctx,mat); 22f4e70085SSatish Balay } 23f4e70085SSatish Balay 24f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat,Vec x,Vec y) 25f4e70085SSatish Balay { 26f4e70085SSatish Balay PetscErrorCode ierr = 0; 27f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[0]))(&mat,&x,&y,&ierr); 28f4e70085SSatish Balay return ierr; 29f4e70085SSatish Balay } 30f4e70085SSatish Balay 31f4e70085SSatish Balay static PetscErrorCode ourmulttranspose(Mat mat,Vec x,Vec y) 32f4e70085SSatish Balay { 33f4e70085SSatish Balay PetscErrorCode ierr = 0; 34f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[2]))(&mat,&x,&y,&ierr); 35f4e70085SSatish Balay return ierr; 36f4e70085SSatish Balay } 37f4e70085SSatish Balay 38f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat,Vec x,Vec y,Vec z) 39f4e70085SSatish Balay { 40f4e70085SSatish Balay PetscErrorCode ierr = 0; 41f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[1]))(&mat,&x,&y,&z,&ierr); 42f4e70085SSatish Balay return ierr; 43f4e70085SSatish Balay } 44f4e70085SSatish Balay 45f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat,Vec x,Vec y,Vec z) 46f4e70085SSatish Balay { 47f4e70085SSatish Balay PetscErrorCode ierr = 0; 48f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[3]))(&mat,&x,&y,&z,&ierr); 49f4e70085SSatish Balay return ierr; 50f4e70085SSatish Balay } 51f4e70085SSatish Balay 5222612f2fSMatthew Knepley static PetscErrorCode ourgetdiagonal(Mat mat,Vec x) 5322612f2fSMatthew Knepley { 5422612f2fSMatthew Knepley PetscErrorCode ierr = 0; 5522612f2fSMatthew Knepley (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[4]))(&mat,&x,&ierr); 5622612f2fSMatthew Knepley return ierr; 5722612f2fSMatthew Knepley } 5822612f2fSMatthew Knepley 59160922c2SBarry Smith static PetscErrorCode ourdiagonalscale(Mat mat,Vec l,Vec r) 60160922c2SBarry Smith { 61160922c2SBarry Smith PetscErrorCode ierr = 0; 62160922c2SBarry Smith if (!l) { 6335b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,(Vec*)PETSC_NULL_OBJECT_Fortran,&r,&ierr); 64160922c2SBarry Smith } else if (!r) { 6535b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,(Vec*)PETSC_NULL_OBJECT_Fortran,&ierr); 66160922c2SBarry Smith } else { 67160922c2SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,&r,&ierr); 68160922c2SBarry Smith } 69160922c2SBarry Smith return ierr; 70160922c2SBarry Smith } 71160922c2SBarry Smith 727911a512SBarry Smith static PetscErrorCode ourgetvecs(Mat mat,Vec *l,Vec *r) 737911a512SBarry Smith { 747911a512SBarry Smith PetscErrorCode ierr = 0; 75501d9185SBarry Smith PetscInt none = -1; 767911a512SBarry Smith if (!l) { 77501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,(Vec*)&none,r,&ierr); 787911a512SBarry Smith } else if (!r) { 79501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,(Vec*)&none,&ierr); 807911a512SBarry Smith } else { 817911a512SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,r,&ierr); 827911a512SBarry Smith } 837911a512SBarry Smith return ierr; 847911a512SBarry Smith } 857911a512SBarry Smith 86f5a4496aSBarry Smith static PetscErrorCode ourdiagonalset(Mat mat,Vec x,InsertMode ins) 87f5a4496aSBarry Smith { 88f5a4496aSBarry Smith PetscErrorCode ierr = 0; 89f5a4496aSBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,InsertMode*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[6]))(&mat,&x,&ins,&ierr); 90f5a4496aSBarry Smith return ierr; 91f5a4496aSBarry Smith } 92f5a4496aSBarry Smith 932950f7e7SBarry Smith static PetscErrorCode ourview(Mat mat,PetscViewer v) 942950f7e7SBarry Smith { 952950f7e7SBarry Smith PetscErrorCode ierr = 0; 962950f7e7SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscViewer*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[8]))(&mat,&v,&ierr); 972950f7e7SBarry Smith return ierr; 982950f7e7SBarry Smith } 992950f7e7SBarry Smith 1003446fae8SBarry Smith static PetscErrorCode oursor(Mat mat,Vec b,PetscReal omega,MatSORType flg,PetscReal shift,PetscInt its,PetscInt lits,Vec x) 1013446fae8SBarry Smith { 1023446fae8SBarry Smith PetscErrorCode ierr = 0; 1033446fae8SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscReal*,MatSORType*,PetscReal*,PetscInt*,PetscInt*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[9]))(&mat,&b,&omega,&flg,&shift,&its,&lits,&x,&ierr); 1043446fae8SBarry Smith return ierr; 1053446fae8SBarry Smith } 1063446fae8SBarry Smith 107cdf26a31SSatish Balay static PetscErrorCode ourshift(Mat mat, PetscScalar a) 108cdf26a31SSatish Balay { 109cdf26a31SSatish Balay PetscErrorCode ierr = 0; 110cdf26a31SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[10]))(&mat,&a,&ierr); 111cdf26a31SSatish Balay return ierr; 112cdf26a31SSatish Balay } 113cdf26a31SSatish Balay 114f4e70085SSatish Balay void PETSC_STDCALL matshellsetoperation_(Mat *mat,MatOperation *op,PetscErrorCode (PETSC_STDCALL *f)(Mat*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode *ierr) 115f4e70085SSatish Balay { 116e32f2f54SBarry Smith MPI_Comm comm; 117e32f2f54SBarry Smith 118e32f2f54SBarry Smith *ierr = PetscObjectGetComm((PetscObject)*mat,&comm);if (*ierr) return; 119cdf26a31SSatish Balay PetscObjectAllocateFortranPointers(*mat,11); 120f4e70085SSatish Balay if (*op == MATOP_MULT) { 121f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmult); 122f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[0] = (PetscVoidFunction)f; 123f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE) { 124f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttranspose); 125f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[2] = (PetscVoidFunction)f; 126f4e70085SSatish Balay } else if (*op == MATOP_MULT_ADD) { 127f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmultadd); 128f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[1] = (PetscVoidFunction)f; 129f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE_ADD) { 130f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttransposeadd); 131f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[3] = (PetscVoidFunction)f; 13222612f2fSMatthew Knepley } else if (*op == MATOP_GET_DIAGONAL) { 13322612f2fSMatthew Knepley *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetdiagonal); 13422612f2fSMatthew Knepley ((PetscObject)*mat)->fortran_func_pointers[4] = (PetscVoidFunction)f; 135160922c2SBarry Smith } else if (*op == MATOP_DIAGONAL_SCALE) { 136160922c2SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalscale); 137160922c2SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[5] = (PetscVoidFunction)f; 13835153367SBarry Smith } else if (*op == MATOP_DIAGONAL_SET) { 139f5a4496aSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalset); 140f5a4496aSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[6] = (PetscVoidFunction)f; 1417911a512SBarry Smith } else if (*op == MATOP_GET_VECS) { 1427911a512SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetvecs); 1437911a512SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[7] = (PetscVoidFunction)f; 1442950f7e7SBarry Smith } else if (*op == MATOP_VIEW) { 1452950f7e7SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourview); 1462950f7e7SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[8] = (PetscVoidFunction)f; 1473446fae8SBarry Smith } else if (*op == MATOP_SOR) { 1483446fae8SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)oursor); 1493446fae8SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[9] = (PetscVoidFunction)f; 150cdf26a31SSatish Balay } else if (*op == MATOP_SHIFT) { 151cdf26a31SSatish Balay *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourshift); 152cdf26a31SSatish Balay ((PetscObject)*mat)->fortran_func_pointers[10] = (PetscVoidFunction)f; 153f4e70085SSatish Balay } else { 154d736bfebSBarry Smith PetscError(comm,__LINE__,"MatShellSetOperation_Fortran",__FILE__,__SDIR__,PETSC_ERR_ARG_WRONG,PETSC_ERROR_INITIAL, 155f4e70085SSatish Balay "Cannot set that matrix operation"); 156f4e70085SSatish Balay *ierr = 1; 157f4e70085SSatish Balay } 158f4e70085SSatish Balay } 159f4e70085SSatish Balay 160f4e70085SSatish Balay EXTERN_C_END 161