xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision 7b6bb2c608b6fc6714ef38fda02c2dbb91c82665)
1 #include <private/fortranimpl.h>
2 #include <petscksp.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define pcasmgetsubksp_            PCASMGETSUBKSP
6 #define pcasmsetlocalsubdomains_   PCASMSETLOCALSUBDOMAINS
7 #define pcasmsetglobalsubdomains_  PCASMSETGLOBALSUBDOMAINS
8 #define pcasmgetlocalsubmatrices_  PCASMGETLOCALSUBMATRICES
9 #define pcasmgetlocalsubdomains_   PCASMGETLOCALSUBDOMAINS
10 #define pcasmcreatesubdomains_     PCASMCREATESUBDOMAINS
11 #define pcasmdestroysubdomains_    PCASMDESTROYSUBDOMAINS
12 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
13 #define pcasmgetsubksp_            pcasmgetsubksp
14 #define pcasmsetlocalsubdomains_   pcasmsetlocalsubdomains
15 #define pcasmsetglobalsubdomains_  pcasmsetglobalsubdomains
16 #define pcasmgetlocalsubmatrices_  pcasmgetlocalsubmatrices
17 #define pcasmgetlocalsubdomains_   pcasmgetlocalsubdomains
18 #define pcasmcreatesubdomains_     pcasmcreatesubdomains
19 #define pcasmdestroysubdomains_    pcasmdestroysubdomains
20 #endif
21 
22 EXTERN_C_BEGIN
23 void PETSC_STDCALL pcasmcreatesubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
24 {
25   PetscInt i;
26   IS       *insubs;
27 
28   *ierr = PCASMCreateSubdomains(*mat,*n,&insubs);if (*ierr) return;
29   for (i=0; i<*n; i++) {
30     subs[i] = insubs[i];
31   }
32   *ierr = PetscFree(insubs);
33 }
34 
35 
36 void PETSC_STDCALL pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
37 {
38   PetscInt i;
39 
40   for (i=0; i<*n; i++) {
41     *ierr = ISDestroy(&subs[i]);if (*ierr) return;
42   }
43 }
44 
45 void PETSC_STDCALL pcasmgetsubksp_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
46 {
47   KSP *tksp;
48   PetscInt  i,nloc;
49   CHKFORTRANNULLINTEGER(n_local);
50   CHKFORTRANNULLINTEGER(first_local);
51   *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp);
52   if (n_local) *n_local = nloc;
53   for (i=0; i<nloc; i++){
54     ksp[i] = tksp[i];
55   }
56 }
57 
58 void PETSC_STDCALL pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
59 {
60   CHKFORTRANNULLOBJECT(is);
61   CHKFORTRANNULLOBJECT(is_local);
62   *ierr = PCASMSetLocalSubdomains(*pc,*n,is,is_local);
63 }
64 
65 void PETSC_STDCALL pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is,IS *is_local, PetscErrorCode *ierr)
66 {
67   CHKFORTRANNULLOBJECT(is);
68   CHKFORTRANNULLOBJECT(is_local);
69   *ierr = PCASMSetTotalSubdomains(*pc,*N,is,is_local);
70 }
71 
72 void PETSC_STDCALL pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr)
73 {
74   PetscInt nloc,i;
75   Mat  *tmat;
76   CHKFORTRANNULLOBJECT(mat);
77   CHKFORTRANNULLINTEGER(n);
78   *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat);
79   if (n) *n = nloc;
80   if (mat) {
81     for (i=0; i<nloc; i++){
82       mat[i] = tmat[i];
83     }
84   }
85 }
86 void PETSC_STDCALL pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
87 {
88   PetscInt nloc,i;
89   IS  *tis, *tis_local;
90   CHKFORTRANNULLOBJECT(is);
91   CHKFORTRANNULLOBJECT(is_local);
92   CHKFORTRANNULLINTEGER(n);
93   *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis,&tis_local);
94   if (n) *n = nloc;
95   if (is) {
96     for (i=0; i<nloc; i++){
97       is[i] = tis[i];
98     }
99   }
100   if (is_local && tis_local) {
101     for (i=0; i<nloc; i++){
102       is[i] = tis_local[i];
103     }
104   }
105 }
106 
107 EXTERN_C_END
108