xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision efa12513287cff49a2b9648ae83199dcbfaad71a)
1 #include <petsc/private/fortranimpl.h>
2 #include <petscksp.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define pcasmgetsubksp1_           PCASMGETSUBKSP1
6 #define pcasmgetsubksp2_           PCASMGETSUBKSP2
7 #define pcasmgetsubksp3_           PCASMGETSUBKSP3
8 #define pcasmgetsubksp4_           PCASMGETSUBKSP4
9 #define pcasmgetsubksp5_           PCASMGETSUBKSP5
10 #define pcasmgetsubksp6_           PCASMGETSUBKSP6
11 #define pcasmgetsubksp7_           PCASMGETSUBKSP7
12 #define pcasmgetsubksp8_           PCASMGETSUBKSP8
13 #define pcasmsetlocalsubdomains_   PCASMSETLOCALSUBDOMAINS
14 #define pcasmsetglobalsubdomains_  PCASMSETGLOBALSUBDOMAINS
15 #define pcasmgetlocalsubmatrices_  PCASMGETLOCALSUBMATRICES
16 #define pcasmgetlocalsubdomains_   PCASMGETLOCALSUBDOMAINS
17 #define pcasmcreatesubdomains_     PCASMCREATESUBDOMAINS
18 #define pcasmdestroysubdomains_    PCASMDESTROYSUBDOMAINS
19 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
20 #define pcasmgetsubksp1_           pcasmgetsubksp1
21 #define pcasmgetsubksp2_           pcasmgetsubksp2
22 #define pcasmgetsubksp3_           pcasmgetsubksp3
23 #define pcasmgetsubksp4_           pcasmgetsubksp4
24 #define pcasmgetsubksp5_           pcasmgetsubksp5
25 #define pcasmgetsubksp6_           pcasmgetsubksp6
26 #define pcasmgetsubksp7_           pcasmgetsubksp7
27 #define pcasmgetsubksp8_           pcasmgetsubksp8
28 #define pcasmsetlocalsubdomains_   pcasmsetlocalsubdomains
29 #define pcasmsetglobalsubdomains_  pcasmsetglobalsubdomains
30 #define pcasmgetlocalsubmatrices_  pcasmgetlocalsubmatrices
31 #define pcasmgetlocalsubdomains_   pcasmgetlocalsubdomains
32 #define pcasmcreatesubdomains_     pcasmcreatesubdomains
33 #define pcasmdestroysubdomains_    pcasmdestroysubdomains
34 #endif
35 
36 PETSC_EXTERN void pcasmcreatesubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
37 {
38   PetscInt i;
39   IS       *insubs;
40 
41   *ierr = PCASMCreateSubdomains(*mat,*n,&insubs);if (*ierr) return;
42   for (i=0; i<*n; i++) subs[i] = insubs[i];
43   *ierr = PetscFree(insubs);
44 }
45 
46 
47 PETSC_EXTERN void pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
48 {
49   PetscInt i;
50 
51   for (i=0; i<*n; i++) {
52     *ierr = ISDestroy(&subs[i]);if (*ierr) return;
53   }
54 }
55 
56 PETSC_EXTERN void pcasmgetsubksp1_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
57 {
58   KSP      *tksp;
59   PetscInt i,nloc;
60   CHKFORTRANNULLINTEGER(n_local);
61   CHKFORTRANNULLINTEGER(first_local);
62   CHKFORTRANNULLOBJECT(ksp);
63   *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp);
64   if (n_local) *n_local = nloc;
65   if (ksp) {
66     for (i=0; i<nloc; i++) ksp[i] = tksp[i];
67   }
68 }
69 
70 PETSC_EXTERN void pcasmgetsubksp2_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
71 {
72   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
73 }
74 
75 PETSC_EXTERN void pcasmgetsubksp3_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
76 {
77   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
78 }
79 
80 PETSC_EXTERN void pcasmgetsubksp4_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
81 {
82   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
83 }
84 
85 PETSC_EXTERN void pcasmgetsubksp5_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
86 {
87   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
88 }
89 
90 PETSC_EXTERN void pcasmgetsubksp6_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
91 {
92   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
93 }
94 
95 PETSC_EXTERN void pcasmgetsubksp7_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
96 {
97   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
98 }
99 
100 PETSC_EXTERN void pcasmgetsubksp8_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
101 {
102   pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr);
103 }
104 
105 PETSC_EXTERN void pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
106 {
107   CHKFORTRANNULLOBJECT(is);
108   CHKFORTRANNULLOBJECT(is_local);
109   *ierr = PCASMSetLocalSubdomains(*pc,*n,is,is_local);
110 }
111 
112 PETSC_EXTERN void pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is,IS *is_local, PetscErrorCode *ierr)
113 {
114   CHKFORTRANNULLOBJECT(is);
115   CHKFORTRANNULLOBJECT(is_local);
116   *ierr = PCASMSetTotalSubdomains(*pc,*N,is,is_local);
117 }
118 
119 PETSC_EXTERN void pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr)
120 {
121   PetscInt nloc,i;
122   Mat      *tmat;
123   CHKFORTRANNULLOBJECT(mat);
124   CHKFORTRANNULLINTEGER(n);
125   *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat);
126   if (n) *n = nloc;
127   if (mat) {
128     for (i=0; i<nloc; i++) mat[i] = tmat[i];
129   }
130 }
131 PETSC_EXTERN void pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
132 {
133   PetscInt nloc,i;
134   IS       *tis, *tis_local;
135   CHKFORTRANNULLOBJECT(is);
136   CHKFORTRANNULLOBJECT(is_local);
137   CHKFORTRANNULLINTEGER(n);
138   *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis,&tis_local);
139   if (n) *n = nloc;
140   if (is) {
141     for (i=0; i<nloc; i++) is[i] = tis[i];
142   }
143   if (is_local && tis_local) {
144     for (i=0; i<nloc; i++) is_local[i] = tis_local[i];
145   }
146 }
147 
148