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