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