1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscsf.h> 3 4 #include <petsc-private/petscfeimpl.h> 5 #include <petsc-private/petscfvimpl.h> 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "DMPlexGetScale" 9 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 10 { 11 DM_Plex *mesh = (DM_Plex*) dm->data; 12 13 PetscFunctionBegin; 14 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15 PetscValidPointer(scale, 3); 16 *scale = mesh->scale[unit]; 17 PetscFunctionReturn(0); 18 } 19 20 #undef __FUNCT__ 21 #define __FUNCT__ "DMPlexSetScale" 22 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 23 { 24 DM_Plex *mesh = (DM_Plex*) dm->data; 25 26 PetscFunctionBegin; 27 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 28 mesh->scale[unit] = scale; 29 PetscFunctionReturn(0); 30 } 31 32 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 33 { 34 switch (i) { 35 case 0: 36 switch (j) { 37 case 0: return 0; 38 case 1: 39 switch (k) { 40 case 0: return 0; 41 case 1: return 0; 42 case 2: return 1; 43 } 44 case 2: 45 switch (k) { 46 case 0: return 0; 47 case 1: return -1; 48 case 2: return 0; 49 } 50 } 51 case 1: 52 switch (j) { 53 case 0: 54 switch (k) { 55 case 0: return 0; 56 case 1: return 0; 57 case 2: return -1; 58 } 59 case 1: return 0; 60 case 2: 61 switch (k) { 62 case 0: return 1; 63 case 1: return 0; 64 case 2: return 0; 65 } 66 } 67 case 2: 68 switch (j) { 69 case 0: 70 switch (k) { 71 case 0: return 0; 72 case 1: return 1; 73 case 2: return 0; 74 } 75 case 1: 76 switch (k) { 77 case 0: return -1; 78 case 1: return 0; 79 case 2: return 0; 80 } 81 case 2: return 0; 82 } 83 } 84 return 0; 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "DMPlexCreateRigidBody" 89 /*@C 90 DMPlexCreateRigidBody - create rigid body modes from coordinates 91 92 Collective on DM 93 94 Input Arguments: 95 + dm - the DM 96 . section - the local section associated with the rigid field, or NULL for the default section 97 - globalSection - the global section associated with the rigid field, or NULL for the default section 98 99 Output Argument: 100 . sp - the null space 101 102 Note: This is necessary to take account of Dirichlet conditions on the displacements 103 104 Level: advanced 105 106 .seealso: MatNullSpaceCreate() 107 @*/ 108 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 109 { 110 MPI_Comm comm; 111 Vec coordinates, localMode, mode[6]; 112 PetscSection coordSection; 113 PetscScalar *coords; 114 PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 115 PetscErrorCode ierr; 116 117 PetscFunctionBegin; 118 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 119 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 120 if (dim == 1) { 121 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 125 if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 126 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 127 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 128 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 129 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 130 m = (dim*(dim+1))/2; 131 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 132 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 133 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 134 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 135 /* Assume P1 */ 136 ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 137 for (d = 0; d < dim; ++d) { 138 PetscScalar values[3] = {0.0, 0.0, 0.0}; 139 140 values[d] = 1.0; 141 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 142 for (v = vStart; v < vEnd; ++v) { 143 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 144 } 145 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 146 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 147 } 148 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 149 for (d = dim; d < dim*(dim+1)/2; ++d) { 150 PetscInt i, j, k = dim > 2 ? d - dim : d; 151 152 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 153 for (v = vStart; v < vEnd; ++v) { 154 PetscScalar values[3] = {0.0, 0.0, 0.0}; 155 PetscInt off; 156 157 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 158 for (i = 0; i < dim; ++i) { 159 for (j = 0; j < dim; ++j) { 160 values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 161 } 162 } 163 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 164 } 165 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 166 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 167 } 168 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 169 ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 170 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 171 /* Orthonormalize system */ 172 for (i = dim; i < m; ++i) { 173 PetscScalar dots[6]; 174 175 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 176 for (j = 0; j < i; ++j) dots[j] *= -1.0; 177 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 178 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 179 } 180 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 181 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "DMPlexSetMaxProjectionHeight" 187 /*@ 188 DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs 189 are computed by associating the basis function with one of the mesh points in its transitively-closed support, and 190 evaluating the dual space basis of that point. A basis function is associated with the point in its 191 transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum 192 projection height, which is set with this function. By default, the maximum projection height is zero, which means 193 that only mesh cells are used to project basis functions. A height of one, for example, evaluates a cell-interior 194 basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face. 195 196 Input Parameters: 197 + dm - the DMPlex object 198 - height - the maximum projection height >= 0 199 200 Level: advanced 201 202 .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal() 203 @*/ 204 PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height) 205 { 206 DM_Plex *plex = (DM_Plex *) dm->data; 207 208 PetscFunctionBegin; 209 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 210 plex->maxProjectionHeight = height; 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "DMPlexGetMaxProjectionHeight" 216 /*@ 217 DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in 218 DMPlexProjectXXXLocal() functions. 219 220 Input Parameters: 221 . dm - the DMPlex object 222 223 Output Parameters: 224 . height - the maximum projection height 225 226 Level: intermediate 227 228 .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal() 229 @*/ 230 PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height) 231 { 232 DM_Plex *plex = (DM_Plex *) dm->data; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 236 *height = plex->maxProjectionHeight; 237 PetscFunctionReturn(0); 238 } 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal" 242 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 243 { 244 PetscDualSpace *sp, *cellsp; 245 PetscInt *numComp; 246 PetscSection section; 247 PetscScalar *values; 248 PetscBool *fieldActive; 249 PetscInt numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, cStart, cEnd, cEndInterior, f, d, v, i, comp, maxHeight, h; 250 PetscErrorCode ierr; 251 252 PetscFunctionBegin; 253 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 254 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 255 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 256 if (cEnd <= cStart) PetscFunctionReturn(0); 257 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 258 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 259 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 260 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 261 ierr = PetscMalloc2(numFields,&sp,numFields,&numComp);CHKERRQ(ierr); 262 ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr); 263 if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);} 264 if (maxHeight > 0) {ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);} 265 else {cellsp = sp;} 266 for (h = 0; h <= maxHeight; h++) { 267 ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 268 if (!h) {pStart = cStart; pEnd = cEnd;} 269 if (pEnd <= pStart) continue; 270 totDim = 0; 271 for (f = 0; f < numFields; ++f) { 272 PetscObject obj; 273 PetscClassId id; 274 275 ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 276 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 277 if (id == PETSCFE_CLASSID) { 278 PetscFE fe = (PetscFE) obj; 279 280 ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 281 if (!h) { 282 ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 283 sp[f] = cellsp[f]; 284 ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr); 285 } else { 286 ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr); 287 if (!sp[f]) continue; 288 } 289 } else if (id == PETSCFV_CLASSID) { 290 PetscFV fv = (PetscFV) obj; 291 PetscQuadrature q; 292 293 ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 294 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 295 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr); 296 ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr); 297 ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 298 ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr); 299 ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr); 300 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 301 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 302 totDim += spDim*numComp[f]; 303 } 304 ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr); 305 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim); 306 if (!totDim) continue; 307 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 308 ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 309 for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 310 for (i = 0; i < numIds; ++i) { 311 IS pointIS; 312 const PetscInt *points; 313 PetscInt n, p; 314 315 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 316 ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 317 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 318 for (p = 0; p < n; ++p) { 319 const PetscInt point = points[p]; 320 PetscFECellGeom geom; 321 322 if ((point < pStart) || (point >= pEnd)) continue; 323 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 324 geom.dim = dim - h; 325 geom.dimEmbed = dimEmbed; 326 for (f = 0, v = 0; f < numFields; ++f) { 327 void * const ctx = ctxs ? ctxs[f] : NULL; 328 329 if (!sp[f]) continue; 330 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 331 for (d = 0; d < spDim; ++d) { 332 if (funcs[f]) { 333 ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 334 } else { 335 for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0; 336 } 337 v += numComp[f]; 338 } 339 } 340 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 341 } 342 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 343 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 344 } 345 if (h) { 346 for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 347 } 348 } 349 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 350 ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 351 for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 352 ierr = PetscFree2(sp, numComp);CHKERRQ(ierr); 353 if (maxHeight > 0) { 354 ierr = PetscFree(cellsp);CHKERRQ(ierr); 355 } 356 PetscFunctionReturn(0); 357 } 358 359 #undef __FUNCT__ 360 #define __FUNCT__ "DMPlexProjectFunctionLocal" 361 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 362 { 363 PetscDualSpace *sp, *cellsp; 364 PetscInt *numComp; 365 PetscSection section; 366 PetscScalar *values; 367 PetscInt numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, f, d, v, comp, h, maxHeight; 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 372 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 373 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 374 if (cEnd <= cStart) PetscFunctionReturn(0); 375 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 376 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 377 ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr); 378 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 379 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 380 ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr); 381 if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);} 382 if (maxHeight > 0) { 383 ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr); 384 } 385 else { 386 cellsp = sp; 387 } 388 for (h = 0; h <= maxHeight; h++) { 389 ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 390 if (!h) {pStart = cStart; pEnd = cEnd;} 391 if (pEnd <= pStart) continue; 392 totDim = 0; 393 for (f = 0; f < numFields; ++f) { 394 PetscObject obj; 395 PetscClassId id; 396 397 ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 398 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 399 if (id == PETSCFE_CLASSID) { 400 PetscFE fe = (PetscFE) obj; 401 402 ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 403 if (!h) { 404 ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 405 sp[f] = cellsp[f]; 406 ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr); 407 } 408 else { 409 ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr); 410 if (!sp[f]) { 411 continue; 412 } 413 } 414 } else if (id == PETSCFV_CLASSID) { 415 PetscFV fv = (PetscFV) obj; 416 PetscQuadrature q; 417 418 if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume"); 419 ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 420 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 421 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr); 422 ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr); 423 ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 424 ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr); 425 ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr); 426 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 427 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 428 totDim += spDim*numComp[f]; 429 } 430 ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr); 431 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim); 432 if (!totDim) continue; 433 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 434 for (p = pStart; p < pEnd; ++p) { 435 PetscFECellGeom geom; 436 437 ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 438 geom.dim = dim - h; 439 geom.dimEmbed = dimEmbed; 440 for (f = 0, v = 0; f < numFields; ++f) { 441 void * const ctx = ctxs ? ctxs[f] : NULL; 442 443 if (!sp[f]) continue; 444 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 445 for (d = 0; d < spDim; ++d) { 446 if (funcs[f]) { 447 ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 448 } else { 449 for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0; 450 } 451 v += numComp[f]; 452 } 453 } 454 ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr); 455 } 456 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 457 if (h || !maxHeight) { 458 for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 459 } 460 } 461 ierr = PetscFree2(sp, numComp);CHKERRQ(ierr); 462 if (maxHeight > 0) { 463 for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);} 464 ierr = PetscFree(cellsp);CHKERRQ(ierr); 465 } 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "DMPlexProjectFunction" 471 /*@C 472 DMPlexProjectFunction - This projects the given function into the function space provided. 473 474 Input Parameters: 475 + dm - The DM 476 . funcs - The coordinate functions to evaluate, one per field 477 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 478 - mode - The insertion mode for values 479 480 Output Parameter: 481 . X - vector 482 483 Level: developer 484 485 .seealso: DMPlexComputeL2Diff() 486 @*/ 487 PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 488 { 489 Vec localX; 490 PetscErrorCode ierr; 491 492 PetscFunctionBegin; 493 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 494 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 495 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr); 496 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 497 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 498 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 499 PetscFunctionReturn(0); 500 } 501 502 #undef __FUNCT__ 503 #define __FUNCT__ "DMPlexProjectFieldLocal" 504 PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX) 505 { 506 DM dmAux; 507 PetscDS prob, probAux = NULL; 508 Vec A; 509 PetscSection section, sectionAux = NULL; 510 PetscScalar *values, *u, *u_x, *a, *a_x; 511 PetscReal *x, *v0, *J, *invJ, detJ; 512 PetscInt Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight; 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 /* TODO: Move out dual space extraction so that we do not create one for every cell (see ProjectFunctionLocal) */ 517 ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr); 518 if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");} 519 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 520 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 521 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 522 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 523 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 524 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 525 ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 526 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 527 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 528 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 529 if (dmAux) { 530 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 531 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 532 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 533 } 534 ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 535 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 536 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 537 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 538 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 539 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 540 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 541 for (c = cStart; c < cEnd; ++c) { 542 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 543 544 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 545 ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 546 if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 547 for (f = 0, v = 0; f < Nf; ++f) { 548 PetscObject obj; 549 PetscClassId id; 550 PetscDualSpace sp; 551 PetscInt Ncf; 552 553 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 554 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 555 if (id == PETSCFE_CLASSID) { 556 PetscFE fe = (PetscFE) obj; 557 558 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 559 ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr); 560 } else if (id == PETSCFV_CLASSID) { 561 PetscFV fv = (PetscFV) obj; 562 PetscQuadrature q; 563 564 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 565 ierr = PetscFVGetNumComponents(fv, &Ncf);CHKERRQ(ierr); 566 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp);CHKERRQ(ierr); 567 ierr = PetscDualSpaceSetDM(sp, dm);CHKERRQ(ierr); 568 ierr = PetscDualSpaceSetType(sp, PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 569 ierr = PetscDualSpaceSimpleSetDimension(sp, 1);CHKERRQ(ierr); 570 ierr = PetscDualSpaceSimpleSetFunctional(sp, 0, q);CHKERRQ(ierr); 571 } 572 ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 573 for (d = 0; d < spDim; ++d) { 574 PetscQuadrature quad; 575 const PetscReal *points, *weights; 576 PetscInt numPoints, q; 577 578 if (funcs[f]) { 579 ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr); 580 ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 581 for (q = 0; q < numPoints; ++q) { 582 CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x); 583 ierr = EvaluateFieldJets(prob, PETSC_FALSE, q, invJ, coefficients, NULL, u, u_x, NULL);CHKERRQ(ierr); 584 ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr); 585 (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]); 586 } 587 } else { 588 for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0; 589 } 590 v += Ncf; 591 } 592 if (id == PETSCFV_CLASSID) {ierr = PetscDualSpaceDestroy(&sp);CHKERRQ(ierr);} 593 } 594 ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 595 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 596 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 597 } 598 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 599 ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal" 605 static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX) 606 { 607 void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 608 void **ctxs; 609 PetscInt numFields; 610 PetscErrorCode ierr; 611 612 PetscFunctionBegin; 613 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 614 ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 615 funcs[field] = func; 616 ctxs[field] = ctx; 617 ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 618 ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 619 PetscFunctionReturn(0); 620 } 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal" 624 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, 625 PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX) 626 { 627 PetscDS prob; 628 PetscSF sf; 629 DM dmFace, dmCell, dmGrad; 630 const PetscScalar *facegeom, *cellgeom, *grad; 631 const PetscInt *leaves; 632 PetscScalar *x, *fx; 633 PetscInt dim, nleaves, loc, fStart, fEnd, pdim, i; 634 PetscErrorCode ierr; 635 636 PetscFunctionBegin; 637 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 638 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 639 nleaves = PetscMax(0, nleaves); 640 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 641 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 642 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 643 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 644 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 645 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 646 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 647 if (Grad) { 648 PetscFV fv; 649 650 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr); 651 ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr); 652 ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 653 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 654 ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 655 } 656 ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 657 for (i = 0; i < numids; ++i) { 658 IS faceIS; 659 const PetscInt *faces; 660 PetscInt numFaces, f; 661 662 ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 663 if (!faceIS) continue; /* No points with that id on this process */ 664 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 665 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 666 for (f = 0; f < numFaces; ++f) { 667 const PetscInt face = faces[f], *cells; 668 const PetscFVFaceGeom *fg; 669 670 if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 671 ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr); 672 if (loc >= 0) continue; 673 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 674 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 675 if (Grad) { 676 const PetscFVCellGeom *cg; 677 const PetscScalar *cx, *cgrad; 678 PetscScalar *xG; 679 PetscReal dx[3]; 680 PetscInt d; 681 682 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 683 ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 684 ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 685 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 686 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 687 for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 688 ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 689 } else { 690 const PetscScalar *xI; 691 PetscScalar *xG; 692 693 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 694 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 695 ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 696 } 697 } 698 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 699 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 700 } 701 ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 702 if (Grad) { 703 ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 704 ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 705 } 706 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 707 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 708 PetscFunctionReturn(0); 709 } 710 711 #undef __FUNCT__ 712 #define __FUNCT__ "DMPlexInsertBoundaryValues" 713 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 714 { 715 PetscInt numBd, b; 716 PetscErrorCode ierr; 717 718 PetscFunctionBegin; 719 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 720 PetscValidHeaderSpecific(locX, VEC_CLASSID, 2); 721 if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);} 722 if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);} 723 if (gradFVM) {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);} 724 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 725 for (b = 0; b < numBd; ++b) { 726 PetscBool isEssential; 727 const char *labelname; 728 DMLabel label; 729 PetscInt field; 730 PetscObject obj; 731 PetscClassId id; 732 void (*func)(); 733 PetscInt numids; 734 const PetscInt *ids; 735 void *ctx; 736 737 ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 738 if (!isEssential) continue; 739 ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 740 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 741 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 742 if (id == PETSCFE_CLASSID) { 743 ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr); 744 } else if (id == PETSCFV_CLASSID) { 745 if (!faceGeomFVM) continue; 746 ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, 747 field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr); 748 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 749 } 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "DMPlexComputeL2Diff" 755 /*@C 756 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 757 758 Input Parameters: 759 + dm - The DM 760 . funcs - The functions to evaluate for each field component 761 . ctxs - Optional array of contexts to pass to each function, or NULL. 762 - X - The coefficient vector u_h 763 764 Output Parameter: 765 . diff - The diff ||u - u_h||_2 766 767 Level: developer 768 769 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff() 770 @*/ 771 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 772 { 773 const PetscInt debug = 0; 774 PetscSection section; 775 PetscQuadrature quad; 776 Vec localX; 777 PetscScalar *funcVal, *interpolant; 778 PetscReal *coords, *v0, *J, *invJ, detJ; 779 PetscReal localDiff = 0.0; 780 const PetscReal *quadPoints, *quadWeights; 781 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset; 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 786 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 787 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 788 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 789 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 790 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 791 for (field = 0; field < numFields; ++field) { 792 PetscObject obj; 793 PetscClassId id; 794 PetscInt Nc; 795 796 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 797 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 798 if (id == PETSCFE_CLASSID) { 799 PetscFE fe = (PetscFE) obj; 800 801 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 802 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 803 } else if (id == PETSCFV_CLASSID) { 804 PetscFV fv = (PetscFV) obj; 805 806 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 807 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 808 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 809 numComponents += Nc; 810 } 811 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 812 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 813 ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 814 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 815 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 816 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 817 for (c = cStart; c < cEnd; ++c) { 818 PetscScalar *x = NULL; 819 PetscReal elemDiff = 0.0; 820 821 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 822 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 823 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 824 825 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 826 PetscObject obj; 827 PetscClassId id; 828 void * const ctx = ctxs ? ctxs[field] : NULL; 829 PetscInt Nb, Nc, q, fc; 830 831 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 832 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 833 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 834 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 835 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 836 if (debug) { 837 char title[1024]; 838 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 839 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 840 } 841 for (q = 0; q < numQuadPoints; ++q) { 842 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 843 (*funcs[field])(coords, funcVal, ctx); 844 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 845 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 846 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 847 for (fc = 0; fc < Nc; ++fc) { 848 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 849 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 850 } 851 } 852 fieldOffset += Nb*Nc; 853 } 854 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 855 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 856 localDiff += elemDiff; 857 } 858 ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 859 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 860 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 861 *diff = PetscSqrtReal(*diff); 862 PetscFunctionReturn(0); 863 } 864 865 #undef __FUNCT__ 866 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 867 /*@C 868 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 869 870 Input Parameters: 871 + dm - The DM 872 . funcs - The gradient functions to evaluate for each field component 873 . ctxs - Optional array of contexts to pass to each function, or NULL. 874 . X - The coefficient vector u_h 875 - n - The vector to project along 876 877 Output Parameter: 878 . diff - The diff ||(grad u - grad u_h) . n||_2 879 880 Level: developer 881 882 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 883 @*/ 884 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 885 { 886 const PetscInt debug = 0; 887 PetscSection section; 888 PetscQuadrature quad; 889 Vec localX; 890 PetscScalar *funcVal, *interpolantVec; 891 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 892 PetscReal localDiff = 0.0; 893 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp; 894 PetscErrorCode ierr; 895 896 PetscFunctionBegin; 897 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 898 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 899 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 900 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 901 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 902 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 903 for (field = 0; field < numFields; ++field) { 904 PetscFE fe; 905 PetscInt Nc; 906 907 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 908 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 909 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 910 numComponents += Nc; 911 } 912 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 913 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 914 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 915 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 916 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 917 for (c = cStart; c < cEnd; ++c) { 918 PetscScalar *x = NULL; 919 PetscReal elemDiff = 0.0; 920 921 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 922 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 923 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 924 925 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 926 PetscFE fe; 927 void * const ctx = ctxs ? ctxs[field] : NULL; 928 const PetscReal *quadPoints, *quadWeights; 929 PetscReal *basisDer; 930 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 931 932 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 933 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 934 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 935 ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr); 936 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 937 if (debug) { 938 char title[1024]; 939 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 940 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 941 } 942 for (q = 0; q < numQuadPoints; ++q) { 943 for (d = 0; d < dim; d++) { 944 coords[d] = v0[d]; 945 for (e = 0; e < dim; e++) { 946 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 947 } 948 } 949 (*funcs[field])(coords, n, funcVal, ctx); 950 for (fc = 0; fc < Ncomp; ++fc) { 951 PetscScalar interpolant = 0.0; 952 953 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 954 for (f = 0; f < Nb; ++f) { 955 const PetscInt fidx = f*Ncomp+fc; 956 957 for (d = 0; d < dim; ++d) { 958 realSpaceDer[d] = 0.0; 959 for (g = 0; g < dim; ++g) { 960 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 961 } 962 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 963 } 964 } 965 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 966 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);} 967 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 968 } 969 } 970 comp += Ncomp; 971 fieldOffset += Nb*Ncomp; 972 } 973 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 974 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 975 localDiff += elemDiff; 976 } 977 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 978 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 979 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 980 *diff = PetscSqrtReal(*diff); 981 PetscFunctionReturn(0); 982 } 983 984 #undef __FUNCT__ 985 #define __FUNCT__ "DMPlexComputeL2FieldDiff" 986 /*@C 987 DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components. 988 989 Input Parameters: 990 + dm - The DM 991 . funcs - The functions to evaluate for each field component 992 . ctxs - Optional array of contexts to pass to each function, or NULL. 993 - X - The coefficient vector u_h 994 995 Output Parameter: 996 . diff - The array of differences, ||u^f - u^f_h||_2 997 998 Level: developer 999 1000 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff() 1001 @*/ 1002 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 1003 { 1004 const PetscInt debug = 0; 1005 PetscSection section; 1006 PetscQuadrature quad; 1007 Vec localX; 1008 PetscScalar *funcVal, *interpolant; 1009 PetscReal *coords, *v0, *J, *invJ, detJ; 1010 PetscReal *localDiff; 1011 const PetscReal *quadPoints, *quadWeights; 1012 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset; 1013 PetscErrorCode ierr; 1014 1015 PetscFunctionBegin; 1016 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1017 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1018 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 1019 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 1020 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1021 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1022 for (field = 0; field < numFields; ++field) { 1023 PetscObject obj; 1024 PetscClassId id; 1025 PetscInt Nc; 1026 1027 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1028 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1029 if (id == PETSCFE_CLASSID) { 1030 PetscFE fe = (PetscFE) obj; 1031 1032 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1033 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1034 } else if (id == PETSCFV_CLASSID) { 1035 PetscFV fv = (PetscFV) obj; 1036 1037 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 1038 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1039 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1040 numComponents += Nc; 1041 } 1042 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 1043 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 1044 ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 1045 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1046 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1047 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1048 for (c = cStart; c < cEnd; ++c) { 1049 PetscScalar *x = NULL; 1050 1051 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 1052 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 1053 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1054 1055 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 1056 PetscObject obj; 1057 PetscClassId id; 1058 void * const ctx = ctxs ? ctxs[field] : NULL; 1059 PetscInt Nb, Nc, q, fc; 1060 1061 PetscReal elemDiff = 0.0; 1062 1063 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1064 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1065 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 1066 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 1067 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1068 if (debug) { 1069 char title[1024]; 1070 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 1071 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 1072 } 1073 for (q = 0; q < numQuadPoints; ++q) { 1074 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 1075 (*funcs[field])(coords, funcVal, ctx); 1076 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1077 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1078 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1079 for (fc = 0; fc < Nc; ++fc) { 1080 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 1081 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 1082 } 1083 } 1084 fieldOffset += Nb*Nc; 1085 localDiff[field] += elemDiff; 1086 } 1087 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1088 } 1089 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1090 ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1091 for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 1092 ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 1093 PetscFunctionReturn(0); 1094 } 1095 1096 #undef __FUNCT__ 1097 #define __FUNCT__ "DMPlexComputeIntegralFEM" 1098 /*@ 1099 DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 1100 1101 Input Parameters: 1102 + dm - The mesh 1103 . X - Local input vector 1104 - user - The user context 1105 1106 Output Parameter: 1107 . integral - Local integral for each field 1108 1109 Level: developer 1110 1111 .seealso: DMPlexComputeResidualFEM() 1112 @*/ 1113 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 1114 { 1115 DM_Plex *mesh = (DM_Plex *) dm->data; 1116 DM dmAux; 1117 Vec localX, A; 1118 PetscDS prob, probAux = NULL; 1119 PetscQuadrature q; 1120 PetscSection section, sectionAux; 1121 PetscFECellGeom *cgeom; 1122 PetscScalar *u, *a = NULL; 1123 PetscReal *lintegral, *vol; 1124 PetscInt dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c; 1125 PetscInt totDim, totDimAux; 1126 PetscErrorCode ierr; 1127 1128 PetscFunctionBegin; 1129 ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1130 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 1131 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1132 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1133 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1134 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1135 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1136 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1137 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1138 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1139 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1140 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1141 numCells = cEnd - cStart; 1142 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1143 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1144 if (dmAux) { 1145 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1146 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1147 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1148 } 1149 ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1150 ierr = PetscMalloc4(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeom,numCells,&vol);CHKERRQ(ierr); 1151 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1152 for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;} 1153 for (c = cStart; c < cEnd; ++c) { 1154 PetscScalar *x = NULL; 1155 PetscInt i; 1156 1157 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr); 1158 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol[c], NULL, NULL);CHKERRQ(ierr); 1159 if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c); 1160 ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 1161 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1162 ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 1163 if (dmAux) { 1164 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1165 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1166 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1167 } 1168 } 1169 for (f = 0; f < Nf; ++f) { 1170 PetscObject obj; 1171 PetscClassId id; 1172 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1173 1174 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1175 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1176 if (id == PETSCFE_CLASSID) { 1177 PetscFE fe = (PetscFE) obj; 1178 PetscQuadrature q; 1179 PetscInt Nq, Nb; 1180 1181 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1182 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1183 ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1184 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1185 blockSize = Nb*Nq; 1186 batchSize = numBlocks * blockSize; 1187 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1188 numChunks = numCells / (numBatches*batchSize); 1189 Ne = numChunks*numBatches*batchSize; 1190 Nr = numCells % (numBatches*batchSize); 1191 offset = numCells - Nr; 1192 ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, lintegral);CHKERRQ(ierr); 1193 ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr); 1194 } else if (id == PETSCFV_CLASSID) { 1195 PetscFV fv = (PetscFV) obj; 1196 PetscInt foff; 1197 void (*obj_func)(const PetscScalar u[], const PetscScalar u_x[], const PetscScalar u_t[], const PetscScalar a[], const PetscScalar a_x[], const PetscScalar a_t[], const PetscReal x[], PetscScalar obj[]); 1198 1199 ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr); 1200 ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1201 if (obj_func) { 1202 for (c = 0; c < numCells; ++c) { 1203 /* TODO: Need full pointwise interpolation and get centroid for x argument */ 1204 obj_func(&u[totDim*c+foff], NULL, NULL, NULL, NULL, NULL, NULL, &lintegral[f]); 1205 lintegral[f] *= vol[c]; 1206 } 1207 } 1208 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1209 } 1210 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1211 if (mesh->printFEM) { 1212 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 1213 for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);} 1214 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 1215 } 1216 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1217 ierr = MPI_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 1218 ierr = PetscFree4(lintegral,u,cgeom,vol);CHKERRQ(ierr); 1219 ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 1223 #undef __FUNCT__ 1224 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1225 /*@ 1226 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1227 1228 Input Parameters: 1229 + dmf - The fine mesh 1230 . dmc - The coarse mesh 1231 - user - The user context 1232 1233 Output Parameter: 1234 . In - The interpolation matrix 1235 1236 Note: 1237 The first member of the user context must be an FEMContext. 1238 1239 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1240 like a GPU, or vectorize on a multicore machine. 1241 1242 Level: developer 1243 1244 .seealso: DMPlexComputeJacobianFEM() 1245 @*/ 1246 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1247 { 1248 DM_Plex *mesh = (DM_Plex *) dmc->data; 1249 const char *name = "Interpolator"; 1250 PetscDS prob; 1251 PetscFE *feRef; 1252 PetscSection fsection, fglobalSection; 1253 PetscSection csection, cglobalSection; 1254 PetscScalar *elemMat; 1255 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c; 1256 PetscInt cTotDim, rTotDim = 0; 1257 PetscErrorCode ierr; 1258 1259 PetscFunctionBegin; 1260 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1261 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1262 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1263 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1264 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1265 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1266 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1267 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1268 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1269 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1270 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1271 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1272 for (f = 0; f < Nf; ++f) { 1273 PetscFE fe; 1274 PetscInt rNb, Nc; 1275 1276 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1277 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1278 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1279 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1280 rTotDim += rNb*Nc; 1281 } 1282 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1283 ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 1284 ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1285 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1286 PetscDualSpace Qref; 1287 PetscQuadrature f; 1288 const PetscReal *qpoints, *qweights; 1289 PetscReal *points; 1290 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1291 1292 /* Compose points from all dual basis functionals */ 1293 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1294 ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1295 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1296 for (i = 0; i < fpdim; ++i) { 1297 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1298 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1299 npoints += Np; 1300 } 1301 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1302 for (i = 0, k = 0; i < fpdim; ++i) { 1303 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1304 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1305 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1306 } 1307 1308 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1309 PetscFE fe; 1310 PetscReal *B; 1311 PetscInt NcJ, cpdim, j; 1312 1313 /* Evaluate basis at points */ 1314 ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr); 1315 ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 1316 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1317 /* For now, fields only interpolate themselves */ 1318 if (fieldI == fieldJ) { 1319 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); 1320 ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1321 for (i = 0, k = 0; i < fpdim; ++i) { 1322 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1323 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1324 for (p = 0; p < Np; ++p, ++k) { 1325 for (j = 0; j < cpdim; ++j) { 1326 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p]; 1327 } 1328 } 1329 } 1330 ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1331 } 1332 offsetJ += cpdim*NcJ; 1333 } 1334 offsetI += fpdim*Nc; 1335 ierr = PetscFree(points);CHKERRQ(ierr); 1336 } 1337 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 1338 /* Preallocate matrix */ 1339 { 1340 PetscHashJK ht; 1341 PetscLayout rLayout; 1342 PetscInt *dnz, *onz, *cellCIndices, *cellFIndices; 1343 PetscInt locRows, rStart, rEnd, cell, r; 1344 1345 ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 1346 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 1347 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 1348 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 1349 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 1350 ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 1351 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 1352 ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 1353 ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr); 1354 for (cell = cStart; cell < cEnd; ++cell) { 1355 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 1356 for (r = 0; r < rTotDim; ++r) { 1357 PetscHashJKKey key; 1358 PetscHashJKIter missing, iter; 1359 1360 key.j = cellFIndices[r]; 1361 if (key.j < 0) continue; 1362 for (c = 0; c < cTotDim; ++c) { 1363 key.k = cellCIndices[c]; 1364 if (key.k < 0) continue; 1365 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr); 1366 if (missing) { 1367 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr); 1368 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart]; 1369 else ++onz[key.j-rStart]; 1370 } 1371 } 1372 } 1373 } 1374 ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr); 1375 ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 1376 ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1377 ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr); 1378 } 1379 /* Fill matrix */ 1380 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1381 for (c = cStart; c < cEnd; ++c) { 1382 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1383 } 1384 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1385 ierr = PetscFree(feRef);CHKERRQ(ierr); 1386 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1387 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1388 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1389 if (mesh->printFEM) { 1390 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1391 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1392 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1393 } 1394 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1395 PetscFunctionReturn(0); 1396 } 1397 1398 #undef __FUNCT__ 1399 #define __FUNCT__ "DMPlexComputeInjectorFEM" 1400 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 1401 { 1402 PetscDS prob; 1403 PetscFE *feRef; 1404 Vec fv, cv; 1405 IS fis, cis; 1406 PetscSection fsection, fglobalSection, csection, cglobalSection; 1407 PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 1408 PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m; 1409 PetscErrorCode ierr; 1410 1411 PetscFunctionBegin; 1412 ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1413 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1414 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1415 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1416 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1417 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1418 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1419 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1420 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1421 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1422 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 1423 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1424 for (f = 0; f < Nf; ++f) { 1425 PetscFE fe; 1426 PetscInt fNb, Nc; 1427 1428 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1429 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1430 ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 1431 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1432 fTotDim += fNb*Nc; 1433 } 1434 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1435 ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 1436 for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 1437 PetscFE feC; 1438 PetscDualSpace QF, QC; 1439 PetscInt NcF, NcC, fpdim, cpdim; 1440 1441 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 1442 ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 1443 ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 1444 if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC); 1445 ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 1446 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 1447 ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 1448 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 1449 for (c = 0; c < cpdim; ++c) { 1450 PetscQuadrature cfunc; 1451 const PetscReal *cqpoints; 1452 PetscInt NpC; 1453 1454 ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 1455 ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 1456 if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 1457 for (f = 0; f < fpdim; ++f) { 1458 PetscQuadrature ffunc; 1459 const PetscReal *fqpoints; 1460 PetscReal sum = 0.0; 1461 PetscInt NpF, comp; 1462 1463 ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 1464 ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 1465 if (NpC != NpF) continue; 1466 for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 1467 if (sum > 1.0e-9) continue; 1468 for (comp = 0; comp < NcC; ++comp) { 1469 cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 1470 } 1471 break; 1472 } 1473 } 1474 offsetC += cpdim*NcC; 1475 offsetF += fpdim*NcF; 1476 } 1477 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1478 ierr = PetscFree(feRef);CHKERRQ(ierr); 1479 1480 ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 1481 ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 1482 ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr); 1483 ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 1484 ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1485 ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1486 ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 1487 for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 1488 for (c = cStart; c < cEnd; ++c) { 1489 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 1490 for (d = 0; d < cTotDim; ++d) { 1491 if (cellCIndices[d] < 0) continue; 1492 if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]); 1493 cindices[cellCIndices[d]-startC] = cellCIndices[d]; 1494 findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 1495 } 1496 } 1497 ierr = PetscFree(cmap);CHKERRQ(ierr); 1498 ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 1499 1500 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 1501 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 1502 ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 1503 ierr = ISDestroy(&cis);CHKERRQ(ierr); 1504 ierr = ISDestroy(&fis);CHKERRQ(ierr); 1505 ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 1506 ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 1507 ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1508 PetscFunctionReturn(0); 1509 } 1510