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