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 NULL 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 Parameter: 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 = NULL; 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 Parameters: 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 field 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 .seealso: DMFieldType, 116 @*/ 117 PetscErrorCode DMFieldSetType(DMField field,DMFieldType type) 118 { 119 PetscErrorCode ierr,(*r)(DMField); 120 PetscBool match; 121 122 PetscFunctionBegin; 123 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 124 PetscValidCharPointer(type,2); 125 126 ierr = PetscObjectTypeCompare((PetscObject)field,type,&match);CHKERRQ(ierr); 127 if (match) PetscFunctionReturn(0); 128 129 ierr = PetscFunctionListFind(DMFieldList,type,&r);CHKERRQ(ierr); 130 PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type); 131 /* Destroy the previous private DMField context */ 132 if (field->ops->destroy) { 133 ierr = (*(field)->ops->destroy)(field);CHKERRQ(ierr); 134 } 135 ierr = PetscMemzero(field->ops,sizeof(*field->ops));CHKERRQ(ierr); 136 ierr = PetscObjectChangeTypeName((PetscObject)field,type);CHKERRQ(ierr); 137 field->ops->create = r; 138 ierr = (*r)(field);CHKERRQ(ierr); 139 PetscFunctionReturn(0); 140 } 141 142 /*@C 143 DMFieldGetType - Gets the DMField type name (as a string) from the DMField. 144 145 Not Collective 146 147 Input Parameter: 148 . field - The DMField context 149 150 Output Parameter: 151 . type - The DMField type name 152 153 Level: advanced 154 155 .seealso: DMFieldSetType() 156 @*/ 157 PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) 158 { 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 PetscValidHeaderSpecific(field, DMFIELD_CLASSID,1); 163 PetscValidPointer(type,2); 164 ierr = DMFieldRegisterAll();CHKERRQ(ierr); 165 *type = ((PetscObject)field)->type_name; 166 PetscFunctionReturn(0); 167 } 168 169 /*@ 170 DMFieldGetNumComponents - Returns the number of components in the field 171 172 Not collective 173 174 Input Parameter: 175 . field - The DMField object 176 177 Output Parameter: 178 . nc - The number of field components 179 180 Level: intermediate 181 182 .seealso: DMFieldEvaluate() 183 @*/ 184 PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) 185 { 186 PetscFunctionBegin; 187 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 188 PetscValidIntPointer(nc,2); 189 *nc = field->numComponents; 190 PetscFunctionReturn(0); 191 } 192 193 /*@ 194 DMFieldGetDM - Returns the DM for the manifold over which the field is defined. 195 196 Not collective 197 198 Input Parameter: 199 . field - The DMField object 200 201 Output Parameter: 202 . dm - The DM object 203 204 Level: intermediate 205 206 .seealso: DMFieldEvaluate() 207 @*/ 208 PetscErrorCode DMFieldGetDM(DMField field, DM *dm) 209 { 210 PetscFunctionBegin; 211 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 212 PetscValidPointer(dm,2); 213 *dm = field->dm; 214 PetscFunctionReturn(0); 215 } 216 217 /*@ 218 DMFieldEvaluate - Evaluate the field and its derivatives on a set of points 219 220 Collective on points 221 222 Input Parameters: 223 + field - The DMField object 224 . points - The points at which to evaluate the field. Should have size d x n, 225 where d is the coordinate dimension of the manifold and n is the number 226 of points 227 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 228 If the field is complex and datatype is PETSC_REAL, the real part of the 229 field is returned. 230 231 Output Parameters: 232 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 233 If B is not NULL, the values of the field are written in this array, varying first by component, 234 then by point. 235 . D - pointer to data of size d * c * n * sizeof(datatype). 236 If D is not NULL, the values of the field's spatial derivatives are written in this array, 237 varying first by the partial derivative component, then by field component, then by point. 238 - H - pointer to data of size d * d * c * n * sizeof(datatype). 239 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 240 varying first by the second partial derivative component, then by field component, then by point. 241 242 Level: intermediate 243 244 .seealso: DMFieldGetDM(), DMFieldGetNumComponents(), DMFieldEvaluateFE(), DMFieldEvaluateFV() 245 @*/ 246 PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) 247 { 248 PetscErrorCode ierr; 249 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 252 PetscValidHeaderSpecific(points,VEC_CLASSID,2); 253 if (B) PetscValidPointer(B,4); 254 if (D) PetscValidPointer(D,5); 255 if (H) PetscValidPointer(H,6); 256 if (field->ops->evaluate) { 257 ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr); 258 } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 259 PetscFunctionReturn(0); 260 } 261 262 /*@ 263 DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from 264 quadrature points on a reference point. The derivatives are taken with respect to the 265 reference coordinates. 266 267 Not collective 268 269 Input Parameters: 270 + field - The DMField object 271 . cellIS - Index set for cells on which to evaluate the field 272 . points - The quadature containing the points in the reference cell at which to evaluate the field. 273 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 274 If the field is complex and datatype is PETSC_REAL, the real part of the 275 field is returned. 276 277 Output Parameters: 278 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 279 If B is not NULL, the values of the field are written in this array, varying first by component, 280 then by point. 281 . D - pointer to data of size d * c * n * sizeof(datatype). 282 If D is not NULL, the values of the field's spatial derivatives are written in this array, 283 varying first by the partial derivative component, then by field component, then by point. 284 - H - pointer to data of size d * d * c * n * sizeof(datatype). 285 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 286 varying first by the second partial derivative component, then by field component, then by point. 287 288 Level: intermediate 289 290 .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFV() 291 @*/ 292 PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) 293 { 294 PetscErrorCode ierr; 295 296 PetscFunctionBegin; 297 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 298 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 299 PetscValidHeader(points,3); 300 if (B) PetscValidPointer(B,5); 301 if (D) PetscValidPointer(D,6); 302 if (H) PetscValidPointer(H,7); 303 if (field->ops->evaluateFE) { 304 ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr); 305 } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 306 PetscFunctionReturn(0); 307 } 308 309 /*@ 310 DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points. 311 312 Not collective 313 314 Input Parameters: 315 + field - The DMField object 316 . cellIS - Index set for cells on which to evaluate the field 317 - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR. 318 If the field is complex and datatype is PETSC_REAL, the real part of the 319 field is returned. 320 321 Output Parameters: 322 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. 323 If B is not NULL, the values of the field are written in this array, varying first by component, 324 then by point. 325 . D - pointer to data of size d * c * n * sizeof(datatype). 326 If D is not NULL, the values of the field's spatial derivatives are written in this array, 327 varying first by the partial derivative component, then by field component, then by point. 328 - H - pointer to data of size d * d * c * n * sizeof(datatype). 329 If H is not NULL, the values of the field's second spatial derivatives are written in this array, 330 varying first by the second partial derivative component, then by field component, then by point. 331 332 Level: intermediate 333 334 .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFE() 335 @*/ 336 PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) 337 { 338 PetscErrorCode ierr; 339 340 PetscFunctionBegin; 341 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 342 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 343 if (B) PetscValidPointer(B,4); 344 if (D) PetscValidPointer(D,5); 345 if (H) PetscValidPointer(H,6); 346 if (field->ops->evaluateFV) { 347 ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr); 348 } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type"); 349 PetscFunctionReturn(0); 350 } 351 352 /*@ 353 DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the 354 reference element 355 356 Not collective 357 358 Input Parameters: 359 + field - the DMField object 360 - cellIS - the index set of points over which we want know the invariance 361 362 Output Parameters: 363 + minDegree - the degree of the largest polynomial space contained in the field on each element 364 - maxDegree - the largest degree of the smallest polynomial space containing the field on any element 365 366 Level: intermediate 367 368 .seealso: DMFieldEvaluateFE() 369 @*/ 370 PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) 371 { 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 376 PetscValidHeaderSpecific(cellIS,IS_CLASSID,2); 377 if (minDegree) PetscValidPointer(minDegree,3); 378 if (maxDegree) PetscValidPointer(maxDegree,4); 379 380 if (minDegree) *minDegree = -1; 381 if (maxDegree) *maxDegree = PETSC_MAX_INT; 382 383 if (field->ops->getDegree) { 384 ierr = (*field->ops->getDegree) (field,cellIS,minDegree,maxDegree);CHKERRQ(ierr); 385 } 386 PetscFunctionReturn(0); 387 } 388 389 /*@ 390 DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected 391 points via pullback onto the reference element 392 393 Not collective 394 395 Input Parameters: 396 + field - the DMField object 397 - pointIS - the index set of points over which we wish to integrate the field 398 399 Output Parameter: 400 . quad - a PetscQuadrature object 401 402 Level: developer 403 404 .seealso: DMFieldEvaluteFE(), DMFieldGetDegree() 405 @*/ 406 PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) 407 { 408 PetscErrorCode ierr; 409 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 412 PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 413 PetscValidPointer(quad,3); 414 415 *quad = NULL; 416 if (field->ops->createDefaultQuadrature) { 417 ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr); 418 } 419 PetscFunctionReturn(0); 420 } 421 422 /*@C 423 DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field 424 425 Not collective 426 427 Input Parameters: 428 + field - the DMField object 429 . pointIS - the index set of points over which we wish to integrate the field 430 . quad - the quadrature points at which to evaluate the geometric factors 431 - faceData - whether additional data for facets (the normal vectors and adjacent cells) should 432 be calculated 433 434 Output Parameter: 435 . geom - the geometric factors 436 437 Level: developer 438 439 .seealso: DMFieldEvaluateFE(), DMFieldCreateDefaulteQuadrature(), DMFieldGetDegree() 440 @*/ 441 PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 442 { 443 PetscInt dim, dE; 444 PetscInt nPoints; 445 PetscInt maxDegree; 446 PetscFEGeom *g; 447 PetscErrorCode ierr; 448 449 PetscFunctionBegin; 450 PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1); 451 PetscValidHeaderSpecific(pointIS,IS_CLASSID,2); 452 PetscValidHeader(quad,3); 453 ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr); 454 dE = field->numComponents; 455 ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr); 456 ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr); 457 dim = g->dim; 458 if (dE > dim) { 459 /* space out J and make square Jacobians */ 460 PetscInt i, j, k, N = g->numPoints * g->numCells; 461 462 for (i = N-1; i >= 0; i--) { 463 PetscReal J[9]; 464 465 for (j = 0; j < dE; j++) { 466 for (k = 0; k < dim; k++) { 467 J[j*dE + k] = g->J[i*dE*dim + j*dim + k]; 468 } 469 } 470 switch (dim) { 471 case 0: 472 for (j = 0; j < dE; j++) { 473 for (k = 0; k < dE; k++) { 474 J[j * dE + k] = (j == k) ? 1. : 0.; 475 } 476 } 477 break; 478 case 1: 479 if (dE == 2) { 480 PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]); 481 482 J[1] = -J[2] / norm; 483 J[3] = J[0] / norm; 484 } else { 485 PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]); 486 PetscReal x = J[0] * inorm; 487 PetscReal y = J[3] * inorm; 488 PetscReal z = J[6] * inorm; 489 490 if (x > 0.) { 491 PetscReal inv1pX = 1./ (1. + x); 492 493 J[1] = -y; J[2] = -z; 494 J[4] = 1. - y*y*inv1pX; J[5] = -y*z*inv1pX; 495 J[7] = -y*z*inv1pX; J[8] = 1. - z*z*inv1pX; 496 } else { 497 PetscReal inv1mX = 1./ (1. - x); 498 499 J[1] = z; J[2] = y; 500 J[4] = -y*z*inv1mX; J[5] = 1. - y*y*inv1mX; 501 J[7] = 1. - z*z*inv1mX; J[8] = -y*z*inv1mX; 502 } 503 } 504 break; 505 case 2: 506 { 507 PetscReal inorm; 508 509 J[2] = J[3] * J[7] - J[6] * J[4]; 510 J[5] = J[6] * J[1] - J[0] * J[7]; 511 J[8] = J[0] * J[4] - J[3] * J[1]; 512 513 inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]); 514 515 J[2] *= inorm; 516 J[5] *= inorm; 517 J[8] *= inorm; 518 } 519 break; 520 } 521 for (j = 0; j < dE*dE; j++) { 522 g->J[i*dE*dE + j] = J[j]; 523 } 524 } 525 } 526 ierr = PetscFEGeomComplete(g);CHKERRQ(ierr); 527 ierr = DMFieldGetDegree(field,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 528 g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE; 529 if (faceData) { 530 PetscCheckFalse(!field->ops->computeFaceData,PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data"); 531 ierr = (*field->ops->computeFaceData) (field, pointIS, quad, g);CHKERRQ(ierr); 532 } 533 *geom = g; 534 PetscFunctionReturn(0); 535 } 536