xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision c36a59bee2e66f2e80e3e9e7603aedea4197893c)
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   #define pcasmcreatesubdomains2d_  PCASMCREATESUBDOMAINS2D
20 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
21   #define pcasmgetsubksp1_          pcasmgetsubksp1
22   #define pcasmgetsubksp2_          pcasmgetsubksp2
23   #define pcasmgetsubksp3_          pcasmgetsubksp3
24   #define pcasmgetsubksp4_          pcasmgetsubksp4
25   #define pcasmgetsubksp5_          pcasmgetsubksp5
26   #define pcasmgetsubksp6_          pcasmgetsubksp6
27   #define pcasmgetsubksp7_          pcasmgetsubksp7
28   #define pcasmgetsubksp8_          pcasmgetsubksp8
29   #define pcasmsetlocalsubdomains_  pcasmsetlocalsubdomains
30   #define pcasmsetglobalsubdomains_ pcasmsetglobalsubdomains
31   #define pcasmgetlocalsubmatrices_ pcasmgetlocalsubmatrices
32   #define pcasmgetlocalsubdomains_  pcasmgetlocalsubdomains
33   #define pcasmcreatesubdomains_    pcasmcreatesubdomains
34   #define pcasmdestroysubdomains_   pcasmdestroysubdomains
35   #define pcasmcreatesubdomains2d_  pcasmcreatesubdomains2d
36 #endif
37 
38 PETSC_EXTERN void pcasmcreatesubdomains2d_(PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N, PetscInt *dof, PetscInt *overlap, PetscInt *Nsub, IS *is, IS *isl, int *ierr)
39 {
40   IS *iis, *iisl;
41   *ierr = PCASMCreateSubdomains2D(*m, *n, *M, *N, *dof, *overlap, Nsub, &iis, &iisl);
42   if (*ierr) return;
43   *ierr = PetscMemcpy(is, iis, *Nsub * sizeof(IS));
44   if (*ierr) return;
45   *ierr = PetscMemcpy(isl, iisl, *Nsub * sizeof(IS));
46   if (*ierr) return;
47   *ierr = PetscFree(iis);
48   if (*ierr) return;
49   *ierr = PetscFree(iisl);
50 }
51 
52 PETSC_EXTERN void pcasmcreatesubdomains_(Mat *mat, PetscInt *n, IS *subs, PetscErrorCode *ierr)
53 {
54   PetscInt i;
55   IS      *insubs;
56 
57   *ierr = PCASMCreateSubdomains(*mat, *n, &insubs);
58   if (*ierr) return;
59   for (i = 0; i < *n; i++) subs[i] = insubs[i];
60   *ierr = PetscFree(insubs);
61 }
62 
63 PETSC_EXTERN void pcasmdestroysubdomains_(Mat *mat, PetscInt *n, IS *subs, PetscErrorCode *ierr)
64 {
65   PetscInt i;
66 
67   for (i = 0; i < *n; i++) {
68     *ierr = ISDestroy(&subs[i]);
69     if (*ierr) return;
70   }
71 }
72 
73 PETSC_EXTERN void pcasmgetsubksp1_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
74 {
75   KSP     *tksp;
76   PetscInt i, nloc;
77   CHKFORTRANNULLINTEGER(n_local);
78   CHKFORTRANNULLINTEGER(first_local);
79   CHKFORTRANNULLOBJECT(ksp);
80   *ierr = PCASMGetSubKSP(*pc, &nloc, first_local, &tksp);
81   if (n_local) *n_local = nloc;
82   if (ksp) {
83     for (i = 0; i < nloc; i++) ksp[i] = tksp[i];
84   }
85 }
86 
87 PETSC_EXTERN void pcasmgetsubksp2_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
88 {
89   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
90 }
91 
92 PETSC_EXTERN void pcasmgetsubksp3_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
93 {
94   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
95 }
96 
97 PETSC_EXTERN void pcasmgetsubksp4_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
98 {
99   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
100 }
101 
102 PETSC_EXTERN void pcasmgetsubksp5_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
103 {
104   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
105 }
106 
107 PETSC_EXTERN void pcasmgetsubksp6_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
108 {
109   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
110 }
111 
112 PETSC_EXTERN void pcasmgetsubksp7_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
113 {
114   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
115 }
116 
117 PETSC_EXTERN void pcasmgetsubksp8_(PC *pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp, PetscErrorCode *ierr)
118 {
119   pcasmgetsubksp1_(pc, n_local, first_local, ksp, ierr);
120 }
121 
122 PETSC_EXTERN void pcasmsetlocalsubdomains_(PC *pc, PetscInt *n, IS *is, IS *is_local, PetscErrorCode *ierr)
123 {
124   CHKFORTRANNULLOBJECT(is);
125   CHKFORTRANNULLOBJECT(is_local);
126   *ierr = PCASMSetLocalSubdomains(*pc, *n, is, is_local);
127 }
128 
129 PETSC_EXTERN void pcasmsettotalsubdomains_(PC *pc, PetscInt *N, IS *is, IS *is_local, PetscErrorCode *ierr)
130 {
131   CHKFORTRANNULLOBJECT(is);
132   CHKFORTRANNULLOBJECT(is_local);
133   *ierr = PCASMSetTotalSubdomains(*pc, *N, is, is_local);
134 }
135 
136 PETSC_EXTERN void pcasmgetlocalsubmatrices_(PC *pc, PetscInt *n, Mat *mat, PetscErrorCode *ierr)
137 {
138   PetscInt nloc, i;
139   Mat     *tmat;
140   CHKFORTRANNULLOBJECT(mat);
141   CHKFORTRANNULLINTEGER(n);
142   *ierr = PCASMGetLocalSubmatrices(*pc, &nloc, &tmat);
143   if (n) *n = nloc;
144   if (mat) {
145     for (i = 0; i < nloc; i++) mat[i] = tmat[i];
146   }
147 }
148 PETSC_EXTERN void pcasmgetlocalsubdomains_(PC *pc, PetscInt *n, IS *is, IS *is_local, PetscErrorCode *ierr)
149 {
150   PetscInt nloc, i;
151   IS      *tis, *tis_local;
152   CHKFORTRANNULLOBJECT(is);
153   CHKFORTRANNULLOBJECT(is_local);
154   CHKFORTRANNULLINTEGER(n);
155   *ierr = PCASMGetLocalSubdomains(*pc, &nloc, &tis, &tis_local);
156   if (n) *n = nloc;
157   if (is) {
158     for (i = 0; i < nloc; i++) is[i] = tis[i];
159   }
160   if (is_local && tis_local) {
161     for (i = 0; i < nloc; i++) is_local[i] = tis_local[i];
162   }
163 }
164