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