xref: /petsc/src/ksp/pc/impls/gasm/ftn-custom/zgasmf.c (revision fa084801f6b15df01ac44a0e53249c011483a183)
1 #include <petsc/private/f90impl.h>
2 #include <petscksp.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define pcgasmdestroysubdomains_  PCGASMDESTROYSUBDOMAINS
6   #define pcgasmgetsubksp           PCGASMGETSUBKSP
7   #define pcgasmcreatesubdomains2d_ PCGASMCREATESUBDOMAINS2D
8 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9   #define pcgasmdestroysubdomains_  pcgasmdestroysubdomains
10   #define pcgasmcreatesubdomains2d_ pcgasmcreatesubdomains2d
11 #endif
12 
13 PETSC_EXTERN void pcgasmdestroysubdomains_(PetscInt *n, F90Array1d *is1, F90Array1d *is2, int *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
14 {
15   IS *isa, *isb;
16 
17   *ierr = F90Array1dAccess(is1, MPIU_FORTRANADDR, (void **)&isa PETSC_F90_2PTR_PARAM(ptrd1));
18   if (*ierr) return;
19   *ierr = F90Array1dAccess(is2, MPIU_FORTRANADDR, (void **)&isb PETSC_F90_2PTR_PARAM(ptrd2));
20   if (*ierr) return;
21   *ierr = F90Array1dDestroy(is1, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd1));
22   if (*ierr) return;
23   *ierr = F90Array1dDestroy(is2, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd2));
24   if (*ierr) return;
25   *ierr = PCGASMDestroySubdomains(*n, &isa, &isb);
26 }
27 
28 PETSC_EXTERN void pcgasmcreatesubdomains2d_(PC *pc, 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))
29 {
30   IS *iis, *iisl;
31   *ierr = PCGASMCreateSubdomains2D(*pc, *m, *n, *M, *N, *dof, *overlap, Nsub, &iis, &iisl);
32   if (*ierr) return;
33   *ierr = F90Array1dCreate(iis, MPIU_FORTRANADDR, 1, *Nsub, is1 PETSC_F90_2PTR_PARAM(ptrd1));
34   if (*ierr) return;
35   *ierr = F90Array1dCreate(iisl, MPIU_FORTRANADDR, 1, *Nsub, is2 PETSC_F90_2PTR_PARAM(ptrd2));
36   if (*ierr) return;
37 }
38 
39 PETSC_EXTERN void pcgasmgetsubksp_(PC *pc, PetscInt *n_local, PetscInt *first_local, F90Array1d *ksp, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
40 {
41   KSP     *tksp;
42   PetscInt nloc, flocal;
43 
44   CHKFORTRANNULLINTEGER(n_local);
45   CHKFORTRANNULLINTEGER(first_local);
46   *ierr = PCGASMGetSubKSP(*pc, &nloc, &flocal, &tksp);
47   if (n_local) *n_local = nloc;
48   if (first_local) *first_local = flocal;
49   *ierr = F90Array1dCreate(tksp, MPIU_FORTRANADDR, 1, nloc, ksp PETSC_F90_2PTR_PARAM(ptrd));
50 }
51