xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
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 dmcompositegetaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
71 {
72   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
73   *ierr = DMCompositeGetAccess(*dm,*v,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
74 }
75 
76 PETSC_EXTERN void PETSC_STDCALL dmcompositescatter4_(DM *dm,Vec *v,void *v1,void *p1,void *v2,void *p2,PetscErrorCode *ierr)
77 {
78   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
79   *ierr = DMCompositeScatter(*dm,*v,*vv1,(PetscScalar*)p1,*vv2,(PetscScalar*)p2);
80 }
81 
82 PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
83 {
84   *ierr = DMCompositeRestoreAccess(*dm,*v,(Vec*)v1,0,(Vec*)v2,0);
85 }
86 
87 PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
88 {
89   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
90   *ierr = DMCompositeGetLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
91 }
92 
93 PETSC_EXTERN void PETSC_STDCALL dmcompositerestorelocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
94 {
95   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
96   *ierr = DMCompositeRestoreLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
97 }
98 
99 PETSC_EXTERN void PETSC_STDCALL dmcompositegetglobaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
100 {
101   IS      *ais;
102   PetscInt i,ndm;
103   *ierr = DMCompositeGetGlobalISs(*dm,&ais); if (*ierr) return;
104   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
105   for (i=0; i<ndm; i++) iss[i] = ais[i];
106   *ierr = PetscFree(ais);
107 }
108 
109 PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
110 {
111   IS      *ais;
112   PetscInt i,ndm;
113   *ierr = DMCompositeGetLocalISs(*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 dmcompositegetaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
120 {
121   CHKFORTRANNULLINTEGER(wanted);
122   *ierr = DMCompositeGetAccessArray(*dm,*gvec,*n,wanted,vecs);
123 }
124 
125 PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
126 {
127   CHKFORTRANNULLINTEGER(wanted);
128   *ierr = DMCompositeRestoreAccessArray(*dm,*gvec,*n,wanted,vecs);
129 }
130 
131 PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
132 {
133   CHKFORTRANNULLINTEGER(wanted);
134   *ierr = DMCompositeGetLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
135 }
136 
137 PETSC_EXTERN void PETSC_STDCALL dmcompositerestorelocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
138 {
139   CHKFORTRANNULLINTEGER(wanted);
140   *ierr = DMCompositeRestoreLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
141 }
142