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. ctxs itself may be 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 PetscInt numQuadPoints = quad.numPoints; 327 const PetscReal *quadPoints = quad.points; 328 const PetscReal *quadWeights = quad.weights; 329 PetscReal *basis; 330 PetscInt numBasisFuncs, numBasisComps, q, d, e, fc, f; 331 332 ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 333 ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 334 ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 335 if (debug) { 336 char title[1024]; 337 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 338 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 339 } 340 for (q = 0; q < numQuadPoints; ++q) { 341 for (d = 0; d < dim; d++) { 342 coords[d] = v0[d]; 343 for (e = 0; e < dim; e++) { 344 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 345 } 346 } 347 (*funcs[field])(coords, funcVal, ctx); 348 for (fc = 0; fc < numBasisComps; ++fc) { 349 PetscScalar interpolant = 0.0; 350 351 for (f = 0; f < numBasisFuncs; ++f) { 352 const PetscInt fidx = f*numBasisComps+fc; 353 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 354 } 355 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);} 356 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 357 } 358 } 359 comp += numBasisComps; 360 fieldOffset += numBasisFuncs*numBasisComps; 361 } 362 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 363 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 364 localDiff += elemDiff; 365 } 366 ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 367 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 368 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 369 *diff = PetscSqrtReal(*diff); 370 PetscFunctionReturn(0); 371 } 372 373 #if 0 374 375 #undef __FUNCT__ 376 #define __FUNCT__ "DMPlexComputeResidualFEM" 377 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 378 { 379 DM_Plex *mesh = (DM_Plex *) dm->data; 380 PetscFEM *fem = (PetscFEM *) user; 381 PetscQuadrature *quad = fem->quad; 382 PetscQuadrature *quadBd = fem->quadBd; 383 PetscSection section; 384 PetscReal *v0, *n, *J, *invJ, *detJ; 385 PetscScalar *elemVec, *u; 386 PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 387 PetscInt cellDof, numComponents; 388 PetscBool has; 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 if (has && quadBd) { 393 DMLabel label; 394 IS pointIS; 395 const PetscInt *points; 396 PetscInt numPoints, p; 397 398 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 399 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 400 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 401 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 402 for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) { 403 cellDof += quadBd[field].numBasisFuncs*quadBd[field].numComponents; 404 numComponents += quadBd[field].numComponents; 405 } 406 ierr = PetscMalloc7(numPoints*cellDof,&u,numPoints*dim,&v0,numPoints*dim,&n,numPoints*dim*dim,&J,numPoints*dim*dim,&invJ,numPoints,&detJ,numPoints*cellDof,&elemVec);CHKERRQ(ierr); 407 for (p = 0; p < numPoints; ++p) { 408 const PetscInt point = points[p]; 409 PetscScalar *x; 410 PetscInt i; 411 412 /* TODO: Add normal determination here */ 413 ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr); 414 if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point); 415 ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 416 417 for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i]; 418 ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 419 } 420 for (field = 0; field < numFields; ++field) { 421 const PetscInt numQuadPoints = quadBd[field].numQuadPoints; 422 const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs; 423 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field]; 424 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field]; 425 /* Conforming batches */ 426 PetscInt blockSize = numBasisFuncs*numQuadPoints; 427 PetscInt numBlocks = 1; 428 PetscInt batchSize = numBlocks * blockSize; 429 PetscInt numBatches = numBatchesTmp; 430 PetscInt numChunks = numPoints / (numBatches*batchSize); 431 /* Remainder */ 432 PetscInt numRemainder = numPoints % (numBatches * batchSize); 433 PetscInt offset = numPoints - numRemainder; 434 435 ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr); 436 ierr = (*mesh->integrateBdResidualFEM)(numRemainder, numFields, field, quadBd, &u[offset*cellDof], &v0[offset*dim], &n[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset], 437 f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 438 } 439 for (p = 0; p < numPoints; ++p) { 440 const PetscInt point = points[p]; 441 442 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);} 443 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr); 444 } 445 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 446 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 447 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 448 } 449 PetscFunctionReturn(0); 450 } 451 452 #else 453 454 #undef __FUNCT__ 455 #define __FUNCT__ "DMPlexComputeResidualFEM" 456 /*@ 457 DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 458 459 Input Parameters: 460 + dm - The mesh 461 . X - Local input vector 462 - user - The user context 463 464 Output Parameter: 465 . F - Local output vector 466 467 Note: 468 The second member of the user context must be an FEMContext. 469 470 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 471 like a GPU, or vectorize on a multicore machine. 472 473 Level: developer 474 475 .seealso: DMPlexComputeJacobianActionFEM() 476 @*/ 477 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 478 { 479 DM_Plex *mesh = (DM_Plex *) dm->data; 480 PetscFEM *fem = (PetscFEM *) user; 481 PetscFE *fe = fem->fe; 482 PetscFE *feAux = fem->feAux; 483 PetscFE *feBd = fem->feBd; 484 const char *name = "Residual"; 485 DM dmAux; 486 Vec A; 487 PetscQuadrature q; 488 PetscCellGeometry geom; 489 PetscSection section, sectionAux; 490 PetscReal *v0, *J, *invJ, *detJ; 491 PetscScalar *elemVec, *u, *a; 492 PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 493 PetscInt cellDof = 0, numComponents = 0; 494 PetscInt cellDofAux = 0, numComponentsAux = 0; 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 499 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 500 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 501 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 502 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 503 numCells = cEnd - cStart; 504 for (f = 0; f < Nf; ++f) { 505 PetscInt Nb, Nc; 506 507 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 508 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 509 cellDof += Nb*Nc; 510 numComponents += Nc; 511 } 512 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 513 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 514 if (dmAux) { 515 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 516 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 517 } 518 for (f = 0; f < NfAux; ++f) { 519 PetscInt Nb, Nc; 520 521 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 522 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 523 cellDofAux += Nb*Nc; 524 numComponentsAux += Nc; 525 } 526 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 527 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 528 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 529 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 530 for (c = cStart; c < cEnd; ++c) { 531 PetscScalar *x = NULL; 532 PetscInt i; 533 534 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 535 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 536 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 537 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 538 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 539 if (dmAux) { 540 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 541 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 542 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 543 } 544 } 545 for (f = 0; f < Nf; ++f) { 546 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 547 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 548 PetscInt Nb; 549 /* Conforming batches */ 550 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 551 /* Remainder */ 552 PetscInt Nr, offset; 553 554 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 555 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 556 ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 557 blockSize = Nb*q.numPoints; 558 batchSize = numBlocks * blockSize; 559 ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 560 numChunks = numCells / (numBatches*batchSize); 561 Ne = numChunks*numBatches*batchSize; 562 Nr = numCells % (numBatches*batchSize); 563 offset = numCells - Nr; 564 geom.v0 = v0; 565 geom.J = J; 566 geom.invJ = invJ; 567 geom.detJ = detJ; 568 ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 569 geom.v0 = &v0[offset*dim]; 570 geom.J = &J[offset*dim*dim]; 571 geom.invJ = &invJ[offset*dim*dim]; 572 geom.detJ = &detJ[offset]; 573 ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 574 } 575 for (c = cStart; c < cEnd; ++c) { 576 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 577 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 578 } 579 ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 580 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 581 if (feBd) { 582 DMLabel label, depth; 583 IS pointIS; 584 const PetscInt *points; 585 PetscInt dep, numPoints, p, numFaces; 586 PetscReal *n; 587 588 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 589 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 590 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 591 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 592 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 593 for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 594 PetscInt Nb, Nc; 595 596 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 597 ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 598 cellDof += Nb*Nc; 599 numComponents += Nc; 600 } 601 for (p = 0, numFaces = 0; p < numPoints; ++p) { 602 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 603 if (dep == dim-1) ++numFaces; 604 } 605 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); 606 for (p = 0, f = 0; p < numPoints; ++p) { 607 const PetscInt point = points[p]; 608 PetscScalar *x = NULL; 609 PetscInt i; 610 611 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 612 if (dep != dim-1) continue; 613 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 614 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 615 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 616 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 617 for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 618 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 619 ++f; 620 } 621 for (f = 0; f < Nf; ++f) { 622 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 623 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 624 PetscInt Nb; 625 /* Conforming batches */ 626 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 627 /* Remainder */ 628 PetscInt Nr, offset; 629 630 ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 631 ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 632 ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 633 blockSize = Nb*q.numPoints; 634 batchSize = numBlocks * blockSize; 635 ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 636 numChunks = numFaces / (numBatches*batchSize); 637 Ne = numChunks*numBatches*batchSize; 638 Nr = numFaces % (numBatches*batchSize); 639 offset = numFaces - Nr; 640 geom.v0 = v0; 641 geom.n = n; 642 geom.J = J; 643 geom.invJ = invJ; 644 geom.detJ = detJ; 645 ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 646 geom.v0 = &v0[offset*dim]; 647 geom.n = &n[offset*dim]; 648 geom.J = &J[offset*dim*dim]; 649 geom.invJ = &invJ[offset*dim*dim]; 650 geom.detJ = &detJ[offset]; 651 ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 652 } 653 for (p = 0, f = 0; p < numPoints; ++p) { 654 const PetscInt point = points[p]; 655 656 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 657 if (dep != dim-1) continue; 658 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 659 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 660 ++f; 661 } 662 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 663 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 664 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 665 } 666 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 667 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 668 PetscFunctionReturn(0); 669 } 670 671 #endif 672 673 #undef __FUNCT__ 674 #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 675 /*@C 676 DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 677 678 Input Parameters: 679 + dm - The mesh 680 . J - The Jacobian shell matrix 681 . X - Local input vector 682 - user - The user context 683 684 Output Parameter: 685 . F - Local output vector 686 687 Note: 688 The second member of the user context must be an FEMContext. 689 690 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 691 like a GPU, or vectorize on a multicore machine. 692 693 Level: developer 694 695 .seealso: DMPlexComputeResidualFEM() 696 @*/ 697 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 698 { 699 DM_Plex *mesh = (DM_Plex *) dm->data; 700 PetscFEM *fem = (PetscFEM *) user; 701 PetscFE *fe = fem->fe; 702 PetscQuadrature quad; 703 PetscCellGeometry geom; 704 PetscSection section; 705 JacActionCtx *jctx; 706 PetscReal *v0, *J, *invJ, *detJ; 707 PetscScalar *elemVec, *u, *a; 708 PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 709 PetscInt cellDof = 0; 710 PetscErrorCode ierr; 711 712 PetscFunctionBegin; 713 /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 714 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 715 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 716 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 717 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 718 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 719 numCells = cEnd - cStart; 720 for (field = 0; field < numFields; ++field) { 721 PetscInt Nb, Nc; 722 723 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 724 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 725 cellDof += Nb*Nc; 726 } 727 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 728 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); 729 for (c = cStart; c < cEnd; ++c) { 730 PetscScalar *x = NULL; 731 PetscInt i; 732 733 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 734 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 735 ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 736 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 737 ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 738 ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 739 for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 740 ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 741 } 742 for (field = 0; field < numFields; ++field) { 743 PetscInt Nb; 744 /* Conforming batches */ 745 PetscInt numBlocks = 1; 746 PetscInt numBatches = 1; 747 PetscInt numChunks, Ne, blockSize, batchSize; 748 /* Remainder */ 749 PetscInt Nr, offset; 750 751 ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 752 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 753 blockSize = Nb*quad.numPoints; 754 batchSize = numBlocks * blockSize; 755 numChunks = numCells / (numBatches*batchSize); 756 Ne = numChunks*numBatches*batchSize; 757 Nr = numCells % (numBatches*batchSize); 758 offset = numCells - Nr; 759 geom.v0 = v0; 760 geom.J = J; 761 geom.invJ = invJ; 762 geom.detJ = detJ; 763 ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 764 geom.v0 = &v0[offset*dim]; 765 geom.J = &J[offset*dim*dim]; 766 geom.invJ = &invJ[offset*dim*dim]; 767 geom.detJ = &detJ[offset]; 768 ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 769 fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 770 } 771 for (c = cStart; c < cEnd; ++c) { 772 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 773 ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 774 } 775 ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 776 if (mesh->printFEM) { 777 PetscMPIInt rank, numProcs; 778 PetscInt p; 779 780 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 781 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 782 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 783 for (p = 0; p < numProcs; ++p) { 784 if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 785 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 786 } 787 } 788 /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 789 PetscFunctionReturn(0); 790 } 791 792 #undef __FUNCT__ 793 #define __FUNCT__ "DMPlexComputeJacobianFEM" 794 /*@ 795 DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 796 797 Input Parameters: 798 + dm - The mesh 799 . X - Local input vector 800 - user - The user context 801 802 Output Parameter: 803 . Jac - Jacobian matrix 804 805 Note: 806 The second member of the user context must be an FEMContext. 807 808 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 809 like a GPU, or vectorize on a multicore machine. 810 811 Level: developer 812 813 .seealso: FormFunctionLocal() 814 @*/ 815 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 816 { 817 DM_Plex *mesh = (DM_Plex *) dm->data; 818 PetscFEM *fem = (PetscFEM *) user; 819 PetscFE *fe = fem->fe; 820 PetscFE *feAux = fem->feAux; 821 const char *name = "Jacobian"; 822 DM dmAux; 823 Vec A; 824 PetscQuadrature quad; 825 PetscCellGeometry geom; 826 PetscSection section, globalSection, sectionAux; 827 PetscReal *v0, *J, *invJ, *detJ; 828 PetscScalar *elemMat, *u, *a; 829 PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 830 PetscInt cellDof = 0, numComponents = 0; 831 PetscInt cellDofAux = 0, numComponentsAux = 0; 832 PetscBool isShell; 833 PetscErrorCode ierr; 834 835 PetscFunctionBegin; 836 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 837 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 838 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 839 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 840 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 841 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 842 numCells = cEnd - cStart; 843 for (f = 0; f < Nf; ++f) { 844 PetscInt Nb, Nc; 845 846 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 847 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 848 cellDof += Nb*Nc; 849 numComponents += Nc; 850 } 851 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 852 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 853 if (dmAux) { 854 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 855 ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 856 } 857 for (f = 0; f < NfAux; ++f) { 858 PetscInt Nb, Nc; 859 860 ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 861 ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 862 cellDofAux += Nb*Nc; 863 numComponentsAux += Nc; 864 } 865 ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 866 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 867 ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 868 if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 869 for (c = cStart; c < cEnd; ++c) { 870 PetscScalar *x = NULL; 871 PetscInt i; 872 873 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 874 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 875 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 876 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 877 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 878 if (dmAux) { 879 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 880 for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 881 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 882 } 883 } 884 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 885 for (fieldI = 0; fieldI < Nf; ++fieldI) { 886 PetscInt Nb; 887 /* Conforming batches */ 888 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 889 /* Remainder */ 890 PetscInt Nr, offset; 891 892 ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 893 ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 894 ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 895 blockSize = Nb*quad.numPoints; 896 batchSize = numBlocks * blockSize; 897 ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 898 numChunks = numCells / (numBatches*batchSize); 899 Ne = numChunks*numBatches*batchSize; 900 Nr = numCells % (numBatches*batchSize); 901 offset = numCells - Nr; 902 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 903 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 904 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 905 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 906 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 907 908 geom.v0 = v0; 909 geom.J = J; 910 geom.invJ = invJ; 911 geom.detJ = detJ; 912 ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 913 geom.v0 = &v0[offset*dim]; 914 geom.J = &J[offset*dim*dim]; 915 geom.invJ = &invJ[offset*dim*dim]; 916 geom.detJ = &detJ[offset]; 917 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); 918 } 919 } 920 for (c = cStart; c < cEnd; ++c) { 921 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 922 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 923 } 924 ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 925 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 926 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 927 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 928 if (mesh->printFEM) { 929 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 930 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 931 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 932 } 933 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 934 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 935 if (isShell) { 936 JacActionCtx *jctx; 937 938 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 939 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 940 } 941 *str = SAME_NONZERO_PATTERN; 942 PetscFunctionReturn(0); 943 } 944