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 PETSC_EXTERN void pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr) 47 { 48 PetscInt i; 49 50 for (i=0; i<*n; i++) { 51 *ierr = ISDestroy(&subs[i]);if (*ierr) return; 52 } 53 } 54 55 PETSC_EXTERN void pcasmgetsubksp1_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 56 { 57 KSP *tksp; 58 PetscInt i,nloc; 59 CHKFORTRANNULLINTEGER(n_local); 60 CHKFORTRANNULLINTEGER(first_local); 61 CHKFORTRANNULLOBJECT(ksp); 62 *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp); 63 if (n_local) *n_local = nloc; 64 if (ksp) { 65 for (i=0; i<nloc; i++) ksp[i] = tksp[i]; 66 } 67 } 68 69 PETSC_EXTERN void pcasmgetsubksp2_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 70 { 71 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 72 } 73 74 PETSC_EXTERN void pcasmgetsubksp3_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 75 { 76 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 77 } 78 79 PETSC_EXTERN void pcasmgetsubksp4_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 80 { 81 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 82 } 83 84 PETSC_EXTERN void pcasmgetsubksp5_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 85 { 86 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 87 } 88 89 PETSC_EXTERN void pcasmgetsubksp6_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 90 { 91 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 92 } 93 94 PETSC_EXTERN void pcasmgetsubksp7_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 95 { 96 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 97 } 98 99 PETSC_EXTERN void pcasmgetsubksp8_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 100 { 101 pcasmgetsubksp1_(pc,n_local,first_local,ksp,ierr); 102 } 103 104 PETSC_EXTERN void pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr) 105 { 106 CHKFORTRANNULLOBJECT(is); 107 CHKFORTRANNULLOBJECT(is_local); 108 *ierr = PCASMSetLocalSubdomains(*pc,*n,is,is_local); 109 } 110 111 PETSC_EXTERN void pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is,IS *is_local, PetscErrorCode *ierr) 112 { 113 CHKFORTRANNULLOBJECT(is); 114 CHKFORTRANNULLOBJECT(is_local); 115 *ierr = PCASMSetTotalSubdomains(*pc,*N,is,is_local); 116 } 117 118 PETSC_EXTERN void pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr) 119 { 120 PetscInt nloc,i; 121 Mat *tmat; 122 CHKFORTRANNULLOBJECT(mat); 123 CHKFORTRANNULLINTEGER(n); 124 *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat); 125 if (n) *n = nloc; 126 if (mat) { 127 for (i=0; i<nloc; i++) mat[i] = tmat[i]; 128 } 129 } 130 PETSC_EXTERN void pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr) 131 { 132 PetscInt nloc,i; 133 IS *tis, *tis_local; 134 CHKFORTRANNULLOBJECT(is); 135 CHKFORTRANNULLOBJECT(is_local); 136 CHKFORTRANNULLINTEGER(n); 137 *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis,&tis_local); 138 if (n) *n = nloc; 139 if (is) { 140 for (i=0; i<nloc; i++) is[i] = tis[i]; 141 } 142 if (is_local && tis_local) { 143 for (i=0; i<nloc; i++) is_local[i] = tis_local[i]; 144 } 145 } 146 147