13da551e6SToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 20145028aSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/ 30145028aSToby Isaac #include <petscdmplex.h> 43da551e6SToby Isaac 59371c9d4SSatish Balay const char *const DMFieldContinuities[] = {"VERTEX", "EDGE", "FACET", "CELL", NULL}; 63da551e6SToby Isaac 7d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm, PetscInt numComponents, DMFieldContinuity continuity, DMField *field) 8d71ae5a4SJacob Faibussowitsch { 93da551e6SToby Isaac DMField b; 103da551e6SToby Isaac 113da551e6SToby Isaac PetscFunctionBegin; 123da551e6SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 133da551e6SToby Isaac PetscValidPointer(field, 2); 149566063dSJacob Faibussowitsch PetscCall(DMFieldInitializePackage()); 153da551e6SToby Isaac 169566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(b, DMFIELD_CLASSID, "DMField", "Field over DM", "DM", PetscObjectComm((PetscObject)dm), DMFieldDestroy, DMFieldView)); 179566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 183da551e6SToby Isaac b->dm = dm; 193da551e6SToby Isaac b->continuity = continuity; 203da551e6SToby Isaac b->numComponents = numComponents; 213da551e6SToby Isaac *field = b; 223da551e6SToby Isaac PetscFunctionReturn(0); 233da551e6SToby Isaac } 243da551e6SToby Isaac 253da551e6SToby Isaac /*@ 26*dce8aebaSBarry Smith DMFieldDestroy - destroy a `DMField` 273da551e6SToby Isaac 283da551e6SToby Isaac Collective 293da551e6SToby Isaac 304165533cSJose E. Roman Input Parameter: 31*dce8aebaSBarry Smith . field - address of `DMField` 323da551e6SToby Isaac 333da551e6SToby Isaac Level: advanced 343da551e6SToby Isaac 35*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldCreate()` 363da551e6SToby Isaac @*/ 37d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldDestroy(DMField *field) 38d71ae5a4SJacob Faibussowitsch { 393da551e6SToby Isaac PetscFunctionBegin; 403da551e6SToby Isaac if (!*field) PetscFunctionReturn(0); 413da551e6SToby Isaac PetscValidHeaderSpecific((*field), DMFIELD_CLASSID, 1); 429371c9d4SSatish Balay if (--((PetscObject)(*field))->refct > 0) { 439371c9d4SSatish Balay *field = NULL; 449371c9d4SSatish Balay PetscFunctionReturn(0); 459371c9d4SSatish Balay } 46dbbe0bcdSBarry Smith PetscTryTypeMethod((*field), destroy); 479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&((*field)->dm))); 489566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(field)); 493da551e6SToby Isaac PetscFunctionReturn(0); 503da551e6SToby Isaac } 513da551e6SToby Isaac 523da551e6SToby Isaac /*@C 533da551e6SToby Isaac DMFieldView - view a DMField 543da551e6SToby Isaac 553da551e6SToby Isaac Collective 563da551e6SToby Isaac 574165533cSJose E. Roman Input Parameters: 583da551e6SToby Isaac + field - DMField 593da551e6SToby Isaac - viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD 603da551e6SToby Isaac 613da551e6SToby Isaac Level: advanced 623da551e6SToby Isaac 63*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldCreate()` 643da551e6SToby Isaac @*/ 65d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldView(DMField field, PetscViewer viewer) 66d71ae5a4SJacob Faibussowitsch { 673da551e6SToby Isaac PetscBool iascii; 683da551e6SToby Isaac 693da551e6SToby Isaac PetscFunctionBegin; 703da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 719566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field), &viewer)); 723da551e6SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 733da551e6SToby Isaac PetscCheckSameComm(field, 1, viewer, 2); 749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 753da551e6SToby Isaac if (iascii) { 769566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)field, viewer)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " components\n", field->numComponents)); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s continuity\n", DMFieldContinuities[field->continuity])); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 819566063dSJacob Faibussowitsch PetscCall(DMView(field->dm, viewer)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 833da551e6SToby Isaac } 84dbbe0bcdSBarry Smith PetscTryTypeMethod(field, view, viewer); 851baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPopTab(viewer)); 863da551e6SToby Isaac PetscFunctionReturn(0); 873da551e6SToby Isaac } 883da551e6SToby Isaac 893da551e6SToby Isaac /*@C 90*dce8aebaSBarry Smith DMFieldSetType - set the `DMField` implementation 913da551e6SToby Isaac 92d083f849SBarry Smith Collective on field 933da551e6SToby Isaac 943da551e6SToby Isaac Input Parameters: 95*dce8aebaSBarry Smith + field - the `DMField` context 963da551e6SToby Isaac - type - a known method 973da551e6SToby Isaac 983da551e6SToby Isaac Notes: 993da551e6SToby Isaac See "include/petscvec.h" for available methods (for instance) 100*dce8aebaSBarry Smith + `DMFIELDDA` - a field defined only by its values at the corners of a `DMDA` 101*dce8aebaSBarry Smith . `DMFIELDDS` - a field defined by a discretization over a mesh set with `DMSetField()` 102*dce8aebaSBarry Smith - `DMFIELDSHELL` - a field defined by arbitrary callbacks 1033da551e6SToby Isaac 1043da551e6SToby Isaac Level: advanced 1053da551e6SToby Isaac 106*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldGetType()`, `DMFieldType`, 1073da551e6SToby Isaac @*/ 108d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldSetType(DMField field, DMFieldType type) 109d71ae5a4SJacob Faibussowitsch { 1103da551e6SToby Isaac PetscBool match; 1115f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(DMField); 1123da551e6SToby Isaac 1133da551e6SToby Isaac PetscFunctionBegin; 1143da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 1153da551e6SToby Isaac PetscValidCharPointer(type, 2); 1163da551e6SToby Isaac 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)field, type, &match)); 1183da551e6SToby Isaac if (match) PetscFunctionReturn(0); 1193da551e6SToby Isaac 1209566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(DMFieldList, type, &r)); 1215f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested DMField type %s", type); 1223da551e6SToby Isaac /* Destroy the previous private DMField context */ 123dbbe0bcdSBarry Smith PetscTryTypeMethod(field, destroy); 1245f80ce2aSJacob Faibussowitsch 1259566063dSJacob Faibussowitsch PetscCall(PetscMemzero(field->ops, sizeof(*field->ops))); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)field, type)); 1273da551e6SToby Isaac field->ops->create = r; 1289566063dSJacob Faibussowitsch PetscCall((*r)(field)); 1293da551e6SToby Isaac PetscFunctionReturn(0); 1303da551e6SToby Isaac } 1313da551e6SToby Isaac 1323da551e6SToby Isaac /*@C 133*dce8aebaSBarry Smith DMFieldGetType - Gets the `DMFieldType` name (as a string) from the `DMField`. 1343da551e6SToby Isaac 1353da551e6SToby Isaac Not Collective 1363da551e6SToby Isaac 1373da551e6SToby Isaac Input Parameter: 138*dce8aebaSBarry Smith . field - The `DMField` context 1393da551e6SToby Isaac 1403da551e6SToby Isaac Output Parameter: 141*dce8aebaSBarry Smith . type - The `DMFieldType` name 1423da551e6SToby Isaac 1433da551e6SToby Isaac Level: advanced 1443da551e6SToby Isaac 145*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldSetType()`, `DMFieldType` 1463da551e6SToby Isaac @*/ 147d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 148d71ae5a4SJacob Faibussowitsch { 1493da551e6SToby Isaac PetscFunctionBegin; 150064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 1513da551e6SToby Isaac PetscValidPointer(type, 2); 1529566063dSJacob Faibussowitsch PetscCall(DMFieldRegisterAll()); 1533da551e6SToby Isaac *type = ((PetscObject)field)->type_name; 1543da551e6SToby Isaac PetscFunctionReturn(0); 1553da551e6SToby Isaac } 1563da551e6SToby Isaac 157da311338SToby Isaac /*@ 158da311338SToby Isaac DMFieldGetNumComponents - Returns the number of components in the field 159da311338SToby Isaac 160da311338SToby Isaac Not collective 161da311338SToby Isaac 162da311338SToby Isaac Input Parameter: 163*dce8aebaSBarry Smith . field - The `DMField` object 164da311338SToby Isaac 165da311338SToby Isaac Output Parameter: 166da311338SToby Isaac . nc - The number of field components 167da311338SToby Isaac 168da311338SToby Isaac Level: intermediate 169da311338SToby Isaac 170*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldEvaluate()` 171da311338SToby Isaac @*/ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 173d71ae5a4SJacob Faibussowitsch { 1743da551e6SToby Isaac PetscFunctionBegin; 1753da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 1763da551e6SToby Isaac PetscValidIntPointer(nc, 2); 1773da551e6SToby Isaac *nc = field->numComponents; 1783da551e6SToby Isaac PetscFunctionReturn(0); 1793da551e6SToby Isaac } 1803da551e6SToby Isaac 181da311338SToby Isaac /*@ 182*dce8aebaSBarry Smith DMFieldGetDM - Returns the `DM` for the manifold over which the field is defined. 183da311338SToby Isaac 184da311338SToby Isaac Not collective 185da311338SToby Isaac 186da311338SToby Isaac Input Parameter: 187*dce8aebaSBarry Smith . field - The `DMField` object 188da311338SToby Isaac 189da311338SToby Isaac Output Parameter: 190*dce8aebaSBarry Smith . dm - The `DM` object 191da311338SToby Isaac 192da311338SToby Isaac Level: intermediate 193da311338SToby Isaac 194*dce8aebaSBarry Smith .seealso: `DMField`, `DM`, `DMFieldEvaluate()` 195da311338SToby Isaac @*/ 196d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 197d71ae5a4SJacob Faibussowitsch { 1983da551e6SToby Isaac PetscFunctionBegin; 1993da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 2003da551e6SToby Isaac PetscValidPointer(dm, 2); 2013da551e6SToby Isaac *dm = field->dm; 2023da551e6SToby Isaac PetscFunctionReturn(0); 2033da551e6SToby Isaac } 2043da551e6SToby Isaac 205da311338SToby Isaac /*@ 206da311338SToby Isaac DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 207da311338SToby Isaac 208d083f849SBarry Smith Collective on points 209da311338SToby Isaac 210d8d19677SJose E. Roman Input Parameters: 211*dce8aebaSBarry Smith + field - The `DMField` object 212da311338SToby Isaac . points - The points at which to evaluate the field. Should have size d x n, 213da311338SToby Isaac where d is the coordinate dimension of the manifold and n is the number 214da311338SToby Isaac of points 215*dce8aebaSBarry Smith - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. 216*dce8aebaSBarry Smith If the field is complex and datatype is `PETSC_REAL`, the real part of the 217da311338SToby Isaac field is returned. 218da311338SToby Isaac 219d8d19677SJose E. Roman Output Parameters: 220da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 221da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 222da311338SToby Isaac then by point. 223da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 224da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 225da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 226da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 227da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 228da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 229da311338SToby Isaac 230da311338SToby Isaac Level: intermediate 231da311338SToby Isaac 232*dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldGetDM()`, `DMFieldGetNumComponents()`, `DMFieldEvaluateFE()`, `DMFieldEvaluateFV()`, `PetscDataType` 233da311338SToby Isaac @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 235d71ae5a4SJacob Faibussowitsch { 2363da551e6SToby Isaac PetscFunctionBegin; 2373da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 2383da551e6SToby Isaac PetscValidHeaderSpecific(points, VEC_CLASSID, 2); 239064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B, 4); 240064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D, 5); 241064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H, 6); 2423da551e6SToby Isaac if (field->ops->evaluate) { 2439566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluate)(field, points, datatype, B, D, H)); 2443da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type"); 2453da551e6SToby Isaac PetscFunctionReturn(0); 2463da551e6SToby Isaac } 2473da551e6SToby Isaac 248da311338SToby Isaac /*@ 249da311338SToby Isaac DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 250da311338SToby Isaac quadrature points on a reference point. The derivatives are taken with respect to the 251da311338SToby Isaac reference coordinates. 252da311338SToby Isaac 253da311338SToby Isaac Not collective 254da311338SToby Isaac 255d8d19677SJose E. Roman Input Parameters: 256*dce8aebaSBarry Smith + field - The `DMField` object 257da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 258da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field. 259da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 260da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 261da311338SToby Isaac field is returned. 262da311338SToby Isaac 263d8d19677SJose E. Roman Output Parameters: 264da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 265da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 266da311338SToby Isaac then by point. 267da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 268da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 269da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 270da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 271da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 272da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 273da311338SToby Isaac 274da311338SToby Isaac Level: intermediate 275da311338SToby Isaac 276*dce8aebaSBarry Smith .seealso: `DMField`, `DM`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFV()` 277da311338SToby Isaac @*/ 278d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 279d71ae5a4SJacob Faibussowitsch { 2803da551e6SToby Isaac PetscFunctionBegin; 2813da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 2823da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 2833da551e6SToby Isaac PetscValidHeader(points, 3); 284064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B, 5); 285064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D, 6); 286064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H, 7); 2873da551e6SToby Isaac if (field->ops->evaluateFE) { 2889566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluateFE)(field, cellIS, points, datatype, B, D, H)); 2893da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type"); 2903da551e6SToby Isaac PetscFunctionReturn(0); 2913da551e6SToby Isaac } 2923da551e6SToby Isaac 293da311338SToby Isaac /*@ 294da311338SToby Isaac DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 295da311338SToby Isaac 296da311338SToby Isaac Not collective 297da311338SToby Isaac 298d8d19677SJose E. Roman Input Parameters: 299*dce8aebaSBarry Smith + field - The `DMField` object 300da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 301*dce8aebaSBarry Smith - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. 302*dce8aebaSBarry Smith If the field is complex and datatype is `PETSC_REAL`, the real part of the 303da311338SToby Isaac field is returned. 304da311338SToby Isaac 305d8d19677SJose E. Roman Output Parameters: 306da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 307da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 308da311338SToby Isaac then by point. 309da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 310da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 311da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 312da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 313da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 314da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 315da311338SToby Isaac 316da311338SToby Isaac Level: intermediate 317da311338SToby Isaac 318*dce8aebaSBarry Smith .seealso: `DMField`, `IS`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFE()`, `PetscDataType` 319da311338SToby Isaac @*/ 320d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 321d71ae5a4SJacob Faibussowitsch { 3223da551e6SToby Isaac PetscFunctionBegin; 3233da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3243da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 325064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B, 4); 326064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D, 5); 327064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H, 6); 3283da551e6SToby Isaac if (field->ops->evaluateFV) { 3299566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluateFV)(field, cellIS, datatype, B, D, H)); 3303da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type"); 3313da551e6SToby Isaac PetscFunctionReturn(0); 3323da551e6SToby Isaac } 3333da551e6SToby Isaac 334da311338SToby Isaac /*@ 335b7260050SToby Isaac DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the 336da311338SToby Isaac reference element 337da311338SToby Isaac 338da311338SToby Isaac Not collective 339da311338SToby Isaac 3404165533cSJose E. Roman Input Parameters: 341*dce8aebaSBarry Smith + field - the `DMField` object 342da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance 343da311338SToby Isaac 3444165533cSJose E. Roman Output Parameters: 345b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element 346b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element 347da311338SToby Isaac 348da311338SToby Isaac Level: intermediate 349da311338SToby Isaac 350*dce8aebaSBarry Smith .seealso: `DMField`, `IS`, `DMFieldEvaluateFE()` 351da311338SToby Isaac @*/ 352d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) 353d71ae5a4SJacob Faibussowitsch { 3543da551e6SToby Isaac PetscFunctionBegin; 3553da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3563da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 357dadcf809SJacob Faibussowitsch if (minDegree) PetscValidIntPointer(minDegree, 3); 358dadcf809SJacob Faibussowitsch if (maxDegree) PetscValidIntPointer(maxDegree, 4); 3593da551e6SToby Isaac 360b7260050SToby Isaac if (minDegree) *minDegree = -1; 361b7260050SToby Isaac if (maxDegree) *maxDegree = PETSC_MAX_INT; 3623da551e6SToby Isaac 3631baa6e33SBarry Smith if (field->ops->getDegree) PetscCall((*field->ops->getDegree)(field, cellIS, minDegree, maxDegree)); 3643da551e6SToby Isaac PetscFunctionReturn(0); 3653da551e6SToby Isaac } 3663da551e6SToby Isaac 367da311338SToby Isaac /*@ 368da311338SToby Isaac DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 369da311338SToby Isaac points via pullback onto the reference element 370da311338SToby Isaac 371da311338SToby Isaac Not collective 372da311338SToby Isaac 3734165533cSJose E. Roman Input Parameters: 374*dce8aebaSBarry Smith + field - the `DMField` object 375da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field 376da311338SToby Isaac 3774165533cSJose E. Roman Output Parameter: 378*dce8aebaSBarry Smith . quad - a `PetscQuadrature` object 379da311338SToby Isaac 380da311338SToby Isaac Level: developer 381da311338SToby Isaac 382*dce8aebaSBarry Smith .seealso: `DMField`, `PetscQuadrature`, `IS`, `DMFieldEvaluteFE()`, `DMFieldGetDegree()` 383da311338SToby Isaac @*/ 384d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 385d71ae5a4SJacob Faibussowitsch { 3863da551e6SToby Isaac PetscFunctionBegin; 3873da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3883da551e6SToby Isaac PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2); 3893da551e6SToby Isaac PetscValidPointer(quad, 3); 3903da551e6SToby Isaac 3913da551e6SToby Isaac *quad = NULL; 392dbbe0bcdSBarry Smith PetscTryTypeMethod(field, createDefaultQuadrature, pointIS, quad); 3933da551e6SToby Isaac PetscFunctionReturn(0); 3943da551e6SToby Isaac } 3953da551e6SToby Isaac 396da311338SToby Isaac /*@C 397da311338SToby Isaac DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 398da311338SToby Isaac 399da311338SToby Isaac Not collective 400da311338SToby Isaac 4014165533cSJose E. Roman Input Parameters: 402*dce8aebaSBarry Smith + field - the `DMField` object 403da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field 404da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors 405da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 406da311338SToby Isaac be calculated 407da311338SToby Isaac 4084165533cSJose E. Roman Output Parameter: 409da311338SToby Isaac . geom - the geometric factors 410da311338SToby Isaac 411da311338SToby Isaac Level: developer 412da311338SToby Isaac 413*dce8aebaSBarry Smith .seealso: `DMField`, `PetscQuadrature`, `IS`, `PetscFEGeom`, `DMFieldEvaluateFE()`, `DMFieldCreateDefaulteQuadrature()`, `DMFieldGetDegree()` 41406dd6b0eSSatish Balay @*/ 415d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 416d71ae5a4SJacob Faibussowitsch { 4172bb1671aSToby Isaac PetscInt dim, dE; 4182bb1671aSToby Isaac PetscInt nPoints; 419b7260050SToby Isaac PetscInt maxDegree; 4202bb1671aSToby Isaac PetscFEGeom *g; 4212bb1671aSToby Isaac 4222bb1671aSToby Isaac PetscFunctionBegin; 4232bb1671aSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 4242bb1671aSToby Isaac PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2); 4252bb1671aSToby Isaac PetscValidHeader(quad, 3); 4269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &nPoints)); 4272bb1671aSToby Isaac dE = field->numComponents; 4289566063dSJacob Faibussowitsch PetscCall(PetscFEGeomCreate(quad, nPoints, dE, faceData, &g)); 4299566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field, pointIS, quad, PETSC_REAL, g->v, g->J, NULL)); 4302bb1671aSToby Isaac dim = g->dim; 4312bb1671aSToby Isaac if (dE > dim) { 4322bb1671aSToby Isaac /* space out J and make square Jacobians */ 4332bb1671aSToby Isaac PetscInt i, j, k, N = g->numPoints * g->numCells; 4342bb1671aSToby Isaac 4352bb1671aSToby Isaac for (i = N - 1; i >= 0; i--) { 4362bb1671aSToby Isaac PetscReal J[9]; 4372bb1671aSToby Isaac 4382bb1671aSToby Isaac for (j = 0; j < dE; j++) { 439ad540459SPierre Jolivet for (k = 0; k < dim; k++) J[j * dE + k] = g->J[i * dE * dim + j * dim + k]; 4402bb1671aSToby Isaac } 4412bb1671aSToby Isaac switch (dim) { 4422bb1671aSToby Isaac case 0: 4432bb1671aSToby Isaac for (j = 0; j < dE; j++) { 444ad540459SPierre Jolivet for (k = 0; k < dE; k++) J[j * dE + k] = (j == k) ? 1. : 0.; 4452bb1671aSToby Isaac } 4462bb1671aSToby Isaac break; 4472bb1671aSToby Isaac case 1: 4482bb1671aSToby Isaac if (dE == 2) { 4492bb1671aSToby Isaac PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 4502bb1671aSToby Isaac 4512bb1671aSToby Isaac J[1] = -J[2] / norm; 4522bb1671aSToby Isaac J[3] = J[0] / norm; 4532bb1671aSToby Isaac } else { 4542bb1671aSToby Isaac PetscReal inorm = 1. / PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 4552bb1671aSToby Isaac PetscReal x = J[0] * inorm; 4562bb1671aSToby Isaac PetscReal y = J[3] * inorm; 4572bb1671aSToby Isaac PetscReal z = J[6] * inorm; 4582bb1671aSToby Isaac 4592bb1671aSToby Isaac if (x > 0.) { 4602bb1671aSToby Isaac PetscReal inv1pX = 1. / (1. + x); 4612bb1671aSToby Isaac 4629371c9d4SSatish Balay J[1] = -y; 4639371c9d4SSatish Balay J[2] = -z; 4649371c9d4SSatish Balay J[4] = 1. - y * y * inv1pX; 4659371c9d4SSatish Balay J[5] = -y * z * inv1pX; 4669371c9d4SSatish Balay J[7] = -y * z * inv1pX; 4679371c9d4SSatish Balay J[8] = 1. - z * z * inv1pX; 4682bb1671aSToby Isaac } else { 4692bb1671aSToby Isaac PetscReal inv1mX = 1. / (1. - x); 4702bb1671aSToby Isaac 4719371c9d4SSatish Balay J[1] = z; 4729371c9d4SSatish Balay J[2] = y; 4739371c9d4SSatish Balay J[4] = -y * z * inv1mX; 4749371c9d4SSatish Balay J[5] = 1. - y * y * inv1mX; 4759371c9d4SSatish Balay J[7] = 1. - z * z * inv1mX; 4769371c9d4SSatish Balay J[8] = -y * z * inv1mX; 4772bb1671aSToby Isaac } 4782bb1671aSToby Isaac } 4792bb1671aSToby Isaac break; 4809371c9d4SSatish Balay case 2: { 4812bb1671aSToby Isaac PetscReal inorm; 4822bb1671aSToby Isaac 4832bb1671aSToby Isaac J[2] = J[3] * J[7] - J[6] * J[4]; 4842bb1671aSToby Isaac J[5] = J[6] * J[1] - J[0] * J[7]; 4852bb1671aSToby Isaac J[8] = J[0] * J[4] - J[3] * J[1]; 4862bb1671aSToby Isaac 4872bb1671aSToby Isaac inorm = 1. / PetscSqrtReal(J[2] * J[2] + J[5] * J[5] + J[8] * J[8]); 4882bb1671aSToby Isaac 4892bb1671aSToby Isaac J[2] *= inorm; 4902bb1671aSToby Isaac J[5] *= inorm; 4912bb1671aSToby Isaac J[8] *= inorm; 4929371c9d4SSatish Balay } break; 4932bb1671aSToby Isaac } 494ad540459SPierre Jolivet for (j = 0; j < dE * dE; j++) g->J[i * dE * dE + j] = J[j]; 4952bb1671aSToby Isaac } 4962bb1671aSToby Isaac } 4979566063dSJacob Faibussowitsch PetscCall(PetscFEGeomComplete(g)); 4989566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(field, pointIS, NULL, &maxDegree)); 499b7260050SToby Isaac g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE; 50048a46eb9SPierre Jolivet if (faceData) PetscCall((*field->ops->computeFaceData)(field, pointIS, quad, g)); 5012bb1671aSToby Isaac *geom = g; 5022bb1671aSToby Isaac PetscFunctionReturn(0); 5032bb1671aSToby Isaac } 504