xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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 PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries1_(DM *dm,DM *da1,PetscErrorCode *ierr)
44 {
45   *ierr = DMCompositeGetEntries(*dm,da1);
46 }
47 
48 PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries2_(DM *dm,DM *da1,DM *da2,PetscErrorCode *ierr)
49 {
50   *ierr = DMCompositeGetEntries(*dm,da1,da2);
51 }
52 
53 PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries3_(DM *dm,DM *da1,DM *da2,DM *da3,PetscErrorCode *ierr)
54 {
55   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3);
56 }
57 
58 PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries4_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,PetscErrorCode *ierr)
59 {
60   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4);
61 }
62 
63 PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries5_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,DM *da5,PetscErrorCode *ierr)
64 {
65   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4,da5);
66 }
67 
68 PETSC_EXTERN void PETSC_STDCALL dmcompositecreate_(MPI_Fint * comm,DM *A, int *ierr)
69 {
70   *ierr = DMCompositeCreate(MPI_Comm_f2c(*(comm)),A);
71 }
72 
73 PETSC_EXTERN void PETSC_STDCALL dmcompositeadddm_(DM *dm,DM *da,PetscErrorCode *ierr)
74 {
75   *ierr = DMCompositeAddDM(*dm,*da);
76 }
77 
78 PETSC_EXTERN void PETSC_STDCALL dmcompositedestroy_(DM *dm,PetscErrorCode *ierr)
79 {
80   *ierr = DMDestroy(dm);
81 }
82 
83 PETSC_EXTERN void PETSC_STDCALL dmcompositegetaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
84 {
85   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
86   *ierr = DMCompositeGetAccess(*dm,*v,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
87 }
88 
89 PETSC_EXTERN void PETSC_STDCALL dmcompositescatter4_(DM *dm,Vec *v,void *v1,void *p1,void *v2,void *p2,PetscErrorCode *ierr)
90 {
91   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
92   *ierr = DMCompositeScatter(*dm,*v,*vv1,(PetscScalar*)p1,*vv2,(PetscScalar*)p2);
93 }
94 
95 PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
96 {
97   *ierr = DMCompositeRestoreAccess(*dm,*v,(Vec*)v1,0,(Vec*)v2,0);
98 }
99 
100 PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
101 {
102   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
103   *ierr = DMCompositeGetLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
104 }
105 
106 PETSC_EXTERN void PETSC_STDCALL dmcompositerestorelocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
107 {
108   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
109   *ierr = DMCompositeRestoreLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
110 }
111 
112 PETSC_EXTERN void PETSC_STDCALL dmcompositegetglobaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
113 {
114   IS      *ais;
115   PetscInt i,ndm;
116   *ierr = DMCompositeGetGlobalISs(*dm,&ais); if (*ierr) return;
117   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
118   for (i=0; i<ndm; i++) iss[i] = ais[i];
119   *ierr = PetscFree(ais);
120 }
121 
122 PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
123 {
124   IS      *ais;
125   PetscInt i,ndm;
126   *ierr = DMCompositeGetLocalISs(*dm,&ais); if (*ierr) return;
127   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
128   for (i=0; i<ndm; i++) iss[i] = ais[i];
129   *ierr = PetscFree(ais);
130 }
131 
132 PETSC_EXTERN void PETSC_STDCALL dmcompositegetaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
133 {
134   CHKFORTRANNULLINTEGER(wanted);
135   *ierr = DMCompositeGetAccessArray(*dm,*gvec,*n,wanted,vecs);
136 }
137 
138 PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
139 {
140   CHKFORTRANNULLINTEGER(wanted);
141   *ierr = DMCompositeRestoreAccessArray(*dm,*gvec,*n,wanted,vecs);
142 }
143