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