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 129371c9d4SSatish Balay static PetscErrorCode DMFieldDestroy_DA(DMField field) { 1351a454edSToby Isaac DMField_DA *dafield; 1451a454edSToby Isaac 1551a454edSToby Isaac PetscFunctionBegin; 1651a454edSToby Isaac dafield = (DMField_DA *)field->data; 179566063dSJacob Faibussowitsch PetscCall(PetscFree3(dafield->cornerVals, dafield->cornerCoeffs, dafield->work)); 189566063dSJacob Faibussowitsch PetscCall(PetscFree(dafield)); 1951a454edSToby Isaac PetscFunctionReturn(0); 2051a454edSToby Isaac } 2151a454edSToby Isaac 229371c9d4SSatish Balay static PetscErrorCode DMFieldView_DA(DMField field, PetscViewer viewer) { 2351a454edSToby Isaac DMField_DA *dafield = (DMField_DA *)field->data; 2451a454edSToby Isaac PetscBool iascii; 2551a454edSToby Isaac 2651a454edSToby Isaac PetscFunctionBegin; 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2851a454edSToby Isaac if (iascii) { 2951a454edSToby Isaac PetscInt i, c, dim; 3051a454edSToby Isaac PetscInt nc; 3151a454edSToby Isaac DM dm = field->dm; 3251a454edSToby Isaac 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Field corner values:\n")); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 359566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3651a454edSToby Isaac nc = field->numComponents; 3751a454edSToby Isaac for (i = 0, c = 0; i < (1 << dim); i++) { 3851a454edSToby Isaac PetscInt j; 3951a454edSToby Isaac 4051a454edSToby Isaac for (j = 0; j < nc; j++, c++) { 4151a454edSToby Isaac PetscScalar val = dafield->cornerVals[nc * i + j]; 4251a454edSToby Isaac 4351a454edSToby Isaac #if !defined(PETSC_USE_COMPLEX) 449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)val)); 4551a454edSToby Isaac #else 469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%g+i%g ", (double)PetscRealPart(val), (double)PetscImaginaryPart(val))); 4751a454edSToby Isaac #endif 4851a454edSToby Isaac } 499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 5051a454edSToby Isaac } 519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 5251a454edSToby Isaac } 5351a454edSToby Isaac PetscFunctionReturn(0); 5451a454edSToby Isaac } 5551a454edSToby Isaac 5651a454edSToby Isaac #define MEdot(y, A, x, m, c, cast) \ 5751a454edSToby Isaac do { \ 5851a454edSToby Isaac PetscInt _k, _l; \ 5951a454edSToby Isaac for (_k = 0; _k < (c); _k++) (y)[_k] = 0.; \ 6051a454edSToby Isaac for (_l = 0; _l < (m); _l++) { \ 619371c9d4SSatish Balay for (_k = 0; _k < (c); _k++) { (y)[_k] += cast((A)[(c)*_l + _k]) * (x)[_l]; } \ 6251a454edSToby Isaac } \ 6351a454edSToby Isaac } while (0) 6451a454edSToby Isaac 6551a454edSToby Isaac #define MEHess(out, cf, etaB, etaD, dim, nc, cast) \ 6651a454edSToby Isaac do { \ 6751a454edSToby Isaac PetscInt _m, _j, _k; \ 6851a454edSToby Isaac for (_m = 0; _m < (nc) * (dim) * (dim); _m++) (out)[_m] = 0.; \ 6951a454edSToby Isaac for (_j = 0; _j < (dim); _j++) { \ 7051a454edSToby Isaac for (_k = _j + 1; _k < (dim); _k++) { \ 7151a454edSToby Isaac PetscInt _ind = (1 << _j) + (1 << _k); \ 7251a454edSToby Isaac for (_m = 0; _m < (nc); _m++) { \ 7351a454edSToby Isaac PetscScalar c = (cf)[_m] * (etaB)[_ind] * (etaD)[_ind]; \ 7451a454edSToby Isaac (out)[(_m * (dim) + _k) * (dim) + _j] += cast(c); \ 7551a454edSToby Isaac (out)[(_m * (dim) + _j) * (dim) + _k] += cast(c); \ 7651a454edSToby Isaac } \ 7751a454edSToby Isaac } \ 7851a454edSToby Isaac } \ 7951a454edSToby Isaac } while (0) 8051a454edSToby Isaac 819371c9d4SSatish Balay 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) { 8251a454edSToby Isaac PetscInt i, j, k, l, m; 8351a454edSToby Isaac PetscInt whol = 1 << dim; 8451a454edSToby Isaac PetscInt half = whol >> 1; 8551a454edSToby Isaac 8651a454edSToby Isaac PetscFunctionBeginHot; 8751a454edSToby Isaac if (!B && !D && !H) PetscFunctionReturnVoid(); 8851a454edSToby Isaac for (i = 0; i < nPoints; i++) { 8951a454edSToby Isaac const PetscScalar *point = &points[dim * i]; 9051a454edSToby Isaac PetscReal deta[3] = {0.}; 9151a454edSToby Isaac PetscReal etaB[8] = {1., 1., 1., 1., 1., 1., 1., 1.}; 9251a454edSToby Isaac PetscReal etaD[8] = {1., 1., 1., 1., 1., 1., 1., 1.}; 9351a454edSToby Isaac PetscReal work[8]; 9451a454edSToby Isaac 9551a454edSToby Isaac for (j = 0; j < dim; j++) { 9651a454edSToby Isaac PetscReal e, d; 9751a454edSToby Isaac 9851a454edSToby Isaac e = (PetscRealPart(point[j]) - coordRange[j][0]) / coordRange[j][1]; 9951a454edSToby Isaac deta[j] = d = 1. / coordRange[j][1]; 10051a454edSToby Isaac for (k = 0; k < whol; k++) { work[k] = etaB[k]; } 10151a454edSToby Isaac for (k = 0; k < half; k++) { 10251a454edSToby Isaac etaB[k] = work[2 * k] * e; 10351a454edSToby Isaac etaB[k + half] = work[2 * k + 1]; 10451a454edSToby Isaac } 10551a454edSToby Isaac if (H) { 10651a454edSToby Isaac for (k = 0; k < whol; k++) { work[k] = etaD[k]; } 10751a454edSToby Isaac for (k = 0; k < half; k++) { 10851a454edSToby Isaac etaD[k + half] = work[2 * k]; 10951a454edSToby Isaac etaD[k] = work[2 * k + 1] * d; 11051a454edSToby Isaac } 11151a454edSToby Isaac } 11251a454edSToby Isaac } 11351a454edSToby Isaac if (B) { 11451a454edSToby Isaac if (datatype == PETSC_SCALAR) { 11551a454edSToby Isaac PetscScalar *out = &((PetscScalar *)B)[nc * i]; 11651a454edSToby Isaac 11751a454edSToby Isaac MEdot(out, cf, etaB, (1 << dim), nc, (PetscScalar)); 11851a454edSToby Isaac } else { 11951a454edSToby Isaac PetscReal *out = &((PetscReal *)B)[nc * i]; 12051a454edSToby Isaac 12151a454edSToby Isaac MEdot(out, cf, etaB, (1 << dim), nc, PetscRealPart); 12251a454edSToby Isaac } 12351a454edSToby Isaac } 12451a454edSToby Isaac if (D) { 12551a454edSToby Isaac if (datatype == PETSC_SCALAR) { 12651a454edSToby Isaac PetscScalar *out = &((PetscScalar *)D)[nc * dim * i]; 12751a454edSToby Isaac 12851a454edSToby Isaac for (m = 0; m < nc * dim; m++) out[m] = 0.; 12951a454edSToby Isaac } else { 13051a454edSToby Isaac PetscReal *out = &((PetscReal *)D)[nc * dim * i]; 13151a454edSToby Isaac 13251a454edSToby Isaac for (m = 0; m < nc * dim; m++) out[m] = 0.; 13351a454edSToby Isaac } 13451a454edSToby Isaac for (j = 0; j < dim; j++) { 13551a454edSToby Isaac PetscReal d = deta[j]; 13651a454edSToby Isaac 13751a454edSToby Isaac for (k = 0; k < whol * nc; k++) { cfWork[k] = cf[k]; } 13851a454edSToby Isaac for (k = 0; k < whol; k++) { work[k] = etaB[k]; } 13951a454edSToby Isaac for (k = 0; k < half; k++) { 14051a454edSToby Isaac PetscReal e; 14151a454edSToby Isaac 14251a454edSToby Isaac etaB[k] = work[2 * k]; 14351a454edSToby Isaac etaB[k + half] = e = work[2 * k + 1]; 14451a454edSToby Isaac 14551a454edSToby Isaac for (l = 0; l < nc; l++) { 14651a454edSToby Isaac cf[k * nc + l] = cfWork[2 * k * nc + l]; 14751a454edSToby Isaac cf[(k + half) * nc + l] = cfWork[(2 * k + 1) * nc + l]; 14851a454edSToby Isaac } 14951a454edSToby Isaac if (datatype == PETSC_SCALAR) { 15051a454edSToby Isaac PetscScalar *out = &((PetscScalar *)D)[nc * dim * i]; 15151a454edSToby Isaac 1529371c9d4SSatish Balay for (l = 0; l < nc; l++) { out[l * dim + j] += d * e * cf[k * nc + l]; } 15351a454edSToby Isaac } else { 15451a454edSToby Isaac PetscReal *out = &((PetscReal *)D)[nc * dim * i]; 15551a454edSToby Isaac 1569371c9d4SSatish Balay for (l = 0; l < nc; l++) { out[l * dim + j] += d * e * PetscRealPart(cf[k * nc + l]); } 15751a454edSToby Isaac } 15851a454edSToby Isaac } 15951a454edSToby Isaac } 16051a454edSToby Isaac } 16151a454edSToby Isaac if (H) { 16251a454edSToby Isaac if (datatype == PETSC_SCALAR) { 16351a454edSToby Isaac PetscScalar *out = &((PetscScalar *)H)[nc * dim * dim * i]; 16451a454edSToby Isaac 16551a454edSToby Isaac MEHess(out, cf, etaB, etaD, dim, nc, (PetscScalar)); 16651a454edSToby Isaac } else { 16751a454edSToby Isaac PetscReal *out = &((PetscReal *)H)[nc * dim * dim * i]; 16851a454edSToby Isaac 16951a454edSToby Isaac MEHess(out, cf, etaB, etaD, dim, nc, PetscRealPart); 17051a454edSToby Isaac } 17151a454edSToby Isaac } 17251a454edSToby Isaac } 17351a454edSToby Isaac PetscFunctionReturnVoid(); 17451a454edSToby Isaac } 17551a454edSToby Isaac 1769371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluate_DA(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) { 17751a454edSToby Isaac DM dm; 17851a454edSToby Isaac DMField_DA *dafield; 17951a454edSToby Isaac PetscInt dim; 18051a454edSToby Isaac PetscInt N, n, nc; 18151a454edSToby Isaac const PetscScalar *array; 18251a454edSToby Isaac PetscReal(*coordRange)[2]; 18351a454edSToby Isaac 18451a454edSToby Isaac PetscFunctionBegin; 18551a454edSToby Isaac dm = field->dm; 18651a454edSToby Isaac nc = field->numComponents; 18751a454edSToby Isaac dafield = (DMField_DA *)field->data; 1889566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1899566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(points, &N)); 19063a3b9bcSJacob 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); 19151a454edSToby Isaac n = N / dim; 19251a454edSToby Isaac coordRange = &(dafield->coordRange[0]); 1939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(points, &array)); 19451a454edSToby Isaac MultilinearEvaluate(dim, coordRange, nc, dafield->cornerCoeffs, dafield->work, n, array, datatype, B, D, H); 1959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(points, &array)); 19651a454edSToby Isaac PetscFunctionReturn(0); 19751a454edSToby Isaac } 19851a454edSToby Isaac 1999371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluateFE_DA(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) { 20051a454edSToby Isaac PetscInt c, i, j, k, dim, cellsPer[3] = {0}, first[3] = {0}, whol, half; 20151a454edSToby Isaac PetscReal stepPer[3] = {0.}; 2029371c9d4SSatish Balay PetscReal cellCoordRange[3][2] = { 2039371c9d4SSatish Balay {0., 1.}, 2049371c9d4SSatish Balay {0., 1.}, 2059371c9d4SSatish Balay {0., 1.} 2069371c9d4SSatish Balay }; 20751a454edSToby Isaac PetscScalar *cellCoeffs, *work; 20851a454edSToby Isaac DM dm; 20951a454edSToby Isaac DMDALocalInfo info; 21051a454edSToby Isaac PetscInt cStart, cEnd; 21151a454edSToby Isaac PetscInt nq, nc; 21251a454edSToby Isaac const PetscReal *q; 21351a454edSToby Isaac #if defined(PETSC_USE_COMPLEX) 21451a454edSToby Isaac PetscScalar *qs; 21551a454edSToby Isaac #else 21651a454edSToby Isaac const PetscScalar *qs; 21751a454edSToby Isaac #endif 21851a454edSToby Isaac DMField_DA *dafield; 21951a454edSToby Isaac PetscBool isStride; 22051a454edSToby Isaac const PetscInt *cells = NULL; 22151a454edSToby Isaac PetscInt sfirst = -1, stride = -1, nCells; 22251a454edSToby Isaac 22351a454edSToby Isaac PetscFunctionBegin; 22451a454edSToby Isaac dafield = (DMField_DA *)field->data; 22551a454edSToby Isaac dm = field->dm; 22651a454edSToby Isaac nc = field->numComponents; 2279566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 22851a454edSToby Isaac dim = info.dim; 22951a454edSToby Isaac work = dafield->work; 23051a454edSToby Isaac stepPer[0] = 1. / info.mx; 23151a454edSToby Isaac stepPer[1] = 1. / info.my; 23251a454edSToby Isaac stepPer[2] = 1. / info.mz; 23351a454edSToby Isaac first[0] = info.gxs; 23451a454edSToby Isaac first[1] = info.gys; 23551a454edSToby Isaac first[2] = info.gzs; 23651a454edSToby Isaac cellsPer[0] = info.gxm; 23751a454edSToby Isaac cellsPer[1] = info.gym; 23851a454edSToby Isaac cellsPer[2] = info.gzm; 23951a454edSToby Isaac /* TODO: probably take components into account */ 2409566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(points, NULL, NULL, &nq, &q, NULL)); 24151a454edSToby Isaac #if defined(PETSC_USE_COMPLEX) 2429566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, nq * dim, MPIU_SCALAR, &qs)); 24351a454edSToby Isaac for (i = 0; i < nq * dim; i++) qs[i] = q[i]; 24451a454edSToby Isaac #else 24551a454edSToby Isaac qs = q; 24651a454edSToby Isaac #endif 2479566063dSJacob Faibussowitsch PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd)); 2489566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs)); 24951a454edSToby Isaac whol = (1 << dim); 25051a454edSToby Isaac half = whol >> 1; 2519566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(cellIS, &nCells)); 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride)); 25351a454edSToby Isaac if (isStride) { 2549566063dSJacob Faibussowitsch PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride)); 2551baa6e33SBarry Smith } else PetscCall(ISGetIndices(cellIS, &cells)); 25651a454edSToby Isaac for (c = 0; c < nCells; c++) { 25751a454edSToby Isaac PetscInt cell = isStride ? (sfirst + c * stride) : cells[c]; 25851a454edSToby Isaac PetscInt rem = cell; 25951a454edSToby Isaac PetscInt ijk[3] = {0}; 26051a454edSToby Isaac void *cB, *cD, *cH; 26151a454edSToby Isaac 26251a454edSToby Isaac if (datatype == PETSC_SCALAR) { 26351a454edSToby Isaac cB = B ? &((PetscScalar *)B)[nc * nq * c] : NULL; 26451a454edSToby Isaac cD = D ? &((PetscScalar *)D)[nc * nq * dim * c] : NULL; 26551a454edSToby Isaac cH = H ? &((PetscScalar *)H)[nc * nq * dim * dim * c] : NULL; 26651a454edSToby Isaac } else { 26751a454edSToby Isaac cB = B ? &((PetscReal *)B)[nc * nq * c] : NULL; 26851a454edSToby Isaac cD = D ? &((PetscReal *)D)[nc * nq * dim * c] : NULL; 26951a454edSToby Isaac cH = H ? &((PetscReal *)H)[nc * nq * dim * dim * c] : NULL; 27051a454edSToby Isaac } 2711dca8a05SBarry 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); 27251a454edSToby Isaac for (i = 0; i < nc * whol; i++) { work[i] = dafield->cornerCoeffs[i]; } 27351a454edSToby Isaac for (j = 0; j < dim; j++) { 27451a454edSToby Isaac PetscReal e, d; 27551a454edSToby Isaac ijk[j] = (rem % cellsPer[j]); 27651a454edSToby Isaac rem /= cellsPer[j]; 27751a454edSToby Isaac 27851a454edSToby Isaac e = 2. * (ijk[j] + first[j] + 0.5) * stepPer[j] - 1.; 27951a454edSToby Isaac d = stepPer[j]; 28051a454edSToby Isaac for (i = 0; i < half; i++) { 28151a454edSToby Isaac for (k = 0; k < nc; k++) { 28251a454edSToby Isaac cellCoeffs[i * nc + k] = work[2 * i * nc + k] * d; 28351a454edSToby Isaac cellCoeffs[(i + half) * nc + k] = work[2 * i * nc + k] * e + work[(2 * i + 1) * nc + k]; 28451a454edSToby Isaac } 28551a454edSToby Isaac } 28651a454edSToby Isaac for (i = 0; i < whol * nc; i++) { work[i] = cellCoeffs[i]; } 28751a454edSToby Isaac } 28851a454edSToby Isaac MultilinearEvaluate(dim, cellCoordRange, nc, cellCoeffs, dafield->work, nq, qs, datatype, cB, cD, cH); 28951a454edSToby Isaac } 290*48a46eb9SPierre Jolivet if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells)); 2919566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs)); 29251a454edSToby Isaac #if defined(PETSC_USE_COMPLEX) 2939566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, nq * dim, MPIU_SCALAR, &qs)); 29451a454edSToby Isaac #endif 29551a454edSToby Isaac PetscFunctionReturn(0); 29651a454edSToby Isaac } 29751a454edSToby Isaac 2989371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluateFV_DA(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) { 29951a454edSToby Isaac PetscInt c, i, dim, cellsPer[3] = {0}, first[3] = {0}; 30051a454edSToby Isaac PetscReal stepPer[3] = {0.}; 30151a454edSToby Isaac DM dm; 30251a454edSToby Isaac DMDALocalInfo info; 30351a454edSToby Isaac PetscInt cStart, cEnd, numCells; 30451a454edSToby Isaac PetscInt nc; 3054d1a973fSToby Isaac PetscScalar *points; 30651a454edSToby Isaac DMField_DA *dafield; 30751a454edSToby Isaac PetscBool isStride; 30851a454edSToby Isaac const PetscInt *cells = NULL; 30951a454edSToby Isaac PetscInt sfirst = -1, stride = -1; 31051a454edSToby Isaac 31151a454edSToby Isaac PetscFunctionBegin; 31251a454edSToby Isaac dafield = (DMField_DA *)field->data; 31351a454edSToby Isaac dm = field->dm; 31451a454edSToby Isaac nc = field->numComponents; 3159566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 31651a454edSToby Isaac dim = info.dim; 31751a454edSToby Isaac stepPer[0] = 1. / info.mx; 31851a454edSToby Isaac stepPer[1] = 1. / info.my; 31951a454edSToby Isaac stepPer[2] = 1. / info.mz; 32051a454edSToby Isaac first[0] = info.gxs; 32151a454edSToby Isaac first[1] = info.gys; 32251a454edSToby Isaac first[2] = info.gzs; 32351a454edSToby Isaac cellsPer[0] = info.gxm; 32451a454edSToby Isaac cellsPer[1] = info.gym; 32551a454edSToby Isaac cellsPer[2] = info.gzm; 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd)); 3279566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(cellIS, &numCells)); 3289566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dim * numCells, MPIU_SCALAR, &points)); 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride)); 33051a454edSToby Isaac if (isStride) { 3319566063dSJacob Faibussowitsch PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride)); 3321baa6e33SBarry Smith } else PetscCall(ISGetIndices(cellIS, &cells)); 33351a454edSToby Isaac for (c = 0; c < numCells; c++) { 33451a454edSToby Isaac PetscInt cell = isStride ? (sfirst + c * stride) : cells[c]; 33551a454edSToby Isaac PetscInt rem = cell; 33651a454edSToby Isaac PetscInt ijk[3] = {0}; 33751a454edSToby Isaac 3381dca8a05SBarry 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); 33951a454edSToby Isaac for (i = 0; i < dim; i++) { 34051a454edSToby Isaac ijk[i] = (rem % cellsPer[i]); 34151a454edSToby Isaac rem /= cellsPer[i]; 34251a454edSToby Isaac points[dim * c + i] = (ijk[i] + first[i] + 0.5) * stepPer[i]; 34351a454edSToby Isaac } 34451a454edSToby Isaac } 345*48a46eb9SPierre Jolivet if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells)); 34651a454edSToby Isaac MultilinearEvaluate(dim, dafield->coordRange, nc, dafield->cornerCoeffs, dafield->work, numCells, points, datatype, B, D, H); 3479566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dim * numCells, MPIU_SCALAR, &points)); 34851a454edSToby Isaac PetscFunctionReturn(0); 34951a454edSToby Isaac } 35051a454edSToby Isaac 3519371c9d4SSatish Balay static PetscErrorCode DMFieldGetDegree_DA(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree) { 35251a454edSToby Isaac DM dm; 35351a454edSToby Isaac PetscInt dim, h, imin; 35451a454edSToby Isaac 35551a454edSToby Isaac PetscFunctionBegin; 35651a454edSToby Isaac dm = field->dm; 3579566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(pointIS, &imin, NULL)); 3589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 35951a454edSToby Isaac for (h = 0; h <= dim; h++) { 36051a454edSToby Isaac PetscInt hEnd; 36151a454edSToby Isaac 3629566063dSJacob Faibussowitsch PetscCall(DMDAGetHeightStratum(dm, h, NULL, &hEnd)); 36351a454edSToby Isaac if (imin < hEnd) break; 36451a454edSToby Isaac } 36551a454edSToby Isaac dim -= h; 366b7260050SToby Isaac if (minDegree) *minDegree = 1; 367b7260050SToby Isaac if (maxDegree) *maxDegree = dim; 36851a454edSToby Isaac PetscFunctionReturn(0); 36951a454edSToby Isaac } 37051a454edSToby Isaac 3719371c9d4SSatish Balay static PetscErrorCode DMFieldCreateDefaultQuadrature_DA(DMField field, IS cellIS, PetscQuadrature *quad) { 37251a454edSToby Isaac PetscInt h, dim, imax, imin; 37351a454edSToby Isaac DM dm; 37451a454edSToby Isaac 37551a454edSToby Isaac PetscFunctionBegin; 37651a454edSToby Isaac dm = field->dm; 3779566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(cellIS, &imin, &imax)); 3789566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 37951a454edSToby Isaac *quad = NULL; 38051a454edSToby Isaac for (h = 0; h <= dim; h++) { 38151a454edSToby Isaac PetscInt hStart, hEnd; 38251a454edSToby Isaac 3839566063dSJacob Faibussowitsch PetscCall(DMDAGetHeightStratum(dm, h, &hStart, &hEnd)); 38451a454edSToby Isaac if (imin >= hStart && imax < hEnd) break; 38551a454edSToby Isaac } 38651a454edSToby Isaac dim -= h; 387*48a46eb9SPierre Jolivet if (dim > 0) PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad)); 38851a454edSToby Isaac 38951a454edSToby Isaac PetscFunctionReturn(0); 39051a454edSToby Isaac } 39151a454edSToby Isaac 3929371c9d4SSatish Balay static PetscErrorCode DMFieldInitialize_DA(DMField field) { 39351a454edSToby Isaac DM dm; 39451a454edSToby Isaac Vec coords = NULL; 39551a454edSToby Isaac PetscInt dim, i, j, k; 3964d1a973fSToby Isaac DMField_DA *dafield = (DMField_DA *)field->data; 39751a454edSToby Isaac 39851a454edSToby Isaac PetscFunctionBegin; 39951a454edSToby Isaac field->ops->destroy = DMFieldDestroy_DA; 40051a454edSToby Isaac field->ops->evaluate = DMFieldEvaluate_DA; 40151a454edSToby Isaac field->ops->evaluateFE = DMFieldEvaluateFE_DA; 40251a454edSToby Isaac field->ops->evaluateFV = DMFieldEvaluateFV_DA; 403b7260050SToby Isaac field->ops->getDegree = DMFieldGetDegree_DA; 40451a454edSToby Isaac field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DA; 40551a454edSToby Isaac field->ops->view = DMFieldView_DA; 40651a454edSToby Isaac dm = field->dm; 4079566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4086858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 40951a454edSToby Isaac if (coords) { 41051a454edSToby Isaac PetscInt n; 41151a454edSToby Isaac const PetscScalar *array; 4129371c9d4SSatish Balay PetscReal mins[3][2] = { 4139371c9d4SSatish Balay {PETSC_MAX_REAL, PETSC_MAX_REAL}, 4149371c9d4SSatish Balay {PETSC_MAX_REAL, PETSC_MAX_REAL}, 4159371c9d4SSatish Balay {PETSC_MAX_REAL, PETSC_MAX_REAL} 4169371c9d4SSatish Balay }; 41751a454edSToby Isaac 4189566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &n)); 41951a454edSToby Isaac n /= dim; 4209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &array)); 42151a454edSToby Isaac for (i = 0, k = 0; i < n; i++) { 42251a454edSToby Isaac for (j = 0; j < dim; j++, k++) { 42351a454edSToby Isaac PetscReal val = PetscRealPart(array[k]); 42451a454edSToby Isaac 42551a454edSToby Isaac mins[j][0] = PetscMin(mins[j][0], val); 42651a454edSToby Isaac mins[j][1] = PetscMin(mins[j][1], -val); 42751a454edSToby Isaac } 42851a454edSToby Isaac } 4299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &array)); 4301c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((PetscReal *)mins, &(dafield->coordRange[0][0]), 2 * dim, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm))); 4319371c9d4SSatish Balay for (j = 0; j < dim; j++) { dafield->coordRange[j][1] = -dafield->coordRange[j][1]; } 43251a454edSToby Isaac } else { 43351a454edSToby Isaac for (j = 0; j < dim; j++) { 43451a454edSToby Isaac dafield->coordRange[j][0] = 0.; 43551a454edSToby Isaac dafield->coordRange[j][1] = 1.; 43651a454edSToby Isaac } 43751a454edSToby Isaac } 43851a454edSToby Isaac for (j = 0; j < dim; j++) { 4394d1a973fSToby Isaac PetscReal avg = 0.5 * (dafield->coordRange[j][1] + dafield->coordRange[j][0]); 4404d1a973fSToby Isaac PetscReal dif = 0.5 * (dafield->coordRange[j][1] - dafield->coordRange[j][0]); 44151a454edSToby Isaac 44251a454edSToby Isaac dafield->coordRange[j][0] = avg; 44351a454edSToby Isaac dafield->coordRange[j][1] = dif; 44451a454edSToby Isaac } 44551a454edSToby Isaac PetscFunctionReturn(0); 44651a454edSToby Isaac } 44751a454edSToby Isaac 4489371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMFieldCreate_DA(DMField field) { 44951a454edSToby Isaac DMField_DA *dafield; 45051a454edSToby Isaac 45151a454edSToby Isaac PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(PetscNewLog(field, &dafield)); 45351a454edSToby Isaac field->data = dafield; 4549566063dSJacob Faibussowitsch PetscCall(DMFieldInitialize_DA(field)); 45551a454edSToby Isaac PetscFunctionReturn(0); 45651a454edSToby Isaac } 45751a454edSToby Isaac 4589371c9d4SSatish Balay PetscErrorCode DMFieldCreateDA(DM dm, PetscInt nc, const PetscScalar *cornerValues, DMField *field) { 45951a454edSToby Isaac DMField b; 46051a454edSToby Isaac DMField_DA *dafield; 46151a454edSToby Isaac PetscInt dim, nv, i, j, k; 46251a454edSToby Isaac PetscInt half; 46351a454edSToby Isaac PetscScalar *cv, *cf, *work; 46451a454edSToby Isaac 46551a454edSToby Isaac PetscFunctionBegin; 4669566063dSJacob Faibussowitsch PetscCall(DMFieldCreate(dm, nc, DMFIELD_VERTEX, &b)); 4679566063dSJacob Faibussowitsch PetscCall(DMFieldSetType(b, DMFIELDDA)); 46851a454edSToby Isaac dafield = (DMField_DA *)b->data; 4699566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 47051a454edSToby Isaac nv = (1 << dim) * nc; 4719566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nv, &cv, nv, &cf, nv, &work)); 47251a454edSToby Isaac for (i = 0; i < nv; i++) cv[i] = cornerValues[i]; 47351a454edSToby Isaac for (i = 0; i < nv; i++) cf[i] = cv[i]; 47451a454edSToby Isaac dafield->cornerVals = cv; 47551a454edSToby Isaac dafield->cornerCoeffs = cf; 47651a454edSToby Isaac dafield->work = work; 47751a454edSToby Isaac half = (1 << (dim - 1)); 47851a454edSToby Isaac for (i = 0; i < dim; i++) { 47951a454edSToby Isaac PetscScalar *w; 48051a454edSToby Isaac 48151a454edSToby Isaac w = work; 48251a454edSToby Isaac for (j = 0; j < half; j++) { 4839371c9d4SSatish Balay for (k = 0; k < nc; k++) { w[j * nc + k] = 0.5 * (cf[(2 * j + 1) * nc + k] - cf[(2 * j) * nc + k]); } 48451a454edSToby Isaac } 48551a454edSToby Isaac w = &work[j * nc]; 48651a454edSToby Isaac for (j = 0; j < half; j++) { 4879371c9d4SSatish Balay for (k = 0; k < nc; k++) { w[j * nc + k] = 0.5 * (cf[(2 * j) * nc + k] + cf[(2 * j + 1) * nc + k]); } 48851a454edSToby Isaac } 48951a454edSToby Isaac for (j = 0; j < nv; j++) { cf[j] = work[j]; } 49051a454edSToby Isaac } 49151a454edSToby Isaac *field = b; 49251a454edSToby Isaac PetscFunctionReturn(0); 49351a454edSToby Isaac } 494