1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petscfe.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DMPlexGetScale" 7 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 8 { 9 DM_Plex *mesh = (DM_Plex*) dm->data; 10 11 PetscFunctionBegin; 12 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13 PetscValidPointer(scale, 3); 14 *scale = mesh->scale[unit]; 15 PetscFunctionReturn(0); 16 } 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "DMPlexSetScale" 20 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 21 { 22 DM_Plex *mesh = (DM_Plex*) dm->data; 23 24 PetscFunctionBegin; 25 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 26 mesh->scale[unit] = scale; 27 PetscFunctionReturn(0); 28 } 29 30 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 31 { 32 switch (i) { 33 case 0: 34 switch (j) { 35 case 0: return 0; 36 case 1: 37 switch (k) { 38 case 0: return 0; 39 case 1: return 0; 40 case 2: return 1; 41 } 42 case 2: 43 switch (k) { 44 case 0: return 0; 45 case 1: return -1; 46 case 2: return 0; 47 } 48 } 49 case 1: 50 switch (j) { 51 case 0: 52 switch (k) { 53 case 0: return 0; 54 case 1: return 0; 55 case 2: return -1; 56 } 57 case 1: return 0; 58 case 2: 59 switch (k) { 60 case 0: return 1; 61 case 1: return 0; 62 case 2: return 0; 63 } 64 } 65 case 2: 66 switch (j) { 67 case 0: 68 switch (k) { 69 case 0: return 0; 70 case 1: return 1; 71 case 2: return 0; 72 } 73 case 1: 74 switch (k) { 75 case 0: return -1; 76 case 1: return 0; 77 case 2: return 0; 78 } 79 case 2: return 0; 80 } 81 } 82 return 0; 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "DMPlexCreateRigidBody" 87 /*@C 88 DMPlexCreateRigidBody - create rigid body modes from coordinates 89 90 Collective on DM 91 92 Input Arguments: 93 + dm - the DM 94 . section - the local section associated with the rigid field, or NULL for the default section 95 - globalSection - the global section associated with the rigid field, or NULL for the default section 96 97 Output Argument: 98 . sp - the null space 99 100 Note: This is necessary to take account of Dirichlet conditions on the displacements 101 102 Level: advanced 103 104 .seealso: MatNullSpaceCreate() 105 @*/ 106 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 107 { 108 MPI_Comm comm; 109 Vec coordinates, localMode, mode[6]; 110 PetscSection coordSection; 111 PetscScalar *coords; 112 PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 113 PetscErrorCode ierr; 114 115 PetscFunctionBegin; 116 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 117 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 118 if (dim == 1) { 119 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 120 PetscFunctionReturn(0); 121 } 122 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 123 if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 124 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 125 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 126 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 127 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 128 m = (dim*(dim+1))/2; 129 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 130 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 131 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 132 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 133 /* Assume P1 */ 134 ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 135 for (d = 0; d < dim; ++d) { 136 PetscScalar values[3] = {0.0, 0.0, 0.0}; 137 138 values[d] = 1.0; 139 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 140 for (v = vStart; v < vEnd; ++v) { 141 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 142 } 143 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 144 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145 } 146 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 147 for (d = dim; d < dim*(dim+1)/2; ++d) { 148 PetscInt i, j, k = dim > 2 ? d - dim : d; 149 150 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 151 for (v = vStart; v < vEnd; ++v) { 152 PetscScalar values[3] = {0.0, 0.0, 0.0}; 153 PetscInt off; 154 155 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 156 for (i = 0; i < dim; ++i) { 157 for (j = 0; j < dim; ++j) { 158 values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 159 } 160 } 161 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 162 } 163 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 164 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165 } 166 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 167 ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 168 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 169 /* Orthonormalize system */ 170 for (i = dim; i < m; ++i) { 171 PetscScalar dots[6]; 172 173 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 174 for (j = 0; j < i; ++j) dots[j] *= -1.0; 175 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 176 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 177 } 178 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 179 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 180 PetscFunctionReturn(0); 181 } 182 183 #undef __FUNCT__ 184 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal" 185 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 186 { 187 PetscDualSpace *sp; 188 PetscSection section; 189 PetscScalar *values; 190 PetscReal *v0, *J, detJ; 191 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 196 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 197 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 198 ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr); 199 for (f = 0; f < numFields; ++f) { 200 ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 201 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 202 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 203 totDim += spDim*numComp; 204 } 205 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 206 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 207 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 208 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 209 for (i = 0; i < numIds; ++i) { 210 IS pointIS; 211 const PetscInt *points; 212 PetscInt n, p; 213 214 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 215 ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 216 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 217 for (p = 0; p < n; ++p) { 218 const PetscInt point = points[p]; 219 PetscCellGeometry geom; 220 221 if ((point < cStart) || (point >= cEnd)) continue; 222 ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr); 223 geom.v0 = v0; 224 geom.J = J; 225 geom.detJ = &detJ; 226 for (f = 0, v = 0; f < numFields; ++f) { 227 void * const ctx = ctxs ? ctxs[f] : NULL; 228 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 229 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 230 for (d = 0; d < spDim; ++d) { 231 if (funcs[f]) { 232 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 233 } else { 234 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 235 } 236 v += numComp; 237 } 238 } 239 ierr = DMPlexVecSetClosure(dm, section, localX, point, values, mode);CHKERRQ(ierr); 240 } 241 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 242 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 243 } 244 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 245 ierr = PetscFree3(sp,v0,J);CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "DMPlexProjectFunctionLocal" 251 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 252 { 253 PetscDualSpace *sp; 254 PetscSection section; 255 PetscScalar *values; 256 PetscReal *v0, *J, detJ; 257 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 262 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 263 ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr); 264 for (f = 0; f < numFields; ++f) { 265 ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 266 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 267 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 268 totDim += spDim*numComp; 269 } 270 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 271 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 272 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 273 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 274 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 275 ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr); 276 for (c = cStart; c < cEnd; ++c) { 277 PetscCellGeometry geom; 278 279 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 280 geom.v0 = v0; 281 geom.J = J; 282 geom.detJ = &detJ; 283 for (f = 0, v = 0; f < numFields; ++f) { 284 void * const ctx = ctxs ? ctxs[f] : NULL; 285 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 286 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 287 for (d = 0; d < spDim; ++d) { 288 if (funcs[f]) { 289 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 290 } else { 291 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 292 } 293 v += numComp; 294 } 295 } 296 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 297 } 298 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 299 ierr = PetscFree2(v0,J);CHKERRQ(ierr); 300 ierr = PetscFree(sp);CHKERRQ(ierr); 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "DMPlexProjectFunction" 306 /*@C 307 DMPlexProjectFunction - This projects the given function into the function space provided. 308 309 Input Parameters: 310 + dm - The DM 311 . fe - The PetscFE associated with the field 312 . funcs - The coordinate functions to evaluate, one per field 313 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 314 - mode - The insertion mode for values 315 316 Output Parameter: 317 . X - vector 318 319 Level: developer 320 321 .seealso: DMPlexComputeL2Diff() 322 @*/ 323 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 324 { 325 Vec localX; 326 PetscErrorCode ierr; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 330 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 331 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr); 332 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 333 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 334 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 335 PetscFunctionReturn(0); 336 } 337 338 #undef __FUNCT__ 339 #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM" 340 PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX) 341 { 342 void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 343 void **ctxs; 344 PetscFE *fe; 345 PetscInt numFields, f, numBd, b; 346 PetscErrorCode ierr; 347 348 PetscFunctionBegin; 349 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 350 PetscValidHeaderSpecific(localX, VEC_CLASSID, 2); 351 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 352 ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 353 for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);} 354 /* OPT: Could attempt to do multiple BCs at once */ 355 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 356 for (b = 0; b < numBd; ++b) { 357 DMLabel label; 358 const PetscInt *ids; 359 const char *name; 360 PetscInt numids, field; 361 PetscBool isEssential; 362 void (*func)(); 363 void *ctx; 364 365 /* TODO: We need to set only the part indicated by the ids */ 366 ierr = DMPlexGetBoundary(dm, b, &isEssential, &name, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 367 ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr); 368 for (f = 0; f < numFields; ++f) { 369 funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL; 370 ctxs[f] = field == f ? ctx : NULL; 371 } 372 ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 373 } 374 ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 /* Assuming dim == 3 */ 379 typedef struct { 380 PetscScalar normal[3]; /* Area-scaled normals */ 381 PetscScalar centroid[3]; /* Location of centroid (quadrature point) */ 382 PetscScalar grad[2][3]; /* Face contribution to gradient in left and right cell */ 383 } FaceGeom; 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM" 387 PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscReal time, Vec locX) 388 { 389 DM dmFace; 390 Vec faceGeometry; 391 DMLabel label; 392 const PetscScalar *facegeom; 393 PetscScalar *x; 394 PetscInt numBd, b, fStart, fEnd; 395 PetscErrorCode ierr; 396 397 PetscFunctionBegin; 398 /* TODO Pull this geometry calculation into the library */ 399 ierr = PetscObjectQuery((PetscObject) dm, "FaceGeometry", (PetscObject *) &faceGeometry);CHKERRQ(ierr); 400 ierr = DMPlexGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 401 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 402 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 403 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 404 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 405 ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 406 for (b = 0; b < numBd; ++b) { 407 PetscErrorCode (*func)(PetscReal,const PetscScalar*,const PetscScalar*,const PetscScalar*,PetscScalar*,void*); 408 const PetscInt *ids; 409 PetscInt numids, i; 410 void *ctx; 411 412 ierr = DMPlexGetBoundary(dm, b, NULL, NULL, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr); 413 for (i = 0; i < numids; ++i) { 414 IS faceIS; 415 const PetscInt *faces; 416 PetscInt numFaces, f; 417 418 ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 419 if (!faceIS) continue; /* No points with that id on this process */ 420 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 421 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 422 for (f = 0; f < numFaces; ++f) { 423 const PetscInt face = faces[f], *cells; 424 const PetscScalar *xI; 425 PetscScalar *xG; 426 const FaceGeom *fg; 427 428 if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 429 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 430 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 431 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 432 ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 433 ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 434 } 435 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 436 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 437 } 438 } 439 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 440 ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "DMPlexComputeL2Diff" 446 /*@C 447 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 448 449 Input Parameters: 450 + dm - The DM 451 . fe - The PetscFE object for each field 452 . funcs - The functions to evaluate for each field component 453 . ctxs - Optional array of contexts to pass to each function, or NULL. 454 - X - The coefficient vector u_h 455 456 Output Parameter: 457 . diff - The diff ||u - u_h||_2 458 459 Level: developer 460 461 .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff() 462 @*/ 463 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 464 { 465 const PetscInt debug = 0; 466 PetscSection section; 467 PetscQuadrature quad; 468 Vec localX; 469 PetscScalar *funcVal; 470 PetscReal *coords, *v0, *J, *invJ, detJ; 471 PetscReal localDiff = 0.0; 472 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 473 PetscErrorCode ierr; 474 475 PetscFunctionBegin; 476 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 477 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 478 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 479 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 480 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 481 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 482 for (field = 0; field < numFields; ++field) { 483 PetscInt Nc; 484 485 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 486 numComponents += Nc; 487 } 488 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 489 ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 490 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 491 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 492 for (c = cStart; c < cEnd; ++c) { 493 PetscScalar *x = NULL; 494 PetscReal elemDiff = 0.0; 495 496 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 497 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 498 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 499 500 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 501 void * const ctx = ctxs ? ctxs[field] : NULL; 502 const PetscReal *quadPoints, *quadWeights; 503 PetscReal *basis; 504 PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 505 506 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 507 ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 508 ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 509 ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 510 if (debug) { 511 char title[1024]; 512 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 513 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 514 } 515 for (q = 0; q < numQuadPoints; ++q) { 516 for (d = 0; d < dim; d++) { 517 coords[d] = v0[d]; 518 for (e = 0; e < dim; e++) { 519 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 520 } 521 } 522 (*funcs[field])(coords, funcVal, ctx); 523 for (fc = 0; fc < numBasisComps; ++fc) { 524 PetscScalar interpolant = 0.0; 525 526 for (f = 0; f < numBasisFuncs; ++f) { 527 const PetscInt fidx = f*numBasisComps+fc; 528 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 529 } 530 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 531 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 532 } 533 } 534 comp += numBasisComps; 535 fieldOffset += numBasisFuncs*numBasisComps; 536 } 537 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 538 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 539 localDiff += elemDiff; 540 } 541 ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 542 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 543 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 544 *diff = PetscSqrtReal(*diff); 545 PetscFunctionReturn(0); 546 } 547 548 #undef __FUNCT__ 549 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 550 /*@C 551 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 552 553 Input Parameters: 554 + dm - The DM 555 . fe - The PetscFE object for each field 556 . funcs - The gradient functions to evaluate for each field component 557 . ctxs - Optional array of contexts to pass to each function, or NULL. 558 . X - The coefficient vector u_h 559 - n - The vector to project along 560 561 Output Parameter: 562 . diff - The diff ||(grad u - grad u_h) . n||_2 563 564 Level: developer 565 566 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 567 @*/ 568 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 569 { 570 const PetscInt debug = 0; 571 PetscSection section; 572 PetscQuadrature quad; 573 Vec localX; 574 PetscScalar *funcVal, *interpolantVec; 575 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 576 PetscReal localDiff = 0.0; 577 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 578 PetscErrorCode ierr; 579 580 PetscFunctionBegin; 581 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 582 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 583 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 584 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 585 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 586 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 587 for (field = 0; field < numFields; ++field) { 588 PetscInt Nc; 589 590 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 591 numComponents += Nc; 592 } 593 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 594 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 595 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 596 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 597 for (c = cStart; c < cEnd; ++c) { 598 PetscScalar *x = NULL; 599 PetscReal elemDiff = 0.0; 600 601 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 602 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 603 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 604 605 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 606 void * const ctx = ctxs ? ctxs[field] : NULL; 607 const PetscReal *quadPoints, *quadWeights; 608 PetscReal *basisDer; 609 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 610 611 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 612 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 613 ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr); 614 ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr); 615 if (debug) { 616 char title[1024]; 617 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 618 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 619 } 620 for (q = 0; q < numQuadPoints; ++q) { 621 for (d = 0; d < dim; d++) { 622 coords[d] = v0[d]; 623 for (e = 0; e < dim; e++) { 624 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 625 } 626 } 627 (*funcs[field])(coords, n, funcVal, ctx); 628 for (fc = 0; fc < Ncomp; ++fc) { 629 PetscScalar interpolant = 0.0; 630 631 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 632 for (f = 0; f < Nb; ++f) { 633 const PetscInt fidx = f*Ncomp+fc; 634 635 for (d = 0; d < dim; ++d) { 636 realSpaceDer[d] = 0.0; 637 for (g = 0; g < dim; ++g) { 638 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 639 } 640 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 641 } 642 } 643 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 644 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 645 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 646 } 647 } 648 comp += Ncomp; 649 fieldOffset += Nb*Ncomp; 650 } 651 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 652 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 653 localDiff += elemDiff; 654 } 655 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 656 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 657 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 658 *diff = PetscSqrtReal(*diff); 659 PetscFunctionReturn(0); 660 } 661 662 #undef __FUNCT__ 663 #define __FUNCT__ "DMPlexComputeResidualFEM" 664 /*@ 665 DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 666 667 Input Parameters: 668 + dm - The mesh 669 . X - Local input vector 670 - user - The user context 671 672 Output Parameter: 673 . F - Local output vector 674 675 Note: 676 The first member of the user context must be an FEMContext. 677 678 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 679 like a GPU, or vectorize on a multicore machine. 680 681 Level: developer 682 683 .seealso: DMPlexComputeJacobianActionFEM() 684 @*/ 685 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 686 { 687 DM_Plex *mesh = (DM_Plex *) dm->data; 688 PetscFEM *fem = (PetscFEM *) user; 689 PetscFE *fe = fem->fe; 690 PetscFE *feAux = fem->feAux; 691 PetscFE *feBd = fem->feBd; 692 const char *name = "Residual"; 693 DM dmAux; 694 Vec A; 695 PetscQuadrature q; 696 PetscCellGeometry geom; 697 PetscSection section, sectionAux; 698 PetscReal *v0, *J, *invJ, *detJ; 699 PetscScalar *elemVec, *u, *a = NULL; 700 PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 701 PetscInt cellDof = 0, numComponents = 0; 702 PetscInt cellDofAux = 0, numComponentsAux = 0; 703 PetscErrorCode ierr; 704 705 PetscFunctionBegin; 706 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 707 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 708 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 709 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 710 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 711 numCells = cEnd - cStart; 712 for (f = 0; f < Nf; ++f) { 713 PetscInt Nb, Nc; 714 715 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 716 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 717 cellDof += Nb*Nc; 718 numComponents += Nc; 719 } 720 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 721 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 722 if (dmAux) { 723 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 724 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 725 } 726 for (f = 0; f < NfAux; ++f) { 727 PetscInt Nb, Nc; 728 729 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 730 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 731 cellDofAux += Nb*Nc; 732 numComponentsAux += Nc; 733 } 734 ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr); 735 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 736 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 737 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 738 for (c = cStart; c < cEnd; ++c) { 739 PetscScalar *x = NULL; 740 PetscInt i; 741 742 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 743 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 744 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 745 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 746 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 747 if (dmAux) { 748 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 749 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 750 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 751 } 752 } 753 for (f = 0; f < Nf; ++f) { 754 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 755 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 756 PetscInt numQuadPoints, Nb; 757 /* Conforming batches */ 758 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 759 /* Remainder */ 760 PetscInt Nr, offset; 761 762 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 763 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 764 ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 765 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 766 blockSize = Nb*numQuadPoints; 767 batchSize = numBlocks * blockSize; 768 ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 769 numChunks = numCells / (numBatches*batchSize); 770 Ne = numChunks*numBatches*batchSize; 771 Nr = numCells % (numBatches*batchSize); 772 offset = numCells - Nr; 773 geom.v0 = v0; 774 geom.J = J; 775 geom.invJ = invJ; 776 geom.detJ = detJ; 777 ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 778 geom.v0 = &v0[offset*dim]; 779 geom.J = &J[offset*dim*dim]; 780 geom.invJ = &invJ[offset*dim*dim]; 781 geom.detJ = &detJ[offset]; 782 ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 783 } 784 for (c = cStart; c < cEnd; ++c) { 785 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 786 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 787 } 788 ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 789 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 790 if (feBd) { 791 DMLabel depth; 792 PetscInt numBd, bd; 793 794 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 795 PetscInt Nb, Nc; 796 797 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 798 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 799 cellDof += Nb*Nc; 800 numComponents += Nc; 801 } 802 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 803 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 804 for (bd = 0; bd < numBd; ++bd) { 805 const char *bdLabel; 806 DMLabel label; 807 IS pointIS; 808 const PetscInt *points; 809 const PetscInt *values; 810 PetscReal *n; 811 PetscInt field, numValues, numPoints, p, dep, numFaces; 812 PetscBool isEssential; 813 814 ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 815 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 816 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 817 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 818 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 819 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 820 for (p = 0, numFaces = 0; p < numPoints; ++p) { 821 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 822 if (dep == dim-1) ++numFaces; 823 } 824 ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr); 825 for (p = 0, f = 0; p < numPoints; ++p) { 826 const PetscInt point = points[p]; 827 PetscScalar *x = NULL; 828 PetscInt i; 829 830 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 831 if (dep != dim-1) continue; 832 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 833 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 834 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 835 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 836 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 837 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 838 ++f; 839 } 840 for (f = 0; f < Nf; ++f) { 841 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 842 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 843 PetscInt numQuadPoints, Nb; 844 /* Conforming batches */ 845 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 846 /* Remainder */ 847 PetscInt Nr, offset; 848 849 ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 850 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 851 ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 852 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 853 blockSize = Nb*numQuadPoints; 854 batchSize = numBlocks * blockSize; 855 ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 856 numChunks = numFaces / (numBatches*batchSize); 857 Ne = numChunks*numBatches*batchSize; 858 Nr = numFaces % (numBatches*batchSize); 859 offset = numFaces - Nr; 860 geom.v0 = v0; 861 geom.n = n; 862 geom.J = J; 863 geom.invJ = invJ; 864 geom.detJ = detJ; 865 ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 866 geom.v0 = &v0[offset*dim]; 867 geom.n = &n[offset*dim]; 868 geom.J = &J[offset*dim*dim]; 869 geom.invJ = &invJ[offset*dim*dim]; 870 geom.detJ = &detJ[offset]; 871 ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 872 } 873 for (p = 0, f = 0; p < numPoints; ++p) { 874 const PetscInt point = points[p]; 875 876 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 877 if (dep != dim-1) continue; 878 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 879 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 880 ++f; 881 } 882 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 883 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 884 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 885 } 886 } 887 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 888 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 889 PetscFunctionReturn(0); 890 } 891 892 #undef __FUNCT__ 893 #define __FUNCT__ "DMPlexComputeIFunctionFEM" 894 /*@ 895 DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user 896 897 Input Parameters: 898 + dm - The mesh 899 . time - The current time 900 . X - Local input vector 901 . X_t - Time derivative of the local input vector 902 - user - The user context 903 904 Output Parameter: 905 . F - Local output vector 906 907 Note: 908 The first member of the user context must be an FEMContext. 909 910 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 911 like a GPU, or vectorize on a multicore machine. 912 913 Level: developer 914 915 .seealso: DMPlexComputeResidualFEM() 916 @*/ 917 PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user) 918 { 919 DM_Plex *mesh = (DM_Plex *) dm->data; 920 PetscFEM *fem = (PetscFEM *) user; 921 PetscFE *fe = fem->fe; 922 PetscFE *feAux = fem->feAux; 923 PetscFE *feBd = fem->feBd; 924 const char *name = "Residual"; 925 DM dmAux; 926 Vec A; 927 PetscQuadrature q; 928 PetscCellGeometry geom; 929 PetscSection section, sectionAux; 930 PetscReal *v0, *J, *invJ, *detJ; 931 PetscScalar *elemVec, *u, *u_t, *a = NULL; 932 PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 933 PetscInt cellDof = 0, numComponents = 0; 934 PetscInt cellDofAux = 0, numComponentsAux = 0; 935 PetscErrorCode ierr; 936 937 PetscFunctionBegin; 938 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 939 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 940 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 941 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 942 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 943 numCells = cEnd - cStart; 944 for (f = 0; f < Nf; ++f) { 945 PetscInt Nb, Nc; 946 947 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 948 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 949 cellDof += Nb*Nc; 950 numComponents += Nc; 951 } 952 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 953 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 954 if (dmAux) { 955 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 956 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 957 } 958 for (f = 0; f < NfAux; ++f) { 959 PetscInt Nb, Nc; 960 961 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 962 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 963 cellDofAux += Nb*Nc; 964 numComponentsAux += Nc; 965 } 966 ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr); 967 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 968 ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 969 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 970 for (c = cStart; c < cEnd; ++c) { 971 PetscScalar *x = NULL, *x_t = NULL; 972 PetscInt i; 973 974 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 975 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 976 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 977 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 978 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 979 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 980 for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i]; 981 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 982 if (dmAux) { 983 PetscScalar *x_a = NULL; 984 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr); 985 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i]; 986 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr); 987 } 988 } 989 for (f = 0; f < Nf; ++f) { 990 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f]; 991 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f]; 992 PetscInt numQuadPoints, Nb; 993 /* Conforming batches */ 994 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 995 /* Remainder */ 996 PetscInt Nr, offset; 997 998 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 999 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 1000 ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1001 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1002 blockSize = Nb*numQuadPoints; 1003 batchSize = numBlocks * blockSize; 1004 ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1005 numChunks = numCells / (numBatches*batchSize); 1006 Ne = numChunks*numBatches*batchSize; 1007 Nr = numCells % (numBatches*batchSize); 1008 offset = numCells - Nr; 1009 geom.v0 = v0; 1010 geom.J = J; 1011 geom.invJ = invJ; 1012 geom.detJ = detJ; 1013 ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 1014 geom.v0 = &v0[offset*dim]; 1015 geom.J = &J[offset*dim*dim]; 1016 geom.invJ = &invJ[offset*dim*dim]; 1017 geom.detJ = &detJ[offset]; 1018 ierr = PetscFEIntegrateIFunction(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 1019 } 1020 for (c = cStart; c < cEnd; ++c) { 1021 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 1022 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 1023 } 1024 ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 1025 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1026 if (feBd) { 1027 DMLabel label, depth; 1028 IS pointIS; 1029 const PetscInt *points; 1030 PetscInt dep, numPoints, p, numFaces; 1031 PetscReal *n; 1032 1033 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 1034 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1035 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1036 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1037 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1038 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 1039 PetscInt Nb, Nc; 1040 1041 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 1042 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 1043 cellDof += Nb*Nc; 1044 numComponents += Nc; 1045 } 1046 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1047 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1048 if (dep == dim-1) ++numFaces; 1049 } 1050 ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr); 1051 ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr); 1052 for (p = 0, f = 0; p < numPoints; ++p) { 1053 const PetscInt point = points[p]; 1054 PetscScalar *x = NULL; 1055 PetscInt i; 1056 1057 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1058 if (dep != dim-1) continue; 1059 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 1060 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 1061 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 1062 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1063 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 1064 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1065 ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1066 for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i]; 1067 ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1068 ++f; 1069 } 1070 for (f = 0; f < Nf; ++f) { 1071 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f]; 1072 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f]; 1073 PetscInt numQuadPoints, Nb; 1074 /* Conforming batches */ 1075 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1076 /* Remainder */ 1077 PetscInt Nr, offset; 1078 1079 ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 1080 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 1081 ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1082 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1083 blockSize = Nb*numQuadPoints; 1084 batchSize = numBlocks * blockSize; 1085 ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1086 numChunks = numFaces / (numBatches*batchSize); 1087 Ne = numChunks*numBatches*batchSize; 1088 Nr = numFaces % (numBatches*batchSize); 1089 offset = numFaces - Nr; 1090 geom.v0 = v0; 1091 geom.n = n; 1092 geom.J = J; 1093 geom.invJ = invJ; 1094 geom.detJ = detJ; 1095 ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 1096 geom.v0 = &v0[offset*dim]; 1097 geom.n = &n[offset*dim]; 1098 geom.J = &J[offset*dim*dim]; 1099 geom.invJ = &invJ[offset*dim*dim]; 1100 geom.detJ = &detJ[offset]; 1101 ierr = PetscFEIntegrateBdIFunction(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 1102 } 1103 for (p = 0, f = 0; p < numPoints; ++p) { 1104 const PetscInt point = points[p]; 1105 1106 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1107 if (dep != dim-1) continue; 1108 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 1109 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 1110 ++f; 1111 } 1112 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1113 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1114 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 1115 ierr = PetscFree(u_t);CHKERRQ(ierr); 1116 } 1117 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 1118 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1119 PetscFunctionReturn(0); 1120 } 1121 1122 #undef __FUNCT__ 1123 #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 1124 /*@C 1125 DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 1126 1127 Input Parameters: 1128 + dm - The mesh 1129 . J - The Jacobian shell matrix 1130 . X - Local input vector 1131 - user - The user context 1132 1133 Output Parameter: 1134 . F - Local output vector 1135 1136 Note: 1137 The first member of the user context must be an FEMContext. 1138 1139 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1140 like a GPU, or vectorize on a multicore machine. 1141 1142 Level: developer 1143 1144 .seealso: DMPlexComputeResidualFEM() 1145 @*/ 1146 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 1147 { 1148 DM_Plex *mesh = (DM_Plex *) dm->data; 1149 PetscFEM *fem = (PetscFEM *) user; 1150 PetscFE *fe = fem->fe; 1151 PetscQuadrature quad; 1152 PetscCellGeometry geom; 1153 PetscSection section; 1154 JacActionCtx *jctx; 1155 PetscReal *v0, *J, *invJ, *detJ; 1156 PetscScalar *elemVec, *u, *a; 1157 PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 1158 PetscInt cellDof = 0; 1159 PetscErrorCode ierr; 1160 1161 PetscFunctionBegin; 1162 /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 1163 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 1164 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1165 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1166 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 1167 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1168 numCells = cEnd - cStart; 1169 for (field = 0; field < numFields; ++field) { 1170 PetscInt Nb, Nc; 1171 1172 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 1173 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 1174 cellDof += Nb*Nc; 1175 } 1176 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 1177 ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&a,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 1178 for (c = cStart; c < cEnd; ++c) { 1179 PetscScalar *x = NULL; 1180 PetscInt i; 1181 1182 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1183 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1184 ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 1185 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 1186 ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 1187 ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 1188 for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 1189 ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 1190 } 1191 for (field = 0; field < numFields; ++field) { 1192 PetscInt numQuadPoints, Nb; 1193 /* Conforming batches */ 1194 PetscInt numBlocks = 1; 1195 PetscInt numBatches = 1; 1196 PetscInt numChunks, Ne, blockSize, batchSize; 1197 /* Remainder */ 1198 PetscInt Nr, offset; 1199 1200 ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 1201 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 1202 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1203 blockSize = Nb*numQuadPoints; 1204 batchSize = numBlocks * blockSize; 1205 numChunks = numCells / (numBatches*batchSize); 1206 Ne = numChunks*numBatches*batchSize; 1207 Nr = numCells % (numBatches*batchSize); 1208 offset = numCells - Nr; 1209 geom.v0 = v0; 1210 geom.J = J; 1211 geom.invJ = invJ; 1212 geom.detJ = detJ; 1213 ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 1214 geom.v0 = &v0[offset*dim]; 1215 geom.J = &J[offset*dim*dim]; 1216 geom.invJ = &invJ[offset*dim*dim]; 1217 geom.detJ = &detJ[offset]; 1218 ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 1219 fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 1220 } 1221 for (c = cStart; c < cEnd; ++c) { 1222 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 1223 ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 1224 } 1225 ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 1226 if (mesh->printFEM) { 1227 PetscMPIInt rank, numProcs; 1228 PetscInt p; 1229 1230 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 1231 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 1232 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 1233 for (p = 0; p < numProcs; ++p) { 1234 if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 1235 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 1236 } 1237 } 1238 /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 1239 PetscFunctionReturn(0); 1240 } 1241 1242 #undef __FUNCT__ 1243 #define __FUNCT__ "DMPlexComputeJacobianFEM" 1244 /*@ 1245 DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 1246 1247 Input Parameters: 1248 + dm - The mesh 1249 . X - Local input vector 1250 - user - The user context 1251 1252 Output Parameter: 1253 . Jac - Jacobian matrix 1254 1255 Note: 1256 The first member of the user context must be an FEMContext. 1257 1258 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1259 like a GPU, or vectorize on a multicore machine. 1260 1261 Level: developer 1262 1263 .seealso: FormFunctionLocal() 1264 @*/ 1265 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 1266 { 1267 DM_Plex *mesh = (DM_Plex *) dm->data; 1268 PetscFEM *fem = (PetscFEM *) user; 1269 PetscFE *fe = fem->fe; 1270 PetscFE *feAux = fem->feAux; 1271 PetscFE *feBd = fem->feBd; 1272 const char *name = "Jacobian"; 1273 DM dmAux; 1274 Vec A; 1275 PetscQuadrature quad; 1276 PetscCellGeometry geom; 1277 PetscSection section, globalSection, sectionAux; 1278 PetscReal *v0, *J, *invJ, *detJ; 1279 PetscScalar *elemMat, *u, *a; 1280 PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 1281 PetscInt cellDof = 0, numComponents = 0; 1282 PetscInt cellDofAux = 0, numComponentsAux = 0; 1283 PetscBool isShell; 1284 PetscErrorCode ierr; 1285 1286 PetscFunctionBegin; 1287 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 1288 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1289 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1290 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 1291 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1292 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1293 numCells = cEnd - cStart; 1294 for (f = 0; f < Nf; ++f) { 1295 PetscInt Nb, Nc; 1296 1297 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 1298 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1299 cellDof += Nb*Nc; 1300 numComponents += Nc; 1301 } 1302 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1303 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1304 if (dmAux) { 1305 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1306 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 1307 } 1308 for (f = 0; f < NfAux; ++f) { 1309 PetscInt Nb, Nc; 1310 1311 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 1312 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 1313 cellDofAux += Nb*Nc; 1314 numComponentsAux += Nc; 1315 } 1316 ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr); 1317 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 1318 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 1319 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 1320 for (c = cStart; c < cEnd; ++c) { 1321 PetscScalar *x = NULL; 1322 PetscInt i; 1323 1324 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1325 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1326 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1327 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 1328 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1329 if (dmAux) { 1330 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1331 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 1332 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1333 } 1334 } 1335 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1336 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1337 PetscInt numQuadPoints, Nb; 1338 /* Conforming batches */ 1339 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1340 /* Remainder */ 1341 PetscInt Nr, offset; 1342 1343 ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 1344 ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 1345 ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1346 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1347 blockSize = Nb*numQuadPoints; 1348 batchSize = numBlocks * blockSize; 1349 ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1350 numChunks = numCells / (numBatches*batchSize); 1351 Ne = numChunks*numBatches*batchSize; 1352 Nr = numCells % (numBatches*batchSize); 1353 offset = numCells - Nr; 1354 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1355 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 1356 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 1357 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 1358 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 1359 1360 geom.v0 = v0; 1361 geom.J = J; 1362 geom.invJ = invJ; 1363 geom.detJ = detJ; 1364 ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 1365 geom.v0 = &v0[offset*dim]; 1366 geom.J = &J[offset*dim*dim]; 1367 geom.invJ = &invJ[offset*dim*dim]; 1368 geom.detJ = &detJ[offset]; 1369 ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, Nf, fe, fieldI, fieldJ, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 1370 } 1371 } 1372 for (c = cStart; c < cEnd; ++c) { 1373 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 1374 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 1375 } 1376 ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1377 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1378 if (feBd) { 1379 DMLabel depth; 1380 PetscInt numBd, bd; 1381 1382 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 1383 PetscInt Nb, Nc; 1384 1385 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 1386 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 1387 cellDof += Nb*Nc; 1388 numComponents += Nc; 1389 } 1390 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1391 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1392 for (bd = 0; bd < numBd; ++bd) { 1393 const char *bdLabel; 1394 DMLabel label; 1395 IS pointIS; 1396 const PetscInt *points; 1397 const PetscInt *values; 1398 PetscReal *n; 1399 PetscInt field, numValues, numPoints, p, dep, numFaces; 1400 PetscBool isEssential; 1401 1402 ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1403 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 1404 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1405 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1406 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1407 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1408 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1409 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1410 if (dep == dim-1) ++numFaces; 1411 } 1412 ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 1413 for (p = 0, f = 0; p < numPoints; ++p) { 1414 const PetscInt point = points[p]; 1415 PetscScalar *x = NULL; 1416 PetscInt i; 1417 1418 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1419 if (dep != dim-1) continue; 1420 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 1421 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 1422 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 1423 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1424 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 1425 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1426 ++f; 1427 } 1428 ierr = PetscMemzero(elemMat, numFaces*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1429 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1430 PetscInt numQuadPoints, Nb; 1431 /* Conforming batches */ 1432 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1433 /* Remainder */ 1434 PetscInt Nr, offset; 1435 1436 ierr = PetscFEGetQuadrature(feBd[fieldI], &quad);CHKERRQ(ierr); 1437 ierr = PetscFEGetDimension(feBd[fieldI], &Nb);CHKERRQ(ierr); 1438 ierr = PetscFEGetTileSizes(feBd[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1439 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1440 blockSize = Nb*numQuadPoints; 1441 batchSize = numBlocks * blockSize; 1442 ierr = PetscFESetTileSizes(feBd[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1443 numChunks = numFaces / (numBatches*batchSize); 1444 Ne = numChunks*numBatches*batchSize; 1445 Nr = numFaces % (numBatches*batchSize); 1446 offset = numFaces - Nr; 1447 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1448 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g0BdFuncs[fieldI*Nf+fieldJ]; 1449 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g1BdFuncs[fieldI*Nf+fieldJ]; 1450 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g2BdFuncs[fieldI*Nf+fieldJ]; 1451 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g3BdFuncs[fieldI*Nf+fieldJ]; 1452 1453 geom.v0 = v0; 1454 geom.n = n; 1455 geom.J = J; 1456 geom.invJ = invJ; 1457 geom.detJ = detJ; 1458 ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Ne, Nf, feBd, fieldI, fieldJ, geom, u, 0, NULL, NULL, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 1459 geom.v0 = &v0[offset*dim]; 1460 geom.n = &n[offset*dim]; 1461 geom.J = &J[offset*dim*dim]; 1462 geom.invJ = &invJ[offset*dim*dim]; 1463 geom.detJ = &detJ[offset]; 1464 ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Nr, Nf, feBd, fieldI, fieldJ, geom, &u[offset*cellDof], 0, NULL, NULL, g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 1465 } 1466 } 1467 for (p = 0, f = 0; p < numPoints; ++p) { 1468 const PetscInt point = points[p]; 1469 1470 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1471 if (dep != dim-1) continue; 1472 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", cellDof, cellDof, &elemMat[f*cellDof*cellDof]);CHKERRQ(ierr);} 1473 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 1474 ++f; 1475 } 1476 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1477 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1478 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1479 } 1480 } 1481 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1482 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1483 if (mesh->printFEM) { 1484 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1485 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 1486 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1487 } 1488 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 1489 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 1490 if (isShell) { 1491 JacActionCtx *jctx; 1492 1493 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 1494 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 1495 } 1496 PetscFunctionReturn(0); 1497 } 1498 1499 #undef __FUNCT__ 1500 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1501 /*@ 1502 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1503 1504 Input Parameters: 1505 + dmf - The fine mesh 1506 . dmc - The coarse mesh 1507 - user - The user context 1508 1509 Output Parameter: 1510 . In - The interpolation matrix 1511 1512 Note: 1513 The first member of the user context must be an FEMContext. 1514 1515 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1516 like a GPU, or vectorize on a multicore machine. 1517 1518 Level: developer 1519 1520 .seealso: DMPlexComputeJacobianFEM() 1521 @*/ 1522 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1523 { 1524 DM_Plex *mesh = (DM_Plex *) dmc->data; 1525 PetscFEM *fem = (PetscFEM *) user; 1526 PetscFE *fe = fem->fe; 1527 const char *name = "Interpolator"; 1528 PetscFE *feRef; 1529 PetscSection fsection, fglobalSection; 1530 PetscSection csection, cglobalSection; 1531 PetscScalar *elemMat; 1532 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 1533 PetscInt rCellDof = 0, cCellDof = 0; 1534 PetscErrorCode ierr; 1535 1536 PetscFunctionBegin; 1537 #if 0 1538 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1539 #endif 1540 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1541 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1542 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1543 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1544 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1545 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1546 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1547 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1548 for (f = 0; f < Nf; ++f) { 1549 PetscInt rNb, cNb, Nc; 1550 1551 ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr); 1552 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1553 ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr); 1554 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1555 rCellDof += rNb*Nc; 1556 cCellDof += cNb*Nc; 1557 } 1558 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1559 ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr); 1560 ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1561 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1562 PetscDualSpace Qref; 1563 PetscQuadrature f; 1564 const PetscReal *qpoints, *qweights; 1565 PetscReal *points; 1566 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1567 1568 /* Compose points from all dual basis functionals */ 1569 ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr); 1570 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1571 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1572 for (i = 0; i < fpdim; ++i) { 1573 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1574 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1575 npoints += Np; 1576 } 1577 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1578 for (i = 0, k = 0; i < fpdim; ++i) { 1579 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1580 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1581 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1582 } 1583 1584 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1585 PetscReal *B; 1586 PetscInt NcJ, cpdim, j; 1587 1588 /* Evaluate basis at points */ 1589 ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr); 1590 if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ); 1591 ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr); 1592 /* For now, fields only interpolate themselves */ 1593 if (fieldI == fieldJ) { 1594 ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1595 for (i = 0, k = 0; i < fpdim; ++i) { 1596 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1597 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1598 for (p = 0; p < Np; ++p, ++k) { 1599 for (j = 0; j < cpdim; ++j) { 1600 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cCellDof + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p]; 1601 } 1602 } 1603 } 1604 ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1605 } 1606 offsetJ += cpdim*NcJ; 1607 } 1608 offsetI += fpdim*Nc; 1609 ierr = PetscFree(points);CHKERRQ(ierr); 1610 } 1611 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);} 1612 for (c = cStart; c < cEnd; ++c) { 1613 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1614 } 1615 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1616 ierr = PetscFree(feRef);CHKERRQ(ierr); 1617 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1618 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1619 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1620 if (mesh->printFEM) { 1621 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1622 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1623 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1624 } 1625 #if 0 1626 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1627 #endif 1628 PetscFunctionReturn(0); 1629 } 1630 1631 #undef __FUNCT__ 1632 #define __FUNCT__ "DMPlexAddBoundary" 1633 /* The ids can be overridden by the command line option -bc_<boundary name> */ 1634 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx) 1635 { 1636 DM_Plex *mesh = (DM_Plex *) dm->data; 1637 DMBoundary b; 1638 PetscErrorCode ierr; 1639 1640 PetscFunctionBegin; 1641 ierr = PetscNew(&b);CHKERRQ(ierr); 1642 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 1643 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 1644 ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr); 1645 b->essential = isEssential; 1646 b->field = field; 1647 b->func = bcFunc; 1648 b->numids = numids; 1649 b->ctx = ctx; 1650 b->next = mesh->boundary; 1651 mesh->boundary = b; 1652 PetscFunctionReturn(0); 1653 } 1654 1655 #undef __FUNCT__ 1656 #define __FUNCT__ "DMPlexGetNumBoundary" 1657 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd) 1658 { 1659 DM_Plex *mesh = (DM_Plex *) dm->data; 1660 DMBoundary b = mesh->boundary; 1661 1662 PetscFunctionBegin; 1663 *numBd = 0; 1664 while (b) {++(*numBd); b = b->next;} 1665 PetscFunctionReturn(0); 1666 } 1667 1668 #undef __FUNCT__ 1669 #define __FUNCT__ "DMPlexGetBoundary" 1670 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx) 1671 { 1672 DM_Plex *mesh = (DM_Plex *) dm->data; 1673 DMBoundary b = mesh->boundary; 1674 PetscInt n = 0; 1675 1676 PetscFunctionBegin; 1677 while (b) { 1678 if (n == bd) break; 1679 b = b->next; 1680 ++n; 1681 } 1682 if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 1683 if (isEssential) { 1684 PetscValidPointer(isEssential, 3); 1685 *isEssential = b->essential; 1686 } 1687 if (name) { 1688 PetscValidPointer(name, 4); 1689 *name = b->name; 1690 } 1691 if (field) { 1692 PetscValidPointer(field, 5); 1693 *field = b->field; 1694 } 1695 if (func) { 1696 PetscValidPointer(func, 6); 1697 *func = b->func; 1698 } 1699 if (numids) { 1700 PetscValidPointer(numids, 7); 1701 *numids = b->numids; 1702 } 1703 if (ids) { 1704 PetscValidPointer(ids, 8); 1705 *ids = b->ids; 1706 } 1707 if (ctx) { 1708 PetscValidPointer(ctx, 9); 1709 *ctx = b->ctx; 1710 } 1711 PetscFunctionReturn(0); 1712 } 1713