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 = DMPlexGetCoordinateSection(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, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX) 186 { 187 Vec coordinates; 188 PetscSection section, cSection; 189 PetscInt dim, vStart, vEnd, v, c, d; 190 PetscScalar *values, *cArray; 191 PetscReal *coords; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 196 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 197 ierr = DMPlexGetCoordinateSection(dm, &cSection);CHKERRQ(ierr); 198 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 199 ierr = PetscMalloc(numComp * sizeof(PetscScalar), &values);CHKERRQ(ierr); 200 ierr = VecGetArray(coordinates, &cArray);CHKERRQ(ierr); 201 ierr = PetscSectionGetDof(cSection, vStart, &dim);CHKERRQ(ierr); 202 ierr = PetscMalloc(dim * sizeof(PetscReal),&coords);CHKERRQ(ierr); 203 for (v = vStart; v < vEnd; ++v) { 204 PetscInt dof, off; 205 206 ierr = PetscSectionGetDof(cSection, v, &dof);CHKERRQ(ierr); 207 ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr); 208 if (dof > dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot have more coordinates %d then dimensions %d", dof, dim); 209 for (d = 0; d < dof; ++d) coords[d] = PetscRealPart(cArray[off+d]); 210 for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]); 211 ierr = VecSetValuesSection(localX, section, v, values, mode);CHKERRQ(ierr); 212 } 213 ierr = VecRestoreArray(coordinates, &cArray);CHKERRQ(ierr); 214 /* Temporary, must be replaced by a projection on the finite element basis */ 215 { 216 PetscInt eStart = 0, eEnd = 0, e, depth; 217 218 ierr = DMPlexGetLabelSize(dm, "depth", &depth);CHKERRQ(ierr); 219 --depth; 220 if (depth > 1) {ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);} 221 for (e = eStart; e < eEnd; ++e) { 222 const PetscInt *cone = NULL; 223 PetscInt coneSize, d; 224 PetscScalar *coordsA, *coordsB; 225 226 ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr); 227 ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr); 228 if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Cone size %d for point %d should be 2", coneSize, e); 229 ierr = VecGetValuesSection(coordinates, cSection, cone[0], &coordsA);CHKERRQ(ierr); 230 ierr = VecGetValuesSection(coordinates, cSection, cone[1], &coordsB);CHKERRQ(ierr); 231 for (d = 0; d < dim; ++d) { 232 coords[d] = 0.5*(PetscRealPart(coordsA[d]) + PetscRealPart(coordsB[d])); 233 } 234 for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]); 235 ierr = VecSetValuesSection(localX, section, e, values, mode);CHKERRQ(ierr); 236 } 237 } 238 239 ierr = PetscFree(coords);CHKERRQ(ierr); 240 ierr = PetscFree(values);CHKERRQ(ierr); 241 #if 0 242 const PetscInt localDof = this->_mesh->sizeWithBC(s, *cells->begin()); 243 PetscReal detJ; 244 245 ierr = PetscMalloc(localDof * sizeof(PetscScalar), &values);CHKERRQ(ierr); 246 ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr); 247 ALE::ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> pV(PetscPowInt(this->_mesh->getSieve()->getMaxConeSize(),dim+1), true); 248 249 for (PetscInt c = cStart; c < cEnd; ++c) { 250 ALE::ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*this->_mesh->getSieve(), c, pV); 251 const PETSC_MESH_TYPE::point_type *oPoints = pV.getPoints(); 252 const int oSize = pV.getSize(); 253 int v = 0; 254 255 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 256 for (PetscInt cl = 0; cl < oSize; ++cl) { 257 const PetscInt fDim; 258 259 ierr = PetscSectionGetDof(oPoints[cl], &fDim);CHKERRQ(ierr); 260 if (pointDim) { 261 for (PetscInt d = 0; d < fDim; ++d, ++v) { 262 values[v] = (*this->_options.integrate)(v0, J, v, initFunc); 263 } 264 } 265 } 266 ierr = DMPlexVecSetClosure(dm, NULL, localX, c, values);CHKERRQ(ierr); 267 pV.clear(); 268 } 269 ierr = PetscFree2(v0,J);CHKERRQ(ierr); 270 ierr = PetscFree(values);CHKERRQ(ierr); 271 #endif 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNCT__ 276 #define __FUNCT__ "DMPlexProjectFunction" 277 /*@C 278 DMPlexProjectFunction - This projects the given function into the function space provided. 279 280 Input Parameters: 281 + dm - The DM 282 . numComp - The number of components (functions) 283 . funcs - The coordinate functions to evaluate 284 - mode - The insertion mode for values 285 286 Output Parameter: 287 . X - vector 288 289 Level: developer 290 291 Note: 292 This currently just calls the function with the coordinates of each vertex and edge midpoint, and stores the result in a vector. 293 We will eventually fix it. 294 295 .seealso: DMPlexComputeL2Diff() 296 @*/ 297 PetscErrorCode DMPlexProjectFunction(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X) 298 { 299 Vec localX; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 304 ierr = DMPlexProjectFunctionLocal(dm, numComp, funcs, mode, localX);CHKERRQ(ierr); 305 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 306 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 307 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "DMPlexComputeL2Diff" 313 /*@C 314 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 315 316 Input Parameters: 317 + dm - The DM 318 . quad - The PetscQuadrature object for each field 319 . funcs - The functions to evaluate for each field component 320 - X - The coefficient vector u_h 321 322 Output Parameter: 323 . diff - The diff ||u - u_h||_2 324 325 Level: developer 326 327 .seealso: DMPlexProjectFunction() 328 @*/ 329 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscQuadrature quad[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff) 330 { 331 const PetscInt debug = 0; 332 PetscSection section; 333 Vec localX; 334 PetscReal *coords, *v0, *J, *invJ, detJ; 335 PetscReal localDiff = 0.0; 336 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 337 PetscErrorCode ierr; 338 339 PetscFunctionBegin; 340 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 341 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 342 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 343 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 344 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 345 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 346 for (field = 0; field < numFields; ++field) { 347 numComponents += quad[field].numComponents; 348 } 349 ierr = DMPlexProjectFunctionLocal(dm, numComponents, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 350 ierr = PetscMalloc4(dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr); 351 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 352 for (c = cStart; c < cEnd; ++c) { 353 PetscScalar *x; 354 PetscReal elemDiff = 0.0; 355 356 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 357 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 358 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 359 360 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 361 const PetscInt numQuadPoints = quad[field].numQuadPoints; 362 const PetscReal *quadPoints = quad[field].quadPoints; 363 const PetscReal *quadWeights = quad[field].quadWeights; 364 const PetscInt numBasisFuncs = quad[field].numBasisFuncs; 365 const PetscInt numBasisComps = quad[field].numComponents; 366 const PetscReal *basis = quad[field].basis; 367 PetscInt q, d, e, fc, f; 368 369 if (debug) { 370 char title[1024]; 371 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 372 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 373 } 374 for (q = 0; q < numQuadPoints; ++q) { 375 for (d = 0; d < dim; d++) { 376 coords[d] = v0[d]; 377 for (e = 0; e < dim; e++) { 378 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 379 } 380 } 381 for (fc = 0; fc < numBasisComps; ++fc) { 382 PetscScalar funcVal; 383 PetscScalar interpolant = 0.0; 384 385 (*funcs[comp+fc])(coords, &funcVal); 386 for (f = 0; f < numBasisFuncs; ++f) { 387 const PetscInt fidx = f*numBasisComps+fc; 388 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 389 } 390 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ);CHKERRQ(ierr);} 391 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ; 392 } 393 } 394 comp += numBasisComps; 395 fieldOffset += numBasisFuncs*numBasisComps; 396 } 397 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 398 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 399 localDiff += elemDiff; 400 } 401 ierr = PetscFree4(coords,v0,J,invJ);CHKERRQ(ierr); 402 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 403 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 404 *diff = PetscSqrtReal(*diff); 405 PetscFunctionReturn(0); 406 } 407 408 #if 0 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "DMPlexComputeResidualFEM" 412 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 413 { 414 DM_Plex *mesh = (DM_Plex*) dm->data; 415 PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 416 PetscQuadrature *quad = fem->quad; 417 PetscQuadrature *quadBd = fem->quadBd; 418 PetscSection section; 419 PetscReal *v0, *n, *J, *invJ, *detJ; 420 PetscScalar *elemVec, *u; 421 PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 422 PetscInt cellDof, numComponents; 423 PetscBool has; 424 PetscErrorCode ierr; 425 426 PetscFunctionBegin; 427 if (has && quadBd) { 428 DMLabel label; 429 IS pointIS; 430 const PetscInt *points; 431 PetscInt numPoints, p; 432 433 ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 434 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 435 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 436 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 437 for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) { 438 cellDof += quadBd[field].numBasisFuncs*quadBd[field].numComponents; 439 numComponents += quadBd[field].numComponents; 440 } 441 ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 442 for (p = 0; p < numPoints; ++p) { 443 const PetscInt point = points[p]; 444 PetscScalar *x; 445 PetscInt i; 446 447 /* TODO: Add normal determination here */ 448 ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr); 449 if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point); 450 ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 451 452 for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i]; 453 ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 454 } 455 for (field = 0; field < numFields; ++field) { 456 const PetscInt numQuadPoints = quadBd[field].numQuadPoints; 457 const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs; 458 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field]; 459 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field]; 460 /* Conforming batches */ 461 PetscInt blockSize = numBasisFuncs*numQuadPoints; 462 PetscInt numBlocks = 1; 463 PetscInt batchSize = numBlocks * blockSize; 464 PetscInt numBatches = numBatchesTmp; 465 PetscInt numChunks = numPoints / (numBatches*batchSize); 466 /* Remainder */ 467 PetscInt numRemainder = numPoints % (numBatches * batchSize); 468 PetscInt offset = numPoints - numRemainder; 469 470 ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr); 471 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], 472 f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 473 } 474 for (p = 0; p < numPoints; ++p) { 475 const PetscInt point = points[p]; 476 477 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);} 478 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr); 479 } 480 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 481 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 482 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 483 } 484 PetscFunctionReturn(0); 485 } 486 487 #else 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "DMPlexComputeResidualFEM" 491 /*@ 492 DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 493 494 Input Parameters: 495 + dm - The mesh 496 . X - Local input vector 497 - user - The user context 498 499 Output Parameter: 500 . F - Local output vector 501 502 Note: 503 The second member of the user context must be an FEMContext. 504 505 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 506 like a GPU, or vectorize on a multicore machine. 507 508 Level: developer 509 510 .seealso: DMPlexComputeJacobianActionFEM() 511 @*/ 512 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 513 { 514 DM_Plex *mesh = (DM_Plex*) dm->data; 515 PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 516 PetscFE *fe = fem->fe; 517 const char *name = "Residual"; 518 PetscQuadrature q; 519 PetscCellGeometry geom; 520 PetscSection section; 521 PetscReal *v0, *J, *invJ, *detJ; 522 PetscScalar *elemVec, *u; 523 PetscInt dim, numFields, f, numCells, cStart, cEnd, c; 524 PetscInt cellDof = 0, numComponents = 0; 525 PetscErrorCode ierr; 526 527 PetscFunctionBegin; 528 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 529 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 530 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 531 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 532 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 533 numCells = cEnd - cStart; 534 for (f = 0; f < numFields; ++f) { 535 PetscInt Nb, Nc; 536 537 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 538 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 539 cellDof += Nb*Nc; 540 numComponents += Nc; 541 } 542 ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 543 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 544 ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 545 for (c = cStart; c < cEnd; ++c) { 546 PetscScalar *x = NULL; 547 PetscInt i; 548 549 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 550 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 551 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 552 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 553 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 554 } 555 for (f = 0; f < numFields; ++f) { 556 void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[]) = fem->f0Funcs[f]; 557 void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[]) = fem->f1Funcs[f]; 558 PetscInt Nb; 559 /* Conforming batches */ 560 PetscInt numBlocks = 1; 561 PetscInt numBatches = 1; 562 PetscInt numChunks, Ne, blockSize, batchSize; 563 /* Remainder */ 564 PetscInt Nr, offset; 565 566 ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 567 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 568 blockSize = Nb*q.numQuadPoints; 569 batchSize = numBlocks * blockSize; 570 numChunks = numCells / (numBatches*batchSize); 571 Ne = numChunks*numBatches*batchSize; 572 Nr = numCells % (numBatches*batchSize); 573 offset = numCells - Nr; 574 geom.v0 = v0; 575 geom.J = J; 576 geom.invJ = invJ; 577 geom.detJ = detJ; 578 ierr = PetscFEIntegrateResidual(fe[f], Ne, numFields, fe, f, geom, u, f0, f1, elemVec);CHKERRQ(ierr); 579 geom.v0 = &v0[offset*dim]; 580 geom.J = &J[offset*dim*dim]; 581 geom.invJ = &invJ[offset*dim*dim]; 582 geom.detJ = &detJ[offset]; 583 ierr = PetscFEIntegrateResidual(fe[f], Nr, numFields, fe, f, geom, &u[offset*cellDof], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 584 } 585 for (c = cStart; c < cEnd; ++c) { 586 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 587 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 588 } 589 ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 590 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, F);CHKERRQ(ierr);} 591 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 592 PetscFunctionReturn(0); 593 } 594 595 #endif 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 599 /*@C 600 DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 601 602 Input Parameters: 603 + dm - The mesh 604 . J - The Jacobian shell matrix 605 . X - Local input vector 606 - user - The user context 607 608 Output Parameter: 609 . F - Local output vector 610 611 Note: 612 The second member of the user context must be an FEMContext. 613 614 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 615 like a GPU, or vectorize on a multicore machine. 616 617 Level: developer 618 619 .seealso: DMPlexComputeResidualFEM() 620 @*/ 621 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 622 { 623 DM_Plex *mesh = (DM_Plex*) dm->data; 624 PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 625 PetscFE *fe = fem->fe; 626 PetscQuadrature quad; 627 PetscCellGeometry geom; 628 PetscSection section; 629 JacActionCtx *jctx; 630 PetscReal *v0, *J, *invJ, *detJ; 631 PetscScalar *elemVec, *u, *a; 632 PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 633 PetscInt cellDof = 0; 634 PetscErrorCode ierr; 635 636 PetscFunctionBegin; 637 /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 638 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 639 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 640 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 641 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 642 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 643 numCells = cEnd - cStart; 644 for (field = 0; field < numFields; ++field) { 645 PetscInt Nb, Nc; 646 647 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 648 ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 649 cellDof += Nb*Nc; 650 } 651 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 652 ierr = PetscMalloc7(numCells*cellDof,PetscScalar,&u,numCells*cellDof,PetscScalar,&a,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 653 for (c = cStart; c < cEnd; ++c) { 654 PetscScalar *x; 655 PetscInt i; 656 657 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 658 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 659 ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 660 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 661 ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 662 ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 663 for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 664 ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 665 } 666 for (field = 0; field < numFields; ++field) { 667 PetscInt Nb; 668 /* Conforming batches */ 669 PetscInt numBlocks = 1; 670 PetscInt numBatches = 1; 671 PetscInt numChunks, Ne, blockSize, batchSize; 672 /* Remainder */ 673 PetscInt Nr, offset; 674 675 ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 676 ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 677 blockSize = Nb*quad.numQuadPoints; 678 batchSize = numBlocks * blockSize; 679 numChunks = numCells / (numBatches*batchSize); 680 Ne = numChunks*numBatches*batchSize; 681 Nr = numCells % (numBatches*batchSize); 682 offset = numCells - Nr; 683 geom.v0 = v0; 684 geom.J = J; 685 geom.invJ = invJ; 686 geom.detJ = detJ; 687 ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 688 geom.v0 = &v0[offset*dim]; 689 geom.J = &J[offset*dim*dim]; 690 geom.invJ = &invJ[offset*dim*dim]; 691 geom.detJ = &detJ[offset]; 692 ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 693 fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 694 } 695 for (c = cStart; c < cEnd; ++c) { 696 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 697 ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 698 } 699 ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 700 if (mesh->printFEM) { 701 PetscMPIInt rank, numProcs; 702 PetscInt p; 703 704 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 705 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 706 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 707 for (p = 0; p < numProcs; ++p) { 708 if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 709 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 710 } 711 } 712 /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 713 PetscFunctionReturn(0); 714 } 715 716 #undef __FUNCT__ 717 #define __FUNCT__ "DMPlexComputeJacobianFEM" 718 /*@ 719 DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 720 721 Input Parameters: 722 + dm - The mesh 723 . X - Local input vector 724 - user - The user context 725 726 Output Parameter: 727 . Jac - Jacobian matrix 728 729 Note: 730 The second member of the user context must be an FEMContext. 731 732 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 733 like a GPU, or vectorize on a multicore machine. 734 735 Level: developer 736 737 .seealso: FormFunctionLocal() 738 @*/ 739 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 740 { 741 DM_Plex *mesh = (DM_Plex*) dm->data; 742 PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 743 PetscFE *fe = fem->fe; 744 const char *name = "Jacobian"; 745 PetscQuadrature quad; 746 PetscCellGeometry geom; 747 PetscSection section, globalSection; 748 PetscReal *v0, *J, *invJ, *detJ; 749 PetscScalar *elemMat, *u; 750 PetscInt dim, numFields, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 751 PetscInt cellDof = 0, numComponents = 0; 752 PetscBool isShell; 753 PetscErrorCode ierr; 754 755 PetscFunctionBegin; 756 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 757 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 758 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 759 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 760 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 761 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 762 numCells = cEnd - cStart; 763 for (f = 0; f < numFields; ++f) { 764 PetscInt Nb, Nc; 765 766 ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 767 ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 768 cellDof += Nb*Nc; 769 numComponents += Nc; 770 } 771 ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 772 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 773 ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof*cellDof,PetscScalar,&elemMat);CHKERRQ(ierr); 774 for (c = cStart; c < cEnd; ++c) { 775 PetscScalar *x; 776 PetscInt i; 777 778 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 779 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 780 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 781 for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 782 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 783 } 784 ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 785 for (fieldI = 0; fieldI < numFields; ++fieldI) { 786 PetscInt Nb; 787 ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 788 ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 789 for (fieldJ = 0; fieldJ < numFields; ++fieldJ) { 790 void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*numFields+fieldJ]; 791 void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*numFields+fieldJ]; 792 void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*numFields+fieldJ]; 793 void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*numFields+fieldJ]; 794 /* Conforming batches */ 795 PetscInt numBlocks = 1; 796 PetscInt numBatches = 1; 797 PetscInt numChunks, Ne, blockSize, batchSize; 798 /* Remainder */ 799 PetscInt Nr, offset; 800 801 blockSize = Nb*quad.numQuadPoints; 802 batchSize = numBlocks * blockSize; 803 numChunks = numCells / (numBatches*batchSize); 804 Ne = numChunks*numBatches*batchSize; 805 Nr = numCells % (numBatches*batchSize); 806 offset = numCells - Nr; 807 geom.v0 = v0; 808 geom.J = J; 809 geom.invJ = invJ; 810 geom.detJ = detJ; 811 ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, numFields, fe, fieldI, fieldJ, geom, u, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 812 geom.v0 = &v0[offset*dim]; 813 geom.J = &J[offset*dim*dim]; 814 geom.invJ = &invJ[offset*dim*dim]; 815 geom.detJ = &detJ[offset]; 816 ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, numFields, fe, fieldI, fieldJ, geom, &u[offset*cellDof], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 817 } 818 } 819 for (c = cStart; c < cEnd; ++c) { 820 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 821 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 822 } 823 ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 824 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 825 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 826 if (mesh->printFEM) { 827 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 828 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 829 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 830 } 831 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 832 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 833 if (isShell) { 834 JacActionCtx *jctx; 835 836 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 837 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 838 } 839 *str = SAME_NONZERO_PATTERN; 840 PetscFunctionReturn(0); 841 } 842