xref: /petsc/src/dm/field/impls/da/dmfieldda.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
251a454edSToby Isaac #include <petsc/private/dmimpl.h>      /*I "petscdm.h" I*/
351a454edSToby Isaac #include <petscdmda.h>
451a454edSToby Isaac 
59371c9d4SSatish Balay typedef struct _n_DMField_DA {
651a454edSToby Isaac   PetscScalar *cornerVals;
751a454edSToby Isaac   PetscScalar *cornerCoeffs;
851a454edSToby Isaac   PetscScalar *work;
951a454edSToby Isaac   PetscReal    coordRange[3][2];
109371c9d4SSatish Balay } DMField_DA;
1151a454edSToby Isaac 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldDestroy_DA(DMField field)
13d71ae5a4SJacob Faibussowitsch {
1451a454edSToby Isaac   DMField_DA *dafield;
1551a454edSToby Isaac 
1651a454edSToby Isaac   PetscFunctionBegin;
1751a454edSToby Isaac   dafield = (DMField_DA *)field->data;
189566063dSJacob Faibussowitsch   PetscCall(PetscFree3(dafield->cornerVals, dafield->cornerCoeffs, dafield->work));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree(dafield));
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2151a454edSToby Isaac }
2251a454edSToby Isaac 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldView_DA(DMField field, PetscViewer viewer)
24d71ae5a4SJacob Faibussowitsch {
2551a454edSToby Isaac   DMField_DA *dafield = (DMField_DA *)field->data;
2651a454edSToby Isaac   PetscBool   iascii;
2751a454edSToby Isaac 
2851a454edSToby Isaac   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3051a454edSToby Isaac   if (iascii) {
3151a454edSToby Isaac     PetscInt i, c, dim;
3251a454edSToby Isaac     PetscInt nc;
3351a454edSToby Isaac     DM       dm = field->dm;
3451a454edSToby Isaac 
359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Field corner values:\n"));
369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
379566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
3851a454edSToby Isaac     nc = field->numComponents;
3951a454edSToby Isaac     for (i = 0, c = 0; i < (1 << dim); i++) {
4051a454edSToby Isaac       PetscInt j;
4151a454edSToby Isaac 
4251a454edSToby Isaac       for (j = 0; j < nc; j++, c++) {
4351a454edSToby Isaac         PetscScalar val = dafield->cornerVals[nc * i + j];
4451a454edSToby Isaac 
4551a454edSToby Isaac #if !defined(PETSC_USE_COMPLEX)
469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)val));
4751a454edSToby Isaac #else
489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%g+i%g ", (double)PetscRealPart(val), (double)PetscImaginaryPart(val)));
4951a454edSToby Isaac #endif
5051a454edSToby Isaac       }
519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
5251a454edSToby Isaac     }
539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
5451a454edSToby Isaac   }
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5651a454edSToby Isaac }
5751a454edSToby Isaac 
5851a454edSToby Isaac #define MEdot(y, A, x, m, c, cast) \
5951a454edSToby Isaac   do { \
6051a454edSToby Isaac     PetscInt _k, _l; \
6151a454edSToby Isaac     for (_k = 0; _k < (c); _k++) (y)[_k] = 0.; \
6251a454edSToby Isaac     for (_l = 0; _l < (m); _l++) { \
63ad540459SPierre Jolivet       for (_k = 0; _k < (c); _k++) (y)[_k] += cast((A)[(c) * _l + _k]) * (x)[_l]; \
6451a454edSToby Isaac     } \
6551a454edSToby Isaac   } while (0)
6651a454edSToby Isaac 
6751a454edSToby Isaac #define MEHess(out, cf, etaB, etaD, dim, nc, cast) \
6851a454edSToby Isaac   do { \
6951a454edSToby Isaac     PetscInt _m, _j, _k; \
7051a454edSToby Isaac     for (_m = 0; _m < (nc) * (dim) * (dim); _m++) (out)[_m] = 0.; \
7151a454edSToby Isaac     for (_j = 0; _j < (dim); _j++) { \
7251a454edSToby Isaac       for (_k = _j + 1; _k < (dim); _k++) { \
7351a454edSToby Isaac         PetscInt _ind = (1 << _j) + (1 << _k); \
7451a454edSToby Isaac         for (_m = 0; _m < (nc); _m++) { \
7551a454edSToby Isaac           PetscScalar c = (cf)[_m] * (etaB)[_ind] * (etaD)[_ind]; \
7651a454edSToby Isaac           (out)[(_m * (dim) + _k) * (dim) + _j] += cast(c); \
7751a454edSToby Isaac           (out)[(_m * (dim) + _j) * (dim) + _k] += cast(c); \
7851a454edSToby Isaac         } \
7951a454edSToby Isaac       } \
8051a454edSToby Isaac     } \
8151a454edSToby Isaac   } while (0)
8251a454edSToby Isaac 
83d71ae5a4SJacob Faibussowitsch static void MultilinearEvaluate(PetscInt dim, PetscReal (*coordRange)[2], PetscInt nc, PetscScalar *cf, PetscScalar *cfWork, PetscInt nPoints, const PetscScalar *points, PetscDataType datatype, void *B, void *D, void *H)
84d71ae5a4SJacob Faibussowitsch {
8551a454edSToby Isaac   PetscInt i, j, k, l, m;
8651a454edSToby Isaac   PetscInt whol = 1 << dim;
8751a454edSToby Isaac   PetscInt half = whol >> 1;
8851a454edSToby Isaac 
8951a454edSToby Isaac   PetscFunctionBeginHot;
9051a454edSToby Isaac   if (!B && !D && !H) PetscFunctionReturnVoid();
9151a454edSToby Isaac   for (i = 0; i < nPoints; i++) {
9251a454edSToby Isaac     const PetscScalar *point   = &points[dim * i];
9351a454edSToby Isaac     PetscReal          deta[3] = {0.};
9451a454edSToby Isaac     PetscReal          etaB[8] = {1., 1., 1., 1., 1., 1., 1., 1.};
9551a454edSToby Isaac     PetscReal          etaD[8] = {1., 1., 1., 1., 1., 1., 1., 1.};
9651a454edSToby Isaac     PetscReal          work[8];
9751a454edSToby Isaac 
9851a454edSToby Isaac     for (j = 0; j < dim; j++) {
9951a454edSToby Isaac       PetscReal e, d;
10051a454edSToby Isaac 
10151a454edSToby Isaac       e       = (PetscRealPart(point[j]) - coordRange[j][0]) / coordRange[j][1];
10251a454edSToby Isaac       deta[j] = d = 1. / coordRange[j][1];
103ad540459SPierre Jolivet       for (k = 0; k < whol; k++) work[k] = etaB[k];
10451a454edSToby Isaac       for (k = 0; k < half; k++) {
10551a454edSToby Isaac         etaB[k]        = work[2 * k] * e;
10651a454edSToby Isaac         etaB[k + half] = work[2 * k + 1];
10751a454edSToby Isaac       }
10851a454edSToby Isaac       if (H) {
109ad540459SPierre Jolivet         for (k = 0; k < whol; k++) work[k] = etaD[k];
11051a454edSToby Isaac         for (k = 0; k < half; k++) {
11151a454edSToby Isaac           etaD[k + half] = work[2 * k];
11251a454edSToby Isaac           etaD[k]        = work[2 * k + 1] * d;
11351a454edSToby Isaac         }
11451a454edSToby Isaac       }
11551a454edSToby Isaac     }
11651a454edSToby Isaac     if (B) {
11751a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
11851a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)B)[nc * i];
11951a454edSToby Isaac 
12051a454edSToby Isaac         MEdot(out, cf, etaB, (1 << dim), nc, (PetscScalar));
12151a454edSToby Isaac       } else {
12251a454edSToby Isaac         PetscReal *out = &((PetscReal *)B)[nc * i];
12351a454edSToby Isaac 
12451a454edSToby Isaac         MEdot(out, cf, etaB, (1 << dim), nc, PetscRealPart);
12551a454edSToby Isaac       }
12651a454edSToby Isaac     }
12751a454edSToby Isaac     if (D) {
12851a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
12951a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
13051a454edSToby Isaac 
13151a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
13251a454edSToby Isaac       } else {
13351a454edSToby Isaac         PetscReal *out = &((PetscReal *)D)[nc * dim * i];
13451a454edSToby Isaac 
13551a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
13651a454edSToby Isaac       }
13751a454edSToby Isaac       for (j = 0; j < dim; j++) {
13851a454edSToby Isaac         PetscReal d = deta[j];
13951a454edSToby Isaac 
140ad540459SPierre Jolivet         for (k = 0; k < whol * nc; k++) cfWork[k] = cf[k];
141ad540459SPierre Jolivet         for (k = 0; k < whol; k++) work[k] = etaB[k];
14251a454edSToby Isaac         for (k = 0; k < half; k++) {
14351a454edSToby Isaac           PetscReal e;
14451a454edSToby Isaac 
14551a454edSToby Isaac           etaB[k]        = work[2 * k];
14651a454edSToby Isaac           etaB[k + half] = e = work[2 * k + 1];
14751a454edSToby Isaac 
14851a454edSToby Isaac           for (l = 0; l < nc; l++) {
14951a454edSToby Isaac             cf[k * nc + l]          = cfWork[2 * k * nc + l];
15051a454edSToby Isaac             cf[(k + half) * nc + l] = cfWork[(2 * k + 1) * nc + l];
15151a454edSToby Isaac           }
15251a454edSToby Isaac           if (datatype == PETSC_SCALAR) {
15351a454edSToby Isaac             PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
15451a454edSToby Isaac 
155ad540459SPierre Jolivet             for (l = 0; l < nc; l++) out[l * dim + j] += d * e * cf[k * nc + l];
15651a454edSToby Isaac           } else {
15751a454edSToby Isaac             PetscReal *out = &((PetscReal *)D)[nc * dim * i];
15851a454edSToby Isaac 
159ad540459SPierre Jolivet             for (l = 0; l < nc; l++) out[l * dim + j] += d * e * PetscRealPart(cf[k * nc + l]);
16051a454edSToby Isaac           }
16151a454edSToby Isaac         }
16251a454edSToby Isaac       }
16351a454edSToby Isaac     }
16451a454edSToby Isaac     if (H) {
16551a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
16651a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)H)[nc * dim * dim * i];
16751a454edSToby Isaac 
16851a454edSToby Isaac         MEHess(out, cf, etaB, etaD, dim, nc, (PetscScalar));
16951a454edSToby Isaac       } else {
17051a454edSToby Isaac         PetscReal *out = &((PetscReal *)H)[nc * dim * dim * i];
17151a454edSToby Isaac 
17251a454edSToby Isaac         MEHess(out, cf, etaB, etaD, dim, nc, PetscRealPart);
17351a454edSToby Isaac       }
17451a454edSToby Isaac     }
17551a454edSToby Isaac   }
17651a454edSToby Isaac   PetscFunctionReturnVoid();
17751a454edSToby Isaac }
17851a454edSToby Isaac 
179d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluate_DA(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
180d71ae5a4SJacob Faibussowitsch {
18151a454edSToby Isaac   DM                 dm;
18251a454edSToby Isaac   DMField_DA        *dafield;
18351a454edSToby Isaac   PetscInt           dim;
18451a454edSToby Isaac   PetscInt           N, n, nc;
18551a454edSToby Isaac   const PetscScalar *array;
18651a454edSToby Isaac   PetscReal(*coordRange)[2];
18751a454edSToby Isaac 
18851a454edSToby Isaac   PetscFunctionBegin;
18951a454edSToby Isaac   dm      = field->dm;
19051a454edSToby Isaac   nc      = field->numComponents;
19151a454edSToby Isaac   dafield = (DMField_DA *)field->data;
1929566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1939566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(points, &N));
19463a3b9bcSJacob Faibussowitsch   PetscCheck(N % dim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point vector size %" PetscInt_FMT " not divisible by coordinate dimension %" PetscInt_FMT, N, dim);
19551a454edSToby Isaac   n          = N / dim;
196f4f49eeaSPierre Jolivet   coordRange = &dafield->coordRange[0];
1979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points, &array));
19851a454edSToby Isaac   MultilinearEvaluate(dim, coordRange, nc, dafield->cornerCoeffs, dafield->work, n, array, datatype, B, D, H);
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points, &array));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20151a454edSToby Isaac }
20251a454edSToby Isaac 
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFE_DA(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
204d71ae5a4SJacob Faibussowitsch {
20551a454edSToby Isaac   PetscInt  c, i, j, k, dim, cellsPer[3] = {0}, first[3] = {0}, whol, half;
20651a454edSToby Isaac   PetscReal stepPer[3]           = {0.};
2079371c9d4SSatish Balay   PetscReal cellCoordRange[3][2] = {
2089371c9d4SSatish Balay     {0., 1.},
2099371c9d4SSatish Balay     {0., 1.},
2109371c9d4SSatish Balay     {0., 1.}
2119371c9d4SSatish Balay   };
21251a454edSToby Isaac   PetscScalar     *cellCoeffs, *work;
21351a454edSToby Isaac   DM               dm;
21451a454edSToby Isaac   DMDALocalInfo    info;
21551a454edSToby Isaac   PetscInt         cStart, cEnd;
21651a454edSToby Isaac   PetscInt         nq, nc;
21751a454edSToby Isaac   const PetscReal *q;
21851a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
21951a454edSToby Isaac   PetscScalar *qs;
22051a454edSToby Isaac #else
22151a454edSToby Isaac   const PetscScalar *qs;
22251a454edSToby Isaac #endif
22351a454edSToby Isaac   DMField_DA     *dafield;
22451a454edSToby Isaac   PetscBool       isStride;
22551a454edSToby Isaac   const PetscInt *cells  = NULL;
22651a454edSToby Isaac   PetscInt        sfirst = -1, stride = -1, nCells;
22751a454edSToby Isaac 
22851a454edSToby Isaac   PetscFunctionBegin;
22951a454edSToby Isaac   dafield = (DMField_DA *)field->data;
23051a454edSToby Isaac   dm      = field->dm;
23151a454edSToby Isaac   nc      = field->numComponents;
2329566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
23351a454edSToby Isaac   dim         = info.dim;
23451a454edSToby Isaac   work        = dafield->work;
23551a454edSToby Isaac   stepPer[0]  = 1. / info.mx;
23651a454edSToby Isaac   stepPer[1]  = 1. / info.my;
23751a454edSToby Isaac   stepPer[2]  = 1. / info.mz;
23851a454edSToby Isaac   first[0]    = info.gxs;
23951a454edSToby Isaac   first[1]    = info.gys;
24051a454edSToby Isaac   first[2]    = info.gzs;
24151a454edSToby Isaac   cellsPer[0] = info.gxm;
24251a454edSToby Isaac   cellsPer[1] = info.gym;
24351a454edSToby Isaac   cellsPer[2] = info.gzm;
24451a454edSToby Isaac   /* TODO: probably take components into account */
2459566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(points, NULL, NULL, &nq, &q, NULL));
24651a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
2479566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, nq * dim, MPIU_SCALAR, &qs));
24851a454edSToby Isaac   for (i = 0; i < nq * dim; i++) qs[i] = q[i];
24951a454edSToby Isaac #else
25051a454edSToby Isaac   qs = q;
25151a454edSToby Isaac #endif
2529566063dSJacob Faibussowitsch   PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd));
2539566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs));
25451a454edSToby Isaac   whol = (1 << dim);
25551a454edSToby Isaac   half = whol >> 1;
2569566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(cellIS, &nCells));
2579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride));
25851a454edSToby Isaac   if (isStride) {
2599566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride));
2601baa6e33SBarry Smith   } else PetscCall(ISGetIndices(cellIS, &cells));
26151a454edSToby Isaac   for (c = 0; c < nCells; c++) {
26251a454edSToby Isaac     PetscInt cell   = isStride ? (sfirst + c * stride) : cells[c];
26351a454edSToby Isaac     PetscInt rem    = cell;
26451a454edSToby Isaac     PetscInt ijk[3] = {0};
26551a454edSToby Isaac     void    *cB, *cD, *cH;
26651a454edSToby Isaac 
26751a454edSToby Isaac     if (datatype == PETSC_SCALAR) {
2688e3a54c0SPierre Jolivet       cB = PetscSafePointerPlusOffset((PetscScalar *)B, nc * nq * c);
2698e3a54c0SPierre Jolivet       cD = PetscSafePointerPlusOffset((PetscScalar *)D, nc * nq * dim * c);
2708e3a54c0SPierre Jolivet       cH = PetscSafePointerPlusOffset((PetscScalar *)H, nc * nq * dim * dim * c);
27151a454edSToby Isaac     } else {
27257508eceSPierre Jolivet       cB = PetscSafePointerPlusOffset((PetscReal *)B, nc * nq * c);
27357508eceSPierre Jolivet       cD = PetscSafePointerPlusOffset((PetscReal *)D, nc * nq * dim * c);
27457508eceSPierre Jolivet       cH = PetscSafePointerPlusOffset((PetscReal *)H, nc * nq * dim * dim * c);
27551a454edSToby Isaac     }
2761dca8a05SBarry Smith     PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " not a cell [%" PetscInt_FMT ",%" PetscInt_FMT "), not implemented yet", cell, cStart, cEnd);
277ad540459SPierre Jolivet     for (i = 0; i < nc * whol; i++) work[i] = dafield->cornerCoeffs[i];
27851a454edSToby Isaac     for (j = 0; j < dim; j++) {
27951a454edSToby Isaac       PetscReal e, d;
28051a454edSToby Isaac       ijk[j] = (rem % cellsPer[j]);
28151a454edSToby Isaac       rem /= cellsPer[j];
28251a454edSToby Isaac 
28351a454edSToby Isaac       e = 2. * (ijk[j] + first[j] + 0.5) * stepPer[j] - 1.;
28451a454edSToby Isaac       d = stepPer[j];
28551a454edSToby Isaac       for (i = 0; i < half; i++) {
28651a454edSToby Isaac         for (k = 0; k < nc; k++) {
28751a454edSToby Isaac           cellCoeffs[i * nc + k]          = work[2 * i * nc + k] * d;
28851a454edSToby Isaac           cellCoeffs[(i + half) * nc + k] = work[2 * i * nc + k] * e + work[(2 * i + 1) * nc + k];
28951a454edSToby Isaac         }
29051a454edSToby Isaac       }
291ad540459SPierre Jolivet       for (i = 0; i < whol * nc; i++) work[i] = cellCoeffs[i];
29251a454edSToby Isaac     }
29351a454edSToby Isaac     MultilinearEvaluate(dim, cellCoordRange, nc, cellCoeffs, dafield->work, nq, qs, datatype, cB, cD, cH);
29451a454edSToby Isaac   }
29548a46eb9SPierre Jolivet   if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells));
2969566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs));
29751a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
2989566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, nq * dim, MPIU_SCALAR, &qs));
29951a454edSToby Isaac #endif
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30151a454edSToby Isaac }
30251a454edSToby Isaac 
303d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFV_DA(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
304d71ae5a4SJacob Faibussowitsch {
30551a454edSToby Isaac   PetscInt        c, i, dim, cellsPer[3] = {0}, first[3] = {0};
30651a454edSToby Isaac   PetscReal       stepPer[3] = {0.};
30751a454edSToby Isaac   DM              dm;
30851a454edSToby Isaac   DMDALocalInfo   info;
30951a454edSToby Isaac   PetscInt        cStart, cEnd, numCells;
31051a454edSToby Isaac   PetscInt        nc;
3114d1a973fSToby Isaac   PetscScalar    *points;
31251a454edSToby Isaac   DMField_DA     *dafield;
31351a454edSToby Isaac   PetscBool       isStride;
31451a454edSToby Isaac   const PetscInt *cells  = NULL;
31551a454edSToby Isaac   PetscInt        sfirst = -1, stride = -1;
31651a454edSToby Isaac 
31751a454edSToby Isaac   PetscFunctionBegin;
31851a454edSToby Isaac   dafield = (DMField_DA *)field->data;
31951a454edSToby Isaac   dm      = field->dm;
32051a454edSToby Isaac   nc      = field->numComponents;
3219566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
32251a454edSToby Isaac   dim         = info.dim;
32351a454edSToby Isaac   stepPer[0]  = 1. / info.mx;
32451a454edSToby Isaac   stepPer[1]  = 1. / info.my;
32551a454edSToby Isaac   stepPer[2]  = 1. / info.mz;
32651a454edSToby Isaac   first[0]    = info.gxs;
32751a454edSToby Isaac   first[1]    = info.gys;
32851a454edSToby Isaac   first[2]    = info.gzs;
32951a454edSToby Isaac   cellsPer[0] = info.gxm;
33051a454edSToby Isaac   cellsPer[1] = info.gym;
33151a454edSToby Isaac   cellsPer[2] = info.gzm;
3329566063dSJacob Faibussowitsch   PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd));
3339566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(cellIS, &numCells));
3349566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, dim * numCells, MPIU_SCALAR, &points));
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride));
33651a454edSToby Isaac   if (isStride) {
3379566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride));
3381baa6e33SBarry Smith   } else PetscCall(ISGetIndices(cellIS, &cells));
33951a454edSToby Isaac   for (c = 0; c < numCells; c++) {
34051a454edSToby Isaac     PetscInt cell   = isStride ? (sfirst + c * stride) : cells[c];
34151a454edSToby Isaac     PetscInt rem    = cell;
34251a454edSToby Isaac     PetscInt ijk[3] = {0};
34351a454edSToby Isaac 
3441dca8a05SBarry Smith     PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " not a cell [%" PetscInt_FMT ",%" PetscInt_FMT "), not implemented yet", cell, cStart, cEnd);
34551a454edSToby Isaac     for (i = 0; i < dim; i++) {
34651a454edSToby Isaac       ijk[i] = (rem % cellsPer[i]);
34751a454edSToby Isaac       rem /= cellsPer[i];
34851a454edSToby Isaac       points[dim * c + i] = (ijk[i] + first[i] + 0.5) * stepPer[i];
34951a454edSToby Isaac     }
35051a454edSToby Isaac   }
35148a46eb9SPierre Jolivet   if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells));
35251a454edSToby Isaac   MultilinearEvaluate(dim, dafield->coordRange, nc, dafield->cornerCoeffs, dafield->work, numCells, points, datatype, B, D, H);
3539566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, dim * numCells, MPIU_SCALAR, &points));
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35551a454edSToby Isaac }
35651a454edSToby Isaac 
357d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldGetDegree_DA(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
358d71ae5a4SJacob Faibussowitsch {
35951a454edSToby Isaac   DM       dm;
36051a454edSToby Isaac   PetscInt dim, h, imin;
36151a454edSToby Isaac 
36251a454edSToby Isaac   PetscFunctionBegin;
36351a454edSToby Isaac   dm = field->dm;
3649566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, NULL));
3659566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
36651a454edSToby Isaac   for (h = 0; h <= dim; h++) {
36751a454edSToby Isaac     PetscInt hEnd;
36851a454edSToby Isaac 
3699566063dSJacob Faibussowitsch     PetscCall(DMDAGetHeightStratum(dm, h, NULL, &hEnd));
37051a454edSToby Isaac     if (imin < hEnd) break;
37151a454edSToby Isaac   }
37251a454edSToby Isaac   dim -= h;
373b7260050SToby Isaac   if (minDegree) *minDegree = 1;
374b7260050SToby Isaac   if (maxDegree) *maxDegree = dim;
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37651a454edSToby Isaac }
37751a454edSToby Isaac 
378d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldCreateDefaultQuadrature_DA(DMField field, IS cellIS, PetscQuadrature *quad)
379d71ae5a4SJacob Faibussowitsch {
38051a454edSToby Isaac   PetscInt h, dim, imax, imin;
38151a454edSToby Isaac   DM       dm;
38251a454edSToby Isaac 
38351a454edSToby Isaac   PetscFunctionBegin;
38451a454edSToby Isaac   dm = field->dm;
3859566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(cellIS, &imin, &imax));
3869566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
38751a454edSToby Isaac   *quad = NULL;
38851a454edSToby Isaac   for (h = 0; h <= dim; h++) {
38951a454edSToby Isaac     PetscInt hStart, hEnd;
39051a454edSToby Isaac 
3919566063dSJacob Faibussowitsch     PetscCall(DMDAGetHeightStratum(dm, h, &hStart, &hEnd));
39251a454edSToby Isaac     if (imin >= hStart && imax < hEnd) break;
39351a454edSToby Isaac   }
39451a454edSToby Isaac   dim -= h;
39548a46eb9SPierre Jolivet   if (dim > 0) PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39751a454edSToby Isaac }
39851a454edSToby Isaac 
399d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldInitialize_DA(DMField field)
400d71ae5a4SJacob Faibussowitsch {
40151a454edSToby Isaac   DM          dm;
40251a454edSToby Isaac   Vec         coords = NULL;
40351a454edSToby Isaac   PetscInt    dim, i, j, k;
4044d1a973fSToby Isaac   DMField_DA *dafield = (DMField_DA *)field->data;
40551a454edSToby Isaac 
40651a454edSToby Isaac   PetscFunctionBegin;
40751a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DA;
40851a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DA;
40951a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DA;
41051a454edSToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DA;
411b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DA;
41251a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DA;
41351a454edSToby Isaac   field->ops->view                    = DMFieldView_DA;
41451a454edSToby Isaac   dm                                  = field->dm;
4159566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4166858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
41751a454edSToby Isaac   if (coords) {
41851a454edSToby Isaac     PetscInt           n;
41951a454edSToby Isaac     const PetscScalar *array;
4209371c9d4SSatish Balay     PetscReal          mins[3][2] = {
4219371c9d4SSatish Balay       {PETSC_MAX_REAL, PETSC_MAX_REAL},
4229371c9d4SSatish Balay       {PETSC_MAX_REAL, PETSC_MAX_REAL},
4239371c9d4SSatish Balay       {PETSC_MAX_REAL, PETSC_MAX_REAL}
4249371c9d4SSatish Balay     };
42551a454edSToby Isaac 
4269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coords, &n));
42751a454edSToby Isaac     n /= dim;
4289566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coords, &array));
42951a454edSToby Isaac     for (i = 0, k = 0; i < n; i++) {
43051a454edSToby Isaac       for (j = 0; j < dim; j++, k++) {
43151a454edSToby Isaac         PetscReal val = PetscRealPart(array[k]);
43251a454edSToby Isaac 
43351a454edSToby Isaac         mins[j][0] = PetscMin(mins[j][0], val);
43451a454edSToby Isaac         mins[j][1] = PetscMin(mins[j][1], -val);
43551a454edSToby Isaac       }
43651a454edSToby Isaac     }
4379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coords, &array));
438*e91c04dfSPierre Jolivet     PetscCallMPI(MPIU_Allreduce(mins, &dafield->coordRange[0][0], 2 * dim, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm)));
439ad540459SPierre Jolivet     for (j = 0; j < dim; j++) dafield->coordRange[j][1] = -dafield->coordRange[j][1];
44051a454edSToby Isaac   } else {
44151a454edSToby Isaac     for (j = 0; j < dim; j++) {
44251a454edSToby Isaac       dafield->coordRange[j][0] = 0.;
44351a454edSToby Isaac       dafield->coordRange[j][1] = 1.;
44451a454edSToby Isaac     }
44551a454edSToby Isaac   }
44651a454edSToby Isaac   for (j = 0; j < dim; j++) {
4474d1a973fSToby Isaac     PetscReal avg = 0.5 * (dafield->coordRange[j][1] + dafield->coordRange[j][0]);
4484d1a973fSToby Isaac     PetscReal dif = 0.5 * (dafield->coordRange[j][1] - dafield->coordRange[j][0]);
44951a454edSToby Isaac 
45051a454edSToby Isaac     dafield->coordRange[j][0] = avg;
45151a454edSToby Isaac     dafield->coordRange[j][1] = dif;
45251a454edSToby Isaac   }
4533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45451a454edSToby Isaac }
45551a454edSToby Isaac 
456d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate_DA(DMField field)
457d71ae5a4SJacob Faibussowitsch {
45851a454edSToby Isaac   DMField_DA *dafield;
45951a454edSToby Isaac 
46051a454edSToby Isaac   PetscFunctionBegin;
4614dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dafield));
46251a454edSToby Isaac   field->data = dafield;
4639566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_DA(field));
4643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46551a454edSToby Isaac }
46651a454edSToby Isaac 
467d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDA(DM dm, PetscInt nc, const PetscScalar *cornerValues, DMField *field)
468d71ae5a4SJacob Faibussowitsch {
46951a454edSToby Isaac   DMField      b;
47051a454edSToby Isaac   DMField_DA  *dafield;
47151a454edSToby Isaac   PetscInt     dim, nv, i, j, k;
47251a454edSToby Isaac   PetscInt     half;
47351a454edSToby Isaac   PetscScalar *cv, *cf, *work;
47451a454edSToby Isaac 
47551a454edSToby Isaac   PetscFunctionBegin;
4769566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm, nc, DMFIELD_VERTEX, &b));
4779566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b, DMFIELDDA));
47851a454edSToby Isaac   dafield = (DMField_DA *)b->data;
4799566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
48051a454edSToby Isaac   nv = (1 << dim) * nc;
4819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nv, &cv, nv, &cf, nv, &work));
48251a454edSToby Isaac   for (i = 0; i < nv; i++) cv[i] = cornerValues[i];
48351a454edSToby Isaac   for (i = 0; i < nv; i++) cf[i] = cv[i];
48451a454edSToby Isaac   dafield->cornerVals   = cv;
48551a454edSToby Isaac   dafield->cornerCoeffs = cf;
48651a454edSToby Isaac   dafield->work         = work;
48751a454edSToby Isaac   half                  = (1 << (dim - 1));
48851a454edSToby Isaac   for (i = 0; i < dim; i++) {
48951a454edSToby Isaac     PetscScalar *w;
49051a454edSToby Isaac 
49151a454edSToby Isaac     w = work;
49251a454edSToby Isaac     for (j = 0; j < half; j++) {
493ad540459SPierre Jolivet       for (k = 0; k < nc; k++) w[j * nc + k] = 0.5 * (cf[(2 * j + 1) * nc + k] - cf[(2 * j) * nc + k]);
49451a454edSToby Isaac     }
49551a454edSToby Isaac     w = &work[j * nc];
49651a454edSToby Isaac     for (j = 0; j < half; j++) {
497ad540459SPierre Jolivet       for (k = 0; k < nc; k++) w[j * nc + k] = 0.5 * (cf[(2 * j) * nc + k] + cf[(2 * j + 1) * nc + k]);
49851a454edSToby Isaac     }
499ad540459SPierre Jolivet     for (j = 0; j < nv; j++) cf[j] = work[j];
50051a454edSToby Isaac   }
50151a454edSToby Isaac   *field = b;
5023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50351a454edSToby Isaac }
504