1c6db04a5SJed Brown #include <petscdmcomposite.h> 2*ce78bad3SBarry Smith #include <petsc/private/f90impl.h> 347c6ae99SBarry Smith 447c6ae99SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 547c6ae99SBarry Smith #define dmcompositegetentries1_ DMCOMPOSITEGETENTRIES1 647c6ae99SBarry Smith #define dmcompositegetentries2_ DMCOMPOSITEGETENTRIES2 747c6ae99SBarry Smith #define dmcompositegetentries3_ DMCOMPOSITEGETENTRIES3 847c6ae99SBarry Smith #define dmcompositegetentries4_ DMCOMPOSITEGETENTRIES4 947c6ae99SBarry Smith #define dmcompositegetentries5_ DMCOMPOSITEGETENTRIES5 1047c6ae99SBarry Smith #define dmcompositegetaccess4_ DMCOMPOSITEGETACCESS4 1147c6ae99SBarry Smith #define dmcompositescatter4_ DMCOMPOSITESCATTER4 1247c6ae99SBarry Smith #define dmcompositerestoreaccess4_ DMCOMPOSITERESTOREACCESS4 1347c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ DMCOMPOSITEGETLOCALVECTORS4 1447c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_ DMCOMPOSITERESTORELOCALVECTORS4 1573e31fe2SJed Brown #define dmcompositegetglobaliss_ DMCOMPOSITEGETGLOBALISS 16*ce78bad3SBarry Smith #define dmcompositerestoreglobaliss_ DMCOMPOSITERESTOREGLOBALISS 1773e31fe2SJed Brown #define dmcompositegetlocaliss_ DMCOMPOSITEGETLOCALISS 18*ce78bad3SBarry Smith #define dmcompositerestorelocaliss_ DMCOMPOSITERESTORELOCALISS 1947c6ae99SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 2047c6ae99SBarry Smith #define dmcompositegetentries1_ dmcompositegetentries1 2147c6ae99SBarry Smith #define dmcompositegetentries2_ dmcompositegetentries2 2247c6ae99SBarry Smith #define dmcompositegetentries3_ dmcompositegetentries3 2347c6ae99SBarry Smith #define dmcompositegetentries4_ dmcompositegetentries4 2447c6ae99SBarry Smith #define dmcompositegetentries5_ dmcompositegetentries5 2547c6ae99SBarry Smith #define dmcompositegetaccess4_ dmcompositegetaccess4 2647c6ae99SBarry Smith #define dmcompositescatter4_ dmcompositescatter4 2747c6ae99SBarry Smith #define dmcompositerestoreaccess4_ dmcompositerestoreaccess4 2847c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ dmcompositegetlocalvectors4 2947c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4 3073e31fe2SJed Brown #define dmcompositegetglobaliss_ dmcompositegetglobaliss 31*ce78bad3SBarry Smith #define dmcompositerestoreglobaliss_ dmcompositerestoreglobaliss 3273e31fe2SJed Brown #define dmcompositegetlocaliss_ dmcompositegetlocaliss 33*ce78bad3SBarry Smith #define dmcompositerestorelocaliss_ dmcompositerestorelocaliss 3447c6ae99SBarry Smith #endif 3547c6ae99SBarry Smith 3619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries1_(DM *dm, DM *da1, PetscErrorCode *ierr) 3747c6ae99SBarry Smith { 3847c6ae99SBarry Smith *ierr = DMCompositeGetEntries(*dm, da1); 3947c6ae99SBarry Smith } 4047c6ae99SBarry Smith 4119caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries2_(DM *dm, DM *da1, DM *da2, PetscErrorCode *ierr) 4247c6ae99SBarry Smith { 4347c6ae99SBarry Smith *ierr = DMCompositeGetEntries(*dm, da1, da2); 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith 4619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries3_(DM *dm, DM *da1, DM *da2, DM *da3, PetscErrorCode *ierr) 4747c6ae99SBarry Smith { 4847c6ae99SBarry Smith *ierr = DMCompositeGetEntries(*dm, da1, da2, da3); 4947c6ae99SBarry Smith } 5047c6ae99SBarry Smith 5119caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries4_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, PetscErrorCode *ierr) 5247c6ae99SBarry Smith { 5347c6ae99SBarry Smith *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4); 5447c6ae99SBarry Smith } 5547c6ae99SBarry Smith 5619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries5_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, DM *da5, PetscErrorCode *ierr) 5747c6ae99SBarry Smith { 5847c6ae99SBarry Smith *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4, da5); 5947c6ae99SBarry Smith } 6047c6ae99SBarry Smith 6119caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetaccess4_(DM *dm, Vec *v, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 6247c6ae99SBarry Smith { 6347c6ae99SBarry Smith Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 6447c6ae99SBarry Smith *ierr = DMCompositeGetAccess(*dm, *v, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith 6719caf8f3SSatish Balay PETSC_EXTERN void dmcompositescatter4_(DM *dm, Vec *v, void *v1, void *p1, void *v2, void *p2, PetscErrorCode *ierr) 6847c6ae99SBarry Smith { 6947c6ae99SBarry Smith Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 7047c6ae99SBarry Smith *ierr = DMCompositeScatter(*dm, *v, *vv1, (PetscScalar *)p1, *vv2, (PetscScalar *)p2); 7147c6ae99SBarry Smith } 7247c6ae99SBarry Smith 7319caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestoreaccess4_(DM *dm, Vec *v, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 7447c6ae99SBarry Smith { 7547c6ae99SBarry Smith *ierr = DMCompositeRestoreAccess(*dm, *v, (Vec *)v1, 0, (Vec *)v2, 0); 7647c6ae99SBarry Smith } 7747c6ae99SBarry Smith 7819caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetlocalvectors4_(DM *dm, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 7947c6ae99SBarry Smith { 8047c6ae99SBarry Smith Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 8147c6ae99SBarry Smith *ierr = DMCompositeGetLocalVectors(*dm, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith 8419caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestorelocalvectors4_(DM *dm, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 8547c6ae99SBarry Smith { 8647c6ae99SBarry Smith Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 8747c6ae99SBarry Smith *ierr = DMCompositeRestoreLocalVectors(*dm, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 8847c6ae99SBarry Smith } 8947c6ae99SBarry Smith 90*ce78bad3SBarry Smith PETSC_EXTERN void dmcompositegetglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd)) 9173e31fe2SJed Brown { 9273e31fe2SJed Brown IS *ais; 93*ce78bad3SBarry Smith PetscInt ndm; 94*ce78bad3SBarry Smith 955975b3b6SBarry Smith *ierr = DMCompositeGetGlobalISs(*dm, &ais); 965975b3b6SBarry Smith if (*ierr) return; 975975b3b6SBarry Smith *ierr = DMCompositeGetNumberDM(*dm, &ndm); 985975b3b6SBarry Smith if (*ierr) return; 99*ce78bad3SBarry Smith *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd)); 10073e31fe2SJed Brown } 10173e31fe2SJed Brown 102*ce78bad3SBarry Smith PETSC_EXTERN void dmcompositerestoreglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd)) 10373e31fe2SJed Brown { 10473e31fe2SJed Brown IS *ais; 105*ce78bad3SBarry Smith PetscInt ndm; 106*ce78bad3SBarry Smith 107*ce78bad3SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd)); 108*ce78bad3SBarry Smith if (*ierr) return; 109*ce78bad3SBarry Smith *ierr = DMCompositeGetNumberDM(*dm, &ndm); 110*ce78bad3SBarry Smith for (PetscInt i = 0; i < ndm; i++) { 111*ce78bad3SBarry Smith *ierr = ISDestroy(&ais[i]); 112*ce78bad3SBarry Smith if (*ierr) return; 113*ce78bad3SBarry Smith } 114*ce78bad3SBarry Smith *ierr = PetscFree(ais); 115*ce78bad3SBarry Smith if (*ierr) return; 116*ce78bad3SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd)); 117*ce78bad3SBarry Smith } 118*ce78bad3SBarry Smith 119*ce78bad3SBarry Smith PETSC_EXTERN void dmcompositegetlocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd)) 120*ce78bad3SBarry Smith { 121*ce78bad3SBarry Smith IS *ais; 122*ce78bad3SBarry Smith PetscInt ndm; 123*ce78bad3SBarry Smith 1245975b3b6SBarry Smith *ierr = DMCompositeGetLocalISs(*dm, &ais); 1255975b3b6SBarry Smith if (*ierr) return; 1265975b3b6SBarry Smith *ierr = DMCompositeGetNumberDM(*dm, &ndm); 1275975b3b6SBarry Smith if (*ierr) return; 128*ce78bad3SBarry Smith *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd)); 129*ce78bad3SBarry Smith } 130*ce78bad3SBarry Smith 131*ce78bad3SBarry Smith PETSC_EXTERN void dmcompositerestorelocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd)) 132*ce78bad3SBarry Smith { 133*ce78bad3SBarry Smith IS *ais; 134*ce78bad3SBarry Smith PetscInt ndm; 135*ce78bad3SBarry Smith 136*ce78bad3SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd)); 137*ce78bad3SBarry Smith if (*ierr) return; 138*ce78bad3SBarry Smith *ierr = DMCompositeGetNumberDM(*dm, &ndm); 139*ce78bad3SBarry Smith for (PetscInt i = 0; i < ndm; i++) { 140*ce78bad3SBarry Smith *ierr = ISDestroy(&ais[i]); 141*ce78bad3SBarry Smith if (*ierr) return; 142*ce78bad3SBarry Smith } 14373e31fe2SJed Brown *ierr = PetscFree(ais); 144*ce78bad3SBarry Smith if (*ierr) return; 145*ce78bad3SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd)); 14673e31fe2SJed Brown } 147