Actual source code: zplexfemf90.c

  1: #include <petsc/private/fortranimpl.h>
  2: #include <petscdmplex.h>
  3: #include <petscds.h>
  4: #include <petsc/private/f90impl.h>

  6: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  7:   #define dmplexgetcellfields_       DMPLEXGETCELLFIELDS
  8:   #define dmplexrestorecellfields_   DMPLEXRESTORECELLFIELDS
  9:   #define dmplexgetfacefields_       DMPLEXGETFACEFIELDS
 10:   #define dmplexrestorefacefields_   DMPLEXRESTOREFACEFIELDS
 11:   #define dmplexgetfacegeometry_     DMPLEXGETFACEGEOMETRY
 12:   #define dmplexrestorefacegeometry_ DMPLEXRESTOREFACEGEOMETRY
 13: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) && !defined(FORTRANDOUBLEUNDERSCORE)
 14:   #define dmplexgetcellfields_       dmplexgetcellfields
 15:   #define dmplexrestorecellfields_   dmplexrestorecellfields
 16:   #define dmplexgetfacefields_       dmplexgetfacefields
 17:   #define dmplexrestorefacefields_   dmplexrestorefacefields
 18:   #define dmplexgetfacegeometry_     dmplexgetfacegeometry
 19:   #define dmplexrestorefacegeometry_ dmplexrestorefacegeometry
 20: #endif

 22: PETSC_EXTERN void dmplexgetcellfields_(DM *dm, IS *cellIS, Vec *locX, Vec *locX_t, Vec *locA, F90Array1d *uPtr, F90Array1d *utPtr, F90Array1d *aPtr, int *ierr PETSC_F90_2PTR_PROTO(uPtrd) PETSC_F90_2PTR_PROTO(utPtrd) PETSC_F90_2PTR_PROTO(aPtrd))
 23: {
 24:   PetscDS      prob;
 25:   PetscScalar *u, *u_t, *a;
 26:   PetscInt     numCells, totDim, totDimAux = 0;

 28:   *ierr = ISGetLocalSize(*cellIS, &numCells);
 29:   if (*ierr) return;
 30:   *ierr = DMPlexGetCellFields(*dm, *cellIS, *locX, *locX_t, *locA, &u, &u_t, &a);
 31:   if (*ierr) return;
 32:   *ierr = DMGetDS(*dm, &prob);
 33:   if (*ierr) return;
 34:   *ierr = PetscDSGetTotalDimension(prob, &totDim);
 35:   if (*ierr) return;
 36:   if (locA) {
 37:     DM      dmAux;
 38:     PetscDS probAux;

 40:     *ierr = VecGetDM(*locA, &dmAux);
 41:     if (*ierr) return;
 42:     *ierr = DMGetDS(dmAux, &probAux);
 43:     if (*ierr) return;
 44:     *ierr = PetscDSGetTotalDimension(probAux, &totDimAux);
 45:     if (*ierr) return;
 46:   }
 47:   *ierr = F90Array1dCreate((void *)u, MPIU_SCALAR, 1, numCells * totDim, uPtr PETSC_F90_2PTR_PARAM(uPtrd));
 48:   if (*ierr) return;
 49:   *ierr = F90Array1dCreate((void *)u_t, MPIU_SCALAR, 1, locX_t ? numCells * totDim : 0, utPtr PETSC_F90_2PTR_PARAM(utPtrd));
 50:   if (*ierr) return;
 51:   *ierr = F90Array1dCreate((void *)a, MPIU_SCALAR, 1, locA ? numCells * totDimAux : 0, aPtr PETSC_F90_2PTR_PARAM(aPtrd));
 52: }

 54: PETSC_EXTERN void dmplexrestorecellfields_(DM *dm, IS *cellIS, Vec *locX, Vec *locX_t, Vec *locA, F90Array1d *uPtr, F90Array1d *utPtr, F90Array1d *aPtr, int *ierr PETSC_F90_2PTR_PROTO(uPtrd) PETSC_F90_2PTR_PROTO(utPtrd) PETSC_F90_2PTR_PROTO(aPtrd))
 55: {
 56:   PetscScalar *u, *u_t, *a;

 58:   *ierr = F90Array1dAccess(uPtr, MPIU_SCALAR, (void **)&u PETSC_F90_2PTR_PARAM(uPtrd));
 59:   if (*ierr) return;
 60:   *ierr = F90Array1dAccess(utPtr, MPIU_SCALAR, (void **)&u_t PETSC_F90_2PTR_PARAM(utPtrd));
 61:   if (*ierr) return;
 62:   *ierr = F90Array1dAccess(aPtr, MPIU_SCALAR, (void **)&a PETSC_F90_2PTR_PARAM(aPtrd));
 63:   if (*ierr) return;
 64:   *ierr = DMPlexRestoreCellFields(*dm, *cellIS, *locX, NULL, NULL, &u, u_t ? &u_t : NULL, a ? &a : NULL);
 65:   if (*ierr) return;
 66:   *ierr = F90Array1dDestroy(uPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(uPtrd));
 67:   if (*ierr) return;
 68:   *ierr = F90Array1dDestroy(utPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(utPtrd));
 69:   if (*ierr) return;
 70:   *ierr = F90Array1dDestroy(aPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(aPtrd));
 71:   if (*ierr) return;
 72: }

 74: PETSC_EXTERN void dmplexgetfacefields_(DM *dm, PetscInt *fStart, PetscInt *fEnd, Vec *locX, Vec *locX_t, Vec *faceGeometry, Vec *cellGeometry, Vec *locGrad, PetscInt *Nface, F90Array1d *uLPtr, F90Array1d *uRPtr, int *ierr PETSC_F90_2PTR_PROTO(uLPtrd) PETSC_F90_2PTR_PROTO(uRPtrd))
 75: {
 76:   PetscDS      prob;
 77:   PetscScalar *uL, *uR;
 78:   PetscInt     numFaces = *fEnd - *fStart, totDim;

 80:   *ierr = DMPlexGetFaceFields(*dm, *fStart, *fEnd, *locX, *locX_t, *faceGeometry, *cellGeometry, *locGrad, Nface, &uL, &uR);
 81:   if (*ierr) return;
 82:   *ierr = DMGetDS(*dm, &prob);
 83:   if (*ierr) return;
 84:   *ierr = PetscDSGetTotalDimension(prob, &totDim);
 85:   if (*ierr) return;
 86:   *ierr = F90Array1dCreate((void *)uL, MPIU_SCALAR, 1, numFaces * totDim, uLPtr PETSC_F90_2PTR_PARAM(uLPtrd));
 87:   if (*ierr) return;
 88:   *ierr = F90Array1dCreate((void *)uR, MPIU_SCALAR, 1, numFaces * totDim, uRPtr PETSC_F90_2PTR_PARAM(uRPtrd));
 89:   if (*ierr) return;
 90: }

 92: PETSC_EXTERN void dmplexrestorefacefields_(DM *dm, PetscInt *fStart, PetscInt *fEnd, Vec *locX, Vec *locX_t, Vec *faceGeometry, Vec *cellGeometry, Vec *locGrad, PetscInt *Nface, F90Array1d *uLPtr, F90Array1d *uRPtr, int *ierr PETSC_F90_2PTR_PROTO(uLPtrd) PETSC_F90_2PTR_PROTO(uRPtrd))
 93: {
 94:   PetscScalar *uL, *uR;

 96:   *ierr = F90Array1dAccess(uLPtr, MPIU_SCALAR, (void **)&uL PETSC_F90_2PTR_PARAM(uLPtrd));
 97:   if (*ierr) return;
 98:   *ierr = F90Array1dAccess(uRPtr, MPIU_SCALAR, (void **)&uR PETSC_F90_2PTR_PARAM(uRPtrd));
 99:   if (*ierr) return;
100:   *ierr = DMPlexRestoreFaceFields(*dm, *fStart, *fEnd, *locX, NULL, *faceGeometry, *cellGeometry, NULL, Nface, &uL, &uR);
101:   if (*ierr) return;
102:   *ierr = F90Array1dDestroy(uLPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(uLPtrd));
103:   if (*ierr) return;
104:   *ierr = F90Array1dDestroy(uRPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(uRPtrd));
105:   if (*ierr) return;
106: }

108: PETSC_EXTERN void dmplexgetfacegeometry_(DM *dm, PetscInt *fStart, PetscInt *fEnd, Vec *faceGeometry, Vec *cellGeometry, PetscInt *Nface, F90Array1d *gPtr, F90Array1d *vPtr, int *ierr PETSC_F90_2PTR_PROTO(gPtrd) PETSC_F90_2PTR_PROTO(vPtrd))
109: {
110:   PetscFVFaceGeom *g;
111:   PetscReal       *v;
112:   PetscInt         numFaces = *fEnd - *fStart, structSize = sizeof(PetscFVFaceGeom) / sizeof(PetscScalar);

114:   *ierr = DMPlexGetFaceGeometry(*dm, *fStart, *fEnd, *faceGeometry, *cellGeometry, Nface, &g, &v);
115:   if (*ierr) return;
116:   *ierr = F90Array1dCreate((void *)g, MPIU_SCALAR, 1, numFaces * structSize, gPtr PETSC_F90_2PTR_PARAM(gPtrd));
117:   if (*ierr) return;
118:   *ierr = F90Array1dCreate((void *)v, MPIU_REAL, 1, numFaces * 2, vPtr PETSC_F90_2PTR_PARAM(vPtrd));
119:   if (*ierr) return;
120: }

122: PETSC_EXTERN void dmplexrestorefacegeometry_(DM *dm, PetscInt *fStart, PetscInt *fEnd, Vec *faceGeometry, Vec *cellGeometry, PetscInt *Nface, F90Array1d *gPtr, F90Array1d *vPtr, int *ierr PETSC_F90_2PTR_PROTO(gPtrd) PETSC_F90_2PTR_PROTO(vPtrd))
123: {
124:   PetscFVFaceGeom *g;
125:   PetscReal       *v;

127:   *ierr = F90Array1dAccess(gPtr, MPIU_SCALAR, (void **)&g PETSC_F90_2PTR_PARAM(gPtrd));
128:   if (*ierr) return;
129:   *ierr = F90Array1dAccess(vPtr, MPIU_REAL, (void **)&v PETSC_F90_2PTR_PARAM(vPtrd));
130:   if (*ierr) return;
131:   *ierr = DMPlexRestoreFaceGeometry(*dm, *fStart, *fEnd, *faceGeometry, *cellGeometry, Nface, &g, &v);
132:   if (*ierr) return;
133:   *ierr = F90Array1dDestroy(gPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(vPtrd));
134:   if (*ierr) return;
135:   *ierr = F90Array1dDestroy(vPtr, MPIU_REAL PETSC_F90_2PTR_PARAM(gPtrd));
136:   if (*ierr) return;
137: }