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