1b45d2f2cSJed 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 /* 13f4e70085SSatish Balay The MatShell Matrix Vector product requires a C routine. 14f4e70085SSatish Balay This C routine then calls the corresponding Fortran routine that was 15f4e70085SSatish Balay set by the user. 16f4e70085SSatish Balay */ 17*8cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matcreateshell_(MPI_Comm *comm,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N,void *ctx,Mat *mat,PetscErrorCode *ierr) 18f4e70085SSatish Balay { 192e843561SJed Brown *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint*)&*comm),*m,*n,*M,*N,ctx,mat); 20f4e70085SSatish Balay } 21f4e70085SSatish Balay 22f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat,Vec x,Vec y) 23f4e70085SSatish Balay { 24f4e70085SSatish Balay PetscErrorCode ierr = 0; 25f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[0]))(&mat,&x,&y,&ierr); 26f4e70085SSatish Balay return ierr; 27f4e70085SSatish Balay } 28f4e70085SSatish Balay 29f4e70085SSatish Balay static PetscErrorCode ourmulttranspose(Mat mat,Vec x,Vec y) 30f4e70085SSatish Balay { 31f4e70085SSatish Balay PetscErrorCode ierr = 0; 32f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[2]))(&mat,&x,&y,&ierr); 33f4e70085SSatish Balay return ierr; 34f4e70085SSatish Balay } 35f4e70085SSatish Balay 36f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat,Vec x,Vec y,Vec z) 37f4e70085SSatish Balay { 38f4e70085SSatish Balay PetscErrorCode ierr = 0; 39f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[1]))(&mat,&x,&y,&z,&ierr); 40f4e70085SSatish Balay return ierr; 41f4e70085SSatish Balay } 42f4e70085SSatish Balay 43f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat,Vec x,Vec y,Vec z) 44f4e70085SSatish Balay { 45f4e70085SSatish Balay PetscErrorCode ierr = 0; 46f4e70085SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[3]))(&mat,&x,&y,&z,&ierr); 47f4e70085SSatish Balay return ierr; 48f4e70085SSatish Balay } 49f4e70085SSatish Balay 5022612f2fSMatthew Knepley static PetscErrorCode ourgetdiagonal(Mat mat,Vec x) 5122612f2fSMatthew Knepley { 5222612f2fSMatthew Knepley PetscErrorCode ierr = 0; 5322612f2fSMatthew Knepley (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[4]))(&mat,&x,&ierr); 5422612f2fSMatthew Knepley return ierr; 5522612f2fSMatthew Knepley } 5622612f2fSMatthew Knepley 57160922c2SBarry Smith static PetscErrorCode ourdiagonalscale(Mat mat,Vec l,Vec r) 58160922c2SBarry Smith { 59160922c2SBarry Smith PetscErrorCode ierr = 0; 60160922c2SBarry Smith if (!l) { 6135b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,(Vec*)PETSC_NULL_OBJECT_Fortran,&r,&ierr); 62160922c2SBarry Smith } else if (!r) { 6335b36911SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,(Vec*)PETSC_NULL_OBJECT_Fortran,&ierr); 64160922c2SBarry Smith } else { 65160922c2SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[5]))(&mat,&l,&r,&ierr); 66160922c2SBarry Smith } 67160922c2SBarry Smith return ierr; 68160922c2SBarry Smith } 69160922c2SBarry Smith 707911a512SBarry Smith static PetscErrorCode ourgetvecs(Mat mat,Vec *l,Vec *r) 717911a512SBarry Smith { 727911a512SBarry Smith PetscErrorCode ierr = 0; 73501d9185SBarry Smith PetscInt none = -1; 747911a512SBarry Smith if (!l) { 75501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,(Vec*)&none,r,&ierr); 767911a512SBarry Smith } else if (!r) { 77501d9185SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,(Vec*)&none,&ierr); 787911a512SBarry Smith } else { 797911a512SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[7]))(&mat,l,r,&ierr); 807911a512SBarry Smith } 817911a512SBarry Smith return ierr; 827911a512SBarry Smith } 837911a512SBarry Smith 84f5a4496aSBarry Smith static PetscErrorCode ourdiagonalset(Mat mat,Vec x,InsertMode ins) 85f5a4496aSBarry Smith { 86f5a4496aSBarry Smith PetscErrorCode ierr = 0; 87f5a4496aSBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,InsertMode*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[6]))(&mat,&x,&ins,&ierr); 88f5a4496aSBarry Smith return ierr; 89f5a4496aSBarry Smith } 90f5a4496aSBarry Smith 912950f7e7SBarry Smith static PetscErrorCode ourview(Mat mat,PetscViewer v) 922950f7e7SBarry Smith { 932950f7e7SBarry Smith PetscErrorCode ierr = 0; 942950f7e7SBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscViewer*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[8]))(&mat,&v,&ierr); 952950f7e7SBarry Smith return ierr; 962950f7e7SBarry Smith } 972950f7e7SBarry Smith 983446fae8SBarry Smith static PetscErrorCode oursor(Mat mat,Vec b,PetscReal omega,MatSORType flg,PetscReal shift,PetscInt its,PetscInt lits,Vec x) 993446fae8SBarry Smith { 1003446fae8SBarry Smith PetscErrorCode ierr = 0; 1013446fae8SBarry 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); 1023446fae8SBarry Smith return ierr; 1033446fae8SBarry Smith } 1043446fae8SBarry Smith 105cdf26a31SSatish Balay static PetscErrorCode ourshift(Mat mat, PetscScalar a) 106cdf26a31SSatish Balay { 107cdf26a31SSatish Balay PetscErrorCode ierr = 0; 108cdf26a31SSatish Balay (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject)mat)->fortran_func_pointers[10]))(&mat,&a,&ierr); 109cdf26a31SSatish Balay return ierr; 110cdf26a31SSatish Balay } 111cdf26a31SSatish Balay 112*8cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matshellsetoperation_(Mat *mat,MatOperation *op,PetscErrorCode (PETSC_STDCALL *f)(Mat*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode *ierr) 113f4e70085SSatish Balay { 114e32f2f54SBarry Smith MPI_Comm comm; 115e32f2f54SBarry Smith 116e32f2f54SBarry Smith *ierr = PetscObjectGetComm((PetscObject)*mat,&comm);if (*ierr) return; 117cdf26a31SSatish Balay PetscObjectAllocateFortranPointers(*mat,11); 118f4e70085SSatish Balay if (*op == MATOP_MULT) { 119f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmult); 120f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[0] = (PetscVoidFunction)f; 121f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE) { 122f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttranspose); 123f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[2] = (PetscVoidFunction)f; 124f4e70085SSatish Balay } else if (*op == MATOP_MULT_ADD) { 125f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmultadd); 126f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[1] = (PetscVoidFunction)f; 127f4e70085SSatish Balay } else if (*op == MATOP_MULT_TRANSPOSE_ADD) { 128f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourmulttransposeadd); 129f68b968cSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[3] = (PetscVoidFunction)f; 13022612f2fSMatthew Knepley } else if (*op == MATOP_GET_DIAGONAL) { 13122612f2fSMatthew Knepley *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetdiagonal); 13222612f2fSMatthew Knepley ((PetscObject)*mat)->fortran_func_pointers[4] = (PetscVoidFunction)f; 133160922c2SBarry Smith } else if (*op == MATOP_DIAGONAL_SCALE) { 134160922c2SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalscale); 135160922c2SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[5] = (PetscVoidFunction)f; 13635153367SBarry Smith } else if (*op == MATOP_DIAGONAL_SET) { 137f5a4496aSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourdiagonalset); 138f5a4496aSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[6] = (PetscVoidFunction)f; 1397911a512SBarry Smith } else if (*op == MATOP_GET_VECS) { 1407911a512SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourgetvecs); 1417911a512SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[7] = (PetscVoidFunction)f; 1422950f7e7SBarry Smith } else if (*op == MATOP_VIEW) { 1432950f7e7SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourview); 1442950f7e7SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[8] = (PetscVoidFunction)f; 1453446fae8SBarry Smith } else if (*op == MATOP_SOR) { 1463446fae8SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)oursor); 1473446fae8SBarry Smith ((PetscObject)*mat)->fortran_func_pointers[9] = (PetscVoidFunction)f; 148cdf26a31SSatish Balay } else if (*op == MATOP_SHIFT) { 149cdf26a31SSatish Balay *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction)ourshift); 150cdf26a31SSatish Balay ((PetscObject)*mat)->fortran_func_pointers[10] = (PetscVoidFunction)f; 151f4e70085SSatish Balay } else { 152d736bfebSBarry Smith PetscError(comm,__LINE__,"MatShellSetOperation_Fortran",__FILE__,__SDIR__,PETSC_ERR_ARG_WRONG,PETSC_ERROR_INITIAL, 153f4e70085SSatish Balay "Cannot set that matrix operation"); 154f4e70085SSatish Balay *ierr = 1; 155f4e70085SSatish Balay } 156f4e70085SSatish Balay } 157f4e70085SSatish Balay 158