xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision 7f296bb328fcd4c99f2da7bfe8ba7ed8a4ebceee)
1 #include <petsc/private/ftnimpl.h>
2 #include <petscksp.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define pcasmgetsubksp_           PCASMGETSUBKSP
6   #define pcasmgetlocalsubmatrices_ PCASMGETLOCALSUBMATRICES
7   #define pcasmgetlocalsubdomains_  PCASMGETLOCALSUBDOMAINS
8   #define pcasmcreatesubdomains_    PCASMCREATESUBDOMAINS
9   #define pcasmdestroysubdomains_   PCASMDESTROYSUBDOMAINS
10   #define pcasmcreatesubdomains2d_  PCASMCREATESUBDOMAINS2D
11 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
12   #define pcasmgetsubksp_           pcasmgetsubksp
13   #define pcasmgetlocalsubmatrices_ pcasmgetlocalsubmatrices
14   #define pcasmgetlocalsubdomains_  pcasmgetlocalsubdomains
15   #define pcasmcreatesubdomains_    pcasmcreatesubdomains
16   #define pcasmdestroysubdomains_   pcasmdestroysubdomains
17   #define pcasmcreatesubdomains2d_  pcasmcreatesubdomains2d
18 #endif
19 
20 PETSC_EXTERN void pcasmcreatesubdomains_(Mat *mat, PetscInt n, F90Array1d *is, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1))
21 {
22   IS *insubs;
23 
24   CHKFORTRANNULLOBJECT(is);
25   *ierr = PCASMCreateSubdomains(*mat, n, &insubs);
26   if (*ierr) return;
27   if (insubs) *ierr = F90Array1dCreate(insubs, MPIU_FORTRANADDR, 1, n, is PETSC_F90_2PTR_PARAM(ptrd1));
28 }
29 
30 PETSC_EXTERN void pcasmgetlocalsubmatrices_(PC *pc, PetscInt *n, F90Array1d *mat, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
31 {
32   PetscInt nloc;
33   Mat     *tmat;
34 
35   CHKFORTRANNULLOBJECT(mat);
36   CHKFORTRANNULLINTEGER(n);
37   *ierr = PCASMGetLocalSubmatrices(*pc, &nloc, &tmat);
38   if (n) *n = nloc;
39   if (mat) *ierr = F90Array1dCreate(tmat, MPIU_FORTRANADDR, 1, nloc, mat PETSC_F90_2PTR_PARAM(ptrd));
40 }
41 
42 PETSC_EXTERN void pcasmgetlocalsubdomains_(PC *pc, PetscInt *n, F90Array1d *is, F90Array1d *is_local, int *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
43 {
44   PetscInt nloc;
45   IS      *tis, *tis_local;
46 
47   CHKFORTRANNULLOBJECT(is);
48   CHKFORTRANNULLOBJECT(is_local);
49   CHKFORTRANNULLINTEGER(n);
50   *ierr = PCASMGetLocalSubdomains(*pc, &nloc, &tis, &tis_local);
51   if (*ierr) return;
52   if (n) *n = nloc;
53   if (is) *ierr = F90Array1dCreate(tis, MPIU_FORTRANADDR, 1, nloc, is PETSC_F90_2PTR_PARAM(ptrd1));
54   if (*ierr) return;
55   if (is_local) *ierr = F90Array1dCreate(tis_local, MPIU_FORTRANADDR, 1, nloc, is_local PETSC_F90_2PTR_PARAM(ptrd2));
56   if (*ierr) return;
57 }
58 
59 PETSC_EXTERN void pcasmdestroysubdomains_(PetscInt *n, F90Array1d *is1, F90Array1d *is2, int *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
60 {
61   IS *isa, *isb;
62 
63   *ierr = F90Array1dAccess(is1, MPIU_FORTRANADDR, (void **)&isa PETSC_F90_2PTR_PARAM(ptrd1));
64   if (*ierr) return;
65   *ierr = F90Array1dAccess(is2, MPIU_FORTRANADDR, (void **)&isb PETSC_F90_2PTR_PARAM(ptrd2));
66   if (*ierr) return;
67   *ierr = PCASMDestroySubdomains(*n, &isa, &isb);
68   if (*ierr) return;
69   *ierr = F90Array1dDestroy(is1, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd1));
70   if (*ierr) return;
71   *ierr = F90Array1dDestroy(is2, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd2));
72   if (*ierr) return;
73 }
74 
75 PETSC_EXTERN void pcasmcreatesubdomains2d_(PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N, PetscInt *dof, PetscInt *overlap, PetscInt *Nsub, F90Array1d *is1, F90Array1d *is2, int *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
76 {
77   IS *iis, *iisl;
78 
79   *ierr = PCASMCreateSubdomains2D(*m, *n, *M, *N, *dof, *overlap, Nsub, &iis, &iisl);
80   if (*ierr) return;
81   *ierr = F90Array1dCreate(iis, MPIU_FORTRANADDR, 1, *Nsub, is1 PETSC_F90_2PTR_PARAM(ptrd1));
82   if (*ierr) return;
83   *ierr = F90Array1dCreate(iisl, MPIU_FORTRANADDR, 1, *Nsub, is2 PETSC_F90_2PTR_PARAM(ptrd2));
84   if (*ierr) return;
85 }
86 
87 PETSC_EXTERN void pcasmgetsubksp_(PC *pc, PetscInt *n_local, PetscInt *first_local, F90Array1d *ksp, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
88 {
89   KSP     *tksp;
90   PetscInt nloc, flocal;
91 
92   CHKFORTRANNULLINTEGER(n_local);
93   CHKFORTRANNULLINTEGER(first_local);
94   *ierr = PCASMGetSubKSP(*pc, &nloc, &flocal, &tksp);
95   if (n_local) *n_local = nloc;
96   if (first_local) *first_local = flocal;
97   *ierr = F90Array1dCreate(tksp, MPIU_FORTRANADDR, 1, nloc, ksp PETSC_F90_2PTR_PARAM(ptrd));
98 }
99