1ce0a2cd1SBarry Smith #include "private/fortranimpl.h" 2e54e4138SSatish Balay #include "petscksp.h" 3e54e4138SSatish Balay 4e54e4138SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS) 5e54e4138SSatish Balay #define pcasmgetsubksp_ PCASMGETSUBKSP 6e54e4138SSatish Balay #define pcasmsetlocalsubdomains_ PCASMSETLOCALSUBDOMAINS 7e54e4138SSatish Balay #define pcasmsetglobalsubdomains_ PCASMSETGLOBALSUBDOMAINS 8e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_ PCASMGETLOCALSUBMATRICES 9e54e4138SSatish Balay #define pcasmgetlocalsubdomains_ PCASMGETLOCALSUBDOMAINS 10*ef356e52SBarry Smith #define pcasmcreatesubdomains_ PCASMCREATESUBDOMAINS 11*ef356e52SBarry Smith #define pcasmdestroysubdomains_ PCASMDESTROYSUBDOMAINS 12e54e4138SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 13e54e4138SSatish Balay #define pcasmgetsubksp_ pcasmgetsubksp 14e54e4138SSatish Balay #define pcasmsetlocalsubdomains_ pcasmsetlocalsubdomains 15e54e4138SSatish Balay #define pcasmsetglobalsubdomains_ pcasmsetglobalsubdomains 16e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_ pcasmgetlocalsubmatrices 17e54e4138SSatish Balay #define pcasmgetlocalsubdomains_ pcasmgetlocalsubdomains 18*ef356e52SBarry Smith #define pcasmcreatesubdomains_ pcasmcreatesubdomains 19*ef356e52SBarry Smith #define pcasmdestroysubdomains_ pcasmdestroysubdomains 20e54e4138SSatish Balay #endif 21e54e4138SSatish Balay 22e54e4138SSatish Balay EXTERN_C_BEGIN 23*ef356e52SBarry Smith void PETSC_STDCALL pcasmcreatesubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr) 24*ef356e52SBarry Smith { 25*ef356e52SBarry Smith PetscInt i; 26*ef356e52SBarry Smith IS *insubs; 27*ef356e52SBarry Smith 28*ef356e52SBarry Smith *ierr = PCASMCreateSubdomains(*mat,*n,&insubs);if (*ierr) return; 29*ef356e52SBarry Smith for (i=0; i<*n; i++) { 30*ef356e52SBarry Smith subs[i] = insubs[i]; 31*ef356e52SBarry Smith } 32*ef356e52SBarry Smith *ierr = PetscFree(insubs); 33*ef356e52SBarry Smith } 34*ef356e52SBarry Smith 35*ef356e52SBarry Smith 36*ef356e52SBarry Smith void PETSC_STDCALL pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr) 37*ef356e52SBarry Smith { 38*ef356e52SBarry Smith PetscInt i; 39*ef356e52SBarry Smith 40*ef356e52SBarry Smith for (i=0; i<*n; i++) { 41*ef356e52SBarry Smith *ierr = ISDestroy(subs[i]);if (*ierr) return; 42*ef356e52SBarry Smith } 43*ef356e52SBarry Smith } 44*ef356e52SBarry Smith 45e54e4138SSatish Balay void PETSC_STDCALL pcasmgetsubksp_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr) 46e54e4138SSatish Balay { 47e54e4138SSatish Balay KSP *tksp; 48e54e4138SSatish Balay PetscInt i,nloc; 49e54e4138SSatish Balay CHKFORTRANNULLINTEGER(n_local); 50e54e4138SSatish Balay CHKFORTRANNULLINTEGER(first_local); 51e54e4138SSatish Balay *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp); 52e54e4138SSatish Balay if (n_local) *n_local = nloc; 53e54e4138SSatish Balay for (i=0; i<nloc; i++){ 54e54e4138SSatish Balay ksp[i] = tksp[i]; 55e54e4138SSatish Balay } 56e54e4138SSatish Balay } 57e54e4138SSatish Balay 58e54e4138SSatish Balay void PETSC_STDCALL pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is, PetscErrorCode *ierr) 59e54e4138SSatish Balay { 60e54e4138SSatish Balay CHKFORTRANNULLOBJECT(is); 61e54e4138SSatish Balay *ierr = PCASMSetLocalSubdomains(*pc,*n,is); 62e54e4138SSatish Balay } 63e54e4138SSatish Balay 64e54e4138SSatish Balay void PETSC_STDCALL pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is, PetscErrorCode *ierr) 65e54e4138SSatish Balay { 66e54e4138SSatish Balay CHKFORTRANNULLOBJECT(is); 67e54e4138SSatish Balay *ierr = PCASMSetTotalSubdomains(*pc,*N,is); 68e54e4138SSatish Balay } 69e54e4138SSatish Balay 70e54e4138SSatish Balay void PETSC_STDCALL pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr) 71e54e4138SSatish Balay { 72e54e4138SSatish Balay PetscInt nloc,i; 73e54e4138SSatish Balay Mat *tmat; 74e54e4138SSatish Balay CHKFORTRANNULLOBJECT(mat); 75e54e4138SSatish Balay CHKFORTRANNULLINTEGER(n); 76e54e4138SSatish Balay *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat); 77e54e4138SSatish Balay if (n) *n = nloc; 78e54e4138SSatish Balay if (mat) { 79e54e4138SSatish Balay for (i=0; i<nloc; i++){ 80e54e4138SSatish Balay mat[i] = tmat[i]; 81e54e4138SSatish Balay } 82e54e4138SSatish Balay } 83e54e4138SSatish Balay } 84e54e4138SSatish Balay void PETSC_STDCALL pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is, PetscErrorCode *ierr) 85e54e4138SSatish Balay { 86e54e4138SSatish Balay PetscInt nloc,i; 87e54e4138SSatish Balay IS *tis; 88e54e4138SSatish Balay CHKFORTRANNULLOBJECT(is); 89e54e4138SSatish Balay CHKFORTRANNULLINTEGER(n); 90e54e4138SSatish Balay *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis); 91e54e4138SSatish Balay if (n) *n = nloc; 92e54e4138SSatish Balay if (is) { 93e54e4138SSatish Balay for (i=0; i<nloc; i++){ 94e54e4138SSatish Balay is[i] = tis[i]; 95e54e4138SSatish Balay } 96e54e4138SSatish Balay } 97e54e4138SSatish Balay } 98e54e4138SSatish Balay 99e54e4138SSatish Balay EXTERN_C_END 100