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", 103da551e6SToby Isaac 0 113da551e6SToby Isaac }; 123da551e6SToby Isaac 133da551e6SToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField *field) 143da551e6SToby Isaac { 153da551e6SToby Isaac PetscErrorCode ierr; 163da551e6SToby Isaac DMField b; 173da551e6SToby Isaac 183da551e6SToby Isaac PetscFunctionBegin; 193da551e6SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 203da551e6SToby Isaac PetscValidPointer(field,2); 213da551e6SToby Isaac ierr = DMFieldInitializePackage();CHKERRQ(ierr); 223da551e6SToby Isaac 233da551e6SToby Isaac ierr = PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView);CHKERRQ(ierr); 243da551e6SToby Isaac ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 253da551e6SToby Isaac b->dm = dm; 263da551e6SToby Isaac b->continuity = continuity; 273da551e6SToby Isaac b->numComponents = numComponents; 283da551e6SToby Isaac *field = b; 293da551e6SToby Isaac PetscFunctionReturn(0); 303da551e6SToby Isaac } 313da551e6SToby Isaac 323da551e6SToby Isaac /*@ 333da551e6SToby Isaac DMFieldDestroy - destroy a DMField 343da551e6SToby Isaac 353da551e6SToby Isaac Collective 363da551e6SToby Isaac 373da551e6SToby Isaac Input Arguments: 383da551e6SToby Isaac . field - address of DMField 393da551e6SToby Isaac 403da551e6SToby Isaac Level: advanced 413da551e6SToby Isaac 423da551e6SToby Isaac .seealso: DMFieldCreate() 433da551e6SToby Isaac @*/ 443da551e6SToby Isaac PetscErrorCode DMFieldDestroy(DMField *field) 453da551e6SToby Isaac { 463da551e6SToby Isaac PetscErrorCode ierr; 473da551e6SToby Isaac 483da551e6SToby Isaac PetscFunctionBegin; 493da551e6SToby Isaac if (!*field) PetscFunctionReturn(0); 503da551e6SToby Isaac PetscValidHeaderSpecific((*field),DMFIELD_CLASSID,1); 513da551e6SToby Isaac if (--((PetscObject)(*field))->refct > 0) {*field = 0; PetscFunctionReturn(0);} 523da551e6SToby Isaac if ((*field)->ops->destroy) {ierr = (*(*field)->ops->destroy)(*field);CHKERRQ(ierr);} 533da551e6SToby Isaac ierr = DMDestroy(&((*field)->dm));CHKERRQ(ierr); 543da551e6SToby Isaac ierr = PetscHeaderDestroy(field);CHKERRQ(ierr); 553da551e6SToby Isaac PetscFunctionReturn(0); 563da551e6SToby Isaac } 573da551e6SToby Isaac 583da551e6SToby Isaac /*@C 593da551e6SToby Isaac DMFieldView - view a DMField 603da551e6SToby Isaac 613da551e6SToby Isaac Collective 623da551e6SToby Isaac 633da551e6SToby Isaac Input Arguments: 643da551e6SToby Isaac + field - DMField 653da551e6SToby Isaac - viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD 663da551e6SToby Isaac 673da551e6SToby Isaac Level: advanced 683da551e6SToby Isaac 693da551e6SToby Isaac .seealso: DMFieldCreate() 703da551e6SToby Isaac @*/ 713da551e6SToby Isaac PetscErrorCode DMFieldView(DMField field,PetscViewer viewer) 723da551e6SToby Isaac { 733da551e6SToby Isaac PetscErrorCode ierr; 743da551e6SToby Isaac PetscBool iascii; 753da551e6SToby Isaac 763da551e6SToby Isaac PetscFunctionBegin; 773da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 783da551e6SToby Isaac if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer);CHKERRQ(ierr);} 793da551e6SToby Isaac PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 803da551e6SToby Isaac PetscCheckSameComm(field,1,viewer,2); 813da551e6SToby Isaac ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 823da551e6SToby Isaac if (iascii) { 833da551e6SToby Isaac ierr = PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer);CHKERRQ(ierr); 843da551e6SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 853da551e6SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"%D components\n",field->numComponents);CHKERRQ(ierr); 863da551e6SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity]);CHKERRQ(ierr); 873da551e6SToby Isaac ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 883da551e6SToby Isaac ierr = DMView(field->dm,viewer);CHKERRQ(ierr); 893da551e6SToby Isaac ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 903da551e6SToby Isaac } 913da551e6SToby Isaac if (field->ops->view) {ierr = (*field->ops->view)(field,viewer);CHKERRQ(ierr);} 923da551e6SToby Isaac if (iascii) { 933da551e6SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 943da551e6SToby Isaac } 953da551e6SToby Isaac PetscFunctionReturn(0); 963da551e6SToby Isaac } 973da551e6SToby Isaac 983da551e6SToby Isaac /*@C 993da551e6SToby Isaac DMFieldSetType - set the DMField implementation 1003da551e6SToby Isaac 1013da551e6SToby Isaac Collective on DMField 1023da551e6SToby Isaac 1033da551e6SToby Isaac Input Parameters: 1043da551e6SToby Isaac + field - the DMField context 1053da551e6SToby Isaac - type - a known method 1063da551e6SToby Isaac 1073da551e6SToby Isaac Notes: 1083da551e6SToby Isaac See "include/petscvec.h" for available methods (for instance) 1093da551e6SToby Isaac + DMFIELDDA - a field defined only by its values at the corners of a DMDA 1103da551e6SToby Isaac . DMFIELDDS - a field defined by a discretization over a mesh set with DMSetField() 1113da551e6SToby Isaac - DMFIELDSHELL - a field defined by arbitrary callbacks 1123da551e6SToby Isaac 1133da551e6SToby Isaac Level: advanced 1143da551e6SToby Isaac 1153da551e6SToby Isaac .keywords: DMField, set, type 1163da551e6SToby Isaac 1173da551e6SToby Isaac .seealso: DMFieldType, 1183da551e6SToby Isaac @*/ 1193da551e6SToby Isaac PetscErrorCode DMFieldSetType(DMField field,DMFieldType type) 1203da551e6SToby Isaac { 1213da551e6SToby Isaac PetscErrorCode ierr,(*r)(DMField); 1223da551e6SToby Isaac PetscBool match; 1233da551e6SToby Isaac 1243da551e6SToby Isaac PetscFunctionBegin; 1253da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1263da551e6SToby Isaac PetscValidCharPointer(type,2); 1273da551e6SToby Isaac 1283da551e6SToby Isaac ierr = PetscObjectTypeCompare((PetscObject)field,type,&match);CHKERRQ(ierr); 1293da551e6SToby Isaac if (match) PetscFunctionReturn(0); 1303da551e6SToby Isaac 1313da551e6SToby Isaac ierr = PetscFunctionListFind(DMFieldList,type,&r);CHKERRQ(ierr); 1323da551e6SToby Isaac if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type); 1333da551e6SToby Isaac /* Destroy the previous private DMField context */ 1343da551e6SToby Isaac if (field->ops->destroy) { 1353da551e6SToby Isaac ierr = (*(field)->ops->destroy)(field);CHKERRQ(ierr); 1363da551e6SToby Isaac } 1373da551e6SToby Isaac ierr = PetscMemzero(field->ops,sizeof(*field->ops));CHKERRQ(ierr); 1383da551e6SToby Isaac ierr = PetscObjectChangeTypeName((PetscObject)field,type);CHKERRQ(ierr); 1393da551e6SToby Isaac field->ops->create = r; 1403da551e6SToby Isaac ierr = (*r)(field);CHKERRQ(ierr); 1413da551e6SToby Isaac PetscFunctionReturn(0); 1423da551e6SToby Isaac } 1433da551e6SToby Isaac 1443da551e6SToby Isaac /*@C 1453da551e6SToby Isaac DMFieldGetType - Gets the DMField type name (as a string) from the DMField. 1463da551e6SToby Isaac 1473da551e6SToby Isaac Not Collective 1483da551e6SToby Isaac 1493da551e6SToby Isaac Input Parameter: 1503da551e6SToby Isaac . field - The DMField context 1513da551e6SToby Isaac 1523da551e6SToby Isaac Output Parameter: 1533da551e6SToby Isaac . type - The DMField type name 1543da551e6SToby Isaac 1553da551e6SToby Isaac Level: advanced 1563da551e6SToby Isaac 1573da551e6SToby Isaac .keywords: DMField, get, type, name 1583da551e6SToby Isaac .seealso: DMFieldSetType() 1593da551e6SToby Isaac @*/ 1603da551e6SToby Isaac PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 1613da551e6SToby Isaac { 1623da551e6SToby Isaac PetscErrorCode ierr; 1633da551e6SToby Isaac 1643da551e6SToby Isaac PetscFunctionBegin; 1653da551e6SToby Isaac PetscValidHeaderSpecific(field, VEC_TAGGER_CLASSID,1); 1663da551e6SToby Isaac PetscValidPointer(type,2); 1673da551e6SToby Isaac ierr = DMFieldRegisterAll();CHKERRQ(ierr); 1683da551e6SToby Isaac *type = ((PetscObject)field)->type_name; 1693da551e6SToby Isaac PetscFunctionReturn(0); 1703da551e6SToby Isaac } 1713da551e6SToby Isaac 172*da311338SToby Isaac /*@ 173*da311338SToby Isaac DMFieldGetNumComponents - Returns the number of components in the field 174*da311338SToby Isaac 175*da311338SToby Isaac Not collective 176*da311338SToby Isaac 177*da311338SToby Isaac Input Parameter: 178*da311338SToby Isaac . field - The DMField object 179*da311338SToby Isaac 180*da311338SToby Isaac Output Parameter: 181*da311338SToby Isaac . nc - The number of field components 182*da311338SToby Isaac 183*da311338SToby Isaac Level: intermediate 184*da311338SToby Isaac 185*da311338SToby Isaac .seealso: DMFieldEvaluate() 186*da311338SToby Isaac @*/ 1873da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 1883da551e6SToby Isaac { 1893da551e6SToby Isaac PetscFunctionBegin; 1903da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 1913da551e6SToby Isaac PetscValidIntPointer(nc,2); 1923da551e6SToby Isaac *nc = field->numComponents; 1933da551e6SToby Isaac PetscFunctionReturn(0); 1943da551e6SToby Isaac } 1953da551e6SToby Isaac 196*da311338SToby Isaac /*@ 197*da311338SToby Isaac DMFieldGetDM - Returns the DM for the manifold over which the field is defined. 198*da311338SToby Isaac 199*da311338SToby Isaac Not collective 200*da311338SToby Isaac 201*da311338SToby Isaac Input Parameter: 202*da311338SToby Isaac . field - The DMField object 203*da311338SToby Isaac 204*da311338SToby Isaac Output Parameter: 205*da311338SToby Isaac . dm - The DM object 206*da311338SToby Isaac 207*da311338SToby Isaac Level: intermediate 208*da311338SToby Isaac 209*da311338SToby Isaac .seealso: DMFieldEvaluate() 210*da311338SToby Isaac @*/ 2113da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 2123da551e6SToby Isaac { 2133da551e6SToby Isaac PetscFunctionBegin; 2143da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2153da551e6SToby Isaac PetscValidPointer(dm,2); 2163da551e6SToby Isaac *dm = field->dm; 2173da551e6SToby Isaac PetscFunctionReturn(0); 2183da551e6SToby Isaac } 2193da551e6SToby Isaac 220*da311338SToby Isaac /*@ 221*da311338SToby Isaac DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 222*da311338SToby Isaac 223*da311338SToby Isaac Collective on Vec 224*da311338SToby Isaac 225*da311338SToby Isaac Input Parameter: 226*da311338SToby Isaac + field - The DMField object 227*da311338SToby Isaac . points - The points at which to evaluate the field. Should have size d x n, 228*da311338SToby Isaac where d is the coordinate dimension of the manifold and n is the number 229*da311338SToby Isaac of points 230*da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 231*da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 232*da311338SToby Isaac field is returned. 233*da311338SToby Isaac 234*da311338SToby Isaac 235*da311338SToby Isaac Output Parameter: 236*da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 237*da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 238*da311338SToby Isaac then by point. 239*da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 240*da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 241*da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 242*da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 243*da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 244*da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 245*da311338SToby Isaac 246*da311338SToby Isaac Level: intermediate 247*da311338SToby Isaac 248*da311338SToby Isaac .seealso: DMFieldGetDM(), DMFieldGetNumComponents(), DMFieldEvaluateFE(), DMFieldEvaluateFV() 249*da311338SToby Isaac @*/ 2503da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 2513da551e6SToby Isaac { 2523da551e6SToby Isaac PetscErrorCode ierr; 2533da551e6SToby Isaac 2543da551e6SToby Isaac PetscFunctionBegin; 2553da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 2563da551e6SToby Isaac PetscValidHeaderSpecific(points,VEC_CLASSID,2); 2573da551e6SToby Isaac if (B) PetscValidPointer(B,3); 2583da551e6SToby Isaac if (D) PetscValidPointer(D,4); 2593da551e6SToby Isaac if (H) PetscValidPointer(H,5); 2603da551e6SToby Isaac if (field->ops->evaluate) { 2613da551e6SToby Isaac ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr); 2623da551e6SToby Isaac } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 2633da551e6SToby Isaac PetscFunctionReturn(0); 2643da551e6SToby Isaac } 2653da551e6SToby Isaac 266*da311338SToby Isaac /*@ 267*da311338SToby Isaac DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 268*da311338SToby Isaac quadrature points on a reference point. The derivatives are taken with respect to the 269*da311338SToby Isaac reference coordinates. 270*da311338SToby Isaac 271*da311338SToby Isaac Not collective 272*da311338SToby Isaac 273*da311338SToby Isaac Input Parameter: 274*da311338SToby Isaac + field - The DMField object 275*da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 276*da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field. 277*da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 278*da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 279*da311338SToby Isaac field is returned. 280*da311338SToby Isaac 281*da311338SToby Isaac 282*da311338SToby Isaac Output Parameter: 283*da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 284*da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 285*da311338SToby Isaac then by point. 286*da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 287*da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 288*da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 289*da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 290*da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 291*da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 292*da311338SToby Isaac 293*da311338SToby Isaac Level: intermediate 294*da311338SToby Isaac 295*da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFV() 296*da311338SToby Isaac @*/ 2973da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 2983da551e6SToby Isaac { 2993da551e6SToby Isaac PetscErrorCode ierr; 3003da551e6SToby Isaac 3013da551e6SToby Isaac PetscFunctionBegin; 3023da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3033da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 3043da551e6SToby Isaac PetscValidHeader(points,3); 3053da551e6SToby Isaac if (B) PetscValidPointer(B,4); 3063da551e6SToby Isaac if (D) PetscValidPointer(D,5); 3073da551e6SToby Isaac if (H) PetscValidPointer(H,6); 3083da551e6SToby Isaac if (field->ops->evaluateFE) { 3093da551e6SToby Isaac ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr); 3103da551e6SToby Isaac } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 3113da551e6SToby Isaac PetscFunctionReturn(0); 3123da551e6SToby Isaac } 3133da551e6SToby Isaac 314*da311338SToby Isaac /*@ 315*da311338SToby Isaac DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 316*da311338SToby Isaac 317*da311338SToby Isaac Not collective 318*da311338SToby Isaac 319*da311338SToby Isaac Input Parameter: 320*da311338SToby Isaac + field - The DMField object 321*da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field 322*da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 323*da311338SToby Isaac If the field is complex and datatype is PETSC_REAL, the real part of the 324*da311338SToby Isaac field is returned. 325*da311338SToby Isaac 326*da311338SToby Isaac 327*da311338SToby Isaac Output Parameter: 328*da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 329*da311338SToby Isaac If B is not NULL, the values of the field are written in this array, varying first by component, 330*da311338SToby Isaac then by point. 331*da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype). 332*da311338SToby Isaac If D is not NULL, the values of the field's spatial derivatives are written in this array, 333*da311338SToby Isaac varying first by the partial derivative component, then by field component, then by point. 334*da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype). 335*da311338SToby Isaac If H is not NULL, the values of the field's second spatial derivatives are written in this array, 336*da311338SToby Isaac varying first by the second partial derivative component, then by field component, then by point. 337*da311338SToby Isaac 338*da311338SToby Isaac Level: intermediate 339*da311338SToby Isaac 340*da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFE() 341*da311338SToby Isaac @*/ 3423da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 3433da551e6SToby Isaac { 3443da551e6SToby Isaac PetscErrorCode ierr; 3453da551e6SToby Isaac 3463da551e6SToby Isaac PetscFunctionBegin; 3473da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3483da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 3493da551e6SToby Isaac if (B) PetscValidPointer(B,3); 3503da551e6SToby Isaac if (D) PetscValidPointer(D,4); 3513da551e6SToby Isaac if (H) PetscValidPointer(H,5); 3523da551e6SToby Isaac if (field->ops->evaluateFV) { 3533da551e6SToby Isaac ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr); 3543da551e6SToby Isaac } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 3553da551e6SToby Isaac PetscFunctionReturn(0); 3563da551e6SToby Isaac } 3573da551e6SToby Isaac 358*da311338SToby Isaac /*@ 359*da311338SToby Isaac DMFieldGetFEInvariance - Get information about the complexity of a field when pulled back onto the 360*da311338SToby Isaac reference element 361*da311338SToby Isaac 362*da311338SToby Isaac Not collective 363*da311338SToby Isaac 364*da311338SToby Isaac Input Arguments: 365*da311338SToby Isaac + field - the DMField object 366*da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance 367*da311338SToby Isaac 368*da311338SToby Isaac Output Arguments: 369*da311338SToby Isaac + isConstant - Returns PETSC_TRUE if the pullback of the field is constant for each point in the index set 370*da311338SToby Isaac . isAffine - Returns PETSC_TRUE if the pullback of the field is affine for each point in the index set 371*da311338SToby Isaac - isQuadratic - Returns PETSC_TRUE if the pullback of the field is quadratic for each point in the index set 372*da311338SToby Isaac 373*da311338SToby Isaac Level: intermediate 374*da311338SToby Isaac 375*da311338SToby Isaac .seealso: DMFieldEvaluateFE() 376*da311338SToby Isaac @*/ 3773da551e6SToby Isaac PetscErrorCode DMFieldGetFEInvariance(DMField field, IS cellIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic) 3783da551e6SToby Isaac { 3793da551e6SToby Isaac PetscErrorCode ierr; 3803da551e6SToby Isaac 3813da551e6SToby Isaac PetscFunctionBegin; 3823da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 3833da551e6SToby Isaac PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 3843da551e6SToby Isaac if (isConstant) PetscValidPointer(isConstant,3); 3853da551e6SToby Isaac if (isAffine) PetscValidPointer(isAffine,4); 3863da551e6SToby Isaac if (isQuadratic) PetscValidPointer(isQuadratic,5); 3873da551e6SToby Isaac 3883da551e6SToby Isaac if (isConstant) *isConstant = PETSC_FALSE; 3893da551e6SToby Isaac if (isAffine) *isAffine = PETSC_FALSE; 3903da551e6SToby Isaac if (isQuadratic) *isQuadratic = PETSC_FALSE; 3913da551e6SToby Isaac 3923da551e6SToby Isaac if (field->ops->getFEInvariance) { 3933da551e6SToby Isaac ierr = (*field->ops->getFEInvariance) (field,cellIS,isConstant,isAffine,isQuadratic);CHKERRQ(ierr); 3943da551e6SToby Isaac } 3953da551e6SToby Isaac PetscFunctionReturn(0); 3963da551e6SToby Isaac } 3973da551e6SToby Isaac 398*da311338SToby Isaac /*@ 399*da311338SToby Isaac DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 400*da311338SToby Isaac points via pullback onto the reference element 401*da311338SToby Isaac 402*da311338SToby Isaac Not collective 403*da311338SToby Isaac 404*da311338SToby Isaac Input Arguments: 405*da311338SToby Isaac + field - the DMField object 406*da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field 407*da311338SToby Isaac 408*da311338SToby Isaac Output Arguments: 409*da311338SToby Isaac . quad - a PetscQuadrature object 410*da311338SToby Isaac 411*da311338SToby Isaac Level: developer 412*da311338SToby Isaac 413*da311338SToby Isaac .seealso: DMFieldEvaluteFE(), DMFIeldGetFEInvariance() 414*da311338SToby Isaac @*/ 4153da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 4163da551e6SToby Isaac { 4173da551e6SToby Isaac PetscErrorCode ierr; 4183da551e6SToby Isaac 4193da551e6SToby Isaac PetscFunctionBegin; 4203da551e6SToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 4213da551e6SToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 4223da551e6SToby Isaac PetscValidPointer(quad,3); 4233da551e6SToby Isaac 4243da551e6SToby Isaac *quad = NULL; 4253da551e6SToby Isaac if (field->ops->createDefaultQuadrature) { 4263da551e6SToby Isaac ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr); 4273da551e6SToby Isaac } 4283da551e6SToby Isaac PetscFunctionReturn(0); 4293da551e6SToby Isaac } 4303da551e6SToby Isaac 431*da311338SToby Isaac /*@C 432*da311338SToby Isaac DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 433*da311338SToby Isaac 434*da311338SToby Isaac Not collective 435*da311338SToby Isaac 436*da311338SToby Isaac Input Arguments: 437*da311338SToby Isaac + field - the DMField object 438*da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field 439*da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors 440*da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 441*da311338SToby Isaac be calculated 442*da311338SToby Isaac 443*da311338SToby Isaac Output Arguments: 444*da311338SToby Isaac . geom - the geometric factors 445*da311338SToby Isaac 446*da311338SToby Isaac Level: developer 447*da311338SToby Isaac 448*da311338SToby Isaac .seealso: DMFieldEvaluateFE(), DMFieldCreateDefaulteQuadrature(), DMFIeldGetFEInvariance() 449*da311338SToby Isaac C@*/ 4502bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 4512bb1671aSToby Isaac { 4522bb1671aSToby Isaac PetscInt dim, dE; 4532bb1671aSToby Isaac PetscInt nPoints; 4542bb1671aSToby Isaac PetscFEGeom *g; 4552bb1671aSToby Isaac PetscErrorCode ierr; 4562bb1671aSToby Isaac 4572bb1671aSToby Isaac PetscFunctionBegin; 4582bb1671aSToby Isaac PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 4592bb1671aSToby Isaac PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 4602bb1671aSToby Isaac PetscValidHeader(quad,3); 4612bb1671aSToby Isaac ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr); 4622bb1671aSToby Isaac dE = field->numComponents; 4632bb1671aSToby Isaac ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr); 4642bb1671aSToby Isaac ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr); 4652bb1671aSToby Isaac dim = g->dim; 4662bb1671aSToby Isaac if (dE > dim) { 4672bb1671aSToby Isaac /* space out J and make square Jacobians */ 4682bb1671aSToby Isaac PetscInt i, j, k, N = g->numPoints * g->numCells; 4692bb1671aSToby Isaac 4702bb1671aSToby Isaac for (i = N-1; i >= 0; i--) { 4712bb1671aSToby Isaac PetscReal J[9]; 4722bb1671aSToby Isaac 4732bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4742bb1671aSToby Isaac for (k = 0; k < dim; k++) { 4752bb1671aSToby Isaac J[j*dE + k] = g->J[i*dE*dim + j*dim + k]; 4762bb1671aSToby Isaac } 4772bb1671aSToby Isaac } 4782bb1671aSToby Isaac switch (dim) { 4792bb1671aSToby Isaac case 0: 4802bb1671aSToby Isaac for (j = 0; j < dE; j++) { 4812bb1671aSToby Isaac for (k = 0; k < dE; k++) { 4822bb1671aSToby Isaac J[j * dE + k] = (j == k) ? 1. : 0.; 4832bb1671aSToby Isaac } 4842bb1671aSToby Isaac } 4852bb1671aSToby Isaac break; 4862bb1671aSToby Isaac case 1: 4872bb1671aSToby Isaac if (dE == 2) { 4882bb1671aSToby Isaac PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 4892bb1671aSToby Isaac 4902bb1671aSToby Isaac J[1] = -J[2] / norm; 4912bb1671aSToby Isaac J[3] = J[0] / norm; 4922bb1671aSToby Isaac } else { 4932bb1671aSToby Isaac PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 4942bb1671aSToby Isaac PetscReal x = J[0] * inorm; 4952bb1671aSToby Isaac PetscReal y = J[3] * inorm; 4962bb1671aSToby Isaac PetscReal z = J[6] * inorm; 4972bb1671aSToby Isaac 4982bb1671aSToby Isaac if (x > 0.) { 4992bb1671aSToby Isaac PetscReal inv1pX = 1./ (1. + x); 5002bb1671aSToby Isaac 5012bb1671aSToby Isaac J[1] = -y; J[2] = -z; 5022bb1671aSToby Isaac J[4] = 1. - y*y*inv1pX; J[5] = -y*z*inv1pX; 5032bb1671aSToby Isaac J[7] = -y*z*inv1pX; J[8] = 1. - z*z*inv1pX; 5042bb1671aSToby Isaac } else { 5052bb1671aSToby Isaac PetscReal inv1mX = 1./ (1. - x); 5062bb1671aSToby Isaac 5072bb1671aSToby Isaac J[1] = z; J[2] = y; 5082bb1671aSToby Isaac J[4] = -y*z*inv1mX; J[5] = 1. - y*y*inv1mX; 5092bb1671aSToby Isaac J[7] = 1. - z*z*inv1mX; J[8] = -y*z*inv1mX; 5102bb1671aSToby Isaac } 5112bb1671aSToby Isaac } 5122bb1671aSToby Isaac break; 5132bb1671aSToby Isaac case 2: 5142bb1671aSToby Isaac { 5152bb1671aSToby Isaac PetscReal inorm; 5162bb1671aSToby Isaac 5172bb1671aSToby Isaac J[2] = J[3] * J[7] - J[6] * J[4]; 5182bb1671aSToby Isaac J[5] = J[6] * J[1] - J[0] * J[7]; 5192bb1671aSToby Isaac J[8] = J[0] * J[4] - J[3] * J[1]; 5202bb1671aSToby Isaac 5212bb1671aSToby Isaac inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]); 5222bb1671aSToby Isaac 5232bb1671aSToby Isaac J[2] *= inorm; 5242bb1671aSToby Isaac J[5] *= inorm; 5252bb1671aSToby Isaac J[8] *= inorm; 5262bb1671aSToby Isaac } 5272bb1671aSToby Isaac break; 5282bb1671aSToby Isaac } 5292bb1671aSToby Isaac for (j = 0; j < dE*dE; j++) { 5302bb1671aSToby Isaac g->J[i*dE*dE + j] = J[j]; 5312bb1671aSToby Isaac } 5322bb1671aSToby Isaac } 5332bb1671aSToby Isaac } 5342bb1671aSToby Isaac ierr = PetscFEGeomComplete(g);CHKERRQ(ierr); 5352bb1671aSToby Isaac ierr = DMFieldGetFEInvariance(field,pointIS,NULL,&g->isAffine,NULL);CHKERRQ(ierr); 5360145028aSToby Isaac if (faceData) { 5370145028aSToby Isaac if (!field->ops->computeFaceData) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data\n"); 5380145028aSToby Isaac ierr = (*field->ops->computeFaceData) (field, pointIS, quad, g);CHKERRQ(ierr); 5390145028aSToby Isaac } 5402bb1671aSToby Isaac *geom = g; 5412bb1671aSToby Isaac PetscFunctionReturn(0); 5422bb1671aSToby Isaac } 543