xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <petsc/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 PETSC_EXTERN void PETSC_STDCALL pcasmcreatesubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
23 {
24   PetscInt i;
25   IS       *insubs;
26 
27   *ierr = PCASMCreateSubdomains(*mat,*n,&insubs);if (*ierr) return;
28   for (i=0; i<*n; i++) subs[i] = insubs[i];
29   *ierr = PetscFree(insubs);
30 }
31 
32 
33 PETSC_EXTERN void PETSC_STDCALL pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
34 {
35   PetscInt i;
36 
37   for (i=0; i<*n; i++) {
38     *ierr = ISDestroy(&subs[i]);if (*ierr) return;
39   }
40 }
41 
42 PETSC_EXTERN void PETSC_STDCALL pcasmgetsubksp_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
43 {
44   KSP      *tksp;
45   PetscInt i,nloc;
46   CHKFORTRANNULLINTEGER(n_local);
47   CHKFORTRANNULLINTEGER(first_local);
48   CHKFORTRANNULLOBJECT(ksp);
49   *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp);
50   if (n_local) *n_local = nloc;
51   if (ksp) {
52     for (i=0; i<nloc; i++) ksp[i] = tksp[i];
53   }
54 }
55 
56 PETSC_EXTERN void PETSC_STDCALL pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
57 {
58   CHKFORTRANNULLOBJECT(is);
59   CHKFORTRANNULLOBJECT(is_local);
60   *ierr = PCASMSetLocalSubdomains(*pc,*n,is,is_local);
61 }
62 
63 PETSC_EXTERN void PETSC_STDCALL pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is,IS *is_local, PetscErrorCode *ierr)
64 {
65   CHKFORTRANNULLOBJECT(is);
66   CHKFORTRANNULLOBJECT(is_local);
67   *ierr = PCASMSetTotalSubdomains(*pc,*N,is,is_local);
68 }
69 
70 PETSC_EXTERN 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++) mat[i] = tmat[i];
80   }
81 }
82 PETSC_EXTERN void PETSC_STDCALL pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
83 {
84   PetscInt nloc,i;
85   IS       *tis, *tis_local;
86   CHKFORTRANNULLOBJECT(is);
87   CHKFORTRANNULLOBJECT(is_local);
88   CHKFORTRANNULLINTEGER(n);
89   *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis,&tis_local);
90   if (n) *n = nloc;
91   if (is) {
92     for (i=0; i<nloc; i++) is[i] = tis[i];
93   }
94   if (is_local && tis_local) {
95     for (i=0; i<nloc; i++) is[i] = tis_local[i];
96   }
97 }
98 
99