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