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 6 /************************** Interpolation *******************************/ 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "DMInterpolationCreate" 10 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 11 { 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 PetscValidPointer(ctx, 2); 16 ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr); 17 18 (*ctx)->comm = comm; 19 (*ctx)->dim = -1; 20 (*ctx)->nInput = 0; 21 (*ctx)->points = NULL; 22 (*ctx)->cells = NULL; 23 (*ctx)->n = -1; 24 (*ctx)->coords = NULL; 25 PetscFunctionReturn(0); 26 } 27 28 #undef __FUNCT__ 29 #define __FUNCT__ "DMInterpolationSetDim" 30 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 31 { 32 PetscFunctionBegin; 33 if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim); 34 ctx->dim = dim; 35 PetscFunctionReturn(0); 36 } 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "DMInterpolationGetDim" 40 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 41 { 42 PetscFunctionBegin; 43 PetscValidIntPointer(dim, 2); 44 *dim = ctx->dim; 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "DMInterpolationSetDof" 50 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 51 { 52 PetscFunctionBegin; 53 if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof); 54 ctx->dof = dof; 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "DMInterpolationGetDof" 60 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 61 { 62 PetscFunctionBegin; 63 PetscValidIntPointer(dof, 2); 64 *dof = ctx->dof; 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "DMInterpolationAddPoints" 70 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 71 { 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 76 if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 77 ctx->nInput = n; 78 79 ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr); 80 ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr); 81 PetscFunctionReturn(0); 82 } 83 84 #undef __FUNCT__ 85 #define __FUNCT__ "DMInterpolationSetUp" 86 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) 87 { 88 MPI_Comm comm = ctx->comm; 89 PetscScalar *a; 90 PetscInt p, q, i; 91 PetscMPIInt rank, size; 92 PetscErrorCode ierr; 93 Vec pointVec; 94 IS cellIS; 95 PetscLayout layout; 96 PetscReal *globalPoints; 97 PetscScalar *globalPointsScalar; 98 const PetscInt *ranges; 99 PetscMPIInt *counts, *displs; 100 const PetscInt *foundCells; 101 PetscMPIInt *foundProcs, *globalProcs; 102 PetscInt n, N; 103 104 PetscFunctionBegin; 105 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 106 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 107 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 108 if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 109 /* Locate points */ 110 n = ctx->nInput; 111 if (!redundantPoints) { 112 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 113 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 114 ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 115 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 116 ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 117 /* Communicate all points to all processes */ 118 ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr); 119 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 120 for (p = 0; p < size; ++p) { 121 counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 122 displs[p] = ranges[p]*ctx->dim; 123 } 124 ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 125 } else { 126 N = n; 127 globalPoints = ctx->points; 128 counts = displs = NULL; 129 layout = NULL; 130 } 131 #if 0 132 ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 133 /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 134 #else 135 #if defined(PETSC_USE_COMPLEX) 136 ierr = PetscMalloc1(N,&globalPointsScalar);CHKERRQ(ierr); 137 for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i]; 138 #else 139 globalPointsScalar = globalPoints; 140 #endif 141 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr); 142 ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 143 ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr); 144 ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr); 145 #endif 146 for (p = 0; p < N; ++p) { 147 if (foundCells[p] >= 0) foundProcs[p] = rank; 148 else foundProcs[p] = size; 149 } 150 /* Let the lowest rank process own each point */ 151 ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 152 ctx->n = 0; 153 for (p = 0; p < N; ++p) { 154 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); 155 else if (globalProcs[p] == rank) ctx->n++; 156 } 157 /* Create coordinates vector and array of owned cells */ 158 ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr); 159 ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 160 ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 161 ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 162 ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr); 163 ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 164 for (p = 0, q = 0, i = 0; p < N; ++p) { 165 if (globalProcs[p] == rank) { 166 PetscInt d; 167 168 for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d]; 169 ctx->cells[q++] = foundCells[p]; 170 } 171 } 172 ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 173 #if 0 174 ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 175 #else 176 ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 177 ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr); 178 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 179 ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 180 #endif 181 if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);} 182 if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);} 183 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "DMInterpolationGetCoordinates" 189 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 190 { 191 PetscFunctionBegin; 192 PetscValidPointer(coordinates, 2); 193 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 194 *coordinates = ctx->coords; 195 PetscFunctionReturn(0); 196 } 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "DMInterpolationGetVector" 200 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 PetscValidPointer(v, 2); 206 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 207 ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 208 ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 209 ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 210 ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr); 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "DMInterpolationRestoreVector" 216 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 217 { 218 PetscErrorCode ierr; 219 220 PetscFunctionBegin; 221 PetscValidPointer(v, 2); 222 if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 223 ierr = VecDestroy(v);CHKERRQ(ierr); 224 PetscFunctionReturn(0); 225 } 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "DMInterpolate_Triangle_Private" 229 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 230 { 231 PetscReal *v0, *J, *invJ, detJ; 232 PetscScalar *a, *coords; 233 PetscInt p; 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 238 ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 239 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 240 for (p = 0; p < ctx->n; ++p) { 241 PetscInt c = ctx->cells[p]; 242 PetscScalar *x = NULL; 243 PetscReal xi[4]; 244 PetscInt d, f, comp; 245 246 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 247 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 248 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 249 for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 250 251 for (d = 0; d < ctx->dim; ++d) { 252 xi[d] = 0.0; 253 for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 254 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]; 255 } 256 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 257 } 258 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 259 ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 260 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "DMInterpolate_Tetrahedron_Private" 266 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 267 { 268 PetscReal *v0, *J, *invJ, detJ; 269 PetscScalar *a, *coords; 270 PetscInt p; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 275 ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 276 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 277 for (p = 0; p < ctx->n; ++p) { 278 PetscInt c = ctx->cells[p]; 279 const PetscInt order[3] = {2, 1, 3}; 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[order[d]*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 = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 298 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "QuadMap_Private" 304 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 305 { 306 const PetscScalar *vertices = (const PetscScalar*) ctx; 307 const PetscScalar x0 = vertices[0]; 308 const PetscScalar y0 = vertices[1]; 309 const PetscScalar x1 = vertices[2]; 310 const PetscScalar y1 = vertices[3]; 311 const PetscScalar x2 = vertices[4]; 312 const PetscScalar y2 = vertices[5]; 313 const PetscScalar x3 = vertices[6]; 314 const PetscScalar y3 = vertices[7]; 315 const PetscScalar f_1 = x1 - x0; 316 const PetscScalar g_1 = y1 - y0; 317 const PetscScalar f_3 = x3 - x0; 318 const PetscScalar g_3 = y3 - y0; 319 const PetscScalar f_01 = x2 - x1 - x3 + x0; 320 const PetscScalar g_01 = y2 - y1 - y3 + y0; 321 PetscScalar *ref, *real; 322 PetscErrorCode ierr; 323 324 PetscFunctionBegin; 325 ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 326 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 327 { 328 const PetscScalar p0 = ref[0]; 329 const PetscScalar p1 = ref[1]; 330 331 real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 332 real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 333 } 334 ierr = PetscLogFlops(28);CHKERRQ(ierr); 335 ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 336 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 337 PetscFunctionReturn(0); 338 } 339 340 #include <petsc-private/dmimpl.h> 341 #undef __FUNCT__ 342 #define __FUNCT__ "QuadJacobian_Private" 343 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, 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_01 = x2 - x1 - x3 + x0; 355 const PetscScalar g_01 = y2 - y1 - y3 + y0; 356 PetscScalar *ref; 357 PetscErrorCode ierr; 358 359 PetscFunctionBegin; 360 ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 361 { 362 const PetscScalar x = ref[0]; 363 const PetscScalar y = ref[1]; 364 const PetscInt rows[2] = {0, 1}; 365 PetscScalar values[4]; 366 367 values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5; 368 values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5; 369 ierr = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 370 } 371 ierr = PetscLogFlops(30);CHKERRQ(ierr); 372 ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 373 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 374 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 #undef __FUNCT__ 379 #define __FUNCT__ "DMInterpolate_Quad_Private" 380 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 381 { 382 DM dmCoord; 383 SNES snes; 384 KSP ksp; 385 PC pc; 386 Vec coordsLocal, r, ref, real; 387 Mat J; 388 PetscScalar *a, *coords; 389 PetscInt p; 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 394 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 395 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 396 ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 397 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 398 ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 399 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 400 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 401 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 402 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 403 ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 404 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 405 ierr = MatSetUp(J);CHKERRQ(ierr); 406 ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr); 407 ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr); 408 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 409 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 410 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 411 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 412 413 ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 414 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 415 for (p = 0; p < ctx->n; ++p) { 416 PetscScalar *x = NULL, *vertices = NULL; 417 PetscScalar *xi; 418 PetscReal xir[2]; 419 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 420 421 /* Can make this do all points at once */ 422 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 423 if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2); 424 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 425 if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof); 426 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 427 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 428 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 429 xi[0] = coords[p*ctx->dim+0]; 430 xi[1] = coords[p*ctx->dim+1]; 431 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 432 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 433 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 434 xir[0] = PetscRealPart(xi[0]); 435 xir[1] = PetscRealPart(xi[1]); 436 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]; 437 438 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 439 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 440 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 441 } 442 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 443 ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 444 445 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 446 ierr = VecDestroy(&r);CHKERRQ(ierr); 447 ierr = VecDestroy(&ref);CHKERRQ(ierr); 448 ierr = VecDestroy(&real);CHKERRQ(ierr); 449 ierr = MatDestroy(&J);CHKERRQ(ierr); 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "HexMap_Private" 455 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 456 { 457 const PetscScalar *vertices = (const PetscScalar*) ctx; 458 const PetscScalar x0 = vertices[0]; 459 const PetscScalar y0 = vertices[1]; 460 const PetscScalar z0 = vertices[2]; 461 const PetscScalar x1 = vertices[9]; 462 const PetscScalar y1 = vertices[10]; 463 const PetscScalar z1 = vertices[11]; 464 const PetscScalar x2 = vertices[6]; 465 const PetscScalar y2 = vertices[7]; 466 const PetscScalar z2 = vertices[8]; 467 const PetscScalar x3 = vertices[3]; 468 const PetscScalar y3 = vertices[4]; 469 const PetscScalar z3 = vertices[5]; 470 const PetscScalar x4 = vertices[12]; 471 const PetscScalar y4 = vertices[13]; 472 const PetscScalar z4 = vertices[14]; 473 const PetscScalar x5 = vertices[15]; 474 const PetscScalar y5 = vertices[16]; 475 const PetscScalar z5 = vertices[17]; 476 const PetscScalar x6 = vertices[18]; 477 const PetscScalar y6 = vertices[19]; 478 const PetscScalar z6 = vertices[20]; 479 const PetscScalar x7 = vertices[21]; 480 const PetscScalar y7 = vertices[22]; 481 const PetscScalar z7 = vertices[23]; 482 const PetscScalar f_1 = x1 - x0; 483 const PetscScalar g_1 = y1 - y0; 484 const PetscScalar h_1 = z1 - z0; 485 const PetscScalar f_3 = x3 - x0; 486 const PetscScalar g_3 = y3 - y0; 487 const PetscScalar h_3 = z3 - z0; 488 const PetscScalar f_4 = x4 - x0; 489 const PetscScalar g_4 = y4 - y0; 490 const PetscScalar h_4 = z4 - z0; 491 const PetscScalar f_01 = x2 - x1 - x3 + x0; 492 const PetscScalar g_01 = y2 - y1 - y3 + y0; 493 const PetscScalar h_01 = z2 - z1 - z3 + z0; 494 const PetscScalar f_12 = x7 - x3 - x4 + x0; 495 const PetscScalar g_12 = y7 - y3 - y4 + y0; 496 const PetscScalar h_12 = z7 - z3 - z4 + z0; 497 const PetscScalar f_02 = x5 - x1 - x4 + x0; 498 const PetscScalar g_02 = y5 - y1 - y4 + y0; 499 const PetscScalar h_02 = z5 - z1 - z4 + z0; 500 const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 501 const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 502 const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 503 PetscScalar *ref, *real; 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; 507 ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 508 ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 509 { 510 const PetscScalar p0 = ref[0]; 511 const PetscScalar p1 = ref[1]; 512 const PetscScalar p2 = ref[2]; 513 514 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; 515 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; 516 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; 517 } 518 ierr = PetscLogFlops(114);CHKERRQ(ierr); 519 ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 520 ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 521 PetscFunctionReturn(0); 522 } 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "HexJacobian_Private" 526 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 527 { 528 const PetscScalar *vertices = (const PetscScalar*) ctx; 529 const PetscScalar x0 = vertices[0]; 530 const PetscScalar y0 = vertices[1]; 531 const PetscScalar z0 = vertices[2]; 532 const PetscScalar x1 = vertices[9]; 533 const PetscScalar y1 = vertices[10]; 534 const PetscScalar z1 = vertices[11]; 535 const PetscScalar x2 = vertices[6]; 536 const PetscScalar y2 = vertices[7]; 537 const PetscScalar z2 = vertices[8]; 538 const PetscScalar x3 = vertices[3]; 539 const PetscScalar y3 = vertices[4]; 540 const PetscScalar z3 = vertices[5]; 541 const PetscScalar x4 = vertices[12]; 542 const PetscScalar y4 = vertices[13]; 543 const PetscScalar z4 = vertices[14]; 544 const PetscScalar x5 = vertices[15]; 545 const PetscScalar y5 = vertices[16]; 546 const PetscScalar z5 = vertices[17]; 547 const PetscScalar x6 = vertices[18]; 548 const PetscScalar y6 = vertices[19]; 549 const PetscScalar z6 = vertices[20]; 550 const PetscScalar x7 = vertices[21]; 551 const PetscScalar y7 = vertices[22]; 552 const PetscScalar z7 = vertices[23]; 553 const PetscScalar f_xy = x2 - x1 - x3 + x0; 554 const PetscScalar g_xy = y2 - y1 - y3 + y0; 555 const PetscScalar h_xy = z2 - z1 - z3 + z0; 556 const PetscScalar f_yz = x7 - x3 - x4 + x0; 557 const PetscScalar g_yz = y7 - y3 - y4 + y0; 558 const PetscScalar h_yz = z7 - z3 - z4 + z0; 559 const PetscScalar f_xz = x5 - x1 - x4 + x0; 560 const PetscScalar g_xz = y5 - y1 - y4 + y0; 561 const PetscScalar h_xz = z5 - z1 - z4 + z0; 562 const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 563 const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 564 const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 565 PetscScalar *ref; 566 PetscErrorCode ierr; 567 568 PetscFunctionBegin; 569 ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 570 { 571 const PetscScalar x = ref[0]; 572 const PetscScalar y = ref[1]; 573 const PetscScalar z = ref[2]; 574 const PetscInt rows[3] = {0, 1, 2}; 575 PetscScalar values[9]; 576 577 values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0; 578 values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0; 579 values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0; 580 values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0; 581 values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0; 582 values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0; 583 values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0; 584 values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0; 585 values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0; 586 587 ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 588 } 589 ierr = PetscLogFlops(152);CHKERRQ(ierr); 590 ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 591 ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 592 ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 #undef __FUNCT__ 597 #define __FUNCT__ "DMInterpolate_Hex_Private" 598 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 599 { 600 DM dmCoord; 601 SNES snes; 602 KSP ksp; 603 PC pc; 604 Vec coordsLocal, r, ref, real; 605 Mat J; 606 PetscScalar *a, *coords; 607 PetscInt p; 608 PetscErrorCode ierr; 609 610 PetscFunctionBegin; 611 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 612 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 613 ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 614 ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 615 ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 616 ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 617 ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 618 ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 619 ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 620 ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 621 ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 622 ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 623 ierr = MatSetUp(J);CHKERRQ(ierr); 624 ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr); 625 ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr); 626 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 627 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 628 ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 629 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 630 631 ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 632 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 633 for (p = 0; p < ctx->n; ++p) { 634 PetscScalar *x = NULL, *vertices = NULL; 635 PetscScalar *xi; 636 PetscReal xir[3]; 637 PetscInt c = ctx->cells[p], comp, coordSize, xSize; 638 639 /* Can make this do all points at once */ 640 ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 641 if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3); 642 ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 643 if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof); 644 ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 645 ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 646 ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 647 xi[0] = coords[p*ctx->dim+0]; 648 xi[1] = coords[p*ctx->dim+1]; 649 xi[2] = coords[p*ctx->dim+2]; 650 ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 651 ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 652 ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 653 xir[0] = PetscRealPart(xi[0]); 654 xir[1] = PetscRealPart(xi[1]); 655 xir[2] = PetscRealPart(xi[2]); 656 for (comp = 0; comp < ctx->dof; ++comp) { 657 a[p*ctx->dof+comp] = 658 x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) + 659 x[3*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) + 660 x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) + 661 x[1*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) + 662 x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] + 663 x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] + 664 x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] + 665 x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2]; 666 } 667 ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 668 ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 669 ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 670 } 671 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 672 ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 673 674 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 675 ierr = VecDestroy(&r);CHKERRQ(ierr); 676 ierr = VecDestroy(&ref);CHKERRQ(ierr); 677 ierr = VecDestroy(&real);CHKERRQ(ierr); 678 ierr = MatDestroy(&J);CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "DMInterpolationEvaluate" 684 /* 685 Input Parameters: 686 + ctx - The DMInterpolationInfo context 687 . dm - The DM 688 - x - The local vector containing the field to be interpolated 689 690 Output Parameters: 691 . v - The vector containing the interpolated values 692 */ 693 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 694 { 695 PetscInt dim, coneSize, n; 696 PetscErrorCode ierr; 697 698 PetscFunctionBegin; 699 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 700 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 701 PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 702 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 703 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); 704 if (n) { 705 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 706 ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 707 if (dim == 2) { 708 if (coneSize == 3) { 709 ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr); 710 } else if (coneSize == 4) { 711 ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 712 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 713 } else if (dim == 3) { 714 if (coneSize == 4) { 715 ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr); 716 } else { 717 ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 718 } 719 } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 720 } 721 PetscFunctionReturn(0); 722 } 723 724 #undef __FUNCT__ 725 #define __FUNCT__ "DMInterpolationDestroy" 726 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 727 { 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 PetscValidPointer(ctx, 2); 732 ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 733 ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 734 ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 735 ierr = PetscFree(*ctx);CHKERRQ(ierr); 736 *ctx = NULL; 737 PetscFunctionReturn(0); 738 } 739 740 #undef __FUNCT__ 741 #define __FUNCT__ "SNESMonitorFields" 742 /*@C 743 SNESMonitorFields - Monitors the residual for each field separately 744 745 Collective on SNES 746 747 Input Parameters: 748 + snes - the SNES context 749 . its - iteration number 750 . fgnorm - 2-norm of residual 751 - dummy - unused context 752 753 Notes: 754 This routine prints the residual norm at each iteration. 755 756 Level: intermediate 757 758 .keywords: SNES, nonlinear, default, monitor, norm 759 .seealso: SNESMonitorSet(), SNESMonitorDefault() 760 @*/ 761 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 762 { 763 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject) snes)); 764 Vec res; 765 DM dm; 766 PetscSection s; 767 const PetscScalar *r; 768 PetscReal *lnorms, *norms; 769 PetscInt numFields, f, pStart, pEnd, p; 770 PetscErrorCode ierr; 771 772 PetscFunctionBegin; 773 ierr = SNESGetFunction(snes, &res, 0, 0);CHKERRQ(ierr); 774 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 775 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 776 ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr); 777 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 778 ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr); 779 ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr); 780 for (p = pStart; p < pEnd; ++p) { 781 for (f = 0; f < numFields; ++f) { 782 PetscInt fdof, foff, d; 783 784 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 785 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 786 for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d])); 787 } 788 } 789 ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr); 790 ierr = MPI_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 791 ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 792 ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr); 793 for (f = 0; f < numFields; ++f) { 794 if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 795 ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr); 796 } 797 ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr); 798 ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr); 799 ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr); 800 PetscFunctionReturn(0); 801 } 802 803 /********************* Residual Computation **************************/ 804 805 #undef __FUNCT__ 806 #define __FUNCT__ "DMPlexSNESGetGeometryFEM" 807 /*@ 808 DMPlexSNESGetGeometryFEM - Return precomputed geometric data 809 810 Input Parameter: 811 . dm - The DM 812 813 Output Parameters: 814 . cellgeom - The values precomputed from cell geometry 815 816 Level: developer 817 818 .seealso: DMPlexSNESSetFunctionLocal() 819 @*/ 820 PetscErrorCode DMPlexSNESGetGeometryFEM(DM dm, Vec *cellgeom) 821 { 822 DMSNES dmsnes; 823 PetscObject obj; 824 PetscErrorCode ierr; 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 828 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 829 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", &obj);CHKERRQ(ierr); 830 if (!obj) { 831 Vec cellgeom; 832 833 ierr = DMPlexComputeGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 834 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject) cellgeom);CHKERRQ(ierr); 835 ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 836 } 837 if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject *) cellgeom);CHKERRQ(ierr);} 838 PetscFunctionReturn(0); 839 } 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "DMPlexSNESGetGeometryFVM" 843 /*@ 844 DMPlexSNESGetGeometryFVM - Return precomputed geometric data 845 846 Input Parameter: 847 . dm - The DM 848 849 Output Parameters: 850 + facegeom - The values precomputed from face geometry 851 . cellgeom - The values precomputed from cell geometry 852 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell 853 854 Level: developer 855 856 .seealso: DMPlexTSSetRHSFunctionLocal() 857 @*/ 858 PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 859 { 860 DMSNES dmsnes; 861 PetscObject obj; 862 PetscErrorCode ierr; 863 864 PetscFunctionBegin; 865 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 866 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 867 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", &obj);CHKERRQ(ierr); 868 if (!obj) { 869 Vec cellgeom, facegeom; 870 871 ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr); 872 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr); 873 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr); 874 ierr = VecDestroy(&facegeom);CHKERRQ(ierr); 875 ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 876 } 877 if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);} 878 if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);} 879 if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);} 880 PetscFunctionReturn(0); 881 } 882 883 #undef __FUNCT__ 884 #define __FUNCT__ "DMPlexSNESGetGradientDM" 885 /*@C 886 DMPlexSNESGetGradientDM - Return gradient data layout 887 888 Input Parameters: 889 + dm - The DM 890 - fv - The PetscFV 891 892 Output Parameter: 893 . dmGrad - The layout for gradient values 894 895 Level: developer 896 897 .seealso: DMPlexSNESGetGeometryFVM() 898 @*/ 899 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad) 900 { 901 DMSNES dmsnes; 902 PetscObject obj; 903 PetscBool computeGradients; 904 PetscErrorCode ierr; 905 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 908 PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2); 909 PetscValidPointer(dmGrad,3); 910 ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr); 911 if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);} 912 ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr); 913 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", &obj);CHKERRQ(ierr); 914 if (!obj) { 915 DM dmGrad; 916 Vec faceGeometry, cellGeometry; 917 918 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 919 ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr); 920 ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr); 921 ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 922 } 923 ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr); 924 PetscFunctionReturn(0); 925 } 926 927 #undef __FUNCT__ 928 #define __FUNCT__ "DMPlexGetCellFields" 929 /*@C 930 DMPlexGetCellFields - Retrieve the field values values for a chunk of cells 931 932 Input Parameters: 933 + dm - The DM 934 . cStart - The first cell to include 935 . cEnd - The first cell to exclude 936 . locX - A local vector with the solution fields 937 . locX_t - A local vector with solution field time derivatives, or NULL 938 - locA - A local vector with auxiliary fields, or NULL 939 940 Output Parameters: 941 + u - The field coefficients 942 . u_t - The fields derivative coefficients 943 - a - The auxiliary field coefficients 944 945 Level: developer 946 947 .seealso: DMPlexGetFaceFields() 948 @*/ 949 PetscErrorCode DMPlexGetCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a) 950 { 951 DM dmAux; 952 PetscSection section, sectionAux; 953 PetscDS prob; 954 PetscInt numCells = cEnd - cStart, totDim, totDimAux, c; 955 PetscErrorCode ierr; 956 957 PetscFunctionBegin; 958 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 959 PetscValidHeaderSpecific(locX, VEC_CLASSID, 4); 960 if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);} 961 if (locA) {PetscValidHeaderSpecific(locA, VEC_CLASSID, 6);} 962 PetscValidPointer(u, 7); 963 PetscValidPointer(u_t, 8); 964 PetscValidPointer(a, 9); 965 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 966 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 967 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 968 if (locA) { 969 PetscDS probAux; 970 971 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 972 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 973 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 974 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 975 } 976 ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u);CHKERRQ(ierr); 977 ierr = DMGetWorkArray(dm, locX_t ? numCells*totDim : 0, PETSC_SCALAR, u_t);CHKERRQ(ierr); 978 ierr = DMGetWorkArray(dm, locA ? numCells*totDimAux : 0, PETSC_SCALAR, a);CHKERRQ(ierr); 979 for (c = cStart; c < cEnd; ++c) { 980 PetscScalar *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a; 981 PetscInt i; 982 983 ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 984 for (i = 0; i < totDim; ++i) ul[c*totDim+i] = x[i]; 985 ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr); 986 if (locX_t) { 987 ierr = DMPlexVecGetClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr); 988 for (i = 0; i < totDim; ++i) ul_t[c*totDim+i] = x_t[i]; 989 ierr = DMPlexVecRestoreClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr); 990 } 991 if (locA) { 992 ierr = DMPlexVecGetClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 993 for (i = 0; i < totDimAux; ++i) al[c*totDimAux+i] = x[i]; 994 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr); 995 } 996 } 997 PetscFunctionReturn(0); 998 } 999 1000 #undef __FUNCT__ 1001 #define __FUNCT__ "DMPlexRestoreCellFields" 1002 /*@C 1003 DMPlexRestoreCellFields - Restore the field values values for a chunk of cells 1004 1005 Input Parameters: 1006 + dm - The DM 1007 . cStart - The first cell to include 1008 . cEnd - The first cell to exclude 1009 . locX - A local vector with the solution fields 1010 . locX_t - A local vector with solution field time derivatives, or NULL 1011 - locA - A local vector with auxiliary fields, or NULL 1012 1013 Output Parameters: 1014 + u - The field coefficients 1015 . u_t - The fields derivative coefficients 1016 - a - The auxiliary field coefficients 1017 1018 Level: developer 1019 1020 .seealso: DMPlexGetFaceFields() 1021 @*/ 1022 PetscErrorCode DMPlexRestoreCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a) 1023 { 1024 PetscErrorCode ierr; 1025 1026 PetscFunctionBegin; 1027 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u);CHKERRQ(ierr); 1028 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u_t);CHKERRQ(ierr); 1029 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, a);CHKERRQ(ierr); 1030 PetscFunctionReturn(0); 1031 } 1032 1033 #undef __FUNCT__ 1034 #define __FUNCT__ "DMPlexGetFaceFields" 1035 /*@C 1036 DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces 1037 1038 Input Parameters: 1039 + dm - The DM 1040 . fStart - The first face to include 1041 . fEnd - The first face to exclude 1042 . locX - A local vector with the solution fields 1043 . locX_t - A local vector with solution field time derivatives, or NULL 1044 . faceGeometry - A local vector with face geometry 1045 . cellGeometry - A local vector with cell geometry 1046 - locaGrad - A local vector with field gradients, or NULL 1047 1048 Output Parameters: 1049 + uL - The field coefficients at the left side of the face 1050 - uR - The field coefficients at the right side of the face 1051 1052 Level: developer 1053 1054 .seealso: DMPlexGetCellFields() 1055 @*/ 1056 PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR) 1057 { 1058 DM dmFace, dmCell, dmGrad = NULL; 1059 PetscDS prob; 1060 DMLabel ghostLabel; 1061 const PetscScalar *facegeom, *cellgeom, *x, *lgrad; 1062 PetscInt dim, totDim, numFaces = fEnd - fStart, iface, face; 1063 PetscErrorCode ierr; 1064 1065 PetscFunctionBegin; 1066 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1067 PetscValidHeaderSpecific(locX, VEC_CLASSID, 4); 1068 if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);} 1069 PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 6); 1070 PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 7); 1071 if (locGrad) {PetscValidHeaderSpecific(locGrad, VEC_CLASSID, 8);} 1072 PetscValidPointer(uL, 9); 1073 PetscValidPointer(uR, 10); 1074 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1075 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1076 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1077 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1078 ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 1079 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1080 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1081 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1082 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1083 if (locGrad) { 1084 ierr = VecGetDM(locGrad, &dmGrad);CHKERRQ(ierr); 1085 ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1086 } 1087 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, uL);CHKERRQ(ierr); 1088 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, uR);CHKERRQ(ierr); 1089 for (face = fStart, iface = 0; face < fEnd; ++face) { 1090 const PetscInt *cells; 1091 const PetscFVFaceGeom *fg; 1092 const PetscFVCellGeom *cgL, *cgR; 1093 const PetscScalar *xL, *xR, *gL, *gR; 1094 PetscScalar *uLl = *uL, *uRl = *uR; 1095 PetscInt ghost, d; 1096 1097 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1098 if (ghost >= 0) continue; 1099 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1100 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1101 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 1102 ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 1103 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr); 1104 ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr); 1105 if (dmGrad) { 1106 PetscReal dxL[3], dxR[3]; 1107 1108 ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); 1109 ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); 1110 DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL); 1111 DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR); 1112 for (d = 0; d < totDim; ++d) { 1113 uLl[iface*totDim+d] = xL[d] + DMPlex_DotD_Internal(dim, &gL[d*dim], dxL); 1114 uRl[iface*totDim+d] = xR[d] + DMPlex_DotD_Internal(dim, &gR[d*dim], dxR); 1115 } 1116 } else { 1117 for (d = 0; d < totDim; ++d) { 1118 uLl[iface*totDim+d] = xL[d]; 1119 uRl[iface*totDim+d] = xR[d]; 1120 } 1121 } 1122 ++iface; 1123 } 1124 ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 1125 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1126 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1127 if (locGrad) { 1128 ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1129 } 1130 PetscFunctionReturn(0); 1131 } 1132 1133 #undef __FUNCT__ 1134 #define __FUNCT__ "DMPlexRestoreFaceFields" 1135 /*@C 1136 DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces 1137 1138 Input Parameters: 1139 + dm - The DM 1140 . fStart - The first face to include 1141 . fEnd - The first face to exclude 1142 . locX - A local vector with the solution fields 1143 . locX_t - A local vector with solution field time derivatives, or NULL 1144 . faceGeometry - A local vector with face geometry 1145 . cellGeometry - A local vector with cell geometry 1146 - locaGrad - A local vector with field gradients, or NULL 1147 1148 Output Parameters: 1149 + uL - The field coefficients at the left side of the face 1150 - uR - The field coefficients at the right side of the face 1151 1152 Level: developer 1153 1154 .seealso: DMPlexGetFaceFields() 1155 @*/ 1156 PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR) 1157 { 1158 PetscErrorCode ierr; 1159 1160 PetscFunctionBegin; 1161 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uL);CHKERRQ(ierr); 1162 ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uR);CHKERRQ(ierr); 1163 PetscFunctionReturn(0); 1164 } 1165 1166 #undef __FUNCT__ 1167 #define __FUNCT__ "DMPlexGetFaceGeometry" 1168 /*@C 1169 DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces 1170 1171 Input Parameters: 1172 + dm - The DM 1173 . fStart - The first face to include 1174 . fEnd - The first face to exclude 1175 . faceGeometry - A local vector with face geometry 1176 - cellGeometry - A local vector with cell geometry 1177 1178 Output Parameters: 1179 + fgeom - The extract the face centroid and normal 1180 - vol - The cell volume 1181 1182 Level: developer 1183 1184 .seealso: DMPlexGetCellFields() 1185 @*/ 1186 PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol) 1187 { 1188 DM dmFace, dmCell; 1189 DMLabel ghostLabel; 1190 const PetscScalar *facegeom, *cellgeom; 1191 PetscInt dim, numFaces = fEnd - fStart, iface, face; 1192 PetscErrorCode ierr; 1193 1194 PetscFunctionBegin; 1195 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1196 PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 4); 1197 PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 5); 1198 PetscValidPointer(fgeom, 6); 1199 PetscValidPointer(vol, 7); 1200 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1201 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1202 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1203 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1204 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1205 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1206 ierr = PetscMalloc1(numFaces, fgeom);CHKERRQ(ierr); 1207 ierr = DMGetWorkArray(dm, numFaces*2, PETSC_SCALAR, vol);CHKERRQ(ierr); 1208 for (face = fStart, iface = 0; face < fEnd; ++face) { 1209 const PetscInt *cells; 1210 const PetscFVFaceGeom *fg; 1211 const PetscFVCellGeom *cgL, *cgR; 1212 PetscFVFaceGeom *fgeoml = *fgeom; 1213 PetscReal *voll = *vol; 1214 PetscInt ghost, d; 1215 1216 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1217 if (ghost >= 0) continue; 1218 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1219 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1220 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 1221 ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 1222 for (d = 0; d < dim; ++d) { 1223 fgeoml[iface].centroid[d] = fg->centroid[d]; 1224 fgeoml[iface].normal[d] = fg->normal[d]; 1225 } 1226 voll[iface*2+0] = cgL->volume; 1227 voll[iface*2+1] = cgR->volume; 1228 ++iface; 1229 } 1230 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1231 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1232 PetscFunctionReturn(0); 1233 } 1234 1235 #undef __FUNCT__ 1236 #define __FUNCT__ "DMPlexRestoreFaceGeometry" 1237 /*@C 1238 DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces 1239 1240 Input Parameters: 1241 + dm - The DM 1242 . fStart - The first face to include 1243 . fEnd - The first face to exclude 1244 . faceGeometry - A local vector with face geometry 1245 - cellGeometry - A local vector with cell geometry 1246 1247 Output Parameters: 1248 + fgeom - The extract the face centroid and normal 1249 - vol - The cell volume 1250 1251 Level: developer 1252 1253 .seealso: DMPlexGetFaceFields() 1254 @*/ 1255 PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol) 1256 { 1257 PetscErrorCode ierr; 1258 1259 PetscFunctionBegin; 1260 ierr = PetscFree(*fgeom);CHKERRQ(ierr); 1261 ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, vol);CHKERRQ(ierr); 1262 PetscFunctionReturn(0); 1263 } 1264 1265 #undef __FUNCT__ 1266 #define __FUNCT__ "DMPlexReconstructGradients_Internal" 1267 PetscErrorCode DMPlexReconstructGradients_Internal(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, Vec locX, Vec grad) 1268 { 1269 DM dmFace, dmCell, dmGrad; 1270 DMLabel ghostLabel; 1271 PetscDS prob; 1272 PetscFV fvm; 1273 PetscLimiter lim; 1274 const PetscScalar *facegeom, *cellgeom, *x; 1275 PetscScalar *gr; 1276 PetscReal *cellPhi; 1277 PetscInt dim, face, cell, totDim, cStart, cEnd, cEndInterior; 1278 PetscErrorCode ierr; 1279 1280 PetscFunctionBegin; 1281 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1282 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1283 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1284 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1285 ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 1286 ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); 1287 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1288 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1289 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1290 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1291 ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 1292 ierr = VecGetDM(grad, &dmGrad);CHKERRQ(ierr); 1293 ierr = VecZeroEntries(grad);CHKERRQ(ierr); 1294 ierr = VecGetArray(grad, &gr);CHKERRQ(ierr); 1295 /* Reconstruct gradients */ 1296 for (face = fStart; face < fEnd; ++face) { 1297 const PetscInt *cells; 1298 const PetscFVFaceGeom *fg; 1299 const PetscScalar *cx[2]; 1300 PetscScalar *cgrad[2]; 1301 PetscBool boundary; 1302 PetscInt ghost, c, pd, d; 1303 1304 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1305 if (ghost >= 0) continue; 1306 ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 1307 if (boundary) continue; 1308 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1309 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 1310 for (c = 0; c < 2; ++c) { 1311 ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); 1312 ierr = DMPlexPointGlobalRef(dmGrad, cells[c], gr, &cgrad[c]);CHKERRQ(ierr); 1313 } 1314 for (pd = 0; pd < totDim; ++pd) { 1315 PetscScalar delta = cx[1][pd] - cx[0][pd]; 1316 1317 for (d = 0; d < dim; ++d) { 1318 if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; 1319 if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; 1320 } 1321 } 1322 } 1323 /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ 1324 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1325 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1326 cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior; 1327 ierr = DMGetWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 1328 for (cell = dmGrad && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { 1329 const PetscInt *faces; 1330 const PetscScalar *cx; 1331 const PetscFVCellGeom *cg; 1332 PetscScalar *cgrad; 1333 PetscInt coneSize, f, pd, d; 1334 1335 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1336 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1337 ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); 1338 ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); 1339 ierr = DMPlexPointGlobalRef(dmGrad, cell, gr, &cgrad);CHKERRQ(ierr); 1340 if (!cgrad) continue; /* Unowned overlap cell, we do not compute */ 1341 /* Limiter will be minimum value over all neighbors */ 1342 for (d = 0; d < totDim; ++d) cellPhi[d] = PETSC_MAX_REAL; 1343 for (f = 0; f < coneSize; ++f) { 1344 const PetscScalar *ncx; 1345 const PetscFVCellGeom *ncg; 1346 const PetscInt *fcells; 1347 PetscInt face = faces[f], ncell, ghost; 1348 PetscReal v[3]; 1349 PetscBool boundary; 1350 1351 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1352 ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 1353 if ((ghost >= 0) || boundary) continue; 1354 ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); 1355 ncell = cell == fcells[0] ? fcells[1] : fcells[0]; 1356 ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); 1357 ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); 1358 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v); 1359 for (d = 0; d < totDim; ++d) { 1360 /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ 1361 PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DMPlex_DotD_Internal(dim, &cgrad[d*dim], v); 1362 1363 ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); 1364 cellPhi[d] = PetscMin(cellPhi[d], phi); 1365 } 1366 } 1367 /* Apply limiter to gradient */ 1368 for (pd = 0; pd < totDim; ++pd) 1369 /* Scalar limiter applied to each component separately */ 1370 for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; 1371 } 1372 ierr = DMRestoreWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 1373 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 1374 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 1375 ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 1376 ierr = VecRestoreArray(grad, &gr);CHKERRQ(ierr); 1377 PetscFunctionReturn(0); 1378 } 1379 1380 #undef __FUNCT__ 1381 #define __FUNCT__ "DMPlexComputeBdResidual_Internal" 1382 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, Vec locF, void *user) 1383 { 1384 DM_Plex *mesh = (DM_Plex *) dm->data; 1385 PetscSection section; 1386 PetscDS prob; 1387 DMLabel depth; 1388 PetscFECellGeom *cgeom; 1389 PetscScalar *u = NULL, *u_t = NULL, *elemVec = NULL; 1390 PetscInt dim, Nf, f, totDimBd, numBd, bd; 1391 PetscErrorCode ierr; 1392 1393 PetscFunctionBegin; 1394 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1395 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1396 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1397 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1398 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 1399 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1400 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1401 for (bd = 0; bd < numBd; ++bd) { 1402 const char *bdLabel; 1403 DMLabel label; 1404 IS pointIS; 1405 const PetscInt *points; 1406 const PetscInt *values; 1407 PetscInt field, numValues, numPoints, p, dep, numFaces; 1408 PetscBool isEssential; 1409 1410 ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1411 if (isEssential) continue; 1412 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 1413 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1414 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1415 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1416 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1417 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1418 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1419 if (dep == dim-1) ++numFaces; 1420 } 1421 ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd,&elemVec);CHKERRQ(ierr); 1422 if (locX_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 1423 for (p = 0, f = 0; p < numPoints; ++p) { 1424 const PetscInt point = points[p]; 1425 PetscScalar *x = NULL; 1426 PetscInt i; 1427 1428 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1429 if (dep != dim-1) continue; 1430 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr); 1431 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n); 1432 if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point); 1433 ierr = DMPlexVecGetClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr); 1434 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 1435 ierr = DMPlexVecRestoreClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr); 1436 if (locX_t) { 1437 ierr = DMPlexVecGetClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr); 1438 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 1439 ierr = DMPlexVecRestoreClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr); 1440 } 1441 ++f; 1442 } 1443 for (f = 0; f < Nf; ++f) { 1444 PetscFE fe; 1445 PetscQuadrature q; 1446 PetscInt numQuadPoints, Nb; 1447 /* Conforming batches */ 1448 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1449 /* Remainder */ 1450 PetscInt Nr, offset; 1451 1452 ierr = PetscDSGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1453 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1454 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1455 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1456 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1457 blockSize = Nb*numQuadPoints; 1458 batchSize = numBlocks * blockSize; 1459 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1460 numChunks = numFaces / (numBatches*batchSize); 1461 Ne = numChunks*numBatches*batchSize; 1462 Nr = numFaces % (numBatches*batchSize); 1463 offset = numFaces - Nr; 1464 ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, cgeom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr); 1465 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); 1466 } 1467 for (p = 0, f = 0; p < numPoints; ++p) { 1468 const PetscInt point = points[p]; 1469 1470 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 1471 if (dep != dim-1) continue; 1472 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);} 1473 ierr = DMPlexVecSetClosure(dm, NULL, locF, point, &elemVec[f*totDimBd], ADD_VALUES);CHKERRQ(ierr); 1474 ++f; 1475 } 1476 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 1477 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1478 ierr = PetscFree3(u,cgeom,elemVec);CHKERRQ(ierr); 1479 if (locX_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 1480 } 1481 PetscFunctionReturn(0); 1482 } 1483 1484 #undef __FUNCT__ 1485 #define __FUNCT__ "DMPlexComputeResidual_Internal" 1486 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) 1487 { 1488 DM_Plex *mesh = (DM_Plex *) dm->data; 1489 const char *name = "Residual"; 1490 DM dmAux = NULL; 1491 DM dmGrad = NULL; 1492 DMLabel ghostLabel = NULL; 1493 PetscDS prob = NULL; 1494 PetscDS probAux = NULL; 1495 PetscSection section = NULL; 1496 PetscBool useFEM = PETSC_FALSE; 1497 PetscBool useFVM = PETSC_FALSE; 1498 PetscFV fvm = NULL; 1499 PetscFECellGeom *cgeomFEM = NULL; 1500 PetscFVCellGeom *cgeomFVM = NULL; 1501 PetscFVFaceGeom *fgeomFVM = NULL; 1502 Vec locA, cellGeometryFEM = NULL, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL; 1503 PetscScalar *u, *u_t, *a, *uL, *uR; 1504 PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, cStart, cEnd, fStart, fEnd; 1505 PetscErrorCode ierr; 1506 1507 PetscFunctionBegin; 1508 ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1509 /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */ 1510 /* TODO The FVM geometry is over-manipulated. Make the preclac functions return exactly what we need */ 1511 /* FEM+FVM */ 1512 /* 1: Get sizes from dm and dmAux */ 1513 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1514 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1515 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1516 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1517 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1518 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr); 1519 if (locA) { 1520 ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr); 1521 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1522 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1523 } 1524 /* 2: Get geometric data */ 1525 for (f = 0; f < Nf; ++f) { 1526 PetscObject obj; 1527 PetscClassId id; 1528 1529 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1530 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1531 if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;} 1532 if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;} 1533 } 1534 if (useFEM) { 1535 ierr = DMPlexSNESGetGeometryFEM(dm, &cellGeometryFEM);CHKERRQ(ierr); 1536 ierr = VecGetArrayRead(cellGeometryFEM, (const PetscScalar **) &cgeomFEM);CHKERRQ(ierr); 1537 } 1538 if (useFVM) { 1539 ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr); 1540 ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 1541 ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1542 /* Reconstruct and limit cell gradients */ 1543 ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr); 1544 if (dmGrad) { 1545 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1546 ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1547 ierr = DMPlexReconstructGradients_Internal(dm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr); 1548 /* Communicate gradient values */ 1549 ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1550 ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1551 ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1552 ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 1553 } 1554 } 1555 /* Handle boundary values */ 1556 ierr = DMPlexInsertBoundaryValues(dm, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr); 1557 /* Loop over chunks */ 1558 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1559 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1560 numChunks = 1; 1561 cellChunkSize = (cEnd - cStart)/numChunks; 1562 faceChunkSize = (fEnd - fStart)/numChunks; 1563 for (chunk = 0; chunk < numChunks; ++chunk) { 1564 PetscScalar *elemVec, *fluxL, *fluxR; 1565 PetscReal *vol; 1566 PetscFVFaceGeom *fgeom; 1567 PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, cell; 1568 PetscInt fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = fE - fS, face; 1569 1570 /* Extract field coefficients */ 1571 if (useFEM) { 1572 ierr = DMPlexGetCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1573 ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr); 1574 } 1575 if (useFVM) { 1576 ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr); 1577 ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr); 1578 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr); 1579 ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr); 1580 } 1581 /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */ 1582 /* Loop over fields */ 1583 for (f = 0; f < Nf; ++f) { 1584 PetscObject obj; 1585 PetscClassId id; 1586 PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 1587 1588 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 1589 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 1590 if (id == PETSCFE_CLASSID) { 1591 PetscFE fe = (PetscFE) obj; 1592 PetscQuadrature q; 1593 PetscInt Nq, Nb; 1594 1595 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1596 1597 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1598 ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1599 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1600 blockSize = Nb*Nq; 1601 batchSize = numBlocks * blockSize; 1602 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1603 numChunks = numCells / (numBatches*batchSize); 1604 Ne = numChunks*numBatches*batchSize; 1605 Nr = numCells % (numBatches*batchSize); 1606 offset = numCells - Nr; 1607 /* Integrate FE residual to get elemVec (need fields at quadrature points) */ 1608 /* 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) */ 1609 ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeomFEM, u, u_t, probAux, a, elemVec);CHKERRQ(ierr); 1610 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); 1611 } else if (id == PETSCFV_CLASSID) { 1612 PetscFV fv = (PetscFV) obj; 1613 1614 Ne = numFaces; 1615 Nr = 0; 1616 /* Riemann solve over faces (need fields at face centroids) */ 1617 /* We need to evaluate FE fields at those coordinates */ 1618 ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr); 1619 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1620 } 1621 /* Loop over domain */ 1622 if (useFEM) { 1623 /* Add elemVec to locX */ 1624 for (cell = cS; cell < cE; ++cell) { 1625 if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cell*totDim]);CHKERRQ(ierr);} 1626 ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cell*totDim], ADD_VALUES);CHKERRQ(ierr); 1627 } 1628 } 1629 if (useFVM) { 1630 PetscScalar *fa; 1631 PetscInt iface; 1632 1633 /* Accumulate fluxes to cells */ 1634 ierr = VecGetArray(locF, &fa);CHKERRQ(ierr); 1635 for (face = fS, iface = 0; face < fE; ++face) { 1636 const PetscInt *cells; 1637 PetscScalar *fL, *fR; 1638 PetscInt ghost, d; 1639 1640 ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 1641 if (ghost >= 0) continue; 1642 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 1643 ierr = DMPlexPointGlobalRef(dm, cells[0], fa, &fL);CHKERRQ(ierr); 1644 ierr = DMPlexPointGlobalRef(dm, cells[1], fa, &fR);CHKERRQ(ierr); 1645 for (d = 0; d < totDim; ++d) { 1646 if (fL) fL[d] -= fluxL[iface*totDim+d]; 1647 if (fR) fR[d] += fluxR[iface*totDim+d]; 1648 } 1649 ++iface; 1650 } 1651 ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr); 1652 } 1653 if (useFEM) { 1654 ierr = DMPlexRestoreCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 1655 ierr = DMRestoreWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr); 1656 } 1657 if (useFVM) { 1658 ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr); 1659 ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr); 1660 ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr); 1661 ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr); 1662 if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);} 1663 } 1664 } 1665 1666 if (useFEM) {ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, locF, user);CHKERRQ(ierr);} 1667 1668 /* FEM */ 1669 /* 1: Get sizes from dm and dmAux */ 1670 /* 2: Get geometric data */ 1671 /* 3: Handle boundary values */ 1672 /* 4: Loop over domain */ 1673 /* Extract coefficients */ 1674 /* Loop over fields */ 1675 /* Set tiling for FE*/ 1676 /* Integrate FE residual to get elemVec */ 1677 /* Loop over subdomain */ 1678 /* Loop over quad points */ 1679 /* Transform coords to real space */ 1680 /* Evaluate field and aux fields at point */ 1681 /* Evaluate residual at point */ 1682 /* Transform residual to real space */ 1683 /* Add residual to elemVec */ 1684 /* Loop over domain */ 1685 /* Add elemVec to locX */ 1686 1687 /* FVM */ 1688 /* Get geometric data */ 1689 /* If using gradients */ 1690 /* Compute gradient data */ 1691 /* Loop over domain faces */ 1692 /* Count computational faces */ 1693 /* Reconstruct cell gradient */ 1694 /* Loop over domain cells */ 1695 /* Limit cell gradients */ 1696 /* Handle boundary values */ 1697 /* Loop over domain faces */ 1698 /* Read out field, centroid, normal, volume for each side of face */ 1699 /* Riemann solve over faces */ 1700 /* Loop over domain faces */ 1701 /* Accumulate fluxes to cells */ 1702 /* TODO Change printFEM to printDisc here */ 1703 if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, locF);CHKERRQ(ierr);} 1704 ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 1705 PetscFunctionReturn(0); 1706 } 1707 1708 #undef __FUNCT__ 1709 #define __FUNCT__ "DMPlexComputeResidualFEM_Check_Internal" 1710 static PetscErrorCode DMPlexComputeResidualFEM_Check_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user) 1711 { 1712 DM dmCh, dmAux; 1713 Vec A, cellgeom; 1714 PetscDS prob, probCh, probAux = NULL; 1715 PetscQuadrature q; 1716 PetscSection section, sectionAux; 1717 PetscFECellGeom *cgeom; 1718 PetscScalar *elemVec, *elemVecCh, *u, *u_t, *a = NULL; 1719 PetscInt dim, Nf, f, numCells, cStart, cEnd, c; 1720 PetscInt totDim, totDimAux, diffCell = 0; 1721 PetscErrorCode ierr; 1722 1723 PetscFunctionBegin; 1724 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1725 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1726 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1727 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1728 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1729 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1730 numCells = cEnd - cStart; 1731 ierr = PetscObjectQuery((PetscObject) dm, "dmCh", (PetscObject *) &dmCh);CHKERRQ(ierr); 1732 ierr = DMGetDS(dmCh, &probCh);CHKERRQ(ierr); 1733 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1734 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1735 if (dmAux) { 1736 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1737 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1738 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1739 } 1740 ierr = DMPlexInsertBoundaryValues(dm, X, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1741 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 1742 ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim,&elemVec);CHKERRQ(ierr); 1743 ierr = PetscMalloc1(numCells*totDim,&elemVecCh);CHKERRQ(ierr); 1744 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1745 ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 1746 ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr); 1747 for (c = cStart; c < cEnd; ++c) { 1748 PetscScalar *x = NULL, *x_t = NULL; 1749 PetscInt i; 1750 1751 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1752 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1753 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1754 if (X_t) { 1755 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1756 for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i]; 1757 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1758 } 1759 if (dmAux) { 1760 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1761 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1762 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1763 } 1764 } 1765 for (f = 0; f < Nf; ++f) { 1766 PetscFE fe, feCh; 1767 PetscInt numQuadPoints, Nb; 1768 /* Conforming batches */ 1769 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1770 /* Remainder */ 1771 PetscInt Nr, offset; 1772 1773 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1774 ierr = PetscDSGetDiscretization(probCh, f, (PetscObject *) &feCh);CHKERRQ(ierr); 1775 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1776 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1777 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1778 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1779 blockSize = Nb*numQuadPoints; 1780 batchSize = numBlocks * blockSize; 1781 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1782 numChunks = numCells / (numBatches*batchSize); 1783 Ne = numChunks*numBatches*batchSize; 1784 Nr = numCells % (numBatches*batchSize); 1785 offset = numCells - Nr; 1786 ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr); 1787 ierr = PetscFEIntegrateResidual(feCh, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVecCh);CHKERRQ(ierr); 1788 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); 1789 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); 1790 } 1791 for (c = cStart; c < cEnd; ++c) { 1792 PetscBool diff = PETSC_FALSE; 1793 PetscInt d; 1794 1795 for (d = 0; d < totDim; ++d) if (PetscAbsScalar(elemVec[c*totDim+d] - elemVecCh[c*totDim+d]) > 1.0e-7) {diff = PETSC_TRUE;break;} 1796 if (diff) { 1797 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Different cell %d\n", c);CHKERRQ(ierr); 1798 ierr = DMPrintCellVector(c, "Residual", totDim, &elemVec[c*totDim]);CHKERRQ(ierr); 1799 ierr = DMPrintCellVector(c, "Check Residual", totDim, &elemVecCh[c*totDim]);CHKERRQ(ierr); 1800 ++diffCell; 1801 } 1802 if (diffCell > 9) break; 1803 ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_VALUES);CHKERRQ(ierr); 1804 } 1805 ierr = VecRestoreArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr); 1806 ierr = PetscFree3(u,u_t,elemVec);CHKERRQ(ierr); 1807 ierr = PetscFree(elemVecCh);CHKERRQ(ierr); 1808 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1809 PetscFunctionReturn(0); 1810 } 1811 1812 #undef __FUNCT__ 1813 #define __FUNCT__ "DMPlexSNESComputeResidualFEM" 1814 /*@ 1815 DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 1816 1817 Input Parameters: 1818 + dm - The mesh 1819 . X - Local solution 1820 - user - The user context 1821 1822 Output Parameter: 1823 . F - Local output vector 1824 1825 Level: developer 1826 1827 .seealso: DMPlexComputeJacobianActionFEM() 1828 @*/ 1829 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 1830 { 1831 PetscObject check; 1832 PetscErrorCode ierr; 1833 1834 PetscFunctionBegin; 1835 /* The dmCh is used to check two mathematically equivalent discretizations for computational equivalence */ 1836 ierr = PetscObjectQuery((PetscObject) dm, "dmCh", &check);CHKERRQ(ierr); 1837 if (check) {ierr = DMPlexComputeResidualFEM_Check_Internal(dm, X, NULL, F, user);CHKERRQ(ierr);} 1838 else {ierr = DMPlexComputeResidual_Internal(dm, 0.0, X, NULL, F, user);CHKERRQ(ierr);} 1839 PetscFunctionReturn(0); 1840 } 1841 1842 #undef __FUNCT__ 1843 #define __FUNCT__ "DMPlexComputeJacobianFEM_Internal" 1844 PetscErrorCode DMPlexComputeJacobianFEM_Internal(DM dm, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user) 1845 { 1846 DM_Plex *mesh = (DM_Plex *) dm->data; 1847 const char *name = "Jacobian"; 1848 DM dmAux; 1849 DMLabel depth; 1850 Vec A, cellgeom; 1851 PetscDS prob, probAux = NULL; 1852 PetscQuadrature quad; 1853 PetscSection section, globalSection, sectionAux; 1854 PetscFECellGeom *cgeom; 1855 PetscScalar *elemMat, *u, *u_t, *a = NULL; 1856 PetscInt dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 1857 PetscInt totDim, totDimBd, totDimAux, numBd, bd; 1858 PetscBool isShell; 1859 PetscErrorCode ierr; 1860 1861 PetscFunctionBegin; 1862 ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 1863 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1864 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1865 ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 1866 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1867 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1868 ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr); 1869 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 1870 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1871 numCells = cEnd - cStart; 1872 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 1873 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 1874 if (dmAux) { 1875 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 1876 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 1877 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 1878 } 1879 ierr = DMPlexInsertBoundaryValues(dm, X, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 1880 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 1881 ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat);CHKERRQ(ierr); 1882 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 1883 ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr); 1884 ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr); 1885 for (c = cStart; c < cEnd; ++c) { 1886 PetscScalar *x = NULL, *x_t = NULL; 1887 PetscInt i; 1888 1889 ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1890 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 1891 ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 1892 if (X_t) { 1893 ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1894 for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i]; 1895 ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr); 1896 } 1897 if (dmAux) { 1898 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1899 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 1900 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 1901 } 1902 } 1903 ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr); 1904 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1905 PetscFE fe; 1906 PetscInt numQuadPoints, Nb; 1907 /* Conforming batches */ 1908 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1909 /* Remainder */ 1910 PetscInt Nr, offset; 1911 1912 ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1913 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1914 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1915 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1916 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1917 blockSize = Nb*numQuadPoints; 1918 batchSize = numBlocks * blockSize; 1919 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1920 numChunks = numCells / (numBatches*batchSize); 1921 Ne = numChunks*numBatches*batchSize; 1922 Nr = numCells % (numBatches*batchSize); 1923 offset = numCells - Nr; 1924 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 1925 ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr); 1926 ierr = PetscFEIntegrateJacobian(fe, prob, 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); 1927 } 1928 } 1929 for (c = cStart; c < cEnd; ++c) { 1930 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[c*totDim*totDim]);CHKERRQ(ierr);} 1931 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*totDim*totDim], ADD_VALUES);CHKERRQ(ierr); 1932 } 1933 ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr); 1934 ierr = PetscFree3(u,u_t,elemMat);CHKERRQ(ierr); 1935 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 1936 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1937 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1938 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1939 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 1940 for (bd = 0; bd < numBd; ++bd) { 1941 const char *bdLabel; 1942 DMLabel label; 1943 IS pointIS; 1944 const PetscInt *points; 1945 const PetscInt *values; 1946 PetscInt field, numValues, numPoints, p, dep, numFaces; 1947 PetscBool isEssential; 1948 1949 ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr); 1950 if (isEssential) continue; 1951 if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this"); 1952 ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr); 1953 ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 1954 ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 1955 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 1956 for (p = 0, numFaces = 0; p < numPoints; ++p) { 1957 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1958 if (dep == dim-1) ++numFaces; 1959 } 1960 ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr); 1961 if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);} 1962 for (p = 0, f = 0; p < numPoints; ++p) { 1963 const PetscInt point = points[p]; 1964 PetscScalar *x = NULL; 1965 PetscInt i; 1966 1967 ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 1968 if (dep != dim-1) continue; 1969 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr); 1970 ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n); 1971 if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point); 1972 ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1973 for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i]; 1974 ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 1975 if (X_t) { 1976 ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1977 for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i]; 1978 ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr); 1979 } 1980 ++f; 1981 } 1982 ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr); 1983 for (fieldI = 0; fieldI < Nf; ++fieldI) { 1984 PetscFE fe; 1985 PetscInt numQuadPoints, Nb; 1986 /* Conforming batches */ 1987 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1988 /* Remainder */ 1989 PetscInt Nr, offset; 1990 1991 ierr = PetscDSGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr); 1992 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1993 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1994 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1995 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1996 blockSize = Nb*numQuadPoints; 1997 batchSize = numBlocks * blockSize; 1998 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1999 numChunks = numFaces / (numBatches*batchSize); 2000 Ne = numChunks*numBatches*batchSize; 2001 Nr = numFaces % (numBatches*batchSize); 2002 offset = numFaces - Nr; 2003 for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 2004 ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr); 2005 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); 2006 } 2007 } 2008 for (p = 0, f = 0; p < numPoints; ++p) { 2009 const PetscInt point = points[p]; 2010 2011 ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 2012 if (dep != dim-1) continue; 2013 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);} 2014 ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr); 2015 ++f; 2016 } 2017 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 2018 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 2019 ierr = PetscFree3(u,cgeom,elemMat);CHKERRQ(ierr); 2020 if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);} 2021 } 2022 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2023 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2024 if (mesh->printFEM) { 2025 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 2026 ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 2027 ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 2028 } 2029 ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 2030 ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 2031 if (isShell) { 2032 JacActionCtx *jctx; 2033 2034 ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 2035 ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 2036 } 2037 PetscFunctionReturn(0); 2038 } 2039 2040 #undef __FUNCT__ 2041 #define __FUNCT__ "DMPlexSNESComputeJacobianFEM" 2042 /*@ 2043 DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 2044 2045 Input Parameters: 2046 + dm - The mesh 2047 . X - Local input vector 2048 - user - The user context 2049 2050 Output Parameter: 2051 . Jac - Jacobian matrix 2052 2053 Note: 2054 The first member of the user context must be an FEMContext. 2055 2056 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 2057 like a GPU, or vectorize on a multicore machine. 2058 2059 Level: developer 2060 2061 .seealso: FormFunctionLocal() 2062 @*/ 2063 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 2064 { 2065 PetscErrorCode ierr; 2066 2067 PetscFunctionBegin; 2068 ierr = DMPlexComputeJacobianFEM_Internal(dm, X, NULL, Jac, JacP, user);CHKERRQ(ierr); 2069 PetscFunctionReturn(0); 2070 } 2071