16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2c6db04a5SJed Brown #include <petscksp.h>
3e54e4138SSatish Balay
4e54e4138SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
5ce78bad3SBarry Smith #define pcasmgetsubksp_ PCASMGETSUBKSP
6*36083efbSBarry Smith #define pcasmrestoresubksp_ PCASMRESTORESUBKSP
7e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_ PCASMGETLOCALSUBMATRICES
8e54e4138SSatish Balay #define pcasmgetlocalsubdomains_ PCASMGETLOCALSUBDOMAINS
9ef356e52SBarry Smith #define pcasmcreatesubdomains_ PCASMCREATESUBDOMAINS
10ef356e52SBarry Smith #define pcasmdestroysubdomains_ PCASMDESTROYSUBDOMAINS
116141accfSBarry Smith #define pcasmcreatesubdomains2d_ PCASMCREATESUBDOMAINS2D
12e54e4138SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
13ce78bad3SBarry Smith #define pcasmgetsubksp_ pcasmgetsubksp
14*36083efbSBarry Smith #define pcasmrestoresubksp_ pcasmrestoresubksp
15e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_ pcasmgetlocalsubmatrices
16e54e4138SSatish Balay #define pcasmgetlocalsubdomains_ pcasmgetlocalsubdomains
17ef356e52SBarry Smith #define pcasmcreatesubdomains_ pcasmcreatesubdomains
18ef356e52SBarry Smith #define pcasmdestroysubdomains_ pcasmdestroysubdomains
196141accfSBarry Smith #define pcasmcreatesubdomains2d_ pcasmcreatesubdomains2d
20e54e4138SSatish Balay #endif
21e54e4138SSatish Balay
pcasmcreatesubdomains_(Mat * mat,PetscInt n,F90Array1d * is,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd1))22ce78bad3SBarry Smith PETSC_EXTERN void pcasmcreatesubdomains_(Mat *mat, PetscInt n, F90Array1d *is, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1))
236141accfSBarry Smith {
24ef356e52SBarry Smith IS *insubs;
25ef356e52SBarry Smith
26ce78bad3SBarry Smith CHKFORTRANNULLOBJECT(is);
27ce78bad3SBarry Smith *ierr = PCASMCreateSubdomains(*mat, n, &insubs);
285975b3b6SBarry Smith if (*ierr) return;
29ce78bad3SBarry Smith if (insubs) *ierr = F90Array1dCreate(insubs, MPIU_FORTRANADDR, 1, n, is PETSC_F90_2PTR_PARAM(ptrd1));
30ef356e52SBarry Smith }
31ef356e52SBarry Smith
pcasmgetlocalsubmatrices_(PC * pc,PetscInt * n,F90Array1d * mat,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))32ce78bad3SBarry Smith PETSC_EXTERN void pcasmgetlocalsubmatrices_(PC *pc, PetscInt *n, F90Array1d *mat, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
33ef356e52SBarry Smith {
34ce78bad3SBarry Smith PetscInt nloc;
35e54e4138SSatish Balay Mat *tmat;
36ce78bad3SBarry Smith
37e54e4138SSatish Balay CHKFORTRANNULLOBJECT(mat);
38e54e4138SSatish Balay CHKFORTRANNULLINTEGER(n);
39e54e4138SSatish Balay *ierr = PCASMGetLocalSubmatrices(*pc, &nloc, &tmat);
40e54e4138SSatish Balay if (n) *n = nloc;
41ce78bad3SBarry Smith if (mat) *ierr = F90Array1dCreate(tmat, MPIU_FORTRANADDR, 1, nloc, mat PETSC_F90_2PTR_PARAM(ptrd));
42e54e4138SSatish Balay }
43ce78bad3SBarry Smith
pcasmgetlocalsubdomains_(PC * pc,PetscInt * n,F90Array1d * is,F90Array1d * is_local,int * ierr PETSC_F90_2PTR_PROTO (ptrd1)PETSC_F90_2PTR_PROTO (ptrd2))44ce78bad3SBarry Smith 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))
45e54e4138SSatish Balay {
46ce78bad3SBarry Smith PetscInt nloc;
472b691e39SMatthew Knepley IS *tis, *tis_local;
48ce78bad3SBarry Smith
49e54e4138SSatish Balay CHKFORTRANNULLOBJECT(is);
502b691e39SMatthew Knepley CHKFORTRANNULLOBJECT(is_local);
51e54e4138SSatish Balay CHKFORTRANNULLINTEGER(n);
522b691e39SMatthew Knepley *ierr = PCASMGetLocalSubdomains(*pc, &nloc, &tis, &tis_local);
53ce78bad3SBarry Smith if (*ierr) return;
54e54e4138SSatish Balay if (n) *n = nloc;
55ce78bad3SBarry Smith if (is) *ierr = F90Array1dCreate(tis, MPIU_FORTRANADDR, 1, nloc, is PETSC_F90_2PTR_PARAM(ptrd1));
56ce78bad3SBarry Smith if (*ierr) return;
57ce78bad3SBarry Smith if (is_local) *ierr = F90Array1dCreate(tis_local, MPIU_FORTRANADDR, 1, nloc, is_local PETSC_F90_2PTR_PARAM(ptrd2));
58ce78bad3SBarry Smith if (*ierr) return;
59e54e4138SSatish Balay }
60ce78bad3SBarry Smith
pcasmdestroysubdomains_(PetscInt * n,F90Array1d * is1,F90Array1d * is2,int * ierr PETSC_F90_2PTR_PROTO (ptrd1)PETSC_F90_2PTR_PROTO (ptrd2))61ce78bad3SBarry Smith PETSC_EXTERN void pcasmdestroysubdomains_(PetscInt *n, F90Array1d *is1, F90Array1d *is2, int *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
62ce78bad3SBarry Smith {
63ce78bad3SBarry Smith IS *isa, *isb;
64ce78bad3SBarry Smith
65ce78bad3SBarry Smith *ierr = F90Array1dAccess(is1, MPIU_FORTRANADDR, (void **)&isa PETSC_F90_2PTR_PARAM(ptrd1));
66ce78bad3SBarry Smith if (*ierr) return;
67ce78bad3SBarry Smith *ierr = F90Array1dAccess(is2, MPIU_FORTRANADDR, (void **)&isb PETSC_F90_2PTR_PARAM(ptrd2));
68ce78bad3SBarry Smith if (*ierr) return;
69ce78bad3SBarry Smith *ierr = PCASMDestroySubdomains(*n, &isa, &isb);
70ce78bad3SBarry Smith if (*ierr) return;
71ce78bad3SBarry Smith *ierr = F90Array1dDestroy(is1, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd1));
72ce78bad3SBarry Smith if (*ierr) return;
73ce78bad3SBarry Smith *ierr = F90Array1dDestroy(is2, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd2));
74ce78bad3SBarry Smith if (*ierr) return;
752b691e39SMatthew Knepley }
76ce78bad3SBarry Smith
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))77ce78bad3SBarry Smith 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))
78ce78bad3SBarry Smith {
79ce78bad3SBarry Smith IS *iis, *iisl;
80ce78bad3SBarry Smith
81ce78bad3SBarry Smith *ierr = PCASMCreateSubdomains2D(*m, *n, *M, *N, *dof, *overlap, Nsub, &iis, &iisl);
82ce78bad3SBarry Smith if (*ierr) return;
83ce78bad3SBarry Smith *ierr = F90Array1dCreate(iis, MPIU_FORTRANADDR, 1, *Nsub, is1 PETSC_F90_2PTR_PARAM(ptrd1));
84ce78bad3SBarry Smith if (*ierr) return;
85ce78bad3SBarry Smith *ierr = F90Array1dCreate(iisl, MPIU_FORTRANADDR, 1, *Nsub, is2 PETSC_F90_2PTR_PARAM(ptrd2));
86ce78bad3SBarry Smith if (*ierr) return;
87ce78bad3SBarry Smith }
88ce78bad3SBarry Smith
pcasmgetsubksp_(PC * pc,PetscInt * n_local,PetscInt * first_local,F90Array1d * ksp,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))89ce78bad3SBarry Smith PETSC_EXTERN void pcasmgetsubksp_(PC *pc, PetscInt *n_local, PetscInt *first_local, F90Array1d *ksp, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
90ce78bad3SBarry Smith {
91ce78bad3SBarry Smith KSP *tksp;
92ce78bad3SBarry Smith PetscInt nloc, flocal;
93ce78bad3SBarry Smith
94ce78bad3SBarry Smith CHKFORTRANNULLINTEGER(n_local);
95ce78bad3SBarry Smith CHKFORTRANNULLINTEGER(first_local);
96ce78bad3SBarry Smith *ierr = PCASMGetSubKSP(*pc, &nloc, &flocal, &tksp);
97ce78bad3SBarry Smith if (n_local) *n_local = nloc;
98ce78bad3SBarry Smith if (first_local) *first_local = flocal;
99ce78bad3SBarry Smith *ierr = F90Array1dCreate(tksp, MPIU_FORTRANADDR, 1, nloc, ksp PETSC_F90_2PTR_PARAM(ptrd));
100e54e4138SSatish Balay }
101*36083efbSBarry Smith
pcasmrestoresubksp_(PC * pc,PetscInt * n_local,PetscInt * first_local,F90Array1d * ksp,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))102*36083efbSBarry Smith PETSC_EXTERN void pcasmrestoresubksp_(PC *pc, PetscInt *n_local, PetscInt *first_local, F90Array1d *ksp, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
103*36083efbSBarry Smith {
104*36083efbSBarry Smith *ierr = F90Array1dDestroy(ksp, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd));
105*36083efbSBarry Smith }
106