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 53da551e6SToby Isaac const char *const DMFieldContinuities[] = { 63da551e6SToby Isaac "VERTEX", 73da551e6SToby Isaac "EDGE", 83da551e6SToby Isaac "FACET", 93da551e6SToby Isaac "CELL", 10ea78f98cSLisandro Dalcin NULL 113da551e6SToby Isaac }; 123da551e6SToby Isaac 133da551e6SToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField *field) 143da551e6SToby Isaac { 153da551e6SToby Isaac DMField b; 163da551e6SToby Isaac 173da551e6SToby Isaac PetscFunctionBegin; 183da551e6SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 193da551e6SToby Isaac PetscValidPointer(field,2); 209566063dSJacob Faibussowitsch PetscCall(DMFieldInitializePackage()); 213da551e6SToby Isaac 229566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 243da551e6SToby Isaac b->dm = dm; 253da551e6SToby Isaac b->continuity = continuity; 263da551e6SToby Isaac b->numComponents = numComponents; 273da551e6SToby Isaac *field = b; 283da551e6SToby Isaac PetscFunctionReturn(0); 293da551e6SToby Isaac } 303da551e6SToby Isaac 313da551e6SToby Isaac /*@ 323da551e6SToby Isaac DMFieldDestroy - destroy a DMField 333da551e6SToby Isaac 343da551e6SToby Isaac Collective 353da551e6SToby Isaac 364165533cSJose E. Roman Input Parameter: 373da551e6SToby Isaac . field - address of DMField 383da551e6SToby Isaac 393da551e6SToby Isaac Level: advanced 403da551e6SToby Isaac 41db781477SPatrick Sanan .seealso: `DMFieldCreate()` 423da551e6SToby Isaac @*/ 433da551e6SToby Isaac PetscErrorCode DMFieldDestroy(DMField *field) 443da551e6SToby Isaac { 453da551e6SToby Isaac PetscFunctionBegin; 463da551e6SToby Isaac if (!*field) PetscFunctionReturn(0); 473da551e6SToby Isaac PetscValidHeaderSpecific((*field),DMFIELD_CLASSID,1); 48ea78f98cSLisandro Dalcin if (--((PetscObject)(*field))->refct > 0) {*field = NULL; PetscFunctionReturn(0);} 49*dbbe0bcdSBarry Smith PetscTryTypeMethod((*field),destroy); 509566063dSJacob Faibussowitsch PetscCall(DMDestroy(&((*field)->dm))); 519566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(field)); 523da551e6SToby Isaac PetscFunctionReturn(0); 533da551e6SToby Isaac } 543da551e6SToby Isaac 553da551e6SToby Isaac /*@C 563da551e6SToby Isaac DMFieldView - view a DMField 573da551e6SToby Isaac 583da551e6SToby Isaac Collective 593da551e6SToby Isaac 604165533cSJose E. Roman Input Parameters: 613da551e6SToby Isaac + field - DMField 623da551e6SToby Isaac - viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD 633da551e6SToby Isaac 643da551e6SToby Isaac Level: advanced 653da551e6SToby Isaac 66db781477SPatrick Sanan .seealso: `DMFieldCreate()` 673da551e6SToby Isaac @*/ 683da551e6SToby Isaac PetscErrorCode DMFieldView(DMField field,PetscViewer viewer) 693da551e6SToby Isaac { 703da551e6SToby Isaac PetscBool iascii; 713da551e6SToby Isaac 723da551e6SToby Isaac PetscFunctionBegin; 733da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 749566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer)); 753da551e6SToby Isaac PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 763da551e6SToby Isaac PetscCheckSameComm(field,1,viewer,2); 779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 783da551e6SToby Isaac if (iascii) { 799566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " components\n",field->numComponents)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity])); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT)); 849566063dSJacob Faibussowitsch PetscCall(DMView(field->dm,viewer)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 863da551e6SToby Isaac } 87*dbbe0bcdSBarry Smith PetscTryTypeMethod(field,view,viewer); 881baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPopTab(viewer)); 893da551e6SToby Isaac PetscFunctionReturn(0); 903da551e6SToby Isaac } 913da551e6SToby Isaac 923da551e6SToby Isaac /*@C 933da551e6SToby Isaac DMFieldSetType - set the DMField implementation 943da551e6SToby Isaac 95d083f849SBarry Smith Collective on field 963da551e6SToby Isaac 973da551e6SToby Isaac Input Parameters: 983da551e6SToby Isaac + field - the DMField context 993da551e6SToby Isaac - type - a known method 1003da551e6SToby Isaac 1013da551e6SToby Isaac Notes: 1023da551e6SToby Isaac See "include/petscvec.h" for available methods (for instance) 1033da551e6SToby Isaac + DMFIELDDA - a field defined only by its values at the corners of a DMDA 1043da551e6SToby Isaac . DMFIELDDS - a field defined by a discretization over a mesh set with DMSetField() 1053da551e6SToby Isaac - DMFIELDSHELL - a field defined by arbitrary callbacks 1063da551e6SToby Isaac 1073da551e6SToby Isaac Level: advanced 1083da551e6SToby Isaac 109db781477SPatrick Sanan .seealso: `DMFieldType`, 1103da551e6SToby Isaac @*/ 1113da551e6SToby Isaac PetscErrorCode DMFieldSetType(DMField field,DMFieldType type) 1123da551e6SToby Isaac { 1133da551e6SToby Isaac PetscBool match; 1145f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(DMField); 1153da551e6SToby Isaac 1163da551e6SToby Isaac PetscFunctionBegin; 1173da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1183da551e6SToby Isaac PetscValidCharPointer(type,2); 1193da551e6SToby Isaac 1209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)field,type,&match)); 1213da551e6SToby Isaac if (match) PetscFunctionReturn(0); 1223da551e6SToby Isaac 1239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(DMFieldList,type,&r)); 1245f80ce2aSJacob Faibussowitsch PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type); 1253da551e6SToby Isaac /* Destroy the previous private DMField context */ 126*dbbe0bcdSBarry Smith PetscTryTypeMethod(field,destroy); 1275f80ce2aSJacob Faibussowitsch 1289566063dSJacob Faibussowitsch PetscCall(PetscMemzero(field->ops,sizeof(*field->ops))); 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)field,type)); 1303da551e6SToby Isaac field->ops->create = r; 1319566063dSJacob Faibussowitsch PetscCall((*r)(field)); 1323da551e6SToby Isaac PetscFunctionReturn(0); 1333da551e6SToby Isaac } 1343da551e6SToby Isaac 1353da551e6SToby Isaac /*@C 1363da551e6SToby Isaac DMFieldGetType - Gets the DMField type name (as a string) from the DMField. 1373da551e6SToby Isaac 1383da551e6SToby Isaac Not Collective 1393da551e6SToby Isaac 1403da551e6SToby Isaac Input Parameter: 1413da551e6SToby Isaac . field - The DMField context 1423da551e6SToby Isaac 1433da551e6SToby Isaac Output Parameter: 1443da551e6SToby Isaac . type - The DMField type name 1453da551e6SToby Isaac 1463da551e6SToby Isaac Level: advanced 1473da551e6SToby Isaac 148db781477SPatrick Sanan .seealso: `DMFieldSetType()` 1493da551e6SToby Isaac @*/ 1503da551e6SToby Isaac PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 1513da551e6SToby Isaac { 1523da551e6SToby Isaac PetscFunctionBegin; 153064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(field, DMFIELD_CLASSID,1); 1543da551e6SToby Isaac PetscValidPointer(type,2); 1559566063dSJacob Faibussowitsch PetscCall(DMFieldRegisterAll()); 1563da551e6SToby Isaac *type = ((PetscObject)field)->type_name; 1573da551e6SToby Isaac PetscFunctionReturn(0); 1583da551e6SToby Isaac } 1593da551e6SToby Isaac 160da311338SToby Isaac /*@ 161da311338SToby Isaac DMFieldGetNumComponents - Returns the number of components in the field 162da311338SToby Isaac 163da311338SToby Isaac Not collective 164da311338SToby Isaac 165da311338SToby Isaac Input Parameter: 166da311338SToby Isaac . field - The DMField object 167da311338SToby Isaac 168da311338SToby Isaac Output Parameter: 169da311338SToby Isaac . nc - The number of field components 170da311338SToby Isaac 171da311338SToby Isaac Level: intermediate 172da311338SToby Isaac 173db781477SPatrick Sanan .seealso: `DMFieldEvaluate()` 174da311338SToby Isaac @*/ 1753da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 1763da551e6SToby Isaac { 1773da551e6SToby Isaac PetscFunctionBegin; 1783da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1793da551e6SToby Isaac PetscValidIntPointer(nc,2); 1803da551e6SToby Isaac *nc = field->numComponents; 1813da551e6SToby Isaac PetscFunctionReturn(0); 1823da551e6SToby Isaac } 1833da551e6SToby Isaac 184da311338SToby Isaac /*@ 185da311338SToby Isaac DMFieldGetDM - Returns the DM for the manifold over which the field is defined. 186da311338SToby Isaac 187da311338SToby Isaac Not collective 188da311338SToby Isaac 189da311338SToby Isaac Input Parameter: 190da311338SToby Isaac . field - The DMField object 191da311338SToby Isaac 192da311338SToby Isaac Output Parameter: 193da311338SToby Isaac . dm - The DM object 194da311338SToby Isaac 195da311338SToby Isaac Level: intermediate 196da311338SToby Isaac 197db781477SPatrick Sanan .seealso: `DMFieldEvaluate()` 198da311338SToby Isaac @*/ 1993da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 2003da551e6SToby Isaac { 2013da551e6SToby Isaac PetscFunctionBegin; 2023da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2033da551e6SToby Isaac PetscValidPointer(dm,2); 2043da551e6SToby Isaac *dm = field->dm; 2053da551e6SToby Isaac PetscFunctionReturn(0); 2063da551e6SToby Isaac } 2073da551e6SToby Isaac 208da311338SToby Isaac /*@ 209da311338SToby Isaac DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 210da311338SToby Isaac 211d083f849SBarry Smith Collective on points 212da311338SToby Isaac 213d8d19677SJose E. Roman Input Parameters: 214da311338SToby Isaac + field - The DMField object 215da311338SToby Isaac . points - The points at which to evaluate the field. Should have size d x n, 216da311338SToby Isaac where d is the coordinate dimension of the manifold and n is the number 217da311338SToby Isaac of points 218da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 219da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 220da311338SToby Isaac field is returned. 221da311338SToby Isaac 222d8d19677SJose E. Roman Output Parameters: 223da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 224da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 225da311338SToby Isaac then by point. 226da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 227da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 228da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 229da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 230da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 231da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 232da311338SToby Isaac 233da311338SToby Isaac Level: intermediate 234da311338SToby Isaac 235db781477SPatrick Sanan .seealso: `DMFieldGetDM()`, `DMFieldGetNumComponents()`, `DMFieldEvaluateFE()`, `DMFieldEvaluateFV()` 236da311338SToby Isaac @*/ 2373da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 2383da551e6SToby Isaac { 2393da551e6SToby Isaac PetscFunctionBegin; 2403da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2413da551e6SToby Isaac PetscValidHeaderSpecific(points,VEC_CLASSID,2); 242064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,4); 243064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,5); 244064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,6); 2453da551e6SToby Isaac if (field->ops->evaluate) { 2469566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluate) (field, points, datatype, B, D, H)); 2473da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 2483da551e6SToby Isaac PetscFunctionReturn(0); 2493da551e6SToby Isaac } 2503da551e6SToby Isaac 251da311338SToby Isaac /*@ 252da311338SToby Isaac DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 253da311338SToby Isaac quadrature points on a reference point. The derivatives are taken with respect to the 254da311338SToby Isaac reference coordinates. 255da311338SToby Isaac 256da311338SToby Isaac Not collective 257da311338SToby Isaac 258d8d19677SJose E. Roman Input Parameters: 259da311338SToby Isaac + field - The DMField object 260da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 261da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field. 262da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 263da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 264da311338SToby Isaac field is returned. 265da311338SToby Isaac 266d8d19677SJose E. Roman Output Parameters: 267da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 268da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 269da311338SToby Isaac then by point. 270da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 271da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 272da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 273da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 274da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 275da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 276da311338SToby Isaac 277da311338SToby Isaac Level: intermediate 278da311338SToby Isaac 279db781477SPatrick Sanan .seealso: `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFV()` 280da311338SToby Isaac @*/ 2813da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 2823da551e6SToby Isaac { 2833da551e6SToby Isaac PetscFunctionBegin; 2843da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2853da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 2863da551e6SToby Isaac PetscValidHeader(points,3); 287064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,5); 288064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,6); 289064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,7); 2903da551e6SToby Isaac if (field->ops->evaluateFE) { 2919566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H)); 2923da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 2933da551e6SToby Isaac PetscFunctionReturn(0); 2943da551e6SToby Isaac } 2953da551e6SToby Isaac 296da311338SToby Isaac /*@ 297da311338SToby Isaac DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 298da311338SToby Isaac 299da311338SToby Isaac Not collective 300da311338SToby Isaac 301d8d19677SJose E. Roman Input Parameters: 302da311338SToby Isaac + field - The DMField object 303da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 304da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 305da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 306da311338SToby Isaac field is returned. 307da311338SToby Isaac 308d8d19677SJose E. Roman Output Parameters: 309da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 310da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 311da311338SToby Isaac then by point. 312da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 313da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 314da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 315da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 316da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 317da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 318da311338SToby Isaac 319da311338SToby Isaac Level: intermediate 320da311338SToby Isaac 321db781477SPatrick Sanan .seealso: `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFE()` 322da311338SToby Isaac @*/ 3233da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 3243da551e6SToby Isaac { 3253da551e6SToby Isaac PetscFunctionBegin; 3263da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3273da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 328064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,4); 329064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,5); 330064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,6); 3313da551e6SToby Isaac if (field->ops->evaluateFV) { 3329566063dSJacob Faibussowitsch PetscCall((*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H)); 3333da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 3343da551e6SToby Isaac PetscFunctionReturn(0); 3353da551e6SToby Isaac } 3363da551e6SToby Isaac 337da311338SToby Isaac /*@ 338b7260050SToby Isaac DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the 339da311338SToby Isaac reference element 340da311338SToby Isaac 341da311338SToby Isaac Not collective 342da311338SToby Isaac 3434165533cSJose E. Roman Input Parameters: 344da311338SToby Isaac + field - the DMField object 345da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance 346da311338SToby Isaac 3474165533cSJose E. Roman Output Parameters: 348b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element 349b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element 350da311338SToby Isaac 351da311338SToby Isaac Level: intermediate 352da311338SToby Isaac 353db781477SPatrick Sanan .seealso: `DMFieldEvaluateFE()` 354da311338SToby Isaac @*/ 355b7260050SToby Isaac PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) 3563da551e6SToby Isaac { 3573da551e6SToby Isaac PetscFunctionBegin; 3583da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3593da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 360dadcf809SJacob Faibussowitsch if (minDegree) PetscValidIntPointer(minDegree,3); 361dadcf809SJacob Faibussowitsch if (maxDegree) PetscValidIntPointer(maxDegree,4); 3623da551e6SToby Isaac 363b7260050SToby Isaac if (minDegree) *minDegree = -1; 364b7260050SToby Isaac if (maxDegree) *maxDegree = PETSC_MAX_INT; 3653da551e6SToby Isaac 3661baa6e33SBarry Smith if (field->ops->getDegree) PetscCall((*field->ops->getDegree) (field,cellIS,minDegree,maxDegree)); 3673da551e6SToby Isaac PetscFunctionReturn(0); 3683da551e6SToby Isaac } 3693da551e6SToby Isaac 370da311338SToby Isaac /*@ 371da311338SToby Isaac DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 372da311338SToby Isaac points via pullback onto the reference element 373da311338SToby Isaac 374da311338SToby Isaac Not collective 375da311338SToby Isaac 3764165533cSJose E. Roman Input Parameters: 377da311338SToby Isaac + field - the DMField object 378da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field 379da311338SToby Isaac 3804165533cSJose E. Roman Output Parameter: 381da311338SToby Isaac . quad - a PetscQuadrature object 382da311338SToby Isaac 383da311338SToby Isaac Level: developer 384da311338SToby Isaac 385db781477SPatrick Sanan .seealso: `DMFieldEvaluteFE()`, `DMFieldGetDegree()` 386da311338SToby Isaac @*/ 3873da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 3883da551e6SToby Isaac { 3893da551e6SToby Isaac PetscFunctionBegin; 3903da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3913da551e6SToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 3923da551e6SToby Isaac PetscValidPointer(quad,3); 3933da551e6SToby Isaac 3943da551e6SToby Isaac *quad = NULL; 395*dbbe0bcdSBarry Smith PetscTryTypeMethod(field,createDefaultQuadrature, pointIS, quad); 3963da551e6SToby Isaac PetscFunctionReturn(0); 3973da551e6SToby Isaac } 3983da551e6SToby Isaac 399da311338SToby Isaac /*@C 400da311338SToby Isaac DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 401da311338SToby Isaac 402da311338SToby Isaac Not collective 403da311338SToby Isaac 4044165533cSJose E. Roman Input Parameters: 405da311338SToby Isaac + field - the DMField object 406da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field 407da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors 408da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 409da311338SToby Isaac be calculated 410da311338SToby Isaac 4114165533cSJose E. Roman Output Parameter: 412da311338SToby Isaac . geom - the geometric factors 413da311338SToby Isaac 414da311338SToby Isaac Level: developer 415da311338SToby Isaac 416db781477SPatrick Sanan .seealso: `DMFieldEvaluateFE()`, `DMFieldCreateDefaulteQuadrature()`, `DMFieldGetDegree()` 41706dd6b0eSSatish Balay @*/ 4182bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 4192bb1671aSToby Isaac { 4202bb1671aSToby Isaac PetscInt dim, dE; 4212bb1671aSToby Isaac PetscInt nPoints; 422b7260050SToby Isaac PetscInt maxDegree; 4232bb1671aSToby Isaac PetscFEGeom *g; 4242bb1671aSToby Isaac 4252bb1671aSToby Isaac PetscFunctionBegin; 4262bb1671aSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 4272bb1671aSToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 4282bb1671aSToby Isaac PetscValidHeader(quad,3); 4299566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS,&nPoints)); 4302bb1671aSToby Isaac dE = field->numComponents; 4319566063dSJacob Faibussowitsch PetscCall(PetscFEGeomCreate(quad,nPoints,dE,faceData,&g)); 4329566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL)); 4332bb1671aSToby Isaac dim = g->dim; 4342bb1671aSToby Isaac if (dE > dim) { 4352bb1671aSToby Isaac /* space out J and make square Jacobians */ 4362bb1671aSToby Isaac PetscInt i, j, k, N = g->numPoints * g->numCells; 4372bb1671aSToby Isaac 4382bb1671aSToby Isaac for (i = N-1; i >= 0; i--) { 4392bb1671aSToby Isaac PetscReal J[9]; 4402bb1671aSToby Isaac 4412bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4422bb1671aSToby Isaac for (k = 0; k < dim; k++) { 4432bb1671aSToby Isaac J[j*dE + k] = g->J[i*dE*dim + j*dim + k]; 4442bb1671aSToby Isaac } 4452bb1671aSToby Isaac } 4462bb1671aSToby Isaac switch (dim) { 4472bb1671aSToby Isaac case 0: 4482bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4492bb1671aSToby Isaac for (k = 0; k < dE; k++) { 4502bb1671aSToby Isaac J[j * dE + k] = (j == k) ? 1. : 0.; 4512bb1671aSToby Isaac } 4522bb1671aSToby Isaac } 4532bb1671aSToby Isaac break; 4542bb1671aSToby Isaac case 1: 4552bb1671aSToby Isaac if (dE == 2) { 4562bb1671aSToby Isaac PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 4572bb1671aSToby Isaac 4582bb1671aSToby Isaac J[1] = -J[2] / norm; 4592bb1671aSToby Isaac J[3] = J[0] / norm; 4602bb1671aSToby Isaac } else { 4612bb1671aSToby Isaac PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 4622bb1671aSToby Isaac PetscReal x = J[0] * inorm; 4632bb1671aSToby Isaac PetscReal y = J[3] * inorm; 4642bb1671aSToby Isaac PetscReal z = J[6] * inorm; 4652bb1671aSToby Isaac 4662bb1671aSToby Isaac if (x > 0.) { 4672bb1671aSToby Isaac PetscReal inv1pX = 1./ (1. + x); 4682bb1671aSToby Isaac 4692bb1671aSToby Isaac J[1] = -y; J[2] = -z; 4702bb1671aSToby Isaac J[4] = 1. - y*y*inv1pX; J[5] = -y*z*inv1pX; 4712bb1671aSToby Isaac J[7] = -y*z*inv1pX; J[8] = 1. - z*z*inv1pX; 4722bb1671aSToby Isaac } else { 4732bb1671aSToby Isaac PetscReal inv1mX = 1./ (1. - x); 4742bb1671aSToby Isaac 4752bb1671aSToby Isaac J[1] = z; J[2] = y; 4762bb1671aSToby Isaac J[4] = -y*z*inv1mX; J[5] = 1. - y*y*inv1mX; 4772bb1671aSToby Isaac J[7] = 1. - z*z*inv1mX; J[8] = -y*z*inv1mX; 4782bb1671aSToby Isaac } 4792bb1671aSToby Isaac } 4802bb1671aSToby Isaac break; 4812bb1671aSToby Isaac case 2: 4822bb1671aSToby Isaac { 4832bb1671aSToby Isaac PetscReal inorm; 4842bb1671aSToby Isaac 4852bb1671aSToby Isaac J[2] = J[3] * J[7] - J[6] * J[4]; 4862bb1671aSToby Isaac J[5] = J[6] * J[1] - J[0] * J[7]; 4872bb1671aSToby Isaac J[8] = J[0] * J[4] - J[3] * J[1]; 4882bb1671aSToby Isaac 4892bb1671aSToby Isaac inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]); 4902bb1671aSToby Isaac 4912bb1671aSToby Isaac J[2] *= inorm; 4922bb1671aSToby Isaac J[5] *= inorm; 4932bb1671aSToby Isaac J[8] *= inorm; 4942bb1671aSToby Isaac } 4952bb1671aSToby Isaac break; 4962bb1671aSToby Isaac } 4972bb1671aSToby Isaac for (j = 0; j < dE*dE; j++) { 4982bb1671aSToby Isaac g->J[i*dE*dE + j] = J[j]; 4992bb1671aSToby Isaac } 5002bb1671aSToby Isaac } 5012bb1671aSToby Isaac } 5029566063dSJacob Faibussowitsch PetscCall(PetscFEGeomComplete(g)); 5039566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(field,pointIS,NULL,&maxDegree)); 504b7260050SToby Isaac g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE; 5050145028aSToby Isaac if (faceData) { 5069566063dSJacob Faibussowitsch PetscCall((*field->ops->computeFaceData) (field, pointIS, quad, g)); 5070145028aSToby Isaac } 5082bb1671aSToby Isaac *geom = g; 5092bb1671aSToby Isaac PetscFunctionReturn(0); 5102bb1671aSToby Isaac } 511