151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 251a454edSToby Isaac 39371c9d4SSatish Balay typedef struct _n_DMField_Shell { 4028944aaSToby Isaac void *ctx; 5028944aaSToby Isaac PetscErrorCode (*destroy)(DMField); 69371c9d4SSatish Balay } DMField_Shell; 7028944aaSToby Isaac 8d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellGetContext(DMField field, void *ctx) 9d71ae5a4SJacob Faibussowitsch { 10028944aaSToby Isaac PetscBool flg; 11028944aaSToby Isaac 1251a454edSToby Isaac PetscFunctionBegin; 13028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 144f572ea9SToby Isaac PetscAssertPointer(ctx, 2); 159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)field, DMFIELDSHELL, &flg)); 16*f4f49eeaSPierre Jolivet if (flg) *(void **)ctx = ((DMField_Shell *)field->data)->ctx; 17028944aaSToby Isaac else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Cannot get context from non-shell shield"); 183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1951a454edSToby Isaac } 2051a454edSToby Isaac 21d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldDestroy_Shell(DMField field) 22d71ae5a4SJacob Faibussowitsch { 23028944aaSToby Isaac DMField_Shell *shell = (DMField_Shell *)field->data; 2451a454edSToby Isaac 25028944aaSToby Isaac PetscFunctionBegin; 26*f4f49eeaSPierre Jolivet if (shell->destroy) PetscCall((*shell->destroy)(field)); 279566063dSJacob Faibussowitsch PetscCall(PetscFree(field->data)); 283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29028944aaSToby Isaac } 30028944aaSToby Isaac 31d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellEvaluateFEDefault(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) 32d71ae5a4SJacob Faibussowitsch { 33028944aaSToby Isaac DM dm = field->dm; 34028944aaSToby Isaac DMField coordField; 35028944aaSToby Isaac PetscFEGeom *geom; 36028944aaSToby Isaac Vec pushforward; 37028944aaSToby Isaac PetscInt dimC, dim, numPoints, Nq, p, Nc; 38028944aaSToby Isaac PetscScalar *pfArray; 39028944aaSToby Isaac 40028944aaSToby Isaac PetscFunctionBegin; 41028944aaSToby Isaac Nc = field->numComponents; 429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 439566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom)); 449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimC)); 459566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &Nq, NULL, NULL)); 469566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numPoints)); 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dimC * Nq * numPoints, &pfArray)); 48028944aaSToby Isaac for (p = 0; p < numPoints * dimC * Nq; p++) pfArray[p] = (PetscScalar)geom->v[p]; 499566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * Nq * numPoints, PETSC_DETERMINE, pfArray, &pushforward)); 509566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluate(field, pushforward, type, B, D, H)); 51028944aaSToby Isaac /* TODO: handle covariant/contravariant pullbacks */ 52028944aaSToby Isaac if (D) { 53028944aaSToby Isaac if (type == PETSC_SCALAR) { 54028944aaSToby Isaac PetscScalar *sD = (PetscScalar *)D; 55028944aaSToby Isaac PetscInt q; 56028944aaSToby Isaac 57028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 58028944aaSToby Isaac for (q = 0; q < Nc; q++) { 59028944aaSToby Isaac PetscScalar d[3]; 60028944aaSToby Isaac 61028944aaSToby Isaac PetscInt i, j; 62028944aaSToby Isaac 63028944aaSToby Isaac for (i = 0; i < dimC; i++) d[i] = sD[(p * Nc + q) * dimC + i]; 64028944aaSToby Isaac for (i = 0; i < dimC; i++) sD[(p * Nc + q) * dimC + i] = 0.; 65028944aaSToby Isaac for (i = 0; i < dimC; i++) { 66ad540459SPierre Jolivet for (j = 0; j < dimC; j++) sD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j]; 67028944aaSToby Isaac } 68028944aaSToby Isaac } 69028944aaSToby Isaac } 70028944aaSToby Isaac } else { 71028944aaSToby Isaac PetscReal *rD = (PetscReal *)D; 72028944aaSToby Isaac PetscInt q; 73028944aaSToby Isaac 74028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 75028944aaSToby Isaac for (q = 0; q < Nc; q++) { 76028944aaSToby Isaac PetscReal d[3]; 77028944aaSToby Isaac 78028944aaSToby Isaac PetscInt i, j; 79028944aaSToby Isaac 80028944aaSToby Isaac for (i = 0; i < dimC; i++) d[i] = rD[(p * Nc + q) * dimC + i]; 81028944aaSToby Isaac for (i = 0; i < dimC; i++) rD[(p * Nc + q) * dimC + i] = 0.; 82028944aaSToby Isaac for (i = 0; i < dimC; i++) { 83ad540459SPierre Jolivet for (j = 0; j < dimC; j++) rD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j]; 84028944aaSToby Isaac } 85028944aaSToby Isaac } 86028944aaSToby Isaac } 87028944aaSToby Isaac } 88028944aaSToby Isaac } 89028944aaSToby Isaac if (H) { 90028944aaSToby Isaac if (type == PETSC_SCALAR) { 91028944aaSToby Isaac PetscScalar *sH = (PetscScalar *)H; 92028944aaSToby Isaac PetscInt q; 93028944aaSToby Isaac 94028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 95028944aaSToby Isaac for (q = 0; q < Nc; q++) { 96028944aaSToby Isaac PetscScalar d[3][3]; 97028944aaSToby Isaac 98028944aaSToby Isaac PetscInt i, j, k, l; 99028944aaSToby Isaac 1009371c9d4SSatish Balay for (i = 0; i < dimC; i++) 1019371c9d4SSatish Balay for (j = 0; j < dimC; j++) d[i][j] = sH[((p * Nc + q) * dimC + i) * dimC + j]; 1029371c9d4SSatish Balay for (i = 0; i < dimC; i++) 1039371c9d4SSatish Balay for (j = 0; j < dimC; j++) sH[((p * Nc + q) * dimC + i) * dimC + j] = 0.; 104028944aaSToby Isaac for (i = 0; i < dimC; i++) { 105028944aaSToby Isaac for (j = 0; j < dimC; j++) { 106028944aaSToby Isaac for (k = 0; k < dimC; k++) { 107ad540459SPierre Jolivet for (l = 0; l < dimC; l++) sH[((p * Nc + q) * dimC + i) * dimC + j] += geom->J[(p * dimC + k) * dimC + i] * geom->J[(p * dimC + l) * dimC + j] * d[k][l]; 108028944aaSToby Isaac } 109028944aaSToby Isaac } 110028944aaSToby Isaac } 111028944aaSToby Isaac } 112028944aaSToby Isaac } 113028944aaSToby Isaac } else { 114028944aaSToby Isaac PetscReal *rH = (PetscReal *)H; 115028944aaSToby Isaac PetscInt q; 116028944aaSToby Isaac 117028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 118028944aaSToby Isaac for (q = 0; q < Nc; q++) { 119028944aaSToby Isaac PetscReal d[3][3]; 120028944aaSToby Isaac 121028944aaSToby Isaac PetscInt i, j, k, l; 122028944aaSToby Isaac 1239371c9d4SSatish Balay for (i = 0; i < dimC; i++) 1249371c9d4SSatish Balay for (j = 0; j < dimC; j++) d[i][j] = rH[((p * Nc + q) * dimC + i) * dimC + j]; 1259371c9d4SSatish Balay for (i = 0; i < dimC; i++) 1269371c9d4SSatish Balay for (j = 0; j < dimC; j++) rH[((p * Nc + q) * dimC + i) * dimC + j] = 0.; 127028944aaSToby Isaac for (i = 0; i < dimC; i++) { 128028944aaSToby Isaac for (j = 0; j < dimC; j++) { 129028944aaSToby Isaac for (k = 0; k < dimC; k++) { 130ad540459SPierre Jolivet for (l = 0; l < dimC; l++) rH[((p * Nc + q) * dimC + i) * dimC + j] += geom->J[(p * dimC + k) * dimC + i] * geom->J[(p * dimC + l) * dimC + j] * d[k][l]; 131028944aaSToby Isaac } 132028944aaSToby Isaac } 133028944aaSToby Isaac } 134028944aaSToby Isaac } 135028944aaSToby Isaac } 136028944aaSToby Isaac } 137028944aaSToby Isaac } 1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pushforward)); 1399566063dSJacob Faibussowitsch PetscCall(PetscFree(pfArray)); 1409566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142028944aaSToby Isaac } 143028944aaSToby Isaac 144d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellEvaluateFVDefault(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) 145d71ae5a4SJacob Faibussowitsch { 146028944aaSToby Isaac DM dm = field->dm; 147028944aaSToby Isaac DMField coordField; 148028944aaSToby Isaac PetscFEGeom *geom; 149028944aaSToby Isaac Vec pushforward; 150028944aaSToby Isaac PetscInt dimC, dim, numPoints, Nq, p; 151028944aaSToby Isaac PetscScalar *pfArray; 152028944aaSToby Isaac PetscQuadrature quad; 153a2aba86cSMatthew G. Knepley MPI_Comm comm; 154028944aaSToby Isaac 155028944aaSToby Isaac PetscFunctionBegin; 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)field, &comm)); 1579566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimC)); 1599566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField)); 1609566063dSJacob Faibussowitsch PetscCall(DMFieldGetFVQuadrature_Internal(coordField, pointIS, &quad)); 161a2aba86cSMatthew G. Knepley PetscCheck(quad, comm, PETSC_ERR_ARG_WRONGSTATE, "coordinate field must have default quadrature for FV computation"); 1629566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 163a2aba86cSMatthew G. Knepley PetscCheck(Nq == 1, comm, PETSC_ERR_ARG_WRONGSTATE, "quadrature must have only one point"); 1649566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom)); 1659566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numPoints)); 1669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dimC * numPoints, &pfArray)); 167028944aaSToby Isaac for (p = 0; p < numPoints * dimC; p++) pfArray[p] = (PetscScalar)geom->v[p]; 1689566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * numPoints, PETSC_DETERMINE, pfArray, &pushforward)); 1699566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluate(field, pushforward, type, B, D, H)); 1709566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 1719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pushforward)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFree(pfArray)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175028944aaSToby Isaac } 176028944aaSToby Isaac 177d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetDestroy(DMField field, PetscErrorCode (*destroy)(DMField)) 178d71ae5a4SJacob Faibussowitsch { 179028944aaSToby Isaac DMField_Shell *shell = (DMField_Shell *)field->data; 180028944aaSToby Isaac 181028944aaSToby Isaac PetscFunctionBegin; 182028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 183028944aaSToby Isaac shell->destroy = destroy; 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185028944aaSToby Isaac } 186028944aaSToby Isaac 187d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetEvaluate(DMField field, PetscErrorCode (*evaluate)(DMField, Vec, PetscDataType, void *, void *, void *)) 188d71ae5a4SJacob Faibussowitsch { 189028944aaSToby Isaac PetscFunctionBegin; 190028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 191028944aaSToby Isaac field->ops->evaluate = evaluate; 1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 193028944aaSToby Isaac } 194028944aaSToby Isaac 195d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetEvaluateFE(DMField field, PetscErrorCode (*evaluateFE)(DMField, IS, PetscQuadrature, PetscDataType, void *, void *, void *)) 196d71ae5a4SJacob Faibussowitsch { 197028944aaSToby Isaac PetscFunctionBegin; 198028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 199028944aaSToby Isaac field->ops->evaluateFE = evaluateFE; 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 201028944aaSToby Isaac } 202028944aaSToby Isaac 203d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetEvaluateFV(DMField field, PetscErrorCode (*evaluateFV)(DMField, IS, PetscDataType, void *, void *, void *)) 204d71ae5a4SJacob Faibussowitsch { 205028944aaSToby Isaac PetscFunctionBegin; 206028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 207028944aaSToby Isaac field->ops->evaluateFV = evaluateFV; 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209028944aaSToby Isaac } 210028944aaSToby Isaac 211d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetGetDegree(DMField field, PetscErrorCode (*getDegree)(DMField, IS, PetscInt *, PetscInt *)) 212d71ae5a4SJacob Faibussowitsch { 213028944aaSToby Isaac PetscFunctionBegin; 214028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 215b7260050SToby Isaac field->ops->getDegree = getDegree; 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217028944aaSToby Isaac } 218028944aaSToby Isaac 219d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldShellSetCreateDefaultQuadrature(DMField field, PetscErrorCode (*createDefaultQuadrature)(DMField, IS, PetscQuadrature *)) 220d71ae5a4SJacob Faibussowitsch { 221028944aaSToby Isaac PetscFunctionBegin; 222028944aaSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 223028944aaSToby Isaac field->ops->createDefaultQuadrature = createDefaultQuadrature; 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225028944aaSToby Isaac } 226028944aaSToby Isaac 227d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldInitialize_Shell(DMField field) 228d71ae5a4SJacob Faibussowitsch { 229028944aaSToby Isaac PetscFunctionBegin; 230028944aaSToby Isaac field->ops->destroy = DMFieldDestroy_Shell; 231028944aaSToby Isaac field->ops->evaluate = NULL; 232028944aaSToby Isaac field->ops->evaluateFE = DMFieldShellEvaluateFEDefault; 233028944aaSToby Isaac field->ops->evaluateFV = DMFieldShellEvaluateFVDefault; 234b7260050SToby Isaac field->ops->getDegree = NULL; 235028944aaSToby Isaac field->ops->createDefaultQuadrature = NULL; 236028944aaSToby Isaac field->ops->view = NULL; 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 238028944aaSToby Isaac } 239028944aaSToby Isaac 240d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate_Shell(DMField field) 241d71ae5a4SJacob Faibussowitsch { 242028944aaSToby Isaac DMField_Shell *shell; 243028944aaSToby Isaac 244028944aaSToby Isaac PetscFunctionBegin; 2454dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&shell)); 246028944aaSToby Isaac field->data = shell; 2479566063dSJacob Faibussowitsch PetscCall(DMFieldInitialize_Shell(field)); 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 249028944aaSToby Isaac } 250028944aaSToby Isaac 251d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateShell(DM dm, PetscInt numComponents, DMFieldContinuity continuity, void *ctx, DMField *field) 252d71ae5a4SJacob Faibussowitsch { 253028944aaSToby Isaac DMField b; 254028944aaSToby Isaac DMField_Shell *shell; 255028944aaSToby Isaac 256028944aaSToby Isaac PetscFunctionBegin; 257028944aaSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2584f572ea9SToby Isaac if (ctx) PetscAssertPointer(ctx, 4); 2594f572ea9SToby Isaac PetscAssertPointer(field, 5); 2609566063dSJacob Faibussowitsch PetscCall(DMFieldCreate(dm, numComponents, continuity, &b)); 2619566063dSJacob Faibussowitsch PetscCall(DMFieldSetType(b, DMFIELDSHELL)); 262028944aaSToby Isaac shell = (DMField_Shell *)b->data; 263028944aaSToby Isaac shell->ctx = ctx; 264028944aaSToby Isaac *field = b; 2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 266028944aaSToby Isaac } 267