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