xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision 5d00a290c924ab1b2c83000ba3a0ae1df9cca2e4)
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, PetscErrorCode *ierr)
59 {
60   CHKFORTRANNULLOBJECT(is);
61   *ierr = PCASMSetLocalSubdomains(*pc,*n,is);
62 }
63 
64 void PETSC_STDCALL pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is, PetscErrorCode *ierr)
65 {
66   CHKFORTRANNULLOBJECT(is);
67   *ierr = PCASMSetTotalSubdomains(*pc,*N,is);
68 }
69 
70 void PETSC_STDCALL pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr)
71 {
72   PetscInt nloc,i;
73   Mat  *tmat;
74   CHKFORTRANNULLOBJECT(mat);
75   CHKFORTRANNULLINTEGER(n);
76   *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat);
77   if (n) *n = nloc;
78   if (mat) {
79     for (i=0; i<nloc; i++){
80       mat[i] = tmat[i];
81     }
82   }
83 }
84 void PETSC_STDCALL pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is, PetscErrorCode *ierr)
85 {
86   PetscInt nloc,i;
87   IS  *tis;
88   CHKFORTRANNULLOBJECT(is);
89   CHKFORTRANNULLINTEGER(n);
90   *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis);
91   if (n) *n = nloc;
92   if (is) {
93     for (i=0; i<nloc; i++){
94       is[i] = tis[i];
95     }
96   }
97 }
98 
99 EXTERN_C_END
100