xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 19caf8f3c08b1f0ca9f5469bde385c134aa76c82)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
3c6db04a5SJed Brown #include <petscdmcomposite.h>
447c6ae99SBarry Smith 
547c6ae99SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
647c6ae99SBarry Smith #define dmcompositegetentries1_      DMCOMPOSITEGETENTRIES1
747c6ae99SBarry Smith #define dmcompositegetentries2_      DMCOMPOSITEGETENTRIES2
847c6ae99SBarry Smith #define dmcompositegetentries3_      DMCOMPOSITEGETENTRIES3
947c6ae99SBarry Smith #define dmcompositegetentries4_      DMCOMPOSITEGETENTRIES4
1047c6ae99SBarry Smith #define dmcompositegetentries5_      DMCOMPOSITEGETENTRIES5
1147c6ae99SBarry Smith #define dmcompositeadddm_            DMCOMPOSITEADDDM
1247c6ae99SBarry Smith #define dmcompositedestroy_          DMCOMPOSITEDESTROY
1347c6ae99SBarry Smith #define dmcompositegetaccess4_       DMCOMPOSITEGETACCESS4
1447c6ae99SBarry Smith #define dmcompositescatter4_         DMCOMPOSITESCATTER4
1547c6ae99SBarry Smith #define dmcompositerestoreaccess4_   DMCOMPOSITERESTOREACCESS4
1647c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ DMCOMPOSITEGETLOCALVECTORS4
1747c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_  DMCOMPOSITERESTORELOCALVECTORS4
1873e31fe2SJed Brown #define dmcompositegetglobaliss_     DMCOMPOSITEGETGLOBALISS
1973e31fe2SJed Brown #define dmcompositegetlocaliss_      DMCOMPOSITEGETLOCALISS
20f73e5cebSJed Brown #define dmcompositegetaccessarray_   DMCOMPOSITEGETACCESSARRAY
21f73e5cebSJed Brown #define dmcompositerestoreaccessarray_  DMCOMPOSITERESTOREACCESSARRAY
227ac2b803SAlex Fikl #define dmcompositegetlocalaccessarray_ DMCOMPOSITEGETLOCALACCESSARRAY
237ac2b803SAlex Fikl #define dmcompositerestorelocalaccessarray_ DMCOMPOSITERESTORELOCALACCESSARRAY
2447c6ae99SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
2547c6ae99SBarry Smith #define dmcompositegetentries1_      dmcompositegetentries1
2647c6ae99SBarry Smith #define dmcompositegetentries2_      dmcompositegetentries2
2747c6ae99SBarry Smith #define dmcompositegetentries3_      dmcompositegetentries3
2847c6ae99SBarry Smith #define dmcompositegetentries4_      dmcompositegetentries4
2947c6ae99SBarry Smith #define dmcompositegetentries5_      dmcompositegetentries5
3047c6ae99SBarry Smith #define dmcompositeadddm_            dmcompositeadddm
3147c6ae99SBarry Smith #define dmcompositedestroy_          dmcompositedestroy
3247c6ae99SBarry Smith #define dmcompositegetaccess4_       dmcompositegetaccess4
3347c6ae99SBarry Smith #define dmcompositescatter4_         dmcompositescatter4
3447c6ae99SBarry Smith #define dmcompositerestoreaccess4_   dmcompositerestoreaccess4
3547c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ dmcompositegetlocalvectors4
3647c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4
3773e31fe2SJed Brown #define dmcompositegetglobaliss_     dmcompositegetglobaliss
3873e31fe2SJed Brown #define dmcompositegetlocaliss_      dmcompositegetlocaliss
39f73e5cebSJed Brown #define dmcompositegetaccessarray_   dmcompositegetaccessarray
40f73e5cebSJed Brown #define dmcompositerestoreaccessarray_  dmcompositerestoreaccessarray
417ac2b803SAlex Fikl #define dmcompositegetlocalaccessarray_ dmcompositegetlocalaccessarray
427ac2b803SAlex Fikl #define dmcompositerestorelocalaccessarray_ dmcompositerestorelocalaccessarray
4347c6ae99SBarry Smith #endif
4447c6ae99SBarry Smith 
45*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries1_(DM *dm,DM *da1,PetscErrorCode *ierr)
4647c6ae99SBarry Smith {
4747c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1);
4847c6ae99SBarry Smith }
4947c6ae99SBarry Smith 
50*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries2_(DM *dm,DM *da1,DM *da2,PetscErrorCode *ierr)
5147c6ae99SBarry Smith {
5247c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2);
5347c6ae99SBarry Smith }
5447c6ae99SBarry Smith 
55*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries3_(DM *dm,DM *da1,DM *da2,DM *da3,PetscErrorCode *ierr)
5647c6ae99SBarry Smith {
5747c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3);
5847c6ae99SBarry Smith }
5947c6ae99SBarry Smith 
60*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries4_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,PetscErrorCode *ierr)
6147c6ae99SBarry Smith {
6247c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4);
6347c6ae99SBarry Smith }
6447c6ae99SBarry Smith 
65*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetentries5_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,DM *da5,PetscErrorCode *ierr)
6647c6ae99SBarry Smith {
6747c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4,da5);
6847c6ae99SBarry Smith }
6947c6ae99SBarry Smith 
70*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
7147c6ae99SBarry Smith {
7247c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
7347c6ae99SBarry Smith   *ierr = DMCompositeGetAccess(*dm,*v,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith 
76*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositescatter4_(DM *dm,Vec *v,void *v1,void *p1,void *v2,void *p2,PetscErrorCode *ierr)
7747c6ae99SBarry Smith {
7847c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
7947c6ae99SBarry Smith   *ierr = DMCompositeScatter(*dm,*v,*vv1,(PetscScalar*)p1,*vv2,(PetscScalar*)p2);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
82*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestoreaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
8347c6ae99SBarry Smith {
8447c6ae99SBarry Smith   *ierr = DMCompositeRestoreAccess(*dm,*v,(Vec*)v1,0,(Vec*)v2,0);
8547c6ae99SBarry Smith }
8647c6ae99SBarry Smith 
87*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetlocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
8847c6ae99SBarry Smith {
8947c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
9047c6ae99SBarry Smith   *ierr = DMCompositeGetLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
9147c6ae99SBarry Smith }
9247c6ae99SBarry Smith 
93*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestorelocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
9447c6ae99SBarry Smith {
9547c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
9647c6ae99SBarry Smith   *ierr = DMCompositeRestoreLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith 
99*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetglobaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
10073e31fe2SJed Brown {
10173e31fe2SJed Brown   IS      *ais;
10273e31fe2SJed Brown   PetscInt i,ndm;
10373e31fe2SJed Brown   *ierr = DMCompositeGetGlobalISs(*dm,&ais); if (*ierr) return;
10473e31fe2SJed Brown   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
10573e31fe2SJed Brown   for (i=0; i<ndm; i++) iss[i] = ais[i];
10673e31fe2SJed Brown   *ierr = PetscFree(ais);
10773e31fe2SJed Brown }
10873e31fe2SJed Brown 
109*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetlocaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
11073e31fe2SJed Brown {
11173e31fe2SJed Brown   IS      *ais;
11273e31fe2SJed Brown   PetscInt i,ndm;
11373e31fe2SJed Brown   *ierr = DMCompositeGetLocalISs(*dm,&ais); if (*ierr) return;
11473e31fe2SJed Brown   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
11573e31fe2SJed Brown   for (i=0; i<ndm; i++) iss[i] = ais[i];
11673e31fe2SJed Brown   *ierr = PetscFree(ais);
11773e31fe2SJed Brown }
118f73e5cebSJed Brown 
119*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
120f73e5cebSJed Brown {
121f73e5cebSJed Brown   CHKFORTRANNULLINTEGER(wanted);
122f73e5cebSJed Brown   *ierr = DMCompositeGetAccessArray(*dm,*gvec,*n,wanted,vecs);
123f73e5cebSJed Brown }
124f73e5cebSJed Brown 
125*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestoreaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
126f73e5cebSJed Brown {
127f73e5cebSJed Brown   CHKFORTRANNULLINTEGER(wanted);
128f73e5cebSJed Brown   *ierr = DMCompositeRestoreAccessArray(*dm,*gvec,*n,wanted,vecs);
129f73e5cebSJed Brown }
1307ac2b803SAlex Fikl 
131*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositegetlocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
1327ac2b803SAlex Fikl {
1337ac2b803SAlex Fikl   CHKFORTRANNULLINTEGER(wanted);
1347ac2b803SAlex Fikl   *ierr = DMCompositeGetLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
1357ac2b803SAlex Fikl }
1367ac2b803SAlex Fikl 
137*19caf8f3SSatish Balay PETSC_EXTERN void dmcompositerestorelocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
1387ac2b803SAlex Fikl {
1397ac2b803SAlex Fikl   CHKFORTRANNULLINTEGER(wanted);
1407ac2b803SAlex Fikl   *ierr = DMCompositeRestoreLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
1417ac2b803SAlex Fikl }
142