151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 251a454edSToby Isaac 3028944aaSToby Isaac typedef struct _n_DMField_Shell 451a454edSToby Isaac { 5028944aaSToby Isaac void *ctx; 6028944aaSToby Isaac PetscErrorCode (*destroy) (DMField); 7028944aaSToby Isaac } 8028944aaSToby Isaac DMField_Shell; 9028944aaSToby Isaac 10028944aaSToby Isaac PetscErrorCode DMFieldShellGetContext(DMField field, void *ctx) 11028944aaSToby Isaac { 12028944aaSToby Isaac PetscBool flg; 13028944aaSToby Isaac PetscErrorCode ierr; 14028944aaSToby Isaac 1551a454edSToby Isaac PetscFunctionBegin; 16028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 17028944aaSToby Isaac PetscValidPointer(ctx,2); 18028944aaSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)field,DMFIELDSHELL,&flg);CHKERRQ(ierr); 19028944aaSToby Isaac if (flg) *(void**)ctx = ((DMField_Shell*)(field->data))->ctx; 20028944aaSToby Isaac else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Cannot get context from non-shell shield"); 2151a454edSToby Isaac PetscFunctionReturn(0); 2251a454edSToby Isaac } 2351a454edSToby Isaac 24028944aaSToby Isaac static PetscErrorCode DMFieldDestroy_Shell(DMField field) 25028944aaSToby Isaac { 26028944aaSToby Isaac DMField_Shell *shell = (DMField_Shell *) field->data; 27028944aaSToby Isaac PetscErrorCode ierr; 2851a454edSToby Isaac 29028944aaSToby Isaac PetscFunctionBegin; 30028944aaSToby Isaac if (shell->destroy) {ierr = (*(shell->destroy)) (field);CHKERRQ(ierr);} 31028944aaSToby Isaac ierr = PetscFree(field->data);CHKERRQ(ierr); 32028944aaSToby Isaac PetscFunctionReturn(0); 33028944aaSToby Isaac } 34028944aaSToby Isaac 35028944aaSToby Isaac PetscErrorCode DMFieldShellEvaluateFEDefault(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) 36028944aaSToby Isaac { 37028944aaSToby Isaac DM dm = field->dm; 38028944aaSToby Isaac DMField coordField; 39028944aaSToby Isaac PetscFEGeom *geom; 40028944aaSToby Isaac Vec pushforward; 41028944aaSToby Isaac PetscInt dimC, dim, numPoints, Nq, p, Nc; 42028944aaSToby Isaac PetscScalar *pfArray; 43028944aaSToby Isaac PetscErrorCode ierr; 44028944aaSToby Isaac 45028944aaSToby Isaac PetscFunctionBegin; 46028944aaSToby Isaac Nc = field->numComponents; 47028944aaSToby Isaac ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 48028944aaSToby Isaac ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom);CHKERRQ(ierr); 49028944aaSToby Isaac ierr = DMGetCoordinateDim(dm, &dimC);CHKERRQ(ierr); 50028944aaSToby Isaac ierr = PetscQuadratureGetData(quad, &dim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 51028944aaSToby Isaac ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 52028944aaSToby Isaac ierr = PetscMalloc1(dimC * Nq * numPoints, &pfArray);CHKERRQ(ierr); 53028944aaSToby Isaac for (p = 0; p < numPoints * dimC * Nq; p++) pfArray[p] = (PetscScalar) geom->v[p]; 54028944aaSToby Isaac ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * Nq * numPoints, PETSC_DETERMINE, pfArray, &pushforward);CHKERRQ(ierr); 55028944aaSToby Isaac ierr = DMFieldEvaluate(field, pushforward, type, B, D, H);CHKERRQ(ierr); 56028944aaSToby Isaac /* TODO: handle covariant/contravariant pullbacks */ 57028944aaSToby Isaac if (D) { 58028944aaSToby Isaac if (type == PETSC_SCALAR) { 59028944aaSToby Isaac PetscScalar *sD = (PetscScalar *) D; 60028944aaSToby Isaac PetscInt q; 61028944aaSToby Isaac 62028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 63028944aaSToby Isaac for (q = 0; q < Nc; q++) { 64028944aaSToby Isaac PetscScalar d[3]; 65028944aaSToby Isaac 66028944aaSToby Isaac PetscInt i, j; 67028944aaSToby Isaac 68028944aaSToby Isaac for (i = 0; i < dimC; i++) d[i] = sD[(p * Nc + q) * dimC + i]; 69028944aaSToby Isaac for (i = 0; i < dimC; i++) sD[(p * Nc + q) * dimC + i] = 0.; 70028944aaSToby Isaac for (i = 0; i < dimC; i++) { 71028944aaSToby Isaac for (j = 0; j < dimC; j++) { 72028944aaSToby Isaac sD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j]; 73028944aaSToby Isaac } 74028944aaSToby Isaac } 75028944aaSToby Isaac } 76028944aaSToby Isaac } 77028944aaSToby Isaac } else { 78028944aaSToby Isaac PetscReal *rD = (PetscReal *) D; 79028944aaSToby Isaac PetscInt q; 80028944aaSToby Isaac 81028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 82028944aaSToby Isaac for (q = 0; q < Nc; q++) { 83028944aaSToby Isaac PetscReal d[3]; 84028944aaSToby Isaac 85028944aaSToby Isaac PetscInt i, j; 86028944aaSToby Isaac 87028944aaSToby Isaac for (i = 0; i < dimC; i++) d[i] = rD[(p * Nc + q) * dimC + i]; 88028944aaSToby Isaac for (i = 0; i < dimC; i++) rD[(p * Nc + q) * dimC + i] = 0.; 89028944aaSToby Isaac for (i = 0; i < dimC; i++) { 90028944aaSToby Isaac for (j = 0; j < dimC; j++) { 91028944aaSToby Isaac rD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j]; 92028944aaSToby Isaac } 93028944aaSToby Isaac } 94028944aaSToby Isaac } 95028944aaSToby Isaac } 96028944aaSToby Isaac } 97028944aaSToby Isaac } 98028944aaSToby Isaac if (H) { 99028944aaSToby Isaac if (type == PETSC_SCALAR) { 100028944aaSToby Isaac PetscScalar *sH = (PetscScalar *) H; 101028944aaSToby Isaac PetscInt q; 102028944aaSToby Isaac 103028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 104028944aaSToby Isaac for (q = 0; q < Nc; q++) { 105028944aaSToby Isaac PetscScalar d[3][3]; 106028944aaSToby Isaac 107028944aaSToby Isaac PetscInt i, j, k, l; 108028944aaSToby Isaac 109028944aaSToby Isaac for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) d[i][j] = sH[((p * Nc + q) * dimC + i) * dimC + j]; 110028944aaSToby Isaac for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) sH[((p * Nc + q) * dimC + i) * dimC + j] = 0.; 111028944aaSToby Isaac for (i = 0; i < dimC; i++) { 112028944aaSToby Isaac for (j = 0; j < dimC; j++) { 113028944aaSToby Isaac for (k = 0; k < dimC; k++) { 114028944aaSToby Isaac for (l = 0; l < dimC; l++) { 115028944aaSToby Isaac 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]; 116028944aaSToby Isaac } 117028944aaSToby Isaac } 118028944aaSToby Isaac } 119028944aaSToby Isaac } 120028944aaSToby Isaac } 121028944aaSToby Isaac } 122028944aaSToby Isaac } else { 123028944aaSToby Isaac PetscReal *rH = (PetscReal *) H; 124028944aaSToby Isaac PetscInt q; 125028944aaSToby Isaac 126028944aaSToby Isaac for (p = 0; p < numPoints * Nq; p++) { 127028944aaSToby Isaac for (q = 0; q < Nc; q++) { 128028944aaSToby Isaac PetscReal d[3][3]; 129028944aaSToby Isaac 130028944aaSToby Isaac PetscInt i, j, k, l; 131028944aaSToby Isaac 132028944aaSToby Isaac for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) d[i][j] = rH[((p * Nc + q) * dimC + i) * dimC + j]; 133028944aaSToby Isaac for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) rH[((p * Nc + q) * dimC + i) * dimC + j] = 0.; 134028944aaSToby Isaac for (i = 0; i < dimC; i++) { 135028944aaSToby Isaac for (j = 0; j < dimC; j++) { 136028944aaSToby Isaac for (k = 0; k < dimC; k++) { 137028944aaSToby Isaac for (l = 0; l < dimC; l++) { 138028944aaSToby Isaac 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]; 139028944aaSToby Isaac } 140028944aaSToby Isaac } 141028944aaSToby Isaac } 142028944aaSToby Isaac } 143028944aaSToby Isaac } 144028944aaSToby Isaac } 145028944aaSToby Isaac } 146028944aaSToby Isaac } 147028944aaSToby Isaac ierr = VecDestroy(&pushforward);CHKERRQ(ierr); 148028944aaSToby Isaac ierr = PetscFree(pfArray);CHKERRQ(ierr); 149028944aaSToby Isaac ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr); 150028944aaSToby Isaac PetscFunctionReturn(0); 151028944aaSToby Isaac } 152028944aaSToby Isaac 153028944aaSToby Isaac PetscErrorCode DMFieldShellEvaluateFVDefault(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) 154028944aaSToby Isaac { 155028944aaSToby Isaac DM dm = field->dm; 156028944aaSToby Isaac DMField coordField; 157028944aaSToby Isaac PetscFEGeom *geom; 158028944aaSToby Isaac Vec pushforward; 159028944aaSToby Isaac PetscInt dimC, dim, numPoints, Nq, p; 160028944aaSToby Isaac PetscScalar *pfArray; 161028944aaSToby Isaac PetscQuadrature quad; 162*a2aba86cSMatthew G. Knepley MPI_Comm comm; 163028944aaSToby Isaac PetscErrorCode ierr; 164028944aaSToby Isaac 165028944aaSToby Isaac PetscFunctionBegin; 166*a2aba86cSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) field, &comm);CHKERRQ(ierr); 167*a2aba86cSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 168028944aaSToby Isaac ierr = DMGetCoordinateDim(dm, &dimC);CHKERRQ(ierr); 169*a2aba86cSMatthew G. Knepley ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 170*a2aba86cSMatthew G. Knepley ierr = DMFieldGetFVQuadrature_Internal(coordField, pointIS, &quad);CHKERRQ(ierr); 171*a2aba86cSMatthew G. Knepley PetscCheck(quad, comm, PETSC_ERR_ARG_WRONGSTATE, "coordinate field must have default quadrature for FV computation"); 172*a2aba86cSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 173*a2aba86cSMatthew G. Knepley PetscCheck(Nq == 1, comm, PETSC_ERR_ARG_WRONGSTATE, "quadrature must have only one point"); 174*a2aba86cSMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom);CHKERRQ(ierr); 175028944aaSToby Isaac ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 176028944aaSToby Isaac ierr = PetscMalloc1(dimC * numPoints, &pfArray);CHKERRQ(ierr); 177028944aaSToby Isaac for (p = 0; p < numPoints * dimC; p++) pfArray[p] = (PetscScalar) geom->v[p]; 178028944aaSToby Isaac ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * numPoints, PETSC_DETERMINE, pfArray, &pushforward);CHKERRQ(ierr); 179028944aaSToby Isaac ierr = DMFieldEvaluate(field, pushforward, type, B, D, H);CHKERRQ(ierr); 180028944aaSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 181028944aaSToby Isaac ierr = VecDestroy(&pushforward);CHKERRQ(ierr); 182028944aaSToby Isaac ierr = PetscFree(pfArray);CHKERRQ(ierr); 183028944aaSToby Isaac ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr); 184028944aaSToby Isaac PetscFunctionReturn(0); 185028944aaSToby Isaac } 186028944aaSToby Isaac 187028944aaSToby Isaac PetscErrorCode DMFieldShellSetDestroy(DMField field, PetscErrorCode (*destroy)(DMField)) 188028944aaSToby Isaac { 189028944aaSToby Isaac DMField_Shell *shell = (DMField_Shell *) field->data; 190028944aaSToby Isaac 191028944aaSToby Isaac PetscFunctionBegin; 192028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 193028944aaSToby Isaac shell->destroy = destroy; 194028944aaSToby Isaac PetscFunctionReturn(0); 195028944aaSToby Isaac } 196028944aaSToby Isaac 197028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluate(DMField field, PetscErrorCode (*evaluate)(DMField,Vec,PetscDataType,void*,void*,void*)) 198028944aaSToby Isaac { 199028944aaSToby Isaac PetscFunctionBegin; 200028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 201028944aaSToby Isaac field->ops->evaluate = evaluate; 202028944aaSToby Isaac PetscFunctionReturn(0); 203028944aaSToby Isaac } 204028944aaSToby Isaac 205028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluateFE(DMField field, PetscErrorCode (*evaluateFE)(DMField,IS,PetscQuadrature,PetscDataType,void*,void*,void*)) 206028944aaSToby Isaac { 207028944aaSToby Isaac PetscFunctionBegin; 208028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 209028944aaSToby Isaac field->ops->evaluateFE = evaluateFE; 210028944aaSToby Isaac PetscFunctionReturn(0); 211028944aaSToby Isaac } 212028944aaSToby Isaac 213028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluateFV(DMField field, PetscErrorCode (*evaluateFV)(DMField,IS,PetscDataType,void*,void*,void*)) 214028944aaSToby Isaac { 215028944aaSToby Isaac PetscFunctionBegin; 216028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 217028944aaSToby Isaac field->ops->evaluateFV = evaluateFV; 218028944aaSToby Isaac PetscFunctionReturn(0); 219028944aaSToby Isaac } 220028944aaSToby Isaac 221b7260050SToby Isaac PetscErrorCode DMFieldShellSetGetDegree(DMField field, PetscErrorCode (*getDegree)(DMField,IS,PetscInt*,PetscInt*)) 222028944aaSToby Isaac { 223028944aaSToby Isaac PetscFunctionBegin; 224028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 225b7260050SToby Isaac field->ops->getDegree = getDegree; 226028944aaSToby Isaac PetscFunctionReturn(0); 227028944aaSToby Isaac } 228028944aaSToby Isaac 229028944aaSToby Isaac PetscErrorCode DMFieldShellSetCreateDefaultQuadrature(DMField field, PetscErrorCode (*createDefaultQuadrature)(DMField,IS,PetscQuadrature*)) 230028944aaSToby Isaac { 231028944aaSToby Isaac PetscFunctionBegin; 232028944aaSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 233028944aaSToby Isaac field->ops->createDefaultQuadrature = createDefaultQuadrature; 234028944aaSToby Isaac PetscFunctionReturn(0); 235028944aaSToby Isaac } 236028944aaSToby Isaac 237028944aaSToby Isaac static PetscErrorCode DMFieldInitialize_Shell(DMField field) 238028944aaSToby Isaac { 239028944aaSToby Isaac PetscFunctionBegin; 240028944aaSToby Isaac field->ops->destroy = DMFieldDestroy_Shell; 241028944aaSToby Isaac field->ops->evaluate = NULL; 242028944aaSToby Isaac field->ops->evaluateFE = DMFieldShellEvaluateFEDefault; 243028944aaSToby Isaac field->ops->evaluateFV = DMFieldShellEvaluateFVDefault; 244b7260050SToby Isaac field->ops->getDegree = NULL; 245028944aaSToby Isaac field->ops->createDefaultQuadrature = NULL; 246028944aaSToby Isaac field->ops->view = NULL; 247028944aaSToby Isaac PetscFunctionReturn(0); 248028944aaSToby Isaac } 249028944aaSToby Isaac 250028944aaSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_Shell(DMField field) 251028944aaSToby Isaac { 252028944aaSToby Isaac DMField_Shell *shell; 253028944aaSToby Isaac PetscErrorCode ierr; 254028944aaSToby Isaac 255028944aaSToby Isaac PetscFunctionBegin; 256028944aaSToby Isaac ierr = PetscNewLog(field,&shell);CHKERRQ(ierr); 257028944aaSToby Isaac field->data = shell; 258028944aaSToby Isaac ierr = DMFieldInitialize_Shell(field);CHKERRQ(ierr); 259028944aaSToby Isaac PetscFunctionReturn(0); 260028944aaSToby Isaac } 261028944aaSToby Isaac 262028944aaSToby Isaac PetscErrorCode DMFieldCreateShell(DM dm, PetscInt numComponents, DMFieldContinuity continuity, void *ctx, DMField *field) 263028944aaSToby Isaac { 264028944aaSToby Isaac DMField b; 265028944aaSToby Isaac DMField_Shell *shell; 266028944aaSToby Isaac PetscErrorCode ierr; 267028944aaSToby Isaac 268028944aaSToby Isaac PetscFunctionBegin; 269028944aaSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 270028944aaSToby Isaac if (ctx) PetscValidPointer(ctx, 4); 271028944aaSToby Isaac PetscValidPointer(field, 5); 272028944aaSToby Isaac ierr = DMFieldCreate(dm,numComponents,continuity,&b);CHKERRQ(ierr); 273028944aaSToby Isaac ierr = DMFieldSetType(b,DMFIELDSHELL);CHKERRQ(ierr); 274028944aaSToby Isaac shell = (DMField_Shell *) b->data; 275028944aaSToby Isaac shell->ctx = ctx; 276028944aaSToby Isaac *field = b; 277028944aaSToby Isaac PetscFunctionReturn(0); 278028944aaSToby Isaac } 279