1c5566c22SJunchao Zhang static char help[] = "Copy of ex5.c\n"; 2c5566c22SJunchao Zhang 3c5566c22SJunchao Zhang /* ------------------------------------------------------------------------ 4c5566c22SJunchao Zhang 5c5566c22SJunchao Zhang Copy of ex5.c. 6c5566c22SJunchao Zhang Once petsc test harness supports conditional linking, we can remove this duplicate. 7c5566c22SJunchao Zhang See https://gitlab.com/petsc/petsc/-/issues/1173 8c5566c22SJunchao Zhang ------------------------------------------------------------------------- */ 9c5566c22SJunchao Zhang 10c5566c22SJunchao Zhang /* 11c5566c22SJunchao Zhang Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 12c5566c22SJunchao Zhang Include "petscsnes.h" so that we can use SNES solvers. Note that this 13c5566c22SJunchao Zhang */ 14c5566c22SJunchao Zhang #include <petscdm.h> 15c5566c22SJunchao Zhang #include <petscdmda.h> 16c5566c22SJunchao Zhang #include <petscsnes.h> 17c5566c22SJunchao Zhang #include <petscmatlab.h> 18c5566c22SJunchao Zhang #include <petsc/private/snesimpl.h> /* For SNES_Solve event */ 19c5566c22SJunchao Zhang #include "ex55.h" 20c5566c22SJunchao Zhang 21c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */ 22c5566c22SJunchao Zhang /* 23c5566c22SJunchao Zhang FormInitialGuess - Forms initial approximation. 24c5566c22SJunchao Zhang 25c5566c22SJunchao Zhang Input Parameters: 26c5566c22SJunchao Zhang da - The DM 27c5566c22SJunchao Zhang user - user-defined application context 28c5566c22SJunchao Zhang 29c5566c22SJunchao Zhang Output Parameter: 30c5566c22SJunchao Zhang X - vector 31c5566c22SJunchao Zhang */ 32d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormInitialGuess(DM da, AppCtx *user, Vec X) 33d71ae5a4SJacob Faibussowitsch { 34c5566c22SJunchao Zhang PetscInt i, j, Mx, My, xs, ys, xm, ym; 35c5566c22SJunchao Zhang PetscReal lambda, temp1, temp, hx, hy; 36c5566c22SJunchao Zhang PetscScalar **x; 37c5566c22SJunchao Zhang 38c5566c22SJunchao Zhang PetscFunctionBeginUser; 39c5566c22SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 40c5566c22SJunchao Zhang 41c5566c22SJunchao Zhang lambda = user->param; 42c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 43c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 44c5566c22SJunchao Zhang temp1 = lambda / (lambda + 1.0); 45c5566c22SJunchao Zhang 46c5566c22SJunchao Zhang /* 47c5566c22SJunchao Zhang Get a pointer to vector data. 48c5566c22SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 49c5566c22SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 50c5566c22SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 51c5566c22SJunchao Zhang the array. 52c5566c22SJunchao Zhang */ 53c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(da, X, &x)); 54c5566c22SJunchao Zhang 55c5566c22SJunchao Zhang /* 56c5566c22SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 57c5566c22SJunchao Zhang xs, ys - starting grid indices (no ghost points) 58c5566c22SJunchao Zhang xm, ym - widths of local grid (no ghost points) 59c5566c22SJunchao Zhang 60c5566c22SJunchao Zhang */ 61c5566c22SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 62c5566c22SJunchao Zhang 63c5566c22SJunchao Zhang /* 64c5566c22SJunchao Zhang Compute initial guess over the locally owned part of the grid 65c5566c22SJunchao Zhang */ 66c5566c22SJunchao Zhang for (j = ys; j < ys + ym; j++) { 67c5566c22SJunchao Zhang temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy; 68c5566c22SJunchao Zhang for (i = xs; i < xs + xm; i++) { 69c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 70c5566c22SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 71c5566c22SJunchao Zhang x[j][i] = 0.0; 72c5566c22SJunchao Zhang } else { 73c5566c22SJunchao Zhang x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp)); 74c5566c22SJunchao Zhang } 75c5566c22SJunchao Zhang } 76c5566c22SJunchao Zhang } 77c5566c22SJunchao Zhang 78c5566c22SJunchao Zhang /* 79c5566c22SJunchao Zhang Restore vector 80c5566c22SJunchao Zhang */ 81c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, X, &x)); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83c5566c22SJunchao Zhang } 84c5566c22SJunchao Zhang 85c5566c22SJunchao Zhang /* 86c5566c22SJunchao Zhang FormExactSolution - Forms MMS solution 87c5566c22SJunchao Zhang 88c5566c22SJunchao Zhang Input Parameters: 89c5566c22SJunchao Zhang da - The DM 90c5566c22SJunchao Zhang user - user-defined application context 91c5566c22SJunchao Zhang 92c5566c22SJunchao Zhang Output Parameter: 93c5566c22SJunchao Zhang X - vector 94c5566c22SJunchao Zhang */ 95d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U) 96d71ae5a4SJacob Faibussowitsch { 97c5566c22SJunchao Zhang DM coordDA; 98c5566c22SJunchao Zhang Vec coordinates; 99c5566c22SJunchao Zhang DMDACoor2d **coords; 100c5566c22SJunchao Zhang PetscScalar **u; 101c5566c22SJunchao Zhang PetscInt xs, ys, xm, ym, i, j; 102c5566c22SJunchao Zhang 103c5566c22SJunchao Zhang PetscFunctionBeginUser; 104c5566c22SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 105c5566c22SJunchao Zhang PetscCall(DMGetCoordinateDM(da, &coordDA)); 106c5566c22SJunchao Zhang PetscCall(DMGetCoordinates(da, &coordinates)); 107c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 108c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(da, U, &u)); 109c5566c22SJunchao Zhang for (j = ys; j < ys + ym; ++j) { 1103ba16761SJacob Faibussowitsch for (i = xs; i < xs + xm; ++i) PetscCall(user->mms_solution(user, &coords[j][i], &u[j][i])); 111c5566c22SJunchao Zhang } 112c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, U, &u)); 113c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115c5566c22SJunchao Zhang } 116c5566c22SJunchao Zhang 117d71ae5a4SJacob Faibussowitsch static PetscErrorCode ZeroBCSolution(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 118d71ae5a4SJacob Faibussowitsch { 119c5566c22SJunchao Zhang u[0] = 0.; 1203ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 121c5566c22SJunchao Zhang } 122c5566c22SJunchao Zhang 123c5566c22SJunchao Zhang /* The functions below evaluate the MMS solution u(x,y) and associated forcing 124c5566c22SJunchao Zhang 125c5566c22SJunchao Zhang f(x,y) = -u_xx - u_yy - lambda exp(u) 126c5566c22SJunchao Zhang 127c5566c22SJunchao Zhang such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term. 128c5566c22SJunchao Zhang */ 129d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution1(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 130d71ae5a4SJacob Faibussowitsch { 131c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 132c5566c22SJunchao Zhang u[0] = x * (1 - x) * y * (1 - y); 1333ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(5)); 1343ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 135c5566c22SJunchao Zhang } 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing1(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) 137d71ae5a4SJacob Faibussowitsch { 138c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 139c5566c22SJunchao Zhang f[0] = 2 * x * (1 - x) + 2 * y * (1 - y) - user->param * PetscExpReal(x * (1 - x) * y * (1 - y)); 1403ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 141c5566c22SJunchao Zhang } 142c5566c22SJunchao Zhang 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution2(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 144d71ae5a4SJacob Faibussowitsch { 145c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 146c5566c22SJunchao Zhang u[0] = PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y); 1473ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(5)); 1483ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 149c5566c22SJunchao Zhang } 150d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing2(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) 151d71ae5a4SJacob Faibussowitsch { 152c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 153c5566c22SJunchao Zhang f[0] = 2 * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y) - user->param * PetscExpReal(PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y)); 1543ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 155c5566c22SJunchao Zhang } 156c5566c22SJunchao Zhang 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution3(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 158d71ae5a4SJacob Faibussowitsch { 159c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 160c5566c22SJunchao Zhang u[0] = PetscSinReal(user->m * PETSC_PI * x * (1 - y)) * PetscSinReal(user->n * PETSC_PI * y * (1 - x)); 1613ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(5)); 1623ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 163c5566c22SJunchao Zhang } 164d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing3(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) 165d71ae5a4SJacob Faibussowitsch { 166c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 167c5566c22SJunchao Zhang PetscReal m = user->m, n = user->n, lambda = user->param; 1689371c9d4SSatish Balay f[0] = (-(PetscExpReal(PetscSinReal(m * PETSC_PI * x * (1 - y)) * PetscSinReal(n * PETSC_PI * (1 - x) * y)) * lambda) + PetscSqr(PETSC_PI) * (-2 * m * n * ((-1 + x) * x + (-1 + y) * y) * PetscCosReal(m * PETSC_PI * x * (-1 + y)) * PetscCosReal(n * PETSC_PI * (-1 + x) * y) + (PetscSqr(m) * (PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n) * (PetscSqr(-1 + x) + PetscSqr(y))) * PetscSinReal(m * PETSC_PI * x * (-1 + y)) * PetscSinReal(n * PETSC_PI * (-1 + x) * y))); 1693ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 170c5566c22SJunchao Zhang } 171c5566c22SJunchao Zhang 172d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution4(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 173d71ae5a4SJacob Faibussowitsch { 174c5566c22SJunchao Zhang const PetscReal Lx = 1., Ly = 1.; 175c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 176c5566c22SJunchao Zhang u[0] = (PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y)); 1773ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(9)); 1783ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 179c5566c22SJunchao Zhang } 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing4(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) 181d71ae5a4SJacob Faibussowitsch { 182c5566c22SJunchao Zhang const PetscReal Lx = 1., Ly = 1.; 183c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 1849371c9d4SSatish Balay f[0] = (2 * PetscSqr(x) * (PetscSqr(x) - PetscSqr(Lx)) * (PetscSqr(Ly) - 6 * PetscSqr(y)) + 2 * PetscSqr(y) * (PetscSqr(Lx) - 6 * PetscSqr(x)) * (PetscSqr(y) - PetscSqr(Ly)) - user->param * PetscExpReal((PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y)))); 1853ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 186c5566c22SJunchao Zhang } 187c5566c22SJunchao Zhang 188c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */ 189c5566c22SJunchao Zhang /* 190c5566c22SJunchao Zhang FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch 191c5566c22SJunchao Zhang 192c5566c22SJunchao Zhang */ 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) 194d71ae5a4SJacob Faibussowitsch { 195c5566c22SJunchao Zhang PetscInt i, j; 196c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx; 197c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing; 198c5566c22SJunchao Zhang DMDACoor2d c; 199c5566c22SJunchao Zhang 200c5566c22SJunchao Zhang PetscFunctionBeginUser; 201c5566c22SJunchao Zhang lambda = user->param; 202c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 203c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 204c5566c22SJunchao Zhang hxdhy = hx / hy; 205c5566c22SJunchao Zhang hydhx = hy / hx; 206c5566c22SJunchao Zhang /* 207c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 208c5566c22SJunchao Zhang */ 209c5566c22SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 210c5566c22SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 211c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 2129371c9d4SSatish Balay c.x = i * hx; 2139371c9d4SSatish Balay c.y = j * hy; 214c5566c22SJunchao Zhang PetscCall(user->mms_solution(user, &c, &mms_solution)); 215c5566c22SJunchao Zhang f[j][i] = 2.0 * (hydhx + hxdhy) * (x[j][i] - mms_solution); 216c5566c22SJunchao Zhang } else { 217c5566c22SJunchao Zhang u = x[j][i]; 218c5566c22SJunchao Zhang uw = x[j][i - 1]; 219c5566c22SJunchao Zhang ue = x[j][i + 1]; 220c5566c22SJunchao Zhang un = x[j - 1][i]; 221c5566c22SJunchao Zhang us = x[j + 1][i]; 222c5566c22SJunchao Zhang 223c5566c22SJunchao Zhang /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 2249371c9d4SSatish Balay if (i - 1 == 0) { 2259371c9d4SSatish Balay c.x = (i - 1) * hx; 2269371c9d4SSatish Balay c.y = j * hy; 2279371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &uw)); 2289371c9d4SSatish Balay } 2299371c9d4SSatish Balay if (i + 1 == info->mx - 1) { 2309371c9d4SSatish Balay c.x = (i + 1) * hx; 2319371c9d4SSatish Balay c.y = j * hy; 2329371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &ue)); 2339371c9d4SSatish Balay } 2349371c9d4SSatish Balay if (j - 1 == 0) { 2359371c9d4SSatish Balay c.x = i * hx; 2369371c9d4SSatish Balay c.y = (j - 1) * hy; 2379371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &un)); 2389371c9d4SSatish Balay } 2399371c9d4SSatish Balay if (j + 1 == info->my - 1) { 2409371c9d4SSatish Balay c.x = i * hx; 2419371c9d4SSatish Balay c.y = (j + 1) * hy; 2429371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &us)); 2439371c9d4SSatish Balay } 244c5566c22SJunchao Zhang 245c5566c22SJunchao Zhang uxx = (2.0 * u - uw - ue) * hydhx; 246c5566c22SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 247c5566c22SJunchao Zhang mms_forcing = 0; 2489371c9d4SSatish Balay c.x = i * hx; 2499371c9d4SSatish Balay c.y = j * hy; 250c5566c22SJunchao Zhang if (user->mms_forcing) PetscCall(user->mms_forcing(user, &c, &mms_forcing)); 251c5566c22SJunchao Zhang f[j][i] = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing); 252c5566c22SJunchao Zhang } 253c5566c22SJunchao Zhang } 254c5566c22SJunchao Zhang } 255c5566c22SJunchao Zhang PetscCall(PetscLogFlops(11.0 * info->ym * info->xm)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257c5566c22SJunchao Zhang } 258c5566c22SJunchao Zhang 259c5566c22SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */ 260d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *obj, AppCtx *user) 261d71ae5a4SJacob Faibussowitsch { 262c5566c22SJunchao Zhang PetscInt i, j; 263c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0; 264c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxux, uyuy; 265c5566c22SJunchao Zhang MPI_Comm comm; 266c5566c22SJunchao Zhang 267c5566c22SJunchao Zhang PetscFunctionBeginUser; 268c5566c22SJunchao Zhang *obj = 0; 269c5566c22SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm)); 270c5566c22SJunchao Zhang lambda = user->param; 271c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 272c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 273c5566c22SJunchao Zhang sc = hx * hy * lambda; 274c5566c22SJunchao Zhang hxdhy = hx / hy; 275c5566c22SJunchao Zhang hydhx = hy / hx; 276c5566c22SJunchao Zhang /* 277c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 278c5566c22SJunchao Zhang */ 279c5566c22SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 280c5566c22SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 281c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 282c5566c22SJunchao Zhang lobj += PetscRealPart((hydhx + hxdhy) * x[j][i] * x[j][i]); 283c5566c22SJunchao Zhang } else { 284c5566c22SJunchao Zhang u = x[j][i]; 285c5566c22SJunchao Zhang uw = x[j][i - 1]; 286c5566c22SJunchao Zhang ue = x[j][i + 1]; 287c5566c22SJunchao Zhang un = x[j - 1][i]; 288c5566c22SJunchao Zhang us = x[j + 1][i]; 289c5566c22SJunchao Zhang 290c5566c22SJunchao Zhang if (i - 1 == 0) uw = 0.; 291c5566c22SJunchao Zhang if (i + 1 == info->mx - 1) ue = 0.; 292c5566c22SJunchao Zhang if (j - 1 == 0) un = 0.; 293c5566c22SJunchao Zhang if (j + 1 == info->my - 1) us = 0.; 294c5566c22SJunchao Zhang 295c5566c22SJunchao Zhang /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 296c5566c22SJunchao Zhang 297c5566c22SJunchao Zhang uxux = u * (2. * u - ue - uw) * hydhx; 298c5566c22SJunchao Zhang uyuy = u * (2. * u - un - us) * hxdhy; 299c5566c22SJunchao Zhang 300c5566c22SJunchao Zhang lobj += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u)); 301c5566c22SJunchao Zhang } 302c5566c22SJunchao Zhang } 303c5566c22SJunchao Zhang } 304c5566c22SJunchao Zhang PetscCall(PetscLogFlops(12.0 * info->ym * info->xm)); 305*0bf52853SStefano Zampini *obj = lobj; 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 307c5566c22SJunchao Zhang } 308c5566c22SJunchao Zhang 309c5566c22SJunchao Zhang /* 310c5566c22SJunchao Zhang FormJacobianLocal - Evaluates Jacobian matrix on local process patch 311c5566c22SJunchao Zhang */ 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat jac, Mat jacpre, AppCtx *user) 313d71ae5a4SJacob Faibussowitsch { 314c5566c22SJunchao Zhang PetscInt i, j, k; 315c5566c22SJunchao Zhang MatStencil col[5], row; 316c5566c22SJunchao Zhang PetscScalar lambda, v[5], hx, hy, hxdhy, hydhx, sc; 317c5566c22SJunchao Zhang DM coordDA; 318c5566c22SJunchao Zhang Vec coordinates; 319c5566c22SJunchao Zhang DMDACoor2d **coords; 320c5566c22SJunchao Zhang 321c5566c22SJunchao Zhang PetscFunctionBeginUser; 322c5566c22SJunchao Zhang lambda = user->param; 323c5566c22SJunchao Zhang /* Extract coordinates */ 324c5566c22SJunchao Zhang PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 325c5566c22SJunchao Zhang PetscCall(DMGetCoordinates(info->da, &coordinates)); 326c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 327c5566c22SJunchao Zhang hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs + 1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0; 328c5566c22SJunchao Zhang hy = info->ym > 1 ? PetscRealPart(coords[info->ys + 1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0; 329c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 330c5566c22SJunchao Zhang hxdhy = hx / hy; 331c5566c22SJunchao Zhang hydhx = hy / hx; 332c5566c22SJunchao Zhang sc = hx * hy * lambda; 333c5566c22SJunchao Zhang 334c5566c22SJunchao Zhang /* 335c5566c22SJunchao Zhang Compute entries for the locally owned part of the Jacobian. 336c5566c22SJunchao Zhang - Currently, all PETSc parallel matrix formats are partitioned by 337c5566c22SJunchao Zhang contiguous chunks of rows across the processors. 338c5566c22SJunchao Zhang - Each processor needs to insert only elements that it owns 339c5566c22SJunchao Zhang locally (but any non-local elements will be sent to the 340c5566c22SJunchao Zhang appropriate processor during matrix assembly). 341c5566c22SJunchao Zhang - Here, we set all entries for a particular row at once. 342c5566c22SJunchao Zhang - We can set matrix entries either using either 343c5566c22SJunchao Zhang MatSetValuesLocal() or MatSetValues(), as discussed above. 344c5566c22SJunchao Zhang */ 345c5566c22SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 346c5566c22SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 3479371c9d4SSatish Balay row.j = j; 3489371c9d4SSatish Balay row.i = i; 349c5566c22SJunchao Zhang /* boundary points */ 350c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 351c5566c22SJunchao Zhang v[0] = 2.0 * (hydhx + hxdhy); 352c5566c22SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre, 1, &row, 1, &row, v, INSERT_VALUES)); 353c5566c22SJunchao Zhang } else { 354c5566c22SJunchao Zhang k = 0; 355c5566c22SJunchao Zhang /* interior grid points */ 356c5566c22SJunchao Zhang if (j - 1 != 0) { 357c5566c22SJunchao Zhang v[k] = -hxdhy; 3589371c9d4SSatish Balay col[k].j = j - 1; 3599371c9d4SSatish Balay col[k].i = i; 360c5566c22SJunchao Zhang k++; 361c5566c22SJunchao Zhang } 362c5566c22SJunchao Zhang if (i - 1 != 0) { 363c5566c22SJunchao Zhang v[k] = -hydhx; 3649371c9d4SSatish Balay col[k].j = j; 3659371c9d4SSatish Balay col[k].i = i - 1; 366c5566c22SJunchao Zhang k++; 367c5566c22SJunchao Zhang } 368c5566c22SJunchao Zhang 3699371c9d4SSatish Balay v[k] = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(x[j][i]); 3709371c9d4SSatish Balay col[k].j = row.j; 3719371c9d4SSatish Balay col[k].i = row.i; 3729371c9d4SSatish Balay k++; 373c5566c22SJunchao Zhang 374c5566c22SJunchao Zhang if (i + 1 != info->mx - 1) { 375c5566c22SJunchao Zhang v[k] = -hydhx; 3769371c9d4SSatish Balay col[k].j = j; 3779371c9d4SSatish Balay col[k].i = i + 1; 378c5566c22SJunchao Zhang k++; 379c5566c22SJunchao Zhang } 380c5566c22SJunchao Zhang if (j + 1 != info->mx - 1) { 381c5566c22SJunchao Zhang v[k] = -hxdhy; 3829371c9d4SSatish Balay col[k].j = j + 1; 3839371c9d4SSatish Balay col[k].i = i; 384c5566c22SJunchao Zhang k++; 385c5566c22SJunchao Zhang } 386c5566c22SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre, 1, &row, k, col, v, INSERT_VALUES)); 387c5566c22SJunchao Zhang } 388c5566c22SJunchao Zhang } 389c5566c22SJunchao Zhang } 390c5566c22SJunchao Zhang 391c5566c22SJunchao Zhang /* 392c5566c22SJunchao Zhang Assemble matrix, using the 2-step process: 393c5566c22SJunchao Zhang MatAssemblyBegin(), MatAssemblyEnd(). 394c5566c22SJunchao Zhang */ 395c5566c22SJunchao Zhang PetscCall(MatAssemblyBegin(jacpre, MAT_FINAL_ASSEMBLY)); 396c5566c22SJunchao Zhang PetscCall(MatAssemblyEnd(jacpre, MAT_FINAL_ASSEMBLY)); 397c5566c22SJunchao Zhang 398c5566c22SJunchao Zhang /* 399c5566c22SJunchao Zhang Tell the matrix we will never add a new nonzero location to the 400c5566c22SJunchao Zhang matrix. If we do, it will generate an error. 401c5566c22SJunchao Zhang */ 402c5566c22SJunchao Zhang PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 404c5566c22SJunchao Zhang } 405c5566c22SJunchao Zhang 406d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionMatlab(SNES snes, Vec X, Vec F, void *ptr) 407d71ae5a4SJacob Faibussowitsch { 408d1e78c4fSBarry Smith #if PetscDefined(HAVE_MATLAB) 409c5566c22SJunchao Zhang AppCtx *user = (AppCtx *)ptr; 410c5566c22SJunchao Zhang PetscInt Mx, My; 411c5566c22SJunchao Zhang PetscReal lambda, hx, hy; 412c5566c22SJunchao Zhang Vec localX, localF; 413c5566c22SJunchao Zhang MPI_Comm comm; 414c5566c22SJunchao Zhang DM da; 415c5566c22SJunchao Zhang 416c5566c22SJunchao Zhang PetscFunctionBeginUser; 417c5566c22SJunchao Zhang PetscCall(SNESGetDM(snes, &da)); 418c5566c22SJunchao Zhang PetscCall(DMGetLocalVector(da, &localX)); 419c5566c22SJunchao Zhang PetscCall(DMGetLocalVector(da, &localF)); 420c5566c22SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localX, "localX")); 421c5566c22SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localF, "localF")); 422c5566c22SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 423c5566c22SJunchao Zhang 424c5566c22SJunchao Zhang lambda = user->param; 425c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 426c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 427c5566c22SJunchao Zhang 428c5566c22SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 429c5566c22SJunchao Zhang /* 430c5566c22SJunchao Zhang Scatter ghost points to local vector,using the 2-step process 431c5566c22SJunchao Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 432c5566c22SJunchao Zhang By placing code between these two statements, computations can be 433c5566c22SJunchao Zhang done while messages are in transition. 434c5566c22SJunchao Zhang */ 435c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 436c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 437c5566c22SJunchao Zhang PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localX)); 438c5566c22SJunchao Zhang PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm), "localF=ex5m(localX,%18.16e,%18.16e,%18.16e)", (double)hx, (double)hy, (double)lambda)); 439c5566c22SJunchao Zhang PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localF)); 440c5566c22SJunchao Zhang 441c5566c22SJunchao Zhang /* 442c5566c22SJunchao Zhang Insert values into global vector 443c5566c22SJunchao Zhang */ 444c5566c22SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da, localF, INSERT_VALUES, F)); 445c5566c22SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da, localF, INSERT_VALUES, F)); 446c5566c22SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localX)); 447c5566c22SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localF)); 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 449c5566c22SJunchao Zhang #else 4503ba16761SJacob Faibussowitsch return PETSC_SUCCESS; /* Never called */ 451c5566c22SJunchao Zhang #endif 452c5566c22SJunchao Zhang } 453c5566c22SJunchao Zhang 454c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */ 455c5566c22SJunchao Zhang /* 456c5566c22SJunchao Zhang Applies some sweeps on nonlinear Gauss-Seidel on each process 457c5566c22SJunchao Zhang 458c5566c22SJunchao Zhang */ 459d71ae5a4SJacob Faibussowitsch static PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx) 460d71ae5a4SJacob Faibussowitsch { 461c5566c22SJunchao Zhang PetscInt i, j, k, Mx, My, xs, ys, xm, ym, its, tot_its, sweeps, l; 462c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc; 463c5566c22SJunchao Zhang PetscScalar **x, **b, bij, F, F0 = 0, J, u, un, us, ue, eu, uw, uxx, uyy, y; 464c5566c22SJunchao Zhang PetscReal atol, rtol, stol; 465c5566c22SJunchao Zhang DM da; 466c5566c22SJunchao Zhang AppCtx *user; 467c5566c22SJunchao Zhang Vec localX, localB; 468c5566c22SJunchao Zhang 469c5566c22SJunchao Zhang PetscFunctionBeginUser; 470c5566c22SJunchao Zhang tot_its = 0; 471c5566c22SJunchao Zhang PetscCall(SNESNGSGetSweeps(snes, &sweeps)); 472c5566c22SJunchao Zhang PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its)); 473c5566c22SJunchao Zhang PetscCall(SNESGetDM(snes, &da)); 474c5566c22SJunchao Zhang PetscCall(DMGetApplicationContext(da, &user)); 475c5566c22SJunchao Zhang 476c5566c22SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 477c5566c22SJunchao Zhang 478c5566c22SJunchao Zhang lambda = user->param; 479c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 480c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 481c5566c22SJunchao Zhang sc = hx * hy * lambda; 482c5566c22SJunchao Zhang hxdhy = hx / hy; 483c5566c22SJunchao Zhang hydhx = hy / hx; 484c5566c22SJunchao Zhang 485c5566c22SJunchao Zhang PetscCall(DMGetLocalVector(da, &localX)); 48648a46eb9SPierre Jolivet if (B) PetscCall(DMGetLocalVector(da, &localB)); 487c5566c22SJunchao Zhang for (l = 0; l < sweeps; l++) { 488c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 489c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 490c5566c22SJunchao Zhang if (B) { 491c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB)); 492c5566c22SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB)); 493c5566c22SJunchao Zhang } 494c5566c22SJunchao Zhang /* 495c5566c22SJunchao Zhang Get a pointer to vector data. 496c5566c22SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 497c5566c22SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 498c5566c22SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 499c5566c22SJunchao Zhang the array. 500c5566c22SJunchao Zhang */ 501c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(da, localX, &x)); 502c5566c22SJunchao Zhang if (B) PetscCall(DMDAVecGetArray(da, localB, &b)); 503c5566c22SJunchao Zhang /* 504c5566c22SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 505c5566c22SJunchao Zhang xs, ys - starting grid indices (no ghost points) 506c5566c22SJunchao Zhang xm, ym - widths of local grid (no ghost points) 507c5566c22SJunchao Zhang */ 508c5566c22SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 509c5566c22SJunchao Zhang 510c5566c22SJunchao Zhang for (j = ys; j < ys + ym; j++) { 511c5566c22SJunchao Zhang for (i = xs; i < xs + xm; i++) { 512c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 513c5566c22SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 514c5566c22SJunchao Zhang x[j][i] = 0.0; 515c5566c22SJunchao Zhang } else { 516c5566c22SJunchao Zhang if (B) bij = b[j][i]; 517c5566c22SJunchao Zhang else bij = 0.; 518c5566c22SJunchao Zhang 519c5566c22SJunchao Zhang u = x[j][i]; 520c5566c22SJunchao Zhang un = x[j - 1][i]; 521c5566c22SJunchao Zhang us = x[j + 1][i]; 522c5566c22SJunchao Zhang ue = x[j][i - 1]; 523c5566c22SJunchao Zhang uw = x[j][i + 1]; 524c5566c22SJunchao Zhang 525c5566c22SJunchao Zhang for (k = 0; k < its; k++) { 526c5566c22SJunchao Zhang eu = PetscExpScalar(u); 527c5566c22SJunchao Zhang uxx = (2.0 * u - ue - uw) * hydhx; 528c5566c22SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 529c5566c22SJunchao Zhang F = uxx + uyy - sc * eu - bij; 530c5566c22SJunchao Zhang if (k == 0) F0 = F; 531c5566c22SJunchao Zhang J = 2.0 * (hydhx + hxdhy) - sc * eu; 532c5566c22SJunchao Zhang y = F / J; 533c5566c22SJunchao Zhang u -= y; 534c5566c22SJunchao Zhang tot_its++; 535c5566c22SJunchao Zhang 536ad540459SPierre Jolivet if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) break; 537c5566c22SJunchao Zhang } 538c5566c22SJunchao Zhang x[j][i] = u; 539c5566c22SJunchao Zhang } 540c5566c22SJunchao Zhang } 541c5566c22SJunchao Zhang } 542c5566c22SJunchao Zhang /* 543c5566c22SJunchao Zhang Restore vector 544c5566c22SJunchao Zhang */ 545c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, localX, &x)); 546c5566c22SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X)); 547c5566c22SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X)); 548c5566c22SJunchao Zhang } 549c5566c22SJunchao Zhang PetscCall(PetscLogFlops(tot_its * (21.0))); 550c5566c22SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localX)); 551c5566c22SJunchao Zhang if (B) { 552c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, localB, &b)); 553c5566c22SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localB)); 554c5566c22SJunchao Zhang } 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 556c5566c22SJunchao Zhang } 557c5566c22SJunchao Zhang 558d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 559d71ae5a4SJacob Faibussowitsch { 560c5566c22SJunchao Zhang SNES snes; /* nonlinear solver */ 561c5566c22SJunchao Zhang Vec x; /* solution vector */ 562c5566c22SJunchao Zhang AppCtx user; /* user-defined work context */ 563c5566c22SJunchao Zhang PetscInt its; /* iterations for convergence */ 564c5566c22SJunchao Zhang PetscReal bratu_lambda_max = 6.81; 565c5566c22SJunchao Zhang PetscReal bratu_lambda_min = 0.; 566c5566c22SJunchao Zhang PetscInt MMS = 1; 5670aeb1f43SStefano Zampini PetscBool flg, setMMS; 568c5566c22SJunchao Zhang DM da; 569c5566c22SJunchao Zhang Vec r = NULL; 570c5566c22SJunchao Zhang KSP ksp; 571c5566c22SJunchao Zhang PetscInt lits, slits; 572c5566c22SJunchao Zhang PetscBool useKokkos = PETSC_FALSE; 573c5566c22SJunchao Zhang 574c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 575c5566c22SJunchao Zhang Initialize program 576c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 577c5566c22SJunchao Zhang 578327415f7SBarry Smith PetscFunctionBeginUser; 579c5566c22SJunchao Zhang PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 580c5566c22SJunchao Zhang 581c5566c22SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_kokkos", &useKokkos, NULL)); 582c5566c22SJunchao Zhang 583c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 584c5566c22SJunchao Zhang Initialize problem parameters 585c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 5860aeb1f43SStefano Zampini user.ncoo = 0; 587c5566c22SJunchao Zhang user.param = 6.0; 588c5566c22SJunchao Zhang PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 589c5566c22SJunchao Zhang PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max); 590c5566c22SJunchao Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mms", &MMS, &setMMS)); 591c5566c22SJunchao Zhang if (MMS == 3) { 592c5566c22SJunchao Zhang PetscInt mPar = 2, nPar = 1; 593c5566c22SJunchao Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_par", &mPar, NULL)); 594c5566c22SJunchao Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_par", &nPar, NULL)); 595c5566c22SJunchao Zhang user.m = PetscPowInt(2, mPar); 596c5566c22SJunchao Zhang user.n = PetscPowInt(2, nPar); 597c5566c22SJunchao Zhang } 598c5566c22SJunchao Zhang 599c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 600c5566c22SJunchao Zhang Create nonlinear solver context 601c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 602c5566c22SJunchao Zhang PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 603c5566c22SJunchao Zhang PetscCall(SNESSetCountersReset(snes, PETSC_FALSE)); 604c5566c22SJunchao Zhang PetscCall(SNESSetNGS(snes, NonlinearGS, NULL)); 605c5566c22SJunchao Zhang 606c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 607c5566c22SJunchao Zhang Create distributed array (DMDA) to manage parallel grid and vectors 608c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 609c5566c22SJunchao Zhang PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 6100aeb1f43SStefano Zampini if (useKokkos) { 6110aeb1f43SStefano Zampini PetscCall(DMSetVecType(da, VECKOKKOS)); 6120aeb1f43SStefano Zampini PetscCall(DMSetMatType(da, MATAIJKOKKOS)); 6130aeb1f43SStefano Zampini } 614c5566c22SJunchao Zhang PetscCall(DMSetFromOptions(da)); 615c5566c22SJunchao Zhang PetscCall(DMSetUp(da)); 616c5566c22SJunchao Zhang PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 617c5566c22SJunchao Zhang PetscCall(DMSetApplicationContext(da, &user)); 618c5566c22SJunchao Zhang PetscCall(SNESSetDM(snes, da)); 619c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 620c5566c22SJunchao Zhang Extract global vectors from DMDA; then duplicate for remaining 621c5566c22SJunchao Zhang vectors that are the same types 622c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 623c5566c22SJunchao Zhang PetscCall(DMCreateGlobalVector(da, &x)); 624c5566c22SJunchao Zhang 625c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 626c5566c22SJunchao Zhang Set local function evaluation routine 627c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 628c5566c22SJunchao Zhang switch (MMS) { 6299371c9d4SSatish Balay case 0: 6309371c9d4SSatish Balay user.mms_solution = ZeroBCSolution; 6319371c9d4SSatish Balay user.mms_forcing = NULL; 6329371c9d4SSatish Balay break; 6339371c9d4SSatish Balay case 1: 6349371c9d4SSatish Balay user.mms_solution = MMSSolution1; 6359371c9d4SSatish Balay user.mms_forcing = MMSForcing1; 6369371c9d4SSatish Balay break; 6379371c9d4SSatish Balay case 2: 6389371c9d4SSatish Balay user.mms_solution = MMSSolution2; 6399371c9d4SSatish Balay user.mms_forcing = MMSForcing2; 6409371c9d4SSatish Balay break; 6419371c9d4SSatish Balay case 3: 6429371c9d4SSatish Balay user.mms_solution = MMSSolution3; 6439371c9d4SSatish Balay user.mms_forcing = MMSForcing3; 6449371c9d4SSatish Balay break; 6459371c9d4SSatish Balay case 4: 6469371c9d4SSatish Balay user.mms_solution = MMSSolution4; 6479371c9d4SSatish Balay user.mms_forcing = MMSForcing4; 6489371c9d4SSatish Balay break; 649d71ae5a4SJacob Faibussowitsch default: 650d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Unknown MMS type %" PetscInt_FMT, MMS); 651c5566c22SJunchao Zhang } 652c5566c22SJunchao Zhang 653c5566c22SJunchao Zhang if (useKokkos) { 654c5566c22SJunchao Zhang PetscCheck(MMS == 1, PETSC_COMM_WORLD, PETSC_ERR_USER, "FormFunctionLocalVec_Kokkos only works with MMS 1"); 655c5566c22SJunchao Zhang PetscCall(DMDASNESSetFunctionLocalVec(da, INSERT_VALUES, (DMDASNESFunctionVec)FormFunctionLocalVec, &user)); 656c5566c22SJunchao Zhang } else { 657c5566c22SJunchao Zhang PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunction)FormFunctionLocal, &user)); 658c5566c22SJunchao Zhang } 659c5566c22SJunchao Zhang 6600aeb1f43SStefano Zampini flg = PETSC_FALSE; 661c5566c22SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-fd", &flg, NULL)); 662c5566c22SJunchao Zhang if (!flg) { 663c5566c22SJunchao Zhang if (useKokkos) PetscCall(DMDASNESSetJacobianLocalVec(da, (DMDASNESJacobianVec)FormJacobianLocalVec, &user)); 664c5566c22SJunchao Zhang else PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, &user)); 665c5566c22SJunchao Zhang } 666c5566c22SJunchao Zhang 6670aeb1f43SStefano Zampini flg = PETSC_FALSE; 668c5566c22SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-obj", &flg, NULL)); 669c5566c22SJunchao Zhang if (flg) { 670c5566c22SJunchao Zhang if (useKokkos) PetscCall(DMDASNESSetObjectiveLocalVec(da, (DMDASNESObjectiveVec)FormObjectiveLocalVec, &user)); 671c5566c22SJunchao Zhang else PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjective)FormObjectiveLocal, &user)); 672c5566c22SJunchao Zhang } 673c5566c22SJunchao Zhang 674d1e78c4fSBarry Smith if (PetscDefined(HAVE_MATLAB)) { 675c5566c22SJunchao Zhang PetscBool matlab_function = PETSC_FALSE; 676c5566c22SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-matlab_function", &matlab_function, 0)); 677c5566c22SJunchao Zhang if (matlab_function) { 678c5566c22SJunchao Zhang PetscCall(VecDuplicate(x, &r)); 679c5566c22SJunchao Zhang PetscCall(SNESSetFunction(snes, r, FormFunctionMatlab, &user)); 680c5566c22SJunchao Zhang } 681c5566c22SJunchao Zhang } 682c5566c22SJunchao Zhang 683c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 684c5566c22SJunchao Zhang Customize nonlinear solver; set runtime options 685c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 686c5566c22SJunchao Zhang PetscCall(SNESSetFromOptions(snes)); 687c5566c22SJunchao Zhang 688c5566c22SJunchao Zhang PetscCall(FormInitialGuess(da, &user, x)); 689c5566c22SJunchao Zhang 690c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 691c5566c22SJunchao Zhang Solve nonlinear system 692c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 693c5566c22SJunchao Zhang PetscCall(SNESSolve(snes, NULL, x)); 694c5566c22SJunchao Zhang PetscCall(SNESGetIterationNumber(snes, &its)); 695c5566c22SJunchao Zhang 696c5566c22SJunchao Zhang PetscCall(SNESGetLinearSolveIterations(snes, &slits)); 697c5566c22SJunchao Zhang PetscCall(SNESGetKSP(snes, &ksp)); 698c5566c22SJunchao Zhang PetscCall(KSPGetTotalIterations(ksp, &lits)); 699c5566c22SJunchao Zhang PetscCheck(lits == slits, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT, slits, lits); 700c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 701c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 702c5566c22SJunchao Zhang 703c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 704c5566c22SJunchao Zhang If using MMS, check the l_2 error 705c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 706c5566c22SJunchao Zhang if (setMMS) { 707c5566c22SJunchao Zhang Vec e; 708c5566c22SJunchao Zhang PetscReal errorl2, errorinf; 709c5566c22SJunchao Zhang PetscInt N; 710c5566c22SJunchao Zhang 711c5566c22SJunchao Zhang PetscCall(VecDuplicate(x, &e)); 712c5566c22SJunchao Zhang PetscCall(PetscObjectViewFromOptions((PetscObject)x, NULL, "-sol_view")); 713c5566c22SJunchao Zhang PetscCall(FormExactSolution(da, &user, e)); 714c5566c22SJunchao Zhang PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-exact_view")); 715c5566c22SJunchao Zhang PetscCall(VecAXPY(e, -1.0, x)); 716c5566c22SJunchao Zhang PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-error_view")); 717c5566c22SJunchao Zhang PetscCall(VecNorm(e, NORM_2, &errorl2)); 718c5566c22SJunchao Zhang PetscCall(VecNorm(e, NORM_INFINITY, &errorinf)); 719c5566c22SJunchao Zhang PetscCall(VecGetSize(e, &N)); 720c5566c22SJunchao Zhang PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2 / PetscSqrtReal((PetscReal)N)), (double)errorinf)); 721c5566c22SJunchao Zhang PetscCall(VecDestroy(&e)); 722c5566c22SJunchao Zhang PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N)); 723c5566c22SJunchao Zhang PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2 / PetscSqrtReal(N))); 724c5566c22SJunchao Zhang } 725c5566c22SJunchao Zhang 726c5566c22SJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 727c5566c22SJunchao Zhang Free work space. All PETSc objects should be destroyed when they 728c5566c22SJunchao Zhang are no longer needed. 729c5566c22SJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 730c5566c22SJunchao Zhang PetscCall(VecDestroy(&r)); 731c5566c22SJunchao Zhang PetscCall(VecDestroy(&x)); 732c5566c22SJunchao Zhang PetscCall(SNESDestroy(&snes)); 733c5566c22SJunchao Zhang PetscCall(DMDestroy(&da)); 734c5566c22SJunchao Zhang PetscCall(PetscFinalize()); 735c5566c22SJunchao Zhang return 0; 736c5566c22SJunchao Zhang } 737c5566c22SJunchao Zhang 738c5566c22SJunchao Zhang /*TEST 739c5566c22SJunchao Zhang build: 7400aeb1f43SStefano Zampini requires: !windows_compilers 741c5566c22SJunchao Zhang depends: ex55k.kokkos.cxx 742c5566c22SJunchao Zhang 743c5566c22SJunchao Zhang testset: 744c5566c22SJunchao Zhang output_file: output/ex55_asm_0.out 745c5566c22SJunchao Zhang requires: !single 746c5566c22SJunchao Zhang args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu 747c5566c22SJunchao Zhang filter: grep -v "type" 748c5566c22SJunchao Zhang 749c5566c22SJunchao Zhang test: 750c5566c22SJunchao Zhang suffix: asm_0 7510aeb1f43SStefano Zampini args: -fd {{0 1}} 7520aeb1f43SStefano Zampini 753c5566c22SJunchao Zhang test: 7540aeb1f43SStefano Zampini requires: kokkos_kernels 755c5566c22SJunchao Zhang suffix: asm_0_kok 7560aeb1f43SStefano Zampini args: -use_kokkos -fd {{0 1}} 7570aeb1f43SStefano Zampini 7580aeb1f43SStefano Zampini testset: 7590aeb1f43SStefano Zampini output_file: output/ex55_1.out 7600aeb1f43SStefano Zampini requires: !single 7610aeb1f43SStefano Zampini args: -snes_monitor 7620aeb1f43SStefano Zampini filter: grep -v "type" 7630aeb1f43SStefano Zampini 7640aeb1f43SStefano Zampini test: 7650aeb1f43SStefano Zampini suffix: 1 7660aeb1f43SStefano Zampini args: -fd {{0 1}} 7670aeb1f43SStefano Zampini 7680aeb1f43SStefano Zampini test: 7690aeb1f43SStefano Zampini requires: kokkos_kernels 7700aeb1f43SStefano Zampini suffix: 1_kok 7710aeb1f43SStefano Zampini args: -use_kokkos -fd {{0 1}} 7720aeb1f43SStefano Zampini 7730aeb1f43SStefano Zampini test: 7740aeb1f43SStefano Zampini requires: h2opus 7750aeb1f43SStefano Zampini suffix: 1_h2opus 7760aeb1f43SStefano Zampini args: -pc_type h2opus -fd {{0 1}} 7770aeb1f43SStefano Zampini 7780aeb1f43SStefano Zampini test: 7790aeb1f43SStefano Zampini requires: h2opus kokkos_kernels 7800aeb1f43SStefano Zampini suffix: 1_h2opus_k 7810aeb1f43SStefano Zampini args: -use_kokkos -pc_type h2opus -fd {{0 1}} 782c5566c22SJunchao Zhang 783c5566c22SJunchao Zhang TEST*/ 784