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