1 #include <petsc/private/fortranimpl.h> 2 #include <petscdmcomposite.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define dmcompositegetentries1_ DMCOMPOSITEGETENTRIES1 6 #define dmcompositegetentries2_ DMCOMPOSITEGETENTRIES2 7 #define dmcompositegetentries3_ DMCOMPOSITEGETENTRIES3 8 #define dmcompositegetentries4_ DMCOMPOSITEGETENTRIES4 9 #define dmcompositegetentries5_ DMCOMPOSITEGETENTRIES5 10 #define dmcompositeadddm_ DMCOMPOSITEADDDM 11 #define dmcompositedestroy_ DMCOMPOSITEDESTROY 12 #define dmcompositegetaccess4_ DMCOMPOSITEGETACCESS4 13 #define dmcompositescatter4_ DMCOMPOSITESCATTER4 14 #define dmcompositerestoreaccess4_ DMCOMPOSITERESTOREACCESS4 15 #define dmcompositegetlocalvectors4_ DMCOMPOSITEGETLOCALVECTORS4 16 #define dmcompositerestorelocalvectors4_ DMCOMPOSITERESTORELOCALVECTORS4 17 #define dmcompositegetglobaliss_ DMCOMPOSITEGETGLOBALISS 18 #define dmcompositegetlocaliss_ DMCOMPOSITEGETLOCALISS 19 #define dmcompositegetaccessarray_ DMCOMPOSITEGETACCESSARRAY 20 #define dmcompositerestoreaccessarray_ DMCOMPOSITERESTOREACCESSARRAY 21 #define dmcompositegetlocalaccessarray_ DMCOMPOSITEGETLOCALACCESSARRAY 22 #define dmcompositerestorelocalaccessarray_ DMCOMPOSITERESTORELOCALACCESSARRAY 23 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 24 #define dmcompositegetentries1_ dmcompositegetentries1 25 #define dmcompositegetentries2_ dmcompositegetentries2 26 #define dmcompositegetentries3_ dmcompositegetentries3 27 #define dmcompositegetentries4_ dmcompositegetentries4 28 #define dmcompositegetentries5_ dmcompositegetentries5 29 #define dmcompositeadddm_ dmcompositeadddm 30 #define dmcompositedestroy_ dmcompositedestroy 31 #define dmcompositegetaccess4_ dmcompositegetaccess4 32 #define dmcompositescatter4_ dmcompositescatter4 33 #define dmcompositerestoreaccess4_ dmcompositerestoreaccess4 34 #define dmcompositegetlocalvectors4_ dmcompositegetlocalvectors4 35 #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4 36 #define dmcompositegetglobaliss_ dmcompositegetglobaliss 37 #define dmcompositegetlocaliss_ dmcompositegetlocaliss 38 #define dmcompositegetaccessarray_ dmcompositegetaccessarray 39 #define dmcompositerestoreaccessarray_ dmcompositerestoreaccessarray 40 #define dmcompositegetlocalaccessarray_ dmcompositegetlocalaccessarray 41 #define dmcompositerestorelocalaccessarray_ dmcompositerestorelocalaccessarray 42 #endif 43 44 PETSC_EXTERN void dmcompositegetentries1_(DM *dm, DM *da1, PetscErrorCode *ierr) 45 { 46 *ierr = DMCompositeGetEntries(*dm, da1); 47 } 48 49 PETSC_EXTERN void dmcompositegetentries2_(DM *dm, DM *da1, DM *da2, PetscErrorCode *ierr) 50 { 51 *ierr = DMCompositeGetEntries(*dm, da1, da2); 52 } 53 54 PETSC_EXTERN void dmcompositegetentries3_(DM *dm, DM *da1, DM *da2, DM *da3, PetscErrorCode *ierr) 55 { 56 *ierr = DMCompositeGetEntries(*dm, da1, da2, da3); 57 } 58 59 PETSC_EXTERN void dmcompositegetentries4_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, PetscErrorCode *ierr) 60 { 61 *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4); 62 } 63 64 PETSC_EXTERN void dmcompositegetentries5_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, DM *da5, PetscErrorCode *ierr) 65 { 66 *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4, da5); 67 } 68 69 PETSC_EXTERN void dmcompositegetaccess4_(DM *dm, Vec *v, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 70 { 71 Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 72 *ierr = DMCompositeGetAccess(*dm, *v, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 73 } 74 75 PETSC_EXTERN void dmcompositescatter4_(DM *dm, Vec *v, void *v1, void *p1, void *v2, void *p2, PetscErrorCode *ierr) 76 { 77 Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 78 *ierr = DMCompositeScatter(*dm, *v, *vv1, (PetscScalar *)p1, *vv2, (PetscScalar *)p2); 79 } 80 81 PETSC_EXTERN void dmcompositerestoreaccess4_(DM *dm, Vec *v, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 82 { 83 *ierr = DMCompositeRestoreAccess(*dm, *v, (Vec *)v1, 0, (Vec *)v2, 0); 84 } 85 86 PETSC_EXTERN void dmcompositegetlocalvectors4_(DM *dm, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 87 { 88 Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 89 *ierr = DMCompositeGetLocalVectors(*dm, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 90 } 91 92 PETSC_EXTERN void dmcompositerestorelocalvectors4_(DM *dm, void **v1, void **p1, void **v2, void **p2, PetscErrorCode *ierr) 93 { 94 Vec *vv1 = (Vec *)v1, *vv2 = (Vec *)v2; 95 *ierr = DMCompositeRestoreLocalVectors(*dm, vv1, (PetscScalar *)p1, vv2, (PetscScalar *)p2); 96 } 97 98 PETSC_EXTERN void dmcompositegetglobaliss_(DM *dm, IS *iss, PetscErrorCode *ierr) 99 { 100 IS *ais; 101 PetscInt i, ndm; 102 *ierr = DMCompositeGetGlobalISs(*dm, &ais); 103 if (*ierr) return; 104 *ierr = DMCompositeGetNumberDM(*dm, &ndm); 105 if (*ierr) return; 106 for (i = 0; i < ndm; i++) iss[i] = ais[i]; 107 *ierr = PetscFree(ais); 108 } 109 110 PETSC_EXTERN void dmcompositegetlocaliss_(DM *dm, IS *iss, PetscErrorCode *ierr) 111 { 112 IS *ais; 113 PetscInt i, ndm; 114 *ierr = DMCompositeGetLocalISs(*dm, &ais); 115 if (*ierr) return; 116 *ierr = DMCompositeGetNumberDM(*dm, &ndm); 117 if (*ierr) return; 118 for (i = 0; i < ndm; i++) iss[i] = ais[i]; 119 *ierr = PetscFree(ais); 120 } 121 122 PETSC_EXTERN void dmcompositegetaccessarray_(DM *dm, Vec *gvec, PetscInt *n, const PetscInt *wanted, Vec *vecs, PetscErrorCode *ierr) 123 { 124 CHKFORTRANNULLINTEGER(wanted); 125 *ierr = DMCompositeGetAccessArray(*dm, *gvec, *n, wanted, vecs); 126 } 127 128 PETSC_EXTERN void dmcompositerestoreaccessarray_(DM *dm, Vec *gvec, PetscInt *n, const PetscInt *wanted, Vec *vecs, PetscErrorCode *ierr) 129 { 130 CHKFORTRANNULLINTEGER(wanted); 131 *ierr = DMCompositeRestoreAccessArray(*dm, *gvec, *n, wanted, vecs); 132 } 133 134 PETSC_EXTERN void dmcompositegetlocalaccessarray_(DM *dm, Vec *lvec, PetscInt *n, const PetscInt *wanted, Vec *vecs, PetscErrorCode *ierr) 135 { 136 CHKFORTRANNULLINTEGER(wanted); 137 *ierr = DMCompositeGetLocalAccessArray(*dm, *lvec, *n, wanted, vecs); 138 } 139 140 PETSC_EXTERN void dmcompositerestorelocalaccessarray_(DM *dm, Vec *lvec, PetscInt *n, const PetscInt *wanted, Vec *vecs, PetscErrorCode *ierr) 141 { 142 CHKFORTRANNULLINTEGER(wanted); 143 *ierr = DMCompositeRestoreLocalAccessArray(*dm, *lvec, *n, wanted, vecs); 144 } 145