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