xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 09b68a49ed2854d1e4985cc2aa6af33c7c4e69b3)
1 #include <petscdmcomposite.h>
2 #include <petsc/private/ftnimpl.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define dmcompositegetentries1_          DMCOMPOSITEGETENTRIES1
6   #define dmcompositegetentries2_          DMCOMPOSITEGETENTRIES2
7   #define dmcompositegetentries3_          DMCOMPOSITEGETENTRIES3
8   #define dmcompositegetentries4_          DMCOMPOSITEGETENTRIES4
9   #define dmcompositegetentries5_          DMCOMPOSITEGETENTRIES5
10   #define dmcompositegetaccess4_           DMCOMPOSITEGETACCESS4
11   #define dmcompositescatter4_             DMCOMPOSITESCATTER4
12   #define dmcompositerestoreaccess4_       DMCOMPOSITERESTOREACCESS4
13   #define dmcompositegetlocalvectors4_     DMCOMPOSITEGETLOCALVECTORS4
14   #define dmcompositerestorelocalvectors4_ DMCOMPOSITERESTORELOCALVECTORS4
15   #define dmcompositegetglobaliss_         DMCOMPOSITEGETGLOBALISS
16   #define dmcompositerestoreglobaliss_     DMCOMPOSITERESTOREGLOBALISS
17   #define dmcompositegetlocaliss_          DMCOMPOSITEGETLOCALISS
18   #define dmcompositerestorelocaliss_      DMCOMPOSITERESTORELOCALISS
19 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
20   #define dmcompositegetentries1_          dmcompositegetentries1
21   #define dmcompositegetentries2_          dmcompositegetentries2
22   #define dmcompositegetentries3_          dmcompositegetentries3
23   #define dmcompositegetentries4_          dmcompositegetentries4
24   #define dmcompositegetentries5_          dmcompositegetentries5
25   #define dmcompositegetaccess4_           dmcompositegetaccess4
26   #define dmcompositescatter4_             dmcompositescatter4
27   #define dmcompositerestoreaccess4_       dmcompositerestoreaccess4
28   #define dmcompositegetlocalvectors4_     dmcompositegetlocalvectors4
29   #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4
30   #define dmcompositegetglobaliss_         dmcompositegetglobaliss
31   #define dmcompositerestoreglobaliss_     dmcompositerestoreglobaliss
32   #define dmcompositegetlocaliss_          dmcompositegetlocaliss
33   #define dmcompositerestorelocaliss_      dmcompositerestorelocaliss
34 #endif
35 
dmcompositegetentries1_(DM * dm,DM * da1,PetscErrorCode * ierr)36 PETSC_EXTERN void dmcompositegetentries1_(DM *dm, DM *da1, PetscErrorCode *ierr)
37 {
38   *ierr = DMCompositeGetEntries(*dm, da1);
39 }
40 
dmcompositegetentries2_(DM * dm,DM * da1,DM * da2,PetscErrorCode * ierr)41 PETSC_EXTERN void dmcompositegetentries2_(DM *dm, DM *da1, DM *da2, PetscErrorCode *ierr)
42 {
43   *ierr = DMCompositeGetEntries(*dm, da1, da2);
44 }
45 
dmcompositegetentries3_(DM * dm,DM * da1,DM * da2,DM * da3,PetscErrorCode * ierr)46 PETSC_EXTERN void dmcompositegetentries3_(DM *dm, DM *da1, DM *da2, DM *da3, PetscErrorCode *ierr)
47 {
48   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3);
49 }
50 
dmcompositegetentries4_(DM * dm,DM * da1,DM * da2,DM * da3,DM * da4,PetscErrorCode * ierr)51 PETSC_EXTERN void dmcompositegetentries4_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, PetscErrorCode *ierr)
52 {
53   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4);
54 }
55 
dmcompositegetentries5_(DM * dm,DM * da1,DM * da2,DM * da3,DM * da4,DM * da5,PetscErrorCode * ierr)56 PETSC_EXTERN void dmcompositegetentries5_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, DM *da5, PetscErrorCode *ierr)
57 {
58   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4, da5);
59 }
60 
dmcompositegetglobaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))61 PETSC_EXTERN void dmcompositegetglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
62 {
63   IS      *ais;
64   PetscInt ndm;
65 
66   *ierr = DMCompositeGetGlobalISs(*dm, &ais);
67   if (*ierr) return;
68   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
69   if (*ierr) return;
70   *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd));
71 }
72 
dmcompositerestoreglobaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))73 PETSC_EXTERN void dmcompositerestoreglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
74 {
75   IS      *ais;
76   PetscInt ndm;
77 
78   *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd));
79   if (*ierr) return;
80   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
81   for (PetscInt i = 0; i < ndm; i++) {
82     *ierr = ISDestroy(&ais[i]);
83     if (*ierr) return;
84   }
85   *ierr = PetscFree(ais);
86   if (*ierr) return;
87   *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd));
88 }
89 
dmcompositegetlocaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))90 PETSC_EXTERN void dmcompositegetlocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
91 {
92   IS      *ais;
93   PetscInt ndm;
94 
95   *ierr = DMCompositeGetLocalISs(*dm, &ais);
96   if (*ierr) return;
97   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
98   if (*ierr) return;
99   *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd));
100 }
101 
dmcompositerestorelocaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))102 PETSC_EXTERN void dmcompositerestorelocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
103 {
104   IS      *ais;
105   PetscInt ndm;
106 
107   *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd));
108   if (*ierr) return;
109   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
110   for (PetscInt i = 0; i < ndm; i++) {
111     *ierr = ISDestroy(&ais[i]);
112     if (*ierr) return;
113   }
114   *ierr = PetscFree(ais);
115   if (*ierr) return;
116   *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd));
117 }
118