xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 09b68a49ed2854d1e4985cc2aa6af33c7c4e69b3)
1c6db04a5SJed Brown #include <petscdmcomposite.h>
2*6dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
547c6ae99SBarry Smith   #define dmcompositegetentries1_          DMCOMPOSITEGETENTRIES1
647c6ae99SBarry Smith   #define dmcompositegetentries2_          DMCOMPOSITEGETENTRIES2
747c6ae99SBarry Smith   #define dmcompositegetentries3_          DMCOMPOSITEGETENTRIES3
847c6ae99SBarry Smith   #define dmcompositegetentries4_          DMCOMPOSITEGETENTRIES4
947c6ae99SBarry Smith   #define dmcompositegetentries5_          DMCOMPOSITEGETENTRIES5
1047c6ae99SBarry Smith   #define dmcompositegetaccess4_           DMCOMPOSITEGETACCESS4
1147c6ae99SBarry Smith   #define dmcompositescatter4_             DMCOMPOSITESCATTER4
1247c6ae99SBarry Smith   #define dmcompositerestoreaccess4_       DMCOMPOSITERESTOREACCESS4
1347c6ae99SBarry Smith   #define dmcompositegetlocalvectors4_     DMCOMPOSITEGETLOCALVECTORS4
1447c6ae99SBarry Smith   #define dmcompositerestorelocalvectors4_ DMCOMPOSITERESTORELOCALVECTORS4
1573e31fe2SJed Brown   #define dmcompositegetglobaliss_         DMCOMPOSITEGETGLOBALISS
16ce78bad3SBarry Smith   #define dmcompositerestoreglobaliss_     DMCOMPOSITERESTOREGLOBALISS
1773e31fe2SJed Brown   #define dmcompositegetlocaliss_          DMCOMPOSITEGETLOCALISS
18ce78bad3SBarry Smith   #define dmcompositerestorelocaliss_      DMCOMPOSITERESTORELOCALISS
1947c6ae99SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
2047c6ae99SBarry Smith   #define dmcompositegetentries1_          dmcompositegetentries1
2147c6ae99SBarry Smith   #define dmcompositegetentries2_          dmcompositegetentries2
2247c6ae99SBarry Smith   #define dmcompositegetentries3_          dmcompositegetentries3
2347c6ae99SBarry Smith   #define dmcompositegetentries4_          dmcompositegetentries4
2447c6ae99SBarry Smith   #define dmcompositegetentries5_          dmcompositegetentries5
2547c6ae99SBarry Smith   #define dmcompositegetaccess4_           dmcompositegetaccess4
2647c6ae99SBarry Smith   #define dmcompositescatter4_             dmcompositescatter4
2747c6ae99SBarry Smith   #define dmcompositerestoreaccess4_       dmcompositerestoreaccess4
2847c6ae99SBarry Smith   #define dmcompositegetlocalvectors4_     dmcompositegetlocalvectors4
2947c6ae99SBarry Smith   #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4
3073e31fe2SJed Brown   #define dmcompositegetglobaliss_         dmcompositegetglobaliss
31ce78bad3SBarry Smith   #define dmcompositerestoreglobaliss_     dmcompositerestoreglobaliss
3273e31fe2SJed Brown   #define dmcompositegetlocaliss_          dmcompositegetlocaliss
33ce78bad3SBarry Smith   #define dmcompositerestorelocaliss_      dmcompositerestorelocaliss
3447c6ae99SBarry Smith #endif
3547c6ae99SBarry Smith 
dmcompositegetentries1_(DM * dm,DM * da1,PetscErrorCode * ierr)3619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries1_(DM *dm, DM *da1, PetscErrorCode *ierr)
3747c6ae99SBarry Smith {
3847c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm, da1);
3947c6ae99SBarry Smith }
4047c6ae99SBarry Smith 
dmcompositegetentries2_(DM * dm,DM * da1,DM * da2,PetscErrorCode * ierr)4119caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries2_(DM *dm, DM *da1, DM *da2, PetscErrorCode *ierr)
4247c6ae99SBarry Smith {
4347c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm, da1, da2);
4447c6ae99SBarry Smith }
4547c6ae99SBarry Smith 
dmcompositegetentries3_(DM * dm,DM * da1,DM * da2,DM * da3,PetscErrorCode * ierr)4619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries3_(DM *dm, DM *da1, DM *da2, DM *da3, PetscErrorCode *ierr)
4747c6ae99SBarry Smith {
4847c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3);
4947c6ae99SBarry Smith }
5047c6ae99SBarry Smith 
dmcompositegetentries4_(DM * dm,DM * da1,DM * da2,DM * da3,DM * da4,PetscErrorCode * ierr)5119caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries4_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, PetscErrorCode *ierr)
5247c6ae99SBarry Smith {
5347c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
dmcompositegetentries5_(DM * dm,DM * da1,DM * da2,DM * da3,DM * da4,DM * da5,PetscErrorCode * ierr)5619caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries5_(DM *dm, DM *da1, DM *da2, DM *da3, DM *da4, DM *da5, PetscErrorCode *ierr)
5747c6ae99SBarry Smith {
5847c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm, da1, da2, da3, da4, da5);
5947c6ae99SBarry Smith }
6047c6ae99SBarry Smith 
dmcompositegetglobaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))61ce78bad3SBarry Smith PETSC_EXTERN void dmcompositegetglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
6273e31fe2SJed Brown {
6373e31fe2SJed Brown   IS      *ais;
64ce78bad3SBarry Smith   PetscInt ndm;
65ce78bad3SBarry Smith 
665975b3b6SBarry Smith   *ierr = DMCompositeGetGlobalISs(*dm, &ais);
675975b3b6SBarry Smith   if (*ierr) return;
685975b3b6SBarry Smith   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
695975b3b6SBarry Smith   if (*ierr) return;
70ce78bad3SBarry Smith   *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd));
7173e31fe2SJed Brown }
7273e31fe2SJed Brown 
dmcompositerestoreglobaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))73ce78bad3SBarry Smith PETSC_EXTERN void dmcompositerestoreglobaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
7473e31fe2SJed Brown {
7573e31fe2SJed Brown   IS      *ais;
76ce78bad3SBarry Smith   PetscInt ndm;
77ce78bad3SBarry Smith 
78ce78bad3SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd));
79ce78bad3SBarry Smith   if (*ierr) return;
80ce78bad3SBarry Smith   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
81ce78bad3SBarry Smith   for (PetscInt i = 0; i < ndm; i++) {
82ce78bad3SBarry Smith     *ierr = ISDestroy(&ais[i]);
83ce78bad3SBarry Smith     if (*ierr) return;
84ce78bad3SBarry Smith   }
85ce78bad3SBarry Smith   *ierr = PetscFree(ais);
86ce78bad3SBarry Smith   if (*ierr) return;
87ce78bad3SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd));
88ce78bad3SBarry Smith }
89ce78bad3SBarry Smith 
dmcompositegetlocaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))90ce78bad3SBarry Smith PETSC_EXTERN void dmcompositegetlocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
91ce78bad3SBarry Smith {
92ce78bad3SBarry Smith   IS      *ais;
93ce78bad3SBarry Smith   PetscInt ndm;
94ce78bad3SBarry Smith 
955975b3b6SBarry Smith   *ierr = DMCompositeGetLocalISs(*dm, &ais);
965975b3b6SBarry Smith   if (*ierr) return;
975975b3b6SBarry Smith   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
985975b3b6SBarry Smith   if (*ierr) return;
99ce78bad3SBarry Smith   *ierr = F90Array1dCreate((void *)ais, MPIU_FORTRANADDR, 1, ndm, ptr PETSC_F90_2PTR_PARAM(ptrd));
100ce78bad3SBarry Smith }
101ce78bad3SBarry Smith 
dmcompositerestorelocaliss_(DM * dm,F90Array1d * ptr,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd))102ce78bad3SBarry Smith PETSC_EXTERN void dmcompositerestorelocaliss_(DM *dm, F90Array1d *ptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
103ce78bad3SBarry Smith {
104ce78bad3SBarry Smith   IS      *ais;
105ce78bad3SBarry Smith   PetscInt ndm;
106ce78bad3SBarry Smith 
107ce78bad3SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_FORTRANADDR, (void **)&ais PETSC_F90_2PTR_PARAM(ptrd));
108ce78bad3SBarry Smith   if (*ierr) return;
109ce78bad3SBarry Smith   *ierr = DMCompositeGetNumberDM(*dm, &ndm);
110ce78bad3SBarry Smith   for (PetscInt i = 0; i < ndm; i++) {
111ce78bad3SBarry Smith     *ierr = ISDestroy(&ais[i]);
112ce78bad3SBarry Smith     if (*ierr) return;
113ce78bad3SBarry Smith   }
11473e31fe2SJed Brown   *ierr = PetscFree(ais);
115ce78bad3SBarry Smith   if (*ierr) return;
116ce78bad3SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd));
11773e31fe2SJed Brown }
118