xref: /petsc/src/mat/impls/shell/ftn-custom/zshellf.c (revision a5b7ff6b211d69bb2e6523989771c5d3a0ae7db2)
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,
31*a5b7ff6bSBarry Smith   FORTRAN_MATOP_GET_DIAGONAL_BLOCK = 15,
32*a5b7ff6bSBarry Smith   FORTRAN_MATOP_SIZE = 16
3386686b9bSAlex Fikl };
3486686b9bSAlex Fikl 
35f4e70085SSatish Balay /*
36f4e70085SSatish Balay   The MatShell Matrix Vector product requires a C routine.
37f4e70085SSatish Balay   This C routine then calls the corresponding Fortran routine that was
38f4e70085SSatish Balay   set by the user.
39f4e70085SSatish Balay */
408cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matcreateshell_(MPI_Comm *comm, PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N, void *ctx, Mat *mat, PetscErrorCode *ierr)
41f4e70085SSatish Balay {
422e843561SJed Brown   *ierr = MatCreateShell(MPI_Comm_f2c(*(MPI_Fint*)&*comm), *m, *n, *M, *N, ctx, mat);
43f4e70085SSatish Balay }
44f4e70085SSatish Balay 
45f4e70085SSatish Balay static PetscErrorCode ourmult(Mat mat, Vec x, Vec y)
46f4e70085SSatish Balay {
47f4e70085SSatish Balay   PetscErrorCode ierr = 0;
48f4e70085SSatish Balay 
4986686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT]))(&mat, &x, &y, &ierr);
50f4e70085SSatish Balay   return ierr;
51f4e70085SSatish Balay }
52f4e70085SSatish Balay 
53f4e70085SSatish Balay static PetscErrorCode ourmultadd(Mat mat, Vec x, Vec y, Vec z)
54f4e70085SSatish Balay {
55f4e70085SSatish Balay   PetscErrorCode ierr = 0;
5686686b9bSAlex Fikl 
5786686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD]))(&mat, &x, &y, &z, &ierr);
5886686b9bSAlex Fikl   return ierr;
5986686b9bSAlex Fikl }
6086686b9bSAlex Fikl 
6186686b9bSAlex Fikl static PetscErrorCode ourmulttranspose(Mat mat, Vec x, Vec y)
6286686b9bSAlex Fikl {
6386686b9bSAlex Fikl   PetscErrorCode ierr = 0;
6486686b9bSAlex Fikl 
6586686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE]))(&mat, &x, &y, &ierr);
66f4e70085SSatish Balay   return ierr;
67f4e70085SSatish Balay }
68f4e70085SSatish Balay 
69f4e70085SSatish Balay static PetscErrorCode ourmulttransposeadd(Mat mat, Vec x, Vec y, Vec z)
70f4e70085SSatish Balay {
71f4e70085SSatish Balay   PetscErrorCode ierr = 0;
72f4e70085SSatish Balay 
7386686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD]))(&mat, &x, &y, &z, &ierr);
742950f7e7SBarry Smith   return ierr;
752950f7e7SBarry Smith }
762950f7e7SBarry Smith 
773446fae8SBarry Smith static PetscErrorCode oursor(Mat mat, Vec b, PetscReal omega, MatSORType flg, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
783446fae8SBarry Smith {
793446fae8SBarry Smith   PetscErrorCode ierr = 0;
8086686b9bSAlex Fikl 
8186686b9bSAlex 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);
8286686b9bSAlex Fikl   return ierr;
8386686b9bSAlex Fikl }
8486686b9bSAlex Fikl 
8586686b9bSAlex Fikl static PetscErrorCode ourtranspose(Mat mat, MatReuse reuse, Mat *B)
8686686b9bSAlex Fikl {
8786686b9bSAlex Fikl   PetscErrorCode ierr = 0;
8886686b9bSAlex Fikl   Mat *b = (!B ? (Mat *) PETSC_NULL_OBJECT_Fortran : B);
8986686b9bSAlex Fikl 
9086686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, MatReuse*, Mat *, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE]))(&mat, &reuse, b, &ierr);
9186686b9bSAlex Fikl   return ierr;
9286686b9bSAlex Fikl }
9386686b9bSAlex Fikl 
9486686b9bSAlex Fikl static PetscErrorCode ourgetdiagonal(Mat mat, Vec x)
9586686b9bSAlex Fikl {
9686686b9bSAlex Fikl   PetscErrorCode ierr = 0;
9786686b9bSAlex Fikl 
9886686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL]))(&mat, &x, &ierr);
9986686b9bSAlex Fikl   return ierr;
10086686b9bSAlex Fikl }
10186686b9bSAlex Fikl 
10286686b9bSAlex Fikl static PetscErrorCode ourdiagonalscale(Mat mat, Vec l, Vec r)
10386686b9bSAlex Fikl {
10486686b9bSAlex Fikl   PetscErrorCode ierr = 0;
10586686b9bSAlex Fikl   Vec *a = (!l ? (Vec*) PETSC_NULL_OBJECT_Fortran : &l);
10686686b9bSAlex Fikl   Vec *b = (!r ? (Vec*) PETSC_NULL_OBJECT_Fortran : &r);
10786686b9bSAlex Fikl 
10886686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE]))(&mat, a, b, &ierr);
10986686b9bSAlex Fikl   return ierr;
11086686b9bSAlex Fikl }
11186686b9bSAlex Fikl 
11286686b9bSAlex Fikl static PetscErrorCode ourzeroentries(Mat mat)
11386686b9bSAlex Fikl {
11486686b9bSAlex Fikl   PetscErrorCode ierr = 0;
11586686b9bSAlex Fikl 
11686686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES]))(&mat, &ierr);
11786686b9bSAlex Fikl   return ierr;
11886686b9bSAlex Fikl }
11986686b9bSAlex Fikl 
12086686b9bSAlex Fikl static PetscErrorCode ouraxpy(Mat mat, PetscScalar a, Mat X, MatStructure str)
12186686b9bSAlex Fikl {
12286686b9bSAlex Fikl   PetscErrorCode ierr = 0;
12386686b9bSAlex Fikl 
12486686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, PetscScalar*, Mat*, MatStructure*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY]))(&mat, &a, &X, &str, &ierr);
1253446fae8SBarry Smith   return ierr;
1263446fae8SBarry Smith }
1273446fae8SBarry Smith 
128cdf26a31SSatish Balay static PetscErrorCode ourshift(Mat mat, PetscScalar a)
129cdf26a31SSatish Balay {
130cdf26a31SSatish Balay   PetscErrorCode ierr = 0;
13186686b9bSAlex Fikl 
13286686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, PetscScalar*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT]))(&mat, &a, &ierr);
13386686b9bSAlex Fikl   return ierr;
13486686b9bSAlex Fikl }
13586686b9bSAlex Fikl 
13686686b9bSAlex Fikl static PetscErrorCode ourdiagonalset(Mat mat, Vec x, InsertMode ins)
13786686b9bSAlex Fikl {
13886686b9bSAlex Fikl   PetscErrorCode ierr = 0;
13986686b9bSAlex Fikl 
14086686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, InsertMode*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET]))(&mat, &x, &ins, &ierr);
14186686b9bSAlex Fikl   return ierr;
14286686b9bSAlex Fikl }
14386686b9bSAlex Fikl 
14486686b9bSAlex Fikl static PetscErrorCode ourdestroy(Mat mat)
14586686b9bSAlex Fikl {
14686686b9bSAlex Fikl   PetscErrorCode ierr = 0;
14786686b9bSAlex Fikl 
14886686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY]))(&mat, &ierr);
14986686b9bSAlex Fikl   return ierr;
15086686b9bSAlex Fikl }
15186686b9bSAlex Fikl 
15286686b9bSAlex Fikl static PetscErrorCode ourview(Mat mat, PetscViewer v)
15386686b9bSAlex Fikl {
15486686b9bSAlex Fikl   PetscErrorCode ierr = 0;
15586686b9bSAlex Fikl 
15686686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, PetscViewer*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW]))(&mat, &v, &ierr);
15786686b9bSAlex Fikl   return ierr;
15886686b9bSAlex Fikl }
15986686b9bSAlex Fikl 
16086686b9bSAlex Fikl static PetscErrorCode ourgetvecs(Mat mat, Vec *l, Vec *r)
16186686b9bSAlex Fikl {
16286686b9bSAlex Fikl   PetscErrorCode ierr = 0;
16386686b9bSAlex Fikl   Vec *a = (!l ? (Vec *) PETSC_NULL_OBJECT_Fortran : l);
16486686b9bSAlex Fikl   Vec *b = (!r ? (Vec *) PETSC_NULL_OBJECT_Fortran : r);
16586686b9bSAlex Fikl 
16686686b9bSAlex Fikl   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Vec*, Vec*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS]))(&mat, a, b, &ierr);
167cdf26a31SSatish Balay   return ierr;
168cdf26a31SSatish Balay }
169cdf26a31SSatish Balay 
170*a5b7ff6bSBarry Smith static PetscErrorCode ourgetdiagonalblock(Mat mat, Mat *l)
171*a5b7ff6bSBarry Smith {
172*a5b7ff6bSBarry Smith   PetscErrorCode ierr = 0;
173*a5b7ff6bSBarry Smith   (*(PetscErrorCode (PETSC_STDCALL *)(Mat*, Mat*, PetscErrorCode*))(((PetscObject) mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK]))(&mat, l, &ierr);
174*a5b7ff6bSBarry Smith   return ierr;
175*a5b7ff6bSBarry Smith }
176*a5b7ff6bSBarry Smith 
1778cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matshellsetoperation_(Mat *mat, MatOperation *op, PetscErrorCode (PETSC_STDCALL *f)(Mat*, Vec*, Vec*, PetscErrorCode*), PetscErrorCode *ierr)
178f4e70085SSatish Balay {
179e32f2f54SBarry Smith   MPI_Comm comm;
180e32f2f54SBarry Smith 
181e32f2f54SBarry Smith   *ierr = PetscObjectGetComm((PetscObject) *mat, &comm);if (*ierr) return;
18286686b9bSAlex Fikl   PetscObjectAllocateFortranPointers(*mat, FORTRAN_MATOP_SIZE);
18386686b9bSAlex Fikl 
18486686b9bSAlex Fikl   switch (*op) {
18586686b9bSAlex Fikl   case MATOP_MULT:
186f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourmult);
18786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT] = (PetscVoidFunction) f;
18886686b9bSAlex Fikl     break;
18986686b9bSAlex Fikl   case MATOP_MULT_ADD:
190f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourmultadd);
19186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_ADD] = (PetscVoidFunction) f;
19286686b9bSAlex Fikl     break;
19386686b9bSAlex Fikl   case MATOP_MULT_TRANSPOSE:
19486686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourmulttranspose);
19586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE] = (PetscVoidFunction) f;
19686686b9bSAlex Fikl     break;
19786686b9bSAlex Fikl   case MATOP_MULT_TRANSPOSE_ADD:
198f68b968cSBarry Smith     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourmulttransposeadd);
19986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_MULT_TRANSPOSE_ADD] = (PetscVoidFunction) f;
20086686b9bSAlex Fikl     break;
20186686b9bSAlex Fikl   case MATOP_SOR:
2023446fae8SBarry Smith     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) oursor);
20386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SOR] = (PetscVoidFunction) f;
20486686b9bSAlex Fikl     break;
20586686b9bSAlex Fikl   case MATOP_TRANSPOSE:
20686686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourtranspose);
20786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_TRANSPOSE] = (PetscVoidFunction) f;
20886686b9bSAlex Fikl     break;
20986686b9bSAlex Fikl   case MATOP_GET_DIAGONAL:
21086686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourgetdiagonal);
21186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL] = (PetscVoidFunction) f;
21286686b9bSAlex Fikl     break;
21386686b9bSAlex Fikl   case MATOP_DIAGONAL_SCALE:
21486686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourdiagonalscale);
21586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SCALE] = (PetscVoidFunction) f;
21686686b9bSAlex Fikl     break;
21786686b9bSAlex Fikl   case MATOP_ZERO_ENTRIES:
21886686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourzeroentries);
21986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_ZERO_ENTRIES] = (PetscVoidFunction) f;
22086686b9bSAlex Fikl     break;
22186686b9bSAlex Fikl   case MATOP_AXPY:
22286686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ouraxpy);
22386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_AXPY] = (PetscVoidFunction) f;
22486686b9bSAlex Fikl     break;
22586686b9bSAlex Fikl   case MATOP_SHIFT:
226cdf26a31SSatish Balay     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourshift);
22786686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_SHIFT] = (PetscVoidFunction) f;
22886686b9bSAlex Fikl     break;
22986686b9bSAlex Fikl   case MATOP_DIAGONAL_SET:
23086686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourdiagonalset);
23186686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DIAGONAL_SET] = (PetscVoidFunction) f;
23286686b9bSAlex Fikl     break;
23386686b9bSAlex Fikl   case MATOP_DESTROY:
23486686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourdestroy);
23586686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_DESTROY] = (PetscVoidFunction) f;
23686686b9bSAlex Fikl     break;
23786686b9bSAlex Fikl   case MATOP_VIEW:
23886686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourview);
23986686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_VIEW] = (PetscVoidFunction) f;
24086686b9bSAlex Fikl     break;
24186686b9bSAlex Fikl   case MATOP_GET_VECS:
24286686b9bSAlex Fikl     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourgetvecs);
24386686b9bSAlex Fikl     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_VECS] = (PetscVoidFunction) f;
24486686b9bSAlex Fikl     break;
245*a5b7ff6bSBarry Smith   case MATOP_GET_DIAGONAL_BLOCK:
246*a5b7ff6bSBarry Smith     *ierr = MatShellSetOperation(*mat, *op, (PetscVoidFunction) ourgetdiagonalblock);
247*a5b7ff6bSBarry Smith     ((PetscObject)*mat)->fortran_func_pointers[FORTRAN_MATOP_GET_DIAGONAL_BLOCK] = (PetscVoidFunction) f;
248*a5b7ff6bSBarry Smith     break;
24986686b9bSAlex Fikl   default:
25086686b9bSAlex Fikl     PetscError(comm, __LINE__, "MatShellSetOperation_Fortran", __FILE__, PETSC_ERR_ARG_WRONG, PETSC_ERROR_INITIAL, "Cannot set that matrix operation");
251f4e70085SSatish Balay     *ierr = 1;
252f4e70085SSatish Balay   }
253f4e70085SSatish Balay }
254