1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petsc-private/petscfeimpl.h> 4 #include <petscfv.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "DMPlexGetScale" 8 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 9 { 10 DM_Plex *mesh = (DM_Plex*) dm->data; 11 12 PetscFunctionBegin; 13 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14 PetscValidPointer(scale, 3); 15 *scale = mesh->scale[unit]; 16 PetscFunctionReturn(0); 17 } 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "DMPlexSetScale" 21 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 22 { 23 DM_Plex *mesh = (DM_Plex*) dm->data; 24 25 PetscFunctionBegin; 26 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 27 mesh->scale[unit] = scale; 28 PetscFunctionReturn(0); 29 } 30 31 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 32 { 33 switch (i) { 34 case 0: 35 switch (j) { 36 case 0: return 0; 37 case 1: 38 switch (k) { 39 case 0: return 0; 40 case 1: return 0; 41 case 2: return 1; 42 } 43 case 2: 44 switch (k) { 45 case 0: return 0; 46 case 1: return -1; 47 case 2: return 0; 48 } 49 } 50 case 1: 51 switch (j) { 52 case 0: 53 switch (k) { 54 case 0: return 0; 55 case 1: return 0; 56 case 2: return -1; 57 } 58 case 1: return 0; 59 case 2: 60 switch (k) { 61 case 0: return 1; 62 case 1: return 0; 63 case 2: return 0; 64 } 65 } 66 case 2: 67 switch (j) { 68 case 0: 69 switch (k) { 70 case 0: return 0; 71 case 1: return 1; 72 case 2: return 0; 73 } 74 case 1: 75 switch (k) { 76 case 0: return -1; 77 case 1: return 0; 78 case 2: return 0; 79 } 80 case 2: return 0; 81 } 82 } 83 return 0; 84 } 85 86 #undef __FUNCT__ 87 #define __FUNCT__ "DMPlexCreateRigidBody" 88 /*@C 89 DMPlexCreateRigidBody - create rigid body modes from coordinates 90 91 Collective on DM 92 93 Input Arguments: 94 + dm - the DM 95 . section - the local section associated with the rigid field, or NULL for the default section 96 - globalSection - the global section associated with the rigid field, or NULL for the default section 97 98 Output Argument: 99 . sp - the null space 100 101 Note: This is necessary to take account of Dirichlet conditions on the displacements 102 103 Level: advanced 104 105 .seealso: MatNullSpaceCreate() 106 @*/ 107 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 108 { 109 MPI_Comm comm; 110 Vec coordinates, localMode, mode[6]; 111 PetscSection coordSection; 112 PetscScalar *coords; 113 PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 118 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 119 if (dim == 1) { 120 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 124 if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 125 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 126 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 127 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 128 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 129 m = (dim*(dim+1))/2; 130 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 131 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 132 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 133 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 134 /* Assume P1 */ 135 ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 136 for (d = 0; d < dim; ++d) { 137 PetscScalar values[3] = {0.0, 0.0, 0.0}; 138 139 values[d] = 1.0; 140 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 141 for (v = vStart; v < vEnd; ++v) { 142 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 143 } 144 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 146 } 147 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 148 for (d = dim; d < dim*(dim+1)/2; ++d) { 149 PetscInt i, j, k = dim > 2 ? d - dim : d; 150 151 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 152 for (v = vStart; v < vEnd; ++v) { 153 PetscScalar values[3] = {0.0, 0.0, 0.0}; 154 PetscInt off; 155 156 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 157 for (i = 0; i < dim; ++i) { 158 for (j = 0; j < dim; ++j) { 159 values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 160 } 161 } 162 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 163 } 164 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 166 } 167 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 168 ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 169 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 170 /* Orthonormalize system */ 171 for (i = dim; i < m; ++i) { 172 PetscScalar dots[6]; 173 174 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 175 for (j = 0; j < i; ++j) dots[j] *= -1.0; 176 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 177 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 178 } 179 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 180 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal" 186 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 187 { 188 PetscDualSpace *sp; 189 PetscSection section; 190 PetscScalar *values; 191 PetscReal *v0, *J, detJ; 192 PetscBool *fieldActive; 193 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 198 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 199 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 200 ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr); 201 for (f = 0; f < numFields; ++f) { 202 ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 203 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 204 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 205 totDim += spDim*numComp; 206 } 207 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 208 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 209 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 210 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 211 ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 212 for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE; 213 for (i = 0; i < numIds; ++i) { 214 IS pointIS; 215 const PetscInt *points; 216 PetscInt n, p; 217 218 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 219 ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 220 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 221 for (p = 0; p < n; ++p) { 222 const PetscInt point = points[p]; 223 PetscCellGeometry geom; 224 225 if ((point < cStart) || (point >= cEnd)) continue; 226 ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr); 227 geom.v0 = v0; 228 geom.J = J; 229 geom.detJ = &detJ; 230 for (f = 0, v = 0; f < numFields; ++f) { 231 void * const ctx = ctxs ? ctxs[f] : NULL; 232 ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 233 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 234 for (d = 0; d < spDim; ++d) { 235 if (funcs[f]) { 236 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 237 } else { 238 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 239 } 240 v += numComp; 241 } 242 } 243 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 244 } 245 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 246 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 247 } 248 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 249 ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 250 ierr = PetscFree3(sp,v0,J);CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "DMPlexProjectFunctionLocal" 256 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 257 { 258 PetscDualSpace *sp; 259 PetscSection section; 260 PetscScalar *values; 261 PetscReal *v0, *J, detJ; 262 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp; 263 PetscErrorCode ierr; 264 265 PetscFunctionBegin; 266 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 267 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 268 ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr); 269 for (f = 0; f < numFields; ++f) { 270 PetscFE fe; 271 272 ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 273 ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr); 274 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 275 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 276 totDim += spDim*numComp; 277 } 278 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 279 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 280 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 281 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 282 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 283 ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr); 284 for (c = cStart; c < cEnd; ++c) { 285 PetscCellGeometry geom; 286 287 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 288 geom.v0 = v0; 289 geom.J = J; 290 geom.detJ = &detJ; 291 for (f = 0, v = 0; f < numFields; ++f) { 292 PetscFE fe; 293 void * const ctx = ctxs ? ctxs[f] : NULL; 294 295 ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 296 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 297 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 298 for (d = 0; d < spDim; ++d) { 299 if (funcs[f]) { 300 ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 301 } else { 302 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 303 } 304 v += numComp; 305 } 306 } 307 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 308 } 309 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 310 ierr = PetscFree2(v0,J);CHKERRQ(ierr); 311 ierr = PetscFree(sp);CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "DMPlexProjectFunction" 317 /*@C 318 DMPlexProjectFunction - This projects the given function into the function space provided. 319 320 Input Parameters: 321 + dm - The DM 322 . funcs - The coordinate functions to evaluate, one per field 323 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 324 - mode - The insertion mode for values 325 326 Output Parameter: 327 . X - vector 328 329 Level: developer 330 331 .seealso: DMPlexComputeL2Diff() 332 @*/ 333 PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 334 { 335 Vec localX; 336 PetscErrorCode ierr; 337 338 PetscFunctionBegin; 339 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 340 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 341 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr); 342 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 343 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 344 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNCT__ 349 #define __FUNCT__ "DMPlexProjectFieldLocal" 350 PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX) 351 { 352 DM dmAux; 353 PetscDS prob, probAux; 354 Vec A; 355 PetscSection section, sectionAux; 356 PetscScalar *values, *u, *u_x, *a, *a_x; 357 PetscReal *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux; 358 PetscInt Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp; 359 PetscErrorCode ierr; 360 361 PetscFunctionBegin; 362 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 363 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 364 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 365 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 366 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 367 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 368 ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr); 369 ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 370 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 371 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 372 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 373 if (dmAux) { 374 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 375 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 376 ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr); 377 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 378 } 379 ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr); 380 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 381 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 382 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 383 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 384 for (c = cStart; c < cEnd; ++c) { 385 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 386 387 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 388 ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 389 if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 390 for (f = 0, v = 0; f < Nf; ++f) { 391 PetscFE fe; 392 PetscDualSpace sp; 393 PetscInt Ncf; 394 395 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 396 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 397 ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr); 398 ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 399 for (d = 0; d < spDim; ++d) { 400 PetscQuadrature quad; 401 const PetscReal *points, *weights; 402 PetscInt numPoints, q; 403 404 if (funcs[f]) { 405 ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr); 406 ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 407 ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 408 for (q = 0; q < numPoints; ++q) { 409 CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x); 410 ierr = EvaluateFieldJets(prob, PETSC_FALSE, q, invJ, coefficients, NULL, u, u_x, NULL);CHKERRQ(ierr); 411 ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr); 412 (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]); 413 } 414 ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 415 } else { 416 for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0; 417 } 418 v += Ncf; 419 } 420 } 421 ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 422 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 423 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 424 } 425 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 426 ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "DMPlexProjectField" 432 /*@C 433 DMPlexProjectField - This projects the given function of the fields into the function space provided. 434 435 Input Parameters: 436 + dm - The DM 437 . U - The input field vector 438 . funcs - The functions to evaluate, one per field 439 - mode - The insertion mode for values 440 441 Output Parameter: 442 . X - The output vector 443 444 Level: developer 445 446 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 447 @*/ 448 PetscErrorCode DMPlexProjectField(DM dm, Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X) 449 { 450 Vec localX, localU; 451 PetscErrorCode ierr; 452 453 PetscFunctionBegin; 454 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 455 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 456 ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr); 457 ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr); 458 ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr); 459 ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr); 460 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 461 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 462 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 463 ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM" 469 PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX) 470 { 471 void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 472 void **ctxs; 473 PetscFE *fe; 474 PetscInt numFields, f, numBd, b; 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 479 PetscValidHeaderSpecific(localX, VEC_CLASSID, 2); 480 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 481 ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 482 for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);} 483 /* OPT: Could attempt to do multiple BCs at once */ 484 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 485 for (b = 0; b < numBd; ++b) { 486 DMLabel label; 487 const PetscInt *ids; 488 const char *labelname; 489 PetscInt numids, field; 490 PetscBool isEssential; 491 void (*func)(); 492 void *ctx; 493 494 ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 495 ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 496 for (f = 0; f < numFields; ++f) { 497 funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL; 498 ctxs[f] = field == f ? ctx : NULL; 499 } 500 ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 501 } 502 ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505 506 #undef __FUNCT__ 507 #define __FUNCT__ "DMPlexComputeL2Diff" 508 /*@C 509 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 510 511 Input Parameters: 512 + dm - The DM 513 . funcs - The functions to evaluate for each field component 514 . ctxs - Optional array of contexts to pass to each function, or NULL. 515 - X - The coefficient vector u_h 516 517 Output Parameter: 518 . diff - The diff ||u - u_h||_2 519 520 Level: developer 521 522 .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff() 523 @*/ 524 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 525 { 526 const PetscInt debug = 0; 527 PetscSection section; 528 PetscQuadrature quad; 529 Vec localX; 530 PetscScalar *funcVal; 531 PetscReal *coords, *v0, *J, *invJ, detJ; 532 PetscReal localDiff = 0.0; 533 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 534 PetscErrorCode ierr; 535 536 PetscFunctionBegin; 537 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 538 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 539 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 540 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 541 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 542 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 543 for (field = 0; field < numFields; ++field) { 544 PetscFE fe; 545 PetscInt Nc; 546 547 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 548 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 549 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 550 numComponents += Nc; 551 } 552 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 553 ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 554 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 555 for (c = cStart; c < cEnd; ++c) { 556 PetscScalar *x = NULL; 557 PetscReal elemDiff = 0.0; 558 559 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 560 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 561 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 562 563 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 564 PetscFE fe; 565 void * const ctx = ctxs ? ctxs[field] : NULL; 566 const PetscReal *quadPoints, *quadWeights; 567 PetscReal *basis; 568 PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 569 570 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 571 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 572 ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr); 573 ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr); 574 ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr); 575 if (debug) { 576 char title[1024]; 577 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 578 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 579 } 580 for (q = 0; q < numQuadPoints; ++q) { 581 for (d = 0; d < dim; d++) { 582 coords[d] = v0[d]; 583 for (e = 0; e < dim; e++) { 584 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 585 } 586 } 587 (*funcs[field])(coords, funcVal, ctx); 588 for (fc = 0; fc < numBasisComps; ++fc) { 589 PetscScalar interpolant = 0.0; 590 591 for (f = 0; f < numBasisFuncs; ++f) { 592 const PetscInt fidx = f*numBasisComps+fc; 593 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 594 } 595 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);} 596 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 597 } 598 } 599 comp += numBasisComps; 600 fieldOffset += numBasisFuncs*numBasisComps; 601 } 602 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 603 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 604 localDiff += elemDiff; 605 } 606 ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 607 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 608 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 609 *diff = PetscSqrtReal(*diff); 610 PetscFunctionReturn(0); 611 } 612 613 #undef __FUNCT__ 614 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 615 /*@C 616 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 617 618 Input Parameters: 619 + dm - The DM 620 . funcs - The gradient functions to evaluate for each field component 621 . ctxs - Optional array of contexts to pass to each function, or NULL. 622 . X - The coefficient vector u_h 623 - n - The vector to project along 624 625 Output Parameter: 626 . diff - The diff ||(grad u - grad u_h) . n||_2 627 628 Level: developer 629 630 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 631 @*/ 632 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 633 { 634 const PetscInt debug = 0; 635 PetscSection section; 636 PetscQuadrature quad; 637 Vec localX; 638 PetscScalar *funcVal, *interpolantVec; 639 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 640 PetscReal localDiff = 0.0; 641 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 642 PetscErrorCode ierr; 643 644 PetscFunctionBegin; 645 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 646 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 647 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 648 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 649 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 650 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 651 for (field = 0; field < numFields; ++field) { 652 PetscFE fe; 653 PetscInt Nc; 654 655 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 656 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 657 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 658 numComponents += Nc; 659 } 660 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 661 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 662 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 663 for (c = cStart; c < cEnd; ++c) { 664 PetscScalar *x = NULL; 665 PetscReal elemDiff = 0.0; 666 667 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 668 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 669 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 670 671 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 672 PetscFE fe; 673 void * const ctx = ctxs ? ctxs[field] : NULL; 674 const PetscReal *quadPoints, *quadWeights; 675 PetscReal *basisDer; 676 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 677 678 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 679 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 680 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 681 ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr); 682 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 683 if (debug) { 684 char title[1024]; 685 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 686 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 687 } 688 for (q = 0; q < numQuadPoints; ++q) { 689 for (d = 0; d < dim; d++) { 690 coords[d] = v0[d]; 691 for (e = 0; e < dim; e++) { 692 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 693 } 694 } 695 (*funcs[field])(coords, n, funcVal, ctx); 696 for (fc = 0; fc < Ncomp; ++fc) { 697 PetscScalar interpolant = 0.0; 698 699 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 700 for (f = 0; f < Nb; ++f) { 701 const PetscInt fidx = f*Ncomp+fc; 702 703 for (d = 0; d < dim; ++d) { 704 realSpaceDer[d] = 0.0; 705 for (g = 0; g < dim; ++g) { 706 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 707 } 708 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 709 } 710 } 711 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 712 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 713 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 714 } 715 } 716 comp += Ncomp; 717 fieldOffset += Nb*Ncomp; 718 } 719 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 720 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 721 localDiff += elemDiff; 722 } 723 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 724 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 725 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 726 *diff = PetscSqrtReal(*diff); 727 PetscFunctionReturn(0); 728 } 729 730 #undef __FUNCT__ 731 #define __FUNCT__ "DMPlexComputeL2FieldDiff" 732 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 733 { 734 const PetscInt debug = 0; 735 PetscSection section; 736 PetscQuadrature quad; 737 Vec localX; 738 PetscScalar *funcVal; 739 PetscReal *coords, *v0, *J, *invJ, detJ; 740 PetscReal *localDiff; 741 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 742 PetscErrorCode ierr; 743 744 PetscFunctionBegin; 745 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 746 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 747 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 748 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 749 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 750 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 751 for (field = 0; field < numFields; ++field) { 752 PetscFE fe; 753 PetscInt Nc; 754 755 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 756 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 757 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 758 numComponents += Nc; 759 } 760 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 761 ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 762 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 763 for (c = cStart; c < cEnd; ++c) { 764 PetscScalar *x = NULL; 765 766 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 767 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 768 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 769 770 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 771 PetscFE fe; 772 void * const ctx = ctxs ? ctxs[field] : NULL; 773 const PetscReal *quadPoints, *quadWeights; 774 PetscReal *basis, elemDiff = 0.0; 775 PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 776 777 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 778 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 779 ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr); 780 ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr); 781 ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr); 782 if (debug) { 783 char title[1024]; 784 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 785 ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 786 } 787 for (q = 0; q < numQuadPoints; ++q) { 788 for (d = 0; d < dim; d++) { 789 coords[d] = v0[d]; 790 for (e = 0; e < dim; e++) { 791 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 792 } 793 } 794 (*funcs[field])(coords, funcVal, ctx); 795 for (fc = 0; fc < numBasisComps; ++fc) { 796 PetscScalar interpolant = 0.0; 797 798 for (f = 0; f < numBasisFuncs; ++f) { 799 const PetscInt fidx = f*numBasisComps+fc; 800 interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 801 } 802 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);} 803 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 804 } 805 } 806 comp += numBasisComps; 807 fieldOffset += numBasisFuncs*numBasisComps; 808 localDiff[field] += elemDiff; 809 } 810 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 811 } 812 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 813 ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 814 for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 815 ierr = PetscFree6(localDiff,funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 816 PetscFunctionReturn(0); 817 } 818 819 #undef __FUNCT__ 820 #define __FUNCT__ "DMPlexComputeIntegralFEM" 821 /*@ 822 DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 823 824 Input Parameters: 825 + dm - The mesh 826 . X - Local input vector 827 - user - The user context 828 829 Output Parameter: 830 . integral - Local integral for each field 831 832 Level: developer 833 834 .seealso: DMPlexComputeResidualFEM() 835 @*/ 836 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 837 { 838 DM_Plex *mesh = (DM_Plex *) dm->data; 839 DM dmAux; 840 Vec localX, A; 841 PetscDS prob, probAux = NULL; 842 PetscQuadrature q; 843 PetscCellGeometry geom; 844 PetscSection section, sectionAux; 845 PetscReal *v0, *J, *invJ, *detJ; 846 PetscScalar *u, *a = NULL; 847 PetscInt dim, Nf, f, numCells, cStart, cEnd, c; 848 PetscInt totDim, totDimAux; 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 853 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 854 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 855 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 856 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 857 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 858 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 859 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 860 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 861 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 862 numCells = cEnd - cStart; 863 for (f = 0; f < Nf; ++f) {integral[f] = 0.0;} 864 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 865 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 866 if (dmAux) { 867 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 868 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 869 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 870 } 871 ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr); 872 ierr = PetscMalloc5(numCells*totDim,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ);CHKERRQ(ierr); 873 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 874 for (c = cStart; c < cEnd; ++c) { 875 PetscScalar *x = NULL; 876 PetscInt i; 877 878 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 879 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 880 ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 881 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 882 ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 883 if (dmAux) { 884 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 885 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 886 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 887 } 888 } 889 for (f = 0; f < Nf; ++f) { 890 PetscFE fe; 891 PetscInt numQuadPoints, Nb; 892 /* Conforming batches */ 893 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 894 /* Remainder */ 895 PetscInt Nr, offset; 896 897 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 898 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 899 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 900 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 901 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 902 blockSize = Nb*numQuadPoints; 903 batchSize = numBlocks * blockSize; 904 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 905 numChunks = numCells / (numBatches*batchSize); 906 Ne = numChunks*numBatches*batchSize; 907 Nr = numCells % (numBatches*batchSize); 908 offset = numCells - Nr; 909 geom.v0 = v0; 910 geom.J = J; 911 geom.invJ = invJ; 912 geom.detJ = detJ; 913 ierr = PetscFEIntegrate(fe, prob, f, Ne, geom, u, probAux, a, integral);CHKERRQ(ierr); 914 geom.v0 = &v0[offset*dim]; 915 geom.J = &J[offset*dim*dim]; 916 geom.invJ = &invJ[offset*dim*dim]; 917 geom.detJ = &detJ[offset]; 918 ierr = PetscFEIntegrate(fe, prob, f, Nr, geom, &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr); 919 } 920 ierr = PetscFree5(u,v0,J,invJ,detJ);CHKERRQ(ierr); 921 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 922 if (mesh->printFEM) { 923 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 924 for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);} 925 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 926 } 927 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 928 /* TODO: Allreduce for integral */ 929 /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "DMPlexComputeResidualFEM_Internal" 935 PetscErrorCode DMPlexComputeResidualFEM_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user) 936 { 937 DM_Plex *mesh = (DM_Plex *) dm->data; 938 const char *name = "Residual"; 939 DM dmAux; 940 DMLabel depth; 941 Vec A; 942 PetscDS prob, probAux = NULL; 943 PetscQuadrature q; 944 PetscCellGeometry geom; 945 PetscSection section, sectionAux; 946 PetscReal *v0, *J, *invJ, *detJ; 947 PetscScalar *elemVec, *u, *u_t, *a = NULL; 948 PetscInt dim, Nf, f, numCells, cStart, cEnd, c, numBd, bd; 949 PetscInt totDim, totDimBd, totDimAux; 950 PetscErrorCode ierr; 951 952 PetscFunctionBegin; 953 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 954 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 955 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 956 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 957 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 958 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 959 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 960 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 961 numCells = cEnd - cStart; 962 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 963 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 964 if (dmAux) { 965 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 966 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 967 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 968 } 969 ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr); 970 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 971 ierr = PetscMalloc7(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*totDim,&elemVec);CHKERRQ(ierr); 972 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 973 for (c = cStart; c < cEnd; ++c) { 974 PetscScalar *x = NULL, *x_t = NULL; 975 PetscInt i; 976 977 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 978 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 979 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 980 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 981 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 982 if (X_t) { 983 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 984 for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i]; 985 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 986 } 987 if (dmAux) { 988 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 989 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 990 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 991 } 992 } 993 for (f = 0; f < Nf; ++f) { 994 PetscFE fe; 995 PetscInt numQuadPoints, Nb; 996 /* Conforming batches */ 997 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 998 /* Remainder */ 999 PetscInt Nr, offset; 1000 1001 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1002 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1003 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1004 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1005 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1006 blockSize = Nb*numQuadPoints; 1007 batchSize = numBlocks * blockSize; 1008 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1009 numChunks = numCells / (numBatches*batchSize); 1010 Ne = numChunks*numBatches*batchSize; 1011 Nr = numCells % (numBatches*batchSize); 1012 offset = numCells - Nr; 1013 geom.v0 = v0; 1014 geom.J = J; 1015 geom.invJ = invJ; 1016 geom.detJ = detJ; 1017 ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, geom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr); 1018 geom.v0 = &v0[offset*dim]; 1019 geom.J = &J[offset*dim*dim]; 1020 geom.invJ = &invJ[offset*dim*dim]; 1021 geom.detJ = &detJ[offset]; 1022 ierr = PetscFEIntegrateResidual(fe, prob, f, Nr, geom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemVec[offset*totDim]);CHKERRQ(ierr); 1023 } 1024 for (c = cStart; c < cEnd; ++c) { 1025 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, totDim, &elemVec[c*totDim]);CHKERRQ(ierr);} 1026 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_VALUES);CHKERRQ(ierr); 1027 } 1028 ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 1029 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1030 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1031 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1032 for (bd = 0; bd < numBd; ++bd) { 1033 const char *bdLabel; 1034 DMLabel label; 1035 IS pointIS; 1036 const PetscInt *points; 1037 const PetscInt *values; 1038 PetscReal *n; 1039 PetscInt field, numValues, numPoints, p, dep, numFaces; 1040 PetscBool isEssential; 1041 1042 ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1043 if (isEssential) continue; 1044 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 1045 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1046 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1047 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1048 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1049 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1050 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1051 if (dep == dim-1) ++numFaces; 1052 } 1053 ierr = PetscMalloc7(numFaces*totDimBd,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*totDimBd,&elemVec);CHKERRQ(ierr); 1054 if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 1055 for (p = 0, f = 0; p < numPoints; ++p) { 1056 const PetscInt point = points[p]; 1057 PetscScalar *x = NULL; 1058 PetscInt i; 1059 1060 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1061 if (dep != dim-1) continue; 1062 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 1063 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 1064 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 1065 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1066 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 1067 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1068 if (X_t) { 1069 ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1070 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 1071 ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1072 } 1073 ++f; 1074 } 1075 for (f = 0; f < Nf; ++f) { 1076 PetscFE fe; 1077 PetscInt numQuadPoints, Nb; 1078 /* Conforming batches */ 1079 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1080 /* Remainder */ 1081 PetscInt Nr, offset; 1082 1083 ierr = PetscDSGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1084 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1085 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1086 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1087 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1088 blockSize = Nb*numQuadPoints; 1089 batchSize = numBlocks * blockSize; 1090 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1091 numChunks = numFaces / (numBatches*batchSize); 1092 Ne = numChunks*numBatches*batchSize; 1093 Nr = numFaces % (numBatches*batchSize); 1094 offset = numFaces - Nr; 1095 geom.v0 = v0; 1096 geom.n = n; 1097 geom.J = J; 1098 geom.invJ = invJ; 1099 geom.detJ = detJ; 1100 ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, geom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr); 1101 geom.v0 = &v0[offset*dim]; 1102 geom.n = &n[offset*dim]; 1103 geom.J = &J[offset*dim*dim]; 1104 geom.invJ = &invJ[offset*dim*dim]; 1105 geom.detJ = &detJ[offset]; 1106 ierr = PetscFEIntegrateBdResidual(fe, prob, f, Nr, geom, &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemVec[offset*totDimBd]);CHKERRQ(ierr); 1107 } 1108 for (p = 0, f = 0; p < numPoints; ++p) { 1109 const PetscInt point = points[p]; 1110 1111 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1112 if (dep != dim-1) continue; 1113 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);} 1114 ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*totDimBd], ADD_VALUES);CHKERRQ(ierr); 1115 ++f; 1116 } 1117 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1118 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1119 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 1120 if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 1121 } 1122 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 1123 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1124 PetscFunctionReturn(0); 1125 } 1126 1127 #undef __FUNCT__ 1128 #define __FUNCT__ "DMPlexSNESComputeResidualFEM" 1129 /*@ 1130 DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 1131 1132 Input Parameters: 1133 + dm - The mesh 1134 . X - Local solution 1135 - user - The user context 1136 1137 Output Parameter: 1138 . F - Local output vector 1139 1140 Level: developer 1141 1142 .seealso: DMPlexComputeJacobianActionFEM() 1143 @*/ 1144 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 1145 { 1146 PetscErrorCode ierr; 1147 1148 PetscFunctionBegin; 1149 ierr = DMPlexComputeResidualFEM_Internal(dm, X, NULL, F, user);CHKERRQ(ierr); 1150 PetscFunctionReturn(0); 1151 } 1152 1153 #undef __FUNCT__ 1154 #define __FUNCT__ "DMPlexTSComputeIFunctionFEM" 1155 /*@ 1156 DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user 1157 1158 Input Parameters: 1159 + dm - The mesh 1160 . t - The time 1161 . X - Local solution 1162 . X_t - Local solution time derivative, or NULL 1163 - user - The user context 1164 1165 Output Parameter: 1166 . F - Local output vector 1167 1168 Level: developer 1169 1170 .seealso: DMPlexComputeJacobianActionFEM() 1171 @*/ 1172 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user) 1173 { 1174 PetscErrorCode ierr; 1175 1176 PetscFunctionBegin; 1177 ierr = DMPlexComputeResidualFEM_Internal(dm, X, X_t, F, user);CHKERRQ(ierr); 1178 PetscFunctionReturn(0); 1179 } 1180 1181 #undef __FUNCT__ 1182 #define __FUNCT__ "DMPlexComputeJacobianFEM_Internal" 1183 PetscErrorCode DMPlexComputeJacobianFEM_Internal(DM dm, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user) 1184 { 1185 DM_Plex *mesh = (DM_Plex *) dm->data; 1186 const char *name = "Jacobian"; 1187 DM dmAux; 1188 DMLabel depth; 1189 Vec A; 1190 PetscDS prob, probAux = NULL; 1191 PetscQuadrature quad; 1192 PetscCellGeometry geom; 1193 PetscSection section, globalSection, sectionAux; 1194 PetscReal *v0, *J, *invJ, *detJ; 1195 PetscScalar *elemMat, *u, *u_t, *a = NULL; 1196 PetscInt dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 1197 PetscInt totDim, totDimBd, totDimAux, numBd, bd; 1198 PetscBool isShell; 1199 PetscErrorCode ierr; 1200 1201 PetscFunctionBegin; 1202 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 1203 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1204 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1205 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 1206 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1207 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1208 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 1209 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1210 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1211 numCells = cEnd - cStart; 1212 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1213 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1214 if (dmAux) { 1215 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1216 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1217 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1218 } 1219 ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr); 1220 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 1221 ierr = PetscMalloc7(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*totDim*totDim,&elemMat);CHKERRQ(ierr); 1222 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1223 for (c = cStart; c < cEnd; ++c) { 1224 PetscScalar *x = NULL, *x_t = NULL; 1225 PetscInt i; 1226 1227 ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1228 if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1229 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1230 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1231 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1232 if (X_t) { 1233 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1234 for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i]; 1235 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1236 } 1237 if (dmAux) { 1238 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1239 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1240 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1241 } 1242 } 1243 ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1244 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1245 PetscFE fe; 1246 PetscInt numQuadPoints, Nb; 1247 /* Conforming batches */ 1248 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1249 /* Remainder */ 1250 PetscInt Nr, offset; 1251 1252 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1253 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1254 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1255 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1256 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1257 blockSize = Nb*numQuadPoints; 1258 batchSize = numBlocks * blockSize; 1259 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1260 numChunks = numCells / (numBatches*batchSize); 1261 Ne = numChunks*numBatches*batchSize; 1262 Nr = numCells % (numBatches*batchSize); 1263 offset = numCells - Nr; 1264 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1265 geom.v0 = v0; 1266 geom.J = J; 1267 geom.invJ = invJ; 1268 geom.detJ = detJ; 1269 ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr); 1270 geom.v0 = &v0[offset*dim]; 1271 geom.J = &J[offset*dim*dim]; 1272 geom.invJ = &invJ[offset*dim*dim]; 1273 geom.detJ = &detJ[offset]; 1274 ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Nr, geom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 1275 } 1276 } 1277 for (c = cStart; c < cEnd; ++c) { 1278 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[c*totDim*totDim]);CHKERRQ(ierr);} 1279 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 1280 } 1281 ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1282 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1283 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1284 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1285 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1286 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1287 for (bd = 0; bd < numBd; ++bd) { 1288 const char *bdLabel; 1289 DMLabel label; 1290 IS pointIS; 1291 const PetscInt *points; 1292 const PetscInt *values; 1293 PetscReal *n; 1294 PetscInt field, numValues, numPoints, p, dep, numFaces; 1295 PetscBool isEssential; 1296 1297 ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1298 if (isEssential) continue; 1299 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 1300 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1301 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1302 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1303 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1304 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1305 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1306 if (dep == dim-1) ++numFaces; 1307 } 1308 ierr = PetscMalloc7(numFaces*totDimBd,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr); 1309 if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 1310 for (p = 0, f = 0; p < numPoints; ++p) { 1311 const PetscInt point = points[p]; 1312 PetscScalar *x = NULL; 1313 PetscInt i; 1314 1315 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1316 if (dep != dim-1) continue; 1317 ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 1318 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 1319 if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 1320 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1321 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 1322 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1323 if (X_t) { 1324 ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1325 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 1326 ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1327 } 1328 ++f; 1329 } 1330 ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr); 1331 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1332 PetscFE fe; 1333 PetscInt numQuadPoints, Nb; 1334 /* Conforming batches */ 1335 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1336 /* Remainder */ 1337 PetscInt Nr, offset; 1338 1339 ierr = PetscDSGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1340 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1341 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1342 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1343 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1344 blockSize = Nb*numQuadPoints; 1345 batchSize = numBlocks * blockSize; 1346 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1347 numChunks = numFaces / (numBatches*batchSize); 1348 Ne = numChunks*numBatches*batchSize; 1349 Nr = numFaces % (numBatches*batchSize); 1350 offset = numFaces - Nr; 1351 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1352 geom.v0 = v0; 1353 geom.n = n; 1354 geom.J = J; 1355 geom.invJ = invJ; 1356 geom.detJ = detJ; 1357 ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr); 1358 geom.v0 = &v0[offset*dim]; 1359 geom.n = &n[offset*dim]; 1360 geom.J = &J[offset*dim*dim]; 1361 geom.invJ = &invJ[offset*dim*dim]; 1362 geom.detJ = &detJ[offset]; 1363 ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Nr, geom, &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemMat[offset*totDimBd*totDimBd]);CHKERRQ(ierr); 1364 } 1365 } 1366 for (p = 0, f = 0; p < numPoints; ++p) { 1367 const PetscInt point = points[p]; 1368 1369 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1370 if (dep != dim-1) continue; 1371 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);} 1372 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr); 1373 ++f; 1374 } 1375 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1376 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1377 ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr); 1378 if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 1379 } 1380 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1381 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1382 if (mesh->printFEM) { 1383 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1384 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 1385 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1386 } 1387 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 1388 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 1389 if (isShell) { 1390 JacActionCtx *jctx; 1391 1392 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 1393 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 1394 } 1395 PetscFunctionReturn(0); 1396 } 1397 1398 #undef __FUNCT__ 1399 #define __FUNCT__ "DMPlexSNESComputeJacobianFEM" 1400 /*@ 1401 DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 1402 1403 Input Parameters: 1404 + dm - The mesh 1405 . X - Local input vector 1406 - user - The user context 1407 1408 Output Parameter: 1409 . Jac - Jacobian matrix 1410 1411 Note: 1412 The first member of the user context must be an FEMContext. 1413 1414 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1415 like a GPU, or vectorize on a multicore machine. 1416 1417 Level: developer 1418 1419 .seealso: FormFunctionLocal() 1420 @*/ 1421 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 1422 { 1423 PetscErrorCode ierr; 1424 1425 PetscFunctionBegin; 1426 ierr = DMPlexComputeJacobianFEM_Internal(dm, X, NULL, Jac, JacP, user);CHKERRQ(ierr); 1427 PetscFunctionReturn(0); 1428 } 1429 1430 #undef __FUNCT__ 1431 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1432 /*@ 1433 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1434 1435 Input Parameters: 1436 + dmf - The fine mesh 1437 . dmc - The coarse mesh 1438 - user - The user context 1439 1440 Output Parameter: 1441 . In - The interpolation matrix 1442 1443 Note: 1444 The first member of the user context must be an FEMContext. 1445 1446 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1447 like a GPU, or vectorize on a multicore machine. 1448 1449 Level: developer 1450 1451 .seealso: DMPlexComputeJacobianFEM() 1452 @*/ 1453 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1454 { 1455 DM_Plex *mesh = (DM_Plex *) dmc->data; 1456 const char *name = "Interpolator"; 1457 PetscDS prob; 1458 PetscFE *feRef; 1459 PetscSection fsection, fglobalSection; 1460 PetscSection csection, cglobalSection; 1461 PetscScalar *elemMat; 1462 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 1463 PetscInt cTotDim, rTotDim = 0; 1464 PetscErrorCode ierr; 1465 1466 PetscFunctionBegin; 1467 #if 0 1468 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1469 #endif 1470 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1471 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1472 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1473 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1474 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1475 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1476 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1477 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1478 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1479 for (f = 0; f < Nf; ++f) { 1480 PetscFE fe; 1481 PetscInt rNb, Nc; 1482 1483 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1484 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1485 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1486 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1487 rTotDim += rNb*Nc; 1488 } 1489 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1490 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1491 ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 1492 ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1493 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1494 PetscDualSpace Qref; 1495 PetscQuadrature f; 1496 const PetscReal *qpoints, *qweights; 1497 PetscReal *points; 1498 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1499 1500 /* Compose points from all dual basis functionals */ 1501 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1502 ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1503 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1504 for (i = 0; i < fpdim; ++i) { 1505 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1506 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1507 npoints += Np; 1508 } 1509 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1510 for (i = 0, k = 0; i < fpdim; ++i) { 1511 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1512 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1513 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1514 } 1515 1516 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1517 PetscFE fe; 1518 PetscReal *B; 1519 PetscInt NcJ, cpdim, j; 1520 1521 /* Evaluate basis at points */ 1522 ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr); 1523 ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 1524 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1525 /* For now, fields only interpolate themselves */ 1526 if (fieldI == fieldJ) { 1527 if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ); 1528 ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1529 for (i = 0, k = 0; i < fpdim; ++i) { 1530 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1531 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1532 for (p = 0; p < Np; ++p, ++k) { 1533 for (j = 0; j < cpdim; ++j) { 1534 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p]; 1535 } 1536 } 1537 } 1538 ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1539 } 1540 offsetJ += cpdim*NcJ; 1541 } 1542 offsetI += fpdim*Nc; 1543 ierr = PetscFree(points);CHKERRQ(ierr); 1544 } 1545 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 1546 for (c = cStart; c < cEnd; ++c) { 1547 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1548 } 1549 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1550 ierr = PetscFree(feRef);CHKERRQ(ierr); 1551 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1552 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1553 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1554 if (mesh->printFEM) { 1555 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1556 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1557 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1558 } 1559 #if 0 1560 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1561 #endif 1562 PetscFunctionReturn(0); 1563 } 1564 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "DMPlexComputeInjectorFEM" 1567 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 1568 { 1569 PetscDS prob; 1570 PetscFE *feRef; 1571 Vec fv, cv; 1572 IS fis, cis; 1573 PetscSection fsection, fglobalSection, csection, cglobalSection; 1574 PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 1575 PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m; 1576 PetscErrorCode ierr; 1577 1578 PetscFunctionBegin; 1579 ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1580 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1581 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1582 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1583 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1584 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1585 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1586 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 1587 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1588 for (f = 0; f < Nf; ++f) { 1589 PetscFE fe; 1590 PetscInt fNb, Nc; 1591 1592 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1593 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1594 ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 1595 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1596 fTotDim += fNb*Nc; 1597 } 1598 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1599 ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 1600 for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 1601 PetscFE feC; 1602 PetscDualSpace QF, QC; 1603 PetscInt NcF, NcC, fpdim, cpdim; 1604 1605 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 1606 ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 1607 ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 1608 if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC); 1609 ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 1610 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 1611 ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 1612 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 1613 for (c = 0; c < cpdim; ++c) { 1614 PetscQuadrature cfunc; 1615 const PetscReal *cqpoints; 1616 PetscInt NpC; 1617 1618 ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 1619 ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 1620 if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 1621 for (f = 0; f < fpdim; ++f) { 1622 PetscQuadrature ffunc; 1623 const PetscReal *fqpoints; 1624 PetscReal sum = 0.0; 1625 PetscInt NpF, comp; 1626 1627 ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 1628 ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 1629 if (NpC != NpF) continue; 1630 for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 1631 if (sum > 1.0e-9) continue; 1632 for (comp = 0; comp < NcC; ++comp) { 1633 cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 1634 } 1635 break; 1636 } 1637 } 1638 offsetC += cpdim*NcC; 1639 offsetF += fpdim*NcF; 1640 } 1641 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1642 ierr = PetscFree(feRef);CHKERRQ(ierr); 1643 1644 ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 1645 ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 1646 ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr); 1647 ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 1648 ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1649 ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1650 ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 1651 for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 1652 for (c = cStart; c < cEnd; ++c) { 1653 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 1654 for (d = 0; d < cTotDim; ++d) { 1655 if (cellCIndices[d] < 0) continue; 1656 if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]); 1657 cindices[cellCIndices[d]-startC] = cellCIndices[d]; 1658 findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 1659 } 1660 } 1661 ierr = PetscFree(cmap);CHKERRQ(ierr); 1662 ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 1663 1664 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 1665 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 1666 ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 1667 ierr = ISDestroy(&cis);CHKERRQ(ierr); 1668 ierr = ISDestroy(&fis);CHKERRQ(ierr); 1669 ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 1670 ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 1671 PetscFunctionReturn(0); 1672 } 1673 1674 #undef __FUNCT__ 1675 #define __FUNCT__ "BoundaryDuplicate" 1676 static PetscErrorCode BoundaryDuplicate(DMBoundary bd, DMBoundary *boundary) 1677 { 1678 DMBoundary b = bd, b2, bold = NULL; 1679 PetscErrorCode ierr; 1680 1681 PetscFunctionBegin; 1682 *boundary = NULL; 1683 for (; b; b = b->next, bold = b2) { 1684 ierr = PetscNew(&b2);CHKERRQ(ierr); 1685 ierr = PetscStrallocpy(b->name, (char **) &b2->name);CHKERRQ(ierr); 1686 ierr = PetscStrallocpy(b->labelname, (char **) &b2->labelname);CHKERRQ(ierr); 1687 ierr = PetscMalloc1(b->numids, &b2->ids);CHKERRQ(ierr); 1688 ierr = PetscMemcpy(b2->ids, b->ids, b->numids*sizeof(PetscInt));CHKERRQ(ierr); 1689 b2->label = NULL; 1690 b2->essential = b->essential; 1691 b2->field = b->field; 1692 b2->func = b->func; 1693 b2->numids = b->numids; 1694 b2->ctx = b->ctx; 1695 b2->next = NULL; 1696 if (!*boundary) *boundary = b2; 1697 if (bold) bold->next = b2; 1698 } 1699 PetscFunctionReturn(0); 1700 } 1701 1702 #undef __FUNCT__ 1703 #define __FUNCT__ "DMPlexCopyBoundary" 1704 PetscErrorCode DMPlexCopyBoundary(DM dm, DM dmNew) 1705 { 1706 DM_Plex *mesh = (DM_Plex *) dm->data; 1707 DM_Plex *meshNew = (DM_Plex *) dmNew->data; 1708 DMBoundary b; 1709 PetscErrorCode ierr; 1710 1711 PetscFunctionBegin; 1712 ierr = BoundaryDuplicate(mesh->boundary, &meshNew->boundary);CHKERRQ(ierr); 1713 for (b = meshNew->boundary; b; b = b->next) { 1714 if (b->labelname) { 1715 ierr = DMPlexGetLabel(dmNew, b->labelname, &b->label);CHKERRQ(ierr); 1716 if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname); 1717 } 1718 } 1719 PetscFunctionReturn(0); 1720 } 1721 1722 #undef __FUNCT__ 1723 #define __FUNCT__ "DMPlexAddBoundary" 1724 /* The ids can be overridden by the command line option -bc_<boundary name> */ 1725 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], const char labelname[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx) 1726 { 1727 DM_Plex *mesh = (DM_Plex *) dm->data; 1728 DMBoundary b; 1729 PetscErrorCode ierr; 1730 1731 PetscFunctionBegin; 1732 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1733 ierr = PetscNew(&b);CHKERRQ(ierr); 1734 ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 1735 ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr); 1736 ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 1737 ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr); 1738 if (b->labelname) { 1739 ierr = DMPlexGetLabel(dm, b->labelname, &b->label);CHKERRQ(ierr); 1740 if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname); 1741 } 1742 b->essential = isEssential; 1743 b->field = field; 1744 b->func = bcFunc; 1745 b->numids = numids; 1746 b->ctx = ctx; 1747 b->next = mesh->boundary; 1748 mesh->boundary = b; 1749 PetscFunctionReturn(0); 1750 } 1751 1752 #undef __FUNCT__ 1753 #define __FUNCT__ "DMPlexGetNumBoundary" 1754 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd) 1755 { 1756 DM_Plex *mesh = (DM_Plex *) dm->data; 1757 DMBoundary b = mesh->boundary; 1758 1759 PetscFunctionBegin; 1760 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1761 PetscValidPointer(numBd, 2); 1762 *numBd = 0; 1763 while (b) {++(*numBd); b = b->next;} 1764 PetscFunctionReturn(0); 1765 } 1766 1767 #undef __FUNCT__ 1768 #define __FUNCT__ "DMPlexGetBoundary" 1769 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, const char **labelname, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx) 1770 { 1771 DM_Plex *mesh = (DM_Plex *) dm->data; 1772 DMBoundary b = mesh->boundary; 1773 PetscInt n = 0; 1774 1775 PetscFunctionBegin; 1776 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1777 while (b) { 1778 if (n == bd) break; 1779 b = b->next; 1780 ++n; 1781 } 1782 if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 1783 if (isEssential) { 1784 PetscValidPointer(isEssential, 3); 1785 *isEssential = b->essential; 1786 } 1787 if (name) { 1788 PetscValidPointer(name, 4); 1789 *name = b->name; 1790 } 1791 if (labelname) { 1792 PetscValidPointer(labelname, 5); 1793 *labelname = b->labelname; 1794 } 1795 if (field) { 1796 PetscValidPointer(field, 6); 1797 *field = b->field; 1798 } 1799 if (func) { 1800 PetscValidPointer(func, 7); 1801 *func = b->func; 1802 } 1803 if (numids) { 1804 PetscValidPointer(numids, 8); 1805 *numids = b->numids; 1806 } 1807 if (ids) { 1808 PetscValidPointer(ids, 9); 1809 *ids = b->ids; 1810 } 1811 if (ctx) { 1812 PetscValidPointer(ctx, 10); 1813 *ctx = b->ctx; 1814 } 1815 PetscFunctionReturn(0); 1816 } 1817 1818 #undef __FUNCT__ 1819 #define __FUNCT__ "DMPlexIsBoundaryPoint" 1820 PetscErrorCode DMPlexIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd) 1821 { 1822 DM_Plex *mesh = (DM_Plex *) dm->data; 1823 DMBoundary b = mesh->boundary; 1824 PetscErrorCode ierr; 1825 1826 PetscFunctionBegin; 1827 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1828 PetscValidPointer(isBd, 3); 1829 *isBd = PETSC_FALSE; 1830 while (b && !(*isBd)) { 1831 if (b->label) { 1832 PetscInt i; 1833 1834 for (i = 0; i < b->numids && !(*isBd); ++i) { 1835 ierr = DMLabelStratumHasPoint(b->label, b->ids[i], point, isBd);CHKERRQ(ierr); 1836 } 1837 } 1838 b = b->next; 1839 } 1840 PetscFunctionReturn(0); 1841 } 1842