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); 134f572ea9SToby Isaac PetscAssertPointer(field, 4); 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; 223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233da551e6SToby Isaac } 243da551e6SToby Isaac 253da551e6SToby Isaac /*@ 26dce8aebaSBarry Smith DMFieldDestroy - destroy a `DMField` 273da551e6SToby Isaac 283da551e6SToby Isaac Collective 293da551e6SToby Isaac 304165533cSJose E. Roman Input Parameter: 31dce8aebaSBarry Smith . field - address of `DMField` 323da551e6SToby Isaac 333da551e6SToby Isaac Level: advanced 343da551e6SToby Isaac 35dce8aebaSBarry Smith .seealso: `DMField`, `DMFieldCreate()` 363da551e6SToby Isaac @*/ 37d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldDestroy(DMField *field) 38d71ae5a4SJacob Faibussowitsch { 393da551e6SToby Isaac PetscFunctionBegin; 403ba16761SJacob Faibussowitsch if (!*field) PetscFunctionReturn(PETSC_SUCCESS); 41f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*field, DMFIELD_CLASSID, 1); 42f4f49eeaSPierre Jolivet if (--((PetscObject)*field)->refct > 0) { 439371c9d4SSatish Balay *field = NULL; 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 459371c9d4SSatish Balay } 46f4f49eeaSPierre Jolivet PetscTryTypeMethod(*field, destroy); 479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&((*field)->dm))); 489566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(field)); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 503da551e6SToby Isaac } 513da551e6SToby Isaac 52*ffeef943SBarry Smith /*@ 5320f4b53cSBarry Smith DMFieldView - view a `DMField` 543da551e6SToby Isaac 553da551e6SToby Isaac Collective 563da551e6SToby Isaac 574165533cSJose E. Roman Input Parameters: 5820f4b53cSBarry Smith + field - `DMField` 5920f4b53cSBarry Smith - viewer - viewer to display field, for example `PETSC_VIEWER_STDOUT_WORLD` 603da551e6SToby Isaac 613da551e6SToby Isaac Level: advanced 623da551e6SToby Isaac 63dce8aebaSBarry 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)); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 873da551e6SToby Isaac } 883da551e6SToby Isaac 89cc4c1da9SBarry Smith /*@ 90dce8aebaSBarry Smith DMFieldSetType - set the `DMField` implementation 913da551e6SToby Isaac 9220f4b53cSBarry Smith Collective 933da551e6SToby Isaac 943da551e6SToby Isaac Input Parameters: 95dce8aebaSBarry Smith + field - the `DMField` context 963da551e6SToby Isaac - type - a known method 973da551e6SToby Isaac 982fe279fdSBarry Smith Level: advanced 992fe279fdSBarry Smith 1003da551e6SToby Isaac Notes: 1013da551e6SToby Isaac See "include/petscvec.h" for available methods (for instance) 102dce8aebaSBarry Smith + `DMFIELDDA` - a field defined only by its values at the corners of a `DMDA` 103dce8aebaSBarry Smith . `DMFIELDDS` - a field defined by a discretization over a mesh set with `DMSetField()` 104dce8aebaSBarry Smith - `DMFIELDSHELL` - a field defined by arbitrary callbacks 1053da551e6SToby Isaac 106dce8aebaSBarry 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); 1154f572ea9SToby Isaac PetscAssertPointer(type, 2); 1163da551e6SToby Isaac 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)field, type, &match)); 1183ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1193da551e6SToby Isaac 1209566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(DMFieldList, type, &r)); 1216adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)field), 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)); 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1303da551e6SToby Isaac } 1313da551e6SToby Isaac 132cc4c1da9SBarry Smith /*@ 133dce8aebaSBarry Smith DMFieldGetType - Gets the `DMFieldType` name (as a string) from the `DMField`. 1343da551e6SToby Isaac 1353da551e6SToby Isaac Not Collective 1363da551e6SToby Isaac 1373da551e6SToby Isaac Input Parameter: 138dce8aebaSBarry Smith . field - The `DMField` context 1393da551e6SToby Isaac 1403da551e6SToby Isaac Output Parameter: 141dce8aebaSBarry Smith . type - The `DMFieldType` name 1423da551e6SToby Isaac 1433da551e6SToby Isaac Level: advanced 1443da551e6SToby Isaac 145dce8aebaSBarry 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); 1514f572ea9SToby Isaac PetscAssertPointer(type, 2); 1529566063dSJacob Faibussowitsch PetscCall(DMFieldRegisterAll()); 1533da551e6SToby Isaac *type = ((PetscObject)field)->type_name; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1553da551e6SToby Isaac } 1563da551e6SToby Isaac 157da311338SToby Isaac /*@ 158da311338SToby Isaac DMFieldGetNumComponents - Returns the number of components in the field 159da311338SToby Isaac 16020f4b53cSBarry Smith Not Collective 161da311338SToby Isaac 162da311338SToby Isaac Input Parameter: 163dce8aebaSBarry 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 170dce8aebaSBarry 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); 1764f572ea9SToby Isaac PetscAssertPointer(nc, 2); 1773da551e6SToby Isaac *nc = field->numComponents; 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1793da551e6SToby Isaac } 1803da551e6SToby Isaac 181da311338SToby Isaac /*@ 182dce8aebaSBarry Smith DMFieldGetDM - Returns the `DM` for the manifold over which the field is defined. 183da311338SToby Isaac 18420f4b53cSBarry Smith Not Collective 185da311338SToby Isaac 186da311338SToby Isaac Input Parameter: 187dce8aebaSBarry Smith . field - The `DMField` object 188da311338SToby Isaac 189da311338SToby Isaac Output Parameter: 190dce8aebaSBarry Smith . dm - The `DM` object 191da311338SToby Isaac 192da311338SToby Isaac Level: intermediate 193da311338SToby Isaac 194dce8aebaSBarry 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); 2004f572ea9SToby Isaac PetscAssertPointer(dm, 2); 2013da551e6SToby Isaac *dm = field->dm; 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2033da551e6SToby Isaac } 2043da551e6SToby Isaac 205da311338SToby Isaac /*@ 206da311338SToby Isaac DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 207da311338SToby Isaac 20820f4b53cSBarry Smith Collective 209da311338SToby Isaac 210d8d19677SJose E. Roman Input Parameters: 211dce8aebaSBarry 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 215dce8aebaSBarry Smith - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. 216dce8aebaSBarry 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 232dce8aebaSBarry 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); 2394f572ea9SToby Isaac if (B) PetscAssertPointer(B, 4); 2404f572ea9SToby Isaac if (D) PetscAssertPointer(D, 5); 2414f572ea9SToby Isaac if (H) PetscAssertPointer(H, 6); 242213acdd3SPierre Jolivet PetscUseTypeMethod(field, evaluate, points, datatype, B, D, H); 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2443da551e6SToby Isaac } 2453da551e6SToby Isaac 246da311338SToby Isaac /*@ 247da311338SToby Isaac DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 248da311338SToby Isaac quadrature points on a reference point. The derivatives are taken with respect to the 249da311338SToby Isaac reference coordinates. 250da311338SToby Isaac 25120f4b53cSBarry Smith Not Collective 252da311338SToby Isaac 253d8d19677SJose E. Roman Input Parameters: 254dce8aebaSBarry Smith + field - The `DMField` object 255da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 256da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field. 25720f4b53cSBarry Smith - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. 25820f4b53cSBarry Smith If the field is complex and datatype is `PETSC_REAL`, the real part of the 259da311338SToby Isaac field is returned. 260da311338SToby Isaac 261d8d19677SJose E. Roman Output Parameters: 262da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 26320f4b53cSBarry Smith If B is not `NULL`, the values of the field are written in this array, varying first by component, 264da311338SToby Isaac then by point. 265da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 26620f4b53cSBarry Smith If D is not `NULL`, the values of the field's spatial derivatives are written in this array, 267da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 268da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 26920f4b53cSBarry Smith If H is not `NULL`, the values of the field's second spatial derivatives are written in this array, 270da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 271da311338SToby Isaac 272da311338SToby Isaac Level: intermediate 273da311338SToby Isaac 274dce8aebaSBarry Smith .seealso: `DMField`, `DM`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFV()` 275da311338SToby Isaac @*/ 276d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 277d71ae5a4SJacob Faibussowitsch { 2783da551e6SToby Isaac PetscFunctionBegin; 2793da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 2803da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 2813da551e6SToby Isaac PetscValidHeader(points, 3); 2824f572ea9SToby Isaac if (B) PetscAssertPointer(B, 5); 2834f572ea9SToby Isaac if (D) PetscAssertPointer(D, 6); 2844f572ea9SToby Isaac if (H) PetscAssertPointer(H, 7); 285213acdd3SPierre Jolivet PetscUseTypeMethod(field, evaluateFE, cellIS, points, datatype, B, D, H); 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2873da551e6SToby Isaac } 2883da551e6SToby Isaac 289da311338SToby Isaac /*@ 290da311338SToby Isaac DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 291da311338SToby Isaac 29220f4b53cSBarry Smith Not Collective 293da311338SToby Isaac 294d8d19677SJose E. Roman Input Parameters: 295dce8aebaSBarry Smith + field - The `DMField` object 296da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 297dce8aebaSBarry Smith - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. 298dce8aebaSBarry Smith If the field is complex and datatype is `PETSC_REAL`, the real part of the 299da311338SToby Isaac field is returned. 300da311338SToby Isaac 301d8d19677SJose E. Roman Output Parameters: 302da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 30320f4b53cSBarry Smith If B is not `NULL`, the values of the field are written in this array, varying first by component, 304da311338SToby Isaac then by point. 305da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 30620f4b53cSBarry Smith If D is not `NULL`, the values of the field's spatial derivatives are written in this array, 307da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 308da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 30920f4b53cSBarry Smith If H is not `NULL`, the values of the field's second spatial derivatives are written in this array, 310da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 311da311338SToby Isaac 312da311338SToby Isaac Level: intermediate 313da311338SToby Isaac 314dce8aebaSBarry Smith .seealso: `DMField`, `IS`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFE()`, `PetscDataType` 315da311338SToby Isaac @*/ 316d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 317d71ae5a4SJacob Faibussowitsch { 3183da551e6SToby Isaac PetscFunctionBegin; 3193da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3203da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 3214f572ea9SToby Isaac if (B) PetscAssertPointer(B, 4); 3224f572ea9SToby Isaac if (D) PetscAssertPointer(D, 5); 3234f572ea9SToby Isaac if (H) PetscAssertPointer(H, 6); 324213acdd3SPierre Jolivet PetscUseTypeMethod(field, evaluateFV, cellIS, datatype, B, D, H); 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3263da551e6SToby Isaac } 3273da551e6SToby Isaac 328da311338SToby Isaac /*@ 329b7260050SToby Isaac DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the 330da311338SToby Isaac reference element 331da311338SToby Isaac 33220f4b53cSBarry Smith Not Collective 333da311338SToby Isaac 3344165533cSJose E. Roman Input Parameters: 335dce8aebaSBarry Smith + field - the `DMField` object 336da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance 337da311338SToby Isaac 3384165533cSJose E. Roman Output Parameters: 339b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element 340b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element 341da311338SToby Isaac 342da311338SToby Isaac Level: intermediate 343da311338SToby Isaac 344dce8aebaSBarry Smith .seealso: `DMField`, `IS`, `DMFieldEvaluateFE()` 345da311338SToby Isaac @*/ 346d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) 347d71ae5a4SJacob Faibussowitsch { 3483da551e6SToby Isaac PetscFunctionBegin; 3493da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3503da551e6SToby Isaac PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2); 3514f572ea9SToby Isaac if (minDegree) PetscAssertPointer(minDegree, 3); 3524f572ea9SToby Isaac if (maxDegree) PetscAssertPointer(maxDegree, 4); 3533da551e6SToby Isaac 354b7260050SToby Isaac if (minDegree) *minDegree = -1; 355b7260050SToby Isaac if (maxDegree) *maxDegree = PETSC_MAX_INT; 3563da551e6SToby Isaac 357213acdd3SPierre Jolivet PetscTryTypeMethod(field, getDegree, cellIS, minDegree, maxDegree); 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3593da551e6SToby Isaac } 3603da551e6SToby Isaac 361da311338SToby Isaac /*@ 362da311338SToby Isaac DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 363da311338SToby Isaac points via pullback onto the reference element 364da311338SToby Isaac 36520f4b53cSBarry Smith Not Collective 366da311338SToby Isaac 3674165533cSJose E. Roman Input Parameters: 368dce8aebaSBarry Smith + field - the `DMField` object 369da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field 370da311338SToby Isaac 3714165533cSJose E. Roman Output Parameter: 372dce8aebaSBarry Smith . quad - a `PetscQuadrature` object 373da311338SToby Isaac 374da311338SToby Isaac Level: developer 375da311338SToby Isaac 376dce8aebaSBarry Smith .seealso: `DMField`, `PetscQuadrature`, `IS`, `DMFieldEvaluteFE()`, `DMFieldGetDegree()` 377da311338SToby Isaac @*/ 378d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 379d71ae5a4SJacob Faibussowitsch { 3803da551e6SToby Isaac PetscFunctionBegin; 3813da551e6SToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 3823da551e6SToby Isaac PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2); 3834f572ea9SToby Isaac PetscAssertPointer(quad, 3); 3843da551e6SToby Isaac 3853da551e6SToby Isaac *quad = NULL; 386dbbe0bcdSBarry Smith PetscTryTypeMethod(field, createDefaultQuadrature, pointIS, quad); 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3883da551e6SToby Isaac } 3893da551e6SToby Isaac 390da311338SToby Isaac /*@C 391da311338SToby Isaac DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 392da311338SToby Isaac 39320f4b53cSBarry Smith Not Collective 394da311338SToby Isaac 3954165533cSJose E. Roman Input Parameters: 396dce8aebaSBarry Smith + field - the `DMField` object 397da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field 398da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors 399da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 400da311338SToby Isaac be calculated 401da311338SToby Isaac 4024165533cSJose E. Roman Output Parameter: 403da311338SToby Isaac . geom - the geometric factors 404da311338SToby Isaac 405da311338SToby Isaac Level: developer 406da311338SToby Isaac 407dce8aebaSBarry Smith .seealso: `DMField`, `PetscQuadrature`, `IS`, `PetscFEGeom`, `DMFieldEvaluateFE()`, `DMFieldCreateDefaulteQuadrature()`, `DMFieldGetDegree()` 40806dd6b0eSSatish Balay @*/ 409d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 410d71ae5a4SJacob Faibussowitsch { 4112bb1671aSToby Isaac PetscInt dim, dE; 4122bb1671aSToby Isaac PetscInt nPoints; 413b7260050SToby Isaac PetscInt maxDegree; 4142bb1671aSToby Isaac PetscFEGeom *g; 4152bb1671aSToby Isaac 4162bb1671aSToby Isaac PetscFunctionBegin; 4172bb1671aSToby Isaac PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1); 4182bb1671aSToby Isaac PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2); 4192bb1671aSToby Isaac PetscValidHeader(quad, 3); 4209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &nPoints)); 4212bb1671aSToby Isaac dE = field->numComponents; 4229566063dSJacob Faibussowitsch PetscCall(PetscFEGeomCreate(quad, nPoints, dE, faceData, &g)); 4239566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field, pointIS, quad, PETSC_REAL, g->v, g->J, NULL)); 4242bb1671aSToby Isaac dim = g->dim; 4252bb1671aSToby Isaac if (dE > dim) { 4262bb1671aSToby Isaac /* space out J and make square Jacobians */ 4272bb1671aSToby Isaac PetscInt i, j, k, N = g->numPoints * g->numCells; 4282bb1671aSToby Isaac 4292bb1671aSToby Isaac for (i = N - 1; i >= 0; i--) { 4309b15cf9aSJacob Faibussowitsch PetscReal J[9] = {0}; 4312bb1671aSToby Isaac 4322bb1671aSToby Isaac for (j = 0; j < dE; j++) { 433ad540459SPierre Jolivet for (k = 0; k < dim; k++) J[j * dE + k] = g->J[i * dE * dim + j * dim + k]; 4342bb1671aSToby Isaac } 4352bb1671aSToby Isaac switch (dim) { 4362bb1671aSToby Isaac case 0: 4372bb1671aSToby Isaac for (j = 0; j < dE; j++) { 438ad540459SPierre Jolivet for (k = 0; k < dE; k++) J[j * dE + k] = (j == k) ? 1. : 0.; 4392bb1671aSToby Isaac } 4402bb1671aSToby Isaac break; 4412bb1671aSToby Isaac case 1: 4422bb1671aSToby Isaac if (dE == 2) { 4432bb1671aSToby Isaac PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 4442bb1671aSToby Isaac 4452bb1671aSToby Isaac J[1] = -J[2] / norm; 4462bb1671aSToby Isaac J[3] = J[0] / norm; 4472bb1671aSToby Isaac } else { 4482bb1671aSToby Isaac PetscReal inorm = 1. / PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 4492bb1671aSToby Isaac PetscReal x = J[0] * inorm; 4502bb1671aSToby Isaac PetscReal y = J[3] * inorm; 4512bb1671aSToby Isaac PetscReal z = J[6] * inorm; 4522bb1671aSToby Isaac 4532bb1671aSToby Isaac if (x > 0.) { 4542bb1671aSToby Isaac PetscReal inv1pX = 1. / (1. + x); 4552bb1671aSToby Isaac 4569371c9d4SSatish Balay J[1] = -y; 4579371c9d4SSatish Balay J[2] = -z; 4589371c9d4SSatish Balay J[4] = 1. - y * y * inv1pX; 4599371c9d4SSatish Balay J[5] = -y * z * inv1pX; 4609371c9d4SSatish Balay J[7] = -y * z * inv1pX; 4619371c9d4SSatish Balay J[8] = 1. - z * z * inv1pX; 4622bb1671aSToby Isaac } else { 4632bb1671aSToby Isaac PetscReal inv1mX = 1. / (1. - x); 4642bb1671aSToby Isaac 4659371c9d4SSatish Balay J[1] = z; 4669371c9d4SSatish Balay J[2] = y; 4679371c9d4SSatish Balay J[4] = -y * z * inv1mX; 4689371c9d4SSatish Balay J[5] = 1. - y * y * inv1mX; 4699371c9d4SSatish Balay J[7] = 1. - z * z * inv1mX; 4709371c9d4SSatish Balay J[8] = -y * z * inv1mX; 4712bb1671aSToby Isaac } 4722bb1671aSToby Isaac } 4732bb1671aSToby Isaac break; 4749371c9d4SSatish Balay case 2: { 4752bb1671aSToby Isaac PetscReal inorm; 4762bb1671aSToby Isaac 4772bb1671aSToby Isaac J[2] = J[3] * J[7] - J[6] * J[4]; 4782bb1671aSToby Isaac J[5] = J[6] * J[1] - J[0] * J[7]; 4792bb1671aSToby Isaac J[8] = J[0] * J[4] - J[3] * J[1]; 4802bb1671aSToby Isaac 4812bb1671aSToby Isaac inorm = 1. / PetscSqrtReal(J[2] * J[2] + J[5] * J[5] + J[8] * J[8]); 4822bb1671aSToby Isaac 4832bb1671aSToby Isaac J[2] *= inorm; 4842bb1671aSToby Isaac J[5] *= inorm; 4852bb1671aSToby Isaac J[8] *= inorm; 4869371c9d4SSatish Balay } break; 4872bb1671aSToby Isaac } 488ad540459SPierre Jolivet for (j = 0; j < dE * dE; j++) g->J[i * dE * dE + j] = J[j]; 4892bb1671aSToby Isaac } 4902bb1671aSToby Isaac } 4919566063dSJacob Faibussowitsch PetscCall(PetscFEGeomComplete(g)); 4929566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(field, pointIS, NULL, &maxDegree)); 493b7260050SToby Isaac g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE; 4949927e4dfSBarry Smith if (faceData) PetscUseTypeMethod(field, computeFaceData, pointIS, quad, g); 4952bb1671aSToby Isaac *geom = g; 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4972bb1671aSToby Isaac } 498