1 #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/ 2 #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/ 3 #include <petscdmplex.h> 4 5 const char *const DMFieldContinuities[] = { 6 "VERTEX", 7 "EDGE", 8 "FACET", 9 "CELL", 10 0 11 }; 12 13 PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField *field) 14 { 15 PetscErrorCode ierr; 16 DMField b; 17 18 PetscFunctionBegin; 19 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 20 PetscValidPointer(field,2); 21 ierr = DMFieldInitializePackage();CHKERRQ(ierr); 22 23 ierr = PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView);CHKERRQ(ierr); 24 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 25 b->dm = dm; 26 b->continuity = continuity; 27 b->numComponents = numComponents; 28 *field = b; 29 PetscFunctionReturn(0); 30 } 31 32 /*@ 33 DMFieldDestroy - destroy a DMField 34 35 Collective 36 37 Input Arguments: 38 . field - address of DMField 39 40 Level: advanced 41 42 .seealso: DMFieldCreate() 43 @*/ 44 PetscErrorCode DMFieldDestroy(DMField *field) 45 { 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 if (!*field) PetscFunctionReturn(0); 50 PetscValidHeaderSpecific((*field),DMFIELD_CLASSID,1); 51 if (--((PetscObject)(*field))->refct > 0) {*field = 0; PetscFunctionReturn(0);} 52 if ((*field)->ops->destroy) {ierr = (*(*field)->ops->destroy)(*field);CHKERRQ(ierr);} 53 ierr = DMDestroy(&((*field)->dm));CHKERRQ(ierr); 54 ierr = PetscHeaderDestroy(field);CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 58 /*@C 59 DMFieldView - view a DMField 60 61 Collective 62 63 Input Arguments: 64 + field - DMField 65 - viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD 66 67 Level: advanced 68 69 .seealso: DMFieldCreate() 70 @*/ 71 PetscErrorCode DMFieldView(DMField field,PetscViewer viewer) 72 { 73 PetscErrorCode ierr; 74 PetscBool iascii; 75 76 PetscFunctionBegin; 77 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 78 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer);CHKERRQ(ierr);} 79 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 80 PetscCheckSameComm(field,1,viewer,2); 81 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 82 if (iascii) { 83 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer);CHKERRQ(ierr); 84 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer,"%D components\n",field->numComponents);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity]);CHKERRQ(ierr); 87 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 88 ierr = DMView(field->dm,viewer);CHKERRQ(ierr); 89 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 90 } 91 if (field->ops->view) {ierr = (*field->ops->view)(field,viewer);CHKERRQ(ierr);} 92 if (iascii) { 93 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 94 } 95 PetscFunctionReturn(0); 96 } 97 98 /*@C 99 DMFieldSetType - set the DMField implementation 100 101 Collective on DMField 102 103 Input Parameters: 104 + field - the DMField context 105 - type - a known method 106 107 Notes: 108 See "include/petscvec.h" for available methods (for instance) 109 + DMFIELDDA - a field defined only by its values at the corners of a DMDA 110 . DMFIELDDS - a field defined by a discretization over a mesh set with DMSetField() 111 - DMFIELDSHELL - a field defined by arbitrary callbacks 112 113 Level: advanced 114 115 .keywords: DMField, set, type 116 117 .seealso: DMFieldType, 118 @*/ 119 PetscErrorCode DMFieldSetType(DMField field,DMFieldType type) 120 { 121 PetscErrorCode ierr,(*r)(DMField); 122 PetscBool match; 123 124 PetscFunctionBegin; 125 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 126 PetscValidCharPointer(type,2); 127 128 ierr = PetscObjectTypeCompare((PetscObject)field,type,&match);CHKERRQ(ierr); 129 if (match) PetscFunctionReturn(0); 130 131 ierr = PetscFunctionListFind(DMFieldList,type,&r);CHKERRQ(ierr); 132 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type); 133 /* Destroy the previous private DMField context */ 134 if (field->ops->destroy) { 135 ierr = (*(field)->ops->destroy)(field);CHKERRQ(ierr); 136 } 137 ierr = PetscMemzero(field->ops,sizeof(*field->ops));CHKERRQ(ierr); 138 ierr = PetscObjectChangeTypeName((PetscObject)field,type);CHKERRQ(ierr); 139 field->ops->create = r; 140 ierr = (*r)(field);CHKERRQ(ierr); 141 PetscFunctionReturn(0); 142 } 143 144 /*@C 145 DMFieldGetType - Gets the DMField type name (as a string) from the DMField. 146 147 Not Collective 148 149 Input Parameter: 150 . field - The DMField context 151 152 Output Parameter: 153 . type - The DMField type name 154 155 Level: advanced 156 157 .keywords: DMField, get, type, name 158 .seealso: DMFieldSetType() 159 @*/ 160 PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 161 { 162 PetscErrorCode ierr; 163 164 PetscFunctionBegin; 165 PetscValidHeaderSpecific(field, VEC_TAGGER_CLASSID,1); 166 PetscValidPointer(type,2); 167 ierr = DMFieldRegisterAll();CHKERRQ(ierr); 168 *type = ((PetscObject)field)->type_name; 169 PetscFunctionReturn(0); 170 } 171 172 /*@ 173 DMFieldGetNumComponents - Returns the number of components in the field 174 175 Not collective 176 177 Input Parameter: 178 . field - The DMField object 179 180 Output Parameter: 181 . nc - The number of field components 182 183 Level: intermediate 184 185 .seealso: DMFieldEvaluate() 186 @*/ 187 PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 188 { 189 PetscFunctionBegin; 190 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 191 PetscValidIntPointer(nc,2); 192 *nc = field->numComponents; 193 PetscFunctionReturn(0); 194 } 195 196 /*@ 197 DMFieldGetDM - Returns the DM for the manifold over which the field is defined. 198 199 Not collective 200 201 Input Parameter: 202 . field - The DMField object 203 204 Output Parameter: 205 . dm - The DM object 206 207 Level: intermediate 208 209 .seealso: DMFieldEvaluate() 210 @*/ 211 PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 212 { 213 PetscFunctionBegin; 214 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 215 PetscValidPointer(dm,2); 216 *dm = field->dm; 217 PetscFunctionReturn(0); 218 } 219 220 /*@ 221 DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 222 223 Collective on Vec 224 225 Input Parameter: 226 + field - The DMField object 227 . points - The points at which to evaluate the field. Should have size d x n, 228 where d is the coordinate dimension of the manifold and n is the number 229 of points 230 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 231 If the field is complex and datatype is PETSC_REAL, the real part of the 232 field is returned. 233 234 235 Output Parameter: 236 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 237 If B is not NULL, the values of the field are written in this array, varying first by component, 238 then by point. 239 . D - pointer to data of size d * c * n * sizeof(datatype). 240 If D is not NULL, the values of the field's spatial derivatives are written in this array, 241 varying first by the partial derivative component, then by field component, then by point. 242 - H - pointer to data of size d * d * c * n * sizeof(datatype). 243 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 244 varying first by the second partial derivative component, then by field component, then by point. 245 246 Level: intermediate 247 248 .seealso: DMFieldGetDM(), DMFieldGetNumComponents(), DMFieldEvaluateFE(), DMFieldEvaluateFV() 249 @*/ 250 PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 251 { 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 256 PetscValidHeaderSpecific(points,VEC_CLASSID,2); 257 if (B) PetscValidPointer(B,3); 258 if (D) PetscValidPointer(D,4); 259 if (H) PetscValidPointer(H,5); 260 if (field->ops->evaluate) { 261 ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr); 262 } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 263 PetscFunctionReturn(0); 264 } 265 266 /*@ 267 DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 268 quadrature points on a reference point. The derivatives are taken with respect to the 269 reference coordinates. 270 271 Not collective 272 273 Input Parameter: 274 + field - The DMField object 275 . cellIS - Index set for cells on which to evaluate the field 276 . points - The quadature containing the points in the reference cell at which to evaluate the field. 277 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 278 If the field is complex and datatype is PETSC_REAL, the real part of the 279 field is returned. 280 281 282 Output Parameter: 283 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 284 If B is not NULL, the values of the field are written in this array, varying first by component, 285 then by point. 286 . D - pointer to data of size d * c * n * sizeof(datatype). 287 If D is not NULL, the values of the field's spatial derivatives are written in this array, 288 varying first by the partial derivative component, then by field component, then by point. 289 - H - pointer to data of size d * d * c * n * sizeof(datatype). 290 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 291 varying first by the second partial derivative component, then by field component, then by point. 292 293 Level: intermediate 294 295 .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFV() 296 @*/ 297 PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 298 { 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 303 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 304 PetscValidHeader(points,3); 305 if (B) PetscValidPointer(B,4); 306 if (D) PetscValidPointer(D,5); 307 if (H) PetscValidPointer(H,6); 308 if (field->ops->evaluateFE) { 309 ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr); 310 } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 311 PetscFunctionReturn(0); 312 } 313 314 /*@ 315 DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 316 317 Not collective 318 319 Input Parameter: 320 + field - The DMField object 321 . cellIS - Index set for cells on which to evaluate the field 322 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 323 If the field is complex and datatype is PETSC_REAL, the real part of the 324 field is returned. 325 326 327 Output Parameter: 328 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 329 If B is not NULL, the values of the field are written in this array, varying first by component, 330 then by point. 331 . D - pointer to data of size d * c * n * sizeof(datatype). 332 If D is not NULL, the values of the field's spatial derivatives are written in this array, 333 varying first by the partial derivative component, then by field component, then by point. 334 - H - pointer to data of size d * d * c * n * sizeof(datatype). 335 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 336 varying first by the second partial derivative component, then by field component, then by point. 337 338 Level: intermediate 339 340 .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFE() 341 @*/ 342 PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 343 { 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 348 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 349 if (B) PetscValidPointer(B,3); 350 if (D) PetscValidPointer(D,4); 351 if (H) PetscValidPointer(H,5); 352 if (field->ops->evaluateFV) { 353 ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr); 354 } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 355 PetscFunctionReturn(0); 356 } 357 358 /*@ 359 DMFieldGetFEInvariance - Get information about the complexity of a field when pulled back onto the 360 reference element 361 362 Not collective 363 364 Input Arguments: 365 + field - the DMField object 366 - cellIS - the index set of points over which we want know the invariance 367 368 Output Arguments: 369 + isConstant - Returns PETSC_TRUE if the pullback of the field is constant for each point in the index set 370 . isAffine - Returns PETSC_TRUE if the pullback of the field is affine for each point in the index set 371 - isQuadratic - Returns PETSC_TRUE if the pullback of the field is quadratic for each point in the index set 372 373 Level: intermediate 374 375 .seealso: DMFieldEvaluateFE() 376 @*/ 377 PetscErrorCode DMFieldGetFEInvariance(DMField field, IS cellIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic) 378 { 379 PetscErrorCode ierr; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 383 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 384 if (isConstant) PetscValidPointer(isConstant,3); 385 if (isAffine) PetscValidPointer(isAffine,4); 386 if (isQuadratic) PetscValidPointer(isQuadratic,5); 387 388 if (isConstant) *isConstant = PETSC_FALSE; 389 if (isAffine) *isAffine = PETSC_FALSE; 390 if (isQuadratic) *isQuadratic = PETSC_FALSE; 391 392 if (field->ops->getFEInvariance) { 393 ierr = (*field->ops->getFEInvariance) (field,cellIS,isConstant,isAffine,isQuadratic);CHKERRQ(ierr); 394 } 395 PetscFunctionReturn(0); 396 } 397 398 /*@ 399 DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 400 points via pullback onto the reference element 401 402 Not collective 403 404 Input Arguments: 405 + field - the DMField object 406 - pointIS - the index set of points over which we wish to integrate the field 407 408 Output Arguments: 409 . quad - a PetscQuadrature object 410 411 Level: developer 412 413 .seealso: DMFieldEvaluteFE(), DMFIeldGetFEInvariance() 414 @*/ 415 PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 416 { 417 PetscErrorCode ierr; 418 419 PetscFunctionBegin; 420 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 421 PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 422 PetscValidPointer(quad,3); 423 424 *quad = NULL; 425 if (field->ops->createDefaultQuadrature) { 426 ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr); 427 } 428 PetscFunctionReturn(0); 429 } 430 431 /*@C 432 DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 433 434 Not collective 435 436 Input Arguments: 437 + field - the DMField object 438 . pointIS - the index set of points over which we wish to integrate the field 439 . quad - the quadrature points at which to evaluate the geometric factors 440 - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 441 be calculated 442 443 Output Arguments: 444 . geom - the geometric factors 445 446 Level: developer 447 448 .seealso: DMFieldEvaluateFE(), DMFieldCreateDefaulteQuadrature(), DMFIeldGetFEInvariance() 449 C@*/ 450 PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 451 { 452 PetscInt dim, dE; 453 PetscInt nPoints; 454 PetscFEGeom *g; 455 PetscErrorCode ierr; 456 457 PetscFunctionBegin; 458 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 459 PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 460 PetscValidHeader(quad,3); 461 ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr); 462 dE = field->numComponents; 463 ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr); 464 ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr); 465 dim = g->dim; 466 if (dE > dim) { 467 /* space out J and make square Jacobians */ 468 PetscInt i, j, k, N = g->numPoints * g->numCells; 469 470 for (i = N-1; i >= 0; i--) { 471 PetscReal J[9]; 472 473 for (j = 0; j < dE; j++) { 474 for (k = 0; k < dim; k++) { 475 J[j*dE + k] = g->J[i*dE*dim + j*dim + k]; 476 } 477 } 478 switch (dim) { 479 case 0: 480 for (j = 0; j < dE; j++) { 481 for (k = 0; k < dE; k++) { 482 J[j * dE + k] = (j == k) ? 1. : 0.; 483 } 484 } 485 break; 486 case 1: 487 if (dE == 2) { 488 PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 489 490 J[1] = -J[2] / norm; 491 J[3] = J[0] / norm; 492 } else { 493 PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 494 PetscReal x = J[0] * inorm; 495 PetscReal y = J[3] * inorm; 496 PetscReal z = J[6] * inorm; 497 498 if (x > 0.) { 499 PetscReal inv1pX = 1./ (1. + x); 500 501 J[1] = -y; J[2] = -z; 502 J[4] = 1. - y*y*inv1pX; J[5] = -y*z*inv1pX; 503 J[7] = -y*z*inv1pX; J[8] = 1. - z*z*inv1pX; 504 } else { 505 PetscReal inv1mX = 1./ (1. - x); 506 507 J[1] = z; J[2] = y; 508 J[4] = -y*z*inv1mX; J[5] = 1. - y*y*inv1mX; 509 J[7] = 1. - z*z*inv1mX; J[8] = -y*z*inv1mX; 510 } 511 } 512 break; 513 case 2: 514 { 515 PetscReal inorm; 516 517 J[2] = J[3] * J[7] - J[6] * J[4]; 518 J[5] = J[6] * J[1] - J[0] * J[7]; 519 J[8] = J[0] * J[4] - J[3] * J[1]; 520 521 inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]); 522 523 J[2] *= inorm; 524 J[5] *= inorm; 525 J[8] *= inorm; 526 } 527 break; 528 } 529 for (j = 0; j < dE*dE; j++) { 530 g->J[i*dE*dE + j] = J[j]; 531 } 532 } 533 } 534 ierr = PetscFEGeomComplete(g);CHKERRQ(ierr); 535 ierr = DMFieldGetFEInvariance(field,pointIS,NULL,&g->isAffine,NULL);CHKERRQ(ierr); 536 if (faceData) { 537 if (!field->ops->computeFaceData) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data\n"); 538 ierr = (*field->ops->computeFaceData) (field, pointIS, quad, g);CHKERRQ(ierr); 539 } 540 *geom = g; 541 PetscFunctionReturn(0); 542 } 543