10aeb1f43SStefano Zampini #include <petscconf.h> 20aeb1f43SStefano Zampini 30aeb1f43SStefano Zampini #if defined(PETSC_HAVE_KOKKOS_KERNELS) 4c5566c22SJunchao Zhang #include <Kokkos_Core.hpp> 5*1693298eSJunchao Zhang #include <petsc_kokkos.hpp> 6c5566c22SJunchao Zhang #include <petscdmda_kokkos.hpp> 7c5566c22SJunchao Zhang 8c5566c22SJunchao Zhang #include <petscdm.h> 9c5566c22SJunchao Zhang #include <petscdmda.h> 10c5566c22SJunchao Zhang #include <petscsnes.h> 11c5566c22SJunchao Zhang #include "ex55.h" 12c5566c22SJunchao Zhang 13c5566c22SJunchao Zhang using DefaultMemorySpace = Kokkos::DefaultExecutionSpace::memory_space; 14c5566c22SJunchao Zhang using ConstPetscScalarKokkosOffsetView2D = Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, DefaultMemorySpace>; 15c5566c22SJunchao Zhang using PetscScalarKokkosOffsetView2D = Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, DefaultMemorySpace>; 16c5566c22SJunchao Zhang 17c5566c22SJunchao Zhang using PetscCountKokkosView = Kokkos::View<PetscCount *, DefaultMemorySpace>; 18c5566c22SJunchao Zhang using PetscIntKokkosView = Kokkos::View<PetscInt *, DefaultMemorySpace>; 19c5566c22SJunchao Zhang using PetscScalarKokkosView = Kokkos::View<PetscScalar *, DefaultMemorySpace>; 2033d966b4SJunchao Zhang using Kokkos::Iterate; 2133d966b4SJunchao Zhang using Kokkos::MDRangePolicy; 229371c9d4SSatish Balay using Kokkos::Rank; 23c5566c22SJunchao Zhang 24d71ae5a4SJacob Faibussowitsch KOKKOS_INLINE_FUNCTION PetscErrorCode MMSSolution1(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) 25d71ae5a4SJacob Faibussowitsch { 26c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 27c5566c22SJunchao Zhang u[0] = x * (1 - x) * y * (1 - y); 283ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 29c5566c22SJunchao Zhang } 30c5566c22SJunchao Zhang 31d71ae5a4SJacob Faibussowitsch KOKKOS_INLINE_FUNCTION PetscErrorCode MMSForcing1(PetscReal user_param, const DMDACoor2d *c, PetscScalar *f) 32d71ae5a4SJacob Faibussowitsch { 33c5566c22SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 34c5566c22SJunchao Zhang f[0] = 2 * x * (1 - x) + 2 * y * (1 - y) - user_param * PetscExpReal(x * (1 - x) * y * (1 - y)); 353ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 36c5566c22SJunchao Zhang } 37c5566c22SJunchao Zhang 38d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocalVec(DMDALocalInfo *info, Vec x, Vec f, AppCtx *user) 39d71ae5a4SJacob Faibussowitsch { 40c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx; 41c5566c22SJunchao Zhang PetscInt xs = info->xs, ys = info->ys, xm = info->xm, ym = info->ym, mx = info->mx, my = info->my; 42c5566c22SJunchao Zhang PetscReal user_param = user->param; 43c5566c22SJunchao Zhang 44c5566c22SJunchao Zhang ConstPetscScalarKokkosOffsetView2D xv; 45c5566c22SJunchao Zhang PetscScalarKokkosOffsetView2D fv; 46*1693298eSJunchao Zhang Kokkos::DefaultExecutionSpace &exec = PetscGetKokkosExecutionSpace(); 47c5566c22SJunchao Zhang 48c5566c22SJunchao Zhang PetscFunctionBeginUser; 49c5566c22SJunchao Zhang lambda = user->param; 50c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 51c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 52c5566c22SJunchao Zhang hxdhy = hx / hy; 53c5566c22SJunchao Zhang hydhx = hy / hx; 54c5566c22SJunchao Zhang /* 55c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 56c5566c22SJunchao Zhang */ 573ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 583ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(info->da, f, &fv)); 59c5566c22SJunchao Zhang 609371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_for( 61*1693298eSJunchao Zhang "FormFunctionLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>(exec, {ys, xs}, {ys + ym, xs + xm}), KOKKOS_LAMBDA(PetscInt j, PetscInt i) { 62c5566c22SJunchao Zhang DMDACoor2d c; 63c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing; 64c5566c22SJunchao Zhang 65c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 669371c9d4SSatish Balay c.x = i * hx; 679371c9d4SSatish Balay c.y = j * hy; 683ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &mms_solution)); 69c5566c22SJunchao Zhang fv(j, i) = 2.0 * (hydhx + hxdhy) * (xv(j, i) - mms_solution); 70c5566c22SJunchao Zhang } else { 71c5566c22SJunchao Zhang u = xv(j, i); 72c5566c22SJunchao Zhang uw = xv(j, i - 1); 73c5566c22SJunchao Zhang ue = xv(j, i + 1); 74c5566c22SJunchao Zhang un = xv(j - 1, i); 75c5566c22SJunchao Zhang us = xv(j + 1, i); 76c5566c22SJunchao Zhang 77c5566c22SJunchao Zhang /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 789371c9d4SSatish Balay if (i - 1 == 0) { 799371c9d4SSatish Balay c.x = (i - 1) * hx; 809371c9d4SSatish Balay c.y = j * hy; 813ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &uw)); 829371c9d4SSatish Balay } 839371c9d4SSatish Balay if (i + 1 == mx - 1) { 849371c9d4SSatish Balay c.x = (i + 1) * hx; 859371c9d4SSatish Balay c.y = j * hy; 863ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &ue)); 879371c9d4SSatish Balay } 889371c9d4SSatish Balay if (j - 1 == 0) { 899371c9d4SSatish Balay c.x = i * hx; 909371c9d4SSatish Balay c.y = (j - 1) * hy; 913ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &un)); 929371c9d4SSatish Balay } 939371c9d4SSatish Balay if (j + 1 == my - 1) { 949371c9d4SSatish Balay c.x = i * hx; 959371c9d4SSatish Balay c.y = (j + 1) * hy; 963ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &us)); 979371c9d4SSatish Balay } 98c5566c22SJunchao Zhang 99c5566c22SJunchao Zhang uxx = (2.0 * u - uw - ue) * hydhx; 100c5566c22SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 101c5566c22SJunchao Zhang mms_forcing = 0; 1029371c9d4SSatish Balay c.x = i * hx; 1039371c9d4SSatish Balay c.y = j * hy; 1043ba16761SJacob Faibussowitsch static_cast<void>(MMSForcing1(user_param, &c, &mms_forcing)); 105c5566c22SJunchao Zhang fv(j, i) = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing); 106c5566c22SJunchao Zhang } 107c5566c22SJunchao Zhang })); 108c5566c22SJunchao Zhang 1093ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 1103ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(info->da, f, &fv)); 111c5566c22SJunchao Zhang 112c5566c22SJunchao Zhang PetscCall(PetscLogFlops(11.0 * info->ym * info->xm)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114c5566c22SJunchao Zhang } 115c5566c22SJunchao Zhang 116d71ae5a4SJacob Faibussowitsch PetscErrorCode FormObjectiveLocalVec(DMDALocalInfo *info, Vec x, PetscReal *obj, AppCtx *user) 117d71ae5a4SJacob Faibussowitsch { 118c5566c22SJunchao Zhang PetscInt xs = info->xs, ys = info->ys, xm = info->xm, ym = info->ym, mx = info->mx, my = info->my; 119c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0; 120c5566c22SJunchao Zhang MPI_Comm comm; 121c5566c22SJunchao Zhang 122c5566c22SJunchao Zhang ConstPetscScalarKokkosOffsetView2D xv; 123*1693298eSJunchao Zhang Kokkos::DefaultExecutionSpace &exec = PetscGetKokkosExecutionSpace(); 124c5566c22SJunchao Zhang 125c5566c22SJunchao Zhang PetscFunctionBeginUser; 126c5566c22SJunchao Zhang *obj = 0; 127c5566c22SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm)); 128c5566c22SJunchao Zhang lambda = user->param; 129c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(mx - 1); 130c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(my - 1); 131c5566c22SJunchao Zhang sc = hx * hy * lambda; 132c5566c22SJunchao Zhang hxdhy = hx / hy; 133c5566c22SJunchao Zhang hydhx = hy / hx; 134c5566c22SJunchao Zhang /* 135c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 136c5566c22SJunchao Zhang */ 1373ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 138c5566c22SJunchao Zhang 1399371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_reduce( 140*1693298eSJunchao Zhang "FormObjectiveLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>(exec, {ys, xs}, {ys + ym, xs + xm}), 1419371c9d4SSatish Balay KOKKOS_LAMBDA(PetscInt j, PetscInt i, PetscReal &update) { 142c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxux, uyuy; 143c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 144c5566c22SJunchao Zhang update += PetscRealPart((hydhx + hxdhy) * xv(j, i) * xv(j, i)); 145c5566c22SJunchao Zhang } else { 146c5566c22SJunchao Zhang u = xv(j, i); 147c5566c22SJunchao Zhang uw = xv(j, i - 1); 148c5566c22SJunchao Zhang ue = xv(j, i + 1); 149c5566c22SJunchao Zhang un = xv(j - 1, i); 150c5566c22SJunchao Zhang us = xv(j + 1, i); 151c5566c22SJunchao Zhang 152c5566c22SJunchao Zhang if (i - 1 == 0) uw = 0.; 153c5566c22SJunchao Zhang if (i + 1 == mx - 1) ue = 0.; 154c5566c22SJunchao Zhang if (j - 1 == 0) un = 0.; 155c5566c22SJunchao Zhang if (j + 1 == my - 1) us = 0.; 156c5566c22SJunchao Zhang 157c5566c22SJunchao Zhang /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 158c5566c22SJunchao Zhang 159c5566c22SJunchao Zhang uxux = u * (2. * u - ue - uw) * hydhx; 160c5566c22SJunchao Zhang uyuy = u * (2. * u - un - us) * hxdhy; 161c5566c22SJunchao Zhang 162c5566c22SJunchao Zhang update += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u)); 163c5566c22SJunchao Zhang } 1649371c9d4SSatish Balay }, 1659371c9d4SSatish Balay lobj)); 166c5566c22SJunchao Zhang 1673ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 168c5566c22SJunchao Zhang PetscCall(PetscLogFlops(12.0 * info->ym * info->xm)); 1690bf52853SStefano Zampini *obj = lobj; 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c5566c22SJunchao Zhang } 172c5566c22SJunchao Zhang 173d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocalVec(DMDALocalInfo *info, Vec x, Mat jac, Mat jacpre, AppCtx *user) 174d71ae5a4SJacob Faibussowitsch { 175c5566c22SJunchao Zhang PetscInt i, j; 176c5566c22SJunchao Zhang PetscInt xs = info->xs, ys = info->ys, xm = info->xm, ym = info->ym, mx = info->mx, my = info->my; 177c5566c22SJunchao Zhang MatStencil col[5], row; 178c5566c22SJunchao Zhang PetscScalar lambda, hx, hy, hxdhy, hydhx, sc; 179c5566c22SJunchao Zhang DM coordDA; 180c5566c22SJunchao Zhang Vec coordinates; 181c5566c22SJunchao Zhang DMDACoor2d **coords; 182c5566c22SJunchao Zhang 183c5566c22SJunchao Zhang PetscFunctionBeginUser; 184c5566c22SJunchao Zhang lambda = user->param; 185c5566c22SJunchao Zhang /* Extract coordinates */ 186c5566c22SJunchao Zhang PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 187c5566c22SJunchao Zhang PetscCall(DMGetCoordinates(info->da, &coordinates)); 188c5566c22SJunchao Zhang 189c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 190c5566c22SJunchao Zhang hx = xm > 1 ? PetscRealPart(coords[ys][xs + 1].x) - PetscRealPart(coords[ys][xs].x) : 1.0; 191c5566c22SJunchao Zhang hy = ym > 1 ? PetscRealPart(coords[ys + 1][xs].y) - PetscRealPart(coords[ys][xs].y) : 1.0; 192c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 193c5566c22SJunchao Zhang 194c5566c22SJunchao Zhang hxdhy = hx / hy; 195c5566c22SJunchao Zhang hydhx = hy / hx; 196c5566c22SJunchao Zhang sc = hx * hy * lambda; 197c5566c22SJunchao Zhang 198c5566c22SJunchao Zhang /* ----------------------------------------- */ 199c5566c22SJunchao Zhang /* MatSetPreallocationCOO() */ 200c5566c22SJunchao Zhang /* ----------------------------------------- */ 2010aeb1f43SStefano Zampini if (!user->ncoo) { 202c5566c22SJunchao Zhang PetscCount ncoo = ((PetscCount)xm) * ((PetscCount)ym) * 5; 203c5566c22SJunchao Zhang PetscInt *coo_i, *coo_j, *ip, *jp; 204c5566c22SJunchao Zhang PetscCall(PetscMalloc2(ncoo, &coo_i, ncoo, &coo_j)); /* 5-point stencil such that each row has at most 5 nonzeros */ 205c5566c22SJunchao Zhang 206c5566c22SJunchao Zhang ip = coo_i; 207c5566c22SJunchao Zhang jp = coo_j; 208c5566c22SJunchao Zhang for (j = ys; j < ys + ym; j++) { 209c5566c22SJunchao Zhang for (i = xs; i < xs + xm; i++) { 2109371c9d4SSatish Balay row.j = j; 2119371c9d4SSatish Balay row.i = i; 212c5566c22SJunchao Zhang /* Initialize neighbors with negative indices */ 213c5566c22SJunchao Zhang col[0].j = col[1].j = col[3].j = col[4].j = -1; 214c5566c22SJunchao Zhang /* boundary points */ 215c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 216c5566c22SJunchao Zhang col[2].j = row.j; 217c5566c22SJunchao Zhang col[2].i = row.i; 218c5566c22SJunchao Zhang } else { 219c5566c22SJunchao Zhang /* interior grid points */ 220c5566c22SJunchao Zhang if (j - 1 != 0) { 221c5566c22SJunchao Zhang col[0].j = j - 1; 222c5566c22SJunchao Zhang col[0].i = i; 223c5566c22SJunchao Zhang } 224c5566c22SJunchao Zhang 225c5566c22SJunchao Zhang if (i - 1 != 0) { 226c5566c22SJunchao Zhang col[1].j = j; 227c5566c22SJunchao Zhang col[1].i = i - 1; 228c5566c22SJunchao Zhang } 229c5566c22SJunchao Zhang 230c5566c22SJunchao Zhang col[2].j = row.j; 231c5566c22SJunchao Zhang col[2].i = row.i; 232c5566c22SJunchao Zhang 233c5566c22SJunchao Zhang if (i + 1 != mx - 1) { 234c5566c22SJunchao Zhang col[3].j = j; 235c5566c22SJunchao Zhang col[3].i = i + 1; 236c5566c22SJunchao Zhang } 237c5566c22SJunchao Zhang 238c5566c22SJunchao Zhang if (j + 1 != mx - 1) { 239c5566c22SJunchao Zhang col[4].j = j + 1; 240c5566c22SJunchao Zhang col[4].i = i; 241c5566c22SJunchao Zhang } 242c5566c22SJunchao Zhang } 243c5566c22SJunchao Zhang PetscCall(DMDAMapMatStencilToGlobal(info->da, 5, col, jp)); 244c5566c22SJunchao Zhang for (PetscInt k = 0; k < 5; k++) ip[k] = jp[2]; 245c5566c22SJunchao Zhang ip += 5; 246c5566c22SJunchao Zhang jp += 5; 247c5566c22SJunchao Zhang } 248c5566c22SJunchao Zhang } 249c5566c22SJunchao Zhang PetscCall(MatSetPreallocationCOO(jacpre, ncoo, coo_i, coo_j)); 250c5566c22SJunchao Zhang PetscCall(PetscFree2(coo_i, coo_j)); 2510aeb1f43SStefano Zampini user->ncoo = ncoo; 2520aeb1f43SStefano Zampini } 253c5566c22SJunchao Zhang 254c5566c22SJunchao Zhang /* ----------------------------------------- */ 255c5566c22SJunchao Zhang /* MatSetValuesCOO() */ 256c5566c22SJunchao Zhang /* ----------------------------------------- */ 2570aeb1f43SStefano Zampini PetscScalarKokkosView coo_v("coo_v", user->ncoo); 258c5566c22SJunchao Zhang ConstPetscScalarKokkosOffsetView2D xv; 259*1693298eSJunchao Zhang Kokkos::DefaultExecutionSpace &exec = PetscGetKokkosExecutionSpace(); 260c5566c22SJunchao Zhang 2613ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 262c5566c22SJunchao Zhang 2639371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_for( 264*1693298eSJunchao Zhang "FormJacobianLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>(exec, {ys, xs}, {ys + ym, xs + xm}), KOKKOS_LAMBDA(PetscCount j, PetscCount i) { 265c5566c22SJunchao Zhang PetscInt p = ((j - ys) * xm + (i - xs)) * 5; 266c5566c22SJunchao Zhang /* boundary points */ 267c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 268c5566c22SJunchao Zhang coo_v(p + 2) = 2.0 * (hydhx + hxdhy); 269c5566c22SJunchao Zhang } else { 270c5566c22SJunchao Zhang /* interior grid points */ 271ad540459SPierre Jolivet if (j - 1 != 0) coo_v(p + 0) = -hxdhy; 272ad540459SPierre Jolivet if (i - 1 != 0) coo_v(p + 1) = -hydhx; 273c5566c22SJunchao Zhang 274c5566c22SJunchao Zhang coo_v(p + 2) = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(xv(j, i)); 275c5566c22SJunchao Zhang 276ad540459SPierre Jolivet if (i + 1 != mx - 1) coo_v(p + 3) = -hydhx; 277ad540459SPierre Jolivet if (j + 1 != mx - 1) coo_v(p + 4) = -hxdhy; 278c5566c22SJunchao Zhang } 279c5566c22SJunchao Zhang })); 280c5566c22SJunchao Zhang PetscCall(MatSetValuesCOO(jacpre, coo_v.data(), INSERT_VALUES)); 2813ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 283c5566c22SJunchao Zhang } 2840aeb1f43SStefano Zampini 2850aeb1f43SStefano Zampini #else 2860aeb1f43SStefano Zampini #include "ex55.h" 2870aeb1f43SStefano Zampini 2880aeb1f43SStefano Zampini PetscErrorCode FormObjectiveLocalVec(DMDALocalInfo *info, Vec x, PetscReal *obj, AppCtx *user) 2890aeb1f43SStefano Zampini { 2900aeb1f43SStefano Zampini PetscFunctionBeginUser; 2910aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 2920aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2930aeb1f43SStefano Zampini } 2940aeb1f43SStefano Zampini 2950aeb1f43SStefano Zampini PetscErrorCode FormFunctionLocalVec(DMDALocalInfo *info, Vec x, Vec f, AppCtx *user) 2960aeb1f43SStefano Zampini { 2970aeb1f43SStefano Zampini PetscFunctionBeginUser; 2980aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 2990aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 3000aeb1f43SStefano Zampini } 3010aeb1f43SStefano Zampini 3020aeb1f43SStefano Zampini PetscErrorCode FormJacobianLocalVec(DMDALocalInfo *info, Vec x, Mat jac, Mat jacpre, AppCtx *user) 3030aeb1f43SStefano Zampini { 3040aeb1f43SStefano Zampini PetscFunctionBeginUser; 3050aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 3060aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 3070aeb1f43SStefano Zampini } 3080aeb1f43SStefano Zampini #endif 30977e5a1f9SBarry Smith 31077e5a1f9SBarry Smith /*TEST 31177e5a1f9SBarry Smith 31277e5a1f9SBarry Smith build: 3130338c944SBarry Smith TODO: 31477e5a1f9SBarry Smith 31577e5a1f9SBarry Smith TEST*/ 316