xref: /petsc/src/dm/impls/composite/ftn-custom/zfddaf.c (revision 7ac2b803e3b92ba6d298600f281e6c25f3c4a7f6)
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 dmcompositecreate_           DMCOMPOSITECREATE
1247c6ae99SBarry Smith #define dmcompositeadddm_            DMCOMPOSITEADDDM
1347c6ae99SBarry Smith #define dmcompositedestroy_          DMCOMPOSITEDESTROY
1447c6ae99SBarry Smith #define dmcompositegetaccess4_       DMCOMPOSITEGETACCESS4
1547c6ae99SBarry Smith #define dmcompositescatter4_         DMCOMPOSITESCATTER4
1647c6ae99SBarry Smith #define dmcompositerestoreaccess4_   DMCOMPOSITERESTOREACCESS4
1747c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ DMCOMPOSITEGETLOCALVECTORS4
1847c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_  DMCOMPOSITERESTORELOCALVECTORS4
1973e31fe2SJed Brown #define dmcompositegetglobaliss_     DMCOMPOSITEGETGLOBALISS
2073e31fe2SJed Brown #define dmcompositegetlocaliss_      DMCOMPOSITEGETLOCALISS
21f73e5cebSJed Brown #define dmcompositegetaccessarray_   DMCOMPOSITEGETACCESSARRAY
22f73e5cebSJed Brown #define dmcompositerestoreaccessarray_  DMCOMPOSITERESTOREACCESSARRAY
23*7ac2b803SAlex Fikl #define dmcompositegetlocalaccessarray_ DMCOMPOSITEGETLOCALACCESSARRAY
24*7ac2b803SAlex Fikl #define dmcompositerestorelocalaccessarray_ DMCOMPOSITERESTORELOCALACCESSARRAY
2547c6ae99SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
2647c6ae99SBarry Smith #define dmcompositegetentries1_      dmcompositegetentries1
2747c6ae99SBarry Smith #define dmcompositegetentries2_      dmcompositegetentries2
2847c6ae99SBarry Smith #define dmcompositegetentries3_      dmcompositegetentries3
2947c6ae99SBarry Smith #define dmcompositegetentries4_      dmcompositegetentries4
3047c6ae99SBarry Smith #define dmcompositegetentries5_      dmcompositegetentries5
3147c6ae99SBarry Smith #define dmcompositecreate_           dmcompositecreate
3247c6ae99SBarry Smith #define dmcompositeadddm_            dmcompositeadddm
3347c6ae99SBarry Smith #define dmcompositedestroy_          dmcompositedestroy
3447c6ae99SBarry Smith #define dmcompositegetaccess4_       dmcompositegetaccess4
3547c6ae99SBarry Smith #define dmcompositescatter4_         dmcompositescatter4
3647c6ae99SBarry Smith #define dmcompositerestoreaccess4_   dmcompositerestoreaccess4
3747c6ae99SBarry Smith #define dmcompositegetlocalvectors4_ dmcompositegetlocalvectors4
3847c6ae99SBarry Smith #define dmcompositerestorelocalvectors4_ dmcompositerestorelocalvectors4
3973e31fe2SJed Brown #define dmcompositegetglobaliss_     dmcompositegetglobaliss
4073e31fe2SJed Brown #define dmcompositegetlocaliss_      dmcompositegetlocaliss
41f73e5cebSJed Brown #define dmcompositegetaccessarray_   dmcompositegetaccessarray
42f73e5cebSJed Brown #define dmcompositerestoreaccessarray_  dmcompositerestoreaccessarray
43*7ac2b803SAlex Fikl #define dmcompositegetlocalaccessarray_ dmcompositegetlocalaccessarray
44*7ac2b803SAlex Fikl #define dmcompositerestorelocalaccessarray_ dmcompositerestorelocalaccessarray
4547c6ae99SBarry Smith #endif
4647c6ae99SBarry Smith 
478cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries1_(DM *dm,DM *da1,PetscErrorCode *ierr)
4847c6ae99SBarry Smith {
4947c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1);
5047c6ae99SBarry Smith }
5147c6ae99SBarry Smith 
528cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries2_(DM *dm,DM *da1,DM *da2,PetscErrorCode *ierr)
5347c6ae99SBarry Smith {
5447c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2);
5547c6ae99SBarry Smith }
5647c6ae99SBarry Smith 
578cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries3_(DM *dm,DM *da1,DM *da2,DM *da3,PetscErrorCode *ierr)
5847c6ae99SBarry Smith {
5947c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3);
6047c6ae99SBarry Smith }
6147c6ae99SBarry Smith 
628cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries4_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,PetscErrorCode *ierr)
6347c6ae99SBarry Smith {
6447c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4);
6547c6ae99SBarry Smith }
6647c6ae99SBarry Smith 
678cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetentries5_(DM *dm,DM *da1,DM *da2,DM *da3,DM *da4,DM *da5,PetscErrorCode *ierr)
6847c6ae99SBarry Smith {
6947c6ae99SBarry Smith   *ierr = DMCompositeGetEntries(*dm,da1,da2,da3,da4,da5);
7047c6ae99SBarry Smith }
7147c6ae99SBarry Smith 
728cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositecreate_(MPI_Fint * comm,DM *A, int *ierr)
73c30958fdSKarl Rupp {
7447c6ae99SBarry Smith   *ierr = DMCompositeCreate(MPI_Comm_f2c(*(comm)),A);
7547c6ae99SBarry Smith }
7647c6ae99SBarry Smith 
778cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositeadddm_(DM *dm,DM *da,PetscErrorCode *ierr)
7847c6ae99SBarry Smith {
7947c6ae99SBarry Smith   *ierr = DMCompositeAddDM(*dm,*da);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
828cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositedestroy_(DM *dm,PetscErrorCode *ierr)
8347c6ae99SBarry Smith {
84fcfd50ebSBarry Smith   *ierr = DMDestroy(dm);
8547c6ae99SBarry Smith }
8647c6ae99SBarry Smith 
878cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
8847c6ae99SBarry Smith {
8947c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
9047c6ae99SBarry Smith   *ierr = DMCompositeGetAccess(*dm,*v,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
9147c6ae99SBarry Smith }
9247c6ae99SBarry Smith 
938cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositescatter4_(DM *dm,Vec *v,void *v1,void *p1,void *v2,void *p2,PetscErrorCode *ierr)
9447c6ae99SBarry Smith {
9547c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
9647c6ae99SBarry Smith   *ierr = DMCompositeScatter(*dm,*v,*vv1,(PetscScalar*)p1,*vv2,(PetscScalar*)p2);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith 
998cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccess4_(DM *dm,Vec *v,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
10047c6ae99SBarry Smith {
10147c6ae99SBarry Smith   *ierr = DMCompositeRestoreAccess(*dm,*v,(Vec*)v1,0,(Vec*)v2,0);
10247c6ae99SBarry Smith }
10347c6ae99SBarry Smith 
1048cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
10547c6ae99SBarry Smith {
10647c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
10747c6ae99SBarry Smith   *ierr = DMCompositeGetLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
10847c6ae99SBarry Smith }
10947c6ae99SBarry Smith 
1108cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositerestorelocalvectors4_(DM *dm,void **v1,void **p1,void **v2,void **p2,PetscErrorCode *ierr)
11147c6ae99SBarry Smith {
11247c6ae99SBarry Smith   Vec *vv1 = (Vec*)v1,*vv2 = (Vec*)v2;
11347c6ae99SBarry Smith   *ierr = DMCompositeRestoreLocalVectors(*dm,vv1,(PetscScalar*)p1,vv2,(PetscScalar*)p2);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
1168cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetglobaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
11773e31fe2SJed Brown {
11873e31fe2SJed Brown   IS      *ais;
11973e31fe2SJed Brown   PetscInt i,ndm;
12073e31fe2SJed Brown   *ierr = DMCompositeGetGlobalISs(*dm,&ais); if (*ierr) return;
12173e31fe2SJed Brown   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
12273e31fe2SJed Brown   for (i=0; i<ndm; i++) iss[i] = ais[i];
12373e31fe2SJed Brown   *ierr = PetscFree(ais);
12473e31fe2SJed Brown }
12573e31fe2SJed Brown 
1268cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocaliss_(DM *dm,IS *iss,PetscErrorCode *ierr)
12773e31fe2SJed Brown {
12873e31fe2SJed Brown   IS      *ais;
12973e31fe2SJed Brown   PetscInt i,ndm;
13073e31fe2SJed Brown   *ierr = DMCompositeGetLocalISs(*dm,&ais); if (*ierr) return;
13173e31fe2SJed Brown   *ierr = DMCompositeGetNumberDM(*dm,&ndm); if (*ierr) return;
13273e31fe2SJed Brown   for (i=0; i<ndm; i++) iss[i] = ais[i];
13373e31fe2SJed Brown   *ierr = PetscFree(ais);
13473e31fe2SJed Brown }
135f73e5cebSJed Brown 
1368cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositegetaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
137f73e5cebSJed Brown {
138f73e5cebSJed Brown   CHKFORTRANNULLINTEGER(wanted);
139f73e5cebSJed Brown   *ierr = DMCompositeGetAccessArray(*dm,*gvec,*n,wanted,vecs);
140f73e5cebSJed Brown }
141f73e5cebSJed Brown 
1428cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccessarray_(DM *dm,Vec *gvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
143f73e5cebSJed Brown {
144f73e5cebSJed Brown   CHKFORTRANNULLINTEGER(wanted);
145f73e5cebSJed Brown   *ierr = DMCompositeRestoreAccessArray(*dm,*gvec,*n,wanted,vecs);
146f73e5cebSJed Brown }
147*7ac2b803SAlex Fikl 
148*7ac2b803SAlex Fikl PETSC_EXTERN void PETSC_STDCALL dmcompositegetlocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
149*7ac2b803SAlex Fikl {
150*7ac2b803SAlex Fikl   CHKFORTRANNULLINTEGER(wanted);
151*7ac2b803SAlex Fikl   *ierr = DMCompositeGetLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
152*7ac2b803SAlex Fikl }
153*7ac2b803SAlex Fikl 
154*7ac2b803SAlex Fikl PETSC_EXTERN void PETSC_STDCALL dmcompositerestorelocalaccessarray_(DM *dm,Vec *lvec,PetscInt *n,const PetscInt *wanted,Vec *vecs,PetscErrorCode *ierr)
155*7ac2b803SAlex Fikl {
156*7ac2b803SAlex Fikl   CHKFORTRANNULLINTEGER(wanted);
157*7ac2b803SAlex Fikl   *ierr = DMCompositeRestoreLocalAccessArray(*dm,*lvec,*n,wanted,vecs);
158*7ac2b803SAlex Fikl }
159