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