1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3 #include <petscds.h> 4 #include <petscblaslapack.h> 5 #include <petsc/private/petscimpl.h> 6 #include <petsc/private/petscfeimpl.h> 7 8 /************************** Interpolation *******************************/ 9 10 static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) 11 { 12 PetscBool isPlex; 13 PetscErrorCode ierr; 14 15 PetscFunctionBegin; 16 ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr); 17 if (isPlex) { 18 *plex = dm; 19 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 20 } else { 21 ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr); 22 if (!*plex) { 23 ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr); 24 ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr); 25 if (copy) { 26 PetscInt i; 27 PetscObject obj; 28 const char *comps[3] = {"A","dmAux","dmCh"}; 29 30 ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr); 31 for (i = 0; i < 3; i++) { 32 ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr); 33 ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr); 34 } 35 } 36 } else { 37 ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr); 38 } 39 } 40 PetscFunctionReturn(0); 41 } 42 43 /*@C 44 DMInterpolationCreate - Creates a DMInterpolationInfo context 45 46 Collective 47 48 Input Parameter: 49 . comm - the communicator 50 51 Output Parameter: 52 . ctx - the context 53 54 Level: beginner 55 56 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy() 57 @*/ 58 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 59 { 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 63 PetscValidPointer(ctx, 2); 64 ierr = PetscNew(ctx);CHKERRQ(ierr); 65 66 (*ctx)->comm = comm; 67 (*ctx)->dim = -1; 68 (*ctx)->nInput = 0; 69 (*ctx)->points = NULL; 70 (*ctx)->cells = NULL; 71 (*ctx)->n = -1; 72 (*ctx)->coords = NULL; 73 PetscFunctionReturn(0); 74 } 75 76 /*@C 77 DMInterpolationSetDim - Sets the spatial dimension for the interpolation context 78 79 Not collective 80 81 Input Parameters: 82 + ctx - the context 83 - dim - the spatial dimension 84 85 Level: intermediate 86 87 .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints() 88 @*/ 89 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 90 { 91 PetscFunctionBegin; 92 if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim); 93 ctx->dim = dim; 94 PetscFunctionReturn(0); 95 } 96 97 /*@C 98 DMInterpolationGetDim - Gets the spatial dimension for the interpolation context 99 100 Not collective 101 102 Input Parameter: 103 . ctx - the context 104 105 Output Parameter: 106 . dim - the spatial dimension 107 108 Level: intermediate 109 110 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints() 111 @*/ 112 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 113 { 114 PetscFunctionBegin; 115 PetscValidIntPointer(dim, 2); 116 *dim = ctx->dim; 117 PetscFunctionReturn(0); 118 } 119 120 /*@C 121 DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context 122 123 Not collective 124 125 Input Parameters: 126 + ctx - the context 127 - dof - the number of fields 128 129 Level: intermediate 130 131 .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints() 132 @*/ 133 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 134 { 135 PetscFunctionBegin; 136 if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof); 137 ctx->dof = dof; 138 PetscFunctionReturn(0); 139 } 140 141 /*@C 142 DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context 143 144 Not collective 145 146 Input Parameter: 147 . ctx - the context 148 149 Output Parameter: 150 . dof - the number of fields 151 152 Level: intermediate 153 154 .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints() 155 @*/ 156 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 157 { 158 PetscFunctionBegin; 159 PetscValidIntPointer(dof, 2); 160 *dof = ctx->dof; 161 PetscFunctionReturn(0); 162 } 163 164 /*@C 165 DMInterpolationAddPoints - Add points at which we will interpolate the fields 166 167 Not collective 168 169 Input Parameters: 170 + ctx - the context 171 . n - the number of points 172 - points - the coordinates for each point, an array of size n * dim 173 174 Note: The coordinate information is copied. 175 176 Level: intermediate 177 178 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate() 179 @*/ 180 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 181 { 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 186 if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 187 ctx->nInput = n; 188 189 ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr); 190 ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 194 /*@C 195 DMInterpolationSetUp - Computea spatial indices that add in point location during interpolation 196 197 Collective on ctx 198 199 Input Parameters: 200 + ctx - the context 201 . dm - the DM for the function space used for interpolation 202 - redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes. 203 204 Level: intermediate 205 206 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate() 207 @*/ 208 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) 209 { 210 MPI_Comm comm = ctx->comm; 211 PetscScalar *a; 212 PetscInt p, q, i; 213 PetscMPIInt rank, size; 214 PetscErrorCode ierr; 215 Vec pointVec; 216 PetscSF cellSF; 217 PetscLayout layout; 218 PetscReal *globalPoints; 219 PetscScalar *globalPointsScalar; 220 const PetscInt *ranges; 221 PetscMPIInt *counts, *displs; 222 const PetscSFNode *foundCells; 223 const PetscInt *foundPoints; 224 PetscMPIInt *foundProcs, *globalProcs; 225 PetscInt n, N, numFound; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 229 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 230 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 231 if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 232 /* Locate points */ 233 n = ctx->nInput; 234 if (!redundantPoints) { 235 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 236 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 237 ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 238 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 239 ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 240 /* Communicate all points to all processes */ 241 ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr); 242 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 243 for (p = 0; p < size; ++p) { 244 counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 245 displs[p] = ranges[p]*ctx->dim; 246 } 247 ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 248 } else { 249 N = n; 250 globalPoints = ctx->points; 251 counts = displs = NULL; 252 layout = NULL; 253 } 254 #if 0 255 ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 256 /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 257 #else 258 #if defined(PETSC_USE_COMPLEX) 259 ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr); 260 for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i]; 261 #else 262 globalPointsScalar = globalPoints; 263 #endif 264 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr); 265 ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 266 for (p = 0; p < N; ++p) {foundProcs[p] = size;} 267 cellSF = NULL; 268 ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr); 269 ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr); 270 #endif 271 for (p = 0; p < numFound; ++p) { 272 if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank; 273 } 274 /* Let the lowest rank process own each point */ 275 ierr = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 276 ctx->n = 0; 277 for (p = 0; p < N; ++p) { 278 if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, (double)globalPoints[p*ctx->dim+0], (double)(ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0), (double)(ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0)); 279 else if (globalProcs[p] == rank) ctx->n++; 280 } 281 /* Create coordinates vector and array of owned cells */ 282 ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr); 283 ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 284 ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 285 ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 286 ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr); 287 ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 288 for (p = 0, q = 0, i = 0; p < N; ++p) { 289 if (globalProcs[p] == rank) { 290 PetscInt d; 291 292 for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d]; 293 ctx->cells[q] = foundCells[q].index; 294 ++q; 295 } 296 } 297 ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 298 #if 0 299 ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 300 #else 301 ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 302 ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr); 303 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 304 #endif 305 if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);} 306 if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);} 307 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 308 PetscFunctionReturn(0); 309 } 310 311 /*@C 312 DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point 313 314 Collective on ctx 315 316 Input Parameter: 317 . ctx - the context 318 319 Output Parameter: 320 . coordinates - the coordinates of interpolation points 321 322 Note: The local vector entries correspond to interpolation points lying on this process, according to the associated DM. This is a borrowed vector that the user should not destroy. 323 324 Level: intermediate 325 326 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate() 327 @*/ 328 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 329 { 330 PetscFunctionBegin; 331 PetscValidPointer(coordinates, 2); 332 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 333 *coordinates = ctx->coords; 334 PetscFunctionReturn(0); 335 } 336 337 /*@C 338 DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values 339 340 Collective on ctx 341 342 Input Parameter: 343 . ctx - the context 344 345 Output Parameter: 346 . v - a vector capable of holding the interpolated field values 347 348 Note: This vector should be returned using DMInterpolationRestoreVector(). 349 350 Level: intermediate 351 352 .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate() 353 @*/ 354 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 355 { 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 PetscValidPointer(v, 2); 360 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 361 ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 362 ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 363 ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 364 ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368 /*@C 369 DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values 370 371 Collective on ctx 372 373 Input Parameters: 374 + ctx - the context 375 - v - a vector capable of holding the interpolated field values 376 377 Level: intermediate 378 379 .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate() 380 @*/ 381 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 382 { 383 PetscErrorCode ierr; 384 385 PetscFunctionBegin; 386 PetscValidPointer(v, 2); 387 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 388 ierr = VecDestroy(v);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 392 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 393 { 394 PetscReal *v0, *J, *invJ, detJ; 395 const PetscScalar *coords; 396 PetscScalar *a; 397 PetscInt p; 398 PetscErrorCode ierr; 399 400 PetscFunctionBegin; 401 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 402 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 403 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 404 for (p = 0; p < ctx->n; ++p) { 405 PetscInt c = ctx->cells[p]; 406 PetscScalar *x = NULL; 407 PetscReal xi[4]; 408 PetscInt d, f, comp; 409 410 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 411 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c); 412 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 413 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 414 415 for (d = 0; d < ctx->dim; ++d) { 416 xi[d] = 0.0; 417 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 418 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d]; 419 } 420 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 421 } 422 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 423 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 424 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 425 PetscFunctionReturn(0); 426 } 427 428 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 429 { 430 PetscReal *v0, *J, *invJ, detJ; 431 const PetscScalar *coords; 432 PetscScalar *a; 433 PetscInt p; 434 PetscErrorCode ierr; 435 436 PetscFunctionBegin; 437 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 438 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 439 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 440 for (p = 0; p < ctx->n; ++p) { 441 PetscInt c = ctx->cells[p]; 442 const PetscInt order[3] = {2, 1, 3}; 443 PetscScalar *x = NULL; 444 PetscReal xi[4]; 445 PetscInt d, f, comp; 446 447 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 448 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c); 449 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 450 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 451 452 for (d = 0; d < ctx->dim; ++d) { 453 xi[d] = 0.0; 454 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 455 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[order[d]*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d]; 456 } 457 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 458 } 459 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 460 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 461 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 462 PetscFunctionReturn(0); 463 } 464 465 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 466 { 467 const PetscScalar *vertices = (const PetscScalar*) ctx; 468 const PetscScalar x0 = vertices[0]; 469 const PetscScalar y0 = vertices[1]; 470 const PetscScalar x1 = vertices[2]; 471 const PetscScalar y1 = vertices[3]; 472 const PetscScalar x2 = vertices[4]; 473 const PetscScalar y2 = vertices[5]; 474 const PetscScalar x3 = vertices[6]; 475 const PetscScalar y3 = vertices[7]; 476 const PetscScalar f_1 = x1 - x0; 477 const PetscScalar g_1 = y1 - y0; 478 const PetscScalar f_3 = x3 - x0; 479 const PetscScalar g_3 = y3 - y0; 480 const PetscScalar f_01 = x2 - x1 - x3 + x0; 481 const PetscScalar g_01 = y2 - y1 - y3 + y0; 482 const PetscScalar *ref; 483 PetscScalar *real; 484 PetscErrorCode ierr; 485 486 PetscFunctionBegin; 487 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 488 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 489 { 490 const PetscScalar p0 = ref[0]; 491 const PetscScalar p1 = ref[1]; 492 493 real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 494 real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 495 } 496 ierr = PetscLogFlops(28);CHKERRQ(ierr); 497 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 498 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 499 PetscFunctionReturn(0); 500 } 501 502 #include <petsc/private/dmimpl.h> 503 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 504 { 505 const PetscScalar *vertices = (const PetscScalar*) ctx; 506 const PetscScalar x0 = vertices[0]; 507 const PetscScalar y0 = vertices[1]; 508 const PetscScalar x1 = vertices[2]; 509 const PetscScalar y1 = vertices[3]; 510 const PetscScalar x2 = vertices[4]; 511 const PetscScalar y2 = vertices[5]; 512 const PetscScalar x3 = vertices[6]; 513 const PetscScalar y3 = vertices[7]; 514 const PetscScalar f_01 = x2 - x1 - x3 + x0; 515 const PetscScalar g_01 = y2 - y1 - y3 + y0; 516 const PetscScalar *ref; 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 521 { 522 const PetscScalar x = ref[0]; 523 const PetscScalar y = ref[1]; 524 const PetscInt rows[2] = {0, 1}; 525 PetscScalar values[4]; 526 527 values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5; 528 values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5; 529 ierr = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 530 } 531 ierr = PetscLogFlops(30);CHKERRQ(ierr); 532 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 533 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 534 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 535 PetscFunctionReturn(0); 536 } 537 538 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 539 { 540 DM dmCoord; 541 PetscFE fem = NULL; 542 SNES snes; 543 KSP ksp; 544 PC pc; 545 Vec coordsLocal, r, ref, real; 546 Mat J; 547 PetscTabulation T; 548 const PetscScalar *coords; 549 PetscScalar *a; 550 PetscReal xir[2]; 551 PetscInt Nf, p; 552 const PetscInt dof = ctx->dof; 553 PetscErrorCode ierr; 554 555 PetscFunctionBegin; 556 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 557 if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);} 558 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 559 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 560 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 561 ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 562 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 563 ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 564 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 565 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 566 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 567 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 568 ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 569 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 570 ierr = MatSetUp(J);CHKERRQ(ierr); 571 ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr); 572 ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr); 573 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 574 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 575 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 576 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 577 578 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 579 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 580 ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr); 581 for (p = 0; p < ctx->n; ++p) { 582 PetscScalar *x = NULL, *vertices = NULL; 583 PetscScalar *xi; 584 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 585 586 /* Can make this do all points at once */ 587 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 588 if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2); 589 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 590 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 591 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 592 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 593 xi[0] = coords[p*ctx->dim+0]; 594 xi[1] = coords[p*ctx->dim+1]; 595 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 596 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 597 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 598 xir[0] = PetscRealPart(xi[0]); 599 xir[1] = PetscRealPart(xi[1]); 600 if (4*dof != xSize) { 601 PetscInt d; 602 603 xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0; 604 ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr); 605 for (comp = 0; comp < dof; ++comp) { 606 a[p*dof+comp] = 0.0; 607 for (d = 0; d < xSize/dof; ++d) { 608 a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp]; 609 } 610 } 611 } else { 612 for (comp = 0; comp < dof; ++comp) 613 a[p*dof+comp] = x[0*dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*dof+comp]*xir[0]*(1 - xir[1]) + x[2*dof+comp]*xir[0]*xir[1] + x[3*dof+comp]*(1 - xir[0])*xir[1]; 614 } 615 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 616 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 617 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 618 } 619 ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr); 620 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 621 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 622 623 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 624 ierr = VecDestroy(&r);CHKERRQ(ierr); 625 ierr = VecDestroy(&ref);CHKERRQ(ierr); 626 ierr = VecDestroy(&real);CHKERRQ(ierr); 627 ierr = MatDestroy(&J);CHKERRQ(ierr); 628 PetscFunctionReturn(0); 629 } 630 631 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 632 { 633 const PetscScalar *vertices = (const PetscScalar*) ctx; 634 const PetscScalar x0 = vertices[0]; 635 const PetscScalar y0 = vertices[1]; 636 const PetscScalar z0 = vertices[2]; 637 const PetscScalar x1 = vertices[9]; 638 const PetscScalar y1 = vertices[10]; 639 const PetscScalar z1 = vertices[11]; 640 const PetscScalar x2 = vertices[6]; 641 const PetscScalar y2 = vertices[7]; 642 const PetscScalar z2 = vertices[8]; 643 const PetscScalar x3 = vertices[3]; 644 const PetscScalar y3 = vertices[4]; 645 const PetscScalar z3 = vertices[5]; 646 const PetscScalar x4 = vertices[12]; 647 const PetscScalar y4 = vertices[13]; 648 const PetscScalar z4 = vertices[14]; 649 const PetscScalar x5 = vertices[15]; 650 const PetscScalar y5 = vertices[16]; 651 const PetscScalar z5 = vertices[17]; 652 const PetscScalar x6 = vertices[18]; 653 const PetscScalar y6 = vertices[19]; 654 const PetscScalar z6 = vertices[20]; 655 const PetscScalar x7 = vertices[21]; 656 const PetscScalar y7 = vertices[22]; 657 const PetscScalar z7 = vertices[23]; 658 const PetscScalar f_1 = x1 - x0; 659 const PetscScalar g_1 = y1 - y0; 660 const PetscScalar h_1 = z1 - z0; 661 const PetscScalar f_3 = x3 - x0; 662 const PetscScalar g_3 = y3 - y0; 663 const PetscScalar h_3 = z3 - z0; 664 const PetscScalar f_4 = x4 - x0; 665 const PetscScalar g_4 = y4 - y0; 666 const PetscScalar h_4 = z4 - z0; 667 const PetscScalar f_01 = x2 - x1 - x3 + x0; 668 const PetscScalar g_01 = y2 - y1 - y3 + y0; 669 const PetscScalar h_01 = z2 - z1 - z3 + z0; 670 const PetscScalar f_12 = x7 - x3 - x4 + x0; 671 const PetscScalar g_12 = y7 - y3 - y4 + y0; 672 const PetscScalar h_12 = z7 - z3 - z4 + z0; 673 const PetscScalar f_02 = x5 - x1 - x4 + x0; 674 const PetscScalar g_02 = y5 - y1 - y4 + y0; 675 const PetscScalar h_02 = z5 - z1 - z4 + z0; 676 const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 677 const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 678 const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 679 const PetscScalar *ref; 680 PetscScalar *real; 681 PetscErrorCode ierr; 682 683 PetscFunctionBegin; 684 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 685 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 686 { 687 const PetscScalar p0 = ref[0]; 688 const PetscScalar p1 = ref[1]; 689 const PetscScalar p2 = ref[2]; 690 691 real[0] = x0 + f_1*p0 + f_3*p1 + f_4*p2 + f_01*p0*p1 + f_12*p1*p2 + f_02*p0*p2 + f_012*p0*p1*p2; 692 real[1] = y0 + g_1*p0 + g_3*p1 + g_4*p2 + g_01*p0*p1 + g_01*p0*p1 + g_12*p1*p2 + g_02*p0*p2 + g_012*p0*p1*p2; 693 real[2] = z0 + h_1*p0 + h_3*p1 + h_4*p2 + h_01*p0*p1 + h_01*p0*p1 + h_12*p1*p2 + h_02*p0*p2 + h_012*p0*p1*p2; 694 } 695 ierr = PetscLogFlops(114);CHKERRQ(ierr); 696 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 697 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 702 { 703 const PetscScalar *vertices = (const PetscScalar*) ctx; 704 const PetscScalar x0 = vertices[0]; 705 const PetscScalar y0 = vertices[1]; 706 const PetscScalar z0 = vertices[2]; 707 const PetscScalar x1 = vertices[9]; 708 const PetscScalar y1 = vertices[10]; 709 const PetscScalar z1 = vertices[11]; 710 const PetscScalar x2 = vertices[6]; 711 const PetscScalar y2 = vertices[7]; 712 const PetscScalar z2 = vertices[8]; 713 const PetscScalar x3 = vertices[3]; 714 const PetscScalar y3 = vertices[4]; 715 const PetscScalar z3 = vertices[5]; 716 const PetscScalar x4 = vertices[12]; 717 const PetscScalar y4 = vertices[13]; 718 const PetscScalar z4 = vertices[14]; 719 const PetscScalar x5 = vertices[15]; 720 const PetscScalar y5 = vertices[16]; 721 const PetscScalar z5 = vertices[17]; 722 const PetscScalar x6 = vertices[18]; 723 const PetscScalar y6 = vertices[19]; 724 const PetscScalar z6 = vertices[20]; 725 const PetscScalar x7 = vertices[21]; 726 const PetscScalar y7 = vertices[22]; 727 const PetscScalar z7 = vertices[23]; 728 const PetscScalar f_xy = x2 - x1 - x3 + x0; 729 const PetscScalar g_xy = y2 - y1 - y3 + y0; 730 const PetscScalar h_xy = z2 - z1 - z3 + z0; 731 const PetscScalar f_yz = x7 - x3 - x4 + x0; 732 const PetscScalar g_yz = y7 - y3 - y4 + y0; 733 const PetscScalar h_yz = z7 - z3 - z4 + z0; 734 const PetscScalar f_xz = x5 - x1 - x4 + x0; 735 const PetscScalar g_xz = y5 - y1 - y4 + y0; 736 const PetscScalar h_xz = z5 - z1 - z4 + z0; 737 const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 738 const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 739 const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 740 const PetscScalar *ref; 741 PetscErrorCode ierr; 742 743 PetscFunctionBegin; 744 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 745 { 746 const PetscScalar x = ref[0]; 747 const PetscScalar y = ref[1]; 748 const PetscScalar z = ref[2]; 749 const PetscInt rows[3] = {0, 1, 2}; 750 PetscScalar values[9]; 751 752 values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0; 753 values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0; 754 values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0; 755 values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0; 756 values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0; 757 values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0; 758 values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0; 759 values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0; 760 values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0; 761 762 ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 763 } 764 ierr = PetscLogFlops(152);CHKERRQ(ierr); 765 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 766 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 767 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770 771 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 772 { 773 DM dmCoord; 774 SNES snes; 775 KSP ksp; 776 PC pc; 777 Vec coordsLocal, r, ref, real; 778 Mat J; 779 const PetscScalar *coords; 780 PetscScalar *a; 781 PetscInt p; 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 786 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 787 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 788 ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 789 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 790 ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 791 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 792 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 793 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 794 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 795 ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 796 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 797 ierr = MatSetUp(J);CHKERRQ(ierr); 798 ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr); 799 ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr); 800 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 801 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 802 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 803 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 804 805 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 806 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 807 for (p = 0; p < ctx->n; ++p) { 808 PetscScalar *x = NULL, *vertices = NULL; 809 PetscScalar *xi; 810 PetscReal xir[3]; 811 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 812 813 /* Can make this do all points at once */ 814 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 815 if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3); 816 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 817 if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof); 818 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 819 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 820 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 821 xi[0] = coords[p*ctx->dim+0]; 822 xi[1] = coords[p*ctx->dim+1]; 823 xi[2] = coords[p*ctx->dim+2]; 824 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 825 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 826 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 827 xir[0] = PetscRealPart(xi[0]); 828 xir[1] = PetscRealPart(xi[1]); 829 xir[2] = PetscRealPart(xi[2]); 830 for (comp = 0; comp < ctx->dof; ++comp) { 831 a[p*ctx->dof+comp] = 832 x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) + 833 x[3*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) + 834 x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) + 835 x[1*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) + 836 x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] + 837 x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] + 838 x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] + 839 x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2]; 840 } 841 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 842 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 843 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 844 } 845 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 846 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 847 848 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 849 ierr = VecDestroy(&r);CHKERRQ(ierr); 850 ierr = VecDestroy(&ref);CHKERRQ(ierr); 851 ierr = VecDestroy(&real);CHKERRQ(ierr); 852 ierr = MatDestroy(&J);CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855 856 /*@C 857 DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points. 858 859 Input Parameters: 860 + ctx - The DMInterpolationInfo context 861 . dm - The DM 862 - x - The local vector containing the field to be interpolated 863 864 Output Parameters: 865 . v - The vector containing the interpolated values 866 867 Note: A suitable v can be obtained using DMInterpolationGetVector(). 868 869 Level: beginner 870 871 .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate() 872 @*/ 873 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 874 { 875 PetscInt dim, coneSize, n; 876 PetscErrorCode ierr; 877 878 PetscFunctionBegin; 879 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 880 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 881 PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 882 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 883 if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %D should be %D", n, ctx->n*ctx->dof); 884 if (n) { 885 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 886 ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 887 if (dim == 2) { 888 if (coneSize == 3) { 889 ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr); 890 } else if (coneSize == 4) { 891 ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 892 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim); 893 } else if (dim == 3) { 894 if (coneSize == 4) { 895 ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr); 896 } else { 897 ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 898 } 899 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim); 900 } 901 PetscFunctionReturn(0); 902 } 903 904 /*@C 905 DMInterpolationDestroy - Destroys a DMInterpolationInfo context 906 907 Collective on ctx 908 909 Input Parameter: 910 . ctx - the context 911 912 Level: beginner 913 914 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate() 915 @*/ 916 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 917 { 918 PetscErrorCode ierr; 919 920 PetscFunctionBegin; 921 PetscValidPointer(ctx, 2); 922 ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 923 ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 924 ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 925 ierr = PetscFree(*ctx);CHKERRQ(ierr); 926 *ctx = NULL; 927 PetscFunctionReturn(0); 928 } 929 930 /*@C 931 SNESMonitorFields - Monitors the residual for each field separately 932 933 Collective on SNES 934 935 Input Parameters: 936 + snes - the SNES context 937 . its - iteration number 938 . fgnorm - 2-norm of residual 939 - vf - PetscViewerAndFormat of type ASCII 940 941 Notes: 942 This routine prints the residual norm at each iteration. 943 944 Level: intermediate 945 946 .seealso: SNESMonitorSet(), SNESMonitorDefault() 947 @*/ 948 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 949 { 950 PetscViewer viewer = vf->viewer; 951 Vec res; 952 DM dm; 953 PetscSection s; 954 const PetscScalar *r; 955 PetscReal *lnorms, *norms; 956 PetscInt numFields, f, pStart, pEnd, p; 957 PetscErrorCode ierr; 958 959 PetscFunctionBegin; 960 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 961 ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr); 962 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 963 ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 964 ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr); 965 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 966 ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr); 967 ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr); 968 for (p = pStart; p < pEnd; ++p) { 969 for (f = 0; f < numFields; ++f) { 970 PetscInt fdof, foff, d; 971 972 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 973 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 974 for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d])); 975 } 976 } 977 ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr); 978 ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 979 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 980 ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 981 ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr); 982 for (f = 0; f < numFields; ++f) { 983 if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 984 ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr); 985 } 986 ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr); 987 ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 988 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 989 ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr); 990 PetscFunctionReturn(0); 991 } 992 993 /********************* Residual Computation **************************/ 994 995 996 /*@ 997 DMPlexSNESGetGeometryFVM - Return precomputed geometric data 998 999 Input Parameter: 1000 . dm - The DM 1001 1002 Output Parameters: 1003 + facegeom - The values precomputed from face geometry 1004 . cellgeom - The values precomputed from cell geometry 1005 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell 1006 1007 Level: developer 1008 1009 .seealso: DMPlexTSSetRHSFunctionLocal() 1010 @*/ 1011 PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 1012 { 1013 DM plex; 1014 PetscErrorCode ierr; 1015 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1018 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 1019 ierr = DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);CHKERRQ(ierr); 1020 if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);} 1021 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1022 PetscFunctionReturn(0); 1023 } 1024 1025 /*@ 1026 DMPlexSNESGetGradientDM - Return gradient data layout 1027 1028 Input Parameters: 1029 + dm - The DM 1030 - fv - The PetscFV 1031 1032 Output Parameter: 1033 . dmGrad - The layout for gradient values 1034 1035 Level: developer 1036 1037 .seealso: DMPlexSNESGetGeometryFVM() 1038 @*/ 1039 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad) 1040 { 1041 DM plex; 1042 PetscBool computeGradients; 1043 PetscErrorCode ierr; 1044 1045 PetscFunctionBegin; 1046 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1047 PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2); 1048 PetscValidPointer(dmGrad,3); 1049 ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr); 1050 if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);} 1051 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 1052 ierr = DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);CHKERRQ(ierr); 1053 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1054 PetscFunctionReturn(0); 1055 } 1056 1057 static PetscErrorCode DMPlexComputeBdResidual_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF, DMField coordField, IS facetIS) 1058 { 1059 DM_Plex *mesh = (DM_Plex *) dm->data; 1060 DM plex = NULL, plexA = NULL; 1061 DMEnclosureType encAux; 1062 PetscDS prob, probAux = NULL; 1063 PetscSection section, sectionAux = NULL; 1064 Vec locA = NULL; 1065 PetscScalar *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL; 1066 PetscInt v; 1067 PetscInt totDim, totDimAux = 0; 1068 PetscErrorCode ierr; 1069 1070 PetscFunctionBegin; 1071 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1072 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 1073 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1074 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1075 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1076 if (locA) { 1077 DM dmAux; 1078 1079 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1080 ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 1081 ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr); 1082 ierr = DMGetDS(plexA, &probAux);CHKERRQ(ierr); 1083 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1084 ierr = DMGetLocalSection(plexA, §ionAux);CHKERRQ(ierr); 1085 } 1086 for (v = 0; v < numValues; ++v) { 1087 PetscFEGeom *fgeom; 1088 PetscInt maxDegree; 1089 PetscQuadrature qGeom = NULL; 1090 IS pointIS; 1091 const PetscInt *points; 1092 PetscInt numFaces, face, Nq; 1093 1094 ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr); 1095 if (!pointIS) continue; /* No points with that id on this process */ 1096 { 1097 IS isectIS; 1098 1099 /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */ 1100 ierr = ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);CHKERRQ(ierr); 1101 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1102 pointIS = isectIS; 1103 } 1104 ierr = ISGetLocalSize(pointIS,&numFaces);CHKERRQ(ierr); 1105 ierr = ISGetIndices(pointIS,&points);CHKERRQ(ierr); 1106 ierr = PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim, &elemVec, locA ? numFaces*totDimAux : 0, &a);CHKERRQ(ierr); 1107 ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 1108 if (maxDegree <= 1) { 1109 ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);CHKERRQ(ierr); 1110 } 1111 if (!qGeom) { 1112 PetscFE fe; 1113 1114 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 1115 ierr = PetscFEGetFaceQuadrature(fe, &qGeom);CHKERRQ(ierr); 1116 ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr); 1117 } 1118 ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1119 ierr = DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr); 1120 for (face = 0; face < numFaces; ++face) { 1121 const PetscInt point = points[face], *support, *cone; 1122 PetscScalar *x = NULL; 1123 PetscInt i, coneSize, faceLoc; 1124 1125 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1126 ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 1127 ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 1128 for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break; 1129 if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]); 1130 fgeom->face[face][0] = faceLoc; 1131 ierr = DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1132 for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i]; 1133 ierr = DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1134 if (locX_t) { 1135 ierr = DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr); 1136 for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i]; 1137 ierr = DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr); 1138 } 1139 if (locA) { 1140 PetscInt subp; 1141 1142 ierr = DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp);CHKERRQ(ierr); 1143 ierr = DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1144 for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i]; 1145 ierr = DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1146 } 1147 } 1148 ierr = PetscArrayzero(elemVec, numFaces*totDim);CHKERRQ(ierr); 1149 { 1150 PetscFE fe; 1151 PetscInt Nb; 1152 PetscFEGeom *chunkGeom = NULL; 1153 /* Conforming batches */ 1154 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1155 /* Remainder */ 1156 PetscInt Nr, offset; 1157 1158 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 1159 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1160 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1161 /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */ 1162 blockSize = Nb; 1163 batchSize = numBlocks * blockSize; 1164 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1165 numChunks = numFaces / (numBatches*batchSize); 1166 Ne = numChunks*numBatches*batchSize; 1167 Nr = numFaces % (numBatches*batchSize); 1168 offset = numFaces - Nr; 1169 ierr = PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);CHKERRQ(ierr); 1170 ierr = PetscFEIntegrateBdResidual(prob, field, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr); 1171 ierr = PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);CHKERRQ(ierr); 1172 ierr = PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr); 1173 ierr = PetscFEIntegrateBdResidual(prob, field, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, &elemVec[offset*totDim]);CHKERRQ(ierr); 1174 ierr = PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr); 1175 } 1176 for (face = 0; face < numFaces; ++face) { 1177 const PetscInt point = points[face], *support; 1178 1179 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDim, &elemVec[face*totDim]);CHKERRQ(ierr);} 1180 ierr = DMPlexGetSupport(plex, point, &support);CHKERRQ(ierr); 1181 ierr = DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face*totDim], ADD_ALL_VALUES);CHKERRQ(ierr); 1182 } 1183 ierr = DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr); 1184 ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr); 1185 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1186 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1187 ierr = PetscFree4(u, u_t, elemVec, a);CHKERRQ(ierr); 1188 } 1189 if (plex) {ierr = DMDestroy(&plex);CHKERRQ(ierr);} 1190 if (plexA) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);} 1191 PetscFunctionReturn(0); 1192 } 1193 1194 PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF) 1195 { 1196 DMField coordField; 1197 DMLabel depthLabel; 1198 IS facetIS; 1199 PetscInt dim; 1200 PetscErrorCode ierr; 1201 1202 PetscFunctionBegin; 1203 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1204 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1205 ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr); 1206 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1207 ierr = DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);CHKERRQ(ierr); 1208 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 1209 PetscFunctionReturn(0); 1210 } 1211 1212 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user) 1213 { 1214 PetscDS prob; 1215 PetscInt numBd, bd; 1216 DMField coordField = NULL; 1217 IS facetIS = NULL; 1218 DMLabel depthLabel; 1219 PetscInt dim; 1220 PetscErrorCode ierr; 1221 1222 PetscFunctionBegin; 1223 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1224 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1225 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1226 ierr = DMLabelGetStratumIS(depthLabel,dim - 1,&facetIS);CHKERRQ(ierr); 1227 ierr = PetscDSGetNumBoundary(prob, &numBd);CHKERRQ(ierr); 1228 for (bd = 0; bd < numBd; ++bd) { 1229 DMBoundaryConditionType type; 1230 const char *bdLabel; 1231 DMLabel label; 1232 const PetscInt *values; 1233 PetscInt field, numValues; 1234 PetscObject obj; 1235 PetscClassId id; 1236 1237 ierr = PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1238 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 1239 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1240 if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue; 1241 if (!facetIS) { 1242 DMLabel depthLabel; 1243 PetscInt dim; 1244 1245 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1246 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1247 ierr = DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS);CHKERRQ(ierr); 1248 } 1249 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1250 ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1251 ierr = DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);CHKERRQ(ierr); 1252 } 1253 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 1254 PetscFunctionReturn(0); 1255 } 1256 1257 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user) 1258 { 1259 DM_Plex *mesh = (DM_Plex *) dm->data; 1260 const char *name = "Residual"; 1261 DM dmAux = NULL; 1262 DM dmGrad = NULL; 1263 DMLabel ghostLabel = NULL; 1264 PetscDS prob = NULL; 1265 PetscDS probAux = NULL; 1266 PetscSection section = NULL; 1267 PetscBool useFEM = PETSC_FALSE; 1268 PetscBool useFVM = PETSC_FALSE; 1269 PetscBool isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE; 1270 PetscFV fvm = NULL; 1271 PetscFVCellGeom *cgeomFVM = NULL; 1272 PetscFVFaceGeom *fgeomFVM = NULL; 1273 DMField coordField = NULL; 1274 Vec locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL; 1275 PetscScalar *u = NULL, *u_t, *a, *uL, *uR; 1276 IS chunkIS; 1277 const PetscInt *cells; 1278 PetscInt cStart, cEnd, numCells; 1279 PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd; 1280 PetscInt maxDegree = PETSC_MAX_INT; 1281 PetscQuadrature affineQuad = NULL, *quads = NULL; 1282 PetscFEGeom *affineGeom = NULL, **geoms = NULL; 1283 PetscErrorCode ierr; 1284 1285 PetscFunctionBegin; 1286 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1287 /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */ 1288 /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */ 1289 /* FEM+FVM */ 1290 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 1291 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1292 /* 1: Get sizes from dm and dmAux */ 1293 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 1294 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1295 ierr = DMGetCellDS(dm, cStart, &prob);CHKERRQ(ierr); 1296 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1297 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1298 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1299 if (locA) { 1300 PetscInt subcell; 1301 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1302 ierr = DMGetEnclosurePoint(dmAux, dm, DM_ENC_UNKNOWN, cStart, &subcell);CHKERRQ(ierr); 1303 ierr = DMGetCellDS(dmAux, subcell, &probAux);CHKERRQ(ierr); 1304 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1305 } 1306 /* 2: Get geometric data */ 1307 for (f = 0; f < Nf; ++f) { 1308 PetscObject obj; 1309 PetscClassId id; 1310 PetscBool fimp; 1311 1312 ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr); 1313 if (isImplicit != fimp) continue; 1314 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1315 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1316 if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;} 1317 if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;} 1318 } 1319 if (useFEM) { 1320 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1321 ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr); 1322 if (maxDegree <= 1) { 1323 ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);CHKERRQ(ierr); 1324 if (affineQuad) { 1325 ierr = DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr); 1326 } 1327 } else { 1328 ierr = PetscCalloc2(Nf,&quads,Nf,&geoms);CHKERRQ(ierr); 1329 for (f = 0; f < Nf; ++f) { 1330 PetscObject obj; 1331 PetscClassId id; 1332 PetscBool fimp; 1333 1334 ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr); 1335 if (isImplicit != fimp) continue; 1336 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1337 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1338 if (id == PETSCFE_CLASSID) { 1339 PetscFE fe = (PetscFE) obj; 1340 1341 ierr = PetscFEGetQuadrature(fe, &quads[f]);CHKERRQ(ierr); 1342 ierr = PetscObjectReference((PetscObject)quads[f]);CHKERRQ(ierr); 1343 ierr = DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);CHKERRQ(ierr); 1344 } 1345 } 1346 } 1347 } 1348 if (useFVM) { 1349 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr); 1350 ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 1351 ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1352 /* Reconstruct and limit cell gradients */ 1353 ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr); 1354 if (dmGrad) { 1355 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1356 ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1357 ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr); 1358 /* Communicate gradient values */ 1359 ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1360 ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1361 ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1362 ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1363 } 1364 /* Handle non-essential (e.g. outflow) boundary values */ 1365 ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr); 1366 } 1367 /* Loop over chunks */ 1368 if (useFEM) {ierr = ISCreate(PETSC_COMM_SELF, &chunkIS);CHKERRQ(ierr);} 1369 numCells = cEnd - cStart; 1370 numChunks = 1; 1371 cellChunkSize = numCells/numChunks; 1372 faceChunkSize = (fEnd - fStart)/numChunks; 1373 numChunks = PetscMin(1,numCells); 1374 for (chunk = 0; chunk < numChunks; ++chunk) { 1375 PetscScalar *elemVec, *fluxL, *fluxR; 1376 PetscReal *vol; 1377 PetscFVFaceGeom *fgeom; 1378 PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c; 1379 PetscInt fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = 0, face; 1380 1381 /* Extract field coefficients */ 1382 if (useFEM) { 1383 ierr = ISGetPointSubrange(chunkIS, cS, cE, cells);CHKERRQ(ierr); 1384 ierr = DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1385 ierr = DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr); 1386 ierr = PetscArrayzero(elemVec, numCells*totDim);CHKERRQ(ierr); 1387 } 1388 if (useFVM) { 1389 ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);CHKERRQ(ierr); 1390 ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);CHKERRQ(ierr); 1391 ierr = DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);CHKERRQ(ierr); 1392 ierr = DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);CHKERRQ(ierr); 1393 ierr = PetscArrayzero(fluxL, numFaces*totDim);CHKERRQ(ierr); 1394 ierr = PetscArrayzero(fluxR, numFaces*totDim);CHKERRQ(ierr); 1395 } 1396 /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */ 1397 /* Loop over fields */ 1398 for (f = 0; f < Nf; ++f) { 1399 PetscObject obj; 1400 PetscClassId id; 1401 PetscBool fimp; 1402 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1403 1404 ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr); 1405 if (isImplicit != fimp) continue; 1406 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1407 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1408 if (id == PETSCFE_CLASSID) { 1409 PetscFE fe = (PetscFE) obj; 1410 PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 1411 PetscFEGeom *chunkGeom = NULL; 1412 PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; 1413 PetscInt Nq, Nb; 1414 1415 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1416 ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1417 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1418 blockSize = Nb; 1419 batchSize = numBlocks * blockSize; 1420 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1421 numChunks = numCells / (numBatches*batchSize); 1422 Ne = numChunks*numBatches*batchSize; 1423 Nr = numCells % (numBatches*batchSize); 1424 offset = numCells - Nr; 1425 /* Integrate FE residual to get elemVec (need fields at quadrature points) */ 1426 /* For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */ 1427 ierr = PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);CHKERRQ(ierr); 1428 ierr = PetscFEIntegrateResidual(prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr); 1429 ierr = PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1430 ierr = PetscFEIntegrateResidual(prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);CHKERRQ(ierr); 1431 ierr = PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1432 } else if (id == PETSCFV_CLASSID) { 1433 PetscFV fv = (PetscFV) obj; 1434 1435 Ne = numFaces; 1436 /* Riemann solve over faces (need fields at face centroids) */ 1437 /* We need to evaluate FE fields at those coordinates */ 1438 ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr); 1439 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 1440 } 1441 /* Loop over domain */ 1442 if (useFEM) { 1443 /* Add elemVec to locX */ 1444 for (c = cS; c < cE; ++c) { 1445 const PetscInt cell = cells ? cells[c] : c; 1446 const PetscInt cind = c - cStart; 1447 1448 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);CHKERRQ(ierr);} 1449 if (ghostLabel) { 1450 PetscInt ghostVal; 1451 1452 ierr = DMLabelGetValue(ghostLabel,cell,&ghostVal);CHKERRQ(ierr); 1453 if (ghostVal > 0) continue; 1454 } 1455 ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);CHKERRQ(ierr); 1456 } 1457 } 1458 if (useFVM) { 1459 PetscScalar *fa; 1460 PetscInt iface; 1461 1462 ierr = VecGetArray(locF, &fa);CHKERRQ(ierr); 1463 for (f = 0; f < Nf; ++f) { 1464 PetscFV fv; 1465 PetscObject obj; 1466 PetscClassId id; 1467 PetscInt foff, pdim; 1468 1469 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1470 ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1471 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1472 if (id != PETSCFV_CLASSID) continue; 1473 fv = (PetscFV) obj; 1474 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 1475 /* Accumulate fluxes to cells */ 1476 for (face = fS, iface = 0; face < fE; ++face) { 1477 const PetscInt *scells; 1478 PetscScalar *fL = NULL, *fR = NULL; 1479 PetscInt ghost, d, nsupp, nchild; 1480 1481 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1482 ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr); 1483 ierr = DMPlexGetTreeChildren(dm, face, &nchild, NULL);CHKERRQ(ierr); 1484 if (ghost >= 0 || nsupp > 2 || nchild > 0) continue; 1485 ierr = DMPlexGetSupport(dm, face, &scells);CHKERRQ(ierr); 1486 ierr = DMLabelGetValue(ghostLabel,scells[0],&ghost);CHKERRQ(ierr); 1487 if (ghost <= 0) {ierr = DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL);CHKERRQ(ierr);} 1488 ierr = DMLabelGetValue(ghostLabel,scells[1],&ghost);CHKERRQ(ierr); 1489 if (ghost <= 0) {ierr = DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR);CHKERRQ(ierr);} 1490 for (d = 0; d < pdim; ++d) { 1491 if (fL) fL[d] -= fluxL[iface*totDim+foff+d]; 1492 if (fR) fR[d] += fluxR[iface*totDim+foff+d]; 1493 } 1494 ++iface; 1495 } 1496 } 1497 ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr); 1498 } 1499 /* Handle time derivative */ 1500 if (locX_t) { 1501 PetscScalar *x_t, *fa; 1502 1503 ierr = VecGetArray(locF, &fa);CHKERRQ(ierr); 1504 ierr = VecGetArray(locX_t, &x_t);CHKERRQ(ierr); 1505 for (f = 0; f < Nf; ++f) { 1506 PetscFV fv; 1507 PetscObject obj; 1508 PetscClassId id; 1509 PetscInt pdim, d; 1510 1511 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1512 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1513 if (id != PETSCFV_CLASSID) continue; 1514 fv = (PetscFV) obj; 1515 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 1516 for (c = cS; c < cE; ++c) { 1517 const PetscInt cell = cells ? cells[c] : c; 1518 PetscScalar *u_t, *r; 1519 1520 if (ghostLabel) { 1521 PetscInt ghostVal; 1522 1523 ierr = DMLabelGetValue(ghostLabel, cell, &ghostVal);CHKERRQ(ierr); 1524 if (ghostVal > 0) continue; 1525 } 1526 ierr = DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);CHKERRQ(ierr); 1527 ierr = DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);CHKERRQ(ierr); 1528 for (d = 0; d < pdim; ++d) r[d] += u_t[d]; 1529 } 1530 } 1531 ierr = VecRestoreArray(locX_t, &x_t);CHKERRQ(ierr); 1532 ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr); 1533 } 1534 if (useFEM) { 1535 ierr = DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1536 ierr = DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr); 1537 } 1538 if (useFVM) { 1539 ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);CHKERRQ(ierr); 1540 ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);CHKERRQ(ierr); 1541 ierr = DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);CHKERRQ(ierr); 1542 ierr = DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);CHKERRQ(ierr); 1543 if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);} 1544 } 1545 } 1546 if (useFEM) {ierr = ISDestroy(&chunkIS);CHKERRQ(ierr);} 1547 ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 1548 1549 if (useFEM) { 1550 ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);CHKERRQ(ierr); 1551 1552 if (maxDegree <= 1) { 1553 ierr = DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr); 1554 ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr); 1555 } else { 1556 for (f = 0; f < Nf; ++f) { 1557 ierr = DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);CHKERRQ(ierr); 1558 ierr = PetscQuadratureDestroy(&quads[f]);CHKERRQ(ierr); 1559 } 1560 ierr = PetscFree2(quads,geoms);CHKERRQ(ierr); 1561 } 1562 } 1563 1564 /* FEM */ 1565 /* 1: Get sizes from dm and dmAux */ 1566 /* 2: Get geometric data */ 1567 /* 3: Handle boundary values */ 1568 /* 4: Loop over domain */ 1569 /* Extract coefficients */ 1570 /* Loop over fields */ 1571 /* Set tiling for FE*/ 1572 /* Integrate FE residual to get elemVec */ 1573 /* Loop over subdomain */ 1574 /* Loop over quad points */ 1575 /* Transform coords to real space */ 1576 /* Evaluate field and aux fields at point */ 1577 /* Evaluate residual at point */ 1578 /* Transform residual to real space */ 1579 /* Add residual to elemVec */ 1580 /* Loop over domain */ 1581 /* Add elemVec to locX */ 1582 1583 /* FVM */ 1584 /* Get geometric data */ 1585 /* If using gradients */ 1586 /* Compute gradient data */ 1587 /* Loop over domain faces */ 1588 /* Count computational faces */ 1589 /* Reconstruct cell gradient */ 1590 /* Loop over domain cells */ 1591 /* Limit cell gradients */ 1592 /* Handle boundary values */ 1593 /* Loop over domain faces */ 1594 /* Read out field, centroid, normal, volume for each side of face */ 1595 /* Riemann solve over faces */ 1596 /* Loop over domain faces */ 1597 /* Accumulate fluxes to cells */ 1598 /* TODO Change printFEM to printDisc here */ 1599 if (mesh->printFEM) { 1600 Vec locFbc; 1601 PetscInt pStart, pEnd, p, maxDof; 1602 PetscScalar *zeroes; 1603 1604 ierr = VecDuplicate(locF,&locFbc);CHKERRQ(ierr); 1605 ierr = VecCopy(locF,locFbc);CHKERRQ(ierr); 1606 ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 1607 ierr = PetscSectionGetMaxDof(section,&maxDof);CHKERRQ(ierr); 1608 ierr = PetscCalloc1(maxDof,&zeroes);CHKERRQ(ierr); 1609 for (p = pStart; p < pEnd; p++) { 1610 ierr = VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);CHKERRQ(ierr); 1611 } 1612 ierr = PetscFree(zeroes);CHKERRQ(ierr); 1613 ierr = DMPrintLocalVec(dm, name, mesh->printTol, locFbc);CHKERRQ(ierr); 1614 ierr = VecDestroy(&locFbc);CHKERRQ(ierr); 1615 } 1616 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1617 PetscFunctionReturn(0); 1618 } 1619 1620 PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM dm, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user) 1621 { 1622 DM_Plex *mesh = (DM_Plex *) dm->data; 1623 const char *name = "Hybrid Residual"; 1624 DM dmAux = NULL; 1625 DMLabel ghostLabel = NULL; 1626 PetscDS prob = NULL; 1627 PetscDS probAux = NULL; 1628 PetscSection section = NULL; 1629 DMField coordField = NULL; 1630 Vec locA; 1631 PetscScalar *u = NULL, *u_t, *a; 1632 PetscScalar *elemVec; 1633 IS chunkIS; 1634 const PetscInt *cells; 1635 PetscInt *faces; 1636 PetscInt cStart, cEnd, numCells; 1637 PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, chunk; 1638 PetscInt maxDegree = PETSC_MAX_INT; 1639 PetscQuadrature affineQuad = NULL, *quads = NULL; 1640 PetscFEGeom *affineGeom = NULL, **geoms = NULL; 1641 PetscErrorCode ierr; 1642 1643 PetscFunctionBegin; 1644 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1645 /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */ 1646 /* FEM */ 1647 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 1648 /* 1: Get sizes from dm and dmAux */ 1649 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 1650 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1651 ierr = DMGetCellDS(dm, cStart, &prob);CHKERRQ(ierr); 1652 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1653 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1654 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1655 if (locA) { 1656 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1657 ierr = DMGetCellDS(dmAux, cStart, &probAux);CHKERRQ(ierr); 1658 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1659 } 1660 /* 2: Setup geometric data */ 1661 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1662 ierr = DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);CHKERRQ(ierr); 1663 if (maxDegree > 1) { 1664 ierr = PetscCalloc2(Nf,&quads,Nf,&geoms);CHKERRQ(ierr); 1665 for (f = 0; f < Nf; ++f) { 1666 PetscFE fe; 1667 1668 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1669 ierr = PetscFEGetQuadrature(fe, &quads[f]);CHKERRQ(ierr); 1670 ierr = PetscObjectReference((PetscObject) quads[f]);CHKERRQ(ierr); 1671 } 1672 } 1673 /* Loop over chunks */ 1674 numCells = cEnd - cStart; 1675 cellChunkSize = numCells; 1676 numChunks = !numCells ? 0 : PetscCeilReal(((PetscReal) numCells)/cellChunkSize); 1677 ierr = PetscCalloc1(2*cellChunkSize, &faces);CHKERRQ(ierr); 1678 ierr = ISCreateGeneral(PETSC_COMM_SELF, cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS);CHKERRQ(ierr); 1679 /* Extract field coefficients */ 1680 /* NOTE This needs the end cap faces to have identical orientations */ 1681 ierr = DMPlexGetCellFields(dm, cellIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1682 ierr = DMGetWorkArray(dm, cellChunkSize*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr); 1683 for (chunk = 0; chunk < numChunks; ++chunk) { 1684 PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c; 1685 1686 ierr = PetscMemzero(elemVec, cellChunkSize*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1687 /* Get faces */ 1688 for (c = cS; c < cE; ++c) { 1689 const PetscInt cell = cells ? cells[c] : c; 1690 const PetscInt *cone; 1691 ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); 1692 faces[(c-cS)*2+0] = cone[0]; 1693 faces[(c-cS)*2+1] = cone[1]; 1694 } 1695 ierr = ISGeneralSetIndices(chunkIS, cellChunkSize, faces, PETSC_USE_POINTER);CHKERRQ(ierr); 1696 /* Get geometric data */ 1697 if (maxDegree <= 1) { 1698 if (!affineQuad) {ierr = DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad);CHKERRQ(ierr);} 1699 if (affineQuad) {ierr = DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom);CHKERRQ(ierr);} 1700 } else { 1701 for (f = 0; f < Nf; ++f) { 1702 ierr = DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]);CHKERRQ(ierr); 1703 } 1704 } 1705 /* Loop over fields */ 1706 for (f = 0; f < Nf; ++f) { 1707 PetscFE fe; 1708 PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 1709 PetscFEGeom *chunkGeom = NULL; 1710 PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; 1711 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb; 1712 1713 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1714 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1715 ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1716 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1717 blockSize = Nb; 1718 batchSize = numBlocks * blockSize; 1719 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1720 numChunks = numCells / (numBatches*batchSize); 1721 Ne = numChunks*numBatches*batchSize; 1722 Nr = numCells % (numBatches*batchSize); 1723 offset = numCells - Nr; 1724 ierr = PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);CHKERRQ(ierr); 1725 ierr = PetscFEIntegrateHybridResidual(prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr); 1726 ierr = PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1727 ierr = PetscFEIntegrateHybridResidual(prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);CHKERRQ(ierr); 1728 ierr = PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr); 1729 } 1730 /* Add elemVec to locX */ 1731 for (c = cS; c < cE; ++c) { 1732 const PetscInt cell = cells ? cells[c] : c; 1733 const PetscInt cind = c - cStart; 1734 1735 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);CHKERRQ(ierr);} 1736 if (ghostLabel) { 1737 PetscInt ghostVal; 1738 1739 ierr = DMLabelGetValue(ghostLabel,cell,&ghostVal);CHKERRQ(ierr); 1740 if (ghostVal > 0) continue; 1741 } 1742 ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);CHKERRQ(ierr); 1743 } 1744 } 1745 ierr = DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1746 ierr = DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr); 1747 ierr = PetscFree(faces);CHKERRQ(ierr); 1748 ierr = ISDestroy(&chunkIS);CHKERRQ(ierr); 1749 ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 1750 if (maxDegree <= 1) { 1751 ierr = DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr); 1752 ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr); 1753 } else { 1754 for (f = 0; f < Nf; ++f) { 1755 if (geoms) {ierr = DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);CHKERRQ(ierr);} 1756 if (quads) {ierr = PetscQuadratureDestroy(&quads[f]);CHKERRQ(ierr);} 1757 } 1758 ierr = PetscFree2(quads,geoms);CHKERRQ(ierr); 1759 } 1760 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1761 PetscFunctionReturn(0); 1762 } 1763 1764 /*@ 1765 DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 1766 1767 Input Parameters: 1768 + dm - The mesh 1769 . X - Local solution 1770 - user - The user context 1771 1772 Output Parameter: 1773 . F - Local output vector 1774 1775 Level: developer 1776 1777 .seealso: DMPlexComputeJacobianAction() 1778 @*/ 1779 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 1780 { 1781 DM plex; 1782 IS cellIS; 1783 PetscInt depth; 1784 PetscErrorCode ierr; 1785 1786 PetscFunctionBegin; 1787 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 1788 ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 1789 ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 1790 if (!cellIS) { 1791 ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr); 1792 } 1793 ierr = DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr); 1794 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 1795 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1796 PetscFunctionReturn(0); 1797 } 1798 1799 /*@ 1800 DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X 1801 1802 Input Parameters: 1803 + dm - The mesh 1804 - user - The user context 1805 1806 Output Parameter: 1807 . X - Local solution 1808 1809 Level: developer 1810 1811 .seealso: DMPlexComputeJacobianAction() 1812 @*/ 1813 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) 1814 { 1815 DM plex; 1816 PetscErrorCode ierr; 1817 1818 PetscFunctionBegin; 1819 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 1820 ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr); 1821 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1822 PetscFunctionReturn(0); 1823 } 1824 1825 PetscErrorCode DMPlexComputeBdJacobian_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, DMField coordField, IS facetIS) 1826 { 1827 DM_Plex *mesh = (DM_Plex *) dm->data; 1828 DM plex = NULL, plexA = NULL, tdm; 1829 DMEnclosureType encAux; 1830 PetscDS prob, probAux = NULL; 1831 PetscSection section, sectionAux = NULL; 1832 PetscSection globalSection, subSection = NULL; 1833 Vec locA = NULL, tv; 1834 PetscScalar *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL; 1835 PetscInt v; 1836 PetscInt Nf, totDim, totDimAux = 0; 1837 PetscBool isMatISP, transform; 1838 PetscErrorCode ierr; 1839 1840 PetscFunctionBegin; 1841 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1842 ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 1843 ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 1844 ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 1845 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 1846 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1847 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1848 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1849 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1850 if (locA) { 1851 DM dmAux; 1852 1853 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1854 ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 1855 ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr); 1856 ierr = DMGetDS(plexA, &probAux);CHKERRQ(ierr); 1857 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1858 ierr = DMGetLocalSection(plexA, §ionAux);CHKERRQ(ierr); 1859 } 1860 1861 ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);CHKERRQ(ierr); 1862 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 1863 if (isMatISP) {ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr);} 1864 for (v = 0; v < numValues; ++v) { 1865 PetscFEGeom *fgeom; 1866 PetscInt maxDegree; 1867 PetscQuadrature qGeom = NULL; 1868 IS pointIS; 1869 const PetscInt *points; 1870 PetscInt numFaces, face, Nq; 1871 1872 ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr); 1873 if (!pointIS) continue; /* No points with that id on this process */ 1874 { 1875 IS isectIS; 1876 1877 /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */ 1878 ierr = ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);CHKERRQ(ierr); 1879 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1880 pointIS = isectIS; 1881 } 1882 ierr = ISGetLocalSize(pointIS, &numFaces);CHKERRQ(ierr); 1883 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1884 ierr = PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim*totDim, &elemMat, locA ? numFaces*totDimAux : 0, &a);CHKERRQ(ierr); 1885 ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 1886 if (maxDegree <= 1) { 1887 ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);CHKERRQ(ierr); 1888 } 1889 if (!qGeom) { 1890 PetscFE fe; 1891 1892 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1893 ierr = PetscFEGetFaceQuadrature(fe, &qGeom);CHKERRQ(ierr); 1894 ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr); 1895 } 1896 ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1897 ierr = DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr); 1898 for (face = 0; face < numFaces; ++face) { 1899 const PetscInt point = points[face], *support, *cone; 1900 PetscScalar *x = NULL; 1901 PetscInt i, coneSize, faceLoc; 1902 1903 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1904 ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); 1905 ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); 1906 for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break; 1907 if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]); 1908 fgeom->face[face][0] = faceLoc; 1909 ierr = DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1910 for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i]; 1911 ierr = DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr); 1912 if (locX_t) { 1913 ierr = DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr); 1914 for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i]; 1915 ierr = DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr); 1916 } 1917 if (locA) { 1918 PetscInt subp; 1919 ierr = DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp);CHKERRQ(ierr); 1920 ierr = DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1921 for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i]; 1922 ierr = DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr); 1923 } 1924 } 1925 ierr = PetscArrayzero(elemMat, numFaces*totDim*totDim);CHKERRQ(ierr); 1926 { 1927 PetscFE fe; 1928 PetscInt Nb; 1929 /* Conforming batches */ 1930 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1931 /* Remainder */ 1932 PetscFEGeom *chunkGeom = NULL; 1933 PetscInt fieldJ, Nr, offset; 1934 1935 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1936 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1937 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1938 blockSize = Nb; 1939 batchSize = numBlocks * blockSize; 1940 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1941 numChunks = numFaces / (numBatches*batchSize); 1942 Ne = numChunks*numBatches*batchSize; 1943 Nr = numFaces % (numBatches*batchSize); 1944 offset = numFaces - Nr; 1945 ierr = PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);CHKERRQ(ierr); 1946 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1947 ierr = PetscFEIntegrateBdJacobian(prob, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr); 1948 } 1949 ierr = PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr); 1950 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1951 ierr = PetscFEIntegrateBdJacobian(prob, fieldI, fieldJ, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 1952 } 1953 ierr = PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr); 1954 } 1955 for (face = 0; face < numFaces; ++face) { 1956 const PetscInt point = points[face], *support; 1957 1958 /* Transform to global basis before insertion in Jacobian */ 1959 ierr = DMPlexGetSupport(plex, point, &support);CHKERRQ(ierr); 1960 if (transform) {ierr = DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face*totDim*totDim]);CHKERRQ(ierr);} 1961 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face*totDim*totDim]);CHKERRQ(ierr);} 1962 if (!isMatISP) { 1963 ierr = DMPlexMatSetClosure(plex, section, globalSection, JacP, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 1964 } else { 1965 Mat lJ; 1966 1967 ierr = MatISGetLocalMat(JacP, &lJ);CHKERRQ(ierr); 1968 ierr = DMPlexMatSetClosure(plex, section, subSection, lJ, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 1969 } 1970 } 1971 ierr = DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr); 1972 ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr); 1973 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1974 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1975 ierr = PetscFree4(u, u_t, elemMat, a);CHKERRQ(ierr); 1976 } 1977 if (plex) {ierr = DMDestroy(&plex);CHKERRQ(ierr);} 1978 if (plexA) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);} 1979 PetscFunctionReturn(0); 1980 } 1981 1982 PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP) 1983 { 1984 DMField coordField; 1985 DMLabel depthLabel; 1986 IS facetIS; 1987 PetscInt dim; 1988 PetscErrorCode ierr; 1989 1990 PetscFunctionBegin; 1991 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1992 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1993 ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr); 1994 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 1995 ierr = DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);CHKERRQ(ierr); 1996 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 1997 PetscFunctionReturn(0); 1998 } 1999 2000 PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user) 2001 { 2002 PetscDS prob; 2003 PetscInt dim, numBd, bd; 2004 DMLabel depthLabel; 2005 DMField coordField = NULL; 2006 IS facetIS; 2007 PetscErrorCode ierr; 2008 2009 PetscFunctionBegin; 2010 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2011 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 2012 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2013 ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr); 2014 ierr = PetscDSGetNumBoundary(prob, &numBd);CHKERRQ(ierr); 2015 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 2016 for (bd = 0; bd < numBd; ++bd) { 2017 DMBoundaryConditionType type; 2018 const char *bdLabel; 2019 DMLabel label; 2020 const PetscInt *values; 2021 PetscInt fieldI, numValues; 2022 PetscObject obj; 2023 PetscClassId id; 2024 2025 ierr = PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &fieldI, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 2026 ierr = PetscDSGetDiscretization(prob, fieldI, &obj);CHKERRQ(ierr); 2027 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2028 if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue; 2029 ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 2030 ierr = DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);CHKERRQ(ierr); 2031 } 2032 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 2033 PetscFunctionReturn(0); 2034 } 2035 2036 PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user) 2037 { 2038 DM_Plex *mesh = (DM_Plex *) dm->data; 2039 const char *name = "Jacobian"; 2040 DM dmAux, plex, tdm; 2041 DMEnclosureType encAux; 2042 Vec A, tv; 2043 DMField coordField; 2044 PetscDS prob, probAux = NULL; 2045 PetscSection section, globalSection, subSection, sectionAux; 2046 PetscScalar *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL; 2047 const PetscInt *cells; 2048 PetscInt Nf, fieldI, fieldJ; 2049 PetscInt totDim, totDimAux, cStart, cEnd, numCells, c; 2050 PetscBool isMatIS, isMatISP, hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE, transform; 2051 PetscErrorCode ierr; 2052 2053 PetscFunctionBegin; 2054 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2055 ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr); 2056 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2057 ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 2058 ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 2059 ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 2060 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 2061 ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);CHKERRQ(ierr); 2062 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 2063 if (isMatISP) {ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr);} 2064 ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr); 2065 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2066 ierr = DMGetCellDS(dm, cStart, &prob);CHKERRQ(ierr); 2067 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 2068 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2069 ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr); 2070 ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr); 2071 /* user passed in the same matrix, avoid double contributions and 2072 only assemble the Jacobian */ 2073 if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE; 2074 ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr); 2075 hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE; 2076 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 2077 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 2078 if (dmAux) { 2079 ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 2080 ierr = DMConvert(dmAux, DMPLEX, &plex);CHKERRQ(ierr); 2081 ierr = DMGetLocalSection(plex, §ionAux);CHKERRQ(ierr); 2082 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 2083 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 2084 } 2085 ierr = PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,hasJac ? numCells*totDim*totDim : 0,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);CHKERRQ(ierr); 2086 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 2087 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 2088 for (c = cStart; c < cEnd; ++c) { 2089 const PetscInt cell = cells ? cells[c] : c; 2090 const PetscInt cind = c - cStart; 2091 PetscScalar *x = NULL, *x_t = NULL; 2092 PetscInt i; 2093 2094 ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr); 2095 for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i]; 2096 ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr); 2097 if (X_t) { 2098 ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr); 2099 for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i]; 2100 ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr); 2101 } 2102 if (dmAux) { 2103 PetscInt subcell; 2104 ierr = DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);CHKERRQ(ierr); 2105 ierr = DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr); 2106 for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i]; 2107 ierr = DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr); 2108 } 2109 } 2110 if (hasJac) {ierr = PetscArrayzero(elemMat, numCells*totDim*totDim);CHKERRQ(ierr);} 2111 if (hasPrec) {ierr = PetscArrayzero(elemMatP, numCells*totDim*totDim);CHKERRQ(ierr);} 2112 if (hasDyn) {ierr = PetscArrayzero(elemMatD, numCells*totDim*totDim);CHKERRQ(ierr);} 2113 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2114 PetscClassId id; 2115 PetscFE fe; 2116 PetscQuadrature qGeom = NULL; 2117 PetscInt Nb; 2118 /* Conforming batches */ 2119 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 2120 /* Remainder */ 2121 PetscInt Nr, offset, Nq; 2122 PetscInt maxDegree; 2123 PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL; 2124 2125 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 2126 ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr); 2127 if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;} 2128 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2129 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2130 ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr); 2131 if (maxDegree <= 1) { 2132 ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr); 2133 } 2134 if (!qGeom) { 2135 ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr); 2136 ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr); 2137 } 2138 ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 2139 ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 2140 blockSize = Nb; 2141 batchSize = numBlocks * blockSize; 2142 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2143 numChunks = numCells / (numBatches*batchSize); 2144 Ne = numChunks*numBatches*batchSize; 2145 Nr = numCells % (numBatches*batchSize); 2146 offset = numCells - Nr; 2147 ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr); 2148 ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr); 2149 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2150 if (hasJac) { 2151 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr); 2152 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 2153 } 2154 if (hasPrec) { 2155 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);CHKERRQ(ierr); 2156 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatP[offset*totDim*totDim]);CHKERRQ(ierr); 2157 } 2158 if (hasDyn) { 2159 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr); 2160 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr); 2161 } 2162 } 2163 ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr); 2164 ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr); 2165 ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 2166 ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr); 2167 } 2168 /* Add contribution from X_t */ 2169 if (hasDyn) {for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];} 2170 if (hasFV) { 2171 PetscClassId id; 2172 PetscFV fv; 2173 PetscInt offsetI, NcI, NbI = 1, fc, f; 2174 2175 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2176 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);CHKERRQ(ierr); 2177 ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr); 2178 ierr = PetscObjectGetClassId((PetscObject) fv, &id);CHKERRQ(ierr); 2179 if (id != PETSCFV_CLASSID) continue; 2180 /* Put in the identity */ 2181 ierr = PetscFVGetNumComponents(fv, &NcI);CHKERRQ(ierr); 2182 for (c = cStart; c < cEnd; ++c) { 2183 const PetscInt cind = c - cStart; 2184 const PetscInt eOffset = cind*totDim*totDim; 2185 for (fc = 0; fc < NcI; ++fc) { 2186 for (f = 0; f < NbI; ++f) { 2187 const PetscInt i = offsetI + f*NcI+fc; 2188 if (hasPrec) { 2189 if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;} 2190 elemMatP[eOffset+i*totDim+i] = 1.0; 2191 } else {elemMat[eOffset+i*totDim+i] = 1.0;} 2192 } 2193 } 2194 } 2195 } 2196 /* No allocated space for FV stuff, so ignore the zero entries */ 2197 ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr); 2198 } 2199 /* Insert values into matrix */ 2200 isMatIS = PETSC_FALSE; 2201 if (hasPrec && hasJac) { 2202 ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatIS);CHKERRQ(ierr); 2203 } 2204 if (isMatIS && !subSection) { 2205 ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr); 2206 } 2207 for (c = cStart; c < cEnd; ++c) { 2208 const PetscInt cell = cells ? cells[c] : c; 2209 const PetscInt cind = c - cStart; 2210 2211 /* Transform to global basis before insertion in Jacobian */ 2212 if (transform) {ierr = DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);} 2213 if (hasPrec) { 2214 if (hasJac) { 2215 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);} 2216 if (!isMatIS) { 2217 ierr = DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2218 } else { 2219 Mat lJ; 2220 2221 ierr = MatISGetLocalMat(Jac,&lJ);CHKERRQ(ierr); 2222 ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2223 } 2224 } 2225 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);CHKERRQ(ierr);} 2226 if (!isMatISP) { 2227 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2228 } else { 2229 Mat lJ; 2230 2231 ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr); 2232 ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2233 } 2234 } else { 2235 if (hasJac) { 2236 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);} 2237 if (!isMatISP) { 2238 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2239 } else { 2240 Mat lJ; 2241 2242 ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr); 2243 ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2244 } 2245 } 2246 } 2247 } 2248 ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2249 if (hasFV) {ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);} 2250 ierr = PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);CHKERRQ(ierr); 2251 if (dmAux) { 2252 ierr = PetscFree(a);CHKERRQ(ierr); 2253 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2254 } 2255 /* Compute boundary integrals */ 2256 ierr = DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user);CHKERRQ(ierr); 2257 /* Assemble matrix */ 2258 if (hasJac && hasPrec) { 2259 ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2260 ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2261 } 2262 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2263 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2264 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2265 PetscFunctionReturn(0); 2266 } 2267 2268 PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user) 2269 { 2270 DM_Plex *mesh = (DM_Plex *) dm->data; 2271 const char *name = "Hybrid Jacobian"; 2272 DM dmAux = NULL; 2273 DM plex = NULL; 2274 DM plexA = NULL; 2275 DMLabel ghostLabel = NULL; 2276 PetscDS prob = NULL; 2277 PetscDS probAux = NULL; 2278 PetscSection section = NULL; 2279 DMField coordField = NULL; 2280 Vec locA; 2281 PetscScalar *u = NULL, *u_t, *a = NULL; 2282 PetscScalar *elemMat, *elemMatP; 2283 PetscSection globalSection, subSection, sectionAux; 2284 IS chunkIS; 2285 const PetscInt *cells; 2286 PetscInt *faces; 2287 PetscInt cStart, cEnd, numCells; 2288 PetscInt Nf, fieldI, fieldJ, totDim, totDimAux, numChunks, cellChunkSize, chunk; 2289 PetscInt maxDegree = PETSC_MAX_INT; 2290 PetscQuadrature affineQuad = NULL, *quads = NULL; 2291 PetscFEGeom *affineGeom = NULL, **geoms = NULL; 2292 PetscBool isMatIS = PETSC_FALSE, isMatISP = PETSC_FALSE, hasBdJac, hasBdPrec; 2293 PetscErrorCode ierr; 2294 2295 PetscFunctionBegin; 2296 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2297 ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr); 2298 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2299 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2300 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 2301 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 2302 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2303 ierr = DMGetCellDS(dm, cStart, &prob);CHKERRQ(ierr); 2304 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 2305 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2306 ierr = PetscDSHasBdJacobian(prob, &hasBdJac);CHKERRQ(ierr); 2307 ierr = PetscDSHasBdJacobianPreconditioner(prob, &hasBdPrec);CHKERRQ(ierr); 2308 ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);CHKERRQ(ierr); 2309 if (isMatISP) {ierr = DMPlexGetSubdomainSection(plex, &subSection);CHKERRQ(ierr);} 2310 if (hasBdPrec && hasBdJac) {ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatIS);CHKERRQ(ierr);} 2311 if (isMatIS && !subSection) {ierr = DMPlexGetSubdomainSection(plex, &subSection);CHKERRQ(ierr);} 2312 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 2313 if (locA) { 2314 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 2315 ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr); 2316 ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 2317 ierr = DMGetCellDS(dmAux, cStart, &probAux);CHKERRQ(ierr); 2318 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 2319 } 2320 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 2321 ierr = DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);CHKERRQ(ierr); 2322 if (maxDegree > 1) { 2323 PetscInt f; 2324 ierr = PetscCalloc2(Nf,&quads,Nf,&geoms);CHKERRQ(ierr); 2325 for (f = 0; f < Nf; ++f) { 2326 PetscFE fe; 2327 2328 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 2329 ierr = PetscFEGetQuadrature(fe, &quads[f]);CHKERRQ(ierr); 2330 ierr = PetscObjectReference((PetscObject) quads[f]);CHKERRQ(ierr); 2331 } 2332 } 2333 cellChunkSize = numCells; 2334 numChunks = !numCells ? 0 : PetscCeilReal(((PetscReal) numCells)/cellChunkSize); 2335 ierr = PetscCalloc1(2*cellChunkSize, &faces);CHKERRQ(ierr); 2336 ierr = ISCreateGeneral(PETSC_COMM_SELF, cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS);CHKERRQ(ierr); 2337 ierr = DMPlexGetCellFields(dm, cellIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 2338 ierr = DMGetWorkArray(dm, hasBdJac ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMat);CHKERRQ(ierr); 2339 ierr = DMGetWorkArray(dm, hasBdPrec ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMatP);CHKERRQ(ierr); 2340 for (chunk = 0; chunk < numChunks; ++chunk) { 2341 PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c; 2342 2343 if (hasBdJac) {ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);} 2344 if (hasBdPrec) {ierr = PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);} 2345 /* Get faces */ 2346 for (c = cS; c < cE; ++c) { 2347 const PetscInt cell = cells ? cells[c] : c; 2348 const PetscInt *cone; 2349 ierr = DMPlexGetCone(plex, cell, &cone);CHKERRQ(ierr); 2350 faces[(c-cS)*2+0] = cone[0]; 2351 faces[(c-cS)*2+1] = cone[1]; 2352 } 2353 ierr = ISGeneralSetIndices(chunkIS, cellChunkSize, faces, PETSC_USE_POINTER);CHKERRQ(ierr); 2354 if (maxDegree <= 1) { 2355 if (!affineQuad) {ierr = DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad);CHKERRQ(ierr);} 2356 if (affineQuad) {ierr = DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom);CHKERRQ(ierr);} 2357 } else { 2358 PetscInt f; 2359 for (f = 0; f < Nf; ++f) { 2360 ierr = DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]);CHKERRQ(ierr); 2361 } 2362 } 2363 2364 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2365 PetscFE fe; 2366 PetscFEGeom *geom = affineGeom ? affineGeom : geoms[fieldI]; 2367 PetscFEGeom *chunkGeom = NULL, *remGeom = NULL; 2368 PetscQuadrature quad = affineQuad ? affineQuad : quads[fieldI]; 2369 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb; 2370 2371 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 2372 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2373 ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 2374 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2375 blockSize = Nb; 2376 batchSize = numBlocks * blockSize; 2377 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2378 numChunks = numCells / (numBatches*batchSize); 2379 Ne = numChunks*numBatches*batchSize; 2380 Nr = numCells % (numBatches*batchSize); 2381 offset = numCells - Nr; 2382 ierr = PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);CHKERRQ(ierr); 2383 ierr = PetscFEGeomGetChunk(geom,offset,numCells,&remGeom);CHKERRQ(ierr); 2384 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2385 if (hasBdJac) { 2386 ierr = PetscFEIntegrateHybridJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr); 2387 ierr = PetscFEIntegrateHybridJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 2388 } 2389 if (hasBdPrec) { 2390 ierr = PetscFEIntegrateHybridJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);CHKERRQ(ierr); 2391 ierr = PetscFEIntegrateHybridJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatP[offset*totDim*totDim]);CHKERRQ(ierr); 2392 } 2393 } 2394 ierr = PetscFEGeomRestoreChunk(geom,offset,numCells,&remGeom);CHKERRQ(ierr); 2395 ierr = PetscFEGeomRestoreChunk(geom,0,offset,&chunkGeom);CHKERRQ(ierr); 2396 } 2397 /* Insert values into matrix */ 2398 for (c = cS; c < cE; ++c) { 2399 const PetscInt cell = cells ? cells[c] : c; 2400 const PetscInt cind = c - cS; 2401 2402 if (hasBdPrec) { 2403 if (hasBdJac) { 2404 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);} 2405 if (!isMatIS) { 2406 ierr = DMPlexMatSetClosure(plex, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2407 } else { 2408 Mat lJ; 2409 2410 ierr = MatISGetLocalMat(Jac,&lJ);CHKERRQ(ierr); 2411 ierr = DMPlexMatSetClosure(plex, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2412 } 2413 } 2414 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);CHKERRQ(ierr);} 2415 if (!isMatISP) { 2416 ierr = DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2417 } else { 2418 Mat lJ; 2419 2420 ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr); 2421 ierr = DMPlexMatSetClosure(plex, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2422 } 2423 } else if (hasBdJac) { 2424 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);} 2425 if (!isMatISP) { 2426 ierr = DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2427 } else { 2428 Mat lJ; 2429 2430 ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr); 2431 ierr = DMPlexMatSetClosure(plex, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2432 } 2433 } 2434 } 2435 } 2436 ierr = DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 2437 ierr = DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMat);CHKERRQ(ierr); 2438 ierr = DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMatP);CHKERRQ(ierr); 2439 ierr = PetscFree(faces);CHKERRQ(ierr); 2440 ierr = ISDestroy(&chunkIS);CHKERRQ(ierr); 2441 ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2442 if (maxDegree <= 1) { 2443 ierr = DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr); 2444 ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr); 2445 } else { 2446 PetscInt f; 2447 for (f = 0; f < Nf; ++f) { 2448 if (geoms) {ierr = DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE, &geoms[f]);CHKERRQ(ierr);} 2449 if (quads) {ierr = PetscQuadratureDestroy(&quads[f]);CHKERRQ(ierr);} 2450 } 2451 ierr = PetscFree2(quads,geoms);CHKERRQ(ierr); 2452 } 2453 if (dmAux) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);} 2454 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2455 /* Assemble matrix */ 2456 if (hasBdJac && hasBdPrec) { 2457 ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2458 ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2459 } 2460 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2461 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2462 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2463 PetscFunctionReturn(0); 2464 } 2465 2466 /*@ 2467 DMPlexComputeJacobianAction - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user. 2468 2469 Input Parameters: 2470 + dm - The mesh 2471 . cellIS - 2472 . t - The time 2473 . X_tShift - The multiplier for the Jacobian with repsect to X_t 2474 . X - Local solution vector 2475 . X_t - Time-derivative of the local solution vector 2476 . Y - Local input vector 2477 - user - The user context 2478 2479 Output Parameter: 2480 . Z - Local output vector 2481 2482 Note: 2483 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 2484 like a GPU, or vectorize on a multicore machine. 2485 2486 Level: developer 2487 2488 .seealso: FormFunctionLocal() 2489 @*/ 2490 PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user) 2491 { 2492 DM_Plex *mesh = (DM_Plex *) dm->data; 2493 const char *name = "Jacobian"; 2494 DM dmAux, plex, plexAux = NULL; 2495 DMEnclosureType encAux; 2496 Vec A; 2497 PetscDS prob, probAux = NULL; 2498 PetscQuadrature quad; 2499 PetscSection section, globalSection, sectionAux; 2500 PetscScalar *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z; 2501 PetscInt Nf, fieldI, fieldJ; 2502 PetscInt totDim, totDimAux = 0; 2503 const PetscInt *cells; 2504 PetscInt cStart, cEnd, numCells, c; 2505 PetscBool hasDyn; 2506 DMField coordField; 2507 PetscErrorCode ierr; 2508 2509 PetscFunctionBegin; 2510 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2511 ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 2512 if (!cellIS) { 2513 PetscInt depth; 2514 2515 ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 2516 ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 2517 if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);} 2518 } else { 2519 ierr = PetscObjectReference((PetscObject) cellIS);CHKERRQ(ierr); 2520 } 2521 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 2522 ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr); 2523 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2524 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2525 ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr); 2526 hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE; 2527 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 2528 ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr); 2529 ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr); 2530 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 2531 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 2532 if (dmAux) { 2533 ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 2534 ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr); 2535 ierr = DMGetLocalSection(plexAux, §ionAux);CHKERRQ(ierr); 2536 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 2537 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 2538 } 2539 ierr = VecSet(Z, 0.0);CHKERRQ(ierr); 2540 ierr = PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);CHKERRQ(ierr); 2541 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 2542 ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 2543 for (c = cStart; c < cEnd; ++c) { 2544 const PetscInt cell = cells ? cells[c] : c; 2545 const PetscInt cind = c - cStart; 2546 PetscScalar *x = NULL, *x_t = NULL; 2547 PetscInt i; 2548 2549 ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr); 2550 for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i]; 2551 ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr); 2552 if (X_t) { 2553 ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr); 2554 for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i]; 2555 ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr); 2556 } 2557 if (dmAux) { 2558 PetscInt subcell; 2559 ierr = DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);CHKERRQ(ierr); 2560 ierr = DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr); 2561 for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i]; 2562 ierr = DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr); 2563 } 2564 ierr = DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr); 2565 for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i]; 2566 ierr = DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr); 2567 } 2568 ierr = PetscArrayzero(elemMat, numCells*totDim*totDim);CHKERRQ(ierr); 2569 if (hasDyn) {ierr = PetscArrayzero(elemMatD, numCells*totDim*totDim);CHKERRQ(ierr);} 2570 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2571 PetscFE fe; 2572 PetscInt Nb; 2573 /* Conforming batches */ 2574 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 2575 /* Remainder */ 2576 PetscInt Nr, offset, Nq; 2577 PetscQuadrature qGeom = NULL; 2578 PetscInt maxDegree; 2579 PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL; 2580 2581 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 2582 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 2583 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2584 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2585 ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr); 2586 if (maxDegree <= 1) {ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);} 2587 if (!qGeom) { 2588 ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr); 2589 ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr); 2590 } 2591 ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 2592 ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 2593 blockSize = Nb; 2594 batchSize = numBlocks * blockSize; 2595 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2596 numChunks = numCells / (numBatches*batchSize); 2597 Ne = numChunks*numBatches*batchSize; 2598 Nr = numCells % (numBatches*batchSize); 2599 offset = numCells - Nr; 2600 ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr); 2601 ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr); 2602 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2603 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr); 2604 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 2605 if (hasDyn) { 2606 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr); 2607 ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr); 2608 } 2609 } 2610 ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr); 2611 ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr); 2612 ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr); 2613 ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr); 2614 } 2615 if (hasDyn) { 2616 for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c]; 2617 } 2618 for (c = cStart; c < cEnd; ++c) { 2619 const PetscInt cell = cells ? cells[c] : c; 2620 const PetscInt cind = c - cStart; 2621 const PetscBLASInt M = totDim, one = 1; 2622 const PetscScalar a = 1.0, b = 0.0; 2623 2624 PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one)); 2625 if (mesh->printFEM > 1) { 2626 ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr); 2627 ierr = DMPrintCellVector(c, "Y", totDim, &y[cind*totDim]);CHKERRQ(ierr); 2628 ierr = DMPrintCellVector(c, "Z", totDim, z);CHKERRQ(ierr); 2629 } 2630 ierr = DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);CHKERRQ(ierr); 2631 } 2632 ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr); 2633 if (mesh->printFEM) { 2634 ierr = PetscPrintf(PETSC_COMM_WORLD, "Z:\n");CHKERRQ(ierr); 2635 ierr = VecView(Z, NULL);CHKERRQ(ierr); 2636 } 2637 ierr = PetscFree(a);CHKERRQ(ierr); 2638 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 2639 ierr = DMDestroy(&plexAux);CHKERRQ(ierr); 2640 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2641 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2642 PetscFunctionReturn(0); 2643 } 2644 2645 /*@ 2646 DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 2647 2648 Input Parameters: 2649 + dm - The mesh 2650 . X - Local input vector 2651 - user - The user context 2652 2653 Output Parameter: 2654 . Jac - Jacobian matrix 2655 2656 Note: 2657 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 2658 like a GPU, or vectorize on a multicore machine. 2659 2660 Level: developer 2661 2662 .seealso: FormFunctionLocal() 2663 @*/ 2664 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 2665 { 2666 DM plex; 2667 PetscDS prob; 2668 IS cellIS; 2669 PetscBool hasJac, hasPrec; 2670 PetscInt depth; 2671 PetscErrorCode ierr; 2672 2673 PetscFunctionBegin; 2674 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 2675 ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 2676 ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 2677 if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);} 2678 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2679 ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr); 2680 ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr); 2681 if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);} 2682 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 2683 ierr = DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr); 2684 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 2685 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2686 PetscFunctionReturn(0); 2687 } 2688 2689 /* 2690 MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context. 2691 2692 Input Parameters: 2693 + X - SNES linearization point 2694 . ovl - index set of overlapping subdomains 2695 2696 Output Parameter: 2697 . J - unassembled (Neumann) local matrix 2698 2699 Level: intermediate 2700 2701 .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat() 2702 */ 2703 static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) 2704 { 2705 SNES snes; 2706 Mat pJ; 2707 DM ovldm,origdm; 2708 DMSNES sdm; 2709 PetscErrorCode (*bfun)(DM,Vec,void*); 2710 PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*); 2711 void *bctx,*jctx; 2712 PetscErrorCode ierr; 2713 2714 PetscFunctionBegin; 2715 ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr); 2716 if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr); 2717 ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr); 2718 if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr); 2719 ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr); 2720 ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr); 2721 ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr); 2722 ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr); 2723 ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr); 2724 ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr); 2725 if (!snes) { 2726 ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr); 2727 ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr); 2728 ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr); 2729 ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr); 2730 } 2731 ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr); 2732 ierr = VecLockReadPush(X);CHKERRQ(ierr); 2733 PetscStackPush("SNES user Jacobian function"); 2734 ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr); 2735 PetscStackPop; 2736 ierr = VecLockReadPop(X);CHKERRQ(ierr); 2737 /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */ 2738 { 2739 Mat locpJ; 2740 2741 ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr); 2742 ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2743 } 2744 PetscFunctionReturn(0); 2745 } 2746 2747 /*@ 2748 DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian. 2749 2750 Input Parameters: 2751 + dm - The DM object 2752 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary()) 2753 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual()) 2754 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian()) 2755 2756 Level: developer 2757 @*/ 2758 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) 2759 { 2760 PetscErrorCode ierr; 2761 2762 PetscFunctionBegin; 2763 ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr); 2764 ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr); 2765 ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr); 2766 ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr); 2767 PetscFunctionReturn(0); 2768 } 2769 2770 /*@C 2771 DMSNESCheckDiscretization - Check the discretization error of the exact solution 2772 2773 Input Parameters: 2774 + snes - the SNES object 2775 . dm - the DM 2776 . u - a DM vector 2777 . exactFuncs - pointwise functions of the exact solution for each field 2778 . ctxs - contexts for the functions 2779 - tol - A tolerance for the check, or -1 to print the results instead 2780 2781 Output Parameters: 2782 . error - An array which holds the discretization error in each field, or NULL 2783 2784 Level: developer 2785 2786 .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian() 2787 @*/ 2788 PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs, PetscReal tol, PetscReal error[]) 2789 { 2790 PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 2791 void **ectxs; 2792 MPI_Comm comm; 2793 PetscDS ds; 2794 PetscReal *err; 2795 PetscInt Nf, f; 2796 PetscErrorCode ierr; 2797 2798 PetscFunctionBegin; 2799 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2800 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 2801 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 2802 if (error) PetscValidRealPointer(error, 6); 2803 ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr); 2804 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 2805 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 2806 ierr = PetscMalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr); 2807 for (f = 0; f < Nf; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);CHKERRQ(ierr);} 2808 ierr = DMProjectFunction(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 2809 ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr); 2810 ierr = PetscObjectSetOptionsPrefix((PetscObject) u, "exact_");CHKERRQ(ierr); 2811 ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr); 2812 if (Nf > 1) { 2813 ierr = DMComputeL2FieldDiff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, u, err);CHKERRQ(ierr); 2814 if (tol >= 0.0) { 2815 for (f = 0; f < Nf; ++f) { 2816 if (err[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) err[f], f, (double) tol); 2817 } 2818 } else if (error) { 2819 for (f = 0; f < Nf; ++f) error[f] = err[f]; 2820 } else { 2821 ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr); 2822 for (f = 0; f < Nf; ++f) { 2823 if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 2824 ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr); 2825 } 2826 ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr); 2827 } 2828 } else { 2829 ierr = DMComputeL2Diff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs , u, &err[0]);CHKERRQ(ierr); 2830 if (tol >= 0.0) { 2831 if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol); 2832 } else if (error) { 2833 error[0] = err[0]; 2834 } else { 2835 ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]);CHKERRQ(ierr); 2836 } 2837 } 2838 ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr); 2839 PetscFunctionReturn(0); 2840 } 2841 2842 /*@C 2843 DMSNESCheckResidual - Check the residual of the exact solution 2844 2845 Input Parameters: 2846 + snes - the SNES object 2847 . dm - the DM 2848 . u - a DM vector 2849 - tol - A tolerance for the check, or -1 to print the results instead 2850 2851 Output Parameters: 2852 . residual - The residual norm of the exact solution, or NULL 2853 2854 Level: developer 2855 2856 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian() 2857 @*/ 2858 PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) 2859 { 2860 PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 2861 void **ectxs; 2862 MPI_Comm comm; 2863 PetscDS ds; 2864 Vec r; 2865 PetscReal res; 2866 PetscInt Nf, f; 2867 PetscBool computeSol = PETSC_FALSE; 2868 PetscErrorCode ierr; 2869 2870 PetscFunctionBegin; 2871 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2872 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 2873 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 2874 if (residual) PetscValidRealPointer(residual, 5); 2875 ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr); 2876 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 2877 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 2878 ierr = PetscMalloc2(Nf, &exacts, Nf, &ectxs);CHKERRQ(ierr); 2879 for (f = 0; f < Nf; ++f) { 2880 ierr = PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);CHKERRQ(ierr); 2881 if (exacts[f]) computeSol = PETSC_TRUE; 2882 } 2883 if (computeSol) {ierr = DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);} 2884 ierr = PetscFree2(exacts, ectxs);CHKERRQ(ierr); 2885 2886 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 2887 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 2888 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 2889 if (tol >= 0.0) { 2890 if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol); 2891 } else if (residual) { 2892 *residual = res; 2893 } else { 2894 ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 2895 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 2896 ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr); 2897 ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr); 2898 ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 2899 } 2900 ierr = VecDestroy(&r);CHKERRQ(ierr); 2901 PetscFunctionReturn(0); 2902 } 2903 2904 /*@C 2905 DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 2906 2907 Input Parameters: 2908 + snes - the SNES object 2909 . dm - the DM 2910 . u - a DM vector 2911 - tol - A tolerance for the check, or -1 to print the results instead 2912 2913 Output Parameters: 2914 + isLinear - Flag indicaing that the function looks linear, or NULL 2915 - convRate - The rate of convergence of the linear model, or NULL 2916 2917 Level: developer 2918 2919 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual() 2920 @*/ 2921 PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 2922 { 2923 PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 2924 void **ectxs; 2925 MPI_Comm comm; 2926 PetscDS ds; 2927 Mat J, M; 2928 MatNullSpace nullspace; 2929 PetscReal slope, intercept; 2930 PetscInt Nf, f; 2931 PetscBool hasJac, hasPrec, isLin = PETSC_FALSE, computeSol = PETSC_FALSE; 2932 PetscErrorCode ierr; 2933 2934 PetscFunctionBegin; 2935 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2936 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 2937 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 2938 if (isLinear) PetscValidBoolPointer(isLinear, 5); 2939 if (convRate) PetscValidRealPointer(convRate, 5); 2940 ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr); 2941 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 2942 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 2943 ierr = PetscMalloc2(Nf, &exacts, Nf, &ectxs);CHKERRQ(ierr); 2944 for (f = 0; f < Nf; ++f) { 2945 ierr = PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);CHKERRQ(ierr); 2946 if (exacts[f]) computeSol = PETSC_TRUE; 2947 } 2948 if (computeSol) {ierr = DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);} 2949 ierr = PetscFree2(exacts, ectxs);CHKERRQ(ierr); 2950 2951 /* Create and view matrices */ 2952 ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 2953 ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 2954 ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 2955 if (hasJac && hasPrec) { 2956 ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 2957 ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr); 2958 ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr); 2959 ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr); 2960 ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr); 2961 ierr = MatDestroy(&M);CHKERRQ(ierr); 2962 } else { 2963 ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr); 2964 } 2965 ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr); 2966 ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr); 2967 ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr); 2968 /* Check nullspace */ 2969 ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 2970 if (nullspace) { 2971 PetscBool isNull; 2972 ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr); 2973 if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 2974 } 2975 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr); 2976 /* Taylor test */ 2977 { 2978 PetscRandom rand; 2979 Vec du, uhat, r, rhat, df; 2980 PetscReal h; 2981 PetscReal *es, *hs, *errors; 2982 PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 2983 PetscInt Nv, v; 2984 2985 /* Choose a perturbation direction */ 2986 ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr); 2987 ierr = VecDuplicate(u, &du);CHKERRQ(ierr); 2988 ierr = VecSetRandom(du, rand); CHKERRQ(ierr); 2989 ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 2990 ierr = VecDuplicate(u, &df);CHKERRQ(ierr); 2991 ierr = MatMult(J, du, df);CHKERRQ(ierr); 2992 /* Evaluate residual at u, F(u), save in vector r */ 2993 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 2994 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 2995 /* Look at the convergence of our Taylor approximation as we approach u */ 2996 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv); 2997 ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr); 2998 ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr); 2999 ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr); 3000 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 3001 ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr); 3002 /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */ 3003 ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr); 3004 ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr); 3005 ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr); 3006 3007 es[Nv] = PetscLog10Real(errors[Nv]); 3008 hs[Nv] = PetscLog10Real(h); 3009 } 3010 ierr = VecDestroy(&uhat);CHKERRQ(ierr); 3011 ierr = VecDestroy(&rhat);CHKERRQ(ierr); 3012 ierr = VecDestroy(&df);CHKERRQ(ierr); 3013 ierr = VecDestroy(&r);CHKERRQ(ierr); 3014 ierr = VecDestroy(&du);CHKERRQ(ierr); 3015 for (v = 0; v < Nv; ++v) { 3016 if ((tol >= 0) && (errors[v] > tol)) break; 3017 else if (errors[v] > PETSC_SMALL) break; 3018 } 3019 if (v == Nv) isLin = PETSC_TRUE; 3020 ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr); 3021 ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr); 3022 /* Slope should be about 2 */ 3023 if (tol >= 0) { 3024 if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope); 3025 } else if (isLinear || convRate) { 3026 if (isLinear) *isLinear = isLin; 3027 if (convRate) *convRate = slope; 3028 } else { 3029 if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);} 3030 else {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);} 3031 } 3032 } 3033 ierr = MatDestroy(&J);CHKERRQ(ierr); 3034 PetscFunctionReturn(0); 3035 } 3036 3037 PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs) 3038 { 3039 PetscErrorCode ierr; 3040 3041 PetscFunctionBegin; 3042 ierr = DMSNESCheckDiscretization(snes, dm, u, exactFuncs, ctxs, -1.0, NULL);CHKERRQ(ierr); 3043 ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr); 3044 ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr); 3045 PetscFunctionReturn(0); 3046 } 3047 3048 /*@C 3049 DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 3050 3051 Input Parameters: 3052 + snes - the SNES object 3053 . u - representative SNES vector 3054 . exactFuncs - pointwise functions of the exact solution for each field 3055 - ctxs - contexts for the functions 3056 3057 Level: developer 3058 @*/ 3059 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs) 3060 { 3061 DM dm; 3062 Vec sol; 3063 PetscBool check; 3064 PetscErrorCode ierr; 3065 3066 PetscFunctionBegin; 3067 ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr); 3068 if (!check) PetscFunctionReturn(0); 3069 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 3070 ierr = VecDuplicate(u, &sol);CHKERRQ(ierr); 3071 ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr); 3072 ierr = DMSNESCheck_Internal(snes, dm, sol, exactFuncs, ctxs);CHKERRQ(ierr); 3073 ierr = VecDestroy(&sol);CHKERRQ(ierr); 3074 PetscFunctionReturn(0); 3075 } 3076