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