1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests recovery from domain errors in MatMult() and PCApply()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown See src/ksp/ksp/tutorials/ex19.c from which this was copied 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscsnes.h> 9c4762a1bSJed Brown #include <petscdm.h> 10c4762a1bSJed Brown #include <petscdmda.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown User-defined routines and data structures 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown PetscScalar u, v, omega, temp; 17c4762a1bSJed Brown } Field; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field **, Field **, void *); 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscReal lidvelocity, prandtl, grashof; /* physical parameters */ 23c4762a1bSJed Brown PetscBool draw_contours; /* flag - 1 indicates drawing contours */ 24c4762a1bSJed Brown PetscBool errorindomain; 25c4762a1bSJed Brown PetscBool errorindomainmf; 26c4762a1bSJed Brown SNES snes; 27c4762a1bSJed Brown } AppCtx; 28c4762a1bSJed Brown 29c4762a1bSJed Brown typedef struct { 30c4762a1bSJed Brown Mat Jmf; 31c4762a1bSJed Brown } MatShellCtx; 32c4762a1bSJed Brown 33c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, DM, Vec); 34c4762a1bSJed Brown extern PetscErrorCode MatMult_MyShell(Mat, Vec, Vec); 35c4762a1bSJed Brown extern PetscErrorCode MatAssemblyEnd_MyShell(Mat, MatAssemblyType); 36c4762a1bSJed Brown extern PetscErrorCode PCApply_MyShell(PC, Vec, Vec); 37c4762a1bSJed Brown extern PetscErrorCode SNESComputeJacobian_MyShell(SNES, Vec, Mat, Mat, void *); 38c4762a1bSJed Brown 39d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 40d71ae5a4SJacob Faibussowitsch { 41c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 42c4762a1bSJed Brown PetscInt mx, my; 43c4762a1bSJed Brown MPI_Comm comm; 44c4762a1bSJed Brown DM da; 45c4762a1bSJed Brown Vec x; 46c4762a1bSJed Brown Mat J = NULL, Jmf = NULL; 47c4762a1bSJed Brown MatShellCtx matshellctx; 48c4762a1bSJed Brown PetscInt mlocal, nlocal; 49c4762a1bSJed Brown PC pc; 50c4762a1bSJed Brown KSP ksp; 51c4762a1bSJed Brown PetscBool errorinmatmult = PETSC_FALSE, errorinpcapply = PETSC_FALSE, errorinpcsetup = PETSC_FALSE; 52c4762a1bSJed Brown 53327415f7SBarry Smith PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_matmult", &errorinmatmult, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_pcapply", &errorinpcapply, NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_pcsetup", &errorinpcsetup, NULL)); 58c4762a1bSJed Brown user.errorindomain = PETSC_FALSE; 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_domain", &user.errorindomain, NULL)); 60c4762a1bSJed Brown user.errorindomainmf = PETSC_FALSE; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_domainmf", &user.errorindomainmf, NULL)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 649566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &user.snes)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown Create distributed array object to manage parallel grid and vectors 68c4762a1bSJed Brown for principal unknowns (x) and governing residuals (f) 69c4762a1bSJed Brown */ 709566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 4, 1, 0, 0, &da)); 719566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 729566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 739566063dSJacob Faibussowitsch PetscCall(SNESSetDM(user.snes, da)); 74c4762a1bSJed Brown 759371c9d4SSatish Balay PetscCall(DMDAGetInfo(da, 0, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Problem parameters (velocity of lid, prandtl, and grashof numbers) 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown user.lidvelocity = 1.0 / (mx * my); 80c4762a1bSJed Brown user.prandtl = 1.0; 81c4762a1bSJed Brown user.grashof = 1.0; 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-lidvelocity", &user.lidvelocity, NULL)); 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-prandtl", &user.prandtl, NULL)); 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-grashof", &user.grashof, NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-contours", &user.draw_contours)); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "x_velocity")); 899566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "y_velocity")); 909566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 2, "Omega")); 919566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 3, "temperature")); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 94c4762a1bSJed Brown Create user context, set problem data, create vector data structures. 95c4762a1bSJed Brown Also, compute the initial guess. 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Create nonlinear solver context 100c4762a1bSJed Brown 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1029566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user)); 1039566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown if (errorinmatmult) { 1069566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(user.snes, &Jmf)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jmf)); 1089566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Jmf, &mlocal, &nlocal)); 109c4762a1bSJed Brown matshellctx.Jmf = Jmf; 1109566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)Jmf), mlocal, nlocal, PETSC_DECIDE, PETSC_DECIDE, &matshellctx, &J)); 1119566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_MyShell)); 1129566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(J, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MyShell)); 1139566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(user.snes, J, J, MatMFFDComputeJacobian, NULL)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(user.snes)); 1179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "lid velocity = %g, prandtl # = %g, grashof # = %g\n", (double)user.lidvelocity, (double)user.prandtl, (double)user.grashof)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown if (errorinpcapply) { 1209566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(user.snes, &ksp)); 1219566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1229566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCSHELL)); 1239566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(pc, PCApply_MyShell)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Solve the nonlinear system 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 1309566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, da, x)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown if (errorinpcsetup) { 1339566063dSJacob Faibussowitsch PetscCall(SNESSetUp(user.snes)); 1349566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(user.snes, NULL, NULL, SNESComputeJacobian_MyShell, NULL)); 135c4762a1bSJed Brown } 1369566063dSJacob Faibussowitsch PetscCall(SNESSolve(user.snes, NULL, x)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 140c4762a1bSJed Brown are no longer needed. 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jmf)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1469566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&user.snes)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 148b122ec5aSJacob Faibussowitsch return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 153c4762a1bSJed Brown 154c4762a1bSJed Brown Input Parameters: 155c4762a1bSJed Brown user - user-defined application context 156c4762a1bSJed Brown X - vector 157c4762a1bSJed Brown 158c4762a1bSJed Brown Output Parameter: 159c4762a1bSJed Brown X - vector 160c4762a1bSJed Brown */ 161d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, DM da, Vec X) 162d71ae5a4SJacob Faibussowitsch { 163c4762a1bSJed Brown PetscInt i, j, mx, xs, ys, xm, ym; 164c4762a1bSJed Brown PetscReal grashof, dx; 165c4762a1bSJed Brown Field **x; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBeginUser; 168c4762a1bSJed Brown grashof = user->grashof; 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 171c4762a1bSJed Brown dx = 1.0 / (mx - 1); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA): 175c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 176c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 177c4762a1bSJed Brown */ 1789566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Get a pointer to vector data. 182c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 183c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 184c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 185c4762a1bSJed Brown the array. 186c4762a1bSJed Brown */ 1879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 191c4762a1bSJed Brown Initial condition is motionless fluid and equilibrium temperature 192c4762a1bSJed Brown */ 193c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 194c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 195c4762a1bSJed Brown x[j][i].u = 0.0; 196c4762a1bSJed Brown x[j][i].v = 0.0; 197c4762a1bSJed Brown x[j][i].omega = 0.0; 198c4762a1bSJed Brown x[j][i].temp = (grashof > 0) * i * dx; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Restore vector 204c4762a1bSJed Brown */ 2059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x)); 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr) 210d71ae5a4SJacob Faibussowitsch { 211c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 212c4762a1bSJed Brown PetscInt xints, xinte, yints, yinte, i, j; 213c4762a1bSJed Brown PetscReal hx, hy, dhx, dhy, hxdhy, hydhx; 214c4762a1bSJed Brown PetscReal grashof, prandtl, lid; 215c4762a1bSJed Brown PetscScalar u, uxx, uyy, vx, vy, avx, avy, vxp, vxm, vyp, vym; 216c4762a1bSJed Brown static PetscInt fail = 0; 217c4762a1bSJed Brown 218c4762a1bSJed Brown PetscFunctionBeginUser; 219c4762a1bSJed Brown if ((fail++ > 7 && user->errorindomainmf) || (fail++ > 36 && user->errorindomain)) { 220c4762a1bSJed Brown PetscMPIInt rank; 2219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)user->snes), &rank)); 22248a46eb9SPierre Jolivet if (rank == 0) PetscCall(SNESSetFunctionDomainError(user->snes)); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown grashof = user->grashof; 225c4762a1bSJed Brown prandtl = user->prandtl; 226c4762a1bSJed Brown lid = user->lidvelocity; 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* 229c4762a1bSJed Brown Define mesh intervals ratios for uniform grid. 230c4762a1bSJed Brown 231c4762a1bSJed Brown Note: FD formulae below are normalized by multiplying through by 232c4762a1bSJed Brown local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions. 233c4762a1bSJed Brown 234c4762a1bSJed Brown */ 2359371c9d4SSatish Balay dhx = (PetscReal)(info->mx - 1); 2369371c9d4SSatish Balay dhy = (PetscReal)(info->my - 1); 2379371c9d4SSatish Balay hx = 1.0 / dhx; 2389371c9d4SSatish Balay hy = 1.0 / dhy; 2399371c9d4SSatish Balay hxdhy = hx * dhy; 2409371c9d4SSatish Balay hydhx = hy * dhx; 241c4762a1bSJed Brown 2429371c9d4SSatish Balay xints = info->xs; 2439371c9d4SSatish Balay xinte = info->xs + info->xm; 2449371c9d4SSatish Balay yints = info->ys; 2459371c9d4SSatish Balay yinte = info->ys + info->ym; 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* Test whether we are on the bottom edge of the global array */ 248c4762a1bSJed Brown if (yints == 0) { 249c4762a1bSJed Brown j = 0; 250c4762a1bSJed Brown yints = yints + 1; 251c4762a1bSJed Brown /* bottom edge */ 252c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 253c4762a1bSJed Brown f[j][i].u = x[j][i].u; 254c4762a1bSJed Brown f[j][i].v = x[j][i].v; 255c4762a1bSJed Brown f[j][i].omega = x[j][i].omega + (x[j + 1][i].u - x[j][i].u) * dhy; 256c4762a1bSJed Brown f[j][i].temp = x[j][i].temp - x[j + 1][i].temp; 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* Test whether we are on the top edge of the global array */ 261c4762a1bSJed Brown if (yinte == info->my) { 262c4762a1bSJed Brown j = info->my - 1; 263c4762a1bSJed Brown yinte = yinte - 1; 264c4762a1bSJed Brown /* top edge */ 265c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 266c4762a1bSJed Brown f[j][i].u = x[j][i].u - lid; 267c4762a1bSJed Brown f[j][i].v = x[j][i].v; 268c4762a1bSJed Brown f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j - 1][i].u) * dhy; 269c4762a1bSJed Brown f[j][i].temp = x[j][i].temp - x[j - 1][i].temp; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* Test whether we are on the left edge of the global array */ 274c4762a1bSJed Brown if (xints == 0) { 275c4762a1bSJed Brown i = 0; 276c4762a1bSJed Brown xints = xints + 1; 277c4762a1bSJed Brown /* left edge */ 278c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 279c4762a1bSJed Brown f[j][i].u = x[j][i].u; 280c4762a1bSJed Brown f[j][i].v = x[j][i].v; 281c4762a1bSJed Brown f[j][i].omega = x[j][i].omega - (x[j][i + 1].v - x[j][i].v) * dhx; 282c4762a1bSJed Brown f[j][i].temp = x[j][i].temp; 283c4762a1bSJed Brown } 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* Test whether we are on the right edge of the global array */ 287c4762a1bSJed Brown if (xinte == info->mx) { 288c4762a1bSJed Brown i = info->mx - 1; 289c4762a1bSJed Brown xinte = xinte - 1; 290c4762a1bSJed Brown /* right edge */ 291c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 292c4762a1bSJed Brown f[j][i].u = x[j][i].u; 293c4762a1bSJed Brown f[j][i].v = x[j][i].v; 294c4762a1bSJed Brown f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i - 1].v) * dhx; 295c4762a1bSJed Brown f[j][i].temp = x[j][i].temp - (PetscReal)(grashof > 0); 296c4762a1bSJed Brown } 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* Compute over the interior points */ 300c4762a1bSJed Brown for (j = yints; j < yinte; j++) { 301c4762a1bSJed Brown for (i = xints; i < xinte; i++) { 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown convective coefficients for upwinding 304c4762a1bSJed Brown */ 3059371c9d4SSatish Balay vx = x[j][i].u; 3069371c9d4SSatish Balay avx = PetscAbsScalar(vx); 3079371c9d4SSatish Balay vxp = .5 * (vx + avx); 3089371c9d4SSatish Balay vxm = .5 * (vx - avx); 3099371c9d4SSatish Balay vy = x[j][i].v; 3109371c9d4SSatish Balay avy = PetscAbsScalar(vy); 3119371c9d4SSatish Balay vyp = .5 * (vy + avy); 3129371c9d4SSatish Balay vym = .5 * (vy - avy); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* U velocity */ 315c4762a1bSJed Brown u = x[j][i].u; 316c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1].u - x[j][i + 1].u) * hydhx; 317c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i].u - x[j + 1][i].u) * hxdhy; 318c4762a1bSJed Brown f[j][i].u = uxx + uyy - .5 * (x[j + 1][i].omega - x[j - 1][i].omega) * hx; 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* V velocity */ 321c4762a1bSJed Brown u = x[j][i].v; 322c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1].v - x[j][i + 1].v) * hydhx; 323c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i].v - x[j + 1][i].v) * hxdhy; 324c4762a1bSJed Brown f[j][i].v = uxx + uyy + .5 * (x[j][i + 1].omega - x[j][i - 1].omega) * hy; 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* Omega */ 327c4762a1bSJed Brown u = x[j][i].omega; 328c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1].omega - x[j][i + 1].omega) * hydhx; 329c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i].omega - x[j + 1][i].omega) * hxdhy; 3309371c9d4SSatish Balay f[j][i].omega = uxx + uyy + (vxp * (u - x[j][i - 1].omega) + vxm * (x[j][i + 1].omega - u)) * hy + (vyp * (u - x[j - 1][i].omega) + vym * (x[j + 1][i].omega - u)) * hx - .5 * grashof * (x[j][i + 1].temp - x[j][i - 1].temp) * hy; 331c4762a1bSJed Brown 332c4762a1bSJed Brown /* Temperature */ 333c4762a1bSJed Brown u = x[j][i].temp; 334c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1].temp - x[j][i + 1].temp) * hydhx; 335c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i].temp - x[j + 1][i].temp) * hxdhy; 3369371c9d4SSatish Balay f[j][i].temp = uxx + uyy + prandtl * ((vxp * (u - x[j][i - 1].temp) + vxm * (x[j][i + 1].temp - u)) * hy + (vyp * (u - x[j - 1][i].temp) + vym * (x[j + 1][i].temp - u)) * hx); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* 341c4762a1bSJed Brown Flop count (multiply-adds are counted as 2 operations) 342c4762a1bSJed Brown */ 3439566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(84.0 * info->ym * info->xm)); 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown 347d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_MyShell(Mat A, Vec x, Vec y) 348d71ae5a4SJacob Faibussowitsch { 349c4762a1bSJed Brown MatShellCtx *matshellctx; 350c4762a1bSJed Brown static PetscInt fail = 0; 351c4762a1bSJed Brown 352c4762a1bSJed Brown PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &matshellctx)); 3549566063dSJacob Faibussowitsch PetscCall(MatMult(matshellctx->Jmf, x, y)); 355c4762a1bSJed Brown if (fail++ > 5) { 356c4762a1bSJed Brown PetscMPIInt rank; 3579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 3589566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecSetInf(y)); 359*3872ee93SStefano Zampini PetscCall(VecAssemblyBegin(y)); 360*3872ee93SStefano Zampini PetscCall(VecAssemblyEnd(y)); 361c4762a1bSJed Brown } 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_MyShell(Mat A, MatAssemblyType tp) 366d71ae5a4SJacob Faibussowitsch { 367c4762a1bSJed Brown MatShellCtx *matshellctx; 368c4762a1bSJed Brown 369c4762a1bSJed Brown PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &matshellctx)); 3719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(matshellctx->Jmf, tp)); 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply_MyShell(PC pc, Vec x, Vec y) 376d71ae5a4SJacob Faibussowitsch { 377c4762a1bSJed Brown static PetscInt fail = 0; 378c4762a1bSJed Brown 379c4762a1bSJed Brown PetscFunctionBegin; 3809566063dSJacob Faibussowitsch PetscCall(VecCopy(x, y)); 381c4762a1bSJed Brown if (fail++ > 3) { 382c4762a1bSJed Brown PetscMPIInt rank; 3839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 3849566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecSetInf(y)); 385*3872ee93SStefano Zampini PetscCall(VecAssemblyBegin(y)); 386*3872ee93SStefano Zampini PetscCall(VecAssemblyEnd(y)); 387c4762a1bSJed Brown } 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES, Vec, Mat, Mat, void *); 392c4762a1bSJed Brown 393d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESComputeJacobian_MyShell(SNES snes, Vec X, Mat A, Mat B, void *ctx) 394d71ae5a4SJacob Faibussowitsch { 395c4762a1bSJed Brown static PetscInt fail = 0; 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscFunctionBegin; 3989566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian_DMDA(snes, X, A, B, ctx)); 39948a46eb9SPierre Jolivet if (fail++ > 0) PetscCall(MatZeroEntries(A)); 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown /*TEST 404c4762a1bSJed Brown 405c4762a1bSJed Brown test: 406c4762a1bSJed Brown args: -snes_converged_reason -ksp_converged_reason 407c4762a1bSJed Brown 408c4762a1bSJed Brown test: 409c4762a1bSJed Brown suffix: 2 410308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_matmult -fp_trap 0 411c4762a1bSJed Brown 412c4762a1bSJed Brown test: 413c4762a1bSJed Brown suffix: 3 414308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_pcapply -fp_trap 0 415c4762a1bSJed Brown 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: 4 418308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -fp_trap 0 419c4762a1bSJed Brown 420c4762a1bSJed Brown test: 421c4762a1bSJed Brown suffix: 5 422308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type bjacobi -fp_trap 0 423c4762a1bSJed Brown 424c4762a1bSJed Brown test: 425c4762a1bSJed Brown suffix: 5_fieldsplit 426308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type fieldsplit -fp_trap 0 427c4762a1bSJed Brown output_file: output/ex69_5.out 428c4762a1bSJed Brown 429c4762a1bSJed Brown test: 430c4762a1bSJed Brown suffix: 6 431308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_domainmf -snes_mf -pc_type none -fp_trap 0 432c4762a1bSJed Brown 433c4762a1bSJed Brown test: 434c4762a1bSJed Brown suffix: 7 435308041e4SStefano Zampini args: -snes_converged_reason -ksp_converged_reason -error_in_domain -fp_trap 0 436c4762a1bSJed Brown 437c4762a1bSJed Brown test: 438c4762a1bSJed Brown suffix: 8 439c4762a1bSJed Brown args: -snes_converged_reason -ksp_converged_reason -error_in_domain -snes_mf -pc_type none 440c4762a1bSJed Brown 441c4762a1bSJed Brown TEST*/ 442