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); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldInitializePackage()); 213da551e6SToby Isaac 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 413da551e6SToby Isaac .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*5f80ce2aSJacob Faibussowitsch if ((*field)->ops->destroy) CHKERRQ((*(*field)->ops->destroy)(*field)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&((*field)->dm))); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 663da551e6SToby Isaac .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); 74*5f80ce2aSJacob Faibussowitsch if (!viewer) CHKERRQ(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer)); 753da551e6SToby Isaac PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 763da551e6SToby Isaac PetscCheckSameComm(field,1,viewer,2); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 783da551e6SToby Isaac if (iascii) { 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D components\n",field->numComponents)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity])); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(field->dm,viewer)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 863da551e6SToby Isaac } 87*5f80ce2aSJacob Faibussowitsch if (field->ops->view) CHKERRQ((*field->ops->view)(field,viewer)); 883da551e6SToby Isaac if (iascii) { 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 903da551e6SToby Isaac } 913da551e6SToby Isaac PetscFunctionReturn(0); 923da551e6SToby Isaac } 933da551e6SToby Isaac 943da551e6SToby Isaac /*@C 953da551e6SToby Isaac DMFieldSetType - set the DMField implementation 963da551e6SToby Isaac 97d083f849SBarry Smith Collective on field 983da551e6SToby Isaac 993da551e6SToby Isaac Input Parameters: 1003da551e6SToby Isaac + field - the DMField context 1013da551e6SToby Isaac - type - a known method 1023da551e6SToby Isaac 1033da551e6SToby Isaac Notes: 1043da551e6SToby Isaac See "include/petscvec.h" for available methods (for instance) 1053da551e6SToby Isaac + DMFIELDDA - a field defined only by its values at the corners of a DMDA 1063da551e6SToby Isaac . DMFIELDDS - a field defined by a discretization over a mesh set with DMSetField() 1073da551e6SToby Isaac - DMFIELDSHELL - a field defined by arbitrary callbacks 1083da551e6SToby Isaac 1093da551e6SToby Isaac Level: advanced 1103da551e6SToby Isaac 1113da551e6SToby Isaac .seealso: DMFieldType, 1123da551e6SToby Isaac @*/ 1133da551e6SToby Isaac PetscErrorCode DMFieldSetType(DMField field,DMFieldType type) 1143da551e6SToby Isaac { 1153da551e6SToby Isaac PetscBool match; 116*5f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(DMField); 1173da551e6SToby Isaac 1183da551e6SToby Isaac PetscFunctionBegin; 1193da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1203da551e6SToby Isaac PetscValidCharPointer(type,2); 1213da551e6SToby Isaac 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)field,type,&match)); 1233da551e6SToby Isaac if (match) PetscFunctionReturn(0); 1243da551e6SToby Isaac 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListFind(DMFieldList,type,&r)); 126*5f80ce2aSJacob Faibussowitsch PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type); 1273da551e6SToby Isaac /* Destroy the previous private DMField context */ 128*5f80ce2aSJacob Faibussowitsch if (field->ops->destroy) CHKERRQ((*(field)->ops->destroy)(field)); 129*5f80ce2aSJacob Faibussowitsch 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(field->ops,sizeof(*field->ops))); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)field,type)); 1323da551e6SToby Isaac field->ops->create = r; 133*5f80ce2aSJacob Faibussowitsch CHKERRQ((*r)(field)); 1343da551e6SToby Isaac PetscFunctionReturn(0); 1353da551e6SToby Isaac } 1363da551e6SToby Isaac 1373da551e6SToby Isaac /*@C 1383da551e6SToby Isaac DMFieldGetType - Gets the DMField type name (as a string) from the DMField. 1393da551e6SToby Isaac 1403da551e6SToby Isaac Not Collective 1413da551e6SToby Isaac 1423da551e6SToby Isaac Input Parameter: 1433da551e6SToby Isaac . field - The DMField context 1443da551e6SToby Isaac 1453da551e6SToby Isaac Output Parameter: 1463da551e6SToby Isaac . type - The DMField type name 1473da551e6SToby Isaac 1483da551e6SToby Isaac Level: advanced 1493da551e6SToby Isaac 1503da551e6SToby Isaac .seealso: DMFieldSetType() 1513da551e6SToby Isaac @*/ 1523da551e6SToby Isaac PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 1533da551e6SToby Isaac { 1543da551e6SToby Isaac PetscFunctionBegin; 155064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(field, DMFIELD_CLASSID,1); 1563da551e6SToby Isaac PetscValidPointer(type,2); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldRegisterAll()); 1583da551e6SToby Isaac *type = ((PetscObject)field)->type_name; 1593da551e6SToby Isaac PetscFunctionReturn(0); 1603da551e6SToby Isaac } 1613da551e6SToby Isaac 162da311338SToby Isaac /*@ 163da311338SToby Isaac DMFieldGetNumComponents - Returns the number of components in the field 164da311338SToby Isaac 165da311338SToby Isaac Not collective 166da311338SToby Isaac 167da311338SToby Isaac Input Parameter: 168da311338SToby Isaac . field - The DMField object 169da311338SToby Isaac 170da311338SToby Isaac Output Parameter: 171da311338SToby Isaac . nc - The number of field components 172da311338SToby Isaac 173da311338SToby Isaac Level: intermediate 174da311338SToby Isaac 175da311338SToby Isaac .seealso: DMFieldEvaluate() 176da311338SToby Isaac @*/ 1773da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 1783da551e6SToby Isaac { 1793da551e6SToby Isaac PetscFunctionBegin; 1803da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1813da551e6SToby Isaac PetscValidIntPointer(nc,2); 1823da551e6SToby Isaac *nc = field->numComponents; 1833da551e6SToby Isaac PetscFunctionReturn(0); 1843da551e6SToby Isaac } 1853da551e6SToby Isaac 186da311338SToby Isaac /*@ 187da311338SToby Isaac DMFieldGetDM - Returns the DM for the manifold over which the field is defined. 188da311338SToby Isaac 189da311338SToby Isaac Not collective 190da311338SToby Isaac 191da311338SToby Isaac Input Parameter: 192da311338SToby Isaac . field - The DMField object 193da311338SToby Isaac 194da311338SToby Isaac Output Parameter: 195da311338SToby Isaac . dm - The DM object 196da311338SToby Isaac 197da311338SToby Isaac Level: intermediate 198da311338SToby Isaac 199da311338SToby Isaac .seealso: DMFieldEvaluate() 200da311338SToby Isaac @*/ 2013da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 2023da551e6SToby Isaac { 2033da551e6SToby Isaac PetscFunctionBegin; 2043da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2053da551e6SToby Isaac PetscValidPointer(dm,2); 2063da551e6SToby Isaac *dm = field->dm; 2073da551e6SToby Isaac PetscFunctionReturn(0); 2083da551e6SToby Isaac } 2093da551e6SToby Isaac 210da311338SToby Isaac /*@ 211da311338SToby Isaac DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 212da311338SToby Isaac 213d083f849SBarry Smith Collective on points 214da311338SToby Isaac 215d8d19677SJose E. Roman Input Parameters: 216da311338SToby Isaac + field - The DMField object 217da311338SToby Isaac . points - The points at which to evaluate the field. Should have size d x n, 218da311338SToby Isaac where d is the coordinate dimension of the manifold and n is the number 219da311338SToby Isaac of points 220da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 221da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 222da311338SToby Isaac field is returned. 223da311338SToby Isaac 224d8d19677SJose E. Roman Output Parameters: 225da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 226da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 227da311338SToby Isaac then by point. 228da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 229da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 230da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 231da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 232da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 233da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 234da311338SToby Isaac 235da311338SToby Isaac Level: intermediate 236da311338SToby Isaac 237da311338SToby Isaac .seealso: DMFieldGetDM(), DMFieldGetNumComponents(), DMFieldEvaluateFE(), DMFieldEvaluateFV() 238da311338SToby Isaac @*/ 2393da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 2403da551e6SToby Isaac { 2413da551e6SToby Isaac PetscFunctionBegin; 2423da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2433da551e6SToby Isaac PetscValidHeaderSpecific(points,VEC_CLASSID,2); 244064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,4); 245064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,5); 246064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,6); 2473da551e6SToby Isaac if (field->ops->evaluate) { 248*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->evaluate) (field, points, datatype, B, D, H)); 2493da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 2503da551e6SToby Isaac PetscFunctionReturn(0); 2513da551e6SToby Isaac } 2523da551e6SToby Isaac 253da311338SToby Isaac /*@ 254da311338SToby Isaac DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 255da311338SToby Isaac quadrature points on a reference point. The derivatives are taken with respect to the 256da311338SToby Isaac reference coordinates. 257da311338SToby Isaac 258da311338SToby Isaac Not collective 259da311338SToby Isaac 260d8d19677SJose E. Roman Input Parameters: 261da311338SToby Isaac + field - The DMField object 262da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 263da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field. 264da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 265da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 266da311338SToby Isaac field is returned. 267da311338SToby Isaac 268d8d19677SJose E. Roman Output Parameters: 269da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 270da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 271da311338SToby Isaac then by point. 272da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 273da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 274da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 275da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 276da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 277da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 278da311338SToby Isaac 279da311338SToby Isaac Level: intermediate 280da311338SToby Isaac 281da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFV() 282da311338SToby Isaac @*/ 2833da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 2843da551e6SToby Isaac { 2853da551e6SToby Isaac PetscFunctionBegin; 2863da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2873da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 2883da551e6SToby Isaac PetscValidHeader(points,3); 289064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,5); 290064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,6); 291064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,7); 2923da551e6SToby Isaac if (field->ops->evaluateFE) { 293*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H)); 2943da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 2953da551e6SToby Isaac PetscFunctionReturn(0); 2963da551e6SToby Isaac } 2973da551e6SToby Isaac 298da311338SToby Isaac /*@ 299da311338SToby Isaac DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 300da311338SToby Isaac 301da311338SToby Isaac Not collective 302da311338SToby Isaac 303d8d19677SJose E. Roman Input Parameters: 304da311338SToby Isaac + field - The DMField object 305da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 306da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 307da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 308da311338SToby Isaac field is returned. 309da311338SToby Isaac 310d8d19677SJose E. Roman Output Parameters: 311da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 312da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 313da311338SToby Isaac then by point. 314da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 315da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 316da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 317da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 318da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 319da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 320da311338SToby Isaac 321da311338SToby Isaac Level: intermediate 322da311338SToby Isaac 323da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFE() 324da311338SToby Isaac @*/ 3253da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 3263da551e6SToby Isaac { 3273da551e6SToby Isaac PetscFunctionBegin; 3283da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3293da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 330064a246eSJacob Faibussowitsch if (B) PetscValidPointer(B,4); 331064a246eSJacob Faibussowitsch if (D) PetscValidPointer(D,5); 332064a246eSJacob Faibussowitsch if (H) PetscValidPointer(H,6); 3333da551e6SToby Isaac if (field->ops->evaluateFV) { 334*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H)); 3353da551e6SToby Isaac } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 3363da551e6SToby Isaac PetscFunctionReturn(0); 3373da551e6SToby Isaac } 3383da551e6SToby Isaac 339da311338SToby Isaac /*@ 340b7260050SToby Isaac DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the 341da311338SToby Isaac reference element 342da311338SToby Isaac 343da311338SToby Isaac Not collective 344da311338SToby Isaac 3454165533cSJose E. Roman Input Parameters: 346da311338SToby Isaac + field - the DMField object 347da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance 348da311338SToby Isaac 3494165533cSJose E. Roman Output Parameters: 350b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element 351b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element 352da311338SToby Isaac 353da311338SToby Isaac Level: intermediate 354da311338SToby Isaac 355da311338SToby Isaac .seealso: DMFieldEvaluateFE() 356da311338SToby Isaac @*/ 357b7260050SToby Isaac PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) 3583da551e6SToby Isaac { 3593da551e6SToby Isaac PetscFunctionBegin; 3603da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3613da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 362b7260050SToby Isaac if (minDegree) PetscValidPointer(minDegree,3); 363b7260050SToby Isaac if (maxDegree) PetscValidPointer(maxDegree,4); 3643da551e6SToby Isaac 365b7260050SToby Isaac if (minDegree) *minDegree = -1; 366b7260050SToby Isaac if (maxDegree) *maxDegree = PETSC_MAX_INT; 3673da551e6SToby Isaac 368b7260050SToby Isaac if (field->ops->getDegree) { 369*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->getDegree) (field,cellIS,minDegree,maxDegree)); 3703da551e6SToby Isaac } 3713da551e6SToby Isaac PetscFunctionReturn(0); 3723da551e6SToby Isaac } 3733da551e6SToby Isaac 374da311338SToby Isaac /*@ 375da311338SToby Isaac DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 376da311338SToby Isaac points via pullback onto the reference element 377da311338SToby Isaac 378da311338SToby Isaac Not collective 379da311338SToby Isaac 3804165533cSJose E. Roman Input Parameters: 381da311338SToby Isaac + field - the DMField object 382da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field 383da311338SToby Isaac 3844165533cSJose E. Roman Output Parameter: 385da311338SToby Isaac . quad - a PetscQuadrature object 386da311338SToby Isaac 387da311338SToby Isaac Level: developer 388da311338SToby Isaac 389b7260050SToby Isaac .seealso: DMFieldEvaluteFE(), DMFieldGetDegree() 390da311338SToby Isaac @*/ 3913da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 3923da551e6SToby Isaac { 3933da551e6SToby Isaac PetscFunctionBegin; 3943da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3953da551e6SToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 3963da551e6SToby Isaac PetscValidPointer(quad,3); 3973da551e6SToby Isaac 3983da551e6SToby Isaac *quad = NULL; 3993da551e6SToby Isaac if (field->ops->createDefaultQuadrature) { 400*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->createDefaultQuadrature)(field, pointIS, quad)); 4013da551e6SToby Isaac } 4023da551e6SToby Isaac PetscFunctionReturn(0); 4033da551e6SToby Isaac } 4043da551e6SToby Isaac 405da311338SToby Isaac /*@C 406da311338SToby Isaac DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 407da311338SToby Isaac 408da311338SToby Isaac Not collective 409da311338SToby Isaac 4104165533cSJose E. Roman Input Parameters: 411da311338SToby Isaac + field - the DMField object 412da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field 413da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors 414da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 415da311338SToby Isaac be calculated 416da311338SToby Isaac 4174165533cSJose E. Roman Output Parameter: 418da311338SToby Isaac . geom - the geometric factors 419da311338SToby Isaac 420da311338SToby Isaac Level: developer 421da311338SToby Isaac 422b7260050SToby Isaac .seealso: DMFieldEvaluateFE(), DMFieldCreateDefaulteQuadrature(), DMFieldGetDegree() 42306dd6b0eSSatish Balay @*/ 4242bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 4252bb1671aSToby Isaac { 4262bb1671aSToby Isaac PetscInt dim, dE; 4272bb1671aSToby Isaac PetscInt nPoints; 428b7260050SToby Isaac PetscInt maxDegree; 4292bb1671aSToby Isaac PetscFEGeom *g; 4302bb1671aSToby Isaac 4312bb1671aSToby Isaac PetscFunctionBegin; 4322bb1671aSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 4332bb1671aSToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 4342bb1671aSToby Isaac PetscValidHeader(quad,3); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(pointIS,&nPoints)); 4362bb1671aSToby Isaac dE = field->numComponents; 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGeomCreate(quad,nPoints,dE,faceData,&g)); 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL)); 4392bb1671aSToby Isaac dim = g->dim; 4402bb1671aSToby Isaac if (dE > dim) { 4412bb1671aSToby Isaac /* space out J and make square Jacobians */ 4422bb1671aSToby Isaac PetscInt i, j, k, N = g->numPoints * g->numCells; 4432bb1671aSToby Isaac 4442bb1671aSToby Isaac for (i = N-1; i >= 0; i--) { 4452bb1671aSToby Isaac PetscReal J[9]; 4462bb1671aSToby Isaac 4472bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4482bb1671aSToby Isaac for (k = 0; k < dim; k++) { 4492bb1671aSToby Isaac J[j*dE + k] = g->J[i*dE*dim + j*dim + k]; 4502bb1671aSToby Isaac } 4512bb1671aSToby Isaac } 4522bb1671aSToby Isaac switch (dim) { 4532bb1671aSToby Isaac case 0: 4542bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4552bb1671aSToby Isaac for (k = 0; k < dE; k++) { 4562bb1671aSToby Isaac J[j * dE + k] = (j == k) ? 1. : 0.; 4572bb1671aSToby Isaac } 4582bb1671aSToby Isaac } 4592bb1671aSToby Isaac break; 4602bb1671aSToby Isaac case 1: 4612bb1671aSToby Isaac if (dE == 2) { 4622bb1671aSToby Isaac PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 4632bb1671aSToby Isaac 4642bb1671aSToby Isaac J[1] = -J[2] / norm; 4652bb1671aSToby Isaac J[3] = J[0] / norm; 4662bb1671aSToby Isaac } else { 4672bb1671aSToby Isaac PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 4682bb1671aSToby Isaac PetscReal x = J[0] * inorm; 4692bb1671aSToby Isaac PetscReal y = J[3] * inorm; 4702bb1671aSToby Isaac PetscReal z = J[6] * inorm; 4712bb1671aSToby Isaac 4722bb1671aSToby Isaac if (x > 0.) { 4732bb1671aSToby Isaac PetscReal inv1pX = 1./ (1. + x); 4742bb1671aSToby Isaac 4752bb1671aSToby Isaac J[1] = -y; J[2] = -z; 4762bb1671aSToby Isaac J[4] = 1. - y*y*inv1pX; J[5] = -y*z*inv1pX; 4772bb1671aSToby Isaac J[7] = -y*z*inv1pX; J[8] = 1. - z*z*inv1pX; 4782bb1671aSToby Isaac } else { 4792bb1671aSToby Isaac PetscReal inv1mX = 1./ (1. - x); 4802bb1671aSToby Isaac 4812bb1671aSToby Isaac J[1] = z; J[2] = y; 4822bb1671aSToby Isaac J[4] = -y*z*inv1mX; J[5] = 1. - y*y*inv1mX; 4832bb1671aSToby Isaac J[7] = 1. - z*z*inv1mX; J[8] = -y*z*inv1mX; 4842bb1671aSToby Isaac } 4852bb1671aSToby Isaac } 4862bb1671aSToby Isaac break; 4872bb1671aSToby Isaac case 2: 4882bb1671aSToby Isaac { 4892bb1671aSToby Isaac PetscReal inorm; 4902bb1671aSToby Isaac 4912bb1671aSToby Isaac J[2] = J[3] * J[7] - J[6] * J[4]; 4922bb1671aSToby Isaac J[5] = J[6] * J[1] - J[0] * J[7]; 4932bb1671aSToby Isaac J[8] = J[0] * J[4] - J[3] * J[1]; 4942bb1671aSToby Isaac 4952bb1671aSToby Isaac inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]); 4962bb1671aSToby Isaac 4972bb1671aSToby Isaac J[2] *= inorm; 4982bb1671aSToby Isaac J[5] *= inorm; 4992bb1671aSToby Isaac J[8] *= inorm; 5002bb1671aSToby Isaac } 5012bb1671aSToby Isaac break; 5022bb1671aSToby Isaac } 5032bb1671aSToby Isaac for (j = 0; j < dE*dE; j++) { 5042bb1671aSToby Isaac g->J[i*dE*dE + j] = J[j]; 5052bb1671aSToby Isaac } 5062bb1671aSToby Isaac } 5072bb1671aSToby Isaac } 508*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGeomComplete(g)); 509*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetDegree(field,pointIS,NULL,&maxDegree)); 510b7260050SToby Isaac g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE; 5110145028aSToby Isaac if (faceData) { 512*5f80ce2aSJacob Faibussowitsch PetscCheck(field->ops->computeFaceData,PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data"); 513*5f80ce2aSJacob Faibussowitsch CHKERRQ((*field->ops->computeFaceData) (field, pointIS, quad, g)); 5140145028aSToby Isaac } 5152bb1671aSToby Isaac *geom = g; 5162bb1671aSToby Isaac PetscFunctionReturn(0); 5172bb1671aSToby Isaac } 518