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 <petsc/private/petscimpl.h> 5 #include <petsc/private/petscfeimpl.h> 6 7 static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[]) { 8 p[0] = u[uOff[1]]; 9 } 10 11 /* 12 SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately. 13 14 Collective on snes 15 16 Input Parameters: 17 + snes - The SNES 18 . pfield - The field number for pressure 19 . nullspace - The pressure nullspace 20 . u - The solution vector 21 - ctx - An optional user context 22 23 Output Parameter: 24 . u - The solution with a continuum pressure integral of zero 25 26 Notes: 27 If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly. 28 29 Level: developer 30 31 .seealso: `SNESConvergedCorrectPressure()` 32 */ 33 static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx) { 34 DM dm; 35 PetscDS ds; 36 const Vec *nullvecs; 37 PetscScalar pintd, *intc, *intn; 38 MPI_Comm comm; 39 PetscInt Nf, Nv; 40 41 PetscFunctionBegin; 42 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 43 PetscCall(SNESGetDM(snes, &dm)); 44 PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM"); 45 PetscCheck(nullspace, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace"); 46 PetscCall(DMGetDS(dm, &ds)); 47 PetscCall(PetscDSSetObjective(ds, pfield, pressure_Private)); 48 PetscCall(MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs)); 49 PetscCheck(Nv == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %" PetscInt_FMT, Nv); 50 PetscCall(VecDot(nullvecs[0], u, &pintd)); 51 PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double)PetscRealPart(pintd)); 52 PetscCall(PetscDSGetNumFields(ds, &Nf)); 53 PetscCall(PetscMalloc2(Nf, &intc, Nf, &intn)); 54 PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx)); 55 PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx)); 56 PetscCall(VecAXPY(u, -intc[pfield] / intn[pfield], nullvecs[0])); 57 #if defined(PETSC_USE_DEBUG) 58 PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx)); 59 PetscCheck(PetscAbsScalar(intc[pfield]) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double)PetscRealPart(intc[pfield])); 60 #endif 61 PetscCall(PetscFree2(intc, intn)); 62 PetscFunctionReturn(0); 63 } 64 65 /*@C 66 SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero. 67 This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics `SNESConvergedDefault()`. 68 69 Logically Collective on snes 70 71 Input Parameters: 72 + snes - the `SNES` context 73 . it - the iteration (0 indicates before any Newton steps) 74 . xnorm - 2-norm of current iterate 75 . snorm - 2-norm of current step 76 . fnorm - 2-norm of function at current iterate 77 - ctx - Optional user context 78 79 Output Parameter: 80 . reason - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FNORM_NAN` 81 82 Notes: 83 In order to use this convergence test, you must set up several PETSc structures. First fields must be added to the `DM`, and a `PetscDS` 84 must be created with discretizations of those fields. We currently assume that the pressure field has index 1. 85 The pressure field must have a nullspace, likely created using the `DMSetNullSpaceConstructor()` interface. 86 Last we must be able to integrate the pressure over the domain, so the `DM` attached to the SNES `must` be a `DMPLEX` at this time. 87 88 Options Database Key: 89 . -snes_convergence_test correct_pressure - see `SNESSetFromOptions()` 90 91 Level: advanced 92 93 .seealso: `SNES`, `DM`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()`, `DMSetNullSpaceConstructor()` 94 @*/ 95 PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx) { 96 PetscBool monitorIntegral = PETSC_FALSE; 97 98 PetscFunctionBegin; 99 PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx)); 100 if (monitorIntegral) { 101 Mat J; 102 Vec u; 103 MatNullSpace nullspace; 104 const Vec *nullvecs; 105 PetscScalar pintd; 106 107 PetscCall(SNESGetSolution(snes, &u)); 108 PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 109 PetscCall(MatGetNullSpace(J, &nullspace)); 110 PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs)); 111 PetscCall(VecDot(nullvecs[0], u, &pintd)); 112 PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double)PetscRealPart(pintd))); 113 } 114 if (*reason > 0) { 115 Mat J; 116 Vec u; 117 MatNullSpace nullspace; 118 PetscInt pfield = 1; 119 120 PetscCall(SNESGetSolution(snes, &u)); 121 PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 122 PetscCall(MatGetNullSpace(J, &nullspace)); 123 PetscCall(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx)); 124 } 125 PetscFunctionReturn(0); 126 } 127 128 /************************** Interpolation *******************************/ 129 130 static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) { 131 PetscBool isPlex; 132 133 PetscFunctionBegin; 134 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex)); 135 if (isPlex) { 136 *plex = dm; 137 PetscCall(PetscObjectReference((PetscObject)dm)); 138 } else { 139 PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex)); 140 if (!*plex) { 141 PetscCall(DMConvert(dm, DMPLEX, plex)); 142 PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex)); 143 if (copy) { 144 PetscCall(DMCopyDMSNES(dm, *plex)); 145 PetscCall(DMCopyAuxiliaryVec(dm, *plex)); 146 } 147 } else { 148 PetscCall(PetscObjectReference((PetscObject)*plex)); 149 } 150 } 151 PetscFunctionReturn(0); 152 } 153 154 /*@C 155 DMInterpolationCreate - Creates a `DMInterpolationInfo` context 156 157 Collective 158 159 Input Parameter: 160 . comm - the communicator 161 162 Output Parameter: 163 . ctx - the context 164 165 Level: beginner 166 167 .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()` 168 @*/ 169 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) { 170 PetscFunctionBegin; 171 PetscValidPointer(ctx, 2); 172 PetscCall(PetscNew(ctx)); 173 174 (*ctx)->comm = comm; 175 (*ctx)->dim = -1; 176 (*ctx)->nInput = 0; 177 (*ctx)->points = NULL; 178 (*ctx)->cells = NULL; 179 (*ctx)->n = -1; 180 (*ctx)->coords = NULL; 181 PetscFunctionReturn(0); 182 } 183 184 /*@C 185 DMInterpolationSetDim - Sets the spatial dimension for the interpolation context 186 187 Not collective 188 189 Input Parameters: 190 + ctx - the context 191 - dim - the spatial dimension 192 193 Level: intermediate 194 195 .seealso: `DMInterpolationInfo`, `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 196 @*/ 197 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) { 198 PetscFunctionBegin; 199 PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim); 200 ctx->dim = dim; 201 PetscFunctionReturn(0); 202 } 203 204 /*@C 205 DMInterpolationGetDim - Gets the spatial dimension for the interpolation context 206 207 Not collective 208 209 Input Parameter: 210 . ctx - the context 211 212 Output Parameter: 213 . dim - the spatial dimension 214 215 Level: intermediate 216 217 .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 218 @*/ 219 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) { 220 PetscFunctionBegin; 221 PetscValidIntPointer(dim, 2); 222 *dim = ctx->dim; 223 PetscFunctionReturn(0); 224 } 225 226 /*@C 227 DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context 228 229 Not collective 230 231 Input Parameters: 232 + ctx - the context 233 - dof - the number of fields 234 235 Level: intermediate 236 237 .seealso: `DMInterpolationInfo`, `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 238 @*/ 239 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) { 240 PetscFunctionBegin; 241 PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof); 242 ctx->dof = dof; 243 PetscFunctionReturn(0); 244 } 245 246 /*@C 247 DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context 248 249 Not collective 250 251 Input Parameter: 252 . ctx - the context 253 254 Output Parameter: 255 . dof - the number of fields 256 257 Level: intermediate 258 259 .seealso: DMInterpolationInfo, `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 260 @*/ 261 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) { 262 PetscFunctionBegin; 263 PetscValidIntPointer(dof, 2); 264 *dof = ctx->dof; 265 PetscFunctionReturn(0); 266 } 267 268 /*@C 269 DMInterpolationAddPoints - Add points at which we will interpolate the fields 270 271 Not collective 272 273 Input Parameters: 274 + ctx - the context 275 . n - the number of points 276 - points - the coordinates for each point, an array of size n * dim 277 278 Note: 279 The coordinate information is copied. 280 281 Level: intermediate 282 283 .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()` 284 @*/ 285 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) { 286 PetscFunctionBegin; 287 PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 288 PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 289 ctx->nInput = n; 290 291 PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points)); 292 PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim)); 293 PetscFunctionReturn(0); 294 } 295 296 /*@C 297 DMInterpolationSetUp - Compute spatial indices for point location during interpolation 298 299 Collective on ctx 300 301 Input Parameters: 302 + ctx - the context 303 . dm - the `DM` for the function space used for interpolation 304 . redundantPoints - If `PETSC_TRUE`, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes. 305 - ignoreOutsideDomain - If `PETSC_TRUE`, ignore points outside the domain, otherwise return an error 306 307 Level: intermediate 308 309 .seealso: DMInterpolationInfo, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 310 @*/ 311 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain) { 312 MPI_Comm comm = ctx->comm; 313 PetscScalar *a; 314 PetscInt p, q, i; 315 PetscMPIInt rank, size; 316 Vec pointVec; 317 PetscSF cellSF; 318 PetscLayout layout; 319 PetscReal *globalPoints; 320 PetscScalar *globalPointsScalar; 321 const PetscInt *ranges; 322 PetscMPIInt *counts, *displs; 323 const PetscSFNode *foundCells; 324 const PetscInt *foundPoints; 325 PetscMPIInt *foundProcs, *globalProcs; 326 PetscInt n, N, numFound; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 330 PetscCallMPI(MPI_Comm_size(comm, &size)); 331 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 332 PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 333 /* Locate points */ 334 n = ctx->nInput; 335 if (!redundantPoints) { 336 PetscCall(PetscLayoutCreate(comm, &layout)); 337 PetscCall(PetscLayoutSetBlockSize(layout, 1)); 338 PetscCall(PetscLayoutSetLocalSize(layout, n)); 339 PetscCall(PetscLayoutSetUp(layout)); 340 PetscCall(PetscLayoutGetSize(layout, &N)); 341 /* Communicate all points to all processes */ 342 PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs)); 343 PetscCall(PetscLayoutGetRanges(layout, &ranges)); 344 for (p = 0; p < size; ++p) { 345 counts[p] = (ranges[p + 1] - ranges[p]) * ctx->dim; 346 displs[p] = ranges[p] * ctx->dim; 347 } 348 PetscCallMPI(MPI_Allgatherv(ctx->points, n * ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm)); 349 } else { 350 N = n; 351 globalPoints = ctx->points; 352 counts = displs = NULL; 353 layout = NULL; 354 } 355 #if 0 356 PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs)); 357 /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 358 #else 359 #if defined(PETSC_USE_COMPLEX) 360 PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar)); 361 for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i]; 362 #else 363 globalPointsScalar = globalPoints; 364 #endif 365 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec)); 366 PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs)); 367 for (p = 0; p < N; ++p) foundProcs[p] = size; 368 cellSF = NULL; 369 PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF)); 370 PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells)); 371 #endif 372 for (p = 0; p < numFound; ++p) { 373 if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank; 374 } 375 /* Let the lowest rank process own each point */ 376 PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm)); 377 ctx->n = 0; 378 for (p = 0; p < N; ++p) { 379 if (globalProcs[p] == size) { 380 PetscCheck(ignoreOutsideDomain, comm, PETSC_ERR_PLIB, "Point %" PetscInt_FMT ": %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), 381 (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0.0)); 382 if (rank == 0) ++ctx->n; 383 } else if (globalProcs[p] == rank) ++ctx->n; 384 } 385 /* Create coordinates vector and array of owned cells */ 386 PetscCall(PetscMalloc1(ctx->n, &ctx->cells)); 387 PetscCall(VecCreate(comm, &ctx->coords)); 388 PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE)); 389 PetscCall(VecSetBlockSize(ctx->coords, ctx->dim)); 390 PetscCall(VecSetType(ctx->coords, VECSTANDARD)); 391 PetscCall(VecGetArray(ctx->coords, &a)); 392 for (p = 0, q = 0, i = 0; p < N; ++p) { 393 if (globalProcs[p] == rank) { 394 PetscInt d; 395 396 for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d]; 397 ctx->cells[q] = foundCells[q].index; 398 ++q; 399 } 400 if (globalProcs[p] == size && rank == 0) { 401 PetscInt d; 402 403 for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.; 404 ctx->cells[q] = -1; 405 ++q; 406 } 407 } 408 PetscCall(VecRestoreArray(ctx->coords, &a)); 409 #if 0 410 PetscCall(PetscFree3(foundCells,foundProcs,globalProcs)); 411 #else 412 PetscCall(PetscFree2(foundProcs, globalProcs)); 413 PetscCall(PetscSFDestroy(&cellSF)); 414 PetscCall(VecDestroy(&pointVec)); 415 #endif 416 if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar)); 417 if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs)); 418 PetscCall(PetscLayoutDestroy(&layout)); 419 PetscFunctionReturn(0); 420 } 421 422 /*@C 423 DMInterpolationGetCoordinates - Gets a `Vec` with the coordinates of each interpolation point 424 425 Collective on ctx 426 427 Input Parameter: 428 . ctx - the context 429 430 Output Parameter: 431 . coordinates - the coordinates of interpolation points 432 433 Note: 434 The local vector entries correspond to interpolation points lying on this process, according to the associated `DM`. 435 This is a borrowed vector that the user should not destroy. 436 437 Level: intermediate 438 439 .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 440 @*/ 441 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) { 442 PetscFunctionBegin; 443 PetscValidPointer(coordinates, 2); 444 PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 445 *coordinates = ctx->coords; 446 PetscFunctionReturn(0); 447 } 448 449 /*@C 450 DMInterpolationGetVector - Gets a `Vec` which can hold all the interpolated field values 451 452 Collective on ctx 453 454 Input Parameter: 455 . ctx - the context 456 457 Output Parameter: 458 . v - a vector capable of holding the interpolated field values 459 460 Note: 461 This vector should be returned using `DMInterpolationRestoreVector()`. 462 463 Level: intermediate 464 465 .seealso: `DMInterpolationInfo`, `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 466 @*/ 467 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) { 468 PetscFunctionBegin; 469 PetscValidPointer(v, 2); 470 PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 471 PetscCall(VecCreate(ctx->comm, v)); 472 PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE)); 473 PetscCall(VecSetBlockSize(*v, ctx->dof)); 474 PetscCall(VecSetType(*v, VECSTANDARD)); 475 PetscFunctionReturn(0); 476 } 477 478 /*@C 479 DMInterpolationRestoreVector - Returns a `Vec` which can hold all the interpolated field values 480 481 Collective on ctx 482 483 Input Parameters: 484 + ctx - the context 485 - v - a vector capable of holding the interpolated field values 486 487 Level: intermediate 488 489 .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 490 @*/ 491 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) { 492 PetscFunctionBegin; 493 PetscValidPointer(v, 2); 494 PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 495 PetscCall(VecDestroy(v)); 496 PetscFunctionReturn(0); 497 } 498 499 static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 500 PetscReal v0, J, invJ, detJ; 501 const PetscInt dof = ctx->dof; 502 const PetscScalar *coords; 503 PetscScalar *a; 504 PetscInt p; 505 506 PetscFunctionBegin; 507 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 508 PetscCall(VecGetArray(v, &a)); 509 for (p = 0; p < ctx->n; ++p) { 510 PetscInt c = ctx->cells[p]; 511 PetscScalar *x = NULL; 512 PetscReal xir[1]; 513 PetscInt xSize, comp; 514 515 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ)); 516 PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 517 xir[0] = invJ * PetscRealPart(coords[p] - v0); 518 PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 519 if (2 * dof == xSize) { 520 for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) + x[1 * dof + comp] * xir[0]; 521 } else if (dof == xSize) { 522 for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 523 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input closure size %" PetscInt_FMT " must be either %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 2 * dof, dof); 524 PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 525 } 526 PetscCall(VecRestoreArray(v, &a)); 527 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 528 PetscFunctionReturn(0); 529 } 530 531 static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 532 PetscReal *v0, *J, *invJ, detJ; 533 const PetscScalar *coords; 534 PetscScalar *a; 535 PetscInt p; 536 537 PetscFunctionBegin; 538 PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 539 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 540 PetscCall(VecGetArray(v, &a)); 541 for (p = 0; p < ctx->n; ++p) { 542 PetscInt c = ctx->cells[p]; 543 PetscScalar *x = NULL; 544 PetscReal xi[4]; 545 PetscInt d, f, comp; 546 547 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 548 PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 549 PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 550 for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 551 552 for (d = 0; d < ctx->dim; ++d) { 553 xi[d] = 0.0; 554 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 555 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]; 556 } 557 PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 558 } 559 PetscCall(VecRestoreArray(v, &a)); 560 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 561 PetscCall(PetscFree3(v0, J, invJ)); 562 PetscFunctionReturn(0); 563 } 564 565 static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 566 PetscReal *v0, *J, *invJ, detJ; 567 const PetscScalar *coords; 568 PetscScalar *a; 569 PetscInt p; 570 571 PetscFunctionBegin; 572 PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 573 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 574 PetscCall(VecGetArray(v, &a)); 575 for (p = 0; p < ctx->n; ++p) { 576 PetscInt c = ctx->cells[p]; 577 const PetscInt order[3] = {2, 1, 3}; 578 PetscScalar *x = NULL; 579 PetscReal xi[4]; 580 PetscInt d, f, comp; 581 582 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 583 PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 584 PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 585 for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 586 587 for (d = 0; d < ctx->dim; ++d) { 588 xi[d] = 0.0; 589 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 590 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]; 591 } 592 PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 593 } 594 PetscCall(VecRestoreArray(v, &a)); 595 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 596 PetscCall(PetscFree3(v0, J, invJ)); 597 PetscFunctionReturn(0); 598 } 599 600 static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) { 601 const PetscScalar *vertices = (const PetscScalar *)ctx; 602 const PetscScalar x0 = vertices[0]; 603 const PetscScalar y0 = vertices[1]; 604 const PetscScalar x1 = vertices[2]; 605 const PetscScalar y1 = vertices[3]; 606 const PetscScalar x2 = vertices[4]; 607 const PetscScalar y2 = vertices[5]; 608 const PetscScalar x3 = vertices[6]; 609 const PetscScalar y3 = vertices[7]; 610 const PetscScalar f_1 = x1 - x0; 611 const PetscScalar g_1 = y1 - y0; 612 const PetscScalar f_3 = x3 - x0; 613 const PetscScalar g_3 = y3 - y0; 614 const PetscScalar f_01 = x2 - x1 - x3 + x0; 615 const PetscScalar g_01 = y2 - y1 - y3 + y0; 616 const PetscScalar *ref; 617 PetscScalar *real; 618 619 PetscFunctionBegin; 620 PetscCall(VecGetArrayRead(Xref, &ref)); 621 PetscCall(VecGetArray(Xreal, &real)); 622 { 623 const PetscScalar p0 = ref[0]; 624 const PetscScalar p1 = ref[1]; 625 626 real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 627 real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 628 } 629 PetscCall(PetscLogFlops(28)); 630 PetscCall(VecRestoreArrayRead(Xref, &ref)); 631 PetscCall(VecRestoreArray(Xreal, &real)); 632 PetscFunctionReturn(0); 633 } 634 635 #include <petsc/private/dmimpl.h> 636 static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) { 637 const PetscScalar *vertices = (const PetscScalar *)ctx; 638 const PetscScalar x0 = vertices[0]; 639 const PetscScalar y0 = vertices[1]; 640 const PetscScalar x1 = vertices[2]; 641 const PetscScalar y1 = vertices[3]; 642 const PetscScalar x2 = vertices[4]; 643 const PetscScalar y2 = vertices[5]; 644 const PetscScalar x3 = vertices[6]; 645 const PetscScalar y3 = vertices[7]; 646 const PetscScalar f_01 = x2 - x1 - x3 + x0; 647 const PetscScalar g_01 = y2 - y1 - y3 + y0; 648 const PetscScalar *ref; 649 650 PetscFunctionBegin; 651 PetscCall(VecGetArrayRead(Xref, &ref)); 652 { 653 const PetscScalar x = ref[0]; 654 const PetscScalar y = ref[1]; 655 const PetscInt rows[2] = {0, 1}; 656 PetscScalar values[4]; 657 658 values[0] = (x1 - x0 + f_01 * y) * 0.5; 659 values[1] = (x3 - x0 + f_01 * x) * 0.5; 660 values[2] = (y1 - y0 + g_01 * y) * 0.5; 661 values[3] = (y3 - y0 + g_01 * x) * 0.5; 662 PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES)); 663 } 664 PetscCall(PetscLogFlops(30)); 665 PetscCall(VecRestoreArrayRead(Xref, &ref)); 666 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 667 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 668 PetscFunctionReturn(0); 669 } 670 671 static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 672 DM dmCoord; 673 PetscFE fem = NULL; 674 SNES snes; 675 KSP ksp; 676 PC pc; 677 Vec coordsLocal, r, ref, real; 678 Mat J; 679 PetscTabulation T = NULL; 680 const PetscScalar *coords; 681 PetscScalar *a; 682 PetscReal xir[2] = {0., 0.}; 683 PetscInt Nf, p; 684 const PetscInt dof = ctx->dof; 685 686 PetscFunctionBegin; 687 PetscCall(DMGetNumFields(dm, &Nf)); 688 if (Nf) { 689 PetscObject obj; 690 PetscClassId id; 691 692 PetscCall(DMGetField(dm, 0, NULL, &obj)); 693 PetscCall(PetscObjectGetClassId(obj, &id)); 694 if (id == PETSCFE_CLASSID) { 695 fem = (PetscFE)obj; 696 PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T)); 697 } 698 } 699 PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 700 PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 701 PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 702 PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_")); 703 PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 704 PetscCall(VecSetSizes(r, 2, 2)); 705 PetscCall(VecSetType(r, dm->vectype)); 706 PetscCall(VecDuplicate(r, &ref)); 707 PetscCall(VecDuplicate(r, &real)); 708 PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 709 PetscCall(MatSetSizes(J, 2, 2, 2, 2)); 710 PetscCall(MatSetType(J, MATSEQDENSE)); 711 PetscCall(MatSetUp(J)); 712 PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL)); 713 PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL)); 714 PetscCall(SNESGetKSP(snes, &ksp)); 715 PetscCall(KSPGetPC(ksp, &pc)); 716 PetscCall(PCSetType(pc, PCLU)); 717 PetscCall(SNESSetFromOptions(snes)); 718 719 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 720 PetscCall(VecGetArray(v, &a)); 721 for (p = 0; p < ctx->n; ++p) { 722 PetscScalar *x = NULL, *vertices = NULL; 723 PetscScalar *xi; 724 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 725 726 /* Can make this do all points at once */ 727 PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 728 PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2); 729 PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 730 PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 731 PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 732 PetscCall(VecGetArray(real, &xi)); 733 xi[0] = coords[p * ctx->dim + 0]; 734 xi[1] = coords[p * ctx->dim + 1]; 735 PetscCall(VecRestoreArray(real, &xi)); 736 PetscCall(SNESSolve(snes, real, ref)); 737 PetscCall(VecGetArray(ref, &xi)); 738 xir[0] = PetscRealPart(xi[0]); 739 xir[1] = PetscRealPart(xi[1]); 740 if (4 * dof == xSize) { 741 for (comp = 0; comp < dof; ++comp) 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]; 742 } else if (dof == xSize) { 743 for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 744 } else { 745 PetscInt d; 746 747 PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE"); 748 xir[0] = 2.0 * xir[0] - 1.0; 749 xir[1] = 2.0 * xir[1] - 1.0; 750 PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T)); 751 for (comp = 0; comp < dof; ++comp) { 752 a[p * dof + comp] = 0.0; 753 for (d = 0; d < xSize / dof; ++d) a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp]; 754 } 755 } 756 PetscCall(VecRestoreArray(ref, &xi)); 757 PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 758 PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 759 } 760 PetscCall(PetscTabulationDestroy(&T)); 761 PetscCall(VecRestoreArray(v, &a)); 762 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 763 764 PetscCall(SNESDestroy(&snes)); 765 PetscCall(VecDestroy(&r)); 766 PetscCall(VecDestroy(&ref)); 767 PetscCall(VecDestroy(&real)); 768 PetscCall(MatDestroy(&J)); 769 PetscFunctionReturn(0); 770 } 771 772 static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) { 773 const PetscScalar *vertices = (const PetscScalar *)ctx; 774 const PetscScalar x0 = vertices[0]; 775 const PetscScalar y0 = vertices[1]; 776 const PetscScalar z0 = vertices[2]; 777 const PetscScalar x1 = vertices[9]; 778 const PetscScalar y1 = vertices[10]; 779 const PetscScalar z1 = vertices[11]; 780 const PetscScalar x2 = vertices[6]; 781 const PetscScalar y2 = vertices[7]; 782 const PetscScalar z2 = vertices[8]; 783 const PetscScalar x3 = vertices[3]; 784 const PetscScalar y3 = vertices[4]; 785 const PetscScalar z3 = vertices[5]; 786 const PetscScalar x4 = vertices[12]; 787 const PetscScalar y4 = vertices[13]; 788 const PetscScalar z4 = vertices[14]; 789 const PetscScalar x5 = vertices[15]; 790 const PetscScalar y5 = vertices[16]; 791 const PetscScalar z5 = vertices[17]; 792 const PetscScalar x6 = vertices[18]; 793 const PetscScalar y6 = vertices[19]; 794 const PetscScalar z6 = vertices[20]; 795 const PetscScalar x7 = vertices[21]; 796 const PetscScalar y7 = vertices[22]; 797 const PetscScalar z7 = vertices[23]; 798 const PetscScalar f_1 = x1 - x0; 799 const PetscScalar g_1 = y1 - y0; 800 const PetscScalar h_1 = z1 - z0; 801 const PetscScalar f_3 = x3 - x0; 802 const PetscScalar g_3 = y3 - y0; 803 const PetscScalar h_3 = z3 - z0; 804 const PetscScalar f_4 = x4 - x0; 805 const PetscScalar g_4 = y4 - y0; 806 const PetscScalar h_4 = z4 - z0; 807 const PetscScalar f_01 = x2 - x1 - x3 + x0; 808 const PetscScalar g_01 = y2 - y1 - y3 + y0; 809 const PetscScalar h_01 = z2 - z1 - z3 + z0; 810 const PetscScalar f_12 = x7 - x3 - x4 + x0; 811 const PetscScalar g_12 = y7 - y3 - y4 + y0; 812 const PetscScalar h_12 = z7 - z3 - z4 + z0; 813 const PetscScalar f_02 = x5 - x1 - x4 + x0; 814 const PetscScalar g_02 = y5 - y1 - y4 + y0; 815 const PetscScalar h_02 = z5 - z1 - z4 + z0; 816 const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 817 const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 818 const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 819 const PetscScalar *ref; 820 PetscScalar *real; 821 822 PetscFunctionBegin; 823 PetscCall(VecGetArrayRead(Xref, &ref)); 824 PetscCall(VecGetArray(Xreal, &real)); 825 { 826 const PetscScalar p0 = ref[0]; 827 const PetscScalar p1 = ref[1]; 828 const PetscScalar p2 = ref[2]; 829 830 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; 831 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; 832 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; 833 } 834 PetscCall(PetscLogFlops(114)); 835 PetscCall(VecRestoreArrayRead(Xref, &ref)); 836 PetscCall(VecRestoreArray(Xreal, &real)); 837 PetscFunctionReturn(0); 838 } 839 840 static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) { 841 const PetscScalar *vertices = (const PetscScalar *)ctx; 842 const PetscScalar x0 = vertices[0]; 843 const PetscScalar y0 = vertices[1]; 844 const PetscScalar z0 = vertices[2]; 845 const PetscScalar x1 = vertices[9]; 846 const PetscScalar y1 = vertices[10]; 847 const PetscScalar z1 = vertices[11]; 848 const PetscScalar x2 = vertices[6]; 849 const PetscScalar y2 = vertices[7]; 850 const PetscScalar z2 = vertices[8]; 851 const PetscScalar x3 = vertices[3]; 852 const PetscScalar y3 = vertices[4]; 853 const PetscScalar z3 = vertices[5]; 854 const PetscScalar x4 = vertices[12]; 855 const PetscScalar y4 = vertices[13]; 856 const PetscScalar z4 = vertices[14]; 857 const PetscScalar x5 = vertices[15]; 858 const PetscScalar y5 = vertices[16]; 859 const PetscScalar z5 = vertices[17]; 860 const PetscScalar x6 = vertices[18]; 861 const PetscScalar y6 = vertices[19]; 862 const PetscScalar z6 = vertices[20]; 863 const PetscScalar x7 = vertices[21]; 864 const PetscScalar y7 = vertices[22]; 865 const PetscScalar z7 = vertices[23]; 866 const PetscScalar f_xy = x2 - x1 - x3 + x0; 867 const PetscScalar g_xy = y2 - y1 - y3 + y0; 868 const PetscScalar h_xy = z2 - z1 - z3 + z0; 869 const PetscScalar f_yz = x7 - x3 - x4 + x0; 870 const PetscScalar g_yz = y7 - y3 - y4 + y0; 871 const PetscScalar h_yz = z7 - z3 - z4 + z0; 872 const PetscScalar f_xz = x5 - x1 - x4 + x0; 873 const PetscScalar g_xz = y5 - y1 - y4 + y0; 874 const PetscScalar h_xz = z5 - z1 - z4 + z0; 875 const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 876 const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 877 const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 878 const PetscScalar *ref; 879 880 PetscFunctionBegin; 881 PetscCall(VecGetArrayRead(Xref, &ref)); 882 { 883 const PetscScalar x = ref[0]; 884 const PetscScalar y = ref[1]; 885 const PetscScalar z = ref[2]; 886 const PetscInt rows[3] = {0, 1, 2}; 887 PetscScalar values[9]; 888 889 values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0; 890 values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0; 891 values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0; 892 values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0; 893 values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0; 894 values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0; 895 values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0; 896 values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0; 897 values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0; 898 899 PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES)); 900 } 901 PetscCall(PetscLogFlops(152)); 902 PetscCall(VecRestoreArrayRead(Xref, &ref)); 903 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 904 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 905 PetscFunctionReturn(0); 906 } 907 908 static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 909 DM dmCoord; 910 SNES snes; 911 KSP ksp; 912 PC pc; 913 Vec coordsLocal, r, ref, real; 914 Mat J; 915 const PetscScalar *coords; 916 PetscScalar *a; 917 PetscInt p; 918 919 PetscFunctionBegin; 920 PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 921 PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 922 PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 923 PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_")); 924 PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 925 PetscCall(VecSetSizes(r, 3, 3)); 926 PetscCall(VecSetType(r, dm->vectype)); 927 PetscCall(VecDuplicate(r, &ref)); 928 PetscCall(VecDuplicate(r, &real)); 929 PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 930 PetscCall(MatSetSizes(J, 3, 3, 3, 3)); 931 PetscCall(MatSetType(J, MATSEQDENSE)); 932 PetscCall(MatSetUp(J)); 933 PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL)); 934 PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL)); 935 PetscCall(SNESGetKSP(snes, &ksp)); 936 PetscCall(KSPGetPC(ksp, &pc)); 937 PetscCall(PCSetType(pc, PCLU)); 938 PetscCall(SNESSetFromOptions(snes)); 939 940 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 941 PetscCall(VecGetArray(v, &a)); 942 for (p = 0; p < ctx->n; ++p) { 943 PetscScalar *x = NULL, *vertices = NULL; 944 PetscScalar *xi; 945 PetscReal xir[3]; 946 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 947 948 /* Can make this do all points at once */ 949 PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 950 PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3); 951 PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 952 PetscCheck((8 * ctx->dof == xSize) || (ctx->dof == xSize), ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input closure size %" PetscInt_FMT " should be %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 8 * ctx->dof, ctx->dof); 953 PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 954 PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 955 PetscCall(VecGetArray(real, &xi)); 956 xi[0] = coords[p * ctx->dim + 0]; 957 xi[1] = coords[p * ctx->dim + 1]; 958 xi[2] = coords[p * ctx->dim + 2]; 959 PetscCall(VecRestoreArray(real, &xi)); 960 PetscCall(SNESSolve(snes, real, ref)); 961 PetscCall(VecGetArray(ref, &xi)); 962 xir[0] = PetscRealPart(xi[0]); 963 xir[1] = PetscRealPart(xi[1]); 964 xir[2] = PetscRealPart(xi[2]); 965 if (8 * ctx->dof == xSize) { 966 for (comp = 0; comp < ctx->dof; ++comp) { 967 a[p * ctx->dof + comp] = x[0 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * (1 - xir[2]) + x[3 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * (1 - xir[2]) + x[2 * ctx->dof + comp] * xir[0] * xir[1] * (1 - xir[2]) + x[1 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * (1 - xir[2]) + 968 x[4 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * xir[2] + x[5 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * xir[2] + x[6 * ctx->dof + comp] * xir[0] * xir[1] * xir[2] + x[7 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * xir[2]; 969 } 970 } else { 971 for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 972 } 973 PetscCall(VecRestoreArray(ref, &xi)); 974 PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 975 PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 976 } 977 PetscCall(VecRestoreArray(v, &a)); 978 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 979 980 PetscCall(SNESDestroy(&snes)); 981 PetscCall(VecDestroy(&r)); 982 PetscCall(VecDestroy(&ref)); 983 PetscCall(VecDestroy(&real)); 984 PetscCall(MatDestroy(&J)); 985 PetscFunctionReturn(0); 986 } 987 988 /*@C 989 DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points. 990 991 Input Parameters: 992 + ctx - The `DMInterpolationInfo` context 993 . dm - The `DM` 994 - x - The local vector containing the field to be interpolated 995 996 Output Parameters: 997 . v - The vector containing the interpolated values 998 999 Note: 1000 A suitable v can be obtained using `DMInterpolationGetVector()`. 1001 1002 Level: beginner 1003 1004 .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 1005 @*/ 1006 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) { 1007 PetscDS ds; 1008 PetscInt n, p, Nf, field; 1009 PetscBool useDS = PETSC_FALSE; 1010 1011 PetscFunctionBegin; 1012 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1013 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 1014 PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 1015 PetscCall(VecGetLocalSize(v, &n)); 1016 PetscCheck(n == ctx->n * ctx->dof, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %" PetscInt_FMT " should be %" PetscInt_FMT, n, ctx->n * ctx->dof); 1017 if (!n) PetscFunctionReturn(0); 1018 PetscCall(DMGetDS(dm, &ds)); 1019 if (ds) { 1020 useDS = PETSC_TRUE; 1021 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1022 for (field = 0; field < Nf; ++field) { 1023 PetscObject obj; 1024 PetscClassId id; 1025 1026 PetscCall(PetscDSGetDiscretization(ds, field, &obj)); 1027 PetscCall(PetscObjectGetClassId(obj, &id)); 1028 if (id != PETSCFE_CLASSID) { 1029 useDS = PETSC_FALSE; 1030 break; 1031 } 1032 } 1033 } 1034 if (useDS) { 1035 const PetscScalar *coords; 1036 PetscScalar *interpolant; 1037 PetscInt cdim, d; 1038 1039 PetscCall(DMGetCoordinateDim(dm, &cdim)); 1040 PetscCall(VecGetArrayRead(ctx->coords, &coords)); 1041 PetscCall(VecGetArrayWrite(v, &interpolant)); 1042 for (p = 0; p < ctx->n; ++p) { 1043 PetscReal pcoords[3], xi[3]; 1044 PetscScalar *xa = NULL; 1045 PetscInt coff = 0, foff = 0, clSize; 1046 1047 if (ctx->cells[p] < 0) continue; 1048 for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]); 1049 PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi)); 1050 PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 1051 for (field = 0; field < Nf; ++field) { 1052 PetscTabulation T; 1053 PetscFE fe; 1054 1055 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 1056 PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T)); 1057 { 1058 const PetscReal *basis = T->T[0]; 1059 const PetscInt Nb = T->Nb; 1060 const PetscInt Nc = T->Nc; 1061 PetscInt f, fc; 1062 1063 for (fc = 0; fc < Nc; ++fc) { 1064 interpolant[p * ctx->dof + coff + fc] = 0.0; 1065 for (f = 0; f < Nb; ++f) interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc]; 1066 } 1067 coff += Nc; 1068 foff += Nb; 1069 } 1070 PetscCall(PetscTabulationDestroy(&T)); 1071 } 1072 PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 1073 PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof); 1074 PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize); 1075 } 1076 PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 1077 PetscCall(VecRestoreArrayWrite(v, &interpolant)); 1078 } else { 1079 DMPolytopeType ct; 1080 1081 /* TODO Check each cell individually */ 1082 PetscCall(DMPlexGetCellType(dm, ctx->cells[0], &ct)); 1083 switch (ct) { 1084 case DM_POLYTOPE_SEGMENT: PetscCall(DMInterpolate_Segment_Private(ctx, dm, x, v)); break; 1085 case DM_POLYTOPE_TRIANGLE: PetscCall(DMInterpolate_Triangle_Private(ctx, dm, x, v)); break; 1086 case DM_POLYTOPE_QUADRILATERAL: PetscCall(DMInterpolate_Quad_Private(ctx, dm, x, v)); break; 1087 case DM_POLYTOPE_TETRAHEDRON: PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, x, v)); break; 1088 case DM_POLYTOPE_HEXAHEDRON: PetscCall(DMInterpolate_Hex_Private(ctx, dm, x, v)); break; 1089 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]); 1090 } 1091 } 1092 PetscFunctionReturn(0); 1093 } 1094 1095 /*@C 1096 DMInterpolationDestroy - Destroys a `DMInterpolationInfo` context 1097 1098 Collective on ctx 1099 1100 Input Parameter: 1101 . ctx - the context 1102 1103 Level: beginner 1104 1105 .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 1106 @*/ 1107 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) { 1108 PetscFunctionBegin; 1109 PetscValidPointer(ctx, 1); 1110 PetscCall(VecDestroy(&(*ctx)->coords)); 1111 PetscCall(PetscFree((*ctx)->points)); 1112 PetscCall(PetscFree((*ctx)->cells)); 1113 PetscCall(PetscFree(*ctx)); 1114 *ctx = NULL; 1115 PetscFunctionReturn(0); 1116 } 1117 1118 /*@C 1119 SNESMonitorFields - Monitors the residual for each field separately 1120 1121 Collective on snes 1122 1123 Input Parameters: 1124 + snes - the `SNES` context 1125 . its - iteration number 1126 . fgnorm - 2-norm of residual 1127 - vf - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII` 1128 1129 Note: 1130 This routine prints the residual norm at each iteration. 1131 1132 Level: intermediate 1133 1134 .seealso: `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()` 1135 @*/ 1136 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) { 1137 PetscViewer viewer = vf->viewer; 1138 Vec res; 1139 DM dm; 1140 PetscSection s; 1141 const PetscScalar *r; 1142 PetscReal *lnorms, *norms; 1143 PetscInt numFields, f, pStart, pEnd, p; 1144 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 1147 PetscCall(SNESGetFunction(snes, &res, NULL, NULL)); 1148 PetscCall(SNESGetDM(snes, &dm)); 1149 PetscCall(DMGetLocalSection(dm, &s)); 1150 PetscCall(PetscSectionGetNumFields(s, &numFields)); 1151 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1152 PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms)); 1153 PetscCall(VecGetArrayRead(res, &r)); 1154 for (p = pStart; p < pEnd; ++p) { 1155 for (f = 0; f < numFields; ++f) { 1156 PetscInt fdof, foff, d; 1157 1158 PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1159 PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff)); 1160 for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d])); 1161 } 1162 } 1163 PetscCall(VecRestoreArrayRead(res, &r)); 1164 PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm))); 1165 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 1166 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 1167 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm)); 1168 for (f = 0; f < numFields; ++f) { 1169 if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 1170 PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f]))); 1171 } 1172 PetscCall(PetscViewerASCIIPrintf(viewer, "]\n")); 1173 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 1174 PetscCall(PetscViewerPopFormat(viewer)); 1175 PetscCall(PetscFree2(lnorms, norms)); 1176 PetscFunctionReturn(0); 1177 } 1178 1179 /********************* Residual Computation **************************/ 1180 1181 PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS) { 1182 PetscInt depth; 1183 1184 PetscFunctionBegin; 1185 PetscCall(DMPlexGetDepth(plex, &depth)); 1186 PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS)); 1187 if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS)); 1188 PetscFunctionReturn(0); 1189 } 1190 1191 /*@ 1192 DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user 1193 1194 Input Parameters: 1195 + dm - The mesh 1196 . X - Local solution 1197 - user - The user context 1198 1199 Output Parameter: 1200 . F - Local output vector 1201 1202 Note: 1203 The residual is summed into F; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in F is intentional. 1204 1205 Level: developer 1206 1207 .seealso: `DM`, `DMPlexComputeJacobianAction()` 1208 @*/ 1209 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) { 1210 DM plex; 1211 IS allcellIS; 1212 PetscInt Nds, s; 1213 1214 PetscFunctionBegin; 1215 PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1216 PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 1217 PetscCall(DMGetNumDS(dm, &Nds)); 1218 for (s = 0; s < Nds; ++s) { 1219 PetscDS ds; 1220 IS cellIS; 1221 PetscFormKey key; 1222 1223 PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds)); 1224 key.value = 0; 1225 key.field = 0; 1226 key.part = 0; 1227 if (!key.label) { 1228 PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1229 cellIS = allcellIS; 1230 } else { 1231 IS pointIS; 1232 1233 key.value = 1; 1234 PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 1235 PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 1236 PetscCall(ISDestroy(&pointIS)); 1237 } 1238 PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 1239 PetscCall(ISDestroy(&cellIS)); 1240 } 1241 PetscCall(ISDestroy(&allcellIS)); 1242 PetscCall(DMDestroy(&plex)); 1243 PetscFunctionReturn(0); 1244 } 1245 1246 PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user) { 1247 DM plex; 1248 IS allcellIS; 1249 PetscInt Nds, s; 1250 1251 PetscFunctionBegin; 1252 PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1253 PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 1254 PetscCall(DMGetNumDS(dm, &Nds)); 1255 for (s = 0; s < Nds; ++s) { 1256 PetscDS ds; 1257 DMLabel label; 1258 IS cellIS; 1259 1260 PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds)); 1261 { 1262 PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1}; 1263 PetscWeakForm wf; 1264 PetscInt Nm = 2, m, Nk = 0, k, kp, off = 0; 1265 PetscFormKey *reskeys; 1266 1267 /* Get unique residual keys */ 1268 for (m = 0; m < Nm; ++m) { 1269 PetscInt Nkm; 1270 PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm)); 1271 Nk += Nkm; 1272 } 1273 PetscCall(PetscMalloc1(Nk, &reskeys)); 1274 for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys)); 1275 PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 1276 PetscCall(PetscFormKeySort(Nk, reskeys)); 1277 for (k = 0, kp = 1; kp < Nk; ++kp) { 1278 if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) { 1279 ++k; 1280 if (kp != k) reskeys[k] = reskeys[kp]; 1281 } 1282 } 1283 Nk = k; 1284 1285 PetscCall(PetscDSGetWeakForm(ds, &wf)); 1286 for (k = 0; k < Nk; ++k) { 1287 DMLabel label = reskeys[k].label; 1288 PetscInt val = reskeys[k].value; 1289 1290 if (!label) { 1291 PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1292 cellIS = allcellIS; 1293 } else { 1294 IS pointIS; 1295 1296 PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 1297 PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 1298 PetscCall(ISDestroy(&pointIS)); 1299 } 1300 PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 1301 PetscCall(ISDestroy(&cellIS)); 1302 } 1303 PetscCall(PetscFree(reskeys)); 1304 } 1305 } 1306 PetscCall(ISDestroy(&allcellIS)); 1307 PetscCall(DMDestroy(&plex)); 1308 PetscFunctionReturn(0); 1309 } 1310 1311 /*@ 1312 DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X 1313 1314 Input Parameters: 1315 + dm - The mesh 1316 - user - The user context 1317 1318 Output Parameter: 1319 . X - Local solution 1320 1321 Level: developer 1322 1323 .seealso: `DMPLEX`, `DMPlexComputeJacobianAction()` 1324 @*/ 1325 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) { 1326 DM plex; 1327 1328 PetscFunctionBegin; 1329 PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1330 PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL)); 1331 PetscCall(DMDestroy(&plex)); 1332 PetscFunctionReturn(0); 1333 } 1334 1335 /*@ 1336 DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y 1337 1338 Input Parameters: 1339 + dm - The `DM` 1340 . X - Local solution vector 1341 . Y - Local input vector 1342 - user - The user context 1343 1344 Output Parameter: 1345 . F - lcoal output vector 1346 1347 Level: developer 1348 1349 Notes: 1350 Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly. 1351 1352 .seealso: `DM`, ``DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()` 1353 @*/ 1354 PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user) { 1355 DM plex; 1356 IS allcellIS; 1357 PetscInt Nds, s; 1358 1359 PetscFunctionBegin; 1360 PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1361 PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 1362 PetscCall(DMGetNumDS(dm, &Nds)); 1363 for (s = 0; s < Nds; ++s) { 1364 PetscDS ds; 1365 DMLabel label; 1366 IS cellIS; 1367 1368 PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds)); 1369 { 1370 PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3}; 1371 PetscWeakForm wf; 1372 PetscInt Nm = 4, m, Nk = 0, k, kp, off = 0; 1373 PetscFormKey *jackeys; 1374 1375 /* Get unique Jacobian keys */ 1376 for (m = 0; m < Nm; ++m) { 1377 PetscInt Nkm; 1378 PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm)); 1379 Nk += Nkm; 1380 } 1381 PetscCall(PetscMalloc1(Nk, &jackeys)); 1382 for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys)); 1383 PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 1384 PetscCall(PetscFormKeySort(Nk, jackeys)); 1385 for (k = 0, kp = 1; kp < Nk; ++kp) { 1386 if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) { 1387 ++k; 1388 if (kp != k) jackeys[k] = jackeys[kp]; 1389 } 1390 } 1391 Nk = k; 1392 1393 PetscCall(PetscDSGetWeakForm(ds, &wf)); 1394 for (k = 0; k < Nk; ++k) { 1395 DMLabel label = jackeys[k].label; 1396 PetscInt val = jackeys[k].value; 1397 1398 if (!label) { 1399 PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1400 cellIS = allcellIS; 1401 } else { 1402 IS pointIS; 1403 1404 PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 1405 PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 1406 PetscCall(ISDestroy(&pointIS)); 1407 } 1408 PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user)); 1409 PetscCall(ISDestroy(&cellIS)); 1410 } 1411 PetscCall(PetscFree(jackeys)); 1412 } 1413 } 1414 PetscCall(ISDestroy(&allcellIS)); 1415 PetscCall(DMDestroy(&plex)); 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /*@ 1420 DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 1421 1422 Input Parameters: 1423 + dm - The mesh 1424 . X - Local input vector 1425 - user - The user context 1426 1427 Output Parameter: 1428 . Jac - Jacobian matrix 1429 1430 Note: 1431 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1432 like a GPU, or vectorize on a multicore machine. 1433 1434 Level: developer 1435 1436 .seealso: `DMPLEX`, `Mat` 1437 @*/ 1438 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user) { 1439 DM plex; 1440 IS allcellIS; 1441 PetscBool hasJac, hasPrec; 1442 PetscInt Nds, s; 1443 1444 PetscFunctionBegin; 1445 PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1446 PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 1447 PetscCall(DMGetNumDS(dm, &Nds)); 1448 for (s = 0; s < Nds; ++s) { 1449 PetscDS ds; 1450 IS cellIS; 1451 PetscFormKey key; 1452 1453 PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds)); 1454 key.value = 0; 1455 key.field = 0; 1456 key.part = 0; 1457 if (!key.label) { 1458 PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1459 cellIS = allcellIS; 1460 } else { 1461 IS pointIS; 1462 1463 key.value = 1; 1464 PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 1465 PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 1466 PetscCall(ISDestroy(&pointIS)); 1467 } 1468 if (!s) { 1469 PetscCall(PetscDSHasJacobian(ds, &hasJac)); 1470 PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 1471 if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac)); 1472 PetscCall(MatZeroEntries(JacP)); 1473 } 1474 PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user)); 1475 PetscCall(ISDestroy(&cellIS)); 1476 } 1477 PetscCall(ISDestroy(&allcellIS)); 1478 PetscCall(DMDestroy(&plex)); 1479 PetscFunctionReturn(0); 1480 } 1481 1482 struct _DMSNESJacobianMFCtx { 1483 DM dm; 1484 Vec X; 1485 void *ctx; 1486 }; 1487 1488 static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A) { 1489 struct _DMSNESJacobianMFCtx *ctx; 1490 1491 PetscFunctionBegin; 1492 PetscCall(MatShellGetContext(A, &ctx)); 1493 PetscCall(MatShellSetContext(A, NULL)); 1494 PetscCall(DMDestroy(&ctx->dm)); 1495 PetscCall(VecDestroy(&ctx->X)); 1496 PetscCall(PetscFree(ctx)); 1497 PetscFunctionReturn(0); 1498 } 1499 1500 static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z) { 1501 struct _DMSNESJacobianMFCtx *ctx; 1502 1503 PetscFunctionBegin; 1504 PetscCall(MatShellGetContext(A, &ctx)); 1505 PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx)); 1506 PetscFunctionReturn(0); 1507 } 1508 1509 /*@ 1510 DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free 1511 1512 Collective on dm 1513 1514 Input Parameters: 1515 + dm - The `DM` 1516 . X - The evaluation point for the Jacobian 1517 - user - A user context, or NULL 1518 1519 Output Parameter: 1520 . J - The `Mat` 1521 1522 Level: advanced 1523 1524 Note: 1525 Vec X is kept in `Mat` J, so updating X then updates the evaluation point. 1526 1527 .seealso: `DM`, `DMSNESComputeJacobianAction()` 1528 @*/ 1529 PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J) { 1530 struct _DMSNESJacobianMFCtx *ctx; 1531 PetscInt n, N; 1532 1533 PetscFunctionBegin; 1534 PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 1535 PetscCall(MatSetType(*J, MATSHELL)); 1536 PetscCall(VecGetLocalSize(X, &n)); 1537 PetscCall(VecGetSize(X, &N)); 1538 PetscCall(MatSetSizes(*J, n, n, N, N)); 1539 PetscCall(PetscObjectReference((PetscObject)dm)); 1540 PetscCall(PetscObjectReference((PetscObject)X)); 1541 PetscCall(PetscMalloc1(1, &ctx)); 1542 ctx->dm = dm; 1543 ctx->X = X; 1544 ctx->ctx = user; 1545 PetscCall(MatShellSetContext(*J, ctx)); 1546 PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private)); 1547 PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private)); 1548 PetscFunctionReturn(0); 1549 } 1550 1551 /* 1552 MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context. 1553 1554 Input Parameters: 1555 + X - SNES linearization point 1556 . ovl - index set of overlapping subdomains 1557 1558 Output Parameter: 1559 . J - unassembled (Neumann) local matrix 1560 1561 Level: intermediate 1562 1563 .seealso: `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()` 1564 */ 1565 static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) { 1566 SNES snes; 1567 Mat pJ; 1568 DM ovldm, origdm; 1569 DMSNES sdm; 1570 PetscErrorCode (*bfun)(DM, Vec, void *); 1571 PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *); 1572 void *bctx, *jctx; 1573 1574 PetscFunctionBegin; 1575 PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ)); 1576 PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat"); 1577 PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm)); 1578 PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM"); 1579 PetscCall(MatGetDM(pJ, &ovldm)); 1580 PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx)); 1581 PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx)); 1582 PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx)); 1583 PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx)); 1584 PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes)); 1585 if (!snes) { 1586 PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes)); 1587 PetscCall(SNESSetDM(snes, ovldm)); 1588 PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes)); 1589 PetscCall(PetscObjectDereference((PetscObject)snes)); 1590 } 1591 PetscCall(DMGetDMSNES(ovldm, &sdm)); 1592 PetscCall(VecLockReadPush(X)); 1593 { 1594 void *ctx; 1595 PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *); 1596 PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx)); 1597 PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx)); 1598 } 1599 PetscCall(VecLockReadPop(X)); 1600 /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */ 1601 { 1602 Mat locpJ; 1603 1604 PetscCall(MatISGetLocalMat(pJ, &locpJ)); 1605 PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN)); 1606 } 1607 PetscFunctionReturn(0); 1608 } 1609 1610 /*@ 1611 DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, residual, and Jacobian. 1612 1613 Input Parameters: 1614 + dm - The `DM` object 1615 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see `PetscDSAddBoundary()`) 1616 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see `PetscDSSetResidual()`) 1617 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see `PetscDSSetJacobian()`) 1618 1619 Level: developer 1620 1621 .seealso: `DMPLEX`, `SNES` 1622 @*/ 1623 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) { 1624 PetscFunctionBegin; 1625 PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, boundaryctx)); 1626 PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, residualctx)); 1627 PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, jacobianctx)); 1628 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex)); 1629 PetscFunctionReturn(0); 1630 } 1631 1632 /*@C 1633 DMSNESCheckDiscretization - Check the discretization error of the exact solution 1634 1635 Input Parameters: 1636 + snes - the `SNES` object 1637 . dm - the `DM` 1638 . t - the time 1639 . u - a `DM` vector 1640 - tol - A tolerance for the check, or -1 to print the results instead 1641 1642 Output Parameters: 1643 . error - An array which holds the discretization error in each field, or NULL 1644 1645 Note: 1646 The user must call `PetscDSSetExactSolution()` beforehand 1647 1648 Level: developer 1649 1650 .seealso: `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`, `PetscDSSetExactSolution()` 1651 @*/ 1652 PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[]) { 1653 PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 1654 void **ectxs; 1655 PetscReal *err; 1656 MPI_Comm comm; 1657 PetscInt Nf, f; 1658 1659 PetscFunctionBegin; 1660 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1661 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1662 PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 1663 if (error) PetscValidRealPointer(error, 6); 1664 1665 PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 1666 PetscCall(VecViewFromOptions(u, NULL, "-vec_view")); 1667 1668 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 1669 PetscCall(DMGetNumFields(dm, &Nf)); 1670 PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err)); 1671 { 1672 PetscInt Nds, s; 1673 1674 PetscCall(DMGetNumDS(dm, &Nds)); 1675 for (s = 0; s < Nds; ++s) { 1676 PetscDS ds; 1677 DMLabel label; 1678 IS fieldIS; 1679 const PetscInt *fields; 1680 PetscInt dsNf, f; 1681 1682 PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds)); 1683 PetscCall(PetscDSGetNumFields(ds, &dsNf)); 1684 PetscCall(ISGetIndices(fieldIS, &fields)); 1685 for (f = 0; f < dsNf; ++f) { 1686 const PetscInt field = fields[f]; 1687 PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field])); 1688 } 1689 PetscCall(ISRestoreIndices(fieldIS, &fields)); 1690 } 1691 } 1692 if (Nf > 1) { 1693 PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err)); 1694 if (tol >= 0.0) { 1695 for (f = 0; f < Nf; ++f) PetscCheck(err[f] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %" PetscInt_FMT " exceeds tolerance %g", (double)err[f], f, (double)tol); 1696 } else if (error) { 1697 for (f = 0; f < Nf; ++f) error[f] = err[f]; 1698 } else { 1699 PetscCall(PetscPrintf(comm, "L_2 Error: [")); 1700 for (f = 0; f < Nf; ++f) { 1701 if (f) PetscCall(PetscPrintf(comm, ", ")); 1702 PetscCall(PetscPrintf(comm, "%g", (double)err[f])); 1703 } 1704 PetscCall(PetscPrintf(comm, "]\n")); 1705 } 1706 } else { 1707 PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0])); 1708 if (tol >= 0.0) { 1709 PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol); 1710 } else if (error) { 1711 error[0] = err[0]; 1712 } else { 1713 PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0])); 1714 } 1715 } 1716 PetscCall(PetscFree3(exacts, ectxs, err)); 1717 PetscFunctionReturn(0); 1718 } 1719 1720 /*@C 1721 DMSNESCheckResidual - Check the residual of the exact solution 1722 1723 Input Parameters: 1724 + snes - the `SNES` object 1725 . dm - the `DM` 1726 . u - a `DM` vector 1727 - tol - A tolerance for the check, or -1 to print the results instead 1728 1729 Output Parameter: 1730 . residual - The residual norm of the exact solution, or NULL 1731 1732 Level: developer 1733 1734 .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()` 1735 @*/ 1736 PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) { 1737 MPI_Comm comm; 1738 Vec r; 1739 PetscReal res; 1740 1741 PetscFunctionBegin; 1742 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1743 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1744 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 1745 if (residual) PetscValidRealPointer(residual, 5); 1746 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 1747 PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 1748 PetscCall(VecDuplicate(u, &r)); 1749 PetscCall(SNESComputeFunction(snes, u, r)); 1750 PetscCall(VecNorm(r, NORM_2, &res)); 1751 if (tol >= 0.0) { 1752 PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol); 1753 } else if (residual) { 1754 *residual = res; 1755 } else { 1756 PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res)); 1757 PetscCall(VecChop(r, 1.0e-10)); 1758 PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual")); 1759 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_")); 1760 PetscCall(VecViewFromOptions(r, NULL, "-vec_view")); 1761 } 1762 PetscCall(VecDestroy(&r)); 1763 PetscFunctionReturn(0); 1764 } 1765 1766 /*@C 1767 DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 1768 1769 Input Parameters: 1770 + snes - the `SNES` object 1771 . dm - the `DM` 1772 . u - a `DM` vector 1773 - tol - A tolerance for the check, or -1 to print the results instead 1774 1775 Output Parameters: 1776 + isLinear - Flag indicaing that the function looks linear, or NULL 1777 - convRate - The rate of convergence of the linear model, or NULL 1778 1779 Level: developer 1780 1781 .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()` 1782 @*/ 1783 PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) { 1784 MPI_Comm comm; 1785 PetscDS ds; 1786 Mat J, M; 1787 MatNullSpace nullspace; 1788 PetscReal slope, intercept; 1789 PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 1790 1791 PetscFunctionBegin; 1792 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1793 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1794 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 1795 if (isLinear) PetscValidBoolPointer(isLinear, 5); 1796 if (convRate) PetscValidRealPointer(convRate, 6); 1797 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 1798 PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 1799 /* Create and view matrices */ 1800 PetscCall(DMCreateMatrix(dm, &J)); 1801 PetscCall(DMGetDS(dm, &ds)); 1802 PetscCall(PetscDSHasJacobian(ds, &hasJac)); 1803 PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 1804 if (hasJac && hasPrec) { 1805 PetscCall(DMCreateMatrix(dm, &M)); 1806 PetscCall(SNESComputeJacobian(snes, u, J, M)); 1807 PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix")); 1808 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_")); 1809 PetscCall(MatViewFromOptions(M, NULL, "-mat_view")); 1810 PetscCall(MatDestroy(&M)); 1811 } else { 1812 PetscCall(SNESComputeJacobian(snes, u, J, J)); 1813 } 1814 PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 1815 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_")); 1816 PetscCall(MatViewFromOptions(J, NULL, "-mat_view")); 1817 /* Check nullspace */ 1818 PetscCall(MatGetNullSpace(J, &nullspace)); 1819 if (nullspace) { 1820 PetscBool isNull; 1821 PetscCall(MatNullSpaceTest(nullspace, J, &isNull)); 1822 PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 1823 } 1824 /* Taylor test */ 1825 { 1826 PetscRandom rand; 1827 Vec du, uhat, r, rhat, df; 1828 PetscReal h; 1829 PetscReal *es, *hs, *errors; 1830 PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 1831 PetscInt Nv, v; 1832 1833 /* Choose a perturbation direction */ 1834 PetscCall(PetscRandomCreate(comm, &rand)); 1835 PetscCall(VecDuplicate(u, &du)); 1836 PetscCall(VecSetRandom(du, rand)); 1837 PetscCall(PetscRandomDestroy(&rand)); 1838 PetscCall(VecDuplicate(u, &df)); 1839 PetscCall(MatMult(J, du, df)); 1840 /* Evaluate residual at u, F(u), save in vector r */ 1841 PetscCall(VecDuplicate(u, &r)); 1842 PetscCall(SNESComputeFunction(snes, u, r)); 1843 /* Look at the convergence of our Taylor approximation as we approach u */ 1844 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) 1845 ; 1846 PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors)); 1847 PetscCall(VecDuplicate(u, &uhat)); 1848 PetscCall(VecDuplicate(u, &rhat)); 1849 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 1850 PetscCall(VecWAXPY(uhat, h, du, u)); 1851 /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */ 1852 PetscCall(SNESComputeFunction(snes, uhat, rhat)); 1853 PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df)); 1854 PetscCall(VecNorm(rhat, NORM_2, &errors[Nv])); 1855 1856 es[Nv] = PetscLog10Real(errors[Nv]); 1857 hs[Nv] = PetscLog10Real(h); 1858 } 1859 PetscCall(VecDestroy(&uhat)); 1860 PetscCall(VecDestroy(&rhat)); 1861 PetscCall(VecDestroy(&df)); 1862 PetscCall(VecDestroy(&r)); 1863 PetscCall(VecDestroy(&du)); 1864 for (v = 0; v < Nv; ++v) { 1865 if ((tol >= 0) && (errors[v] > tol)) break; 1866 else if (errors[v] > PETSC_SMALL) break; 1867 } 1868 if (v == Nv) isLin = PETSC_TRUE; 1869 PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept)); 1870 PetscCall(PetscFree3(es, hs, errors)); 1871 /* Slope should be about 2 */ 1872 if (tol >= 0) { 1873 PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope); 1874 } else if (isLinear || convRate) { 1875 if (isLinear) *isLinear = isLin; 1876 if (convRate) *convRate = slope; 1877 } else { 1878 if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope)); 1879 else PetscCall(PetscPrintf(comm, "Function appears to be linear\n")); 1880 } 1881 } 1882 PetscCall(MatDestroy(&J)); 1883 PetscFunctionReturn(0); 1884 } 1885 1886 PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u) { 1887 PetscFunctionBegin; 1888 PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL)); 1889 PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL)); 1890 PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL)); 1891 PetscFunctionReturn(0); 1892 } 1893 1894 /*@C 1895 DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 1896 1897 Input Parameters: 1898 + snes - the `SNES` object 1899 - u - representative `SNES` vector 1900 1901 Note: 1902 The user must call `PetscDSSetExactSolution()` beforehand 1903 1904 Level: developer 1905 @*/ 1906 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u) { 1907 DM dm; 1908 Vec sol; 1909 PetscBool check; 1910 1911 PetscFunctionBegin; 1912 PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check)); 1913 if (!check) PetscFunctionReturn(0); 1914 PetscCall(SNESGetDM(snes, &dm)); 1915 PetscCall(VecDuplicate(u, &sol)); 1916 PetscCall(SNESSetSolution(snes, sol)); 1917 PetscCall(DMSNESCheck_Internal(snes, dm, sol)); 1918 PetscCall(VecDestroy(&sol)); 1919 PetscFunctionReturn(0); 1920 } 1921