Actual source code: zfddaf90.c
1: #include <petsc/private/f90impl.h>
2: #include <petscdmcomposite.h>
4: #if defined(PETSC_HAVE_FORTRAN_CAPS)
5: #define dmcompositegetaccessvpvp_ DMCOMPOSITEGETACCESSVPVP
6: #define dmcompositerestoreaccessvpvp_ DMCOMPOSITERESTOREACCESSVPVP
7: #define dmcompositegetentriesarray_ DMCOMPOSITEGETENTRIESARRAY
8: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9: #define dmcompositegetaccessvpvp_ dmcompositegetaccessvpvp
10: #define dmcompositerestoreaccessvpvp_ dmcompositerestoreaccessvpvp
11: #define dmcompositegetentriesarray_ dmcompositegetentriesarray
12: #endif
14: PETSC_EXTERN void dmcompositegetaccessvpvp_(DM *dm, Vec *v, Vec *v1, F90Array1d *p1, Vec *v2, F90Array1d *p2, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
15: {
16: PetscScalar *pp1, *pp2;
17: PetscInt np1, np2;
18: *ierr = DMCompositeGetEntries(*dm, 0, &np1, 0, &np2);
19: *ierr = DMCompositeGetAccess(*dm, *v, v1, &pp1, v2, &pp2);
20: *ierr = F90Array1dCreate(pp1, MPIU_SCALAR, 0, np1 - 1, p1 PETSC_F90_2PTR_PARAM(ptrd1));
21: *ierr = F90Array1dCreate(pp2, MPIU_SCALAR, 0, np2 - 1, p2 PETSC_F90_2PTR_PARAM(ptrd2));
22: }
24: PETSC_EXTERN void dmcompositerestoreaccessvpvp_(DM *dm, Vec *v, Vec *v1, F90Array1d *p1, Vec *v2, F90Array1d *p2, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
25: {
26: *ierr = DMCompositeRestoreAccess(*dm, *v, v1, 0, v2, 0);
27: *ierr = F90Array1dDestroy(p1, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd1));
28: *ierr = F90Array1dDestroy(p2, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd2));
29: }
31: PETSC_EXTERN void dmcompositegetentriesarray_(DM *dm, DM *dmarray, PetscErrorCode *ierr)
32: {
33: *ierr = DMCompositeGetEntriesArray(*dm, dmarray);
34: }