1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petscfe.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DMPlexGetScale" 7 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 8 { 9 DM_Plex *mesh = (DM_Plex*) dm->data; 10 11 PetscFunctionBegin; 12 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13 PetscValidPointer(scale, 3); 14 *scale = mesh->scale[unit]; 15 PetscFunctionReturn(0); 16 } 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "DMPlexSetScale" 20 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 21 { 22 DM_Plex *mesh = (DM_Plex*) dm->data; 23 24 PetscFunctionBegin; 25 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 26 mesh->scale[unit] = scale; 27 PetscFunctionReturn(0); 28 } 29 30 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 31 { 32 switch (i) { 33 case 0: 34 switch (j) { 35 case 0: return 0; 36 case 1: 37 switch (k) { 38 case 0: return 0; 39 case 1: return 0; 40 case 2: return 1; 41 } 42 case 2: 43 switch (k) { 44 case 0: return 0; 45 case 1: return -1; 46 case 2: return 0; 47 } 48 } 49 case 1: 50 switch (j) { 51 case 0: 52 switch (k) { 53 case 0: return 0; 54 case 1: return 0; 55 case 2: return -1; 56 } 57 case 1: return 0; 58 case 2: 59 switch (k) { 60 case 0: return 1; 61 case 1: return 0; 62 case 2: return 0; 63 } 64 } 65 case 2: 66 switch (j) { 67 case 0: 68 switch (k) { 69 case 0: return 0; 70 case 1: return 1; 71 case 2: return 0; 72 } 73 case 1: 74 switch (k) { 75 case 0: return -1; 76 case 1: return 0; 77 case 2: return 0; 78 } 79 case 2: return 0; 80 } 81 } 82 return 0; 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "DMPlexCreateRigidBody" 87 /*@C 88 DMPlexCreateRigidBody - create rigid body modes from coordinates 89 90 Collective on DM 91 92 Input Arguments: 93 + dm - the DM 94 . section - the local section associated with the rigid field, or NULL for the default section 95 - globalSection - the global section associated with the rigid field, or NULL for the default section 96 97 Output Argument: 98 . sp - the null space 99 100 Note: This is necessary to take account of Dirichlet conditions on the displacements 101 102 Level: advanced 103 104 .seealso: MatNullSpaceCreate() 105 @*/ 106 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 107 { 108 MPI_Comm comm; 109 Vec coordinates, localMode, mode[6]; 110 PetscSection coordSection; 111 PetscScalar *coords; 112 PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 113 PetscErrorCode ierr; 114 115 PetscFunctionBegin; 116 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 117 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 118 if (dim == 1) { 119 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 120 PetscFunctionReturn(0); 121 } 122 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 123 if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 124 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 125 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 126 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 127 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 128 m = (dim*(dim+1))/2; 129 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 130 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 131 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 132 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 133 /* Assume P1 */ 134 ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 135 for (d = 0; d < dim; ++d) { 136 PetscScalar values[3] = {0.0, 0.0, 0.0}; 137 138 values[d] = 1.0; 139 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 140 for (v = vStart; v < vEnd; ++v) { 141 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 142 } 143 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 144 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145 } 146 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 147 for (d = dim; d < dim*(dim+1)/2; ++d) { 148 PetscInt i, j, k = dim > 2 ? d - dim : d; 149 150 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 151 for (v = vStart; v < vEnd; ++v) { 152 PetscScalar values[3] = {0.0, 0.0, 0.0}; 153 PetscInt off; 154 155 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 156 for (i = 0; i < dim; ++i) { 157 for (j = 0; j < dim; ++j) { 158 values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 159 } 160 } 161 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 162 } 163 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 164 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165 } 166 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 167 ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 168 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 169 /* Orthonormalize system */ 170 for (i = dim; i < m; ++i) { 171 PetscScalar dots[6]; 172 173 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 174 for (j = 0; j < i; ++j) dots[j] *= -1.0; 175 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 176 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 177 } 178 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 179 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 180 PetscFunctionReturn(0); 181 } 182 183 #undef __FUNCT__ 184 #define __FUNCT__ "DMPlexProjectFunctionLocal" 185 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 186 { 187 PetscDualSpace *sp; 188 PetscSection section; 189 PetscScalar *values; 190 PetscReal *v0, *J, detJ; 191 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 196 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 197 ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr); 198 for (f = 0; f < numFields; ++f) { 199 ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 200 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 201 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 202 totDim += spDim*numComp; 203 } 204 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 205 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 206 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 207 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 208 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 209 ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr); 210 for (c = cStart; c < cEnd; ++c) { 211 PetscCellGeometry geom; 212 213 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 214 geom.v0 = v0; 215 geom.J = J; 216 geom.detJ = &detJ; 217 for (f = 0, v = 0; f < numFields; ++f) { 218 void * const ctx = ctxs ? ctxs[f] : NULL; 219 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 220 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 221 for (d = 0; d < spDim; ++d) { 222 if (funcs[f]) { 223 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 224 } else { 225 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 226 } 227 v += numComp; 228 } 229 } 230 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 231 } 232 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 233 ierr = PetscFree2(v0,J);CHKERRQ(ierr); 234 ierr = PetscFree(sp);CHKERRQ(ierr); 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "DMPlexProjectFunction" 240 /*@C 241 DMPlexProjectFunction - This projects the given function into the function space provided. 242 243 Input Parameters: 244 + dm - The DM 245 . fe - The PetscFE associated with the field 246 . funcs - The coordinate functions to evaluate, one per field 247 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 248 - mode - The insertion mode for values 249 250 Output Parameter: 251 . X - vector 252 253 Level: developer 254 255 .seealso: DMPlexComputeL2Diff() 256 @*/ 257 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 258 { 259 Vec localX; 260 PetscErrorCode ierr; 261 262 PetscFunctionBegin; 263 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 264 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 265 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr); 266 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 267 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 268 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "DMPlexComputeL2Diff" 274 /*@C 275 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 276 277 Input Parameters: 278 + dm - The DM 279 . fe - The PetscFE object for each field 280 . funcs - The functions to evaluate for each field component 281 . ctxs - Optional array of contexts to pass to each function, or NULL. 282 - X - The coefficient vector u_h 283 284 Output Parameter: 285 . diff - The diff ||u - u_h||_2 286 287 Level: developer 288 289 .seealso: DMPlexProjectFunction() 290 @*/ 291 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 292 { 293 const PetscInt debug = 0; 294 PetscSection section; 295 PetscQuadrature quad; 296 Vec localX; 297 PetscScalar *funcVal; 298 PetscReal *coords, *v0, *J, *invJ, detJ; 299 PetscReal localDiff = 0.0; 300 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 301 PetscErrorCode ierr; 302 303 PetscFunctionBegin; 304 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 305 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 306 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 307 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 308 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 309 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 310 for (field = 0; field < numFields; ++field) { 311 PetscInt Nc; 312 313 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 314 numComponents += Nc; 315 } 316 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 317 ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 318 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 319 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 320 for (c = cStart; c < cEnd; ++c) { 321 PetscScalar *x = NULL; 322 PetscReal elemDiff = 0.0; 323 324 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 325 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 326 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 327 328 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 329 void * const ctx = ctxs ? ctxs[field] : NULL; 330 const PetscReal *quadPoints, *quadWeights; 331 PetscReal *basis; 332 PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 333 334 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 335 ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 336 ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 337 ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 338 if (debug) { 339 char title[1024]; 340 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 341 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 342 } 343 for (q = 0; q < numQuadPoints; ++q) { 344 for (d = 0; d < dim; d++) { 345 coords[d] = v0[d]; 346 for (e = 0; e < dim; e++) { 347 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 348 } 349 } 350 (*funcs[field])(coords, funcVal, ctx); 351 for (fc = 0; fc < numBasisComps; ++fc) { 352 PetscScalar interpolant = 0.0; 353 354 for (f = 0; f < numBasisFuncs; ++f) { 355 const PetscInt fidx = f*numBasisComps+fc; 356 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 357 } 358 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 359 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 360 } 361 } 362 comp += numBasisComps; 363 fieldOffset += numBasisFuncs*numBasisComps; 364 } 365 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 366 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 367 localDiff += elemDiff; 368 } 369 ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 370 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 371 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 372 *diff = PetscSqrtReal(*diff); 373 PetscFunctionReturn(0); 374 } 375 376 #undef __FUNCT__ 377 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 378 /*@C 379 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 380 381 Input Parameters: 382 + dm - The DM 383 . fe - The PetscFE object for each field 384 . funcs - The gradient functions to evaluate for each field component 385 . ctxs - Optional array of contexts to pass to each function, or NULL. 386 . X - The coefficient vector u_h 387 - n - The vector to project along 388 389 Output Parameter: 390 . diff - The diff ||(grad u - grad u_h) . n||_2 391 392 Level: developer 393 394 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 395 @*/ 396 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 397 { 398 const PetscInt debug = 0; 399 PetscSection section; 400 PetscQuadrature quad; 401 Vec localX; 402 PetscScalar *funcVal, *interpolantVec; 403 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 404 PetscReal localDiff = 0.0; 405 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 406 PetscErrorCode ierr; 407 408 PetscFunctionBegin; 409 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 410 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 411 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 412 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 413 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 414 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 415 for (field = 0; field < numFields; ++field) { 416 PetscInt Nc; 417 418 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 419 numComponents += Nc; 420 } 421 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 422 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 423 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 424 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 425 for (c = cStart; c < cEnd; ++c) { 426 PetscScalar *x = NULL; 427 PetscReal elemDiff = 0.0; 428 429 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 430 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 431 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 432 433 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 434 void * const ctx = ctxs ? ctxs[field] : NULL; 435 const PetscReal *quadPoints, *quadWeights; 436 PetscReal *basisDer; 437 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 438 439 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 440 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 441 ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr); 442 ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr); 443 if (debug) { 444 char title[1024]; 445 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 446 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 447 } 448 for (q = 0; q < numQuadPoints; ++q) { 449 for (d = 0; d < dim; d++) { 450 coords[d] = v0[d]; 451 for (e = 0; e < dim; e++) { 452 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 453 } 454 } 455 (*funcs[field])(coords, n, funcVal, ctx); 456 for (fc = 0; fc < Ncomp; ++fc) { 457 PetscScalar interpolant = 0.0; 458 459 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 460 for (f = 0; f < Nb; ++f) { 461 const PetscInt fidx = f*Ncomp+fc; 462 463 for (d = 0; d < dim; ++d) { 464 realSpaceDer[d] = 0.0; 465 for (g = 0; g < dim; ++g) { 466 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 467 } 468 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 469 } 470 } 471 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 472 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);} 473 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 474 } 475 } 476 comp += Ncomp; 477 fieldOffset += Nb*Ncomp; 478 } 479 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 480 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 481 localDiff += elemDiff; 482 } 483 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 484 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 485 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 486 *diff = PetscSqrtReal(*diff); 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "DMPlexComputeResidualFEM" 492 /*@ 493 DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 494 495 Input Parameters: 496 + dm - The mesh 497 . X - Local input vector 498 - user - The user context 499 500 Output Parameter: 501 . F - Local output vector 502 503 Note: 504 The first member of the user context must be an FEMContext. 505 506 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 507 like a GPU, or vectorize on a multicore machine. 508 509 Level: developer 510 511 .seealso: DMPlexComputeJacobianActionFEM() 512 @*/ 513 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 514 { 515 DM_Plex *mesh = (DM_Plex *) dm->data; 516 PetscFEM *fem = (PetscFEM *) user; 517 PetscFE *fe = fem->fe; 518 PetscFE *feAux = fem->feAux; 519 PetscFE *feBd = fem->feBd; 520 const char *name = "Residual"; 521 DM dmAux; 522 Vec A; 523 PetscQuadrature q; 524 PetscCellGeometry geom; 525 PetscSection section, sectionAux; 526 PetscReal *v0, *J, *invJ, *detJ; 527 PetscScalar *elemVec, *u, *a = NULL; 528 PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 529 PetscInt cellDof = 0, numComponents = 0; 530 PetscInt cellDofAux = 0, numComponentsAux = 0; 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 535 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 536 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 537 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 538 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 539 numCells = cEnd - cStart; 540 for (f = 0; f < Nf; ++f) { 541 PetscInt Nb, Nc; 542 543 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 544 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 545 cellDof += Nb*Nc; 546 numComponents += Nc; 547 } 548 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 549 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 550 if (dmAux) { 551 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 552 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 553 } 554 for (f = 0; f < NfAux; ++f) { 555 PetscInt Nb, Nc; 556 557 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 558 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 559 cellDofAux += Nb*Nc; 560 numComponentsAux += Nc; 561 } 562 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 563 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 564 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 565 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 566 for (c = cStart; c < cEnd; ++c) { 567 PetscScalar *x = NULL; 568 PetscInt i; 569 570 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 571 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 572 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 573 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 574 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 575 if (dmAux) { 576 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 577 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 578 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 579 } 580 } 581 for (f = 0; f < Nf; ++f) { 582 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 583 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 584 PetscInt numQuadPoints, Nb; 585 /* Conforming batches */ 586 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 587 /* Remainder */ 588 PetscInt Nr, offset; 589 590 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 591 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 592 ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 593 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 594 blockSize = Nb*numQuadPoints; 595 batchSize = numBlocks * blockSize; 596 ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 597 numChunks = numCells / (numBatches*batchSize); 598 Ne = numChunks*numBatches*batchSize; 599 Nr = numCells % (numBatches*batchSize); 600 offset = numCells - Nr; 601 geom.v0 = v0; 602 geom.J = J; 603 geom.invJ = invJ; 604 geom.detJ = detJ; 605 ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 606 geom.v0 = &v0[offset*dim]; 607 geom.J = &J[offset*dim*dim]; 608 geom.invJ = &invJ[offset*dim*dim]; 609 geom.detJ = &detJ[offset]; 610 ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 611 } 612 for (c = cStart; c < cEnd; ++c) { 613 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 614 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 615 } 616 ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 617 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 618 if (feBd) { 619 DMLabel label, depth; 620 IS pointIS; 621 const PetscInt *points; 622 PetscInt dep, numPoints, p, numFaces; 623 PetscReal *n; 624 625 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 626 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 627 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 628 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 629 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 630 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 631 PetscInt Nb, Nc; 632 633 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 634 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 635 cellDof += Nb*Nc; 636 numComponents += Nc; 637 } 638 for (p = 0, numFaces = 0; p < numPoints; ++p) { 639 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 640 if (dep == dim-1) ++numFaces; 641 } 642 ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr); 643 for (p = 0, f = 0; p < numPoints; ++p) { 644 const PetscInt point = points[p]; 645 PetscScalar *x = NULL; 646 PetscInt i; 647 648 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 649 if (dep != dim-1) continue; 650 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 651 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 652 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 653 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 654 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 655 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 656 ++f; 657 } 658 for (f = 0; f < Nf; ++f) { 659 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 660 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 661 PetscInt numQuadPoints, Nb; 662 /* Conforming batches */ 663 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 664 /* Remainder */ 665 PetscInt Nr, offset; 666 667 ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 668 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 669 ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 670 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 671 blockSize = Nb*numQuadPoints; 672 batchSize = numBlocks * blockSize; 673 ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 674 numChunks = numFaces / (numBatches*batchSize); 675 Ne = numChunks*numBatches*batchSize; 676 Nr = numFaces % (numBatches*batchSize); 677 offset = numFaces - Nr; 678 geom.v0 = v0; 679 geom.n = n; 680 geom.J = J; 681 geom.invJ = invJ; 682 geom.detJ = detJ; 683 ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 684 geom.v0 = &v0[offset*dim]; 685 geom.n = &n[offset*dim]; 686 geom.J = &J[offset*dim*dim]; 687 geom.invJ = &invJ[offset*dim*dim]; 688 geom.detJ = &detJ[offset]; 689 ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 690 } 691 for (p = 0, f = 0; p < numPoints; ++p) { 692 const PetscInt point = points[p]; 693 694 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 695 if (dep != dim-1) continue; 696 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 697 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 698 ++f; 699 } 700 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 701 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 702 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 703 } 704 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 705 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 711 /*@C 712 DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 713 714 Input Parameters: 715 + dm - The mesh 716 . J - The Jacobian shell matrix 717 . X - Local input vector 718 - user - The user context 719 720 Output Parameter: 721 . F - Local output vector 722 723 Note: 724 The first member of the user context must be an FEMContext. 725 726 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 727 like a GPU, or vectorize on a multicore machine. 728 729 Level: developer 730 731 .seealso: DMPlexComputeResidualFEM() 732 @*/ 733 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 734 { 735 DM_Plex *mesh = (DM_Plex *) dm->data; 736 PetscFEM *fem = (PetscFEM *) user; 737 PetscFE *fe = fem->fe; 738 PetscQuadrature quad; 739 PetscCellGeometry geom; 740 PetscSection section; 741 JacActionCtx *jctx; 742 PetscReal *v0, *J, *invJ, *detJ; 743 PetscScalar *elemVec, *u, *a; 744 PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 745 PetscInt cellDof = 0; 746 PetscErrorCode ierr; 747 748 PetscFunctionBegin; 749 /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 750 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 751 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 752 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 753 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 754 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 755 numCells = cEnd - cStart; 756 for (field = 0; field < numFields; ++field) { 757 PetscInt Nb, Nc; 758 759 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 760 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 761 cellDof += Nb*Nc; 762 } 763 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 764 ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&a,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 765 for (c = cStart; c < cEnd; ++c) { 766 PetscScalar *x = NULL; 767 PetscInt i; 768 769 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 770 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 771 ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 772 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 773 ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 774 ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 775 for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 776 ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 777 } 778 for (field = 0; field < numFields; ++field) { 779 PetscInt numQuadPoints, Nb; 780 /* Conforming batches */ 781 PetscInt numBlocks = 1; 782 PetscInt numBatches = 1; 783 PetscInt numChunks, Ne, blockSize, batchSize; 784 /* Remainder */ 785 PetscInt Nr, offset; 786 787 ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 788 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 789 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 790 blockSize = Nb*numQuadPoints; 791 batchSize = numBlocks * blockSize; 792 numChunks = numCells / (numBatches*batchSize); 793 Ne = numChunks*numBatches*batchSize; 794 Nr = numCells % (numBatches*batchSize); 795 offset = numCells - Nr; 796 geom.v0 = v0; 797 geom.J = J; 798 geom.invJ = invJ; 799 geom.detJ = detJ; 800 ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 801 geom.v0 = &v0[offset*dim]; 802 geom.J = &J[offset*dim*dim]; 803 geom.invJ = &invJ[offset*dim*dim]; 804 geom.detJ = &detJ[offset]; 805 ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 806 fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 807 } 808 for (c = cStart; c < cEnd; ++c) { 809 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 810 ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 811 } 812 ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 813 if (mesh->printFEM) { 814 PetscMPIInt rank, numProcs; 815 PetscInt p; 816 817 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 818 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 819 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 820 for (p = 0; p < numProcs; ++p) { 821 if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 822 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 823 } 824 } 825 /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 826 PetscFunctionReturn(0); 827 } 828 829 #undef __FUNCT__ 830 #define __FUNCT__ "DMPlexComputeJacobianFEM" 831 /*@ 832 DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 833 834 Input Parameters: 835 + dm - The mesh 836 . X - Local input vector 837 - user - The user context 838 839 Output Parameter: 840 . Jac - Jacobian matrix 841 842 Note: 843 The first member of the user context must be an FEMContext. 844 845 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 846 like a GPU, or vectorize on a multicore machine. 847 848 Level: developer 849 850 .seealso: FormFunctionLocal() 851 @*/ 852 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 853 { 854 DM_Plex *mesh = (DM_Plex *) dm->data; 855 PetscFEM *fem = (PetscFEM *) user; 856 PetscFE *fe = fem->fe; 857 PetscFE *feAux = fem->feAux; 858 const char *name = "Jacobian"; 859 DM dmAux; 860 Vec A; 861 PetscQuadrature quad; 862 PetscCellGeometry geom; 863 PetscSection section, globalSection, sectionAux; 864 PetscReal *v0, *J, *invJ, *detJ; 865 PetscScalar *elemMat, *u, *a; 866 PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 867 PetscInt cellDof = 0, numComponents = 0; 868 PetscInt cellDofAux = 0, numComponentsAux = 0; 869 PetscBool isShell; 870 PetscErrorCode ierr; 871 872 PetscFunctionBegin; 873 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 874 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 875 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 876 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 877 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 878 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 879 numCells = cEnd - cStart; 880 for (f = 0; f < Nf; ++f) { 881 PetscInt Nb, Nc; 882 883 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 884 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 885 cellDof += Nb*Nc; 886 numComponents += Nc; 887 } 888 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 889 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 890 if (dmAux) { 891 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 892 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 893 } 894 for (f = 0; f < NfAux; ++f) { 895 PetscInt Nb, Nc; 896 897 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 898 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 899 cellDofAux += Nb*Nc; 900 numComponentsAux += Nc; 901 } 902 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 903 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 904 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 905 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 906 for (c = cStart; c < cEnd; ++c) { 907 PetscScalar *x = NULL; 908 PetscInt i; 909 910 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 911 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 912 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 913 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 914 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 915 if (dmAux) { 916 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 917 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 918 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 919 } 920 } 921 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 922 for (fieldI = 0; fieldI < Nf; ++fieldI) { 923 PetscInt numQuadPoints, Nb; 924 /* Conforming batches */ 925 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 926 /* Remainder */ 927 PetscInt Nr, offset; 928 929 ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 930 ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 931 ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 932 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 933 blockSize = Nb*numQuadPoints; 934 batchSize = numBlocks * blockSize; 935 ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 936 numChunks = numCells / (numBatches*batchSize); 937 Ne = numChunks*numBatches*batchSize; 938 Nr = numCells % (numBatches*batchSize); 939 offset = numCells - Nr; 940 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 941 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 942 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 943 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 944 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 945 946 geom.v0 = v0; 947 geom.J = J; 948 geom.invJ = invJ; 949 geom.detJ = detJ; 950 ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 951 geom.v0 = &v0[offset*dim]; 952 geom.J = &J[offset*dim*dim]; 953 geom.invJ = &invJ[offset*dim*dim]; 954 geom.detJ = &detJ[offset]; 955 ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, Nf, fe, fieldI, fieldJ, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 956 } 957 } 958 for (c = cStart; c < cEnd; ++c) { 959 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 960 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 961 } 962 ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 963 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 964 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 965 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 966 if (mesh->printFEM) { 967 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 968 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 969 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 970 } 971 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 972 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 973 if (isShell) { 974 JacActionCtx *jctx; 975 976 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 977 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 978 } 979 PetscFunctionReturn(0); 980 } 981 982 #if 0 983 984 static void g0_identity_1d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 985 { 986 g3[0] = 1.0; 987 } 988 989 static void g0_identity_2d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 990 { 991 g3[0] = g3[3] = 1.0; 992 } 993 994 static void g0_identity_3d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 995 { 996 g3[0] = g3[4] = g3[8] = 1.0; 997 } 998 999 #undef __FUNCT__ 1000 #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken" 1001 PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user) 1002 { 1003 DM_Plex *mesh = (DM_Plex *) dmc->data; 1004 PetscFEM *fem = (PetscFEM *) user; 1005 PetscFE *fe = fem->fe; 1006 const char *name = "Interpolator"; 1007 PetscFE *feRef; 1008 PetscQuadrature quad, quadOld; 1009 PetscCellGeometry geom; 1010 PetscSection fsection, fglobalSection; 1011 PetscSection csection, cglobalSection; 1012 PetscReal *v0, *J, *invJ, *detJ; 1013 PetscScalar *elemMat; 1014 PetscInt dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 1015 PetscInt rCellDof = 0, cCellDof = 0, numComponents = 0; 1016 PetscErrorCode ierr; 1017 1018 PetscFunctionBegin; 1019 #if 0 1020 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1021 #endif 1022 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1023 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1024 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1025 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1026 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1027 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1028 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1029 numCells = cEnd - cStart; 1030 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1031 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1032 ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr); 1033 } 1034 for (f = 0; f < Nf; ++f) { 1035 PetscInt rNb, cNb, Nc; 1036 1037 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1038 ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr); 1039 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1040 numComponents += Nc; 1041 rCellDof += rNb*Nc; 1042 cCellDof += cNb*Nc; 1043 } 1044 ierr = MatZeroEntries(I);CHKERRQ(ierr); 1045 ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr); 1046 for (c = cStart; c < cEnd; ++c) { 1047 ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1048 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1049 } 1050 ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1051 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1052 PetscInt numQuadPoints, Nb; 1053 /* Conforming batches */ 1054 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1055 /* Remainder */ 1056 PetscInt Nr, offset; 1057 1058 /* Make new fine FE which refines the ref cell and the quadrature rule */ 1059 ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr); 1060 ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr); 1061 ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1062 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1063 blockSize = Nb*numQuadPoints; 1064 batchSize = numBlocks * blockSize; 1065 ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1066 numChunks = numCells / (numBatches*batchSize); 1067 Ne = numChunks*numBatches*batchSize; 1068 Nr = numCells % (numBatches*batchSize); 1069 offset = numCells - Nr; 1070 1071 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1072 /* void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */ 1073 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static; 1074 1075 /* Replace quadrature in coarse FE with refined quadrature */ 1076 ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr); 1077 ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr); 1078 ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr); 1079 geom.v0 = v0; 1080 geom.J = J; 1081 geom.invJ = invJ; 1082 geom.detJ = detJ; 1083 ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr); 1084 geom.v0 = &v0[offset*dim]; 1085 geom.J = &J[offset*dim*dim]; 1086 geom.invJ = &invJ[offset*dim*dim]; 1087 geom.detJ = &detJ[offset]; 1088 ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr); 1089 ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr); 1090 } 1091 } 1092 for (c = cStart; c < cEnd; ++c) { 1093 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);} 1094 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr); 1095 } 1096 ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1097 ierr = PetscFree(feRef);CHKERRQ(ierr); 1098 ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1099 ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1100 if (mesh->printFEM) { 1101 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1102 ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr); 1103 ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1104 } 1105 #if 0 1106 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1107 #endif 1108 PetscFunctionReturn(0); 1109 } 1110 #endif 1111 1112 #undef __FUNCT__ 1113 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1114 /*@ 1115 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1116 1117 Input Parameters: 1118 + dmf - The fine mesh 1119 . dmc - The coarse mesh 1120 - user - The user context 1121 1122 Output Parameter: 1123 . In - The interpolation matrix 1124 1125 Note: 1126 The first member of the user context must be an FEMContext. 1127 1128 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1129 like a GPU, or vectorize on a multicore machine. 1130 1131 Level: developer 1132 1133 .seealso: DMPlexComputeJacobianFEM() 1134 @*/ 1135 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1136 { 1137 DM_Plex *mesh = (DM_Plex *) dmc->data; 1138 PetscFEM *fem = (PetscFEM *) user; 1139 PetscFE *fe = fem->fe; 1140 const char *name = "Interpolator"; 1141 PetscFE *feRef; 1142 PetscSection fsection, fglobalSection; 1143 PetscSection csection, cglobalSection; 1144 PetscScalar *elemMat; 1145 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 1146 PetscInt rCellDof = 0, cCellDof = 0; 1147 PetscErrorCode ierr; 1148 1149 PetscFunctionBegin; 1150 #if 0 1151 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1152 #endif 1153 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1154 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1155 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1156 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1157 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1158 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1159 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1160 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1161 for (f = 0; f < Nf; ++f) { 1162 PetscInt rNb, cNb, Nc; 1163 1164 ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr); 1165 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1166 ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr); 1167 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1168 rCellDof += rNb*Nc; 1169 cCellDof += cNb*Nc; 1170 } 1171 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1172 ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr); 1173 ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1174 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1175 PetscDualSpace Qref; 1176 PetscQuadrature f; 1177 const PetscReal *qpoints, *qweights; 1178 PetscReal *points; 1179 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1180 1181 /* Compose points from all dual basis functionals */ 1182 ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr); 1183 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1184 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1185 for (i = 0; i < fpdim; ++i) { 1186 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1187 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1188 npoints += Np; 1189 } 1190 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1191 for (i = 0, k = 0; i < fpdim; ++i) { 1192 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1193 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1194 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1195 } 1196 1197 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1198 PetscSpace P; 1199 PetscReal *B; 1200 PetscInt NcJ, cpdim, j; 1201 1202 /* For now, fields only interpolate themselves */ 1203 if (fieldI != fieldJ) {offsetJ += cpdim; continue;} 1204 /* Evaluate basis at points */ 1205 ierr = PetscFEGetBasisSpace(fe[fieldJ], &P);CHKERRQ(ierr); 1206 ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr); 1207 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); 1208 ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr); 1209 ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1210 for (i = 0, k = 0; i < fpdim; ++i) { 1211 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1212 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1213 for (p = 0; p < Np; ++p, ++k) { 1214 for (j = 0; j < cpdim; ++j) { 1215 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cCellDof + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p]; 1216 } 1217 } 1218 } 1219 ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1220 offsetJ += cpdim*NcJ; 1221 } 1222 offsetI += fpdim*Nc; 1223 ierr = PetscFree(points);CHKERRQ(ierr); 1224 } 1225 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);} 1226 for (c = cStart; c < cEnd; ++c) { 1227 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1228 } 1229 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1230 ierr = PetscFree(feRef);CHKERRQ(ierr); 1231 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1232 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1233 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1234 if (mesh->printFEM) { 1235 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1236 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1237 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1238 } 1239 #if 0 1240 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1241 #endif 1242 PetscFunctionReturn(0); 1243 } 1244 1245 #undef __FUNCT__ 1246 #define __FUNCT__ "DMPlexAddBoundary" 1247 /* The ids can be overridden by the command line option -bc_<boundary name> */ 1248 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx) 1249 { 1250 DM_Plex *mesh = (DM_Plex *) dm->data; 1251 DMBoundary b; 1252 PetscErrorCode ierr; 1253 1254 PetscFunctionBegin; 1255 ierr = PetscNew(&b);CHKERRQ(ierr); 1256 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 1257 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 1258 ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr); 1259 b->essential = isEssential; 1260 b->field = field; 1261 b->func = bcFunc; 1262 b->numids = numids; 1263 b->ctx = ctx; 1264 b->next = mesh->boundary; 1265 mesh->boundary = b; 1266 PetscFunctionReturn(0); 1267 } 1268 1269 #undef __FUNCT__ 1270 #define __FUNCT__ "DMPlexGetNumBoundary" 1271 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd) 1272 { 1273 DM_Plex *mesh = (DM_Plex *) dm->data; 1274 DMBoundary b = mesh->boundary; 1275 1276 PetscFunctionBegin; 1277 *numBd = 0; 1278 while (b) {++(*numBd); b = b->next;} 1279 PetscFunctionReturn(0); 1280 } 1281 1282 #undef __FUNCT__ 1283 #define __FUNCT__ "DMPlexGetBoundary" 1284 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx) 1285 { 1286 DM_Plex *mesh = (DM_Plex *) dm->data; 1287 DMBoundary b = mesh->boundary; 1288 PetscInt n = 0; 1289 1290 PetscFunctionBegin; 1291 while (b) { 1292 if (n == bd) break; 1293 b = b->next; 1294 ++n; 1295 } 1296 if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 1297 if (isEssential) { 1298 PetscValidPointer(isEssential, 3); 1299 *isEssential = b->essential; 1300 } 1301 if (name) { 1302 PetscValidPointer(name, 4); 1303 *name = b->name; 1304 } 1305 if (field) { 1306 PetscValidPointer(field, 5); 1307 *field = b->field; 1308 } 1309 if (func) { 1310 PetscValidPointer(func, 6); 1311 *func = b->func; 1312 } 1313 if (numids) { 1314 PetscValidPointer(numids, 7); 1315 *numids = b->numids; 1316 } 1317 if (ids) { 1318 PetscValidPointer(ids, 8); 1319 *ids = b->ids; 1320 } 1321 if (ctx) { 1322 PetscValidPointer(ctx, 9); 1323 *ctx = b->ctx; 1324 } 1325 PetscFunctionReturn(0); 1326 } 1327