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