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; 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 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 223 v += numComp; 224 } 225 } 226 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 227 } 228 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 229 ierr = PetscFree2(v0,J);CHKERRQ(ierr); 230 ierr = PetscFree(sp);CHKERRQ(ierr); 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "DMPlexProjectFunction" 236 /*@C 237 DMPlexProjectFunction - This projects the given function into the function space provided. 238 239 Input Parameters: 240 + dm - The DM 241 . fe - The PetscFE associated with the field 242 . funcs - The coordinate functions to evaluate, one per field 243 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 244 - mode - The insertion mode for values 245 246 Output Parameter: 247 . X - vector 248 249 Level: developer 250 251 .seealso: DMPlexComputeL2Diff() 252 @*/ 253 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 254 { 255 Vec localX; 256 PetscErrorCode ierr; 257 258 PetscFunctionBegin; 259 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 260 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 261 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr); 262 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 263 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 264 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 265 PetscFunctionReturn(0); 266 } 267 268 #undef __FUNCT__ 269 #define __FUNCT__ "DMPlexComputeL2Diff" 270 /*@C 271 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 272 273 Input Parameters: 274 + dm - The DM 275 . fe - The PetscFE object for each field 276 . funcs - The functions to evaluate for each field component 277 . ctxs - Optional array of contexts to pass to each function, or NULL. 278 - X - The coefficient vector u_h 279 280 Output Parameter: 281 . diff - The diff ||u - u_h||_2 282 283 Level: developer 284 285 .seealso: DMPlexProjectFunction() 286 @*/ 287 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 288 { 289 const PetscInt debug = 0; 290 PetscSection section; 291 PetscQuadrature quad; 292 Vec localX; 293 PetscScalar *funcVal; 294 PetscReal *coords, *v0, *J, *invJ, detJ; 295 PetscReal localDiff = 0.0; 296 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 301 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 302 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 303 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 304 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 305 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 306 for (field = 0; field < numFields; ++field) { 307 PetscInt Nc; 308 309 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 310 numComponents += Nc; 311 } 312 ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 313 ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 314 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 315 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 316 for (c = cStart; c < cEnd; ++c) { 317 PetscScalar *x = NULL; 318 PetscReal elemDiff = 0.0; 319 320 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 321 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 322 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 323 324 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 325 void * const ctx = ctxs ? ctxs[field] : NULL; 326 const PetscReal *quadPoints, *quadWeights; 327 PetscReal *basis; 328 PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 329 330 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 331 ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 332 ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 333 ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 334 if (debug) { 335 char title[1024]; 336 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 337 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 338 } 339 for (q = 0; q < numQuadPoints; ++q) { 340 for (d = 0; d < dim; d++) { 341 coords[d] = v0[d]; 342 for (e = 0; e < dim; e++) { 343 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 344 } 345 } 346 (*funcs[field])(coords, funcVal, ctx); 347 for (fc = 0; fc < numBasisComps; ++fc) { 348 PetscScalar interpolant = 0.0; 349 350 for (f = 0; f < numBasisFuncs; ++f) { 351 const PetscInt fidx = f*numBasisComps+fc; 352 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 353 } 354 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);} 355 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 356 } 357 } 358 comp += numBasisComps; 359 fieldOffset += numBasisFuncs*numBasisComps; 360 } 361 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 362 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 363 localDiff += elemDiff; 364 } 365 ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 366 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 367 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 368 *diff = PetscSqrtReal(*diff); 369 PetscFunctionReturn(0); 370 } 371 372 #undef __FUNCT__ 373 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 374 /*@C 375 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 376 377 Input Parameters: 378 + dm - The DM 379 . fe - The PetscFE object for each field 380 . funcs - The gradient functions to evaluate for each field component 381 . ctxs - Optional array of contexts to pass to each function, or NULL. 382 . X - The coefficient vector u_h 383 - n - The vector to project along 384 385 Output Parameter: 386 . diff - The diff ||(grad u - grad u_h) . n||_2 387 388 Level: developer 389 390 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 391 @*/ 392 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 393 { 394 const PetscInt debug = 0; 395 PetscSection section; 396 PetscQuadrature quad; 397 Vec localX; 398 PetscScalar *funcVal, *interpolantVec; 399 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 400 PetscReal localDiff = 0.0; 401 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 402 PetscErrorCode ierr; 403 404 PetscFunctionBegin; 405 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 406 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 407 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 408 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 409 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 410 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 411 for (field = 0; field < numFields; ++field) { 412 PetscInt Nc; 413 414 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 415 numComponents += Nc; 416 } 417 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 418 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 419 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 420 ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 421 for (c = cStart; c < cEnd; ++c) { 422 PetscScalar *x = NULL; 423 PetscReal elemDiff = 0.0; 424 425 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 426 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 427 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 428 429 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 430 void * const ctx = ctxs ? ctxs[field] : NULL; 431 const PetscReal *quadPoints, *quadWeights; 432 PetscReal *basisDer; 433 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 434 435 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 436 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 437 ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr); 438 ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr); 439 if (debug) { 440 char title[1024]; 441 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 442 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 443 } 444 for (q = 0; q < numQuadPoints; ++q) { 445 for (d = 0; d < dim; d++) { 446 coords[d] = v0[d]; 447 for (e = 0; e < dim; e++) { 448 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 449 } 450 } 451 (*funcs[field])(coords, n, funcVal, ctx); 452 for (fc = 0; fc < Ncomp; ++fc) { 453 PetscScalar interpolant = 0.0; 454 455 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 456 for (f = 0; f < Nb; ++f) { 457 const PetscInt fidx = f*Ncomp+fc; 458 459 for (d = 0; d < dim; ++d) { 460 realSpaceDer[d] = 0.0; 461 for (g = 0; g < dim; ++g) { 462 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 463 } 464 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 465 } 466 } 467 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 468 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);} 469 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 470 } 471 } 472 comp += Ncomp; 473 fieldOffset += Nb*Ncomp; 474 } 475 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 476 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 477 localDiff += elemDiff; 478 } 479 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 480 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 481 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 482 *diff = PetscSqrtReal(*diff); 483 PetscFunctionReturn(0); 484 } 485 486 #undef __FUNCT__ 487 #define __FUNCT__ "DMPlexComputeResidualFEM" 488 /*@ 489 DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 490 491 Input Parameters: 492 + dm - The mesh 493 . X - Local input vector 494 - user - The user context 495 496 Output Parameter: 497 . F - Local output vector 498 499 Note: 500 The first member of the user context must be an FEMContext. 501 502 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 503 like a GPU, or vectorize on a multicore machine. 504 505 Level: developer 506 507 .seealso: DMPlexComputeJacobianActionFEM() 508 @*/ 509 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 510 { 511 DM_Plex *mesh = (DM_Plex *) dm->data; 512 PetscFEM *fem = (PetscFEM *) user; 513 PetscFE *fe = fem->fe; 514 PetscFE *feAux = fem->feAux; 515 PetscFE *feBd = fem->feBd; 516 const char *name = "Residual"; 517 DM dmAux; 518 Vec A; 519 PetscQuadrature q; 520 PetscCellGeometry geom; 521 PetscSection section, sectionAux; 522 PetscReal *v0, *J, *invJ, *detJ; 523 PetscScalar *elemVec, *u, *a = NULL; 524 PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 525 PetscInt cellDof = 0, numComponents = 0; 526 PetscInt cellDofAux = 0, numComponentsAux = 0; 527 PetscErrorCode ierr; 528 529 PetscFunctionBegin; 530 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 531 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 532 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 533 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 534 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 535 numCells = cEnd - cStart; 536 for (f = 0; f < Nf; ++f) { 537 PetscInt Nb, Nc; 538 539 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 540 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 541 cellDof += Nb*Nc; 542 numComponents += Nc; 543 } 544 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 545 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 546 if (dmAux) { 547 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 548 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 549 } 550 for (f = 0; f < NfAux; ++f) { 551 PetscInt Nb, Nc; 552 553 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 554 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 555 cellDofAux += Nb*Nc; 556 numComponentsAux += Nc; 557 } 558 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 559 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 560 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 561 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 562 for (c = cStart; c < cEnd; ++c) { 563 PetscScalar *x = NULL; 564 PetscInt i; 565 566 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 567 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 568 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 569 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 570 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 571 if (dmAux) { 572 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 573 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 574 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 575 } 576 } 577 for (f = 0; f < Nf; ++f) { 578 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 579 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 580 PetscInt numQuadPoints, Nb; 581 /* Conforming batches */ 582 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 583 /* Remainder */ 584 PetscInt Nr, offset; 585 586 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 587 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 588 ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 589 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 590 blockSize = Nb*numQuadPoints; 591 batchSize = numBlocks * blockSize; 592 ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 593 numChunks = numCells / (numBatches*batchSize); 594 Ne = numChunks*numBatches*batchSize; 595 Nr = numCells % (numBatches*batchSize); 596 offset = numCells - Nr; 597 geom.v0 = v0; 598 geom.J = J; 599 geom.invJ = invJ; 600 geom.detJ = detJ; 601 ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 602 geom.v0 = &v0[offset*dim]; 603 geom.J = &J[offset*dim*dim]; 604 geom.invJ = &invJ[offset*dim*dim]; 605 geom.detJ = &detJ[offset]; 606 ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 607 } 608 for (c = cStart; c < cEnd; ++c) { 609 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 610 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 611 } 612 ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 613 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 614 if (feBd) { 615 DMLabel label, depth; 616 IS pointIS; 617 const PetscInt *points; 618 PetscInt dep, numPoints, p, numFaces; 619 PetscReal *n; 620 621 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 622 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 623 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 624 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 625 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 626 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 627 PetscInt Nb, Nc; 628 629 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 630 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 631 cellDof += Nb*Nc; 632 numComponents += Nc; 633 } 634 for (p = 0, numFaces = 0; p < numPoints; ++p) { 635 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 636 if (dep == dim-1) ++numFaces; 637 } 638 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); 639 for (p = 0, f = 0; p < numPoints; ++p) { 640 const PetscInt point = points[p]; 641 PetscScalar *x = NULL; 642 PetscInt i; 643 644 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 645 if (dep != dim-1) continue; 646 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 647 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 648 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 649 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 650 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 651 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 652 ++f; 653 } 654 for (f = 0; f < Nf; ++f) { 655 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 656 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 657 PetscInt numQuadPoints, Nb; 658 /* Conforming batches */ 659 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 660 /* Remainder */ 661 PetscInt Nr, offset; 662 663 ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 664 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 665 ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 666 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 667 blockSize = Nb*numQuadPoints; 668 batchSize = numBlocks * blockSize; 669 ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 670 numChunks = numFaces / (numBatches*batchSize); 671 Ne = numChunks*numBatches*batchSize; 672 Nr = numFaces % (numBatches*batchSize); 673 offset = numFaces - Nr; 674 geom.v0 = v0; 675 geom.n = n; 676 geom.J = J; 677 geom.invJ = invJ; 678 geom.detJ = detJ; 679 ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 680 geom.v0 = &v0[offset*dim]; 681 geom.n = &n[offset*dim]; 682 geom.J = &J[offset*dim*dim]; 683 geom.invJ = &invJ[offset*dim*dim]; 684 geom.detJ = &detJ[offset]; 685 ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 686 } 687 for (p = 0, f = 0; p < numPoints; ++p) { 688 const PetscInt point = points[p]; 689 690 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 691 if (dep != dim-1) continue; 692 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 693 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 694 ++f; 695 } 696 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 697 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 698 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 699 } 700 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 701 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 #undef __FUNCT__ 706 #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 707 /*@C 708 DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 709 710 Input Parameters: 711 + dm - The mesh 712 . J - The Jacobian shell matrix 713 . X - Local input vector 714 - user - The user context 715 716 Output Parameter: 717 . F - Local output vector 718 719 Note: 720 The first member of the user context must be an FEMContext. 721 722 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 723 like a GPU, or vectorize on a multicore machine. 724 725 Level: developer 726 727 .seealso: DMPlexComputeResidualFEM() 728 @*/ 729 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 730 { 731 DM_Plex *mesh = (DM_Plex *) dm->data; 732 PetscFEM *fem = (PetscFEM *) user; 733 PetscFE *fe = fem->fe; 734 PetscQuadrature quad; 735 PetscCellGeometry geom; 736 PetscSection section; 737 JacActionCtx *jctx; 738 PetscReal *v0, *J, *invJ, *detJ; 739 PetscScalar *elemVec, *u, *a; 740 PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 741 PetscInt cellDof = 0; 742 PetscErrorCode ierr; 743 744 PetscFunctionBegin; 745 /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 746 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 747 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 748 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 749 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 750 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 751 numCells = cEnd - cStart; 752 for (field = 0; field < numFields; ++field) { 753 PetscInt Nb, Nc; 754 755 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 756 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 757 cellDof += Nb*Nc; 758 } 759 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 760 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); 761 for (c = cStart; c < cEnd; ++c) { 762 PetscScalar *x = NULL; 763 PetscInt i; 764 765 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 766 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 767 ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 768 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 769 ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 770 ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 771 for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 772 ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 773 } 774 for (field = 0; field < numFields; ++field) { 775 PetscInt numQuadPoints, Nb; 776 /* Conforming batches */ 777 PetscInt numBlocks = 1; 778 PetscInt numBatches = 1; 779 PetscInt numChunks, Ne, blockSize, batchSize; 780 /* Remainder */ 781 PetscInt Nr, offset; 782 783 ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 784 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 785 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 786 blockSize = Nb*numQuadPoints; 787 batchSize = numBlocks * blockSize; 788 numChunks = numCells / (numBatches*batchSize); 789 Ne = numChunks*numBatches*batchSize; 790 Nr = numCells % (numBatches*batchSize); 791 offset = numCells - Nr; 792 geom.v0 = v0; 793 geom.J = J; 794 geom.invJ = invJ; 795 geom.detJ = detJ; 796 ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 797 geom.v0 = &v0[offset*dim]; 798 geom.J = &J[offset*dim*dim]; 799 geom.invJ = &invJ[offset*dim*dim]; 800 geom.detJ = &detJ[offset]; 801 ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 802 fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 803 } 804 for (c = cStart; c < cEnd; ++c) { 805 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 806 ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 807 } 808 ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 809 if (mesh->printFEM) { 810 PetscMPIInt rank, numProcs; 811 PetscInt p; 812 813 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 814 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 815 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 816 for (p = 0; p < numProcs; ++p) { 817 if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 818 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 819 } 820 } 821 /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 822 PetscFunctionReturn(0); 823 } 824 825 #undef __FUNCT__ 826 #define __FUNCT__ "DMPlexComputeJacobianFEM" 827 /*@ 828 DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 829 830 Input Parameters: 831 + dm - The mesh 832 . X - Local input vector 833 - user - The user context 834 835 Output Parameter: 836 . Jac - Jacobian matrix 837 838 Note: 839 The first member of the user context must be an FEMContext. 840 841 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 842 like a GPU, or vectorize on a multicore machine. 843 844 Level: developer 845 846 .seealso: FormFunctionLocal() 847 @*/ 848 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 849 { 850 DM_Plex *mesh = (DM_Plex *) dm->data; 851 PetscFEM *fem = (PetscFEM *) user; 852 PetscFE *fe = fem->fe; 853 PetscFE *feAux = fem->feAux; 854 const char *name = "Jacobian"; 855 DM dmAux; 856 Vec A; 857 PetscQuadrature quad; 858 PetscCellGeometry geom; 859 PetscSection section, globalSection, sectionAux; 860 PetscReal *v0, *J, *invJ, *detJ; 861 PetscScalar *elemMat, *u, *a; 862 PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 863 PetscInt cellDof = 0, numComponents = 0; 864 PetscInt cellDofAux = 0, numComponentsAux = 0; 865 PetscBool isShell; 866 PetscErrorCode ierr; 867 868 PetscFunctionBegin; 869 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 870 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 871 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 872 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 873 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 874 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 875 numCells = cEnd - cStart; 876 for (f = 0; f < Nf; ++f) { 877 PetscInt Nb, Nc; 878 879 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 880 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 881 cellDof += Nb*Nc; 882 numComponents += Nc; 883 } 884 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 885 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 886 if (dmAux) { 887 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 888 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 889 } 890 for (f = 0; f < NfAux; ++f) { 891 PetscInt Nb, Nc; 892 893 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 894 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 895 cellDofAux += Nb*Nc; 896 numComponentsAux += Nc; 897 } 898 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 899 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 900 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 901 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 902 for (c = cStart; c < cEnd; ++c) { 903 PetscScalar *x = NULL; 904 PetscInt i; 905 906 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 907 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 908 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 909 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 910 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 911 if (dmAux) { 912 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 913 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 914 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 915 } 916 } 917 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 918 for (fieldI = 0; fieldI < Nf; ++fieldI) { 919 PetscInt numQuadPoints, Nb; 920 /* Conforming batches */ 921 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 922 /* Remainder */ 923 PetscInt Nr, offset; 924 925 ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 926 ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 927 ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 928 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 929 blockSize = Nb*numQuadPoints; 930 batchSize = numBlocks * blockSize; 931 ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 932 numChunks = numCells / (numBatches*batchSize); 933 Ne = numChunks*numBatches*batchSize; 934 Nr = numCells % (numBatches*batchSize); 935 offset = numCells - Nr; 936 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 937 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 938 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 939 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 940 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 941 942 geom.v0 = v0; 943 geom.J = J; 944 geom.invJ = invJ; 945 geom.detJ = detJ; 946 ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 947 geom.v0 = &v0[offset*dim]; 948 geom.J = &J[offset*dim*dim]; 949 geom.invJ = &invJ[offset*dim*dim]; 950 geom.detJ = &detJ[offset]; 951 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); 952 } 953 } 954 for (c = cStart; c < cEnd; ++c) { 955 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 956 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 957 } 958 ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 959 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 960 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 961 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 962 if (mesh->printFEM) { 963 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 964 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 965 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 966 } 967 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 968 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 969 if (isShell) { 970 JacActionCtx *jctx; 971 972 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 973 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 974 } 975 *str = SAME_NONZERO_PATTERN; 976 PetscFunctionReturn(0); 977 } 978 979 #undef __FUNCT__ 980 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 981 /*@ 982 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 983 984 Input Parameters: 985 + dmf - The fine mesh 986 . dmc - The coarse mesh 987 - user - The user context 988 989 Output Parameter: 990 . I - The interpolation matrix 991 992 Note: 993 The first member of the user context must be an FEMContext. 994 995 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 996 like a GPU, or vectorize on a multicore machine. 997 998 Level: developer 999 1000 .seealso: DMPlexComputeJacobianFEM() 1001 @*/ 1002 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat I, void *user) 1003 { 1004 DM_Plex *mesh = (DM_Plex *) dmc->data; 1005 PetscFEM *fem = (PetscFEM *) user; 1006 PetscFE *fe = fem->fe; 1007 const char *name = "Interpolator"; 1008 PetscQuadrature quad, quadOld; 1009 PetscCellGeometry geom; 1010 PetscSection section, globalSection; 1011 PetscReal *v0, *J, *invJ, *detJ; 1012 PetscScalar *elemMat; 1013 PetscInt dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 1014 PetscInt cellDof = 0, numComponents = 0; 1015 PetscErrorCode ierr; 1016 1017 PetscFunctionBegin; 1018 #if 0 1019 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1020 #endif 1021 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1022 ierr = DMGetDefaultSection(dmf, §ion);CHKERRQ(ierr); 1023 ierr = DMGetDefaultGlobalSection(dmf, &globalSection);CHKERRQ(ierr); 1024 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1025 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1026 numCells = cEnd - cStart; 1027 for (f = 0; f < Nf; ++f) { 1028 PetscInt Nb, Nc; 1029 1030 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 1031 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1032 cellDof += Nb*Nc; 1033 numComponents += Nc; 1034 } 1035 ierr = MatZeroEntries(I);CHKERRQ(ierr); 1036 ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 1037 for (c = cStart; c < cEnd; ++c) { 1038 ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1039 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1040 } 1041 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1042 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1043 PetscFE feI; 1044 PetscInt numQuadPoints, Nb; 1045 /* Conforming batches */ 1046 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1047 /* Remainder */ 1048 PetscInt Nr, offset; 1049 1050 /* Make new fine FE which refines the ref cell and the quadrature rule */ 1051 ierr = PetscFERefine(fe[fieldI], &feI);CHKERRQ(ierr); 1052 ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 1053 ierr = PetscFEGetDimension(feI, &Nb);CHKERRQ(ierr); 1054 ierr = PetscFEGetTileSizes(feI, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1055 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1056 blockSize = Nb*numQuadPoints; 1057 batchSize = numBlocks * blockSize; 1058 ierr = PetscFESetTileSizes(feI, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1059 numChunks = numCells / (numBatches*batchSize); 1060 Ne = numChunks*numBatches*batchSize; 1061 Nr = numCells % (numBatches*batchSize); 1062 offset = numCells - Nr; 1063 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1064 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 1065 1066 /* Replace quadrature in coarse FE with refined quadrature */ 1067 ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr); 1068 ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr); 1069 ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr); 1070 geom.v0 = v0; 1071 geom.J = J; 1072 geom.invJ = invJ; 1073 geom.detJ = detJ; 1074 ierr = PetscFEIntegrateInterpolator_Basic(feI, Ne, Nf, fe, fieldI, fieldJ, geom, g0, elemMat);CHKERRQ(ierr); 1075 geom.v0 = &v0[offset*dim]; 1076 geom.J = &J[offset*dim*dim]; 1077 geom.invJ = &invJ[offset*dim*dim]; 1078 geom.detJ = &detJ[offset]; 1079 ierr = PetscFEIntegrateInterpolator_Basic(feI, Nr, Nf, fe, fieldI, fieldJ, geom, g0, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 1080 ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr); 1081 } 1082 } 1083 for (c = cStart; c < cEnd; ++c) { 1084 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 1085 ierr = DMPlexMatSetClosure(dmc, section, globalSection, I, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 1086 } 1087 ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1088 ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1089 ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1090 if (mesh->printFEM) { 1091 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1092 ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr); 1093 ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1094 } 1095 #if 0 1096 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1097 #endif 1098 PetscFunctionReturn(0); 1099 } 1100 1101 #undef __FUNCT__ 1102 #define __FUNCT__ "DMPlexAddBoundary" 1103 /* The ids can be overridden by the command line option -bc_<boundary name> */ 1104 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx) 1105 { 1106 DM_Plex *mesh = (DM_Plex *) dm->data; 1107 DMBoundary b; 1108 PetscErrorCode ierr; 1109 1110 PetscFunctionBegin; 1111 ierr = PetscNew(&b);CHKERRQ(ierr); 1112 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 1113 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 1114 ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr); 1115 b->essential = isEssential; 1116 b->field = field; 1117 b->func = bcFunc; 1118 b->numids = numids; 1119 b->ctx = ctx; 1120 b->next = mesh->boundary; 1121 mesh->boundary = b; 1122 PetscFunctionReturn(0); 1123 } 1124 1125 #undef __FUNCT__ 1126 #define __FUNCT__ "DMPlexGetNumBoundary" 1127 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd) 1128 { 1129 DM_Plex *mesh = (DM_Plex *) dm->data; 1130 DMBoundary b = mesh->boundary; 1131 1132 PetscFunctionBegin; 1133 *numBd = 0; 1134 while (b) {++(*numBd); b = b->next;} 1135 PetscFunctionReturn(0); 1136 } 1137 1138 #undef __FUNCT__ 1139 #define __FUNCT__ "DMPlexGetBoundary" 1140 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx) 1141 { 1142 DM_Plex *mesh = (DM_Plex *) dm->data; 1143 DMBoundary b = mesh->boundary; 1144 PetscInt n = 0; 1145 1146 PetscFunctionBegin; 1147 while (b) { 1148 if (n == bd) break; 1149 b = b->next; 1150 ++n; 1151 } 1152 if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 1153 if (isEssential) { 1154 PetscValidPointer(isEssential, 3); 1155 *isEssential = b->essential; 1156 } 1157 if (name) { 1158 PetscValidPointer(name, 4); 1159 *name = b->name; 1160 } 1161 if (field) { 1162 PetscValidPointer(field, 5); 1163 *field = b->field; 1164 } 1165 if (func) { 1166 PetscValidPointer(func, 6); 1167 *func = b->func; 1168 } 1169 if (numids) { 1170 PetscValidPointer(numids, 7); 1171 *numids = b->numids; 1172 } 1173 if (ids) { 1174 PetscValidPointer(ids, 8); 1175 *ids = b->ids; 1176 } 1177 if (ctx) { 1178 PetscValidPointer(ctx, 9); 1179 *ctx = b->ctx; 1180 } 1181 PetscFunctionReturn(0); 1182 } 1183