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 /************************** Interpolation *******************************/ 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "DMSNESConvertPlex" 11 static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) 12 { 13 PetscBool isPlex; 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr); 18 if (isPlex) { 19 *plex = dm; 20 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 21 } else { 22 ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr); 23 if (!*plex) { 24 ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr); 25 ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr); 26 if (copy) { 27 PetscInt i; 28 PetscObject obj; 29 const char *comps[3] = {"A","dmAux","dmCh"}; 30 31 ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr); 32 for (i = 0; i < 3; i++) { 33 ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr); 34 ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr); 35 } 36 } 37 } else { 38 ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr); 39 } 40 } 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "DMInterpolationCreate" 46 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 47 { 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 PetscValidPointer(ctx, 2); 52 ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr); 53 54 (*ctx)->comm = comm; 55 (*ctx)->dim = -1; 56 (*ctx)->nInput = 0; 57 (*ctx)->points = NULL; 58 (*ctx)->cells = NULL; 59 (*ctx)->n = -1; 60 (*ctx)->coords = NULL; 61 PetscFunctionReturn(0); 62 } 63 64 #undef __FUNCT__ 65 #define __FUNCT__ "DMInterpolationSetDim" 66 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 67 { 68 PetscFunctionBegin; 69 if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim); 70 ctx->dim = dim; 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "DMInterpolationGetDim" 76 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 77 { 78 PetscFunctionBegin; 79 PetscValidIntPointer(dim, 2); 80 *dim = ctx->dim; 81 PetscFunctionReturn(0); 82 } 83 84 #undef __FUNCT__ 85 #define __FUNCT__ "DMInterpolationSetDof" 86 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 87 { 88 PetscFunctionBegin; 89 if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof); 90 ctx->dof = dof; 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "DMInterpolationGetDof" 96 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 97 { 98 PetscFunctionBegin; 99 PetscValidIntPointer(dof, 2); 100 *dof = ctx->dof; 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "DMInterpolationAddPoints" 106 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 107 { 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 112 if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 113 ctx->nInput = n; 114 115 ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr); 116 ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 #undef __FUNCT__ 121 #define __FUNCT__ "DMInterpolationSetUp" 122 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) 123 { 124 MPI_Comm comm = ctx->comm; 125 PetscScalar *a; 126 PetscInt p, q, i; 127 PetscMPIInt rank, size; 128 PetscErrorCode ierr; 129 Vec pointVec; 130 PetscSF cellSF; 131 PetscLayout layout; 132 PetscReal *globalPoints; 133 PetscScalar *globalPointsScalar; 134 const PetscInt *ranges; 135 PetscMPIInt *counts, *displs; 136 const PetscSFNode *foundCells; 137 const PetscInt *foundPoints; 138 PetscMPIInt *foundProcs, *globalProcs; 139 PetscInt n, N, numFound; 140 141 PetscFunctionBegin; 142 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 143 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 144 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 145 if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 146 /* Locate points */ 147 n = ctx->nInput; 148 if (!redundantPoints) { 149 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 150 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 151 ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 152 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 153 ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 154 /* Communicate all points to all processes */ 155 ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr); 156 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 157 for (p = 0; p < size; ++p) { 158 counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 159 displs[p] = ranges[p]*ctx->dim; 160 } 161 ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 162 } else { 163 N = n; 164 globalPoints = ctx->points; 165 counts = displs = NULL; 166 layout = NULL; 167 } 168 #if 0 169 ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 170 /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 171 #else 172 #if defined(PETSC_USE_COMPLEX) 173 ierr = PetscMalloc1(N,&globalPointsScalar);CHKERRQ(ierr); 174 for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i]; 175 #else 176 globalPointsScalar = globalPoints; 177 #endif 178 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr); 179 ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 180 cellSF = NULL; 181 ierr = DMLocatePoints(dm, pointVec, &cellSF);CHKERRQ(ierr); 182 ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr); 183 #endif 184 for (p = 0; p < numFound; ++p) { 185 if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank; 186 else foundProcs[foundPoints ? foundPoints[p] : p] = size; 187 } 188 /* Let the lowest rank process own each point */ 189 ierr = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 190 ctx->n = 0; 191 for (p = 0; p < N; ++p) { 192 if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, globalPoints[p*ctx->dim+0], ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0, ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0); 193 else if (globalProcs[p] == rank) ctx->n++; 194 } 195 /* Create coordinates vector and array of owned cells */ 196 ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr); 197 ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 198 ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 199 ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 200 ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr); 201 ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 202 for (p = 0, q = 0, i = 0; p < N; ++p) { 203 if (globalProcs[p] == rank) { 204 PetscInt d; 205 206 for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d]; 207 ctx->cells[q++] = foundCells[p].index; 208 } 209 } 210 ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 211 #if 0 212 ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 213 #else 214 ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 215 ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr); 216 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 217 #endif 218 if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);} 219 if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);} 220 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 221 PetscFunctionReturn(0); 222 } 223 224 #undef __FUNCT__ 225 #define __FUNCT__ "DMInterpolationGetCoordinates" 226 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 227 { 228 PetscFunctionBegin; 229 PetscValidPointer(coordinates, 2); 230 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 231 *coordinates = ctx->coords; 232 PetscFunctionReturn(0); 233 } 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "DMInterpolationGetVector" 237 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 238 { 239 PetscErrorCode ierr; 240 241 PetscFunctionBegin; 242 PetscValidPointer(v, 2); 243 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 244 ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 245 ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 246 ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 247 ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr); 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "DMInterpolationRestoreVector" 253 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 254 { 255 PetscErrorCode ierr; 256 257 PetscFunctionBegin; 258 PetscValidPointer(v, 2); 259 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 260 ierr = VecDestroy(v);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "DMInterpolate_Triangle_Private" 266 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 267 { 268 PetscReal *v0, *J, *invJ, detJ; 269 const PetscScalar *coords; 270 PetscScalar *a; 271 PetscInt p; 272 PetscErrorCode ierr; 273 274 PetscFunctionBegin; 275 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 276 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 277 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 278 for (p = 0; p < ctx->n; ++p) { 279 PetscInt c = ctx->cells[p]; 280 PetscScalar *x = NULL; 281 PetscReal xi[4]; 282 PetscInt d, f, comp; 283 284 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 285 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 286 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 287 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 288 289 for (d = 0; d < ctx->dim; ++d) { 290 xi[d] = 0.0; 291 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 292 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]; 293 } 294 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 295 } 296 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 297 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 298 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "DMInterpolate_Tetrahedron_Private" 304 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 305 { 306 PetscReal *v0, *J, *invJ, detJ; 307 const PetscScalar *coords; 308 PetscScalar *a; 309 PetscInt p; 310 PetscErrorCode ierr; 311 312 PetscFunctionBegin; 313 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 314 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 315 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 316 for (p = 0; p < ctx->n; ++p) { 317 PetscInt c = ctx->cells[p]; 318 const PetscInt order[3] = {2, 1, 3}; 319 PetscScalar *x = NULL; 320 PetscReal xi[4]; 321 PetscInt d, f, comp; 322 323 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 324 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 325 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 326 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 327 328 for (d = 0; d < ctx->dim; ++d) { 329 xi[d] = 0.0; 330 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 331 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]; 332 } 333 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 334 } 335 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 336 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 337 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "QuadMap_Private" 343 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 344 { 345 const PetscScalar *vertices = (const PetscScalar*) ctx; 346 const PetscScalar x0 = vertices[0]; 347 const PetscScalar y0 = vertices[1]; 348 const PetscScalar x1 = vertices[2]; 349 const PetscScalar y1 = vertices[3]; 350 const PetscScalar x2 = vertices[4]; 351 const PetscScalar y2 = vertices[5]; 352 const PetscScalar x3 = vertices[6]; 353 const PetscScalar y3 = vertices[7]; 354 const PetscScalar f_1 = x1 - x0; 355 const PetscScalar g_1 = y1 - y0; 356 const PetscScalar f_3 = x3 - x0; 357 const PetscScalar g_3 = y3 - y0; 358 const PetscScalar f_01 = x2 - x1 - x3 + x0; 359 const PetscScalar g_01 = y2 - y1 - y3 + y0; 360 const PetscScalar *ref; 361 PetscScalar *real; 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 366 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 367 { 368 const PetscScalar p0 = ref[0]; 369 const PetscScalar p1 = ref[1]; 370 371 real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 372 real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 373 } 374 ierr = PetscLogFlops(28);CHKERRQ(ierr); 375 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 376 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 380 #include <petsc/private/dmimpl.h> 381 #undef __FUNCT__ 382 #define __FUNCT__ "QuadJacobian_Private" 383 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 384 { 385 const PetscScalar *vertices = (const PetscScalar*) ctx; 386 const PetscScalar x0 = vertices[0]; 387 const PetscScalar y0 = vertices[1]; 388 const PetscScalar x1 = vertices[2]; 389 const PetscScalar y1 = vertices[3]; 390 const PetscScalar x2 = vertices[4]; 391 const PetscScalar y2 = vertices[5]; 392 const PetscScalar x3 = vertices[6]; 393 const PetscScalar y3 = vertices[7]; 394 const PetscScalar f_01 = x2 - x1 - x3 + x0; 395 const PetscScalar g_01 = y2 - y1 - y3 + y0; 396 const PetscScalar *ref; 397 PetscErrorCode ierr; 398 399 PetscFunctionBegin; 400 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 401 { 402 const PetscScalar x = ref[0]; 403 const PetscScalar y = ref[1]; 404 const PetscInt rows[2] = {0, 1}; 405 PetscScalar values[4]; 406 407 values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5; 408 values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5; 409 ierr = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 410 } 411 ierr = PetscLogFlops(30);CHKERRQ(ierr); 412 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 413 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 414 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "DMInterpolate_Quad_Private" 420 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 421 { 422 DM dmCoord; 423 SNES snes; 424 KSP ksp; 425 PC pc; 426 Vec coordsLocal, r, ref, real; 427 Mat J; 428 const PetscScalar *coords; 429 PetscScalar *a; 430 PetscInt p; 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 435 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 436 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 437 ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 438 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 439 ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 440 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 441 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 442 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 443 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 444 ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 445 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 446 ierr = MatSetUp(J);CHKERRQ(ierr); 447 ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr); 448 ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr); 449 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 450 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 451 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 452 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 453 454 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 455 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 456 for (p = 0; p < ctx->n; ++p) { 457 PetscScalar *x = NULL, *vertices = NULL; 458 PetscScalar *xi; 459 PetscReal xir[2]; 460 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 461 462 /* Can make this do all points at once */ 463 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 464 if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2); 465 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 466 if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof); 467 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 468 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 469 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 470 xi[0] = coords[p*ctx->dim+0]; 471 xi[1] = coords[p*ctx->dim+1]; 472 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 473 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 474 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 475 xir[0] = PetscRealPart(xi[0]); 476 xir[1] = PetscRealPart(xi[1]); 477 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*ctx->dof+comp]*xir[0]*(1 - xir[1]) + x[2*ctx->dof+comp]*xir[0]*xir[1] + x[3*ctx->dof+comp]*(1 - xir[0])*xir[1]; 478 479 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 480 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 481 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 482 } 483 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 484 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 485 486 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 487 ierr = VecDestroy(&r);CHKERRQ(ierr); 488 ierr = VecDestroy(&ref);CHKERRQ(ierr); 489 ierr = VecDestroy(&real);CHKERRQ(ierr); 490 ierr = MatDestroy(&J);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "HexMap_Private" 496 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 497 { 498 const PetscScalar *vertices = (const PetscScalar*) ctx; 499 const PetscScalar x0 = vertices[0]; 500 const PetscScalar y0 = vertices[1]; 501 const PetscScalar z0 = vertices[2]; 502 const PetscScalar x1 = vertices[9]; 503 const PetscScalar y1 = vertices[10]; 504 const PetscScalar z1 = vertices[11]; 505 const PetscScalar x2 = vertices[6]; 506 const PetscScalar y2 = vertices[7]; 507 const PetscScalar z2 = vertices[8]; 508 const PetscScalar x3 = vertices[3]; 509 const PetscScalar y3 = vertices[4]; 510 const PetscScalar z3 = vertices[5]; 511 const PetscScalar x4 = vertices[12]; 512 const PetscScalar y4 = vertices[13]; 513 const PetscScalar z4 = vertices[14]; 514 const PetscScalar x5 = vertices[15]; 515 const PetscScalar y5 = vertices[16]; 516 const PetscScalar z5 = vertices[17]; 517 const PetscScalar x6 = vertices[18]; 518 const PetscScalar y6 = vertices[19]; 519 const PetscScalar z6 = vertices[20]; 520 const PetscScalar x7 = vertices[21]; 521 const PetscScalar y7 = vertices[22]; 522 const PetscScalar z7 = vertices[23]; 523 const PetscScalar f_1 = x1 - x0; 524 const PetscScalar g_1 = y1 - y0; 525 const PetscScalar h_1 = z1 - z0; 526 const PetscScalar f_3 = x3 - x0; 527 const PetscScalar g_3 = y3 - y0; 528 const PetscScalar h_3 = z3 - z0; 529 const PetscScalar f_4 = x4 - x0; 530 const PetscScalar g_4 = y4 - y0; 531 const PetscScalar h_4 = z4 - z0; 532 const PetscScalar f_01 = x2 - x1 - x3 + x0; 533 const PetscScalar g_01 = y2 - y1 - y3 + y0; 534 const PetscScalar h_01 = z2 - z1 - z3 + z0; 535 const PetscScalar f_12 = x7 - x3 - x4 + x0; 536 const PetscScalar g_12 = y7 - y3 - y4 + y0; 537 const PetscScalar h_12 = z7 - z3 - z4 + z0; 538 const PetscScalar f_02 = x5 - x1 - x4 + x0; 539 const PetscScalar g_02 = y5 - y1 - y4 + y0; 540 const PetscScalar h_02 = z5 - z1 - z4 + z0; 541 const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 542 const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 543 const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 544 const PetscScalar *ref; 545 PetscScalar *real; 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 550 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 551 { 552 const PetscScalar p0 = ref[0]; 553 const PetscScalar p1 = ref[1]; 554 const PetscScalar p2 = ref[2]; 555 556 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; 557 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; 558 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; 559 } 560 ierr = PetscLogFlops(114);CHKERRQ(ierr); 561 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 562 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 #undef __FUNCT__ 567 #define __FUNCT__ "HexJacobian_Private" 568 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 569 { 570 const PetscScalar *vertices = (const PetscScalar*) ctx; 571 const PetscScalar x0 = vertices[0]; 572 const PetscScalar y0 = vertices[1]; 573 const PetscScalar z0 = vertices[2]; 574 const PetscScalar x1 = vertices[9]; 575 const PetscScalar y1 = vertices[10]; 576 const PetscScalar z1 = vertices[11]; 577 const PetscScalar x2 = vertices[6]; 578 const PetscScalar y2 = vertices[7]; 579 const PetscScalar z2 = vertices[8]; 580 const PetscScalar x3 = vertices[3]; 581 const PetscScalar y3 = vertices[4]; 582 const PetscScalar z3 = vertices[5]; 583 const PetscScalar x4 = vertices[12]; 584 const PetscScalar y4 = vertices[13]; 585 const PetscScalar z4 = vertices[14]; 586 const PetscScalar x5 = vertices[15]; 587 const PetscScalar y5 = vertices[16]; 588 const PetscScalar z5 = vertices[17]; 589 const PetscScalar x6 = vertices[18]; 590 const PetscScalar y6 = vertices[19]; 591 const PetscScalar z6 = vertices[20]; 592 const PetscScalar x7 = vertices[21]; 593 const PetscScalar y7 = vertices[22]; 594 const PetscScalar z7 = vertices[23]; 595 const PetscScalar f_xy = x2 - x1 - x3 + x0; 596 const PetscScalar g_xy = y2 - y1 - y3 + y0; 597 const PetscScalar h_xy = z2 - z1 - z3 + z0; 598 const PetscScalar f_yz = x7 - x3 - x4 + x0; 599 const PetscScalar g_yz = y7 - y3 - y4 + y0; 600 const PetscScalar h_yz = z7 - z3 - z4 + z0; 601 const PetscScalar f_xz = x5 - x1 - x4 + x0; 602 const PetscScalar g_xz = y5 - y1 - y4 + y0; 603 const PetscScalar h_xz = z5 - z1 - z4 + z0; 604 const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 605 const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 606 const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 607 const PetscScalar *ref; 608 PetscErrorCode ierr; 609 610 PetscFunctionBegin; 611 ierr = VecGetArrayRead(Xref, &ref);CHKERRQ(ierr); 612 { 613 const PetscScalar x = ref[0]; 614 const PetscScalar y = ref[1]; 615 const PetscScalar z = ref[2]; 616 const PetscInt rows[3] = {0, 1, 2}; 617 PetscScalar values[9]; 618 619 values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0; 620 values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0; 621 values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0; 622 values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0; 623 values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0; 624 values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0; 625 values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0; 626 values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0; 627 values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0; 628 629 ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 630 } 631 ierr = PetscLogFlops(152);CHKERRQ(ierr); 632 ierr = VecRestoreArrayRead(Xref, &ref);CHKERRQ(ierr); 633 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 634 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 635 PetscFunctionReturn(0); 636 } 637 638 #undef __FUNCT__ 639 #define __FUNCT__ "DMInterpolate_Hex_Private" 640 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 641 { 642 DM dmCoord; 643 SNES snes; 644 KSP ksp; 645 PC pc; 646 Vec coordsLocal, r, ref, real; 647 Mat J; 648 const PetscScalar *coords; 649 PetscScalar *a; 650 PetscInt p; 651 PetscErrorCode ierr; 652 653 PetscFunctionBegin; 654 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 655 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 656 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 657 ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 658 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 659 ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 660 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 661 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 662 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 663 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 664 ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 665 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 666 ierr = MatSetUp(J);CHKERRQ(ierr); 667 ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr); 668 ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr); 669 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 670 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 671 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 672 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 673 674 ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 675 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 676 for (p = 0; p < ctx->n; ++p) { 677 PetscScalar *x = NULL, *vertices = NULL; 678 PetscScalar *xi; 679 PetscReal xir[3]; 680 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 681 682 /* Can make this do all points at once */ 683 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 684 if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3); 685 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 686 if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof); 687 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 688 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 689 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 690 xi[0] = coords[p*ctx->dim+0]; 691 xi[1] = coords[p*ctx->dim+1]; 692 xi[2] = coords[p*ctx->dim+2]; 693 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 694 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 695 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 696 xir[0] = PetscRealPart(xi[0]); 697 xir[1] = PetscRealPart(xi[1]); 698 xir[2] = PetscRealPart(xi[2]); 699 for (comp = 0; comp < ctx->dof; ++comp) { 700 a[p*ctx->dof+comp] = 701 x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) + 702 x[3*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) + 703 x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) + 704 x[1*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) + 705 x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] + 706 x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] + 707 x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] + 708 x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2]; 709 } 710 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 711 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 712 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 713 } 714 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 715 ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr); 716 717 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 718 ierr = VecDestroy(&r);CHKERRQ(ierr); 719 ierr = VecDestroy(&ref);CHKERRQ(ierr); 720 ierr = VecDestroy(&real);CHKERRQ(ierr); 721 ierr = MatDestroy(&J);CHKERRQ(ierr); 722 PetscFunctionReturn(0); 723 } 724 725 #undef __FUNCT__ 726 #define __FUNCT__ "DMInterpolationEvaluate" 727 /* 728 Input Parameters: 729 + ctx - The DMInterpolationInfo context 730 . dm - The DM 731 - x - The local vector containing the field to be interpolated 732 733 Output Parameters: 734 . v - The vector containing the interpolated values 735 */ 736 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 737 { 738 PetscInt dim, coneSize, n; 739 PetscErrorCode ierr; 740 741 PetscFunctionBegin; 742 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 743 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 744 PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 745 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 746 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); 747 if (n) { 748 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 749 ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 750 if (dim == 2) { 751 if (coneSize == 3) { 752 ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr); 753 } else if (coneSize == 4) { 754 ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 755 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 756 } else if (dim == 3) { 757 if (coneSize == 4) { 758 ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr); 759 } else { 760 ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 761 } 762 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 763 } 764 PetscFunctionReturn(0); 765 } 766 767 #undef __FUNCT__ 768 #define __FUNCT__ "DMInterpolationDestroy" 769 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 770 { 771 PetscErrorCode ierr; 772 773 PetscFunctionBegin; 774 PetscValidPointer(ctx, 2); 775 ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 776 ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 777 ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 778 ierr = PetscFree(*ctx);CHKERRQ(ierr); 779 *ctx = NULL; 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "SNESMonitorFields" 785 /*@C 786 SNESMonitorFields - Monitors the residual for each field separately 787 788 Collective on SNES 789 790 Input Parameters: 791 + snes - the SNES context 792 . its - iteration number 793 . fgnorm - 2-norm of residual 794 - vf - PetscViewerAndFormat of type ASCII 795 796 Notes: 797 This routine prints the residual norm at each iteration. 798 799 Level: intermediate 800 801 .keywords: SNES, nonlinear, default, monitor, norm 802 .seealso: SNESMonitorSet(), SNESMonitorDefault() 803 @*/ 804 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 805 { 806 PetscViewer viewer = vf->viewer; 807 Vec res; 808 DM dm; 809 PetscSection s; 810 const PetscScalar *r; 811 PetscReal *lnorms, *norms; 812 PetscInt numFields, f, pStart, pEnd, p; 813 PetscErrorCode ierr; 814 815 PetscFunctionBegin; 816 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 817 ierr = SNESGetFunction(snes, &res, 0, 0);CHKERRQ(ierr); 818 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 819 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 820 ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr); 821 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 822 ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr); 823 ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr); 824 for (p = pStart; p < pEnd; ++p) { 825 for (f = 0; f < numFields; ++f) { 826 PetscInt fdof, foff, d; 827 828 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 829 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 830 for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d])); 831 } 832 } 833 ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr); 834 ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 835 ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 836 ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 837 ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr); 838 for (f = 0; f < numFields; ++f) { 839 if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 840 ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr); 841 } 842 ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr); 843 ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 844 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 845 ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr); 846 PetscFunctionReturn(0); 847 } 848 849 /********************* Residual Computation **************************/ 850 851 #undef __FUNCT__ 852 #define __FUNCT__ "DMPlexSNESGetGeometryFEM" 853 /*@ 854 DMPlexSNESGetGeometryFEM - Return precomputed geometric data 855 856 Input Parameter: 857 . dm - The DM 858 859 Output Parameters: 860 . cellgeom - The values precomputed from cell geometry 861 862 Level: developer 863 864 .seealso: DMPlexSNESSetFunctionLocal() 865 @*/ 866 PetscErrorCode DMPlexSNESGetGeometryFEM(DM dm, Vec *cellgeom) 867 { 868 DMSNES dmsnes; 869 PetscObject obj; 870 PetscErrorCode ierr; 871 872 PetscFunctionBegin; 873 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 874 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 875 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", &obj);CHKERRQ(ierr); 876 if (!obj) { 877 Vec cellgeom; 878 879 ierr = DMPlexComputeGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 880 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject) cellgeom);CHKERRQ(ierr); 881 ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 882 } 883 if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject *) cellgeom);CHKERRQ(ierr);} 884 PetscFunctionReturn(0); 885 } 886 887 #undef __FUNCT__ 888 #define __FUNCT__ "DMPlexSNESGetGeometryFVM" 889 /*@ 890 DMPlexSNESGetGeometryFVM - Return precomputed geometric data 891 892 Input Parameter: 893 . dm - The DM 894 895 Output Parameters: 896 + facegeom - The values precomputed from face geometry 897 . cellgeom - The values precomputed from cell geometry 898 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell 899 900 Level: developer 901 902 .seealso: DMPlexTSSetRHSFunctionLocal() 903 @*/ 904 PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 905 { 906 DMSNES dmsnes; 907 PetscObject obj; 908 PetscErrorCode ierr; 909 910 PetscFunctionBegin; 911 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 912 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 913 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", &obj);CHKERRQ(ierr); 914 if (!obj) { 915 Vec cellgeom, facegeom; 916 917 ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr); 918 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr); 919 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr); 920 ierr = VecDestroy(&facegeom);CHKERRQ(ierr); 921 ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 922 } 923 if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);} 924 if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);} 925 if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);} 926 PetscFunctionReturn(0); 927 } 928 929 #undef __FUNCT__ 930 #define __FUNCT__ "DMPlexSNESGetGradientDM" 931 /*@ 932 DMPlexSNESGetGradientDM - Return gradient data layout 933 934 Input Parameters: 935 + dm - The DM 936 - fv - The PetscFV 937 938 Output Parameter: 939 . dmGrad - The layout for gradient values 940 941 Level: developer 942 943 .seealso: DMPlexSNESGetGeometryFVM() 944 @*/ 945 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad) 946 { 947 DMSNES dmsnes; 948 PetscObject obj; 949 PetscBool computeGradients; 950 PetscErrorCode ierr; 951 952 PetscFunctionBegin; 953 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 954 PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2); 955 PetscValidPointer(dmGrad,3); 956 ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr); 957 if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);} 958 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 959 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", &obj);CHKERRQ(ierr); 960 if (!obj) { 961 DM dmGrad; 962 Vec faceGeometry, cellGeometry; 963 964 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 965 ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr); 966 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr); 967 ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 968 } 969 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr); 970 PetscFunctionReturn(0); 971 } 972 973 #undef __FUNCT__ 974 #define __FUNCT__ "DMPlexGetCellFields" 975 /*@C 976 DMPlexGetCellFields - Retrieve the field values values for a chunk of cells 977 978 Input Parameters: 979 + dm - The DM 980 . cStart - The first cell to include 981 . cEnd - The first cell to exclude 982 . locX - A local vector with the solution fields 983 . locX_t - A local vector with solution field time derivatives, or NULL 984 - locA - A local vector with auxiliary fields, or NULL 985 986 Output Parameters: 987 + u - The field coefficients 988 . u_t - The fields derivative coefficients 989 - a - The auxiliary field coefficients 990 991 Level: developer 992 993 .seealso: DMPlexGetFaceFields() 994 @*/ 995 PetscErrorCode DMPlexGetCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a) 996 { 997 DM dmAux; 998 PetscSection section, sectionAux; 999 PetscDS prob; 1000 PetscInt numCells = cEnd - cStart, totDim, totDimAux, c; 1001 PetscErrorCode ierr; 1002 1003 PetscFunctionBegin; 1004 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1005 PetscValidHeaderSpecific(locX, VEC_CLASSID, 4); 1006 if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);} 1007 if (locA) {PetscValidHeaderSpecific(locA, VEC_CLASSID, 6);} 1008 PetscValidPointer(u, 7); 1009 PetscValidPointer(u_t, 8); 1010 PetscValidPointer(a, 9); 1011 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1012 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1013 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1014 if (locA) { 1015 PetscDS probAux; 1016 1017 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1018 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1019 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1020 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1021 } 1022 ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u);CHKERRQ(ierr); 1023 if (locX_t) {ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u_t);CHKERRQ(ierr);} else {*u_t = NULL;} 1024 if (locA) {ierr = DMGetWorkArray(dm, numCells*totDimAux, PETSC_SCALAR, a);CHKERRQ(ierr);} else {*a = NULL;} 1025 for (c = cStart; c < cEnd; ++c) { 1026 PetscScalar *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a; 1027 PetscInt i; 1028 1029 ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 1030 for (i = 0; i < totDim; ++i) ul[(c-cStart)*totDim+i] = x[i]; 1031 ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 1032 if (locX_t) { 1033 ierr = DMPlexVecGetClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr); 1034 for (i = 0; i < totDim; ++i) ul_t[(c-cStart)*totDim+i] = x_t[i]; 1035 ierr = DMPlexVecRestoreClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr); 1036 } 1037 if (locA) { 1038 DM dmAuxPlex; 1039 1040 ierr = DMSNESConvertPlex(dmAux, &dmAuxPlex, PETSC_FALSE);CHKERRQ(ierr); 1041 ierr = DMPlexVecGetClosure(dmAuxPlex, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 1042 for (i = 0; i < totDimAux; ++i) al[(c-cStart)*totDimAux+i] = x[i]; 1043 ierr = DMPlexVecRestoreClosure(dmAuxPlex, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 1044 ierr = DMDestroy(&dmAuxPlex);CHKERRQ(ierr); 1045 } 1046 } 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "DMPlexRestoreCellFields" 1052 /*@C 1053 DMPlexRestoreCellFields - Restore the field values values for a chunk of cells 1054 1055 Input Parameters: 1056 + dm - The DM 1057 . cStart - The first cell to include 1058 . cEnd - The first cell to exclude 1059 . locX - A local vector with the solution fields 1060 . locX_t - A local vector with solution field time derivatives, or NULL 1061 - locA - A local vector with auxiliary fields, or NULL 1062 1063 Output Parameters: 1064 + u - The field coefficients 1065 . u_t - The fields derivative coefficients 1066 - a - The auxiliary field coefficients 1067 1068 Level: developer 1069 1070 .seealso: DMPlexGetFaceFields() 1071 @*/ 1072 PetscErrorCode DMPlexRestoreCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a) 1073 { 1074 PetscErrorCode ierr; 1075 1076 PetscFunctionBegin; 1077 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u);CHKERRQ(ierr); 1078 if (*u_t) {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u_t);CHKERRQ(ierr);} 1079 if (*a) {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, a);CHKERRQ(ierr);} 1080 PetscFunctionReturn(0); 1081 } 1082 1083 #undef __FUNCT__ 1084 #define __FUNCT__ "DMPlexGetFaceFields" 1085 /*@C 1086 DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces 1087 1088 Input Parameters: 1089 + dm - The DM 1090 . fStart - The first face to include 1091 . fEnd - The first face to exclude 1092 . locX - A local vector with the solution fields 1093 . locX_t - A local vector with solution field time derivatives, or NULL 1094 . faceGeometry - A local vector with face geometry 1095 . cellGeometry - A local vector with cell geometry 1096 - locaGrad - A local vector with field gradients, or NULL 1097 1098 Output Parameters: 1099 + uL - The field values at the left side of the face 1100 - uR - The field values at the right side of the face 1101 1102 Level: developer 1103 1104 .seealso: DMPlexGetCellFields() 1105 @*/ 1106 PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR) 1107 { 1108 DM dmFace, dmCell, dmGrad = NULL; 1109 PetscSection section; 1110 PetscDS prob; 1111 DMLabel ghostLabel; 1112 const PetscScalar *facegeom, *cellgeom, *x, *lgrad; 1113 PetscBool *isFE; 1114 PetscInt dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face; 1115 PetscErrorCode ierr; 1116 1117 PetscFunctionBegin; 1118 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1119 PetscValidHeaderSpecific(locX, VEC_CLASSID, 4); 1120 if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);} 1121 PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 6); 1122 PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 7); 1123 if (locGrad) {PetscValidHeaderSpecific(locGrad, VEC_CLASSID, 8);} 1124 PetscValidPointer(uL, 9); 1125 PetscValidPointer(uR, 10); 1126 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1127 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1128 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1129 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1130 ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 1131 ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr); 1132 for (f = 0; f < Nf; ++f) { 1133 PetscObject obj; 1134 PetscClassId id; 1135 1136 ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 1137 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1138 if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 1139 else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 1140 else {isFE[f] = PETSC_FALSE;} 1141 } 1142 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1143 ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 1144 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1145 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1146 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1147 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1148 if (locGrad) { 1149 ierr = VecGetDM(locGrad, &dmGrad);CHKERRQ(ierr); 1150 ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1151 } 1152 ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uL);CHKERRQ(ierr); 1153 ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uR);CHKERRQ(ierr); 1154 /* Right now just eat the extra work for FE (could make a cell loop) */ 1155 for (face = fStart, iface = 0; face < fEnd; ++face) { 1156 const PetscInt *cells; 1157 PetscFVFaceGeom *fg; 1158 PetscFVCellGeom *cgL, *cgR; 1159 PetscScalar *xL, *xR, *gL, *gR; 1160 PetscScalar *uLl = *uL, *uRl = *uR; 1161 PetscInt ghost, nsupp; 1162 1163 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1164 ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr); 1165 if (ghost >= 0) continue; 1166 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1167 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1168 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 1169 ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 1170 for (f = 0; f < Nf; ++f) { 1171 PetscInt off; 1172 1173 ierr = PetscDSGetComponentOffset(prob, f, &off);CHKERRQ(ierr); 1174 if (isFE[f]) { 1175 const PetscInt *cone; 1176 PetscInt comp, coneSize, faceLocL, faceLocR, ldof, rdof, d; 1177 1178 xL = xR = NULL; 1179 ierr = DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);CHKERRQ(ierr); 1180 ierr = DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);CHKERRQ(ierr); 1181 ierr = DMPlexGetCone(dm, cells[0], &cone);CHKERRQ(ierr); 1182 ierr = DMPlexGetConeSize(dm, cells[0], &coneSize);CHKERRQ(ierr); 1183 for (faceLocL = 0; faceLocL < coneSize; ++faceLocL) if (cone[faceLocL] == face) break; 1184 if (faceLocL == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of cell %d", face, cells[0]); 1185 ierr = DMPlexGetCone(dm, cells[1], &cone);CHKERRQ(ierr); 1186 ierr = DMPlexGetConeSize(dm, cells[1], &coneSize);CHKERRQ(ierr); 1187 for (faceLocR = 0; faceLocR < coneSize; ++faceLocR) if (cone[faceLocR] == face) break; 1188 if (faceLocR == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of cell %d", face, cells[1]); 1189 /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */ 1190 ierr = EvaluateFaceFields(prob, f, faceLocL, xL, &uLl[iface*Nc+off]);CHKERRQ(ierr); 1191 if (rdof == ldof) {ierr = EvaluateFaceFields(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);CHKERRQ(ierr);} 1192 else {ierr = PetscSectionGetFieldComponents(section, f, &comp);CHKERRQ(ierr); for(d = 0; d < comp; ++d) uRl[iface*Nc+off+d] = uLl[iface*Nc+off+d];} 1193 ierr = DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);CHKERRQ(ierr); 1194 ierr = DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);CHKERRQ(ierr); 1195 } else { 1196 PetscFV fv; 1197 PetscInt numComp, c; 1198 1199 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fv);CHKERRQ(ierr); 1200 ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr); 1201 if (nsupp > 2) { 1202 for (f = 0; f < Nf; ++f) { 1203 PetscInt off; 1204 1205 ierr = PetscDSGetComponentOffset(prob, f, &off);CHKERRQ(ierr); 1206 ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr); 1207 for (c = 0; c < numComp; ++c) { 1208 uLl[iface*Nc+off+c] = 0.; 1209 uRl[iface*Nc+off+c] = 0.; 1210 } 1211 } 1212 continue; 1213 } 1214 ierr = DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL);CHKERRQ(ierr); 1215 ierr = DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR);CHKERRQ(ierr); 1216 if (dmGrad) { 1217 PetscReal dxL[3], dxR[3]; 1218 1219 ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); 1220 ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); 1221 DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL); 1222 DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR); 1223 for (c = 0; c < numComp; ++c) { 1224 uLl[iface*Nc+off+c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c*dim], dxL); 1225 uRl[iface*Nc+off+c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c*dim], dxR); 1226 } 1227 } else { 1228 for (c = 0; c < numComp; ++c) { 1229 uLl[iface*Nc+off+c] = xL[c]; 1230 uRl[iface*Nc+off+c] = xR[c]; 1231 } 1232 } 1233 } 1234 } 1235 ++iface; 1236 } 1237 ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 1238 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1239 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1240 if (locGrad) { 1241 ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1242 } 1243 ierr = PetscFree(isFE);CHKERRQ(ierr); 1244 PetscFunctionReturn(0); 1245 } 1246 1247 #undef __FUNCT__ 1248 #define __FUNCT__ "DMPlexRestoreFaceFields" 1249 /*@C 1250 DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces 1251 1252 Input Parameters: 1253 + dm - The DM 1254 . fStart - The first face to include 1255 . fEnd - The first face to exclude 1256 . locX - A local vector with the solution fields 1257 . locX_t - A local vector with solution field time derivatives, or NULL 1258 . faceGeometry - A local vector with face geometry 1259 . cellGeometry - A local vector with cell geometry 1260 - locaGrad - A local vector with field gradients, or NULL 1261 1262 Output Parameters: 1263 + uL - The field values at the left side of the face 1264 - uR - The field values at the right side of the face 1265 1266 Level: developer 1267 1268 .seealso: DMPlexGetFaceFields() 1269 @*/ 1270 PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR) 1271 { 1272 PetscErrorCode ierr; 1273 1274 PetscFunctionBegin; 1275 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uL);CHKERRQ(ierr); 1276 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uR);CHKERRQ(ierr); 1277 PetscFunctionReturn(0); 1278 } 1279 1280 #undef __FUNCT__ 1281 #define __FUNCT__ "DMPlexGetFaceGeometry" 1282 /*@C 1283 DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces 1284 1285 Input Parameters: 1286 + dm - The DM 1287 . fStart - The first face to include 1288 . fEnd - The first face to exclude 1289 . faceGeometry - A local vector with face geometry 1290 - cellGeometry - A local vector with cell geometry 1291 1292 Output Parameters: 1293 + fgeom - The extract the face centroid and normal 1294 - vol - The cell volume 1295 1296 Level: developer 1297 1298 .seealso: DMPlexGetCellFields() 1299 @*/ 1300 PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol) 1301 { 1302 DM dmFace, dmCell; 1303 DMLabel ghostLabel; 1304 const PetscScalar *facegeom, *cellgeom; 1305 PetscInt dim, numFaces = fEnd - fStart, iface, face; 1306 PetscErrorCode ierr; 1307 1308 PetscFunctionBegin; 1309 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1310 PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 4); 1311 PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 5); 1312 PetscValidPointer(fgeom, 6); 1313 PetscValidPointer(vol, 7); 1314 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1315 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1316 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1317 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1318 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1319 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1320 ierr = PetscMalloc1(numFaces, fgeom);CHKERRQ(ierr); 1321 ierr = DMGetWorkArray(dm, numFaces*2, PETSC_SCALAR, vol);CHKERRQ(ierr); 1322 for (face = fStart, iface = 0; face < fEnd; ++face) { 1323 const PetscInt *cells; 1324 PetscFVFaceGeom *fg; 1325 PetscFVCellGeom *cgL, *cgR; 1326 PetscFVFaceGeom *fgeoml = *fgeom; 1327 PetscReal *voll = *vol; 1328 PetscInt ghost, d; 1329 1330 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1331 if (ghost >= 0) continue; 1332 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1333 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1334 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 1335 ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 1336 for (d = 0; d < dim; ++d) { 1337 fgeoml[iface].centroid[d] = fg->centroid[d]; 1338 fgeoml[iface].normal[d] = fg->normal[d]; 1339 } 1340 voll[iface*2+0] = cgL->volume; 1341 voll[iface*2+1] = cgR->volume; 1342 ++iface; 1343 } 1344 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1345 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1346 PetscFunctionReturn(0); 1347 } 1348 1349 #undef __FUNCT__ 1350 #define __FUNCT__ "DMPlexRestoreFaceGeometry" 1351 /*@C 1352 DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces 1353 1354 Input Parameters: 1355 + dm - The DM 1356 . fStart - The first face to include 1357 . fEnd - The first face to exclude 1358 . faceGeometry - A local vector with face geometry 1359 - cellGeometry - A local vector with cell geometry 1360 1361 Output Parameters: 1362 + fgeom - The extract the face centroid and normal 1363 - vol - The cell volume 1364 1365 Level: developer 1366 1367 .seealso: DMPlexGetFaceFields() 1368 @*/ 1369 PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol) 1370 { 1371 PetscErrorCode ierr; 1372 1373 PetscFunctionBegin; 1374 ierr = PetscFree(*fgeom);CHKERRQ(ierr); 1375 ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, vol);CHKERRQ(ierr); 1376 PetscFunctionReturn(0); 1377 } 1378 1379 #undef __FUNCT__ 1380 #define __FUNCT__ "DMPlexApplyLimiter_Internal" 1381 static PetscErrorCode DMPlexApplyLimiter_Internal (DM dm, DM dmCell, PetscLimiter lim, PetscInt dim, PetscInt totDim, PetscInt cell, PetscInt face, PetscInt fStart, PetscInt fEnd, PetscReal *cellPhi, const PetscScalar *x, 1382 const PetscScalar *cellgeom, const PetscFVCellGeom *cg, const PetscScalar *cx, const PetscScalar *cgrad) 1383 { 1384 const PetscInt *children; 1385 PetscInt numChildren; 1386 PetscErrorCode ierr; 1387 1388 PetscFunctionBegin; 1389 ierr = DMPlexGetTreeChildren(dm,face,&numChildren,&children);CHKERRQ(ierr); 1390 if (numChildren) { 1391 PetscInt c; 1392 1393 for (c = 0; c < numChildren; c++) { 1394 PetscInt childFace = children[c]; 1395 1396 if (childFace >= fStart && childFace < fEnd) { 1397 ierr = DMPlexApplyLimiter_Internal(dm,dmCell,lim,dim,totDim,cell,childFace,fStart,fEnd,cellPhi,x,cellgeom,cg,cx,cgrad);CHKERRQ(ierr); 1398 } 1399 } 1400 } 1401 else { 1402 PetscScalar *ncx; 1403 PetscFVCellGeom *ncg; 1404 const PetscInt *fcells; 1405 PetscInt ncell, d; 1406 PetscReal v[3]; 1407 1408 ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); 1409 ncell = cell == fcells[0] ? fcells[1] : fcells[0]; 1410 ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); 1411 ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); 1412 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v); 1413 for (d = 0; d < totDim; ++d) { 1414 /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ 1415 PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DMPlex_DotD_Internal(dim, &cgrad[d*dim], v); 1416 1417 ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); 1418 cellPhi[d] = PetscMin(cellPhi[d], phi); 1419 } 1420 } 1421 PetscFunctionReturn(0); 1422 } 1423 1424 #undef __FUNCT__ 1425 #define __FUNCT__ "DMPlexReconstructGradients_Internal" 1426 PetscErrorCode DMPlexReconstructGradients_Internal(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, Vec locX, Vec grad) 1427 { 1428 DM dmFace, dmCell, dmGrad; 1429 DMLabel ghostLabel; 1430 PetscDS prob; 1431 PetscFV fvm; 1432 PetscLimiter lim; 1433 const PetscScalar *facegeom, *cellgeom, *x; 1434 PetscScalar *gr; 1435 PetscReal *cellPhi; 1436 PetscInt dim, face, cell, totDim, cStart, cEnd, cEndInterior; 1437 PetscErrorCode ierr; 1438 1439 PetscFunctionBegin; 1440 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1441 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1442 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1443 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1444 ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 1445 ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); 1446 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1447 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1448 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1449 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1450 ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 1451 ierr = VecGetDM(grad, &dmGrad);CHKERRQ(ierr); 1452 ierr = VecZeroEntries(grad);CHKERRQ(ierr); 1453 ierr = VecGetArray(grad, &gr);CHKERRQ(ierr); 1454 /* Reconstruct gradients */ 1455 for (face = fStart; face < fEnd; ++face) { 1456 const PetscInt *cells; 1457 PetscFVFaceGeom *fg; 1458 PetscScalar *cx[2]; 1459 PetscScalar *cgrad[2]; 1460 PetscBool boundary; 1461 PetscInt ghost, c, pd, d, numChildren, numCells; 1462 1463 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1464 ierr = DMIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 1465 ierr = DMPlexGetTreeChildren(dm, face, &numChildren, NULL);CHKERRQ(ierr); 1466 if (ghost >= 0 || boundary || numChildren) continue; 1467 ierr = DMPlexGetSupportSize(dm, face, &numCells);CHKERRQ(ierr); 1468 if (numCells != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "facet %d has %d support points: expected 2",face,numCells); 1469 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1470 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1471 for (c = 0; c < 2; ++c) { 1472 ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); 1473 ierr = DMPlexPointGlobalRef(dmGrad, cells[c], gr, &cgrad[c]);CHKERRQ(ierr); 1474 } 1475 for (pd = 0; pd < totDim; ++pd) { 1476 PetscScalar delta = cx[1][pd] - cx[0][pd]; 1477 1478 for (d = 0; d < dim; ++d) { 1479 if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; 1480 if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; 1481 } 1482 } 1483 } 1484 /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ 1485 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1486 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1487 cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior; 1488 ierr = DMGetWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 1489 for (cell = dmGrad && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { 1490 const PetscInt *faces; 1491 PetscScalar *cx; 1492 PetscFVCellGeom *cg; 1493 PetscScalar *cgrad; 1494 PetscInt coneSize, f, pd, d; 1495 1496 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1497 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1498 ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); 1499 ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); 1500 ierr = DMPlexPointGlobalRef(dmGrad, cell, gr, &cgrad);CHKERRQ(ierr); 1501 if (!cgrad) continue; /* Unowned overlap cell, we do not compute */ 1502 /* Limiter will be minimum value over all neighbors */ 1503 for (d = 0; d < totDim; ++d) cellPhi[d] = PETSC_MAX_REAL; 1504 for (f = 0; f < coneSize; ++f) { 1505 ierr = DMPlexApplyLimiter_Internal(dm,dmCell,lim,dim,totDim,cell,faces[f],fStart,fEnd,cellPhi,x,cellgeom,cg,cx,cgrad);CHKERRQ(ierr); 1506 } 1507 /* Apply limiter to gradient */ 1508 for (pd = 0; pd < totDim; ++pd) 1509 /* Scalar limiter applied to each component separately */ 1510 for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; 1511 } 1512 ierr = DMRestoreWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 1513 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1514 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1515 ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 1516 ierr = VecRestoreArray(grad, &gr);CHKERRQ(ierr); 1517 PetscFunctionReturn(0); 1518 } 1519 1520 #undef __FUNCT__ 1521 #define __FUNCT__ "DMPlexComputeBdResidual_Internal" 1522 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, Vec locF, void *user) 1523 { 1524 DM_Plex *mesh = (DM_Plex *) dm->data; 1525 PetscSection section; 1526 PetscDS prob; 1527 DMLabel depth; 1528 PetscFECellGeom *cgeom; 1529 PetscScalar *u = NULL, *u_t = NULL, *elemVec = NULL; 1530 PetscInt dim, Nf, f, totDimBd, numBd, bd; 1531 PetscErrorCode ierr; 1532 1533 PetscFunctionBegin; 1534 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1535 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1536 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1537 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1538 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 1539 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1540 ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1541 for (bd = 0; bd < numBd; ++bd) { 1542 const char *bdLabel; 1543 DMLabel label; 1544 IS pointIS; 1545 const PetscInt *points; 1546 const PetscInt *values; 1547 PetscInt field, numValues, v, numPoints, p, dep, numFaces; 1548 PetscBool isEssential; 1549 PetscObject obj; 1550 PetscClassId id; 1551 1552 ierr = DMGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1553 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 1554 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1555 if ((id != PETSCFE_CLASSID) || isEssential) continue; 1556 ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1557 for (v = 0; v < numValues; ++v) { 1558 ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr); 1559 ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr); 1560 if (!pointIS) continue; /* No points with that id on this process */ 1561 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1562 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1563 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1564 if (dep == dim-1) ++numFaces; 1565 } 1566 ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd,&elemVec);CHKERRQ(ierr); 1567 if (locX_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 1568 for (p = 0, f = 0; p < numPoints; ++p) { 1569 const PetscInt point = points[p]; 1570 PetscScalar *x = NULL; 1571 PetscInt i; 1572 1573 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1574 if (dep != dim-1) continue; 1575 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr); 1576 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr); 1577 if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point); 1578 /* TODO: Matt, this is wrong if feBd does not match fe: i.e., if the order differs. */ 1579 ierr = DMPlexVecGetClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr); 1580 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 1581 ierr = DMPlexVecRestoreClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr); 1582 if (locX_t) { 1583 ierr = DMPlexVecGetClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr); 1584 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 1585 ierr = DMPlexVecRestoreClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr); 1586 } 1587 ++f; 1588 } 1589 for (f = 0; f < Nf; ++f) { 1590 PetscFE fe; 1591 PetscQuadrature q; 1592 PetscInt numQuadPoints, Nb; 1593 /* Conforming batches */ 1594 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1595 /* Remainder */ 1596 PetscInt Nr, offset; 1597 1598 ierr = PetscDSGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1599 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1600 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1601 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1602 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1603 blockSize = Nb*numQuadPoints; 1604 batchSize = numBlocks * blockSize; 1605 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1606 numChunks = numFaces / (numBatches*batchSize); 1607 Ne = numChunks*numBatches*batchSize; 1608 Nr = numFaces % (numBatches*batchSize); 1609 offset = numFaces - Nr; 1610 ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, cgeom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr); 1611 ierr = PetscFEIntegrateBdResidual(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemVec[offset*totDimBd]);CHKERRQ(ierr); 1612 } 1613 for (p = 0, f = 0; p < numPoints; ++p) { 1614 const PetscInt point = points[p]; 1615 1616 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1617 if (dep != dim-1) continue; 1618 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);} 1619 ierr = DMPlexVecSetClosure(dm, NULL, locF, point, &elemVec[f*totDimBd], ADD_ALL_VALUES);CHKERRQ(ierr); 1620 ++f; 1621 } 1622 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1623 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1624 ierr = PetscFree3(u,cgeom,elemVec);CHKERRQ(ierr); 1625 if (locX_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 1626 } 1627 } 1628 PetscFunctionReturn(0); 1629 } 1630 1631 #undef __FUNCT__ 1632 #define __FUNCT__ "DMPlexReconstructGradientsFVM" 1633 /*@ 1634 DMPlexReconstructGradientsFVM - reconstruct the gradient of a vector using a finite volume method. 1635 1636 Input Parameters: 1637 + dm - the mesh 1638 - locX - the local representation of the vector 1639 1640 Output Parameter: 1641 . grad - the global representation of the gradient 1642 1643 Level: developer 1644 1645 .seealso: DMPlexSNESGetGradientDM() 1646 @*/ 1647 PetscErrorCode DMPlexReconstructGradientsFVM(DM dm, Vec locX, Vec grad) 1648 { 1649 PetscDS prob; 1650 PetscInt Nf, f, fStart, fEnd; 1651 PetscBool useFVM = PETSC_FALSE; 1652 PetscFV fvm = NULL; 1653 Vec faceGeometryFVM, cellGeometryFVM; 1654 PetscFVCellGeom *cgeomFVM = NULL; 1655 PetscFVFaceGeom *fgeomFVM = NULL; 1656 DM dmGrad = NULL; 1657 PetscErrorCode ierr; 1658 1659 PetscFunctionBegin; 1660 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1661 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1662 for (f = 0; f < Nf; ++f) { 1663 PetscObject obj; 1664 PetscClassId id; 1665 1666 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1667 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1668 if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;} 1669 } 1670 if (!useFVM) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"This dm does not have a finite volume discretization"); 1671 ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr); 1672 if (!dmGrad) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"This dm's finite volume discretization does not reconstruct gradients"); 1673 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr); 1674 ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 1675 ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1676 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1677 ierr = DMPlexReconstructGradients_Internal(dm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr); 1678 PetscFunctionReturn(0); 1679 } 1680 1681 #undef __FUNCT__ 1682 #define __FUNCT__ "DMPlexComputeResidual_Internal" 1683 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscInt cStart, PetscInt cEnd, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) 1684 { 1685 DM_Plex *mesh = (DM_Plex *) dm->data; 1686 const char *name = "Residual"; 1687 DM dmAux = NULL; 1688 DM dmGrad = NULL; 1689 DMLabel ghostLabel = NULL; 1690 PetscDS prob = NULL; 1691 PetscDS probAux = NULL; 1692 PetscSection section = NULL; 1693 PetscBool useFEM = PETSC_FALSE; 1694 PetscBool useFVM = PETSC_FALSE; 1695 PetscBool isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE; 1696 PetscFV fvm = NULL; 1697 PetscFECellGeom *cgeomFEM = NULL; 1698 PetscScalar *cgeomScal; 1699 PetscFVCellGeom *cgeomFVM = NULL; 1700 PetscFVFaceGeom *fgeomFVM = NULL; 1701 Vec locA, cellGeometryFEM = NULL, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL; 1702 PetscScalar *u, *u_t, *a, *uL, *uR; 1703 PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd; 1704 PetscErrorCode ierr; 1705 1706 PetscFunctionBegin; 1707 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1708 /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */ 1709 /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */ 1710 /* FEM+FVM */ 1711 /* 1: Get sizes from dm and dmAux */ 1712 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1713 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1714 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1715 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1716 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1717 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1718 if (locA) { 1719 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1720 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1721 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1722 } 1723 /* 2: Get geometric data */ 1724 for (f = 0; f < Nf; ++f) { 1725 PetscObject obj; 1726 PetscClassId id; 1727 PetscBool fimp; 1728 1729 ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr); 1730 if (isImplicit != fimp) continue; 1731 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1732 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1733 if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;} 1734 if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;} 1735 } 1736 if (useFEM) { 1737 ierr = DMPlexSNESGetGeometryFEM(dm, &cellGeometryFEM);CHKERRQ(ierr); 1738 ierr = VecGetArray(cellGeometryFEM, &cgeomScal);CHKERRQ(ierr); 1739 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 1740 DM dmCell; 1741 PetscInt c; 1742 1743 ierr = VecGetDM(cellGeometryFEM,&dmCell);CHKERRQ(ierr); 1744 ierr = PetscMalloc1(cEnd-cStart,&cgeomFEM);CHKERRQ(ierr); 1745 for (c = 0; c < cEnd - cStart; c++) { 1746 PetscScalar *thisgeom; 1747 1748 ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr); 1749 cgeomFEM[c] = *((PetscFECellGeom *) thisgeom); 1750 } 1751 } 1752 else { 1753 cgeomFEM = (PetscFECellGeom *) cgeomScal; 1754 } 1755 } 1756 if (useFVM) { 1757 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr); 1758 ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 1759 ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1760 /* Reconstruct and limit cell gradients */ 1761 ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr); 1762 if (dmGrad) { 1763 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1764 ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1765 ierr = DMPlexReconstructGradients_Internal(dm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr); 1766 /* Communicate gradient values */ 1767 ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1768 ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1769 ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1770 ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1771 } 1772 /* Handle non-essential (e.g. outflow) boundary values */ 1773 ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr); 1774 } 1775 /* Loop over chunks */ 1776 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1777 numChunks = 1; 1778 cellChunkSize = (cEnd - cStart)/numChunks; 1779 faceChunkSize = (fEnd - fStart)/numChunks; 1780 for (chunk = 0; chunk < numChunks; ++chunk) { 1781 PetscScalar *elemVec, *fluxL, *fluxR; 1782 PetscReal *vol; 1783 PetscFVFaceGeom *fgeom; 1784 PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, cell; 1785 PetscInt fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = fE - fS, face; 1786 1787 /* Extract field coefficients */ 1788 if (useFEM) { 1789 ierr = DMPlexGetCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1790 ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr); 1791 ierr = PetscMemzero(elemVec, numCells*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1792 } 1793 if (useFVM) { 1794 ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr); 1795 ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr); 1796 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr); 1797 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr); 1798 ierr = PetscMemzero(fluxL, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1799 ierr = PetscMemzero(fluxR, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1800 } 1801 /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */ 1802 /* Loop over fields */ 1803 for (f = 0; f < Nf; ++f) { 1804 PetscObject obj; 1805 PetscClassId id; 1806 PetscBool fimp; 1807 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1808 1809 ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr); 1810 if (isImplicit != fimp) continue; 1811 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1812 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1813 if (id == PETSCFE_CLASSID) { 1814 PetscFE fe = (PetscFE) obj; 1815 PetscQuadrature q; 1816 PetscInt Nq, Nb; 1817 1818 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1819 1820 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1821 ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1822 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1823 blockSize = Nb*Nq; 1824 batchSize = numBlocks * blockSize; 1825 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1826 numChunks = numCells / (numBatches*batchSize); 1827 Ne = numChunks*numBatches*batchSize; 1828 Nr = numCells % (numBatches*batchSize); 1829 offset = numCells - Nr; 1830 /* Integrate FE residual to get elemVec (need fields at quadrature points) */ 1831 /* For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */ 1832 ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeomFEM, u, u_t, probAux, a, elemVec);CHKERRQ(ierr); 1833 ierr = PetscFEIntegrateResidual(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemVec[offset*totDim]);CHKERRQ(ierr); 1834 } else if (id == PETSCFV_CLASSID) { 1835 PetscFV fv = (PetscFV) obj; 1836 1837 Ne = numFaces; 1838 /* Riemann solve over faces (need fields at face centroids) */ 1839 /* We need to evaluate FE fields at those coordinates */ 1840 ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr); 1841 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1842 } 1843 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 1844 ierr = PetscFree(cgeomFEM);CHKERRQ(ierr); 1845 } 1846 else { 1847 cgeomFEM = NULL; 1848 } 1849 if (cellGeometryFEM) {ierr = VecRestoreArray(cellGeometryFEM, &cgeomScal);CHKERRQ(ierr);} 1850 /* Loop over domain */ 1851 if (useFEM) { 1852 /* Add elemVec to locX */ 1853 for (cell = cS; cell < cE; ++cell) { 1854 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cell*totDim]);CHKERRQ(ierr);} 1855 ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cell*totDim], ADD_ALL_VALUES);CHKERRQ(ierr); 1856 } 1857 } 1858 if (useFVM) { 1859 PetscScalar *fa; 1860 PetscInt iface; 1861 1862 ierr = VecGetArray(locF, &fa);CHKERRQ(ierr); 1863 for (f = 0; f < Nf; ++f) { 1864 PetscFV fv; 1865 PetscObject obj; 1866 PetscClassId id; 1867 PetscInt foff, pdim; 1868 1869 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1870 ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1871 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1872 if (id != PETSCFV_CLASSID) continue; 1873 fv = (PetscFV) obj; 1874 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 1875 /* Accumulate fluxes to cells */ 1876 for (face = fS, iface = 0; face < fE; ++face) { 1877 const PetscInt *cells; 1878 PetscScalar *fL, *fR; 1879 PetscInt ghost, d, nsupp; 1880 1881 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1882 ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr); 1883 if (ghost >= 0) continue; 1884 if (nsupp > 2) { /* noop */ 1885 ++iface; 1886 continue; 1887 } 1888 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1889 ierr = DMPlexPointGlobalFieldRef(dm, cells[0], f, fa, &fL);CHKERRQ(ierr); 1890 ierr = DMPlexPointGlobalFieldRef(dm, cells[1], f, fa, &fR);CHKERRQ(ierr); 1891 for (d = 0; d < pdim; ++d) { 1892 if (fL) fL[d] -= fluxL[iface*totDim+foff+d]; 1893 if (fR) fR[d] += fluxR[iface*totDim+foff+d]; 1894 } 1895 ++iface; 1896 } 1897 } 1898 ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr); 1899 } 1900 /* Handle time derivative */ 1901 if (locX_t) { 1902 PetscScalar *x_t, *fa; 1903 1904 ierr = VecGetArray(locF, &fa);CHKERRQ(ierr); 1905 ierr = VecGetArray(locX_t, &x_t);CHKERRQ(ierr); 1906 for (f = 0; f < Nf; ++f) { 1907 PetscFV fv; 1908 PetscObject obj; 1909 PetscClassId id; 1910 PetscInt pdim, d; 1911 1912 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1913 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1914 if (id != PETSCFV_CLASSID) continue; 1915 fv = (PetscFV) obj; 1916 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 1917 for (cell = cS; cell < cE; ++cell) { 1918 PetscScalar *u_t, *r; 1919 1920 ierr = DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);CHKERRQ(ierr); 1921 ierr = DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);CHKERRQ(ierr); 1922 for (d = 0; d < pdim; ++d) r[d] += u_t[d]; 1923 } 1924 } 1925 ierr = VecRestoreArray(locX_t, &x_t);CHKERRQ(ierr); 1926 ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr); 1927 } 1928 if (useFEM) { 1929 ierr = DMPlexRestoreCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1930 ierr = DMRestoreWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr); 1931 } 1932 if (useFVM) { 1933 ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr); 1934 ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr); 1935 ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr); 1936 ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr); 1937 if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);} 1938 } 1939 } 1940 1941 if (useFEM) {ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, locF, user);CHKERRQ(ierr);} 1942 1943 /* FEM */ 1944 /* 1: Get sizes from dm and dmAux */ 1945 /* 2: Get geometric data */ 1946 /* 3: Handle boundary values */ 1947 /* 4: Loop over domain */ 1948 /* Extract coefficients */ 1949 /* Loop over fields */ 1950 /* Set tiling for FE*/ 1951 /* Integrate FE residual to get elemVec */ 1952 /* Loop over subdomain */ 1953 /* Loop over quad points */ 1954 /* Transform coords to real space */ 1955 /* Evaluate field and aux fields at point */ 1956 /* Evaluate residual at point */ 1957 /* Transform residual to real space */ 1958 /* Add residual to elemVec */ 1959 /* Loop over domain */ 1960 /* Add elemVec to locX */ 1961 1962 /* FVM */ 1963 /* Get geometric data */ 1964 /* If using gradients */ 1965 /* Compute gradient data */ 1966 /* Loop over domain faces */ 1967 /* Count computational faces */ 1968 /* Reconstruct cell gradient */ 1969 /* Loop over domain cells */ 1970 /* Limit cell gradients */ 1971 /* Handle boundary values */ 1972 /* Loop over domain faces */ 1973 /* Read out field, centroid, normal, volume for each side of face */ 1974 /* Riemann solve over faces */ 1975 /* Loop over domain faces */ 1976 /* Accumulate fluxes to cells */ 1977 /* TODO Change printFEM to printDisc here */ 1978 if (mesh->printFEM) { 1979 Vec locFbc; 1980 PetscInt pStart, pEnd, p, maxDof; 1981 PetscScalar *zeroes; 1982 1983 ierr = VecDuplicate(locF,&locFbc);CHKERRQ(ierr); 1984 ierr = VecCopy(locF,locFbc);CHKERRQ(ierr); 1985 ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 1986 ierr = PetscSectionGetMaxDof(section,&maxDof);CHKERRQ(ierr); 1987 ierr = PetscCalloc1(maxDof,&zeroes);CHKERRQ(ierr); 1988 for (p = pStart; p < pEnd; p++) { 1989 ierr = VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);CHKERRQ(ierr); 1990 } 1991 ierr = PetscFree(zeroes);CHKERRQ(ierr); 1992 ierr = DMPrintLocalVec(dm, name, mesh->printTol, locFbc);CHKERRQ(ierr); 1993 ierr = VecDestroy(&locFbc);CHKERRQ(ierr); 1994 } 1995 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1996 PetscFunctionReturn(0); 1997 } 1998 1999 #undef __FUNCT__ 2000 #define __FUNCT__ "DMPlexComputeResidualFEM_Check_Internal" 2001 static PetscErrorCode DMPlexComputeResidualFEM_Check_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user) 2002 { 2003 DM dmCh, dmAux; 2004 Vec A, cellgeom; 2005 PetscDS prob, probCh, probAux = NULL; 2006 PetscQuadrature q; 2007 PetscSection section, sectionAux; 2008 PetscFECellGeom *cgeom = NULL; 2009 PetscScalar *cgeomScal; 2010 PetscScalar *elemVec, *elemVecCh, *u, *u_t, *a = NULL; 2011 PetscInt dim, Nf, f, numCells, cStart, cEnd, c; 2012 PetscInt totDim, totDimAux, diffCell = 0; 2013 PetscErrorCode ierr; 2014 2015 PetscFunctionBegin; 2016 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2017 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 2018 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2019 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2020 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 2021 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2022 numCells = cEnd - cStart; 2023 ierr = PetscObjectQuery((PetscObject) dm, "dmCh", (PetscObject *) &dmCh);CHKERRQ(ierr); 2024 ierr = DMGetDS(dmCh, &probCh);CHKERRQ(ierr); 2025 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 2026 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 2027 if (dmAux) { 2028 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 2029 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 2030 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 2031 } 2032 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 2033 ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim,&elemVec);CHKERRQ(ierr); 2034 ierr = PetscMalloc1(numCells*totDim,&elemVecCh);CHKERRQ(ierr); 2035 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 2036 ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 2037 ierr = VecGetArray(cellgeom, &cgeomScal);CHKERRQ(ierr); 2038 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 2039 DM dmCell; 2040 2041 ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr); 2042 ierr = PetscMalloc1(cEnd-cStart,&cgeom);CHKERRQ(ierr); 2043 for (c = 0; c < cEnd - cStart; c++) { 2044 PetscScalar *thisgeom; 2045 2046 ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr); 2047 cgeom[c] = *((PetscFECellGeom *) thisgeom); 2048 } 2049 } 2050 else { 2051 cgeom = (PetscFECellGeom *) cgeomScal; 2052 } 2053 for (c = cStart; c < cEnd; ++c) { 2054 PetscScalar *x = NULL, *x_t = NULL; 2055 PetscInt i; 2056 2057 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 2058 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 2059 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 2060 if (X_t) { 2061 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 2062 for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i]; 2063 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 2064 } 2065 if (dmAux) { 2066 DM dmAuxPlex; 2067 2068 ierr = DMSNESConvertPlex(dmAux,&dmAuxPlex, PETSC_FALSE);CHKERRQ(ierr); 2069 ierr = DMPlexVecGetClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 2070 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 2071 ierr = DMPlexVecRestoreClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 2072 ierr = DMDestroy(&dmAuxPlex);CHKERRQ(ierr); 2073 } 2074 } 2075 for (f = 0; f < Nf; ++f) { 2076 PetscFE fe, feCh; 2077 PetscInt numQuadPoints, Nb; 2078 /* Conforming batches */ 2079 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 2080 /* Remainder */ 2081 PetscInt Nr, offset; 2082 2083 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 2084 ierr = PetscDSGetDiscretization(probCh, f, (PetscObject *) &feCh);CHKERRQ(ierr); 2085 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 2086 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2087 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2088 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 2089 blockSize = Nb*numQuadPoints; 2090 batchSize = numBlocks * blockSize; 2091 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2092 numChunks = numCells / (numBatches*batchSize); 2093 Ne = numChunks*numBatches*batchSize; 2094 Nr = numCells % (numBatches*batchSize); 2095 offset = numCells - Nr; 2096 ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr); 2097 ierr = PetscFEIntegrateResidual(feCh, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVecCh);CHKERRQ(ierr); 2098 ierr = PetscFEIntegrateResidual(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemVec[offset*totDim]);CHKERRQ(ierr); 2099 ierr = PetscFEIntegrateResidual(feCh, prob, f, Nr, &cgeom[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemVecCh[offset*totDim]);CHKERRQ(ierr); 2100 } 2101 for (c = cStart; c < cEnd; ++c) { 2102 PetscBool diff = PETSC_FALSE; 2103 PetscInt d; 2104 2105 for (d = 0; d < totDim; ++d) if (PetscAbsScalar(elemVec[c*totDim+d] - elemVecCh[c*totDim+d]) > 1.0e-7) {diff = PETSC_TRUE;break;} 2106 if (diff) { 2107 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Different cell %d\n", c);CHKERRQ(ierr); 2108 ierr = DMPrintCellVector(c, "Residual", totDim, &elemVec[c*totDim]);CHKERRQ(ierr); 2109 ierr = DMPrintCellVector(c, "Check Residual", totDim, &elemVecCh[c*totDim]);CHKERRQ(ierr); 2110 ++diffCell; 2111 } 2112 if (diffCell > 9) break; 2113 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_ALL_VALUES);CHKERRQ(ierr); 2114 } 2115 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 2116 ierr = PetscFree(cgeom);CHKERRQ(ierr); 2117 } 2118 else { 2119 cgeom = NULL; 2120 } 2121 ierr = VecRestoreArray(cellgeom, &cgeomScal);CHKERRQ(ierr); 2122 ierr = PetscFree3(u,u_t,elemVec);CHKERRQ(ierr); 2123 ierr = PetscFree(elemVecCh);CHKERRQ(ierr); 2124 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 2125 PetscFunctionReturn(0); 2126 } 2127 2128 #undef __FUNCT__ 2129 #define __FUNCT__ "DMPlexSNESComputeResidualFEM" 2130 /*@ 2131 DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 2132 2133 Input Parameters: 2134 + dm - The mesh 2135 . X - Local solution 2136 - user - The user context 2137 2138 Output Parameter: 2139 . F - Local output vector 2140 2141 Level: developer 2142 2143 .seealso: DMPlexComputeJacobianActionFEM() 2144 @*/ 2145 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 2146 { 2147 PetscObject check; 2148 PetscInt cStart, cEnd, cEndInterior; 2149 DM plex; 2150 PetscErrorCode ierr; 2151 2152 PetscFunctionBegin; 2153 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 2154 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 2155 ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2156 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 2157 /* The dmCh is used to check two mathematically equivalent discretizations for computational equivalence */ 2158 ierr = PetscObjectQuery((PetscObject) plex, "dmCh", &check);CHKERRQ(ierr); 2159 if (check) {ierr = DMPlexComputeResidualFEM_Check_Internal(plex, X, NULL, F, user);CHKERRQ(ierr);} 2160 else {ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, PETSC_MIN_REAL, X, NULL, F, user);CHKERRQ(ierr);} 2161 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2162 PetscFunctionReturn(0); 2163 } 2164 2165 #undef __FUNCT__ 2166 #define __FUNCT__ "DMPlexSNESComputeBoundaryFEM" 2167 /*@ 2168 DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X 2169 2170 Input Parameters: 2171 + dm - The mesh 2172 - user - The user context 2173 2174 Output Parameter: 2175 . X - Local solution 2176 2177 Level: developer 2178 2179 .seealso: DMPlexComputeJacobianActionFEM() 2180 @*/ 2181 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) 2182 { 2183 DM plex; 2184 PetscErrorCode ierr; 2185 2186 PetscFunctionBegin; 2187 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 2188 ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr); 2189 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2190 PetscFunctionReturn(0); 2191 } 2192 2193 #undef __FUNCT__ 2194 #define __FUNCT__ "DMPlexComputeJacobian_Internal" 2195 PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, PetscInt cStart, PetscInt cEnd, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user) 2196 { 2197 DM_Plex *mesh = (DM_Plex *) dm->data; 2198 const char *name = "Jacobian"; 2199 DM dmAux, plex; 2200 DMLabel depth; 2201 Vec A, cellgeom; 2202 PetscDS prob, probAux = NULL; 2203 PetscQuadrature quad; 2204 PetscSection section, globalSection, sectionAux; 2205 PetscFECellGeom *cgeom = NULL; 2206 PetscScalar *cgeomScal; 2207 PetscScalar *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL; 2208 PetscInt dim, Nf, f, fieldI, fieldJ, numCells, c; 2209 PetscInt totDim, totDimBd, totDimAux, numBd, bd; 2210 PetscBool isShell, hasPrec, hasDyn; 2211 PetscErrorCode ierr; 2212 2213 PetscFunctionBegin; 2214 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2215 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2216 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 2217 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 2218 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2219 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 2220 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 2221 ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr); 2222 ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr); 2223 hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE; 2224 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 2225 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2226 numCells = cEnd - cStart; 2227 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 2228 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 2229 if (dmAux) { 2230 ierr = DMConvert(dmAux, DMPLEX, &plex);CHKERRQ(ierr); 2231 ierr = DMGetDefaultSection(plex, §ionAux);CHKERRQ(ierr); 2232 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 2233 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 2234 } 2235 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 2236 ierr = PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);CHKERRQ(ierr); 2237 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 2238 ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 2239 ierr = VecGetArray(cellgeom, &cgeomScal);CHKERRQ(ierr); 2240 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 2241 DM dmCell; 2242 2243 ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr); 2244 ierr = PetscMalloc1(cEnd-cStart,&cgeom);CHKERRQ(ierr); 2245 for (c = 0; c < cEnd - cStart; c++) { 2246 PetscScalar *thisgeom; 2247 2248 ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr); 2249 cgeom[c] = *((PetscFECellGeom *) thisgeom); 2250 } 2251 } 2252 else { 2253 cgeom = (PetscFECellGeom *) cgeomScal; 2254 } 2255 for (c = cStart; c < cEnd; ++c) { 2256 PetscScalar *x = NULL, *x_t = NULL; 2257 PetscInt i; 2258 2259 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 2260 for (i = 0; i < totDim; ++i) u[(c-cStart)*totDim+i] = x[i]; 2261 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 2262 if (X_t) { 2263 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 2264 for (i = 0; i < totDim; ++i) u_t[(c-cStart)*totDim+i] = x_t[i]; 2265 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 2266 } 2267 if (dmAux) { 2268 ierr = DMPlexVecGetClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 2269 for (i = 0; i < totDimAux; ++i) a[(c-cStart)*totDimAux+i] = x[i]; 2270 ierr = DMPlexVecRestoreClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 2271 } 2272 } 2273 ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 2274 if (hasPrec) {ierr = PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);} 2275 if (hasDyn) {ierr = PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);} 2276 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2277 PetscFE fe; 2278 PetscInt numQuadPoints, Nb; 2279 /* Conforming batches */ 2280 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 2281 /* Remainder */ 2282 PetscInt Nr, offset; 2283 2284 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 2285 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 2286 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2287 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2288 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 2289 blockSize = Nb*numQuadPoints; 2290 batchSize = numBlocks * blockSize; 2291 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2292 numChunks = numCells / (numBatches*batchSize); 2293 Ne = numChunks*numBatches*batchSize; 2294 Nr = numCells % (numBatches*batchSize); 2295 offset = numCells - Nr; 2296 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2297 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr); 2298 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, &cgeom[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemMat[offset*totDim*totDim]);CHKERRQ(ierr); 2299 if (hasPrec) { 2300 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, elemMatP);CHKERRQ(ierr); 2301 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Nr, &cgeom[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemMatP[offset*totDim*totDim]);CHKERRQ(ierr); 2302 } 2303 if (hasDyn) { 2304 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, elemMatD);CHKERRQ(ierr); 2305 ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, &cgeom[offset], &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr); 2306 } 2307 } 2308 } 2309 if (hasDyn) { 2310 for (c = 0; c < (cEnd - cStart)*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c]; 2311 } 2312 for (c = cStart; c < cEnd; ++c) { 2313 if (hasPrec) { 2314 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);} 2315 ierr = DMPlexMatSetClosure(dm, section, globalSection, Jac, c, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2316 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMatP[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);} 2317 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMatP[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2318 } else { 2319 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);} 2320 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 2321 } 2322 } 2323 if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) { 2324 ierr = PetscFree(cgeom);CHKERRQ(ierr); 2325 } 2326 else { 2327 cgeom = NULL; 2328 } 2329 ierr = VecRestoreArray(cellgeom, &cgeomScal);CHKERRQ(ierr); 2330 ierr = PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);CHKERRQ(ierr); 2331 if (dmAux) { 2332 ierr = PetscFree(a);CHKERRQ(ierr); 2333 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2334 } 2335 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 2336 ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 2337 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 2338 ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 2339 for (bd = 0; bd < numBd; ++bd) { 2340 const char *bdLabel; 2341 DMLabel label; 2342 IS pointIS; 2343 const PetscInt *points; 2344 const PetscInt *values; 2345 PetscInt field, numValues, v, numPoints, p, dep, numFaces; 2346 PetscBool isEssential; 2347 PetscObject obj; 2348 PetscClassId id; 2349 2350 ierr = DMGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 2351 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 2352 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 2353 if ((id != PETSCFE_CLASSID) || isEssential) continue; 2354 ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 2355 for (v = 0; v < numValues; ++v) { 2356 ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr); 2357 ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr); 2358 if (!pointIS) continue; /* No points with that id on this process */ 2359 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 2360 for (p = 0, numFaces = 0; p < numPoints; ++p) { 2361 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 2362 if (dep == dim-1) ++numFaces; 2363 } 2364 ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr); 2365 if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 2366 for (p = 0, f = 0; p < numPoints; ++p) { 2367 const PetscInt point = points[p]; 2368 PetscScalar *x = NULL; 2369 PetscInt i; 2370 2371 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 2372 if (dep != dim-1) continue; 2373 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr); 2374 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr); 2375 if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point); 2376 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 2377 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 2378 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 2379 if (X_t) { 2380 ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 2381 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 2382 ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 2383 } 2384 ++f; 2385 } 2386 ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr); 2387 for (fieldI = 0; fieldI < Nf; ++fieldI) { 2388 PetscFE fe; 2389 PetscInt numQuadPoints, Nb; 2390 /* Conforming batches */ 2391 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 2392 /* Remainder */ 2393 PetscInt Nr, offset; 2394 2395 ierr = PetscDSGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 2396 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 2397 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 2398 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 2399 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 2400 blockSize = Nb*numQuadPoints; 2401 batchSize = numBlocks * blockSize; 2402 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 2403 numChunks = numFaces / (numBatches*batchSize); 2404 Ne = numChunks*numBatches*batchSize; 2405 Nr = numFaces % (numBatches*batchSize); 2406 offset = numFaces - Nr; 2407 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2408 ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr); 2409 ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Nr, &cgeom[offset], &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemMat[offset*totDimBd*totDimBd]);CHKERRQ(ierr); 2410 } 2411 } 2412 for (p = 0, f = 0; p < numPoints; ++p) { 2413 const PetscInt point = points[p]; 2414 2415 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 2416 if (dep != dim-1) continue; 2417 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);} 2418 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr); 2419 ++f; 2420 } 2421 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 2422 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 2423 ierr = PetscFree3(u,cgeom,elemMat);CHKERRQ(ierr); 2424 if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 2425 } 2426 } 2427 if (hasPrec) { 2428 ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2429 ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2430 } 2431 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2432 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2433 if (mesh->printFEM) { 2434 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 2435 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 2436 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 2437 } 2438 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2439 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 2440 if (isShell) { 2441 JacActionCtx *jctx; 2442 2443 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 2444 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 2445 } 2446 PetscFunctionReturn(0); 2447 } 2448 2449 #undef __FUNCT__ 2450 #define __FUNCT__ "DMPlexSNESComputeJacobianFEM" 2451 /*@ 2452 DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 2453 2454 Input Parameters: 2455 + dm - The mesh 2456 . X - Local input vector 2457 - user - The user context 2458 2459 Output Parameter: 2460 . Jac - Jacobian matrix 2461 2462 Note: 2463 The first member of the user context must be an FEMContext. 2464 2465 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 2466 like a GPU, or vectorize on a multicore machine. 2467 2468 Level: developer 2469 2470 .seealso: FormFunctionLocal() 2471 @*/ 2472 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 2473 { 2474 PetscInt cStart, cEnd, cEndInterior; 2475 DM plex; 2476 PetscErrorCode ierr; 2477 2478 PetscFunctionBegin; 2479 ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 2480 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 2481 ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2482 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 2483 ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr); 2484 ierr = DMDestroy(&plex);CHKERRQ(ierr); 2485 PetscFunctionReturn(0); 2486 } 2487 2488 #undef __FUNCT__ 2489 #define __FUNCT__ "DMPlexSetSNESLocalFEM" 2490 /*@ 2491 DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian. 2492 2493 Input Parameters: 2494 + dm - The DM object 2495 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see DMAddBoundary()) 2496 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual()) 2497 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian()) 2498 2499 Level: developer 2500 @*/ 2501 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) 2502 { 2503 PetscErrorCode ierr; 2504 2505 PetscFunctionBegin; 2506 ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr); 2507 ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr); 2508 ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr); 2509 PetscFunctionReturn(0); 2510 } 2511 2512 #undef __FUNCT__ 2513 #define __FUNCT__ "DMSNESCheckFromOptions_Internal" 2514 PetscErrorCode DMSNESCheckFromOptions_Internal(SNES snes, DM dm, Vec u, Vec sol, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs) 2515 { 2516 Mat J, M; 2517 Vec r, b; 2518 MatNullSpace nullSpace; 2519 PetscReal *error, res = 0.0; 2520 PetscInt numFields; 2521 PetscErrorCode ierr; 2522 2523 PetscFunctionBegin; 2524 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 2525 ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 2526 M = J; 2527 /* TODO Null space for J */ 2528 /* Check discretization error */ 2529 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 2530 ierr = PetscMalloc1(PetscMax(1, numFields), &error);CHKERRQ(ierr); 2531 if (numFields > 1) { 2532 PetscInt f; 2533 2534 ierr = DMComputeL2FieldDiff(dm, 0.0, exactFuncs, ctxs, u, error);CHKERRQ(ierr); 2535 ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: [");CHKERRQ(ierr); 2536 for (f = 0; f < numFields; ++f) { 2537 if (f) {ierr = PetscPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);} 2538 if (error[f] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "%g", error[f]);CHKERRQ(ierr);} 2539 else {ierr = PetscPrintf(PETSC_COMM_WORLD, "< 1.0e-11");CHKERRQ(ierr);} 2540 } 2541 ierr = PetscPrintf(PETSC_COMM_WORLD, "]\n");CHKERRQ(ierr); 2542 } else { 2543 ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, ctxs, u, &error[0]);CHKERRQ(ierr); 2544 if (error[0] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error[0]);CHKERRQ(ierr);} 2545 else {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 2546 } 2547 ierr = PetscFree(error);CHKERRQ(ierr); 2548 /* Check residual */ 2549 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 2550 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 2551 ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);CHKERRQ(ierr); 2552 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 2553 ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr); 2554 ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr); 2555 ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 2556 /* Check Jacobian */ 2557 ierr = SNESComputeJacobian(snes, u, M, M);CHKERRQ(ierr); 2558 ierr = MatGetNullSpace(J, &nullSpace);CHKERRQ(ierr); 2559 if (nullSpace) { 2560 PetscBool isNull; 2561 ierr = MatNullSpaceTest(nullSpace, J, &isNull);CHKERRQ(ierr); 2562 if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 2563 } 2564 ierr = VecDuplicate(u, &b);CHKERRQ(ierr); 2565 ierr = VecSet(r, 0.0);CHKERRQ(ierr); 2566 ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr); 2567 ierr = MatMult(M, u, r);CHKERRQ(ierr); 2568 ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr); 2569 ierr = VecDestroy(&b);CHKERRQ(ierr); 2570 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 2571 ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);CHKERRQ(ierr); 2572 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 2573 ierr = PetscObjectSetName((PetscObject) r, "Au - b = Au + F(0)");CHKERRQ(ierr); 2574 ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"linear_res_");CHKERRQ(ierr); 2575 ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 2576 ierr = VecDestroy(&r);CHKERRQ(ierr); 2577 ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr); 2578 ierr = MatDestroy(&J);CHKERRQ(ierr); 2579 PetscFunctionReturn(0); 2580 } 2581 2582 #undef __FUNCT__ 2583 #define __FUNCT__ "DMSNESCheckFromOptions" 2584 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs) 2585 { 2586 DM dm; 2587 Vec sol; 2588 PetscBool check; 2589 PetscErrorCode ierr; 2590 2591 PetscFunctionBegin; 2592 ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr); 2593 if (!check) PetscFunctionReturn(0); 2594 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 2595 ierr = VecDuplicate(u, &sol);CHKERRQ(ierr); 2596 ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr); 2597 ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr); 2598 ierr = VecDestroy(&sol);CHKERRQ(ierr); 2599 PetscFunctionReturn(0); 2600 } 2601