1af0996ceSBarry Smith #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 1286686b9bSAlex Fikl /** 1386686b9bSAlex Fikl * Subset of MatOperation that is supported by the Fortran wrappers. 1486686b9bSAlex Fikl */ 1586686b9bSAlex Fikl enum FortranMatOperation { 1686686b9bSAlex Fikl FORTRAN_MATOP_MULT = 0, 1786686b9bSAlex Fikl FORTRAN_MATOP_MULT_ADD = 1, 1886686b9bSAlex Fikl FORTRAN_MATOP_MULT_TRANSPOSE = 2, 1986686b9bSAlex Fikl FORTRAN_MATOP_MULT_TRANSPOSE_ADD = 3, 2086686b9bSAlex Fikl FORTRAN_MATOP_SOR = 4, 2186686b9bSAlex Fikl FORTRAN_MATOP_TRANSPOSE = 5, 2286686b9bSAlex Fikl FORTRAN_MATOP_GET_DIAGONAL = 6, 2386686b9bSAlex Fikl FORTRAN_MATOP_DIAGONAL_SCALE = 7, 2486686b9bSAlex Fikl FORTRAN_MATOP_ZERO_ENTRIES = 8, 2586686b9bSAlex Fikl FORTRAN_MATOP_AXPY = 9, 2686686b9bSAlex Fikl FORTRAN_MATOP_SHIFT = 10, 2786686b9bSAlex Fikl FORTRAN_MATOP_DIAGONAL_SET = 11, 2886686b9bSAlex Fikl FORTRAN_MATOP_DESTROY = 12, 2986686b9bSAlex Fikl FORTRAN_MATOP_VIEW = 13, 3086686b9bSAlex Fikl FORTRAN_MATOP_GET_VECS = 14, 31a5b7ff6bSBarry Smith FORTRAN_MATOP_GET_DIAGONAL_BLOCK = 15, 32*626206dbSAlex Fikl FORTRAN_MATOP_COPY = 16, 33*626206dbSAlex Fikl FORTRAN_MATOP_SCALE = 17, 34*626206dbSAlex Fikl FORTRAN_MATOP_SIZE = 18 3586686b9bSAlex Fikl }; 3686686b9bSAlex Fikl 37f4e70085SSatish Balay /* 38f4e70085SSatish Balay The MatShell Matrix Vector product requires a C routine. 39f4e70085SSatish Balay This C routine then calls the corresponding Fortran routine that was 40f4e70085SSatish Balay set by the user. 41f4e70085SSatish Balay */ 428cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matcreateshell_(MPI_Comm *comm,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N,void *ctx,Mat *mat,PetscErrorCode *ierr) 43f4e70085SSatish Balay { 442e843561SJed Brown *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint*)&*comm),*m,*n,*M,*N,ctx,mat); 45f4e70085SSatish Balay } 46f4e70085SSatish Balay 47f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat,Vec x,Vec y) 48f4e70085SSatish Balay { 49f4e70085SSatish Balay PetscErrorCode ierr = 0; 50f4e70085SSatish Balay 5186686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT]))(&mat,&x,&y,&ierr); 52f4e70085SSatish Balay return ierr; 53f4e70085SSatish Balay } 54f4e70085SSatish Balay 55f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat,Vec x,Vec y,Vec z) 56f4e70085SSatish Balay { 57f4e70085SSatish Balay PetscErrorCode ierr = 0; 5886686b9bSAlex Fikl 5986686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD]))(&mat,&x,&y,&z,&ierr); 6086686b9bSAlex Fikl return ierr; 6186686b9bSAlex Fikl } 6286686b9bSAlex Fikl 6386686b9bSAlex Fikl static PetscErrorCode ourmulttranspose(Mat mat,Vec x,Vec y) 6486686b9bSAlex Fikl { 6586686b9bSAlex Fikl PetscErrorCode ierr = 0; 6686686b9bSAlex Fikl 6786686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE]))(&mat,&x,&y,&ierr); 68f4e70085SSatish Balay return ierr; 69f4e70085SSatish Balay } 70f4e70085SSatish Balay 71f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat,Vec x,Vec y,Vec z) 72f4e70085SSatish Balay { 73f4e70085SSatish Balay PetscErrorCode ierr = 0; 74f4e70085SSatish Balay 7586686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD]))(&mat,&x,&y,&z,&ierr); 762950f7e7SBarry Smith return ierr; 772950f7e7SBarry Smith } 782950f7e7SBarry Smith 793446fae8SBarry Smith static PetscErrorCode oursor(Mat mat,Vec b,PetscReal omega,MatSORType flg,PetscReal shift,PetscInt its,PetscInt lits,Vec x) 803446fae8SBarry Smith { 813446fae8SBarry Smith PetscErrorCode ierr = 0; 8286686b9bSAlex Fikl 8386686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscReal*,MatSORType*,PetscReal*,PetscInt*,PetscInt*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SOR]))(&mat,&b,&omega,&flg,&shift,&its,&lits,&x,&ierr); 8486686b9bSAlex Fikl return ierr; 8586686b9bSAlex Fikl } 8686686b9bSAlex Fikl 8786686b9bSAlex Fikl static PetscErrorCode ourtranspose(Mat mat,MatReuse reuse,Mat *B) 8886686b9bSAlex Fikl { 8986686b9bSAlex Fikl PetscErrorCode ierr = 0; 9086686b9bSAlex Fikl Mat *b = (!B ? (Mat *) PETSC_NULL_OBJECT_Fortran : B); 9186686b9bSAlex Fikl 9286686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,MatReuse*,Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE]))(&mat,&reuse,b,&ierr); 9386686b9bSAlex Fikl return ierr; 9486686b9bSAlex Fikl } 9586686b9bSAlex Fikl 9686686b9bSAlex Fikl static PetscErrorCode ourgetdiagonal(Mat mat,Vec x) 9786686b9bSAlex Fikl { 9886686b9bSAlex Fikl PetscErrorCode ierr = 0; 9986686b9bSAlex Fikl 10086686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL]))(&mat,&x,&ierr); 10186686b9bSAlex Fikl return ierr; 10286686b9bSAlex Fikl } 10386686b9bSAlex Fikl 10486686b9bSAlex Fikl static PetscErrorCode ourdiagonalscale(Mat mat,Vec l,Vec r) 10586686b9bSAlex Fikl { 10686686b9bSAlex Fikl PetscErrorCode ierr = 0; 10786686b9bSAlex Fikl Vec *a = (!l ? (Vec*) PETSC_NULL_OBJECT_Fortran : &l); 10886686b9bSAlex Fikl Vec *b = (!r ? (Vec*) PETSC_NULL_OBJECT_Fortran : &r); 10986686b9bSAlex Fikl 11086686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE]))(&mat,a,b,&ierr); 11186686b9bSAlex Fikl return ierr; 11286686b9bSAlex Fikl } 11386686b9bSAlex Fikl 11486686b9bSAlex Fikl static PetscErrorCode ourzeroentries(Mat mat) 11586686b9bSAlex Fikl { 11686686b9bSAlex Fikl PetscErrorCode ierr = 0; 11786686b9bSAlex Fikl 11886686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES]))(&mat,&ierr); 11986686b9bSAlex Fikl return ierr; 12086686b9bSAlex Fikl } 12186686b9bSAlex Fikl 12286686b9bSAlex Fikl static PetscErrorCode ouraxpy(Mat mat,PetscScalar a,Mat X,MatStructure str) 12386686b9bSAlex Fikl { 12486686b9bSAlex Fikl PetscErrorCode ierr = 0; 12586686b9bSAlex Fikl 12686686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,Mat*,MatStructure*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY]))(&mat,&a,&X,&str,&ierr); 1273446fae8SBarry Smith return ierr; 1283446fae8SBarry Smith } 1293446fae8SBarry Smith 130cdf26a31SSatish Balay static PetscErrorCode ourshift(Mat mat,PetscScalar a) 131cdf26a31SSatish Balay { 132cdf26a31SSatish Balay PetscErrorCode ierr = 0; 13386686b9bSAlex Fikl 13486686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT]))(&mat,&a,&ierr); 13586686b9bSAlex Fikl return ierr; 13686686b9bSAlex Fikl } 13786686b9bSAlex Fikl 13886686b9bSAlex Fikl static PetscErrorCode ourdiagonalset(Mat mat,Vec x,InsertMode ins) 13986686b9bSAlex Fikl { 14086686b9bSAlex Fikl PetscErrorCode ierr = 0; 14186686b9bSAlex Fikl 14286686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,InsertMode*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET]))(&mat,&x,&ins,&ierr); 14386686b9bSAlex Fikl return ierr; 14486686b9bSAlex Fikl } 14586686b9bSAlex Fikl 14686686b9bSAlex Fikl static PetscErrorCode ourdestroy(Mat mat) 14786686b9bSAlex Fikl { 14886686b9bSAlex Fikl PetscErrorCode ierr = 0; 14986686b9bSAlex Fikl 15086686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY]))(&mat,&ierr); 15186686b9bSAlex Fikl return ierr; 15286686b9bSAlex Fikl } 15386686b9bSAlex Fikl 15486686b9bSAlex Fikl static PetscErrorCode ourview(Mat mat,PetscViewer v) 15586686b9bSAlex Fikl { 15686686b9bSAlex Fikl PetscErrorCode ierr = 0; 15786686b9bSAlex Fikl 15886686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscViewer*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW]))(&mat,&v,&ierr); 15986686b9bSAlex Fikl return ierr; 16086686b9bSAlex Fikl } 16186686b9bSAlex Fikl 16286686b9bSAlex Fikl static PetscErrorCode ourgetvecs(Mat mat,Vec *l,Vec *r) 16386686b9bSAlex Fikl { 16486686b9bSAlex Fikl PetscErrorCode ierr = 0; 16586686b9bSAlex Fikl Vec *a = (!l ? (Vec *) PETSC_NULL_OBJECT_Fortran : l); 16686686b9bSAlex Fikl Vec *b = (!r ? (Vec *) PETSC_NULL_OBJECT_Fortran : r); 16786686b9bSAlex Fikl 16886686b9bSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS]))(&mat,a,b,&ierr); 169cdf26a31SSatish Balay return ierr; 170cdf26a31SSatish Balay } 171cdf26a31SSatish Balay 172a5b7ff6bSBarry Smith static PetscErrorCode ourgetdiagonalblock(Mat mat,Mat *l) 173a5b7ff6bSBarry Smith { 174a5b7ff6bSBarry Smith PetscErrorCode ierr = 0; 175a5b7ff6bSBarry Smith (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK]))(&mat,l,&ierr); 176a5b7ff6bSBarry Smith return ierr; 177a5b7ff6bSBarry Smith } 178a5b7ff6bSBarry Smith 179*626206dbSAlex Fikl static PetscErrorCode ourcopy(Mat mat,Mat B,MatStructure str) 180*626206dbSAlex Fikl { 181*626206dbSAlex Fikl PetscErrorCode ierr = 0; 182*626206dbSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Mat*,MatStructure*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_COPY]))(&mat,&B,&str,&ierr); 183*626206dbSAlex Fikl return ierr; 184*626206dbSAlex Fikl } 185*626206dbSAlex Fikl 186*626206dbSAlex Fikl static PetscErrorCode ourscale(Mat mat,PetscScalar a) 187*626206dbSAlex Fikl { 188*626206dbSAlex Fikl PetscErrorCode ierr = 0; 189*626206dbSAlex Fikl (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SCALE]))(&mat,&a,&ierr); 190*626206dbSAlex Fikl return ierr; 191*626206dbSAlex Fikl } 192*626206dbSAlex Fikl 1938cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matshellsetoperation_(Mat *mat,MatOperation *op,PetscErrorCode (PETSC_STDCALL *f)(Mat*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode *ierr) 194f4e70085SSatish Balay { 195e32f2f54SBarry Smith MPI_Comm comm; 196e32f2f54SBarry Smith 197e32f2f54SBarry Smith *ierr = PetscObjectGetComm((PetscObject) *mat,&comm);if (*ierr) return; 19886686b9bSAlex Fikl PetscObjectAllocateFortranPointers(*mat,FORTRAN_MATOP_SIZE); 19986686b9bSAlex Fikl 20086686b9bSAlex Fikl switch (*op) { 20186686b9bSAlex Fikl case MATOP_MULT: 202f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmult); 20386686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT] = (PetscVoidFunction) f; 20486686b9bSAlex Fikl break; 20586686b9bSAlex Fikl case MATOP_MULT_ADD: 206f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmultadd); 20786686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD] = (PetscVoidFunction) f; 20886686b9bSAlex Fikl break; 20986686b9bSAlex Fikl case MATOP_MULT_TRANSPOSE: 21086686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmulttranspose); 21186686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE] = (PetscVoidFunction) f; 21286686b9bSAlex Fikl break; 21386686b9bSAlex Fikl case MATOP_MULT_TRANSPOSE_ADD: 214f68b968cSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmulttransposeadd); 21586686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD] = (PetscVoidFunction) f; 21686686b9bSAlex Fikl break; 21786686b9bSAlex Fikl case MATOP_SOR: 2183446fae8SBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) oursor); 21986686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SOR] = (PetscVoidFunction) f; 22086686b9bSAlex Fikl break; 22186686b9bSAlex Fikl case MATOP_TRANSPOSE: 22286686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourtranspose); 22386686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE] = (PetscVoidFunction) f; 22486686b9bSAlex Fikl break; 22586686b9bSAlex Fikl case MATOP_GET_DIAGONAL: 22686686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetdiagonal); 22786686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL] = (PetscVoidFunction) f; 22886686b9bSAlex Fikl break; 22986686b9bSAlex Fikl case MATOP_DIAGONAL_SCALE: 23086686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdiagonalscale); 23186686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE] = (PetscVoidFunction) f; 23286686b9bSAlex Fikl break; 23386686b9bSAlex Fikl case MATOP_ZERO_ENTRIES: 23486686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourzeroentries); 23586686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES] = (PetscVoidFunction) f; 23686686b9bSAlex Fikl break; 23786686b9bSAlex Fikl case MATOP_AXPY: 23886686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ouraxpy); 23986686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY] = (PetscVoidFunction) f; 24086686b9bSAlex Fikl break; 24186686b9bSAlex Fikl case MATOP_SHIFT: 242cdf26a31SSatish Balay *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourshift); 24386686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT] = (PetscVoidFunction) f; 24486686b9bSAlex Fikl break; 24586686b9bSAlex Fikl case MATOP_DIAGONAL_SET: 24686686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdiagonalset); 24786686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET] = (PetscVoidFunction) f; 24886686b9bSAlex Fikl break; 24986686b9bSAlex Fikl case MATOP_DESTROY: 25086686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdestroy); 25186686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY] = (PetscVoidFunction) f; 25286686b9bSAlex Fikl break; 25386686b9bSAlex Fikl case MATOP_VIEW: 25486686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourview); 25586686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW] = (PetscVoidFunction) f; 25686686b9bSAlex Fikl break; 25786686b9bSAlex Fikl case MATOP_GET_VECS: 25886686b9bSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetvecs); 25986686b9bSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS] = (PetscVoidFunction) f; 26086686b9bSAlex Fikl break; 261a5b7ff6bSBarry Smith case MATOP_GET_DIAGONAL_BLOCK: 262a5b7ff6bSBarry Smith *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetdiagonalblock); 263a5b7ff6bSBarry Smith ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK] = (PetscVoidFunction) f; 264a5b7ff6bSBarry Smith break; 265*626206dbSAlex Fikl case MATOP_COPY: 266*626206dbSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourcopy); 267*626206dbSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_COPY] = (PetscVoidFunction) f; 268*626206dbSAlex Fikl break; 269*626206dbSAlex Fikl case MATOP_SCALE: 270*626206dbSAlex Fikl *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourscale); 271*626206dbSAlex Fikl ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SCALE] = (PetscVoidFunction) f; 272*626206dbSAlex Fikl break; 27386686b9bSAlex Fikl default: 27486686b9bSAlex Fikl PetscError(comm,__LINE__,"MatShellSetOperation_Fortran",__FILE__,PETSC_ERR_ARG_WRONG,PETSC_ERROR_INITIAL,"Cannot set that matrix operation"); 275f4e70085SSatish Balay *ierr = 1; 276f4e70085SSatish Balay } 277f4e70085SSatish Balay } 278