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