xref: /petsc/src/mat/impls/shell/ftn-custom/zshellf.c (revision d18b9058a302a5932c30eaa4e142cc3ba0992349)
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,
32626206dbSAlex Fikl   FORTRAN_MATOP_COPY = 16,
33626206dbSAlex Fikl   FORTRAN_MATOP_SCALE = 17,
34*d18b9058SVincent Le Chenadec   FORTRAN_MATOP_SIZE = 18,
35*d18b9058SVincent Le Chenadec   FORTRAN_MATOP_SET_RANDOM = 19
3686686b9bSAlex Fikl };
3786686b9bSAlex Fikl 
38f4e70085SSatish Balay /*
39f4e70085SSatish Balay   The MatShell Matrix Vector product requires a C routine.
40f4e70085SSatish Balay   This C routine then calls the corresponding Fortran routine that was
41f4e70085SSatish Balay   set by the user.
42f4e70085SSatish Balay */
438cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matcreateshell_(MPI_Comm *comm,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N,void *ctx,Mat *mat,PetscErrorCode *ierr)
44f4e70085SSatish Balay {
452e843561SJed Brown   *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint*)&*comm),*m,*n,*M,*N,ctx,mat);
46f4e70085SSatish Balay }
47f4e70085SSatish Balay 
48f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat,Vec x,Vec y)
49f4e70085SSatish Balay {
50f4e70085SSatish Balay   PetscErrorCode ierr = 0;
51f4e70085SSatish Balay 
5286686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT]))(&mat,&x,&y,&ierr);
53f4e70085SSatish Balay   return ierr;
54f4e70085SSatish Balay }
55f4e70085SSatish Balay 
56f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat,Vec x,Vec y,Vec z)
57f4e70085SSatish Balay {
58f4e70085SSatish Balay   PetscErrorCode ierr = 0;
5986686b9bSAlex Fikl 
6086686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD]))(&mat,&x,&y,&z,&ierr);
6186686b9bSAlex Fikl   return ierr;
6286686b9bSAlex Fikl }
6386686b9bSAlex Fikl 
6486686b9bSAlex Fikl static PetscErrorCode ourmulttranspose(Mat mat,Vec x,Vec y)
6586686b9bSAlex Fikl {
6686686b9bSAlex Fikl   PetscErrorCode ierr = 0;
6786686b9bSAlex Fikl 
6886686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE]))(&mat,&x,&y,&ierr);
69f4e70085SSatish Balay   return ierr;
70f4e70085SSatish Balay }
71f4e70085SSatish Balay 
72f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat,Vec x,Vec y,Vec z)
73f4e70085SSatish Balay {
74f4e70085SSatish Balay   PetscErrorCode ierr = 0;
75f4e70085SSatish Balay 
7686686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD]))(&mat,&x,&y,&z,&ierr);
772950f7e7SBarry Smith   return ierr;
782950f7e7SBarry Smith }
792950f7e7SBarry Smith 
803446fae8SBarry Smith static PetscErrorCode oursor(Mat mat,Vec b,PetscReal omega,MatSORType flg,PetscReal shift,PetscInt its,PetscInt lits,Vec x)
813446fae8SBarry Smith {
823446fae8SBarry Smith   PetscErrorCode ierr = 0;
8386686b9bSAlex Fikl 
8486686b9bSAlex 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);
8586686b9bSAlex Fikl   return ierr;
8686686b9bSAlex Fikl }
8786686b9bSAlex Fikl 
8886686b9bSAlex Fikl static PetscErrorCode ourtranspose(Mat mat,MatReuse reuse,Mat *B)
8986686b9bSAlex Fikl {
9086686b9bSAlex Fikl   PetscErrorCode ierr = 0;
9186686b9bSAlex Fikl   Mat *b = (!B ? (Mat *) PETSC_NULL_OBJECT_Fortran : B);
9286686b9bSAlex Fikl 
9386686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,MatReuse*,Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE]))(&mat,&reuse,b,&ierr);
9486686b9bSAlex Fikl   return ierr;
9586686b9bSAlex Fikl }
9686686b9bSAlex Fikl 
9786686b9bSAlex Fikl static PetscErrorCode ourgetdiagonal(Mat mat,Vec x)
9886686b9bSAlex Fikl {
9986686b9bSAlex Fikl   PetscErrorCode ierr = 0;
10086686b9bSAlex Fikl 
10186686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL]))(&mat,&x,&ierr);
10286686b9bSAlex Fikl   return ierr;
10386686b9bSAlex Fikl }
10486686b9bSAlex Fikl 
10586686b9bSAlex Fikl static PetscErrorCode ourdiagonalscale(Mat mat,Vec l,Vec r)
10686686b9bSAlex Fikl {
10786686b9bSAlex Fikl   PetscErrorCode ierr = 0;
10886686b9bSAlex Fikl   Vec *a = (!l ? (Vec*) PETSC_NULL_OBJECT_Fortran : &l);
10986686b9bSAlex Fikl   Vec *b = (!r ? (Vec*) PETSC_NULL_OBJECT_Fortran : &r);
11086686b9bSAlex Fikl 
11186686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE]))(&mat,a,b,&ierr);
11286686b9bSAlex Fikl   return ierr;
11386686b9bSAlex Fikl }
11486686b9bSAlex Fikl 
11586686b9bSAlex Fikl static PetscErrorCode ourzeroentries(Mat mat)
11686686b9bSAlex Fikl {
11786686b9bSAlex Fikl   PetscErrorCode ierr = 0;
11886686b9bSAlex Fikl 
11986686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES]))(&mat,&ierr);
12086686b9bSAlex Fikl   return ierr;
12186686b9bSAlex Fikl }
12286686b9bSAlex Fikl 
12386686b9bSAlex Fikl static PetscErrorCode ouraxpy(Mat mat,PetscScalar a,Mat X,MatStructure str)
12486686b9bSAlex Fikl {
12586686b9bSAlex Fikl   PetscErrorCode ierr = 0;
12686686b9bSAlex Fikl 
12786686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,Mat*,MatStructure*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY]))(&mat,&a,&X,&str,&ierr);
1283446fae8SBarry Smith   return ierr;
1293446fae8SBarry Smith }
1303446fae8SBarry Smith 
131cdf26a31SSatish Balay static PetscErrorCode ourshift(Mat mat,PetscScalar a)
132cdf26a31SSatish Balay {
133cdf26a31SSatish Balay   PetscErrorCode ierr = 0;
13486686b9bSAlex Fikl 
13586686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT]))(&mat,&a,&ierr);
13686686b9bSAlex Fikl   return ierr;
13786686b9bSAlex Fikl }
13886686b9bSAlex Fikl 
13986686b9bSAlex Fikl static PetscErrorCode ourdiagonalset(Mat mat,Vec x,InsertMode ins)
14086686b9bSAlex Fikl {
14186686b9bSAlex Fikl   PetscErrorCode ierr = 0;
14286686b9bSAlex Fikl 
14386686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,InsertMode*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET]))(&mat,&x,&ins,&ierr);
14486686b9bSAlex Fikl   return ierr;
14586686b9bSAlex Fikl }
14686686b9bSAlex Fikl 
14786686b9bSAlex Fikl static PetscErrorCode ourdestroy(Mat mat)
14886686b9bSAlex Fikl {
14986686b9bSAlex Fikl   PetscErrorCode ierr = 0;
15086686b9bSAlex Fikl 
15186686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY]))(&mat,&ierr);
15286686b9bSAlex Fikl   return ierr;
15386686b9bSAlex Fikl }
15486686b9bSAlex Fikl 
15586686b9bSAlex Fikl static PetscErrorCode ourview(Mat mat,PetscViewer v)
15686686b9bSAlex Fikl {
15786686b9bSAlex Fikl   PetscErrorCode ierr = 0;
15886686b9bSAlex Fikl 
15986686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscViewer*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW]))(&mat,&v,&ierr);
16086686b9bSAlex Fikl   return ierr;
16186686b9bSAlex Fikl }
16286686b9bSAlex Fikl 
16386686b9bSAlex Fikl static PetscErrorCode ourgetvecs(Mat mat,Vec *l,Vec *r)
16486686b9bSAlex Fikl {
16586686b9bSAlex Fikl   PetscErrorCode ierr = 0;
16686686b9bSAlex Fikl   Vec *a = (!l ? (Vec *) PETSC_NULL_OBJECT_Fortran : l);
16786686b9bSAlex Fikl   Vec *b = (!r ? (Vec *) PETSC_NULL_OBJECT_Fortran : r);
16886686b9bSAlex Fikl 
16986686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Vec*,Vec*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS]))(&mat,a,b,&ierr);
170cdf26a31SSatish Balay   return ierr;
171cdf26a31SSatish Balay }
172cdf26a31SSatish Balay 
173a5b7ff6bSBarry Smith static PetscErrorCode ourgetdiagonalblock(Mat mat,Mat *l)
174a5b7ff6bSBarry Smith {
175a5b7ff6bSBarry Smith   PetscErrorCode ierr = 0;
176a5b7ff6bSBarry Smith   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Mat*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK]))(&mat,l,&ierr);
177a5b7ff6bSBarry Smith   return ierr;
178a5b7ff6bSBarry Smith }
179a5b7ff6bSBarry Smith 
180626206dbSAlex Fikl static PetscErrorCode ourcopy(Mat mat,Mat B,MatStructure str)
181626206dbSAlex Fikl {
182626206dbSAlex Fikl   PetscErrorCode ierr = 0;
183626206dbSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,Mat*,MatStructure*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_COPY]))(&mat,&B,&str,&ierr);
184626206dbSAlex Fikl   return ierr;
185626206dbSAlex Fikl }
186626206dbSAlex Fikl 
187626206dbSAlex Fikl static PetscErrorCode ourscale(Mat mat,PetscScalar a)
188626206dbSAlex Fikl {
189626206dbSAlex Fikl   PetscErrorCode ierr = 0;
190626206dbSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscScalar*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SCALE]))(&mat,&a,&ierr);
191626206dbSAlex Fikl   return ierr;
192626206dbSAlex Fikl }
193626206dbSAlex Fikl 
194*d18b9058SVincent Le Chenadec static PetscErrorCode oursetrandom(Mat mat,PetscRandom ctx)
195*d18b9058SVincent Le Chenadec {
196*d18b9058SVincent Le Chenadec   PetscErrorCode ierr = 0;
197*d18b9058SVincent Le Chenadec   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*,PetscRandom*,PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SET_RANDOM]))(&mat,&ctx,&ierr);
198*d18b9058SVincent Le Chenadec   return ierr;
199*d18b9058SVincent Le Chenadec }
200*d18b9058SVincent Le Chenadec 
2018cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matshellsetoperation_(Mat *mat,MatOperation *op,PetscErrorCode (PETSC_STDCALL *f)(Mat*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode *ierr)
202f4e70085SSatish Balay {
203e32f2f54SBarry Smith   MPI_Comm comm;
204e32f2f54SBarry Smith 
205e32f2f54SBarry Smith   *ierr = PetscObjectGetComm((PetscObject) *mat,&comm);if (*ierr) return;
20686686b9bSAlex Fikl   PetscObjectAllocateFortranPointers(*mat,FORTRAN_MATOP_SIZE);
20786686b9bSAlex Fikl 
20886686b9bSAlex Fikl   switch (*op) {
20986686b9bSAlex Fikl   case MATOP_MULT:
210f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmult);
21186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT] = (PetscVoidFunction) f;
21286686b9bSAlex Fikl     break;
21386686b9bSAlex Fikl   case MATOP_MULT_ADD:
214f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmultadd);
21586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD] = (PetscVoidFunction) f;
21686686b9bSAlex Fikl     break;
21786686b9bSAlex Fikl   case MATOP_MULT_TRANSPOSE:
21886686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmulttranspose);
21986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE] = (PetscVoidFunction) f;
22086686b9bSAlex Fikl     break;
22186686b9bSAlex Fikl   case MATOP_MULT_TRANSPOSE_ADD:
222f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourmulttransposeadd);
22386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD] = (PetscVoidFunction) f;
22486686b9bSAlex Fikl     break;
22586686b9bSAlex Fikl   case MATOP_SOR:
2263446fae8SBarry Smith     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) oursor);
22786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SOR] = (PetscVoidFunction) f;
22886686b9bSAlex Fikl     break;
22986686b9bSAlex Fikl   case MATOP_TRANSPOSE:
23086686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourtranspose);
23186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE] = (PetscVoidFunction) f;
23286686b9bSAlex Fikl     break;
23386686b9bSAlex Fikl   case MATOP_GET_DIAGONAL:
23486686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetdiagonal);
23586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL] = (PetscVoidFunction) f;
23686686b9bSAlex Fikl     break;
23786686b9bSAlex Fikl   case MATOP_DIAGONAL_SCALE:
23886686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdiagonalscale);
23986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE] = (PetscVoidFunction) f;
24086686b9bSAlex Fikl     break;
24186686b9bSAlex Fikl   case MATOP_ZERO_ENTRIES:
24286686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourzeroentries);
24386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES] = (PetscVoidFunction) f;
24486686b9bSAlex Fikl     break;
24586686b9bSAlex Fikl   case MATOP_AXPY:
24686686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ouraxpy);
24786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY] = (PetscVoidFunction) f;
24886686b9bSAlex Fikl     break;
24986686b9bSAlex Fikl   case MATOP_SHIFT:
250cdf26a31SSatish Balay     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourshift);
25186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT] = (PetscVoidFunction) f;
25286686b9bSAlex Fikl     break;
25386686b9bSAlex Fikl   case MATOP_DIAGONAL_SET:
25486686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdiagonalset);
25586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET] = (PetscVoidFunction) f;
25686686b9bSAlex Fikl     break;
25786686b9bSAlex Fikl   case MATOP_DESTROY:
25886686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourdestroy);
25986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY] = (PetscVoidFunction) f;
26086686b9bSAlex Fikl     break;
26186686b9bSAlex Fikl   case MATOP_VIEW:
26286686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourview);
26386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW] = (PetscVoidFunction) f;
26486686b9bSAlex Fikl     break;
26586686b9bSAlex Fikl   case MATOP_GET_VECS:
26686686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetvecs);
26786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS] = (PetscVoidFunction) f;
26886686b9bSAlex Fikl     break;
269a5b7ff6bSBarry Smith   case MATOP_GET_DIAGONAL_BLOCK:
270a5b7ff6bSBarry Smith     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourgetdiagonalblock);
271a5b7ff6bSBarry Smith     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK] = (PetscVoidFunction) f;
272a5b7ff6bSBarry Smith     break;
273626206dbSAlex Fikl   case MATOP_COPY:
274626206dbSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourcopy);
275626206dbSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_COPY] = (PetscVoidFunction) f;
276626206dbSAlex Fikl     break;
277626206dbSAlex Fikl   case MATOP_SCALE:
278626206dbSAlex Fikl     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) ourscale);
279626206dbSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SCALE] = (PetscVoidFunction) f;
280626206dbSAlex Fikl     break;
281*d18b9058SVincent Le Chenadec   case MATOP_SET_RANDOM:
282*d18b9058SVincent Le Chenadec     *ierr = MatShellSetOperation(*mat,*op,(PetscVoidFunction) oursetrandom);
283*d18b9058SVincent Le Chenadec     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SET_RANDOM] = (PetscVoidFunction) f;
284*d18b9058SVincent Le Chenadec     break;
28586686b9bSAlex Fikl   default:
28686686b9bSAlex Fikl     PetscError(comm,__LINE__,"MatShellSetOperation_Fortran",__FILE__,PETSC_ERR_ARG_WRONG,PETSC_ERROR_INITIAL,"Cannot set that matrix operation");
287f4e70085SSatish Balay     *ierr = 1;
288f4e70085SSatish Balay   }
289f4e70085SSatish Balay }
290