xref: /petsc/src/dm/field/impls/da/dmfieldda.c (revision 51fff3d9040290fe99bbbf7774a3dbfc89bc981c)
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   PetscErrorCode ierr;
1851a454edSToby Isaac 
1951a454edSToby Isaac   PetscFunctionBegin;
2051a454edSToby Isaac   dafield = (DMField_DA *) field->data;
2151a454edSToby Isaac   ierr = PetscFree3(dafield->cornerVals,dafield->cornerCoeffs,dafield->work);CHKERRQ(ierr);
2251a454edSToby Isaac   ierr = PetscFree(dafield);CHKERRQ(ierr);
2351a454edSToby Isaac   PetscFunctionReturn(0);
2451a454edSToby Isaac }
2551a454edSToby Isaac 
2651a454edSToby Isaac static PetscErrorCode DMFieldView_DA(DMField field,PetscViewer viewer)
2751a454edSToby Isaac {
2851a454edSToby Isaac   DMField_DA     *dafield = (DMField_DA *) field->data;
2951a454edSToby Isaac   PetscBool      iascii;
3051a454edSToby Isaac   PetscErrorCode ierr;
3151a454edSToby Isaac 
3251a454edSToby Isaac   PetscFunctionBegin;
3351a454edSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3451a454edSToby Isaac   if (iascii) {
3551a454edSToby Isaac     PetscInt i, c, dim;
3651a454edSToby Isaac     PetscInt nc;
3751a454edSToby Isaac     DM       dm = field->dm;
3851a454edSToby Isaac 
3951a454edSToby Isaac     PetscViewerASCIIPrintf(viewer, "Field corner values:\n");CHKERRQ(ierr);
4051a454edSToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4151a454edSToby Isaac     ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
4251a454edSToby Isaac     nc = field->numComponents;
4351a454edSToby Isaac     for (i = 0, c = 0; i < (1 << dim); i++) {
4451a454edSToby Isaac       PetscInt j;
4551a454edSToby Isaac 
4651a454edSToby Isaac       for (j = 0; j < nc; j++, c++) {
4751a454edSToby Isaac         PetscScalar val = dafield->cornerVals[nc * i + j];
4851a454edSToby Isaac 
4951a454edSToby Isaac #if !defined(PETSC_USE_COMPLEX)
5051a454edSToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"%g ",(double) val);CHKERRQ(ierr);
5151a454edSToby Isaac #else
5251a454edSToby Isaac         ierr = PetscViewerASCIIPrintf(viewer,"%g+i%g ",(double) PetscRealPart(val),(double) PetscImaginaryPart(val));CHKERRQ(ierr);
5351a454edSToby Isaac #endif
5451a454edSToby Isaac       }
5551a454edSToby Isaac       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5651a454edSToby Isaac     }
5751a454edSToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
5851a454edSToby Isaac   }
5951a454edSToby Isaac   PetscFunctionReturn(0);
6051a454edSToby Isaac }
6151a454edSToby Isaac 
6251a454edSToby Isaac #define MEdot(y,A,x,m,c,cast)                          \
6351a454edSToby Isaac   do {                                                 \
6451a454edSToby Isaac     PetscInt _k, _l;                                   \
6551a454edSToby Isaac     for (_k = 0; _k < (c); _k++) (y)[_k] = 0.;         \
6651a454edSToby Isaac     for (_l = 0; _l < (m); _l++) {                     \
6751a454edSToby Isaac       for (_k = 0; _k < (c); _k++) {                   \
6851a454edSToby Isaac         (y)[_k] += cast((A)[(c) * _l + _k]) * (x)[_l]; \
6951a454edSToby Isaac       }                                                \
7051a454edSToby Isaac     }                                                  \
7151a454edSToby Isaac   } while (0)
7251a454edSToby Isaac 
7351a454edSToby Isaac #define MEHess(out,cf,etaB,etaD,dim,nc,cast)                      \
7451a454edSToby Isaac   do {                                                            \
7551a454edSToby Isaac     PetscInt _m, _j, _k;                                          \
7651a454edSToby Isaac     for (_m = 0; _m < (nc) * (dim) * (dim); _m++) (out)[_m] = 0.; \
7751a454edSToby Isaac     for (_j = 0; _j < (dim); _j++) {                              \
7851a454edSToby Isaac       for (_k = _j + 1; _k < (dim); _k++) {                       \
7951a454edSToby Isaac         PetscInt _ind = (1 << _j) + (1 << _k);                    \
8051a454edSToby Isaac         for (_m = 0; _m < (nc); _m++) {                           \
8151a454edSToby Isaac           PetscScalar c = (cf)[_m] * (etaB)[_ind] * (etaD)[_ind];   \
8251a454edSToby Isaac           (out)[(_m * (dim) + _k) * (dim) + _j] += cast(c);       \
8351a454edSToby Isaac           (out)[(_m * (dim) + _j) * (dim) + _k] += cast(c);       \
8451a454edSToby Isaac         }                                                         \
8551a454edSToby Isaac       }                                                           \
8651a454edSToby Isaac     }                                                             \
8751a454edSToby Isaac   } while (0)
8851a454edSToby Isaac 
8951a454edSToby 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)
9051a454edSToby Isaac {
9151a454edSToby Isaac   PetscInt i, j, k, l, m;
9251a454edSToby Isaac   PetscInt  whol = 1 << dim;
9351a454edSToby Isaac   PetscInt  half = whol >> 1;
9451a454edSToby Isaac 
9551a454edSToby Isaac   PetscFunctionBeginHot;
9651a454edSToby Isaac   if (!B && !D && !H) PetscFunctionReturnVoid();
9751a454edSToby Isaac   for (i = 0; i < nPoints; i++) {
9851a454edSToby Isaac     const PetscScalar *point = &points[dim * i];
9951a454edSToby Isaac     PetscReal deta[3] = {0.};
10051a454edSToby Isaac     PetscReal etaB[8] = {1.,1.,1.,1.,1.,1.,1.,1.};
10151a454edSToby Isaac     PetscReal etaD[8] = {1.,1.,1.,1.,1.,1.,1.,1.};
10251a454edSToby Isaac     PetscReal work[8];
10351a454edSToby Isaac 
10451a454edSToby Isaac     for (j = 0; j < dim; j++) {
10551a454edSToby Isaac       PetscReal e, d;
10651a454edSToby Isaac 
10751a454edSToby Isaac       e = (PetscRealPart(point[j]) - coordRange[j][0]) / coordRange[j][1];
10851a454edSToby Isaac       deta[j] = d = 1. / coordRange[j][1];
10951a454edSToby Isaac       for (k = 0; k < whol; k++) {work[k] = etaB[k];}
11051a454edSToby Isaac       for (k = 0; k < half; k++) {
11151a454edSToby Isaac         etaB[k]        = work[2 * k] * e;
11251a454edSToby Isaac         etaB[k + half] = work[2 * k + 1];
11351a454edSToby Isaac       }
11451a454edSToby Isaac       if (H) {
11551a454edSToby Isaac         for (k = 0; k < whol; k++) {work[k] = etaD[k];}
11651a454edSToby Isaac         for (k = 0; k < half; k++) {
11751a454edSToby Isaac           etaD[k + half] = work[2 * k];
11851a454edSToby Isaac           etaD[k       ] = work[2 * k + 1] * d;
11951a454edSToby Isaac         }
12051a454edSToby Isaac       }
12151a454edSToby Isaac     }
12251a454edSToby Isaac     if (B) {
12351a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
12451a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)B)[nc * i];
12551a454edSToby Isaac 
12651a454edSToby Isaac         MEdot(out,cf,etaB,(1 << dim),nc,(PetscScalar));
12751a454edSToby Isaac       } else {
12851a454edSToby Isaac         PetscReal *out = &((PetscReal *)B)[nc * i];
12951a454edSToby Isaac 
13051a454edSToby Isaac         MEdot(out,cf,etaB,(1 << dim),nc,PetscRealPart);
13151a454edSToby Isaac       }
13251a454edSToby Isaac     }
13351a454edSToby Isaac     if (D) {
13451a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
13551a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
13651a454edSToby Isaac 
13751a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
13851a454edSToby Isaac       } else {
13951a454edSToby Isaac         PetscReal *out = &((PetscReal *)D)[nc * dim * i];
14051a454edSToby Isaac 
14151a454edSToby Isaac         for (m = 0; m < nc * dim; m++) out[m] = 0.;
14251a454edSToby Isaac       }
14351a454edSToby Isaac       for (j = 0; j < dim; j++) {
14451a454edSToby Isaac         PetscReal d = deta[j];
14551a454edSToby Isaac 
14651a454edSToby Isaac         for (k = 0; k < whol * nc; k++) {cfWork[k] = cf[k];}
14751a454edSToby Isaac         for (k = 0; k < whol; k++) {work[k] = etaB[k];}
14851a454edSToby Isaac         for (k = 0; k < half; k++) {
14951a454edSToby Isaac           PetscReal e;
15051a454edSToby Isaac 
15151a454edSToby Isaac           etaB[k]        =     work[2 * k];
15251a454edSToby Isaac           etaB[k + half] = e = work[2 * k + 1];
15351a454edSToby Isaac 
15451a454edSToby Isaac           for (l = 0; l < nc; l++) {
15551a454edSToby Isaac             cf[ k         * nc + l] = cfWork[ 2 * k      * nc + l];
15651a454edSToby Isaac             cf[(k + half) * nc + l] = cfWork[(2 * k + 1) * nc + l];
15751a454edSToby Isaac           }
15851a454edSToby Isaac           if (datatype == PETSC_SCALAR) {
15951a454edSToby Isaac             PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];
16051a454edSToby Isaac 
16151a454edSToby Isaac             for (l = 0; l < nc; l++) {
16251a454edSToby Isaac               out[l * dim + j] += d * e * cf[k * nc + l];
16351a454edSToby Isaac             }
16451a454edSToby Isaac           } else {
16551a454edSToby Isaac             PetscReal *out = &((PetscReal *)D)[nc * dim * i];
16651a454edSToby Isaac 
16751a454edSToby Isaac             for (l = 0; l < nc; l++) {
16851a454edSToby Isaac               out[l * dim + j] += d * e * PetscRealPart(cf[k * nc + l]);
16951a454edSToby Isaac             }
17051a454edSToby Isaac           }
17151a454edSToby Isaac         }
17251a454edSToby Isaac       }
17351a454edSToby Isaac     }
17451a454edSToby Isaac     if (H) {
17551a454edSToby Isaac       if (datatype == PETSC_SCALAR) {
17651a454edSToby Isaac         PetscScalar *out = &((PetscScalar *)H)[nc * dim * dim * i];
17751a454edSToby Isaac 
17851a454edSToby Isaac         MEHess(out,cf,etaB,etaD,dim,nc,(PetscScalar));
17951a454edSToby Isaac       } else {
18051a454edSToby Isaac         PetscReal *out = &((PetscReal *)H)[nc * dim * dim * i];
18151a454edSToby Isaac 
18251a454edSToby Isaac         MEHess(out,cf,etaB,etaD,dim,nc,PetscRealPart);
18351a454edSToby Isaac       }
18451a454edSToby Isaac     }
18551a454edSToby Isaac   }
18651a454edSToby Isaac   PetscFunctionReturnVoid();
18751a454edSToby Isaac }
18851a454edSToby Isaac 
18951a454edSToby Isaac static PetscErrorCode DMFieldEvaluate_DA(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
19051a454edSToby Isaac {
19151a454edSToby Isaac   DM             dm;
19251a454edSToby Isaac   DMField_DA     *dafield;
19351a454edSToby Isaac   PetscInt       dim;
19451a454edSToby Isaac   PetscInt       N, n, nc;
19551a454edSToby Isaac   const PetscScalar *array;
19651a454edSToby Isaac   PetscReal (*coordRange)[2];
19751a454edSToby Isaac   PetscErrorCode ierr;
19851a454edSToby Isaac 
19951a454edSToby Isaac   PetscFunctionBegin;
20051a454edSToby Isaac   dm      = field->dm;
20151a454edSToby Isaac   nc      = field->numComponents;
20251a454edSToby Isaac   dafield = (DMField_DA *) field->data;
20351a454edSToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
20451a454edSToby Isaac   ierr = VecGetLocalSize(points,&N);CHKERRQ(ierr);
20551a454edSToby Isaac   if (N % dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Point vector size %D not divisible by coordinate dimension %D\n",N,dim);
20651a454edSToby Isaac   n = N / dim;
20751a454edSToby Isaac   coordRange = &(dafield->coordRange[0]);
20851a454edSToby Isaac   ierr = VecGetArrayRead(points,&array);CHKERRQ(ierr);
20951a454edSToby Isaac   MultilinearEvaluate(dim,coordRange,nc,dafield->cornerCoeffs,dafield->work,n,array,datatype,B,D,H);
21051a454edSToby Isaac   ierr = VecRestoreArrayRead(points,&array);CHKERRQ(ierr);
21151a454edSToby Isaac   PetscFunctionReturn(0);
21251a454edSToby Isaac }
21351a454edSToby Isaac 
21451a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFE_DA(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
21551a454edSToby Isaac {
21651a454edSToby Isaac   PetscInt       c, i, j, k, dim, cellsPer[3] = {0}, first[3] = {0}, whol, half;
21751a454edSToby Isaac   PetscReal      stepPer[3] = {0.};
21851a454edSToby Isaac   PetscReal      cellCoordRange[3][2] = {{0.,1.},{0.,1.},{0.,1.}};
21951a454edSToby Isaac   PetscScalar    *cellCoeffs, *work;
22051a454edSToby Isaac   DM             dm;
22151a454edSToby Isaac   DMDALocalInfo  info;
22251a454edSToby Isaac   PetscInt       cStart, cEnd;
22351a454edSToby Isaac   PetscInt       nq, nc;
22451a454edSToby Isaac   const PetscReal *q;
22551a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
22651a454edSToby Isaac   PetscScalar    *qs;
22751a454edSToby Isaac #else
22851a454edSToby Isaac   const PetscScalar *qs;
22951a454edSToby Isaac #endif
23051a454edSToby Isaac   DMField_DA     *dafield;
23151a454edSToby Isaac   PetscBool      isStride;
23251a454edSToby Isaac   const PetscInt *cells = NULL;
23351a454edSToby Isaac   PetscInt       sfirst = -1, stride = -1, nCells;
23451a454edSToby Isaac   PetscErrorCode ierr;
23551a454edSToby Isaac 
23651a454edSToby Isaac   PetscFunctionBegin;
23751a454edSToby Isaac   dafield = (DMField_DA *) field->data;
23851a454edSToby Isaac   dm = field->dm;
23951a454edSToby Isaac   nc = field->numComponents;
24051a454edSToby Isaac   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
24151a454edSToby Isaac   dim = info.dim;
24251a454edSToby Isaac   work = dafield->work;
24351a454edSToby Isaac   stepPer[0] = 1./ info.mx;
24451a454edSToby Isaac   stepPer[1] = 1./ info.my;
24551a454edSToby Isaac   stepPer[2] = 1./ info.mz;
24651a454edSToby Isaac   first[0] = info.gxs;
24751a454edSToby Isaac   first[1] = info.gys;
24851a454edSToby Isaac   first[2] = info.gzs;
24951a454edSToby Isaac   cellsPer[0] = info.gxm;
25051a454edSToby Isaac   cellsPer[1] = info.gym;
25151a454edSToby Isaac   cellsPer[2] = info.gzm;
25251a454edSToby Isaac   /* TODO: probably take components into account */
25351a454edSToby Isaac   ierr = PetscQuadratureGetData(points, NULL, NULL, &nq, &q, NULL);CHKERRQ(ierr);
25451a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
255*51fff3d9SToby Isaac   ierr = DMGetWorkArray(dm,nq * dim,MPIU_SCALAR,&qs);CHKERRQ(ierr);
25651a454edSToby Isaac   for (i = 0; i < nq * dim; i++) qs[i] = q[i];
25751a454edSToby Isaac #else
25851a454edSToby Isaac   qs = q;
25951a454edSToby Isaac #endif
26051a454edSToby Isaac   ierr = DMDAGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
261*51fff3d9SToby Isaac   ierr = DMGetWorkArray(dm,(1 << dim) * nc,MPIU_SCALAR,&cellCoeffs);CHKERRQ(ierr);
26251a454edSToby Isaac   whol = (1 << dim);
26351a454edSToby Isaac   half = whol >> 1;
26451a454edSToby Isaac   ierr = ISGetLocalSize(cellIS,&nCells);CHKERRQ(ierr);
26551a454edSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)cellIS,ISSTRIDE,&isStride);CHKERRQ(ierr);
26651a454edSToby Isaac   if (isStride) {
26751a454edSToby Isaac     ierr = ISStrideGetInfo(cellIS,&sfirst,&stride);CHKERRQ(ierr);
26851a454edSToby Isaac   } else {
26951a454edSToby Isaac     ierr = ISGetIndices(cellIS,&cells);CHKERRQ(ierr);
27051a454edSToby Isaac   }
27151a454edSToby Isaac   for (c = 0; c < nCells; c++) {
27251a454edSToby Isaac     PetscInt  cell = isStride ? (sfirst + c * stride) : cells[c];
27351a454edSToby Isaac     PetscInt  rem  = cell;
27451a454edSToby Isaac     PetscInt  ijk[3] = {0};
27551a454edSToby Isaac     void *cB, *cD, *cH;
27651a454edSToby Isaac 
27751a454edSToby Isaac     if (datatype == PETSC_SCALAR) {
27851a454edSToby Isaac       cB = B ? &((PetscScalar *)B)[nc * nq * c] : NULL;
27951a454edSToby Isaac       cD = D ? &((PetscScalar *)D)[nc * nq * dim * c] : NULL;
28051a454edSToby Isaac       cH = H ? &((PetscScalar *)H)[nc * nq * dim * dim * c] : NULL;
28151a454edSToby Isaac     } else {
28251a454edSToby Isaac       cB = B ? &((PetscReal *)B)[nc * nq * c] : NULL;
28351a454edSToby Isaac       cD = D ? &((PetscReal *)D)[nc * nq * dim * c] : NULL;
28451a454edSToby Isaac       cH = H ? &((PetscReal *)H)[nc * nq * dim * dim * c] : NULL;
28551a454edSToby Isaac     }
28651a454edSToby Isaac     if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Point %D not a cell [%D,%D), not implemented yet",cell,cStart,cEnd);
28751a454edSToby Isaac     for (i = 0; i < nc * whol; i++) {work[i] = dafield->cornerCoeffs[i];}
28851a454edSToby Isaac     for (j = 0; j < dim; j++) {
28951a454edSToby Isaac       PetscReal e, d;
29051a454edSToby Isaac       ijk[j] = (rem % cellsPer[j]);
29151a454edSToby Isaac       rem /= cellsPer[j];
29251a454edSToby Isaac 
29351a454edSToby Isaac       e = 2. * (ijk[j] + first[j] + 0.5) * stepPer[j] - 1.;
29451a454edSToby Isaac       d = stepPer[j];
29551a454edSToby Isaac       for (i = 0; i < half; i++) {
29651a454edSToby Isaac         for (k = 0; k < nc; k++) {
29751a454edSToby Isaac           cellCoeffs[ i         * nc + k] = work[ 2 * i * nc + k] * d;
29851a454edSToby Isaac           cellCoeffs[(i + half) * nc + k] = work[ 2 * i * nc + k] * e + work[(2 * i + 1) * nc + k];
29951a454edSToby Isaac         }
30051a454edSToby Isaac       }
30151a454edSToby Isaac       for (i = 0; i < whol * nc; i++) {work[i] = cellCoeffs[i];}
30251a454edSToby Isaac     }
30351a454edSToby Isaac     MultilinearEvaluate(dim,cellCoordRange,nc,cellCoeffs,dafield->work,nq,qs,datatype,cB,cD,cH);
30451a454edSToby Isaac   }
30551a454edSToby Isaac   if (!isStride) {
30651a454edSToby Isaac     ierr = ISRestoreIndices(cellIS,&cells);CHKERRQ(ierr);
30751a454edSToby Isaac   }
308*51fff3d9SToby Isaac   ierr = DMRestoreWorkArray(dm,(1 << dim) * nc,MPIU_SCALAR,&cellCoeffs);CHKERRQ(ierr);
30951a454edSToby Isaac #if defined(PETSC_USE_COMPLEX)
310*51fff3d9SToby Isaac   ierr = DMRestoreWorkArray(dm,nq * dim,MPIU_SCALAR,&qs);CHKERRQ(ierr);
31151a454edSToby Isaac #endif
31251a454edSToby Isaac   PetscFunctionReturn(0);
31351a454edSToby Isaac }
31451a454edSToby Isaac 
31551a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFV_DA(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
31651a454edSToby Isaac {
31751a454edSToby Isaac   PetscInt       c, i, dim, cellsPer[3] = {0}, first[3] = {0};
31851a454edSToby Isaac   PetscReal      stepPer[3] = {0.};
31951a454edSToby Isaac   DM             dm;
32051a454edSToby Isaac   DMDALocalInfo  info;
32151a454edSToby Isaac   PetscInt       cStart, cEnd, numCells;
32251a454edSToby Isaac   PetscInt       nc;
32351a454edSToby Isaac   PetscReal      *points;
32451a454edSToby Isaac   DMField_DA     *dafield;
32551a454edSToby Isaac   PetscBool      isStride;
32651a454edSToby Isaac   const PetscInt *cells = NULL;
32751a454edSToby Isaac   PetscInt       sfirst = -1, stride = -1;
32851a454edSToby Isaac   PetscErrorCode ierr;
32951a454edSToby Isaac 
33051a454edSToby Isaac   PetscFunctionBegin;
33151a454edSToby Isaac   dafield = (DMField_DA *) field->data;
33251a454edSToby Isaac   dm = field->dm;
33351a454edSToby Isaac   nc = field->numComponents;
33451a454edSToby Isaac   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
33551a454edSToby Isaac   dim = info.dim;
33651a454edSToby Isaac   stepPer[0] = 1./ info.mx;
33751a454edSToby Isaac   stepPer[1] = 1./ info.my;
33851a454edSToby Isaac   stepPer[2] = 1./ info.mz;
33951a454edSToby Isaac   first[0] = info.gxs;
34051a454edSToby Isaac   first[1] = info.gys;
34151a454edSToby Isaac   first[2] = info.gzs;
34251a454edSToby Isaac   cellsPer[0] = info.gxm;
34351a454edSToby Isaac   cellsPer[1] = info.gym;
34451a454edSToby Isaac   cellsPer[2] = info.gzm;
34551a454edSToby Isaac   ierr = DMDAGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
34651a454edSToby Isaac   ierr = ISGetLocalSize(cellIS,&numCells);CHKERRQ(ierr);
347*51fff3d9SToby Isaac   ierr = DMGetWorkArray(dm,dim * numCells,MPIU_SCALAR,&points);CHKERRQ(ierr);
34851a454edSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)cellIS,ISSTRIDE,&isStride);CHKERRQ(ierr);
34951a454edSToby Isaac   if (isStride) {
35051a454edSToby Isaac     ierr = ISStrideGetInfo(cellIS,&sfirst,&stride);CHKERRQ(ierr);
35151a454edSToby Isaac   } else {
35251a454edSToby Isaac     ierr = ISGetIndices(cellIS,&cells);CHKERRQ(ierr);
35351a454edSToby Isaac   }
35451a454edSToby Isaac   for (c = 0; c < numCells; c++) {
35551a454edSToby Isaac     PetscInt  cell = isStride ? (sfirst + c * stride) : cells[c];
35651a454edSToby Isaac     PetscInt  rem  = cell;
35751a454edSToby Isaac     PetscInt  ijk[3] = {0};
35851a454edSToby Isaac 
35951a454edSToby Isaac     if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Point %D not a cell [%D,%D), not implemented yet",cell,cStart,cEnd);
36051a454edSToby Isaac     for (i = 0; i < dim; i++) {
36151a454edSToby Isaac       ijk[i] = (rem % cellsPer[i]);
36251a454edSToby Isaac       rem /= cellsPer[i];
36351a454edSToby Isaac       points[dim * c + i] = (ijk[i] + first[i] + 0.5) * stepPer[i];
36451a454edSToby Isaac     }
36551a454edSToby Isaac   }
36651a454edSToby Isaac   if (!isStride) {
36751a454edSToby Isaac     ierr = ISRestoreIndices(cellIS,&cells);CHKERRQ(ierr);
36851a454edSToby Isaac   }
36951a454edSToby Isaac   MultilinearEvaluate(dim,dafield->coordRange,nc,dafield->cornerCoeffs,dafield->work,numCells,points,datatype,B,D,H);
370*51fff3d9SToby Isaac   ierr = DMRestoreWorkArray(dm,dim * numCells,MPIU_SCALAR,&points);CHKERRQ(ierr);
37151a454edSToby Isaac   PetscFunctionReturn(0);
37251a454edSToby Isaac }
37351a454edSToby Isaac 
37451a454edSToby Isaac static PetscErrorCode DMFieldGetFEInvariance_DA(DMField field, IS pointIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic)
37551a454edSToby Isaac {
37651a454edSToby Isaac   DM             dm;
37751a454edSToby Isaac   PetscInt       dim, h, imin;
37851a454edSToby Isaac   PetscErrorCode ierr;
37951a454edSToby Isaac 
38051a454edSToby Isaac   PetscFunctionBegin;
38151a454edSToby Isaac   dm = field->dm;
38251a454edSToby Isaac   ierr = ISGetMinMax(pointIS,&imin,NULL);CHKERRQ(ierr);
38351a454edSToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
38451a454edSToby Isaac   for (h = 0; h <= dim; h++) {
38551a454edSToby Isaac     PetscInt hEnd;
38651a454edSToby Isaac 
38751a454edSToby Isaac     ierr = DMDAGetHeightStratum(dm,h,NULL,&hEnd);CHKERRQ(ierr);
38851a454edSToby Isaac     if (imin < hEnd) break;
38951a454edSToby Isaac   }
39051a454edSToby Isaac   dim -= h;
39151a454edSToby Isaac   if (isConstant)  *isConstant  = (dim < 1) ? PETSC_TRUE : PETSC_FALSE;
39251a454edSToby Isaac   if (isAffine)    *isAffine    = (dim < 2) ? PETSC_TRUE : PETSC_FALSE;
39351a454edSToby Isaac   if (isQuadratic) *isQuadratic = (dim < 3) ? PETSC_TRUE : PETSC_FALSE;
39451a454edSToby Isaac   PetscFunctionReturn(0);
39551a454edSToby Isaac }
39651a454edSToby Isaac 
39751a454edSToby Isaac static PetscErrorCode DMFieldCreateDefaultQuadrature_DA(DMField field, IS cellIS, PetscQuadrature *quad)
39851a454edSToby Isaac {
39951a454edSToby Isaac   PetscInt       h, dim, imax, imin;
40051a454edSToby Isaac   DM             dm;
40151a454edSToby Isaac   PetscErrorCode ierr;
40251a454edSToby Isaac 
40351a454edSToby Isaac   PetscFunctionBegin;
40451a454edSToby Isaac   dm = field->dm;
40551a454edSToby Isaac   ierr = ISGetMinMax(cellIS,&imax,&imin);CHKERRQ(ierr);
40651a454edSToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
40751a454edSToby Isaac   *quad = NULL;
40851a454edSToby Isaac   for (h = 0; h <= dim; h++) {
40951a454edSToby Isaac     PetscInt hStart, hEnd;
41051a454edSToby Isaac 
41151a454edSToby Isaac     ierr = DMDAGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
41251a454edSToby Isaac     if (imin >= hStart && imax < hEnd) break;
41351a454edSToby Isaac   }
41451a454edSToby Isaac   dim -= h;
41551a454edSToby Isaac   if (dim > 0) {
41651a454edSToby Isaac     ierr = PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad);CHKERRQ(ierr);
41751a454edSToby Isaac   }
41851a454edSToby Isaac 
41951a454edSToby Isaac   PetscFunctionReturn(0);
42051a454edSToby Isaac }
42151a454edSToby Isaac 
42251a454edSToby Isaac static PetscErrorCode DMFieldInitialize_DA(DMField field)
42351a454edSToby Isaac {
42451a454edSToby Isaac   DM             dm;
42551a454edSToby Isaac   Vec            coords = NULL;
42651a454edSToby Isaac   PetscInt       dim, i, j, k;
42751a454edSToby Isaac   DMField_DA     *dafield = field->data;
42851a454edSToby Isaac   PetscErrorCode ierr;
42951a454edSToby Isaac 
43051a454edSToby Isaac   PetscFunctionBegin;
43151a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DA;
43251a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DA;
43351a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DA;
43451a454edSToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DA;
43551a454edSToby Isaac   field->ops->getFEInvariance         = DMFieldGetFEInvariance_DA;
43651a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DA;
43751a454edSToby Isaac   field->ops->view                    = DMFieldView_DA;
43851a454edSToby Isaac   dm = field->dm;
43951a454edSToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
44051a454edSToby Isaac   if (dm->coordinates) coords = dm->coordinates;
44151a454edSToby Isaac   else if (dm->coordinatesLocal) coords = dm->coordinatesLocal;
44251a454edSToby Isaac   if (coords) {
44351a454edSToby Isaac     PetscInt          n;
44451a454edSToby Isaac     const PetscScalar *array;
44551a454edSToby Isaac     PetscReal         mins[3][2] = {{PETSC_MAX_REAL,PETSC_MAX_REAL},{PETSC_MAX_REAL,PETSC_MAX_REAL},{PETSC_MAX_REAL,PETSC_MAX_REAL}};
44651a454edSToby Isaac 
44751a454edSToby Isaac     ierr = VecGetLocalSize(coords,&n);CHKERRQ(ierr);
44851a454edSToby Isaac     n /= dim;
44951a454edSToby Isaac     ierr = VecGetArrayRead(coords,&array);CHKERRQ(ierr);
45051a454edSToby Isaac     for (i = 0, k = 0; i < n; i++) {
45151a454edSToby Isaac       for (j = 0; j < dim; j++, k++) {
45251a454edSToby Isaac         PetscReal val = PetscRealPart(array[k]);
45351a454edSToby Isaac 
45451a454edSToby Isaac         mins[j][0] = PetscMin(mins[j][0],val);
45551a454edSToby Isaac         mins[j][1] = PetscMin(mins[j][1],-val);
45651a454edSToby Isaac       }
45751a454edSToby Isaac     }
45851a454edSToby Isaac     ierr = VecRestoreArrayRead(coords,&array);CHKERRQ(ierr);
45951a454edSToby Isaac     ierr = MPIU_Allreduce(mins,&(dafield->coordRange[0][0]),2*dim,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
46051a454edSToby Isaac     for (j = 0; j < dim; j++) {
46151a454edSToby Isaac       dafield->coordRange[j][1] = -dafield->coordRange[j][1];
46251a454edSToby Isaac     }
46351a454edSToby Isaac   } else {
46451a454edSToby Isaac     for (j = 0; j < dim; j++) {
46551a454edSToby Isaac       dafield->coordRange[j][0] = 0.;
46651a454edSToby Isaac       dafield->coordRange[j][1] = 1.;
46751a454edSToby Isaac     }
46851a454edSToby Isaac   }
46951a454edSToby Isaac   for (j = 0; j < dim; j++) {
47051a454edSToby Isaac     PetscScalar avg = 0.5 * (dafield->coordRange[j][1] + dafield->coordRange[j][0]);
47151a454edSToby Isaac     PetscScalar dif = 0.5 * (dafield->coordRange[j][1] - dafield->coordRange[j][0]);
47251a454edSToby Isaac 
47351a454edSToby Isaac     dafield->coordRange[j][0] = avg;
47451a454edSToby Isaac     dafield->coordRange[j][1] = dif;
47551a454edSToby Isaac   }
47651a454edSToby Isaac   PetscFunctionReturn(0);
47751a454edSToby Isaac }
47851a454edSToby Isaac 
47951a454edSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_DA(DMField field)
48051a454edSToby Isaac {
48151a454edSToby Isaac   DMField_DA     *dafield;
48251a454edSToby Isaac   PetscErrorCode ierr;
48351a454edSToby Isaac 
48451a454edSToby Isaac   PetscFunctionBegin;
48551a454edSToby Isaac   ierr = PetscNewLog(field,&dafield);CHKERRQ(ierr);
48651a454edSToby Isaac   field->data = dafield;
48751a454edSToby Isaac   ierr = DMFieldInitialize_DA(field);CHKERRQ(ierr);
48851a454edSToby Isaac   PetscFunctionReturn(0);
48951a454edSToby Isaac }
49051a454edSToby Isaac 
49151a454edSToby Isaac PetscErrorCode DMFieldCreateDA(DM dm, PetscInt nc, const PetscScalar *cornerValues,DMField *field)
49251a454edSToby Isaac {
49351a454edSToby Isaac   DMField        b;
49451a454edSToby Isaac   DMField_DA     *dafield;
49551a454edSToby Isaac   PetscInt       dim, nv, i, j, k;
49651a454edSToby Isaac   PetscInt       half;
49751a454edSToby Isaac   PetscScalar    *cv, *cf, *work;
49851a454edSToby Isaac   PetscErrorCode ierr;
49951a454edSToby Isaac 
50051a454edSToby Isaac   PetscFunctionBegin;
50151a454edSToby Isaac   ierr = DMFieldCreate(dm,nc,DMFIELD_VERTEX,&b);CHKERRQ(ierr);
50251a454edSToby Isaac   ierr = DMFieldSetType(b,DMFIELDDA);CHKERRQ(ierr);
50351a454edSToby Isaac   dafield = (DMField_DA *) b->data;
50451a454edSToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
50551a454edSToby Isaac   nv = (1 << dim) * nc;
50651a454edSToby Isaac   ierr = PetscMalloc3(nv,&cv,nv,&cf,nv,&work);CHKERRQ(ierr);
50751a454edSToby Isaac   for (i = 0; i < nv; i++) cv[i] = cornerValues[i];
50851a454edSToby Isaac   for (i = 0; i < nv; i++) cf[i] = cv[i];
50951a454edSToby Isaac   dafield->cornerVals = cv;
51051a454edSToby Isaac   dafield->cornerCoeffs = cf;
51151a454edSToby Isaac   dafield->work = work;
51251a454edSToby Isaac   half = (1 << (dim - 1));
51351a454edSToby Isaac   for (i = 0; i < dim; i++) {
51451a454edSToby Isaac     PetscScalar *w;
51551a454edSToby Isaac 
51651a454edSToby Isaac     w = work;
51751a454edSToby Isaac     for (j = 0; j < half; j++) {
51851a454edSToby Isaac       for (k = 0; k < nc; k++) {
51951a454edSToby Isaac         w[j * nc + k] = 0.5 * (cf[(2 * j + 1) * nc + k] - cf[(2 * j) * nc + k]);
52051a454edSToby Isaac       }
52151a454edSToby Isaac     }
52251a454edSToby Isaac     w = &work[j * nc];
52351a454edSToby Isaac     for (j = 0; j < half; j++) {
52451a454edSToby Isaac       for (k = 0; k < nc; k++) {
52551a454edSToby Isaac         w[j * nc + k] = 0.5 * (cf[(2 * j) * nc + k] + cf[(2 * j + 1) * nc + k]);
52651a454edSToby Isaac       }
52751a454edSToby Isaac     }
52851a454edSToby Isaac     for (j = 0; j < nv; j++) {cf[j] = work[j];}
52951a454edSToby Isaac   }
53051a454edSToby Isaac   *field = b;
53151a454edSToby Isaac   PetscFunctionReturn(0);
53251a454edSToby Isaac }
533