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