1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petscsf.h> 3 4 #include <petsc/private/hashsetij.h> 5 #include <petsc/private/petscfeimpl.h> 6 #include <petsc/private/petscfvimpl.h> 7 8 static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom (void *ctx) 9 { 10 PetscFEGeom *geom = (PetscFEGeom *) ctx; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr); 15 PetscFunctionReturn(0); 16 } 17 18 static PetscErrorCode DMPlexGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 19 { 20 char composeStr[33] = {0}; 21 PetscObjectId id; 22 PetscContainer container; 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 ierr = PetscObjectGetId((PetscObject)quad,&id);CHKERRQ(ierr); 27 ierr = PetscSNPrintf(composeStr, 32, "DMPlexGetFEGeom_%x\n", id);CHKERRQ(ierr); 28 ierr = PetscObjectQuery((PetscObject) pointIS, composeStr, (PetscObject *) &container);CHKERRQ(ierr); 29 if (container) { 30 ierr = PetscContainerGetPointer(container, (void **) geom);CHKERRQ(ierr); 31 } else { 32 ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, faceData, geom);CHKERRQ(ierr); 33 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 34 ierr = PetscContainerSetPointer(container, (void *) *geom);CHKERRQ(ierr); 35 ierr = PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom);CHKERRQ(ierr); 36 ierr = PetscObjectCompose((PetscObject) pointIS, composeStr, (PetscObject) container);CHKERRQ(ierr); 37 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 38 } 39 PetscFunctionReturn(0); 40 } 41 42 static PetscErrorCode DMPlexRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 43 { 44 PetscFunctionBegin; 45 *geom = NULL; 46 PetscFunctionReturn(0); 47 } 48 49 /*@ 50 DMPlexGetScale - Get the scale for the specified fundamental unit 51 52 Not collective 53 54 Input Arguments: 55 + dm - the DM 56 - unit - The SI unit 57 58 Output Argument: 59 . scale - The value used to scale all quantities with this unit 60 61 Level: advanced 62 63 .seealso: DMPlexSetScale(), PetscUnit 64 @*/ 65 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 66 { 67 DM_Plex *mesh = (DM_Plex*) dm->data; 68 69 PetscFunctionBegin; 70 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 71 PetscValidPointer(scale, 3); 72 *scale = mesh->scale[unit]; 73 PetscFunctionReturn(0); 74 } 75 76 /*@ 77 DMPlexSetScale - Set the scale for the specified fundamental unit 78 79 Not collective 80 81 Input Arguments: 82 + dm - the DM 83 . unit - The SI unit 84 - scale - The value used to scale all quantities with this unit 85 86 Level: advanced 87 88 .seealso: DMPlexGetScale(), PetscUnit 89 @*/ 90 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 91 { 92 DM_Plex *mesh = (DM_Plex*) dm->data; 93 94 PetscFunctionBegin; 95 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 96 mesh->scale[unit] = scale; 97 PetscFunctionReturn(0); 98 } 99 100 static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx) 101 { 102 const PetscInt eps[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}}; 103 PetscInt *ctxInt = (PetscInt *) ctx; 104 PetscInt dim2 = ctxInt[0]; 105 PetscInt d = ctxInt[1]; 106 PetscInt i, j, k = dim > 2 ? d - dim : d; 107 108 PetscFunctionBegin; 109 if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2); 110 for (i = 0; i < dim; i++) mode[i] = 0.; 111 if (d < dim) { 112 mode[d] = 1.; /* Translation along axis d */ 113 } else { 114 for (i = 0; i < dim; i++) { 115 for (j = 0; j < dim; j++) { 116 mode[j] += eps[i][j][k]*X[i]; /* Rotation about axis d */ 117 } 118 } 119 } 120 PetscFunctionReturn(0); 121 } 122 123 /*@C 124 DMPlexCreateRigidBody - For the default global section, create rigid body modes by function space interpolation 125 126 Collective on DM 127 128 Input Arguments: 129 . dm - the DM 130 131 Output Argument: 132 . sp - the null space 133 134 Note: This is necessary to provide a suitable coarse space for algebraic multigrid 135 136 Level: advanced 137 138 .seealso: MatNullSpaceCreate(), PCGAMG 139 @*/ 140 PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp) 141 { 142 MPI_Comm comm; 143 Vec mode[6]; 144 PetscSection section, globalSection; 145 PetscInt dim, dimEmbed, n, m, d, i, j; 146 PetscErrorCode ierr; 147 148 PetscFunctionBegin; 149 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 150 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 151 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 152 if (dim == 1) { 153 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 154 PetscFunctionReturn(0); 155 } 156 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 157 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 158 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 159 m = (dim*(dim+1))/2; 160 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 161 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 162 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 163 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 164 for (d = 0; d < m; d++) { 165 PetscInt ctx[2]; 166 PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private; 167 void *voidctx = (void *) (&ctx[0]); 168 169 ctx[0] = dimEmbed; 170 ctx[1] = d; 171 ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 172 } 173 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 174 /* Orthonormalize system */ 175 for (i = dim; i < m; ++i) { 176 PetscScalar dots[6]; 177 178 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 179 for (j = 0; j < i; ++j) dots[j] *= -1.0; 180 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 181 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 182 } 183 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 184 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 185 PetscFunctionReturn(0); 186 } 187 188 /*@C 189 DMPlexCreateRigidBodies - For the default global section, create rigid body modes by function space interpolation 190 191 Collective on DM 192 193 Input Arguments: 194 + dm - the DM 195 . nb - The number of bodies 196 . label - The DMLabel marking each domain 197 . nids - The number of ids per body 198 - ids - An array of the label ids in sequence for each domain 199 200 Output Argument: 201 . sp - the null space 202 203 Note: This is necessary to provide a suitable coarse space for algebraic multigrid 204 205 Level: advanced 206 207 .seealso: MatNullSpaceCreate() 208 @*/ 209 PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp) 210 { 211 MPI_Comm comm; 212 PetscSection section, globalSection; 213 Vec *mode; 214 PetscScalar *dots; 215 PetscInt dim, dimEmbed, n, m, b, d, i, j, off; 216 PetscErrorCode ierr; 217 218 PetscFunctionBegin; 219 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 220 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 221 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 222 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 223 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 224 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 225 m = nb * (dim*(dim+1))/2; 226 ierr = PetscMalloc2(m, &mode, m, &dots);CHKERRQ(ierr); 227 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 228 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 229 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 230 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 231 for (b = 0, off = 0; b < nb; ++b) { 232 for (d = 0; d < m/nb; ++d) { 233 PetscInt ctx[2]; 234 PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private; 235 void *voidctx = (void *) (&ctx[0]); 236 237 ctx[0] = dimEmbed; 238 ctx[1] = d; 239 ierr = DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 240 off += nids[b]; 241 } 242 } 243 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 244 /* Orthonormalize system */ 245 for (i = 0; i < m; ++i) { 246 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 247 for (j = 0; j < i; ++j) dots[j] *= -1.0; 248 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 249 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 250 } 251 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 252 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 253 ierr = PetscFree2(mode, dots);CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 /*@ 258 DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs 259 are computed by associating the basis function with one of the mesh points in its transitively-closed support, and 260 evaluating the dual space basis of that point. A basis function is associated with the point in its 261 transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum 262 projection height, which is set with this function. By default, the maximum projection height is zero, which means 263 that only mesh cells are used to project basis functions. A height of one, for example, evaluates a cell-interior 264 basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face. 265 266 Input Parameters: 267 + dm - the DMPlex object 268 - height - the maximum projection height >= 0 269 270 Level: advanced 271 272 .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal() 273 @*/ 274 PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height) 275 { 276 DM_Plex *plex = (DM_Plex *) dm->data; 277 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 280 plex->maxProjectionHeight = height; 281 PetscFunctionReturn(0); 282 } 283 284 /*@ 285 DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in 286 DMPlexProjectXXXLocal() functions. 287 288 Input Parameters: 289 . dm - the DMPlex object 290 291 Output Parameters: 292 . height - the maximum projection height 293 294 Level: intermediate 295 296 .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal() 297 @*/ 298 PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height) 299 { 300 DM_Plex *plex = (DM_Plex *) dm->data; 301 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 304 *height = plex->maxProjectionHeight; 305 PetscFunctionReturn(0); 306 } 307 308 /*@C 309 DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector 310 311 Input Parameters: 312 + dm - The DM, with a PetscDS that matches the problem being constrained 313 . time - The time 314 . field - The field to constrain 315 . Nc - The number of constrained field components, or 0 for all components 316 . comps - An array of constrained component numbers, or NULL for all components 317 . label - The DMLabel defining constrained points 318 . numids - The number of DMLabel ids for constrained points 319 . ids - An array of ids for constrained points 320 . func - A pointwise function giving boundary values 321 - ctx - An optional user context for bcFunc 322 323 Output Parameter: 324 . locX - A local vector to receives the boundary values 325 326 Level: developer 327 328 .seealso: DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary() 329 @*/ 330 PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM dm, PetscReal time, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX) 331 { 332 PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 333 void **ctxs; 334 PetscInt numFields; 335 PetscErrorCode ierr; 336 337 PetscFunctionBegin; 338 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 339 ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 340 funcs[field] = func; 341 ctxs[field] = ctx; 342 ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 343 ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 344 PetscFunctionReturn(0); 345 } 346 347 /*@C 348 DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector 349 350 Input Parameters: 351 + dm - The DM, with a PetscDS that matches the problem being constrained 352 . time - The time 353 . locU - A local vector with the input solution values 354 . field - The field to constrain 355 . Nc - The number of constrained field components, or 0 for all components 356 . comps - An array of constrained component numbers, or NULL for all components 357 . label - The DMLabel defining constrained points 358 . numids - The number of DMLabel ids for constrained points 359 . ids - An array of ids for constrained points 360 . func - A pointwise function giving boundary values 361 - ctx - An optional user context for bcFunc 362 363 Output Parameter: 364 . locX - A local vector to receives the boundary values 365 366 Level: developer 367 368 .seealso: DMPlexInsertBoundaryValuesEssential(), DMAddBoundary() 369 @*/ 370 PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], 371 void (*func)(PetscInt, PetscInt, PetscInt, 372 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 373 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 374 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], 375 PetscScalar[]), 376 void *ctx, Vec locX) 377 { 378 void (**funcs)(PetscInt, PetscInt, PetscInt, 379 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 380 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 381 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 382 void **ctxs; 383 PetscInt numFields; 384 PetscErrorCode ierr; 385 386 PetscFunctionBegin; 387 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 388 ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 389 funcs[field] = func; 390 ctxs[field] = ctx; 391 ierr = DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 392 ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 /*@C 397 DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector 398 399 Input Parameters: 400 + dm - The DM, with a PetscDS that matches the problem being constrained 401 . time - The time 402 . faceGeometry - A vector with the FVM face geometry information 403 . cellGeometry - A vector with the FVM cell geometry information 404 . Grad - A vector with the FVM cell gradient information 405 . field - The field to constrain 406 . Nc - The number of constrained field components, or 0 for all components 407 . comps - An array of constrained component numbers, or NULL for all components 408 . label - The DMLabel defining constrained points 409 . numids - The number of DMLabel ids for constrained points 410 . ids - An array of ids for constrained points 411 . func - A pointwise function giving boundary values 412 - ctx - An optional user context for bcFunc 413 414 Output Parameter: 415 . locX - A local vector to receives the boundary values 416 417 Note: This implementation currently ignores the numcomps/comps argument from DMAddBoundary() 418 419 Level: developer 420 421 .seealso: DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary() 422 @*/ 423 PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], 424 PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX) 425 { 426 PetscDS prob; 427 PetscSF sf; 428 DM dmFace, dmCell, dmGrad; 429 const PetscScalar *facegeom, *cellgeom = NULL, *grad; 430 const PetscInt *leaves; 431 PetscScalar *x, *fx; 432 PetscInt dim, nleaves, loc, fStart, fEnd, pdim, i; 433 PetscErrorCode ierr, ierru = 0; 434 435 PetscFunctionBegin; 436 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 437 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 438 nleaves = PetscMax(0, nleaves); 439 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 440 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 441 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 442 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 443 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 444 if (cellGeometry) { 445 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 446 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 447 } 448 if (Grad) { 449 PetscFV fv; 450 451 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr); 452 ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr); 453 ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 454 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 455 ierr = DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx);CHKERRQ(ierr); 456 } 457 ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 458 for (i = 0; i < numids; ++i) { 459 IS faceIS; 460 const PetscInt *faces; 461 PetscInt numFaces, f; 462 463 ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 464 if (!faceIS) continue; /* No points with that id on this process */ 465 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 466 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 467 for (f = 0; f < numFaces; ++f) { 468 const PetscInt face = faces[f], *cells; 469 PetscFVFaceGeom *fg; 470 471 if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 472 ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr); 473 if (loc >= 0) continue; 474 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 475 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 476 if (Grad) { 477 PetscFVCellGeom *cg; 478 PetscScalar *cx, *cgrad; 479 PetscScalar *xG; 480 PetscReal dx[3]; 481 PetscInt d; 482 483 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 484 ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 485 ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 486 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 487 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 488 for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 489 ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx); 490 if (ierru) { 491 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 492 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 493 goto cleanup; 494 } 495 } else { 496 PetscScalar *xI; 497 PetscScalar *xG; 498 499 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 500 ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 501 ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx); 502 if (ierru) { 503 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 504 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 505 goto cleanup; 506 } 507 } 508 } 509 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 510 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 511 } 512 cleanup: 513 ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 514 if (Grad) { 515 ierr = DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx);CHKERRQ(ierr); 516 ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 517 } 518 if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);} 519 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 520 CHKERRQ(ierru); 521 PetscFunctionReturn(0); 522 } 523 524 PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 525 { 526 PetscInt numBd, b; 527 PetscErrorCode ierr; 528 529 PetscFunctionBegin; 530 ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr); 531 for (b = 0; b < numBd; ++b) { 532 DMBoundaryConditionType type; 533 const char *labelname; 534 DMLabel label; 535 PetscInt field, Nc; 536 const PetscInt *comps; 537 PetscObject obj; 538 PetscClassId id; 539 void (*func)(void); 540 PetscInt numids; 541 const PetscInt *ids; 542 void *ctx; 543 544 ierr = DMGetBoundary(dm, b, &type, NULL, &labelname, &field, &Nc, &comps, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 545 if (insertEssential != (type & DM_BC_ESSENTIAL)) continue; 546 ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr); 547 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 548 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 549 if (id == PETSCFE_CLASSID) { 550 switch (type) { 551 /* for FEM, there is no insertion to be done for non-essential boundary conditions */ 552 case DM_BC_ESSENTIAL: 553 ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr); 554 ierr = DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr); 555 ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr); 556 break; 557 case DM_BC_ESSENTIAL_FIELD: 558 ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr); 559 ierr = DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, 560 (void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 561 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 562 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) func, ctx, locX);CHKERRQ(ierr); 563 ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr); 564 break; 565 default: break; 566 } 567 } else if (id == PETSCFV_CLASSID) { 568 if (!faceGeomFVM) continue; 569 ierr = DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids, 570 (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr); 571 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 572 } 573 PetscFunctionReturn(0); 574 } 575 576 /*@ 577 DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector 578 579 Input Parameters: 580 + dm - The DM 581 . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions 582 . time - The time 583 . faceGeomFVM - Face geometry data for FV discretizations 584 . cellGeomFVM - Cell geometry data for FV discretizations 585 - gradFVM - Gradient reconstruction data for FV discretizations 586 587 Output Parameters: 588 . locX - Solution updated with boundary values 589 590 Level: developer 591 592 .seealso: DMProjectFunctionLabelLocal() 593 @*/ 594 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 595 { 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 600 PetscValidHeaderSpecific(locX, VEC_CLASSID, 2); 601 if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);} 602 if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);} 603 if (gradFVM) {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);} 604 ierr = PetscTryMethod(dm,"DMPlexInsertBoundaryValues_C",(DM,PetscBool,Vec,PetscReal,Vec,Vec,Vec),(dm,insertEssential,locX,time,faceGeomFVM,cellGeomFVM,gradFVM));CHKERRQ(ierr); 605 PetscFunctionReturn(0); 606 } 607 608 PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 609 { 610 Vec localX; 611 PetscErrorCode ierr; 612 613 PetscFunctionBegin; 614 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 615 ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL);CHKERRQ(ierr); 616 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 617 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 618 ierr = DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff);CHKERRQ(ierr); 619 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 620 PetscFunctionReturn(0); 621 } 622 623 /*@C 624 DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 625 626 Input Parameters: 627 + dm - The DM 628 . time - The time 629 . funcs - The functions to evaluate for each field component 630 . ctxs - Optional array of contexts to pass to each function, or NULL. 631 - localX - The coefficient vector u_h, a local vector 632 633 Output Parameter: 634 . diff - The diff ||u - u_h||_2 635 636 Level: developer 637 638 .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff() 639 @*/ 640 PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff) 641 { 642 const PetscInt debug = ((DM_Plex*)dm->data)->printL2; 643 PetscSection section; 644 PetscQuadrature quad; 645 PetscScalar *funcVal, *interpolant; 646 PetscReal *coords, *detJ, *J; 647 PetscReal localDiff = 0.0; 648 const PetscReal *quadWeights; 649 PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, cEndInterior, c, field, fieldOffset; 650 PetscErrorCode ierr; 651 652 PetscFunctionBegin; 653 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 654 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 655 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 656 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 657 for (field = 0; field < numFields; ++field) { 658 PetscObject obj; 659 PetscClassId id; 660 PetscInt Nc; 661 662 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 663 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 664 if (id == PETSCFE_CLASSID) { 665 PetscFE fe = (PetscFE) obj; 666 667 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 668 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 669 } else if (id == PETSCFV_CLASSID) { 670 PetscFV fv = (PetscFV) obj; 671 672 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 673 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 674 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 675 numComponents += Nc; 676 } 677 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);CHKERRQ(ierr); 678 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 679 ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 680 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 681 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 682 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 683 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 684 for (c = cStart; c < cEnd; ++c) { 685 PetscScalar *x = NULL; 686 PetscReal elemDiff = 0.0; 687 PetscInt qc = 0; 688 689 ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 690 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 691 692 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 693 PetscObject obj; 694 PetscClassId id; 695 void * const ctx = ctxs ? ctxs[field] : NULL; 696 PetscInt Nb, Nc, q, fc; 697 698 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 699 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 700 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 701 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 702 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 703 if (debug) { 704 char title[1024]; 705 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 706 ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr); 707 } 708 for (q = 0; q < Nq; ++q) { 709 if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, point %D", detJ[q], c, q); 710 ierr = (*funcs[field])(coordDim, time, &coords[coordDim * q], Nc, funcVal, ctx); 711 if (ierr) { 712 PetscErrorCode ierr2; 713 ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 714 ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 715 ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2); 716 CHKERRQ(ierr); 717 } 718 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 719 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 720 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 721 for (fc = 0; fc < Nc; ++fc) { 722 const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 723 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 724 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 725 } 726 } 727 fieldOffset += Nb; 728 qc += Nc; 729 } 730 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 731 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 732 localDiff += elemDiff; 733 } 734 ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 735 ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 736 *diff = PetscSqrtReal(*diff); 737 PetscFunctionReturn(0); 738 } 739 740 PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 741 { 742 const PetscInt debug = ((DM_Plex*)dm->data)->printL2; 743 PetscSection section; 744 PetscQuadrature quad; 745 Vec localX; 746 PetscScalar *funcVal, *interpolant; 747 const PetscReal *quadPoints, *quadWeights; 748 PetscReal *coords, *realSpaceDer, *J, *invJ, *detJ; 749 PetscReal localDiff = 0.0; 750 PetscInt dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset; 751 PetscErrorCode ierr; 752 753 PetscFunctionBegin; 754 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 755 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 756 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 757 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 758 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 759 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 760 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 761 for (field = 0; field < numFields; ++field) { 762 PetscFE fe; 763 PetscInt Nc; 764 765 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 766 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 767 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 768 numComponents += Nc; 769 } 770 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 771 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 772 /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 773 ierr = PetscMalloc7(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*Nq,&realSpaceDer,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ,numComponents,&interpolant,Nq,&detJ);CHKERRQ(ierr); 774 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 775 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 776 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 777 for (c = cStart; c < cEnd; ++c) { 778 PetscScalar *x = NULL; 779 PetscReal elemDiff = 0.0; 780 PetscInt qc = 0; 781 782 ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, invJ, detJ);CHKERRQ(ierr); 783 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 784 785 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 786 PetscFE fe; 787 void * const ctx = ctxs ? ctxs[field] : NULL; 788 PetscReal *basisDer; 789 PetscInt Nb, Nc, q, fc; 790 791 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 792 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 793 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 794 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 795 if (debug) { 796 char title[1024]; 797 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 798 ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr); 799 } 800 for (q = 0; q < Nq; ++q) { 801 if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q); 802 ierr = (*funcs[field])(coordDim, time, &coords[q*coordDim], n, numFields, funcVal, ctx); 803 if (ierr) { 804 PetscErrorCode ierr2; 805 ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 806 ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 807 ierr2 = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr2); 808 CHKERRQ(ierr); 809 } 810 ierr = PetscFEInterpolateGradient_Static(fe, &x[fieldOffset], coordDim, invJ, n, q, interpolant);CHKERRQ(ierr); 811 for (fc = 0; fc < Nc; ++fc) { 812 const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 813 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 814 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 815 } 816 } 817 fieldOffset += Nb; 818 qc += Nc; 819 } 820 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 821 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 822 localDiff += elemDiff; 823 } 824 ierr = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr); 825 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 826 ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 827 *diff = PetscSqrtReal(*diff); 828 PetscFunctionReturn(0); 829 } 830 831 PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 832 { 833 const PetscInt debug = ((DM_Plex*)dm->data)->printL2; 834 PetscSection section; 835 PetscQuadrature quad; 836 Vec localX; 837 PetscScalar *funcVal, *interpolant; 838 PetscReal *coords, *detJ, *J; 839 PetscReal *localDiff; 840 const PetscReal *quadPoints, *quadWeights; 841 PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset; 842 PetscErrorCode ierr; 843 844 PetscFunctionBegin; 845 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 846 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 847 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 848 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 849 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 850 ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 851 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 852 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 853 for (field = 0; field < numFields; ++field) { 854 PetscObject obj; 855 PetscClassId id; 856 PetscInt Nc; 857 858 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 859 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 860 if (id == PETSCFE_CLASSID) { 861 PetscFE fe = (PetscFE) obj; 862 863 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 864 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 865 } else if (id == PETSCFV_CLASSID) { 866 PetscFV fv = (PetscFV) obj; 867 868 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 869 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 870 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 871 numComponents += Nc; 872 } 873 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 874 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 875 ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 876 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 877 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 878 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 879 for (c = cStart; c < cEnd; ++c) { 880 PetscScalar *x = NULL; 881 PetscInt qc = 0; 882 883 ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 884 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 885 886 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 887 PetscObject obj; 888 PetscClassId id; 889 void * const ctx = ctxs ? ctxs[field] : NULL; 890 PetscInt Nb, Nc, q, fc; 891 892 PetscReal elemDiff = 0.0; 893 894 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 895 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 896 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 897 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 898 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 899 if (debug) { 900 char title[1024]; 901 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 902 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 903 } 904 for (q = 0; q < Nq; ++q) { 905 if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature point %D", detJ, c, q); 906 ierr = (*funcs[field])(coordDim, time, &coords[coordDim*q], numFields, funcVal, ctx); 907 if (ierr) { 908 PetscErrorCode ierr2; 909 ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 910 ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 911 ierr2 = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2); 912 CHKERRQ(ierr); 913 } 914 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 915 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 916 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 917 for (fc = 0; fc < Nc; ++fc) { 918 const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 919 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d point %g %g %g diff %g\n", c, field, coordDim > 0 ? coords[0] : 0., coordDim > 1 ? coords[1] : 0., coordDim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 920 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 921 } 922 } 923 fieldOffset += Nb; 924 qc += Nc; 925 localDiff[field] += elemDiff; 926 } 927 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 928 } 929 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 930 ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 931 for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 932 ierr = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 933 PetscFunctionReturn(0); 934 } 935 936 /*@C 937 DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec. 938 939 Input Parameters: 940 + dm - The DM 941 . time - The time 942 . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation 943 . ctxs - Optional array of contexts to pass to each function, or NULL. 944 - X - The coefficient vector u_h 945 946 Output Parameter: 947 . D - A Vec which holds the difference ||u - u_h||_2 for each cell 948 949 Level: developer 950 951 .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff() 952 @*/ 953 PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D) 954 { 955 PetscSection section; 956 PetscQuadrature quad; 957 Vec localX; 958 PetscScalar *funcVal, *interpolant; 959 PetscReal *coords, *detJ, *J; 960 const PetscReal *quadPoints, *quadWeights; 961 PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset; 962 PetscErrorCode ierr; 963 964 PetscFunctionBegin; 965 ierr = VecSet(D, 0.0);CHKERRQ(ierr); 966 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 967 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 968 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 969 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 970 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 971 ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 972 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 973 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 974 for (field = 0; field < numFields; ++field) { 975 PetscObject obj; 976 PetscClassId id; 977 PetscInt Nc; 978 979 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 980 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 981 if (id == PETSCFE_CLASSID) { 982 PetscFE fe = (PetscFE) obj; 983 984 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 985 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 986 } else if (id == PETSCFV_CLASSID) { 987 PetscFV fv = (PetscFV) obj; 988 989 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 990 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 991 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 992 numComponents += Nc; 993 } 994 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 995 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 996 ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 997 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 998 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 999 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1000 for (c = cStart; c < cEnd; ++c) { 1001 PetscScalar *x = NULL; 1002 PetscScalar elemDiff = 0.0; 1003 PetscInt qc = 0; 1004 1005 ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 1006 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1007 1008 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 1009 PetscObject obj; 1010 PetscClassId id; 1011 void * const ctx = ctxs ? ctxs[field] : NULL; 1012 PetscInt Nb, Nc, q, fc; 1013 1014 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1015 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1016 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 1017 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 1018 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1019 if (funcs[field]) { 1020 for (q = 0; q < Nq; ++q) { 1021 if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], c, q); 1022 ierr = (*funcs[field])(coordDim, time, &coords[q*coordDim], Nc, funcVal, ctx); 1023 if (ierr) { 1024 PetscErrorCode ierr2; 1025 ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 1026 ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2); 1027 ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 1028 CHKERRQ(ierr); 1029 } 1030 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1031 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 1032 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1033 for (fc = 0; fc < Nc; ++fc) { 1034 const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 1035 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 1036 } 1037 } 1038 } 1039 fieldOffset += Nb; 1040 qc += Nc; 1041 } 1042 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 1043 ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr); 1044 } 1045 ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 1046 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1047 ierr = VecSqrtAbs(D);CHKERRQ(ierr); 1048 PetscFunctionReturn(0); 1049 } 1050 1051 /*@C 1052 DMPlexComputeGradientClementInterpolant - This function computes the L2 projection of the cellwise gradient of a function u onto P1, and stores it in a Vec. 1053 1054 Input Parameters: 1055 + dm - The DM 1056 - LocX - The coefficient vector u_h 1057 1058 Output Parameter: 1059 . locC - A Vec which holds the Clement interpolant of the gradient 1060 1061 Notes: 1062 Add citation to (Clement, 1975) and definition of the interpolant 1063 \nabla u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| \nabla u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| where |T_i| is the cell volume 1064 1065 Level: developer 1066 1067 .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff() 1068 @*/ 1069 PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC) 1070 { 1071 DM_Plex *mesh = (DM_Plex *) dm->data; 1072 PetscInt debug = mesh->printFEM; 1073 DM dmC; 1074 PetscSection section; 1075 PetscQuadrature quad; 1076 PetscScalar *interpolant, *gradsum; 1077 PetscReal *coords, *detJ, *J, *invJ; 1078 const PetscReal *quadPoints, *quadWeights; 1079 PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, vStart, vEnd, v, field, fieldOffset; 1080 PetscErrorCode ierr; 1081 1082 PetscFunctionBegin; 1083 ierr = VecGetDM(locC, &dmC);CHKERRQ(ierr); 1084 ierr = VecSet(locC, 0.0);CHKERRQ(ierr); 1085 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1086 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 1087 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 1088 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 1089 for (field = 0; field < numFields; ++field) { 1090 PetscObject obj; 1091 PetscClassId id; 1092 PetscInt Nc; 1093 1094 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1095 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1096 if (id == PETSCFE_CLASSID) { 1097 PetscFE fe = (PetscFE) obj; 1098 1099 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1100 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1101 } else if (id == PETSCFV_CLASSID) { 1102 PetscFV fv = (PetscFV) obj; 1103 1104 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 1105 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1106 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1107 numComponents += Nc; 1108 } 1109 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1110 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 1111 ierr = PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ);CHKERRQ(ierr); 1112 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1113 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1114 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1115 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1116 for (v = vStart; v < vEnd; ++v) { 1117 PetscScalar volsum = 0.0; 1118 PetscInt *star = NULL; 1119 PetscInt starSize, st, d, fc; 1120 1121 ierr = PetscMemzero(gradsum, coordDim*numComponents * sizeof(PetscScalar));CHKERRQ(ierr); 1122 ierr = DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1123 for (st = 0; st < starSize*2; st += 2) { 1124 const PetscInt cell = star[st]; 1125 PetscScalar *grad = &gradsum[coordDim*numComponents]; 1126 PetscScalar *x = NULL; 1127 PetscReal vol = 0.0; 1128 1129 if ((cell < cStart) || (cell >= cEnd)) continue; 1130 ierr = DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, J, invJ, detJ);CHKERRQ(ierr); 1131 ierr = DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 1132 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 1133 PetscObject obj; 1134 PetscClassId id; 1135 PetscInt Nb, Nc, q, qc = 0; 1136 1137 ierr = PetscMemzero(grad, coordDim*numComponents * sizeof(PetscScalar));CHKERRQ(ierr); 1138 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1139 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1140 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 1141 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 1142 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1143 for (q = 0; q < Nq; ++q) { 1144 if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], cell, q); 1145 if (ierr) { 1146 PetscErrorCode ierr2; 1147 ierr2 = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr2); 1148 ierr2 = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr2); 1149 ierr2 = PetscFree4(interpolant,coords,detJ,J);CHKERRQ(ierr2); 1150 CHKERRQ(ierr); 1151 } 1152 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolateGradient_Static((PetscFE) obj, &x[fieldOffset], coordDim, invJ, NULL, q, interpolant);CHKERRQ(ierr);} 1153 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 1154 for (fc = 0; fc < Nc; ++fc) { 1155 const PetscReal wt = quadWeights[q*qNc+qc+fc]; 1156 1157 for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*detJ[q]; 1158 } 1159 vol += quadWeights[q*qNc]*detJ[q]; 1160 } 1161 fieldOffset += Nb; 1162 qc += Nc; 1163 } 1164 ierr = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 1165 for (fc = 0; fc < numComponents; ++fc) { 1166 for (d = 0; d < coordDim; ++d) { 1167 gradsum[fc*coordDim+d] += grad[fc*coordDim+d]; 1168 } 1169 } 1170 volsum += vol; 1171 if (debug) { 1172 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell);CHKERRQ(ierr); 1173 for (fc = 0; fc < numComponents; ++fc) { 1174 for (d = 0; d < coordDim; ++d) { 1175 if (fc || d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);} 1176 ierr = PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));CHKERRQ(ierr); 1177 } 1178 } 1179 ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr); 1180 } 1181 } 1182 for (fc = 0; fc < numComponents; ++fc) { 1183 for (d = 0; d < coordDim; ++d) gradsum[fc*coordDim+d] /= volsum; 1184 } 1185 ierr = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 1186 ierr = DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);CHKERRQ(ierr); 1187 } 1188 ierr = PetscFree6(gradsum,interpolant,coords,detJ,J,invJ);CHKERRQ(ierr); 1189 PetscFunctionReturn(0); 1190 } 1191 1192 static PetscErrorCode DMPlexComputeIntegral_Internal(DM dm, Vec X, PetscInt cStart, PetscInt cEnd, PetscScalar *cintegral, void *user) 1193 { 1194 DM dmAux = NULL; 1195 PetscDS prob, probAux = NULL; 1196 PetscSection section, sectionAux; 1197 Vec locX, locA; 1198 PetscInt dim, numCells = cEnd - cStart, c, f; 1199 PetscBool useFVM = PETSC_FALSE; 1200 /* DS */ 1201 PetscInt Nf, totDim, *uOff, *uOff_x, numConstants; 1202 PetscInt NfAux, totDimAux, *aOff; 1203 PetscScalar *u, *a; 1204 const PetscScalar *constants; 1205 /* Geometry */ 1206 PetscFEGeom *cgeomFEM; 1207 DM dmGrad; 1208 PetscQuadrature affineQuad = NULL; 1209 Vec cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL; 1210 PetscFVCellGeom *cgeomFVM; 1211 const PetscScalar *lgrad; 1212 PetscBool isAffine; 1213 DMField coordField; 1214 IS cellIS; 1215 PetscErrorCode ierr; 1216 1217 PetscFunctionBegin; 1218 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1219 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1220 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 1221 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1222 /* Determine which discretizations we have */ 1223 for (f = 0; f < Nf; ++f) { 1224 PetscObject obj; 1225 PetscClassId id; 1226 1227 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1228 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1229 if (id == PETSCFV_CLASSID) useFVM = PETSC_TRUE; 1230 } 1231 /* Get local solution with boundary values */ 1232 ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 1233 ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1234 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 1235 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 1236 /* Read DS information */ 1237 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1238 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 1239 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 1240 ierr = ISCreateStride(PETSC_COMM_SELF,numCells,cStart,1,&cellIS);CHKERRQ(ierr); 1241 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 1242 /* Read Auxiliary DS information */ 1243 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1244 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1245 if (dmAux) { 1246 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1247 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 1248 ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 1249 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1250 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 1251 } 1252 /* Allocate data arrays */ 1253 ierr = PetscCalloc1(numCells*totDim, &u);CHKERRQ(ierr); 1254 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1255 /* Read out geometry */ 1256 ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 1257 ierr = DMFieldGetFEInvariance(coordField,cellIS,NULL,&isAffine,NULL);CHKERRQ(ierr); 1258 if (isAffine) { 1259 ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);CHKERRQ(ierr); 1260 if (affineQuad) { 1261 ierr = DMFieldCreateFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 1262 } 1263 } 1264 if (useFVM) { 1265 PetscFV fv = NULL; 1266 Vec grad; 1267 PetscInt fStart, fEnd; 1268 PetscBool compGrad; 1269 1270 for (f = 0; f < Nf; ++f) { 1271 PetscObject obj; 1272 PetscClassId id; 1273 1274 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1275 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1276 if (id == PETSCFV_CLASSID) {fv = (PetscFV) obj; break;} 1277 } 1278 ierr = PetscFVGetComputeGradients(fv, &compGrad);CHKERRQ(ierr); 1279 ierr = PetscFVSetComputeGradients(fv, PETSC_TRUE);CHKERRQ(ierr); 1280 ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr); 1281 ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr); 1282 ierr = PetscFVSetComputeGradients(fv, compGrad);CHKERRQ(ierr); 1283 ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1284 /* Reconstruct and limit cell gradients */ 1285 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1286 ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1287 ierr = DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr); 1288 /* Communicate gradient values */ 1289 ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1290 ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1291 ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1292 ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1293 /* Handle non-essential (e.g. outflow) boundary values */ 1294 ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr); 1295 ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1296 } 1297 /* Read out data from inputs */ 1298 for (c = cStart; c < cEnd; ++c) { 1299 PetscScalar *x = NULL; 1300 PetscInt i; 1301 1302 ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 1303 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1304 ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 1305 if (dmAux) { 1306 ierr = DMPlexVecGetClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 1307 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1308 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 1309 } 1310 } 1311 /* Do integration for each field */ 1312 for (f = 0; f < Nf; ++f) { 1313 PetscObject obj; 1314 PetscClassId id; 1315 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1316 1317 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1318 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1319 if (id == PETSCFE_CLASSID) { 1320 PetscFE fe = (PetscFE) obj; 1321 PetscQuadrature q; 1322 PetscFEGeom *chunkGeom = NULL; 1323 PetscInt Nq, Nb; 1324 1325 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1326 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1327 ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1328 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1329 blockSize = Nb*Nq; 1330 batchSize = numBlocks * blockSize; 1331 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1332 numChunks = numCells / (numBatches*batchSize); 1333 Ne = numChunks*numBatches*batchSize; 1334 Nr = numCells % (numBatches*batchSize); 1335 offset = numCells - Nr; 1336 if (!affineQuad) { 1337 ierr = DMFieldCreateFEGeom(coordField,cellIS,q,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 1338 } 1339 ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr); 1340 ierr = PetscFEIntegrate(fe, prob, f, Ne, chunkGeom, u, probAux, a, cintegral);CHKERRQ(ierr); 1341 ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1342 ierr = PetscFEIntegrate(fe, prob, f, Nr, chunkGeom, &u[offset*totDim], probAux, &a[offset*totDimAux], &cintegral[offset*Nf]);CHKERRQ(ierr); 1343 ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1344 if (!affineQuad) { 1345 ierr = PetscFEGeomDestroy(&cgeomFEM);CHKERRQ(ierr); 1346 } 1347 } else if (id == PETSCFV_CLASSID) { 1348 PetscInt foff; 1349 PetscPointFunc obj_func; 1350 PetscScalar lint; 1351 1352 ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr); 1353 ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1354 if (obj_func) { 1355 for (c = 0; c < numCells; ++c) { 1356 PetscScalar *u_x; 1357 1358 ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr); 1359 obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, &a[totDimAux*c], NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint); 1360 cintegral[c*Nf+f] += PetscRealPart(lint)*cgeomFVM[c].volume; 1361 } 1362 } 1363 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1364 } 1365 /* Cleanup data arrays */ 1366 if (useFVM) { 1367 ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1368 ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1369 ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1370 ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr); 1371 ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr); 1372 ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 1373 } 1374 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1375 ierr = PetscFree(u);CHKERRQ(ierr); 1376 /* Cleanup */ 1377 if (affineQuad) { 1378 ierr = PetscFEGeomDestroy(&cgeomFEM);CHKERRQ(ierr); 1379 } 1380 ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr); 1381 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 1382 ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 1383 PetscFunctionReturn(0); 1384 } 1385 1386 /*@ 1387 DMPlexComputeIntegralFEM - Form the integral over the domain from the global input X using pointwise functions specified by the user 1388 1389 Input Parameters: 1390 + dm - The mesh 1391 . X - Global input vector 1392 - user - The user context 1393 1394 Output Parameter: 1395 . integral - Integral for each field 1396 1397 Level: developer 1398 1399 .seealso: DMPlexComputeResidualFEM() 1400 @*/ 1401 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscScalar *integral, void *user) 1402 { 1403 DM_Plex *mesh = (DM_Plex *) dm->data; 1404 PetscScalar *cintegral, *lintegral; 1405 PetscInt Nf, f, cellHeight, cStart, cEnd, cEndInterior[4], cell; 1406 PetscErrorCode ierr; 1407 1408 PetscFunctionBegin; 1409 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1410 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1411 PetscValidPointer(integral, 3); 1412 ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1413 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 1414 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 1415 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 1416 ierr = DMPlexGetHybridBounds(dm, &cEndInterior[0], &cEndInterior[1], &cEndInterior[2], &cEndInterior[3]);CHKERRQ(ierr); 1417 cEnd = cEndInterior[cellHeight] < 0 ? cEnd : cEndInterior[cellHeight]; 1418 /* TODO Introduce a loop over large chunks (right now this is a single chunk) */ 1419 ierr = PetscCalloc2(Nf, &lintegral, (cEnd-cStart)*Nf, &cintegral);CHKERRQ(ierr); 1420 ierr = DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);CHKERRQ(ierr); 1421 /* Sum up values */ 1422 for (cell = cStart; cell < cEnd; ++cell) { 1423 const PetscInt c = cell - cStart; 1424 1425 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);CHKERRQ(ierr);} 1426 for (f = 0; f < Nf; ++f) lintegral[f] += cintegral[c*Nf+f]; 1427 } 1428 ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 1429 if (mesh->printFEM) { 1430 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Integral:");CHKERRQ(ierr); 1431 for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", (double) PetscRealPart(integral[f]));CHKERRQ(ierr);} 1432 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 1433 } 1434 ierr = PetscFree2(lintegral, cintegral);CHKERRQ(ierr); 1435 ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1436 PetscFunctionReturn(0); 1437 } 1438 1439 /*@ 1440 DMPlexComputeCellwiseIntegralFEM - Form the vector of cellwise integrals F from the global input X using pointwise functions specified by the user 1441 1442 Input Parameters: 1443 + dm - The mesh 1444 . X - Global input vector 1445 - user - The user context 1446 1447 Output Parameter: 1448 . integral - Cellwise integrals for each field 1449 1450 Level: developer 1451 1452 .seealso: DMPlexComputeResidualFEM() 1453 @*/ 1454 PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM dm, Vec X, Vec F, void *user) 1455 { 1456 DM_Plex *mesh = (DM_Plex *) dm->data; 1457 DM dmF; 1458 PetscSection sectionF; 1459 PetscScalar *cintegral, *af; 1460 PetscInt Nf, f, cellHeight, cStart, cEnd, cEndInterior[4], cell; 1461 PetscErrorCode ierr; 1462 1463 PetscFunctionBegin; 1464 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1465 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1466 PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 1467 ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1468 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 1469 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 1470 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 1471 ierr = DMPlexGetHybridBounds(dm, &cEndInterior[0], &cEndInterior[1], &cEndInterior[2], &cEndInterior[3]);CHKERRQ(ierr); 1472 cEnd = cEndInterior[cellHeight] < 0 ? cEnd : cEndInterior[cellHeight]; 1473 /* TODO Introduce a loop over large chunks (right now this is a single chunk) */ 1474 ierr = PetscCalloc1((cEnd-cStart)*Nf, &cintegral);CHKERRQ(ierr); 1475 ierr = DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);CHKERRQ(ierr); 1476 /* Put values in F*/ 1477 ierr = VecGetDM(F, &dmF);CHKERRQ(ierr); 1478 ierr = DMGetSection(dmF, §ionF);CHKERRQ(ierr); 1479 ierr = VecGetArray(F, &af);CHKERRQ(ierr); 1480 for (cell = cStart; cell < cEnd; ++cell) { 1481 const PetscInt c = cell - cStart; 1482 PetscInt dof, off; 1483 1484 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);CHKERRQ(ierr);} 1485 ierr = PetscSectionGetDof(sectionF, cell, &dof);CHKERRQ(ierr); 1486 ierr = PetscSectionGetOffset(sectionF, cell, &off);CHKERRQ(ierr); 1487 if (dof != Nf) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cell dofs %D != %D", dof, Nf); 1488 for (f = 0; f < Nf; ++f) af[off+f] = cintegral[c*Nf+f]; 1489 } 1490 ierr = VecRestoreArray(F, &af);CHKERRQ(ierr); 1491 ierr = PetscFree(cintegral);CHKERRQ(ierr); 1492 ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1493 PetscFunctionReturn(0); 1494 } 1495 1496 static PetscErrorCode DMPlexComputeBdIntegral_Internal(DM dm, Vec locX, IS pointIS, 1497 void (*func)(PetscInt, PetscInt, PetscInt, 1498 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1499 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1500 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 1501 PetscScalar *fintegral, void *user) 1502 { 1503 DM plex = NULL, plexA = NULL; 1504 PetscDS prob, probAux = NULL; 1505 PetscSection section, sectionAux = NULL; 1506 Vec locA = NULL; 1507 DMField coordField; 1508 PetscInt Nf, totDim, *uOff, *uOff_x; 1509 PetscInt NfAux = 0, totDimAux = 0, *aOff = NULL; 1510 PetscScalar *u, *a = NULL; 1511 const PetscScalar *constants; 1512 PetscInt numConstants, f; 1513 PetscErrorCode ierr; 1514 1515 PetscFunctionBegin; 1516 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1517 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1518 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1519 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1520 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1521 /* Determine which discretizations we have */ 1522 for (f = 0; f < Nf; ++f) { 1523 PetscObject obj; 1524 PetscClassId id; 1525 1526 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1527 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1528 if (id == PETSCFV_CLASSID) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not supported for FVM (field %D)", f); 1529 } 1530 /* Read DS information */ 1531 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1532 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 1533 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 1534 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 1535 /* Read Auxiliary DS information */ 1536 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1537 if (locA) { 1538 DM dmAux; 1539 1540 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1541 ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr); 1542 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1543 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 1544 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1545 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1546 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 1547 } 1548 /* Integrate over points */ 1549 { 1550 PetscFEGeom *fgeom, *chunkGeom = NULL; 1551 PetscBool isAffine; 1552 PetscQuadrature qGeom = NULL; 1553 const PetscInt *points; 1554 PetscInt numFaces, face, Nq, field; 1555 PetscInt numChunks, chunkSize, chunk, Nr, offset; 1556 1557 ierr = ISGetLocalSize(pointIS, &numFaces);CHKERRQ(ierr); 1558 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1559 ierr = PetscCalloc2(numFaces*totDim, &u, locA ? numFaces*totDimAux : 0, &a);CHKERRQ(ierr); 1560 ierr = DMFieldGetFEInvariance(coordField, pointIS, NULL, &isAffine, NULL);CHKERRQ(ierr); 1561 for (field = 0; field < Nf; ++field) { 1562 PetscFE fe; 1563 1564 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 1565 if (isAffine) {ierr = DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom);CHKERRQ(ierr);} 1566 if (!qGeom) { 1567 ierr = PetscFEGetFaceQuadrature(fe, &qGeom);CHKERRQ(ierr); 1568 ierr = PetscObjectReference((PetscObject) qGeom);CHKERRQ(ierr); 1569 } 1570 ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1571 ierr = DMPlexGetFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom);CHKERRQ(ierr); 1572 for (face = 0; face < numFaces; ++face) { 1573 const PetscInt point = points[face], *support, *cone; 1574 PetscScalar *x = NULL; 1575 PetscInt i, coneSize, faceLoc; 1576 1577 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1578 ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 1579 ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 1580 for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break; 1581 if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", face, support[0]); 1582 fgeom->face[face][0] = faceLoc; 1583 ierr = DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1584 for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i]; 1585 ierr = DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1586 if (locA) { 1587 PetscInt subp; 1588 ierr = DMPlexGetSubpoint(plexA, support[0], &subp);CHKERRQ(ierr); 1589 ierr = DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1590 for (i = 0; i < totDimAux; ++i) a[f*totDimAux+i] = x[i]; 1591 ierr = DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1592 } 1593 } 1594 /* Get blocking */ 1595 { 1596 PetscQuadrature q; 1597 PetscInt numBatches, batchSize, numBlocks, blockSize; 1598 PetscInt Nq, Nb; 1599 1600 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1601 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1602 ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1603 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1604 blockSize = Nb*Nq; 1605 batchSize = numBlocks * blockSize; 1606 chunkSize = numBatches*batchSize; 1607 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1608 numChunks = numFaces / chunkSize; 1609 Nr = numFaces % chunkSize; 1610 offset = numFaces - Nr; 1611 } 1612 /* Do integration for each field */ 1613 for (chunk = 0; chunk < numChunks; ++chunk) { 1614 ierr = PetscFEGeomGetChunk(fgeom, chunk*chunkSize, (chunk+1)*chunkSize, &chunkGeom);CHKERRQ(ierr); 1615 ierr = PetscFEIntegrateBd(fe, prob, field, func, chunkSize, chunkGeom, u, probAux, a, fintegral);CHKERRQ(ierr); 1616 ierr = PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);CHKERRQ(ierr); 1617 } 1618 ierr = PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom);CHKERRQ(ierr); 1619 ierr = PetscFEIntegrateBd(fe, prob, field, func, Nr, chunkGeom, &u[offset*totDim], probAux, a ? &a[offset*totDimAux] : NULL, &fintegral[offset*Nf]);CHKERRQ(ierr); 1620 ierr = PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom);CHKERRQ(ierr); 1621 /* Cleanup data arrays */ 1622 ierr = DMPlexRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom);CHKERRQ(ierr); 1623 ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr); 1624 ierr = PetscFree2(u, a);CHKERRQ(ierr); 1625 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1626 } 1627 } 1628 if (plex) {ierr = DMDestroy(&plex);CHKERRQ(ierr);} 1629 if (plexA) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);} 1630 PetscFunctionReturn(0); 1631 } 1632 1633 /*@ 1634 DMPlexComputeBdIntegral - Form the integral over the specified boundary from the global input X using pointwise functions specified by the user 1635 1636 Input Parameters: 1637 + dm - The mesh 1638 . X - Global input vector 1639 . label - The boundary DMLabel 1640 . numVals - The number of label values to use, or PETSC_DETERMINE for all values 1641 . vals - The label values to use, or PETSC_NULL for all values 1642 . func = The function to integrate along the boundary 1643 - user - The user context 1644 1645 Output Parameter: 1646 . integral - Integral for each field 1647 1648 Level: developer 1649 1650 .seealso: DMPlexComputeIntegralFEM(), DMPlexComputeBdResidualFEM() 1651 @*/ 1652 PetscErrorCode DMPlexComputeBdIntegral(DM dm, Vec X, DMLabel label, PetscInt numVals, const PetscInt vals[], 1653 void (*func)(PetscInt, PetscInt, PetscInt, 1654 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1655 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1656 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 1657 PetscScalar *integral, void *user) 1658 { 1659 Vec locX; 1660 PetscSection section; 1661 DMLabel depthLabel; 1662 IS facetIS; 1663 PetscInt dim, Nf, f, v; 1664 PetscErrorCode ierr; 1665 1666 PetscFunctionBegin; 1667 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1668 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1669 PetscValidPointer(label, 3); 1670 if (vals) PetscValidPointer(vals, 5); 1671 PetscValidPointer(integral, 6); 1672 ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1673 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1674 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1675 ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr); 1676 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1677 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1678 /* Get local solution with boundary values */ 1679 ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 1680 ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1681 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 1682 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 1683 /* Loop over label values */ 1684 ierr = PetscMemzero(integral, Nf * sizeof(PetscScalar));CHKERRQ(ierr); 1685 for (v = 0; v < numVals; ++v) { 1686 IS pointIS; 1687 PetscInt numFaces, face; 1688 PetscScalar *fintegral; 1689 1690 ierr = DMLabelGetStratumIS(label, vals[v], &pointIS);CHKERRQ(ierr); 1691 if (!pointIS) continue; /* No points with that id on this process */ 1692 { 1693 IS isectIS; 1694 1695 /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */ 1696 ierr = ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS);CHKERRQ(ierr); 1697 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1698 pointIS = isectIS; 1699 } 1700 ierr = ISGetLocalSize(pointIS, &numFaces);CHKERRQ(ierr); 1701 ierr = PetscCalloc1(numFaces*Nf, &fintegral);CHKERRQ(ierr); 1702 ierr = DMPlexComputeBdIntegral_Internal(dm, locX, pointIS, func, fintegral, user);CHKERRQ(ierr); 1703 /* Sum point contributions into integral */ 1704 for (f = 0; f < Nf; ++f) for (face = 0; face < numFaces; ++face) integral[f] += fintegral[face*Nf+f]; 1705 ierr = PetscFree(fintegral);CHKERRQ(ierr); 1706 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1707 } 1708 ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 1709 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 1710 ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 1711 PetscFunctionReturn(0); 1712 } 1713 1714 /*@ 1715 DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1716 1717 Input Parameters: 1718 + dmf - The fine mesh 1719 . dmc - The coarse mesh 1720 - user - The user context 1721 1722 Output Parameter: 1723 . In - The interpolation matrix 1724 1725 Level: developer 1726 1727 .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM() 1728 @*/ 1729 PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user) 1730 { 1731 DM_Plex *mesh = (DM_Plex *) dmc->data; 1732 const char *name = "Interpolator"; 1733 PetscDS prob; 1734 PetscFE *feRef; 1735 PetscFV *fvRef; 1736 PetscSection fsection, fglobalSection; 1737 PetscSection csection, cglobalSection; 1738 PetscScalar *elemMat; 1739 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c; 1740 PetscInt cTotDim, rTotDim = 0; 1741 PetscErrorCode ierr; 1742 1743 PetscFunctionBegin; 1744 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1745 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1746 ierr = DMGetSection(dmf, &fsection);CHKERRQ(ierr); 1747 ierr = DMGetGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1748 ierr = DMGetSection(dmc, &csection);CHKERRQ(ierr); 1749 ierr = DMGetGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1750 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1751 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1752 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1753 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1754 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1755 ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr); 1756 for (f = 0; f < Nf; ++f) { 1757 PetscObject obj; 1758 PetscClassId id; 1759 PetscInt rNb = 0, Nc = 0; 1760 1761 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1762 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1763 if (id == PETSCFE_CLASSID) { 1764 PetscFE fe = (PetscFE) obj; 1765 1766 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1767 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1768 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1769 } else if (id == PETSCFV_CLASSID) { 1770 PetscFV fv = (PetscFV) obj; 1771 PetscDualSpace Q; 1772 1773 ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr); 1774 ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr); 1775 ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr); 1776 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 1777 } 1778 rTotDim += rNb; 1779 } 1780 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1781 ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 1782 ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1783 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1784 PetscDualSpace Qref; 1785 PetscQuadrature f; 1786 const PetscReal *qpoints, *qweights; 1787 PetscReal *points; 1788 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1789 1790 /* Compose points from all dual basis functionals */ 1791 if (feRef[fieldI]) { 1792 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1793 ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1794 } else { 1795 ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr); 1796 ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr); 1797 } 1798 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1799 for (i = 0; i < fpdim; ++i) { 1800 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1801 ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1802 npoints += Np; 1803 } 1804 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1805 for (i = 0, k = 0; i < fpdim; ++i) { 1806 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1807 ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1808 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1809 } 1810 1811 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1812 PetscObject obj; 1813 PetscClassId id; 1814 PetscReal *B; 1815 PetscInt NcJ = 0, cpdim = 0, j, qNc; 1816 1817 ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr); 1818 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1819 if (id == PETSCFE_CLASSID) { 1820 PetscFE fe = (PetscFE) obj; 1821 1822 /* Evaluate basis at points */ 1823 ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 1824 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1825 /* For now, fields only interpolate themselves */ 1826 if (fieldI == fieldJ) { 1827 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); 1828 ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1829 for (i = 0, k = 0; i < fpdim; ++i) { 1830 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1831 ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr); 1832 if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ); 1833 for (p = 0; p < Np; ++p, ++k) { 1834 for (j = 0; j < cpdim; ++j) { 1835 /* 1836 cTotDim: Total columns in element interpolation matrix, sum of number of dual basis functionals in each field 1837 offsetI, offsetJ: Offsets into the larger element interpolation matrix for different fields 1838 fpdim, i, cpdim, j: Dofs for fine and coarse grids, correspond to dual space basis functionals 1839 qNC, Nc, Ncj, c: Number of components in this field 1840 Np, p: Number of quad points in the fine grid functional i 1841 k: i*Np + p, overall point number for the interpolation 1842 */ 1843 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c]; 1844 } 1845 } 1846 } 1847 ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1848 } 1849 } else if (id == PETSCFV_CLASSID) { 1850 PetscFV fv = (PetscFV) obj; 1851 1852 /* Evaluate constant function at points */ 1853 ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr); 1854 cpdim = 1; 1855 /* For now, fields only interpolate themselves */ 1856 if (fieldI == fieldJ) { 1857 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); 1858 for (i = 0, k = 0; i < fpdim; ++i) { 1859 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1860 ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr); 1861 if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ); 1862 for (p = 0; p < Np; ++p, ++k) { 1863 for (j = 0; j < cpdim; ++j) { 1864 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p*qNc+c]; 1865 } 1866 } 1867 } 1868 } 1869 } 1870 offsetJ += cpdim; 1871 } 1872 offsetI += fpdim; 1873 ierr = PetscFree(points);CHKERRQ(ierr); 1874 } 1875 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 1876 /* Preallocate matrix */ 1877 { 1878 Mat preallocator; 1879 PetscScalar *vals; 1880 PetscInt *cellCIndices, *cellFIndices; 1881 PetscInt locRows, locCols, cell; 1882 1883 ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr); 1884 ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr); 1885 ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr); 1886 ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1887 ierr = MatSetUp(preallocator);CHKERRQ(ierr); 1888 ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 1889 for (cell = cStart; cell < cEnd; ++cell) { 1890 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 1891 ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr); 1892 } 1893 ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr); 1894 ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1895 ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1896 ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr); 1897 ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 1898 } 1899 /* Fill matrix */ 1900 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1901 for (c = cStart; c < cEnd; ++c) { 1902 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1903 } 1904 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1905 ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr); 1906 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1907 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1908 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1909 if (mesh->printFEM) { 1910 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1911 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1912 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1913 } 1914 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1915 PetscFunctionReturn(0); 1916 } 1917 1918 PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user) 1919 { 1920 SETERRQ(PetscObjectComm((PetscObject) dmc), PETSC_ERR_SUP, "Laziness"); 1921 } 1922 1923 /*@ 1924 DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM. 1925 1926 Input Parameters: 1927 + dmf - The fine mesh 1928 . dmc - The coarse mesh 1929 - user - The user context 1930 1931 Output Parameter: 1932 . In - The interpolation matrix 1933 1934 Level: developer 1935 1936 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM() 1937 @*/ 1938 PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user) 1939 { 1940 DM_Plex *mesh = (DM_Plex *) dmf->data; 1941 const char *name = "Interpolator"; 1942 PetscDS prob; 1943 PetscSection fsection, csection, globalFSection, globalCSection; 1944 PetscHSetIJ ht; 1945 PetscLayout rLayout; 1946 PetscInt *dnz, *onz; 1947 PetscInt locRows, rStart, rEnd; 1948 PetscReal *x, *v0, *J, *invJ, detJ; 1949 PetscReal *v0c, *Jc, *invJc, detJc; 1950 PetscScalar *elemMat; 1951 PetscInt dim, Nf, field, totDim, cStart, cEnd, cell, ccell; 1952 PetscErrorCode ierr; 1953 1954 PetscFunctionBegin; 1955 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1956 ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr); 1957 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 1958 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 1959 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1960 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 1961 ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr); 1962 ierr = DMGetSection(dmf, &fsection);CHKERRQ(ierr); 1963 ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 1964 ierr = DMGetSection(dmc, &csection);CHKERRQ(ierr); 1965 ierr = DMGetGlobalSection(dmc, &globalCSection);CHKERRQ(ierr); 1966 ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 1967 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1968 ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr); 1969 1970 ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 1971 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 1972 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 1973 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 1974 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 1975 ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 1976 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 1977 ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr); 1978 ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr); 1979 for (field = 0; field < Nf; ++field) { 1980 PetscObject obj; 1981 PetscClassId id; 1982 PetscDualSpace Q = NULL; 1983 PetscQuadrature f; 1984 const PetscReal *qpoints; 1985 PetscInt Nc, Np, fpdim, i, d; 1986 1987 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 1988 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1989 if (id == PETSCFE_CLASSID) { 1990 PetscFE fe = (PetscFE) obj; 1991 1992 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1993 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1994 } else if (id == PETSCFV_CLASSID) { 1995 PetscFV fv = (PetscFV) obj; 1996 1997 ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr); 1998 Nc = 1; 1999 } 2000 ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr); 2001 /* For each fine grid cell */ 2002 for (cell = cStart; cell < cEnd; ++cell) { 2003 PetscInt *findices, *cindices; 2004 PetscInt numFIndices, numCIndices; 2005 2006 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2007 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 2008 if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim); 2009 for (i = 0; i < fpdim; ++i) { 2010 Vec pointVec; 2011 PetscScalar *pV; 2012 PetscSF coarseCellSF = NULL; 2013 const PetscSFNode *coarseCells; 2014 PetscInt numCoarseCells, q, c; 2015 2016 /* Get points from the dual basis functional quadrature */ 2017 ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr); 2018 ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 2019 ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr); 2020 ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 2021 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2022 for (q = 0; q < Np; ++q) { 2023 const PetscReal xi0[3] = {-1., -1., -1.}; 2024 2025 /* Transform point to real space */ 2026 CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x); 2027 for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 2028 } 2029 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2030 /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 2031 /* OPT: Pack all quad points from fine cell */ 2032 ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 2033 ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr); 2034 /* Update preallocation info */ 2035 ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 2036 if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 2037 { 2038 PetscHashIJKey key; 2039 PetscBool missing; 2040 2041 key.i = findices[i]; 2042 if (key.i >= 0) { 2043 /* Get indices for coarse elements */ 2044 for (ccell = 0; ccell < numCoarseCells; ++ccell) { 2045 ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2046 for (c = 0; c < numCIndices; ++c) { 2047 key.j = cindices[c]; 2048 if (key.j < 0) continue; 2049 ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr); 2050 if (missing) { 2051 if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i-rStart]; 2052 else ++onz[key.i-rStart]; 2053 } 2054 } 2055 ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2056 } 2057 } 2058 } 2059 ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 2060 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 2061 } 2062 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2063 } 2064 } 2065 ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr); 2066 ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 2067 ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 2068 ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 2069 for (field = 0; field < Nf; ++field) { 2070 PetscObject obj; 2071 PetscClassId id; 2072 PetscDualSpace Q = NULL; 2073 PetscQuadrature f; 2074 const PetscReal *qpoints, *qweights; 2075 PetscInt Nc, qNc, Np, fpdim, i, d; 2076 2077 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 2078 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2079 if (id == PETSCFE_CLASSID) { 2080 PetscFE fe = (PetscFE) obj; 2081 2082 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 2083 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 2084 } else if (id == PETSCFV_CLASSID) { 2085 PetscFV fv = (PetscFV) obj; 2086 2087 ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr); 2088 Nc = 1; 2089 } else SETERRQ1(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %d",field); 2090 ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr); 2091 /* For each fine grid cell */ 2092 for (cell = cStart; cell < cEnd; ++cell) { 2093 PetscInt *findices, *cindices; 2094 PetscInt numFIndices, numCIndices; 2095 2096 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2097 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 2098 if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim); 2099 for (i = 0; i < fpdim; ++i) { 2100 Vec pointVec; 2101 PetscScalar *pV; 2102 PetscSF coarseCellSF = NULL; 2103 const PetscSFNode *coarseCells; 2104 PetscInt numCoarseCells, cpdim, q, c, j; 2105 2106 /* Get points from the dual basis functional quadrature */ 2107 ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr); 2108 ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);CHKERRQ(ierr); 2109 if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, Nc); 2110 ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr); 2111 ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 2112 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2113 for (q = 0; q < Np; ++q) { 2114 const PetscReal xi0[3] = {-1., -1., -1.}; 2115 2116 /* Transform point to real space */ 2117 CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x); 2118 for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 2119 } 2120 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2121 /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 2122 /* OPT: Read this out from preallocation information */ 2123 ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 2124 /* Update preallocation info */ 2125 ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 2126 if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 2127 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2128 for (ccell = 0; ccell < numCoarseCells; ++ccell) { 2129 PetscReal pVReal[3]; 2130 const PetscReal xi0[3] = {-1., -1., -1.}; 2131 2132 ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2133 /* Transform points from real space to coarse reference space */ 2134 ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr); 2135 for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]); 2136 CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x); 2137 2138 if (id == PETSCFE_CLASSID) { 2139 PetscFE fe = (PetscFE) obj; 2140 PetscReal *B; 2141 2142 /* Evaluate coarse basis on contained point */ 2143 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 2144 ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr); 2145 ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr); 2146 /* Get elemMat entries by multiplying by weight */ 2147 for (j = 0; j < cpdim; ++j) { 2148 for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*qweights[ccell*qNc + c]; 2149 } 2150 ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 2151 } else { 2152 cpdim = 1; 2153 for (j = 0; j < cpdim; ++j) { 2154 for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c]; 2155 } 2156 } 2157 /* Update interpolator */ 2158 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 2159 if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim); 2160 ierr = MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr); 2161 ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2162 } 2163 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2164 ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 2165 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 2166 } 2167 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2168 } 2169 } 2170 ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 2171 ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr); 2172 ierr = PetscFree(elemMat);CHKERRQ(ierr); 2173 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2174 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2175 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 2176 PetscFunctionReturn(0); 2177 } 2178 2179 /*@ 2180 DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix M from the coarse DM to a non-nested fine DM. 2181 2182 Input Parameters: 2183 + dmf - The fine mesh 2184 . dmc - The coarse mesh 2185 - user - The user context 2186 2187 Output Parameter: 2188 . mass - The mass matrix 2189 2190 Level: developer 2191 2192 .seealso: DMPlexComputeMassMatrixNested(), DMPlexComputeInterpolatorNested(), DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM() 2193 @*/ 2194 PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user) 2195 { 2196 DM_Plex *mesh = (DM_Plex *) dmf->data; 2197 const char *name = "Mass Matrix"; 2198 PetscDS prob; 2199 PetscSection fsection, csection, globalFSection, globalCSection; 2200 PetscHSetIJ ht; 2201 PetscLayout rLayout; 2202 PetscInt *dnz, *onz; 2203 PetscInt locRows, rStart, rEnd; 2204 PetscReal *x, *v0, *J, *invJ, detJ; 2205 PetscReal *v0c, *Jc, *invJc, detJc; 2206 PetscScalar *elemMat; 2207 PetscInt dim, Nf, field, totDim, cStart, cEnd, cell, ccell; 2208 PetscErrorCode ierr; 2209 2210 PetscFunctionBegin; 2211 ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr); 2212 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 2213 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 2214 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 2215 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 2216 ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr); 2217 ierr = DMGetSection(dmf, &fsection);CHKERRQ(ierr); 2218 ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 2219 ierr = DMGetSection(dmc, &csection);CHKERRQ(ierr); 2220 ierr = DMGetGlobalSection(dmc, &globalCSection);CHKERRQ(ierr); 2221 ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 2222 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2223 ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr); 2224 2225 ierr = MatGetLocalSize(mass, &locRows, NULL);CHKERRQ(ierr); 2226 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) mass), &rLayout);CHKERRQ(ierr); 2227 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 2228 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 2229 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 2230 ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 2231 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 2232 ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr); 2233 ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr); 2234 for (field = 0; field < Nf; ++field) { 2235 PetscObject obj; 2236 PetscClassId id; 2237 PetscQuadrature quad; 2238 const PetscReal *qpoints; 2239 PetscInt Nq, Nc, i, d; 2240 2241 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 2242 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2243 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);} 2244 else {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);} 2245 ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL);CHKERRQ(ierr); 2246 /* For each fine grid cell */ 2247 for (cell = cStart; cell < cEnd; ++cell) { 2248 Vec pointVec; 2249 PetscScalar *pV; 2250 PetscSF coarseCellSF = NULL; 2251 const PetscSFNode *coarseCells; 2252 PetscInt numCoarseCells, q, c; 2253 PetscInt *findices, *cindices; 2254 PetscInt numFIndices, numCIndices; 2255 2256 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2257 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 2258 /* Get points from the quadrature */ 2259 ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr); 2260 ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 2261 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2262 for (q = 0; q < Nq; ++q) { 2263 const PetscReal xi0[3] = {-1., -1., -1.}; 2264 2265 /* Transform point to real space */ 2266 CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x); 2267 for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 2268 } 2269 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2270 /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 2271 ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 2272 ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr); 2273 /* Update preallocation info */ 2274 ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 2275 if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 2276 { 2277 PetscHashIJKey key; 2278 PetscBool missing; 2279 2280 for (i = 0; i < numFIndices; ++i) { 2281 key.i = findices[i]; 2282 if (key.i >= 0) { 2283 /* Get indices for coarse elements */ 2284 for (ccell = 0; ccell < numCoarseCells; ++ccell) { 2285 ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2286 for (c = 0; c < numCIndices; ++c) { 2287 key.j = cindices[c]; 2288 if (key.j < 0) continue; 2289 ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr); 2290 if (missing) { 2291 if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i-rStart]; 2292 else ++onz[key.i-rStart]; 2293 } 2294 } 2295 ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2296 } 2297 } 2298 } 2299 } 2300 ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 2301 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 2302 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2303 } 2304 } 2305 ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr); 2306 ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 2307 ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 2308 ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 2309 for (field = 0; field < Nf; ++field) { 2310 PetscObject obj; 2311 PetscClassId id; 2312 PetscQuadrature quad; 2313 PetscReal *Bfine; 2314 const PetscReal *qpoints, *qweights; 2315 PetscInt Nq, Nc, i, d; 2316 2317 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 2318 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2319 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);ierr = PetscFEGetDefaultTabulation((PetscFE) obj, &Bfine, NULL, NULL);CHKERRQ(ierr);} 2320 else {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);} 2321 ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights);CHKERRQ(ierr); 2322 /* For each fine grid cell */ 2323 for (cell = cStart; cell < cEnd; ++cell) { 2324 Vec pointVec; 2325 PetscScalar *pV; 2326 PetscSF coarseCellSF = NULL; 2327 const PetscSFNode *coarseCells; 2328 PetscInt numCoarseCells, cpdim, q, c, j; 2329 PetscInt *findices, *cindices; 2330 PetscInt numFIndices, numCIndices; 2331 2332 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2333 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 2334 /* Get points from the quadrature */ 2335 ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr); 2336 ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 2337 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2338 for (q = 0; q < Nq; ++q) { 2339 const PetscReal xi0[3] = {-1., -1., -1.}; 2340 2341 /* Transform point to real space */ 2342 CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x); 2343 for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 2344 } 2345 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2346 /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 2347 ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 2348 /* Update matrix */ 2349 ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 2350 if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 2351 ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 2352 for (ccell = 0; ccell < numCoarseCells; ++ccell) { 2353 PetscReal pVReal[3]; 2354 const PetscReal xi0[3] = {-1., -1., -1.}; 2355 2356 2357 ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2358 /* Transform points from real space to coarse reference space */ 2359 ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr); 2360 for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]); 2361 CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x); 2362 2363 if (id == PETSCFE_CLASSID) { 2364 PetscFE fe = (PetscFE) obj; 2365 PetscReal *B; 2366 2367 /* Evaluate coarse basis on contained point */ 2368 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 2369 ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr); 2370 /* Get elemMat entries by multiplying by weight */ 2371 for (i = 0; i < numFIndices; ++i) { 2372 ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr); 2373 for (j = 0; j < cpdim; ++j) { 2374 for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*Bfine[(ccell*numFIndices + i)*Nc + c]*qweights[ccell*Nc + c]*detJ; 2375 } 2376 /* Update interpolator */ 2377 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 2378 if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim); 2379 ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr); 2380 } 2381 ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 2382 } else { 2383 cpdim = 1; 2384 for (i = 0; i < numFIndices; ++i) { 2385 ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr); 2386 for (j = 0; j < cpdim; ++j) { 2387 for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*1.0*qweights[ccell*Nc + c]*detJ; 2388 } 2389 /* Update interpolator */ 2390 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 2391 ierr = PetscPrintf(PETSC_COMM_SELF, "Nq: %d %d Nf: %d %d Nc: %d %d\n", ccell, Nq, i, numFIndices, j, numCIndices);CHKERRQ(ierr); 2392 if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim); 2393 ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr); 2394 } 2395 } 2396 ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 2397 } 2398 ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 2399 ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 2400 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 2401 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 2402 } 2403 } 2404 ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 2405 ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr); 2406 ierr = PetscFree(elemMat);CHKERRQ(ierr); 2407 ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2408 ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2409 PetscFunctionReturn(0); 2410 } 2411 2412 /*@ 2413 DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns 2414 2415 Input Parameters: 2416 + dmc - The coarse mesh 2417 - dmf - The fine mesh 2418 - user - The user context 2419 2420 Output Parameter: 2421 . sc - The mapping 2422 2423 Level: developer 2424 2425 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM() 2426 @*/ 2427 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 2428 { 2429 PetscDS prob; 2430 PetscFE *feRef; 2431 PetscFV *fvRef; 2432 Vec fv, cv; 2433 IS fis, cis; 2434 PetscSection fsection, fglobalSection, csection, cglobalSection; 2435 PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 2436 PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m; 2437 PetscErrorCode ierr; 2438 2439 PetscFunctionBegin; 2440 ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 2441 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 2442 ierr = DMGetSection(dmf, &fsection);CHKERRQ(ierr); 2443 ierr = DMGetGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 2444 ierr = DMGetSection(dmc, &csection);CHKERRQ(ierr); 2445 ierr = DMGetGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 2446 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 2447 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 2448 ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2449 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 2450 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 2451 ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr); 2452 for (f = 0; f < Nf; ++f) { 2453 PetscObject obj; 2454 PetscClassId id; 2455 PetscInt fNb = 0, Nc = 0; 2456 2457 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 2458 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2459 if (id == PETSCFE_CLASSID) { 2460 PetscFE fe = (PetscFE) obj; 2461 2462 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 2463 ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 2464 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 2465 } else if (id == PETSCFV_CLASSID) { 2466 PetscFV fv = (PetscFV) obj; 2467 PetscDualSpace Q; 2468 2469 ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr); 2470 ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr); 2471 ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr); 2472 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 2473 } 2474 fTotDim += fNb; 2475 } 2476 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 2477 ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 2478 for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 2479 PetscFE feC; 2480 PetscFV fvC; 2481 PetscDualSpace QF, QC; 2482 PetscInt order = -1, NcF, NcC, fpdim, cpdim; 2483 2484 if (feRef[field]) { 2485 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 2486 ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 2487 ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 2488 ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 2489 ierr = PetscDualSpaceGetOrder(QF, &order);CHKERRQ(ierr); 2490 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 2491 ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 2492 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 2493 } else { 2494 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr); 2495 ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr); 2496 ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr); 2497 ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr); 2498 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 2499 ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr); 2500 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 2501 } 2502 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); 2503 for (c = 0; c < cpdim; ++c) { 2504 PetscQuadrature cfunc; 2505 const PetscReal *cqpoints, *cqweights; 2506 PetscInt NqcC, NpC; 2507 PetscBool found = PETSC_FALSE; 2508 2509 ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 2510 ierr = PetscQuadratureGetData(cfunc, NULL, &NqcC, &NpC, &cqpoints, &cqweights);CHKERRQ(ierr); 2511 if (NqcC != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %D must match number of field components", NqcC, NcC); 2512 if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 2513 for (f = 0; f < fpdim; ++f) { 2514 PetscQuadrature ffunc; 2515 const PetscReal *fqpoints, *fqweights; 2516 PetscReal sum = 0.0; 2517 PetscInt NqcF, NpF; 2518 2519 ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 2520 ierr = PetscQuadratureGetData(ffunc, NULL, &NqcF, &NpF, &fqpoints, &fqweights);CHKERRQ(ierr); 2521 if (NqcF != NcF) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %D must match number of field components", NqcF, NcF); 2522 if (NpC != NpF) continue; 2523 for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 2524 if (sum > 1.0e-9) continue; 2525 for (d = 0; d < NcC; ++d) sum += PetscAbsReal(cqweights[d]*fqweights[d]); 2526 if (sum < 1.0e-9) continue; 2527 cmap[offsetC+c] = offsetF+f; 2528 found = PETSC_TRUE; 2529 break; 2530 } 2531 if (!found) { 2532 /* TODO We really want the average here, but some asshole put VecScatter in the interface */ 2533 if (fvRef[field] || (feRef[field] && order == 0)) { 2534 cmap[offsetC+c] = offsetF+0; 2535 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection"); 2536 } 2537 } 2538 offsetC += cpdim; 2539 offsetF += fpdim; 2540 } 2541 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);} 2542 ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr); 2543 2544 ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 2545 ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 2546 ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr); 2547 ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 2548 ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 2549 ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 2550 ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 2551 for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 2552 for (c = cStart; c < cEnd; ++c) { 2553 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 2554 for (d = 0; d < cTotDim; ++d) { 2555 if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue; 2556 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]]); 2557 cindices[cellCIndices[d]-startC] = cellCIndices[d]; 2558 findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 2559 } 2560 } 2561 ierr = PetscFree(cmap);CHKERRQ(ierr); 2562 ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 2563 2564 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 2565 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 2566 ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 2567 ierr = ISDestroy(&cis);CHKERRQ(ierr); 2568 ierr = ISDestroy(&fis);CHKERRQ(ierr); 2569 ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 2570 ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 2571 ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 2572 PetscFunctionReturn(0); 2573 } 2574