xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 047240e14af00aad1ef65e96f6fface8924f7f7e)
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 #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 dmcompositecreate_           dmcompositecreate
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 #endif
42 
43 EXTERN_C_BEGIN
44 
45 void PETSC_STDCALL dmcompositegetentries1_(DM *dm,DM *da1,PetscErrorCode *ierr)
46 {
47   *ierr = DMCompositeGetEntries(*dm,da1);
48 }
49 
50 void PETSC_STDCALL dmcompositegetentries2_(DM *dm,DM *da1,DM *da2,PetscErrorCode *ierr)
51 {
52   *ierr = DMCompositeGetEntries(*dm,da1,da2);
53 }
54 
55 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 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 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 void PETSC_STDCALL  dmcompositecreate_(MPI_Fint * comm,DM *A, int *ierr)
71 {
72   *ierr = DMCompositeCreate(MPI_Comm_f2c(*(comm)),A);
73 }
74 
75 void PETSC_STDCALL dmcompositeadddm_(DM *dm,DM *da,PetscErrorCode *ierr)
76 {
77   *ierr = DMCompositeAddDM(*dm,*da);
78 }
79 
80 void PETSC_STDCALL dmcompositedestroy_(DM *dm,PetscErrorCode *ierr)
81 {
82   *ierr = DMDestroy(dm);
83 }
84 
85 void PETSC_STDCALL dmcompositegetaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
86 {
87   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
88   *ierr = DMCompositeGetAccess(*dm,*v,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
89 }
90 
91 void PETSC_STDCALL dmcompositescatter4_(DM *dm,Vec *v,void *v1,void *p1,void *v2,void *p2,PetscErrorCode *ierr)
92 {
93   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
94   *ierr = DMCompositeScatter(*dm,*v,*vv1,(PetscScalar*)p1,*vv2,(PetscScalar*)p2);
95 }
96 
97 void PETSC_STDCALL dmcompositerestoreaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
98 {
99   *ierr = DMCompositeRestoreAccess(*dm,*v,(Vec*)v1,0,(Vec*)v2,0);
100 }
101 
102 void PETSC_STDCALL dmcompositegetlocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
103 {
104   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
105   *ierr = DMCompositeGetLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
106 }
107 
108 void PETSC_STDCALL dmcompositerestorelocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
109 {
110   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
111   *ierr = DMCompositeRestoreLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
112 }
113 
114 EXTERN_C_END
115 
116 PETSC_EXTERN_C 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_C 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_C 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_C 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