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