xref: /petsc/src/dm/interface/ftn-custom/zdmf.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1 #include <petsc/private/ftnimpl.h>
2 #include <petscdm.h>
3 #include <petscviewer.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6   #define dmcreatesuperdm_                       DMCREATESUPERDM
7   #define dmcreatefielddecompositiongetname_     DMCREATEFIELDDECOMPOSITIONGETNAME
8   #define dmcreatefielddecompositiongetisdm_     DMCREATEFIELDDECOMPOSITIONGETISDM
9   #define dmcreatefielddecompositionrestoreisdm_ DMCREATEFIELDDECOMPOSITIONRESTOREISDM
10 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
11   #define dmcreatesuperdm_                       dmreatesuperdm
12   #define dmcreatefielddecompositiongetname_     dmcreatefielddecompositiongetname
13   #define dmcreatefielddecompositiongetisdm_     dmcreatefielddecompositiongetisdm
14   #define dmcreatefielddecompositionrestoreisdm_ dmcreatefielddecompositionrestoreisdm
15 #endif
16 
dmcreatesuperdm_(DM dms[],PetscInt * len,IS *** is,DM * superdm,PetscErrorCode * ierr)17 PETSC_EXTERN void dmcreatesuperdm_(DM dms[], PetscInt *len, IS ***is, DM *superdm, PetscErrorCode *ierr)
18 {
19   *ierr = DMCreateSuperDM(dms, *len, *is, superdm);
20 }
21 
dmcreatefielddecompositiongetname_(DM * dm,PetscInt * i,char * name,PetscErrorCode * ierr,PETSC_FORTRAN_CHARLEN_T l_b)22 PETSC_EXTERN void dmcreatefielddecompositiongetname_(DM *dm, PetscInt *i, char *name, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T l_b)
23 {
24   PetscInt n;
25   char   **names;
26   *ierr = DMCreateFieldDecomposition(*dm, &n, &names, NULL, NULL);
27   if (*ierr) return;
28   *ierr = PetscStrncpy((char *)name, names[*i - 1], l_b);
29   if (*ierr) return;
30   FIXRETURNCHAR(PETSC_TRUE, name, l_b);
31   for (PetscInt j = 0; j < n; j++) *ierr = PetscFree(names[j]);
32   *ierr = PetscFree(names);
33 }
34 
dmcreatefielddecompositiongetisdm_(DM * dm,F90Array1d * iss,F90Array1d * dms,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd1)PETSC_F90_2PTR_PROTO (ptrd2))35 PETSC_EXTERN void dmcreatefielddecompositiongetisdm_(DM *dm, F90Array1d *iss, F90Array1d *dms, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
36 {
37   PetscInt n;
38   IS      *tis;
39   DM      *tdm;
40 
41   if (iss && dms) {
42     *ierr = DMCreateFieldDecomposition(*dm, &n, NULL, &tis, &tdm);
43   } else if (iss) {
44     *ierr = DMCreateFieldDecomposition(*dm, &n, NULL, &tis, NULL);
45   } else if (dms) {
46     *ierr = DMCreateFieldDecomposition(*dm, &n, NULL, NULL, &tdm);
47   }
48   if (*ierr) return;
49   if (iss) *ierr = F90Array1dCreate(tis, MPIU_FORTRANADDR, 1, n, iss PETSC_F90_2PTR_PARAM(ptrd1));
50   if (*ierr) return;
51   if (dms) *ierr = F90Array1dCreate(tdm, MPIU_FORTRANADDR, 1, n, dms PETSC_F90_2PTR_PARAM(ptrd2));
52 }
53 
dmcreatefielddecompositionrestoreisdm_(DM * dm,F90Array1d * iss,F90Array1d * dms,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (ptrd1)PETSC_F90_2PTR_PROTO (ptrd2))54 PETSC_EXTERN void dmcreatefielddecompositionrestoreisdm_(DM *dm, F90Array1d *iss, F90Array1d *dms, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
55 {
56   PetscInt n;
57 
58   *ierr = DMGetNumFields(*dm, &n);
59   if (*ierr) return;
60   if (iss) {
61     IS *tis;
62     *ierr = F90Array1dAccess(iss, MPIU_FORTRANADDR, (void **)&tis PETSC_F90_2PTR_PARAM(ptrd1));
63     if (*ierr) return;
64     *ierr = F90Array1dDestroy(iss, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd1));
65     if (*ierr) return;
66     for (PetscInt i = 0; i < n; i++) *ierr = ISDestroy(&tis[i]);
67     *ierr = PetscFree(tis);
68     if (*ierr) return;
69   }
70   if (dms) {
71     DM *tdm;
72     *ierr = F90Array1dAccess(dms, MPIU_FORTRANADDR, (void **)&tdm PETSC_F90_2PTR_PARAM(ptrd2));
73     if (*ierr) return;
74     *ierr = F90Array1dDestroy(dms, MPIU_FORTRANADDR PETSC_F90_2PTR_PARAM(ptrd2));
75     if (*ierr) return;
76     for (PetscInt i = 0; i < n; i++) *ierr = DMDestroy(&tdm[i]);
77     *ierr = PetscFree(tdm);
78     if (*ierr) return;
79   }
80 }
81