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 Nf, 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, &Nf);CHKERRQ(ierr); 377 ierr = PetscMalloc2(Nf, &sp, Nf, &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(Nf,&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 < Nf; ++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 < Nf; ++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 < Nf; 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 < Nf; 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 PetscDualSpace *sp; 511 PetscInt *Ncf; 512 PetscScalar *values, *u, *u_x, *a, *a_x; 513 PetscReal *x, *v0, *J, *invJ, detJ; 514 PetscInt Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight; 515 PetscErrorCode ierr; 516 517 PetscFunctionBegin; 518 ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr); 519 if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");} 520 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 521 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 522 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 523 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 524 ierr = PetscMalloc2(Nf, &sp, Nf, &Ncf);CHKERRQ(ierr); 525 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 526 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 527 ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 528 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 529 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 530 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 531 if (dmAux) { 532 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 533 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 534 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 535 } 536 ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 537 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 538 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 539 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 540 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 541 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 542 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 543 for (c = cStart; c < cEnd; ++c) { 544 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 545 546 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 547 ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 548 if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 549 for (f = 0, v = 0; f < Nf; ++f) { 550 PetscObject obj; 551 PetscClassId id; 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[f]);CHKERRQ(ierr); 559 ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr); 560 ierr = PetscFEGetNumComponents(fe, &Ncf[f]);CHKERRQ(ierr); 561 } else if (id == PETSCFV_CLASSID) { 562 PetscFV fv = (PetscFV) obj; 563 PetscQuadrature q; 564 565 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 566 ierr = PetscFVGetNumComponents(fv, &Ncf[f]);CHKERRQ(ierr); 567 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr); 568 ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr); 569 ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 570 ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr); 571 ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr); 572 } 573 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 574 for (d = 0; d < spDim; ++d) { 575 PetscQuadrature quad; 576 const PetscReal *points, *weights; 577 PetscInt numPoints, q; 578 579 if (funcs[f]) { 580 ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr); 581 ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 582 for (q = 0; q < numPoints; ++q) { 583 CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x); 584 ierr = EvaluateFieldJets(prob, PETSC_FALSE, q, invJ, coefficients, NULL, u, u_x, NULL);CHKERRQ(ierr); 585 ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr); 586 (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]); 587 } 588 } else { 589 for (comp = 0; comp < Ncf[f]; ++comp) values[v+comp] = 0.0; 590 } 591 v += Ncf[f]; 592 } 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 for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 601 ierr = PetscFree2(sp, Ncf);CHKERRQ(ierr); 602 PetscFunctionReturn(0); 603 } 604 605 #undef __FUNCT__ 606 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal" 607 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) 608 { 609 void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 610 void **ctxs; 611 PetscInt numFields; 612 PetscErrorCode ierr; 613 614 PetscFunctionBegin; 615 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 616 ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 617 funcs[field] = func; 618 ctxs[field] = ctx; 619 ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 620 ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 621 PetscFunctionReturn(0); 622 } 623 624 #undef __FUNCT__ 625 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal" 626 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, 627 PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX) 628 { 629 PetscDS prob; 630 PetscSF sf; 631 DM dmFace, dmCell, dmGrad; 632 const PetscScalar *facegeom, *cellgeom, *grad; 633 const PetscInt *leaves; 634 PetscScalar *x, *fx; 635 PetscInt dim, nleaves, loc, fStart, fEnd, pdim, i; 636 PetscErrorCode ierr; 637 638 PetscFunctionBegin; 639 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 640 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 641 nleaves = PetscMax(0, nleaves); 642 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 643 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 644 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 645 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 646 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 647 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 648 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 649 if (Grad) { 650 PetscFV fv; 651 652 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr); 653 ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr); 654 ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 655 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 656 ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 657 } 658 ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 659 for (i = 0; i < numids; ++i) { 660 IS faceIS; 661 const PetscInt *faces; 662 PetscInt numFaces, f; 663 664 ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 665 if (!faceIS) continue; /* No points with that id on this process */ 666 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 667 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 668 for (f = 0; f < numFaces; ++f) { 669 const PetscInt face = faces[f], *cells; 670 const PetscFVFaceGeom *fg; 671 672 if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 673 ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr); 674 if (loc >= 0) continue; 675 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 676 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 677 if (Grad) { 678 const PetscFVCellGeom *cg; 679 const PetscScalar *cx, *cgrad; 680 PetscScalar *xG; 681 PetscReal dx[3]; 682 PetscInt d; 683 684 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 685 ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 686 ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 687 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 688 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 689 for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 690 ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 691 } else { 692 const PetscScalar *xI; 693 PetscScalar *xG; 694 695 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 696 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 697 ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 698 } 699 } 700 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 701 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 702 } 703 ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 704 if (Grad) { 705 ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 706 ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 707 } 708 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 709 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 710 PetscFunctionReturn(0); 711 } 712 713 #undef __FUNCT__ 714 #define __FUNCT__ "DMPlexInsertBoundaryValues" 715 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 716 { 717 PetscInt numBd, b; 718 PetscErrorCode ierr; 719 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 722 PetscValidHeaderSpecific(locX, VEC_CLASSID, 2); 723 if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);} 724 if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);} 725 if (gradFVM) {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);} 726 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 727 for (b = 0; b < numBd; ++b) { 728 PetscBool isEssential; 729 const char *labelname; 730 DMLabel label; 731 PetscInt field; 732 PetscObject obj; 733 PetscClassId id; 734 void (*func)(); 735 PetscInt numids; 736 const PetscInt *ids; 737 void *ctx; 738 739 ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 740 if (!isEssential) continue; 741 ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 742 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 743 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 744 if (id == PETSCFE_CLASSID) { 745 ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr); 746 } else if (id == PETSCFV_CLASSID) { 747 if (!faceGeomFVM) continue; 748 ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, 749 field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr); 750 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 751 } 752 PetscFunctionReturn(0); 753 } 754 755 #undef __FUNCT__ 756 #define __FUNCT__ "DMPlexComputeL2Diff" 757 /*@C 758 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 759 760 Input Parameters: 761 + dm - The DM 762 . funcs - The functions to evaluate for each field component 763 . ctxs - Optional array of contexts to pass to each function, or NULL. 764 - X - The coefficient vector u_h 765 766 Output Parameter: 767 . diff - The diff ||u - u_h||_2 768 769 Level: developer 770 771 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff() 772 @*/ 773 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 774 { 775 const PetscInt debug = 0; 776 PetscSection section; 777 PetscQuadrature quad; 778 Vec localX; 779 PetscScalar *funcVal, *interpolant; 780 PetscReal *coords, *v0, *J, *invJ, detJ; 781 PetscReal localDiff = 0.0; 782 const PetscReal *quadPoints, *quadWeights; 783 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset; 784 PetscErrorCode ierr; 785 786 PetscFunctionBegin; 787 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 788 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 789 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 790 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 791 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 792 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 793 for (field = 0; field < numFields; ++field) { 794 PetscObject obj; 795 PetscClassId id; 796 PetscInt Nc; 797 798 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 799 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 800 if (id == PETSCFE_CLASSID) { 801 PetscFE fe = (PetscFE) obj; 802 803 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 804 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 805 } else if (id == PETSCFV_CLASSID) { 806 PetscFV fv = (PetscFV) obj; 807 808 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 809 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 810 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 811 numComponents += Nc; 812 } 813 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 814 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 815 ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 816 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 817 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 818 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 819 for (c = cStart; c < cEnd; ++c) { 820 PetscScalar *x = NULL; 821 PetscReal elemDiff = 0.0; 822 823 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 824 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 825 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 826 827 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 828 PetscObject obj; 829 PetscClassId id; 830 void * const ctx = ctxs ? ctxs[field] : NULL; 831 PetscInt Nb, Nc, q, fc; 832 833 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 834 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 835 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 836 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 837 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 838 if (debug) { 839 char title[1024]; 840 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 841 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 842 } 843 for (q = 0; q < numQuadPoints; ++q) { 844 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 845 (*funcs[field])(coords, funcVal, ctx); 846 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 847 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 848 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 849 for (fc = 0; fc < Nc; ++fc) { 850 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);} 851 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 852 } 853 } 854 fieldOffset += Nb*Nc; 855 } 856 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 857 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 858 localDiff += elemDiff; 859 } 860 ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 861 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 862 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 863 *diff = PetscSqrtReal(*diff); 864 PetscFunctionReturn(0); 865 } 866 867 #undef __FUNCT__ 868 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 869 /*@C 870 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 871 872 Input Parameters: 873 + dm - The DM 874 . funcs - The gradient functions to evaluate for each field component 875 . ctxs - Optional array of contexts to pass to each function, or NULL. 876 . X - The coefficient vector u_h 877 - n - The vector to project along 878 879 Output Parameter: 880 . diff - The diff ||(grad u - grad u_h) . n||_2 881 882 Level: developer 883 884 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 885 @*/ 886 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 887 { 888 const PetscInt debug = 0; 889 PetscSection section; 890 PetscQuadrature quad; 891 Vec localX; 892 PetscScalar *funcVal, *interpolantVec; 893 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 894 PetscReal localDiff = 0.0; 895 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp; 896 PetscErrorCode ierr; 897 898 PetscFunctionBegin; 899 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 900 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 901 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 902 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 903 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 904 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 905 for (field = 0; field < numFields; ++field) { 906 PetscFE fe; 907 PetscInt Nc; 908 909 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 910 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 911 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 912 numComponents += Nc; 913 } 914 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 915 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 916 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 917 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 918 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 919 for (c = cStart; c < cEnd; ++c) { 920 PetscScalar *x = NULL; 921 PetscReal elemDiff = 0.0; 922 923 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 924 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 925 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 926 927 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 928 PetscFE fe; 929 void * const ctx = ctxs ? ctxs[field] : NULL; 930 const PetscReal *quadPoints, *quadWeights; 931 PetscReal *basisDer; 932 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 933 934 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 935 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 936 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 937 ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr); 938 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 939 if (debug) { 940 char title[1024]; 941 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 942 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 943 } 944 for (q = 0; q < numQuadPoints; ++q) { 945 for (d = 0; d < dim; d++) { 946 coords[d] = v0[d]; 947 for (e = 0; e < dim; e++) { 948 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 949 } 950 } 951 (*funcs[field])(coords, n, funcVal, ctx); 952 for (fc = 0; fc < Ncomp; ++fc) { 953 PetscScalar interpolant = 0.0; 954 955 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 956 for (f = 0; f < Nb; ++f) { 957 const PetscInt fidx = f*Ncomp+fc; 958 959 for (d = 0; d < dim; ++d) { 960 realSpaceDer[d] = 0.0; 961 for (g = 0; g < dim; ++g) { 962 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 963 } 964 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 965 } 966 } 967 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 968 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);} 969 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 970 } 971 } 972 comp += Ncomp; 973 fieldOffset += Nb*Ncomp; 974 } 975 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 976 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 977 localDiff += elemDiff; 978 } 979 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 980 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 981 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 982 *diff = PetscSqrtReal(*diff); 983 PetscFunctionReturn(0); 984 } 985 986 #undef __FUNCT__ 987 #define __FUNCT__ "DMPlexComputeL2FieldDiff" 988 /*@C 989 DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components. 990 991 Input Parameters: 992 + dm - The DM 993 . funcs - The functions to evaluate for each field component 994 . ctxs - Optional array of contexts to pass to each function, or NULL. 995 - X - The coefficient vector u_h 996 997 Output Parameter: 998 . diff - The array of differences, ||u^f - u^f_h||_2 999 1000 Level: developer 1001 1002 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff() 1003 @*/ 1004 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 1005 { 1006 const PetscInt debug = 0; 1007 PetscSection section; 1008 PetscQuadrature quad; 1009 Vec localX; 1010 PetscScalar *funcVal, *interpolant; 1011 PetscReal *coords, *v0, *J, *invJ, detJ; 1012 PetscReal *localDiff; 1013 const PetscReal *quadPoints, *quadWeights; 1014 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset; 1015 PetscErrorCode ierr; 1016 1017 PetscFunctionBegin; 1018 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1019 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1020 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 1021 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 1022 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1023 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1024 for (field = 0; field < numFields; ++field) { 1025 PetscObject obj; 1026 PetscClassId id; 1027 PetscInt Nc; 1028 1029 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1030 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1031 if (id == PETSCFE_CLASSID) { 1032 PetscFE fe = (PetscFE) obj; 1033 1034 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1035 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1036 } else if (id == PETSCFV_CLASSID) { 1037 PetscFV fv = (PetscFV) obj; 1038 1039 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 1040 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1041 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1042 numComponents += Nc; 1043 } 1044 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 1045 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 1046 ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 1047 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1048 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1049 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1050 for (c = cStart; c < cEnd; ++c) { 1051 PetscScalar *x = NULL; 1052 1053 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 1054 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 1055 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1056 1057 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 1058 PetscObject obj; 1059 PetscClassId id; 1060 void * const ctx = ctxs ? ctxs[field] : NULL; 1061 PetscInt Nb, Nc, q, fc; 1062 1063 PetscReal elemDiff = 0.0; 1064 1065 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1066 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1067 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 1068 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 1069 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1070 if (debug) { 1071 char title[1024]; 1072 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 1073 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 1074 } 1075 for (q = 0; q < numQuadPoints; ++q) { 1076 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 1077 (*funcs[field])(coords, funcVal, ctx); 1078 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1079 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1080 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1081 for (fc = 0; fc < Nc; ++fc) { 1082 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);} 1083 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 1084 } 1085 } 1086 fieldOffset += Nb*Nc; 1087 localDiff[field] += elemDiff; 1088 } 1089 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1090 } 1091 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1092 ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1093 for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 1094 ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 1095 PetscFunctionReturn(0); 1096 } 1097 1098 #undef __FUNCT__ 1099 #define __FUNCT__ "DMPlexComputeIntegralFEM" 1100 /*@ 1101 DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 1102 1103 Input Parameters: 1104 + dm - The mesh 1105 . X - Local input vector 1106 - user - The user context 1107 1108 Output Parameter: 1109 . integral - Local integral for each field 1110 1111 Level: developer 1112 1113 .seealso: DMPlexComputeResidualFEM() 1114 @*/ 1115 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 1116 { 1117 DM_Plex *mesh = (DM_Plex *) dm->data; 1118 DM dmAux; 1119 Vec localX, A; 1120 PetscDS prob, probAux = NULL; 1121 PetscQuadrature q; 1122 PetscSection section, sectionAux; 1123 PetscFECellGeom *cgeom; 1124 PetscScalar *u, *a = NULL; 1125 PetscReal *lintegral, *vol; 1126 PetscInt dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c; 1127 PetscInt totDim, totDimAux; 1128 PetscErrorCode ierr; 1129 1130 PetscFunctionBegin; 1131 ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1132 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 1133 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1134 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 1135 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1136 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1137 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1138 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1139 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1140 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1141 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1142 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1143 numCells = cEnd - cStart; 1144 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1145 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1146 if (dmAux) { 1147 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1148 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1149 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1150 } 1151 ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1152 ierr = PetscMalloc4(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeom,numCells,&vol);CHKERRQ(ierr); 1153 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1154 for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;} 1155 for (c = cStart; c < cEnd; ++c) { 1156 PetscScalar *x = NULL; 1157 PetscInt i; 1158 1159 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr); 1160 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol[c], NULL, NULL);CHKERRQ(ierr); 1161 if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c); 1162 ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 1163 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1164 ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 1165 if (dmAux) { 1166 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1167 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1168 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1169 } 1170 } 1171 for (f = 0; f < Nf; ++f) { 1172 PetscObject obj; 1173 PetscClassId id; 1174 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1175 1176 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1177 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1178 if (id == PETSCFE_CLASSID) { 1179 PetscFE fe = (PetscFE) obj; 1180 PetscQuadrature q; 1181 PetscInt Nq, Nb; 1182 1183 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1184 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1185 ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1186 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1187 blockSize = Nb*Nq; 1188 batchSize = numBlocks * blockSize; 1189 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1190 numChunks = numCells / (numBatches*batchSize); 1191 Ne = numChunks*numBatches*batchSize; 1192 Nr = numCells % (numBatches*batchSize); 1193 offset = numCells - Nr; 1194 ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, lintegral);CHKERRQ(ierr); 1195 ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr); 1196 } else if (id == PETSCFV_CLASSID) { 1197 PetscFV fv = (PetscFV) obj; 1198 PetscInt foff; 1199 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[]); 1200 1201 ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr); 1202 ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1203 if (obj_func) { 1204 for (c = 0; c < numCells; ++c) { 1205 /* TODO: Need full pointwise interpolation and get centroid for x argument */ 1206 obj_func(&u[totDim*c+foff], NULL, NULL, NULL, NULL, NULL, NULL, &lintegral[f]); 1207 lintegral[f] *= vol[c]; 1208 } 1209 } 1210 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1211 } 1212 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1213 if (mesh->printFEM) { 1214 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 1215 for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);} 1216 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 1217 } 1218 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1219 ierr = MPI_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 1220 ierr = PetscFree4(lintegral,u,cgeom,vol);CHKERRQ(ierr); 1221 ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1222 PetscFunctionReturn(0); 1223 } 1224 1225 #undef __FUNCT__ 1226 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1227 /*@ 1228 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1229 1230 Input Parameters: 1231 + dmf - The fine mesh 1232 . dmc - The coarse mesh 1233 - user - The user context 1234 1235 Output Parameter: 1236 . In - The interpolation matrix 1237 1238 Note: 1239 The first member of the user context must be an FEMContext. 1240 1241 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1242 like a GPU, or vectorize on a multicore machine. 1243 1244 Level: developer 1245 1246 .seealso: DMPlexComputeJacobianFEM() 1247 @*/ 1248 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1249 { 1250 DM_Plex *mesh = (DM_Plex *) dmc->data; 1251 const char *name = "Interpolator"; 1252 PetscDS prob; 1253 PetscFE *feRef; 1254 PetscSection fsection, fglobalSection; 1255 PetscSection csection, cglobalSection; 1256 PetscScalar *elemMat; 1257 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c; 1258 PetscInt cTotDim, rTotDim = 0; 1259 PetscErrorCode ierr; 1260 1261 PetscFunctionBegin; 1262 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1263 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1264 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1265 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1266 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1267 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1268 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1269 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1270 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1271 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1272 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1273 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1274 for (f = 0; f < Nf; ++f) { 1275 PetscFE fe; 1276 PetscInt rNb, Nc; 1277 1278 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1279 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1280 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1281 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1282 rTotDim += rNb*Nc; 1283 } 1284 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1285 ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 1286 ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1287 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1288 PetscDualSpace Qref; 1289 PetscQuadrature f; 1290 const PetscReal *qpoints, *qweights; 1291 PetscReal *points; 1292 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1293 1294 /* Compose points from all dual basis functionals */ 1295 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1296 ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1297 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1298 for (i = 0; i < fpdim; ++i) { 1299 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1300 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1301 npoints += Np; 1302 } 1303 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1304 for (i = 0, k = 0; i < fpdim; ++i) { 1305 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1306 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1307 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1308 } 1309 1310 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1311 PetscFE fe; 1312 PetscReal *B; 1313 PetscInt NcJ, cpdim, j; 1314 1315 /* Evaluate basis at points */ 1316 ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr); 1317 ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 1318 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1319 /* For now, fields only interpolate themselves */ 1320 if (fieldI == fieldJ) { 1321 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); 1322 ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1323 for (i = 0, k = 0; i < fpdim; ++i) { 1324 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1325 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1326 for (p = 0; p < Np; ++p, ++k) { 1327 for (j = 0; j < cpdim; ++j) { 1328 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]; 1329 } 1330 } 1331 } 1332 ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1333 } 1334 offsetJ += cpdim*NcJ; 1335 } 1336 offsetI += fpdim*Nc; 1337 ierr = PetscFree(points);CHKERRQ(ierr); 1338 } 1339 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 1340 /* Preallocate matrix */ 1341 { 1342 PetscHashJK ht; 1343 PetscLayout rLayout; 1344 PetscInt *dnz, *onz, *cellCIndices, *cellFIndices; 1345 PetscInt locRows, rStart, rEnd, cell, r; 1346 1347 ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 1348 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 1349 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 1350 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 1351 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 1352 ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 1353 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 1354 ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 1355 ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr); 1356 for (cell = cStart; cell < cEnd; ++cell) { 1357 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 1358 for (r = 0; r < rTotDim; ++r) { 1359 PetscHashJKKey key; 1360 PetscHashJKIter missing, iter; 1361 1362 key.j = cellFIndices[r]; 1363 if (key.j < 0) continue; 1364 for (c = 0; c < cTotDim; ++c) { 1365 key.k = cellCIndices[c]; 1366 if (key.k < 0) continue; 1367 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr); 1368 if (missing) { 1369 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr); 1370 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart]; 1371 else ++onz[key.j-rStart]; 1372 } 1373 } 1374 } 1375 } 1376 ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr); 1377 ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 1378 ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1379 ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr); 1380 } 1381 /* Fill matrix */ 1382 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1383 for (c = cStart; c < cEnd; ++c) { 1384 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1385 } 1386 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1387 ierr = PetscFree(feRef);CHKERRQ(ierr); 1388 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1389 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1390 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1391 if (mesh->printFEM) { 1392 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1393 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1394 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1395 } 1396 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1397 PetscFunctionReturn(0); 1398 } 1399 1400 #undef __FUNCT__ 1401 #define __FUNCT__ "DMPlexComputeInjectorFEM" 1402 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 1403 { 1404 PetscDS prob; 1405 PetscFE *feRef; 1406 Vec fv, cv; 1407 IS fis, cis; 1408 PetscSection fsection, fglobalSection, csection, cglobalSection; 1409 PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 1410 PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m; 1411 PetscErrorCode ierr; 1412 1413 PetscFunctionBegin; 1414 ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1415 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1416 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1417 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1418 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1419 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1420 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1421 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1422 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1423 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1424 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 1425 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1426 for (f = 0; f < Nf; ++f) { 1427 PetscFE fe; 1428 PetscInt fNb, Nc; 1429 1430 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1431 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1432 ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 1433 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1434 fTotDim += fNb*Nc; 1435 } 1436 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1437 ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 1438 for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 1439 PetscFE feC; 1440 PetscDualSpace QF, QC; 1441 PetscInt NcF, NcC, fpdim, cpdim; 1442 1443 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 1444 ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 1445 ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 1446 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); 1447 ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 1448 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 1449 ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 1450 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 1451 for (c = 0; c < cpdim; ++c) { 1452 PetscQuadrature cfunc; 1453 const PetscReal *cqpoints; 1454 PetscInt NpC; 1455 1456 ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 1457 ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 1458 if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 1459 for (f = 0; f < fpdim; ++f) { 1460 PetscQuadrature ffunc; 1461 const PetscReal *fqpoints; 1462 PetscReal sum = 0.0; 1463 PetscInt NpF, comp; 1464 1465 ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 1466 ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 1467 if (NpC != NpF) continue; 1468 for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 1469 if (sum > 1.0e-9) continue; 1470 for (comp = 0; comp < NcC; ++comp) { 1471 cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 1472 } 1473 break; 1474 } 1475 } 1476 offsetC += cpdim*NcC; 1477 offsetF += fpdim*NcF; 1478 } 1479 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1480 ierr = PetscFree(feRef);CHKERRQ(ierr); 1481 1482 ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 1483 ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 1484 ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr); 1485 ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 1486 ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1487 ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1488 ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 1489 for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 1490 for (c = cStart; c < cEnd; ++c) { 1491 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 1492 for (d = 0; d < cTotDim; ++d) { 1493 if (cellCIndices[d] < 0) continue; 1494 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]]); 1495 cindices[cellCIndices[d]-startC] = cellCIndices[d]; 1496 findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 1497 } 1498 } 1499 ierr = PetscFree(cmap);CHKERRQ(ierr); 1500 ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 1501 1502 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 1503 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 1504 ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 1505 ierr = ISDestroy(&cis);CHKERRQ(ierr); 1506 ierr = ISDestroy(&fis);CHKERRQ(ierr); 1507 ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 1508 ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 1509 ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1510 PetscFunctionReturn(0); 1511 } 1512