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