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