xref: /petsc/src/dm/field/impls/da/dmfieldda.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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 
551a454edSToby Isaac typedef struct _n_DMField_DA
651a454edSToby Isaac {
751a454edSToby Isaac   PetscScalar     *cornerVals;
851a454edSToby Isaac   PetscScalar     *cornerCoeffs;
951a454edSToby Isaac   PetscScalar     *work;
1051a454edSToby Isaac   PetscReal       coordRange[3][2];
1151a454edSToby Isaac }
1251a454edSToby Isaac DMField_DA;
1351a454edSToby Isaac 
1451a454edSToby Isaac static PetscErrorCode DMFieldDestroy_DA(DMField field)
1551a454edSToby Isaac {
1651a454edSToby Isaac   DMField_DA     *dafield;
1751a454edSToby Isaac 
1851a454edSToby Isaac   PetscFunctionBegin;
1951a454edSToby Isaac   dafield = (DMField_DA *) field->data;
209566063dSJacob Faibussowitsch   PetscCall(PetscFree3(dafield->cornerVals,dafield->cornerCoeffs,dafield->work));
219566063dSJacob Faibussowitsch   PetscCall(PetscFree(dafield));
2251a454edSToby Isaac   PetscFunctionReturn(0);
2351a454edSToby Isaac }
2451a454edSToby Isaac 
2551a454edSToby Isaac static PetscErrorCode DMFieldView_DA(DMField field,PetscViewer viewer)
2651a454edSToby Isaac {
2751a454edSToby Isaac   DMField_DA     *dafield = (DMField_DA *) field->data;
2851a454edSToby Isaac   PetscBool      iascii;
2951a454edSToby Isaac 
3051a454edSToby Isaac   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
3251a454edSToby Isaac   if (iascii) {
3351a454edSToby Isaac     PetscInt i, c, dim;
3451a454edSToby Isaac     PetscInt nc;
3551a454edSToby Isaac     DM       dm = field->dm;
3651a454edSToby Isaac 
379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Field corner values:\n"));
389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
399566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm,&dim));
4051a454edSToby Isaac     nc = field->numComponents;
4151a454edSToby Isaac     for (i = 0, c = 0; i < (1 << dim); i++) {
4251a454edSToby Isaac       PetscInt j;
4351a454edSToby Isaac 
4451a454edSToby Isaac       for (j = 0; j < nc; j++, c++) {
4551a454edSToby Isaac         PetscScalar val = dafield->cornerVals[nc * i + j];
4651a454edSToby Isaac 
4751a454edSToby Isaac #if !defined(PETSC_USE_COMPLEX)
489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"%g ",(double) val));
4951a454edSToby Isaac #else
509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"%g+i%g ",(double) PetscRealPart(val),(double) PetscImaginaryPart(val)));
5151a454edSToby Isaac #endif
5251a454edSToby Isaac       }
539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
5451a454edSToby Isaac     }
559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
5651a454edSToby Isaac   }
5751a454edSToby Isaac   PetscFunctionReturn(0);
5851a454edSToby Isaac }
5951a454edSToby Isaac 
6051a454edSToby Isaac #define MEdot(y,A,x,m,c,cast)                          \
6151a454edSToby Isaac   do {                                                 \
6251a454edSToby Isaac     PetscInt _k, _l;                                   \
6351a454edSToby Isaac     for (_k = 0; _k < (c); _k++) (y)[_k] = 0.;         \
6451a454edSToby Isaac     for (_l = 0; _l < (m); _l++) {                     \
6551a454edSToby Isaac       for (_k = 0; _k < (c); _k++) {                   \
6651a454edSToby Isaac         (y)[_k] += cast((A)[(c) * _l + _k]) * (x)[_l]; \
6751a454edSToby Isaac       }                                                \
6851a454edSToby Isaac     }                                                  \
6951a454edSToby Isaac   } while (0)
7051a454edSToby Isaac 
7151a454edSToby Isaac #define MEHess(out,cf,etaB,etaD,dim,nc,cast)                      \
7251a454edSToby Isaac   do {                                                            \
7351a454edSToby Isaac     PetscInt _m, _j, _k;                                          \
7451a454edSToby Isaac     for (_m = 0; _m < (nc) * (dim) * (dim); _m++) (out)[_m] = 0.; \
7551a454edSToby Isaac     for (_j = 0; _j < (dim); _j++) {                              \
7651a454edSToby Isaac       for (_k = _j + 1; _k < (dim); _k++) {                       \
7751a454edSToby Isaac         PetscInt _ind = (1 << _j) + (1 << _k);                    \
7851a454edSToby Isaac         for (_m = 0; _m < (nc); _m++) {                           \
7951a454edSToby Isaac           PetscScalar c = (cf)[_m] * (etaB)[_ind] * (etaD)[_ind];   \
8051a454edSToby Isaac           (out)[(_m * (dim) + _k) * (dim) + _j] += cast(c);       \
8151a454edSToby Isaac           (out)[(_m * (dim) + _j) * (dim) + _k] += cast(c);       \
8251a454edSToby Isaac         }                                                         \
8351a454edSToby Isaac       }                                                           \
8451a454edSToby Isaac     }                                                             \
8551a454edSToby Isaac   } while (0)
8651a454edSToby Isaac 
8751a454edSToby Isaac 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)
8851a454edSToby Isaac {
8951a454edSToby Isaac   PetscInt i, j, k, l, m;
9051a454edSToby Isaac   PetscInt  whol = 1 << dim;
9151a454edSToby Isaac   PetscInt  half = whol >> 1;
9251a454edSToby Isaac 
9351a454edSToby Isaac   PetscFunctionBeginHot;
9451a454edSToby Isaac   if (!B && !D && !H) PetscFunctionReturnVoid();
9551a454edSToby Isaac   for (i = 0; i < nPoints; i++) {
9651a454edSToby Isaac     const PetscScalar *point = &points[dim * i];
9751a454edSToby Isaac     PetscReal deta[3] = {0.};
9851a454edSToby Isaac     PetscReal etaB[8] = {1.,1.,1.,1.,1.,1.,1.,1.};
9951a454edSToby Isaac     PetscReal etaD[8] = {1.,1.,1.,1.,1.,1.,1.,1.};
10051a454edSToby Isaac     PetscReal work[8];
10151a454edSToby Isaac 
10251a454edSToby Isaac     for (j = 0; j < dim; j++) {
10351a454edSToby Isaac       PetscReal e, d;
10451a454edSToby Isaac 
10551a454edSToby Isaac       e = (PetscRealPart(point[j]) - coordRange[j][0]) / coordRange[j][1];
10651a454edSToby Isaac       deta[j] = d = 1. / coordRange[j][1];
10751a454edSToby Isaac       for (k = 0; k < whol; k++) {work[k] = etaB[k];}
10851a454edSToby Isaac       for (k = 0; k < half; k++) {
10951a454edSToby Isaac         etaB[k]        = work[2 * k] * e;
11051a454edSToby Isaac         etaB[k + half] = work[2 * k + 1];
11151a454edSToby Isaac       }
11251a454edSToby Isaac       if (H) {
11351a454edSToby Isaac         for (k = 0; k < whol; k++) {work[k] = etaD[k];}
11451a454edSToby Isaac         for (k = 0; k < half; k++) {
11551a454edSToby Isaac           etaD[k + half] = work[2 * k];
11651a454edSToby Isaac           etaD[k       ] = work[2 * k + 1] * d;
11751a454edSToby Isaac         }
11851a454edSToby Isaac       }
11951a454edSToby Isaac     }
12051a454edSToby Isaac     if (B) {
12151a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
12251a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)B)[nc * i];
12351a454edSToby Isaac 
12451a454edSToby Isaac         MEdot(out,cf,etaB,(1 << dim),nc,(PetscScalar));
12551a454edSToby Isaac       } else {
12651a454edSToby Isaac         PetscReal *out = &((PetscReal *)B)[nc * i];
12751a454edSToby Isaac 
12851a454edSToby Isaac         MEdot(out,cf,etaB,(1 << dim),nc,PetscRealPart);
12951a454edSToby Isaac       }
13051a454edSToby Isaac     }
13151a454edSToby Isaac     if (D) {
13251a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
13351a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
13451a454edSToby Isaac 
13551a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
13651a454edSToby Isaac       } else {
13751a454edSToby Isaac         PetscReal *out = &((PetscReal *)D)[nc * dim * i];
13851a454edSToby Isaac 
13951a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
14051a454edSToby Isaac       }
14151a454edSToby Isaac       for (j = 0; j < dim; j++) {
14251a454edSToby Isaac         PetscReal d = deta[j];
14351a454edSToby Isaac 
14451a454edSToby Isaac         for (k = 0; k < whol * nc; k++) {cfWork[k] = cf[k];}
14551a454edSToby Isaac         for (k = 0; k < whol; k++) {work[k] = etaB[k];}
14651a454edSToby Isaac         for (k = 0; k < half; k++) {
14751a454edSToby Isaac           PetscReal e;
14851a454edSToby Isaac 
14951a454edSToby Isaac           etaB[k]        =     work[2 * k];
15051a454edSToby Isaac           etaB[k + half] = e = work[2 * k + 1];
15151a454edSToby Isaac 
15251a454edSToby Isaac           for (l = 0; l < nc; l++) {
15351a454edSToby Isaac             cf[ k         * nc + l] = cfWork[ 2 * k      * nc + l];
15451a454edSToby Isaac             cf[(k + half) * nc + l] = cfWork[(2 * k + 1) * nc + l];
15551a454edSToby Isaac           }
15651a454edSToby Isaac           if (datatype == PETSC_SCALAR) {
15751a454edSToby Isaac             PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
15851a454edSToby Isaac 
15951a454edSToby Isaac             for (l = 0; l < nc; l++) {
16051a454edSToby Isaac               out[l * dim + j] += d * e * cf[k * nc + l];
16151a454edSToby Isaac             }
16251a454edSToby Isaac           } else {
16351a454edSToby Isaac             PetscReal *out = &((PetscReal *)D)[nc * dim * i];
16451a454edSToby Isaac 
16551a454edSToby Isaac             for (l = 0; l < nc; l++) {
16651a454edSToby Isaac               out[l * dim + j] += d * e * PetscRealPart(cf[k * nc + l]);
16751a454edSToby Isaac             }
16851a454edSToby Isaac           }
16951a454edSToby Isaac         }
17051a454edSToby Isaac       }
17151a454edSToby Isaac     }
17251a454edSToby Isaac     if (H) {
17351a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
17451a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)H)[nc * dim * dim * i];
17551a454edSToby Isaac 
17651a454edSToby Isaac         MEHess(out,cf,etaB,etaD,dim,nc,(PetscScalar));
17751a454edSToby Isaac       } else {
17851a454edSToby Isaac         PetscReal *out = &((PetscReal *)H)[nc * dim * dim * i];
17951a454edSToby Isaac 
18051a454edSToby Isaac         MEHess(out,cf,etaB,etaD,dim,nc,PetscRealPart);
18151a454edSToby Isaac       }
18251a454edSToby Isaac     }
18351a454edSToby Isaac   }
18451a454edSToby Isaac   PetscFunctionReturnVoid();
18551a454edSToby Isaac }
18651a454edSToby Isaac 
18751a454edSToby Isaac static PetscErrorCode DMFieldEvaluate_DA(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
18851a454edSToby Isaac {
18951a454edSToby Isaac   DM             dm;
19051a454edSToby Isaac   DMField_DA     *dafield;
19151a454edSToby Isaac   PetscInt       dim;
19251a454edSToby Isaac   PetscInt       N, n, nc;
19351a454edSToby Isaac   const PetscScalar *array;
19451a454edSToby Isaac   PetscReal (*coordRange)[2];
19551a454edSToby Isaac 
19651a454edSToby Isaac   PetscFunctionBegin;
19751a454edSToby Isaac   dm      = field->dm;
19851a454edSToby Isaac   nc      = field->numComponents;
19951a454edSToby Isaac   dafield = (DMField_DA *) field->data;
2009566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(points,&N));
20263a3b9bcSJacob 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);
20351a454edSToby Isaac   n = N / dim;
20451a454edSToby Isaac   coordRange = &(dafield->coordRange[0]);
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points,&array));
20651a454edSToby Isaac   MultilinearEvaluate(dim,coordRange,nc,dafield->cornerCoeffs,dafield->work,n,array,datatype,B,D,H);
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points,&array));
20851a454edSToby Isaac   PetscFunctionReturn(0);
20951a454edSToby Isaac }
21051a454edSToby Isaac 
21151a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFE_DA(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
21251a454edSToby Isaac {
21351a454edSToby Isaac   PetscInt       c, i, j, k, dim, cellsPer[3] = {0}, first[3] = {0}, whol, half;
21451a454edSToby Isaac   PetscReal      stepPer[3] = {0.};
21551a454edSToby Isaac   PetscReal      cellCoordRange[3][2] = {{0.,1.},{0.,1.},{0.,1.}};
21651a454edSToby Isaac   PetscScalar    *cellCoeffs, *work;
21751a454edSToby Isaac   DM             dm;
21851a454edSToby Isaac   DMDALocalInfo  info;
21951a454edSToby Isaac   PetscInt       cStart, cEnd;
22051a454edSToby Isaac   PetscInt       nq, nc;
22151a454edSToby Isaac   const PetscReal *q;
22251a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
22351a454edSToby Isaac   PetscScalar    *qs;
22451a454edSToby Isaac #else
22551a454edSToby Isaac   const PetscScalar *qs;
22651a454edSToby Isaac #endif
22751a454edSToby Isaac   DMField_DA     *dafield;
22851a454edSToby Isaac   PetscBool      isStride;
22951a454edSToby Isaac   const PetscInt *cells = NULL;
23051a454edSToby Isaac   PetscInt       sfirst = -1, stride = -1, nCells;
23151a454edSToby Isaac 
23251a454edSToby Isaac   PetscFunctionBegin;
23351a454edSToby Isaac   dafield = (DMField_DA *) field->data;
23451a454edSToby Isaac   dm = field->dm;
23551a454edSToby Isaac   nc = field->numComponents;
2369566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
23751a454edSToby Isaac   dim = info.dim;
23851a454edSToby Isaac   work = dafield->work;
23951a454edSToby Isaac   stepPer[0] = 1./ info.mx;
24051a454edSToby Isaac   stepPer[1] = 1./ info.my;
24151a454edSToby Isaac   stepPer[2] = 1./ info.mz;
24251a454edSToby Isaac   first[0] = info.gxs;
24351a454edSToby Isaac   first[1] = info.gys;
24451a454edSToby Isaac   first[2] = info.gzs;
24551a454edSToby Isaac   cellsPer[0] = info.gxm;
24651a454edSToby Isaac   cellsPer[1] = info.gym;
24751a454edSToby Isaac   cellsPer[2] = info.gzm;
24851a454edSToby Isaac   /* TODO: probably take components into account */
2499566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(points, NULL, NULL, &nq, &q, NULL));
25051a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
2519566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm,nq * dim,MPIU_SCALAR,&qs));
25251a454edSToby Isaac   for (i = 0; i < nq * dim; i++) qs[i] = q[i];
25351a454edSToby Isaac #else
25451a454edSToby Isaac   qs = q;
25551a454edSToby Isaac #endif
2569566063dSJacob Faibussowitsch   PetscCall(DMDAGetHeightStratum(dm,0,&cStart,&cEnd));
2579566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm,(1 << dim) * nc,MPIU_SCALAR,&cellCoeffs));
25851a454edSToby Isaac   whol = (1 << dim);
25951a454edSToby Isaac   half = whol >> 1;
2609566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(cellIS,&nCells));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS,ISSTRIDE,&isStride));
26251a454edSToby Isaac   if (isStride) {
2639566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(cellIS,&sfirst,&stride));
26451a454edSToby Isaac   } else {
2659566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(cellIS,&cells));
26651a454edSToby Isaac   }
26751a454edSToby Isaac   for (c = 0; c < nCells; c++) {
26851a454edSToby Isaac     PetscInt  cell = isStride ? (sfirst + c * stride) : cells[c];
26951a454edSToby Isaac     PetscInt  rem  = cell;
27051a454edSToby Isaac     PetscInt  ijk[3] = {0};
27151a454edSToby Isaac     void *cB, *cD, *cH;
27251a454edSToby Isaac 
27351a454edSToby Isaac     if (datatype == PETSC_SCALAR) {
27451a454edSToby Isaac       cB = B ? &((PetscScalar *)B)[nc * nq * c] : NULL;
27551a454edSToby Isaac       cD = D ? &((PetscScalar *)D)[nc * nq * dim * c] : NULL;
27651a454edSToby Isaac       cH = H ? &((PetscScalar *)H)[nc * nq * dim * dim * c] : NULL;
27751a454edSToby Isaac     } else {
27851a454edSToby Isaac       cB = B ? &((PetscReal *)B)[nc * nq * c] : NULL;
27951a454edSToby Isaac       cD = D ? &((PetscReal *)D)[nc * nq * dim * c] : NULL;
28051a454edSToby Isaac       cH = H ? &((PetscReal *)H)[nc * nq * dim * dim * c] : NULL;
28151a454edSToby Isaac     }
282*1dca8a05SBarry 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);
28351a454edSToby Isaac     for (i = 0; i < nc * whol; i++) {work[i] = dafield->cornerCoeffs[i];}
28451a454edSToby Isaac     for (j = 0; j < dim; j++) {
28551a454edSToby Isaac       PetscReal e, d;
28651a454edSToby Isaac       ijk[j] = (rem % cellsPer[j]);
28751a454edSToby Isaac       rem /= cellsPer[j];
28851a454edSToby Isaac 
28951a454edSToby Isaac       e = 2. * (ijk[j] + first[j] + 0.5) * stepPer[j] - 1.;
29051a454edSToby Isaac       d = stepPer[j];
29151a454edSToby Isaac       for (i = 0; i < half; i++) {
29251a454edSToby Isaac         for (k = 0; k < nc; k++) {
29351a454edSToby Isaac           cellCoeffs[ i         * nc + k] = work[ 2 * i * nc + k] * d;
29451a454edSToby Isaac           cellCoeffs[(i + half) * nc + k] = work[ 2 * i * nc + k] * e + work[(2 * i + 1) * nc + k];
29551a454edSToby Isaac         }
29651a454edSToby Isaac       }
29751a454edSToby Isaac       for (i = 0; i < whol * nc; i++) {work[i] = cellCoeffs[i];}
29851a454edSToby Isaac     }
29951a454edSToby Isaac     MultilinearEvaluate(dim,cellCoordRange,nc,cellCoeffs,dafield->work,nq,qs,datatype,cB,cD,cH);
30051a454edSToby Isaac   }
30151a454edSToby Isaac   if (!isStride) {
3029566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(cellIS,&cells));
30351a454edSToby Isaac   }
3049566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm,(1 << dim) * nc,MPIU_SCALAR,&cellCoeffs));
30551a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
3069566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm,nq * dim,MPIU_SCALAR,&qs));
30751a454edSToby Isaac #endif
30851a454edSToby Isaac   PetscFunctionReturn(0);
30951a454edSToby Isaac }
31051a454edSToby Isaac 
31151a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFV_DA(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
31251a454edSToby Isaac {
31351a454edSToby Isaac   PetscInt       c, i, dim, cellsPer[3] = {0}, first[3] = {0};
31451a454edSToby Isaac   PetscReal      stepPer[3] = {0.};
31551a454edSToby Isaac   DM             dm;
31651a454edSToby Isaac   DMDALocalInfo  info;
31751a454edSToby Isaac   PetscInt       cStart, cEnd, numCells;
31851a454edSToby Isaac   PetscInt       nc;
3194d1a973fSToby Isaac   PetscScalar    *points;
32051a454edSToby Isaac   DMField_DA     *dafield;
32151a454edSToby Isaac   PetscBool      isStride;
32251a454edSToby Isaac   const PetscInt *cells = NULL;
32351a454edSToby Isaac   PetscInt       sfirst = -1, stride = -1;
32451a454edSToby Isaac 
32551a454edSToby Isaac   PetscFunctionBegin;
32651a454edSToby Isaac   dafield = (DMField_DA *) field->data;
32751a454edSToby Isaac   dm = field->dm;
32851a454edSToby Isaac   nc = field->numComponents;
3299566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
33051a454edSToby Isaac   dim = info.dim;
33151a454edSToby Isaac   stepPer[0] = 1./ info.mx;
33251a454edSToby Isaac   stepPer[1] = 1./ info.my;
33351a454edSToby Isaac   stepPer[2] = 1./ info.mz;
33451a454edSToby Isaac   first[0] = info.gxs;
33551a454edSToby Isaac   first[1] = info.gys;
33651a454edSToby Isaac   first[2] = info.gzs;
33751a454edSToby Isaac   cellsPer[0] = info.gxm;
33851a454edSToby Isaac   cellsPer[1] = info.gym;
33951a454edSToby Isaac   cellsPer[2] = info.gzm;
3409566063dSJacob Faibussowitsch   PetscCall(DMDAGetHeightStratum(dm,0,&cStart,&cEnd));
3419566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(cellIS,&numCells));
3429566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm,dim * numCells,MPIU_SCALAR,&points));
3439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS,ISSTRIDE,&isStride));
34451a454edSToby Isaac   if (isStride) {
3459566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(cellIS,&sfirst,&stride));
34651a454edSToby Isaac   } else {
3479566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(cellIS,&cells));
34851a454edSToby Isaac   }
34951a454edSToby Isaac   for (c = 0; c < numCells; c++) {
35051a454edSToby Isaac     PetscInt  cell = isStride ? (sfirst + c * stride) : cells[c];
35151a454edSToby Isaac     PetscInt  rem  = cell;
35251a454edSToby Isaac     PetscInt  ijk[3] = {0};
35351a454edSToby Isaac 
354*1dca8a05SBarry 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);
35551a454edSToby Isaac     for (i = 0; i < dim; i++) {
35651a454edSToby Isaac       ijk[i] = (rem % cellsPer[i]);
35751a454edSToby Isaac       rem /= cellsPer[i];
35851a454edSToby Isaac       points[dim * c + i] = (ijk[i] + first[i] + 0.5) * stepPer[i];
35951a454edSToby Isaac     }
36051a454edSToby Isaac   }
36151a454edSToby Isaac   if (!isStride) {
3629566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(cellIS,&cells));
36351a454edSToby Isaac   }
36451a454edSToby Isaac   MultilinearEvaluate(dim,dafield->coordRange,nc,dafield->cornerCoeffs,dafield->work,numCells,points,datatype,B,D,H);
3659566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm,dim * numCells,MPIU_SCALAR,&points));
36651a454edSToby Isaac   PetscFunctionReturn(0);
36751a454edSToby Isaac }
36851a454edSToby Isaac 
369b7260050SToby Isaac static PetscErrorCode DMFieldGetDegree_DA(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
37051a454edSToby Isaac {
37151a454edSToby Isaac   DM             dm;
37251a454edSToby Isaac   PetscInt       dim, h, imin;
37351a454edSToby Isaac 
37451a454edSToby Isaac   PetscFunctionBegin;
37551a454edSToby Isaac   dm = field->dm;
3769566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS,&imin,NULL));
3779566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
37851a454edSToby Isaac   for (h = 0; h <= dim; h++) {
37951a454edSToby Isaac     PetscInt hEnd;
38051a454edSToby Isaac 
3819566063dSJacob Faibussowitsch     PetscCall(DMDAGetHeightStratum(dm,h,NULL,&hEnd));
38251a454edSToby Isaac     if (imin < hEnd) break;
38351a454edSToby Isaac   }
38451a454edSToby Isaac   dim -= h;
385b7260050SToby Isaac   if (minDegree) *minDegree = 1;
386b7260050SToby Isaac   if (maxDegree) *maxDegree = dim;
38751a454edSToby Isaac   PetscFunctionReturn(0);
38851a454edSToby Isaac }
38951a454edSToby Isaac 
39051a454edSToby Isaac static PetscErrorCode DMFieldCreateDefaultQuadrature_DA(DMField field, IS cellIS, PetscQuadrature *quad)
39151a454edSToby Isaac {
39251a454edSToby Isaac   PetscInt       h, dim, imax, imin;
39351a454edSToby Isaac   DM             dm;
39451a454edSToby Isaac 
39551a454edSToby Isaac   PetscFunctionBegin;
39651a454edSToby Isaac   dm = field->dm;
3979566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(cellIS,&imin,&imax));
3989566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
39951a454edSToby Isaac   *quad = NULL;
40051a454edSToby Isaac   for (h = 0; h <= dim; h++) {
40151a454edSToby Isaac     PetscInt hStart, hEnd;
40251a454edSToby Isaac 
4039566063dSJacob Faibussowitsch     PetscCall(DMDAGetHeightStratum(dm,h,&hStart,&hEnd));
40451a454edSToby Isaac     if (imin >= hStart && imax < hEnd) break;
40551a454edSToby Isaac   }
40651a454edSToby Isaac   dim -= h;
40751a454edSToby Isaac   if (dim > 0) {
4089566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
40951a454edSToby Isaac   }
41051a454edSToby Isaac 
41151a454edSToby Isaac   PetscFunctionReturn(0);
41251a454edSToby Isaac }
41351a454edSToby Isaac 
41451a454edSToby Isaac static PetscErrorCode DMFieldInitialize_DA(DMField field)
41551a454edSToby Isaac {
41651a454edSToby Isaac   DM             dm;
41751a454edSToby Isaac   Vec            coords = NULL;
41851a454edSToby Isaac   PetscInt       dim, i, j, k;
4194d1a973fSToby Isaac   DMField_DA     *dafield = (DMField_DA *) field->data;
42051a454edSToby Isaac 
42151a454edSToby Isaac   PetscFunctionBegin;
42251a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DA;
42351a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DA;
42451a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DA;
42551a454edSToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DA;
426b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DA;
42751a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DA;
42851a454edSToby Isaac   field->ops->view                    = DMFieldView_DA;
42951a454edSToby Isaac   dm = field->dm;
4309566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
43151a454edSToby Isaac   if (dm->coordinates) coords = dm->coordinates;
43251a454edSToby Isaac   else if (dm->coordinatesLocal) coords = dm->coordinatesLocal;
43351a454edSToby Isaac   if (coords) {
43451a454edSToby Isaac     PetscInt          n;
43551a454edSToby Isaac     const PetscScalar *array;
43651a454edSToby Isaac     PetscReal         mins[3][2] = {{PETSC_MAX_REAL,PETSC_MAX_REAL},{PETSC_MAX_REAL,PETSC_MAX_REAL},{PETSC_MAX_REAL,PETSC_MAX_REAL}};
43751a454edSToby Isaac 
4389566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coords,&n));
43951a454edSToby Isaac     n /= dim;
4409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coords,&array));
44151a454edSToby Isaac     for (i = 0, k = 0; i < n; i++) {
44251a454edSToby Isaac       for (j = 0; j < dim; j++, k++) {
44351a454edSToby Isaac         PetscReal val = PetscRealPart(array[k]);
44451a454edSToby Isaac 
44551a454edSToby Isaac         mins[j][0] = PetscMin(mins[j][0],val);
44651a454edSToby Isaac         mins[j][1] = PetscMin(mins[j][1],-val);
44751a454edSToby Isaac       }
44851a454edSToby Isaac     }
4499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coords,&array));
4501c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce((PetscReal *) mins,&(dafield->coordRange[0][0]),2*dim,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm)));
45151a454edSToby Isaac     for (j = 0; j < dim; j++) {
45251a454edSToby Isaac       dafield->coordRange[j][1] = -dafield->coordRange[j][1];
45351a454edSToby Isaac     }
45451a454edSToby Isaac   } else {
45551a454edSToby Isaac     for (j = 0; j < dim; j++) {
45651a454edSToby Isaac       dafield->coordRange[j][0] = 0.;
45751a454edSToby Isaac       dafield->coordRange[j][1] = 1.;
45851a454edSToby Isaac     }
45951a454edSToby Isaac   }
46051a454edSToby Isaac   for (j = 0; j < dim; j++) {
4614d1a973fSToby Isaac     PetscReal avg = 0.5 * (dafield->coordRange[j][1] + dafield->coordRange[j][0]);
4624d1a973fSToby Isaac     PetscReal dif = 0.5 * (dafield->coordRange[j][1] - dafield->coordRange[j][0]);
46351a454edSToby Isaac 
46451a454edSToby Isaac     dafield->coordRange[j][0] = avg;
46551a454edSToby Isaac     dafield->coordRange[j][1] = dif;
46651a454edSToby Isaac   }
46751a454edSToby Isaac   PetscFunctionReturn(0);
46851a454edSToby Isaac }
46951a454edSToby Isaac 
47051a454edSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_DA(DMField field)
47151a454edSToby Isaac {
47251a454edSToby Isaac   DMField_DA     *dafield;
47351a454edSToby Isaac 
47451a454edSToby Isaac   PetscFunctionBegin;
4759566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(field,&dafield));
47651a454edSToby Isaac   field->data = dafield;
4779566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_DA(field));
47851a454edSToby Isaac   PetscFunctionReturn(0);
47951a454edSToby Isaac }
48051a454edSToby Isaac 
48151a454edSToby Isaac PetscErrorCode DMFieldCreateDA(DM dm, PetscInt nc, const PetscScalar *cornerValues,DMField *field)
48251a454edSToby Isaac {
48351a454edSToby Isaac   DMField        b;
48451a454edSToby Isaac   DMField_DA     *dafield;
48551a454edSToby Isaac   PetscInt       dim, nv, i, j, k;
48651a454edSToby Isaac   PetscInt       half;
48751a454edSToby Isaac   PetscScalar    *cv, *cf, *work;
48851a454edSToby Isaac 
48951a454edSToby Isaac   PetscFunctionBegin;
4909566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm,nc,DMFIELD_VERTEX,&b));
4919566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b,DMFIELDDA));
49251a454edSToby Isaac   dafield = (DMField_DA *) b->data;
4939566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
49451a454edSToby Isaac   nv = (1 << dim) * nc;
4959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nv,&cv,nv,&cf,nv,&work));
49651a454edSToby Isaac   for (i = 0; i < nv; i++) cv[i] = cornerValues[i];
49751a454edSToby Isaac   for (i = 0; i < nv; i++) cf[i] = cv[i];
49851a454edSToby Isaac   dafield->cornerVals = cv;
49951a454edSToby Isaac   dafield->cornerCoeffs = cf;
50051a454edSToby Isaac   dafield->work = work;
50151a454edSToby Isaac   half = (1 << (dim - 1));
50251a454edSToby Isaac   for (i = 0; i < dim; i++) {
50351a454edSToby Isaac     PetscScalar *w;
50451a454edSToby Isaac 
50551a454edSToby Isaac     w = work;
50651a454edSToby Isaac     for (j = 0; j < half; j++) {
50751a454edSToby Isaac       for (k = 0; k < nc; k++) {
50851a454edSToby Isaac         w[j * nc + k] = 0.5 * (cf[(2 * j + 1) * nc + k] - cf[(2 * j) * nc + k]);
50951a454edSToby Isaac       }
51051a454edSToby Isaac     }
51151a454edSToby Isaac     w = &work[j * nc];
51251a454edSToby Isaac     for (j = 0; j < half; j++) {
51351a454edSToby Isaac       for (k = 0; k < nc; k++) {
51451a454edSToby Isaac         w[j * nc + k] = 0.5 * (cf[(2 * j) * nc + k] + cf[(2 * j + 1) * nc + k]);
51551a454edSToby Isaac       }
51651a454edSToby Isaac     }
51751a454edSToby Isaac     for (j = 0; j < nv; j++) {cf[j] = work[j];}
51851a454edSToby Isaac   }
51951a454edSToby Isaac   *field = b;
52051a454edSToby Isaac   PetscFunctionReturn(0);
52151a454edSToby Isaac }
522