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