10aeb1f43SStefano Zampini #include <petscconf.h> 20aeb1f43SStefano Zampini 30aeb1f43SStefano Zampini #if defined(PETSC_HAVE_KOKKOS_KERNELS) 4c5566c22SJunchao Zhang #include <Kokkos_Core.hpp> 5c5566c22SJunchao Zhang #include <petscdmda_kokkos.hpp> 6c5566c22SJunchao Zhang 7c5566c22SJunchao Zhang #include <petscdm.h> 8c5566c22SJunchao Zhang #include <petscdmda.h> 9c5566c22SJunchao Zhang #include <petscsnes.h> 10c5566c22SJunchao Zhang #include "ex55.h" 11c5566c22SJunchao Zhang 12c5566c22SJunchao Zhang using DefaultMemorySpace = Kokkos::DefaultExecutionSpace::memory_space; 13c5566c22SJunchao Zhang using ConstPetscScalarKokkosOffsetView2D = Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, DefaultMemorySpace>; 14c5566c22SJunchao Zhang using PetscScalarKokkosOffsetView2D = Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, DefaultMemorySpace>; 15c5566c22SJunchao Zhang 16c5566c22SJunchao Zhang using PetscCountKokkosView = Kokkos::View<PetscCount *, DefaultMemorySpace>; 17c5566c22SJunchao Zhang using PetscIntKokkosView = Kokkos::View<PetscInt *, DefaultMemorySpace>; 18c5566c22SJunchao Zhang using PetscCountKokkosViewHost = Kokkos::View<PetscCount *, Kokkos::HostSpace>; 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; 46c5566c22SJunchao Zhang 47c5566c22SJunchao Zhang PetscFunctionBeginUser; 48c5566c22SJunchao Zhang lambda = user->param; 49c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 50c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 51c5566c22SJunchao Zhang hxdhy = hx / hy; 52c5566c22SJunchao Zhang hydhx = hy / hx; 53c5566c22SJunchao Zhang /* 54c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 55c5566c22SJunchao Zhang */ 563ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 573ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(info->da, f, &fv)); 58c5566c22SJunchao Zhang 599371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_for( 609371c9d4SSatish Balay "FormFunctionLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>({ys, xs}, {ys + ym, xs + xm}), KOKKOS_LAMBDA(PetscInt j, PetscInt i) { 61c5566c22SJunchao Zhang DMDACoor2d c; 62c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing; 63c5566c22SJunchao Zhang 64c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 659371c9d4SSatish Balay c.x = i * hx; 669371c9d4SSatish Balay c.y = j * hy; 673ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &mms_solution)); 68c5566c22SJunchao Zhang fv(j, i) = 2.0 * (hydhx + hxdhy) * (xv(j, i) - mms_solution); 69c5566c22SJunchao Zhang } else { 70c5566c22SJunchao Zhang u = xv(j, i); 71c5566c22SJunchao Zhang uw = xv(j, i - 1); 72c5566c22SJunchao Zhang ue = xv(j, i + 1); 73c5566c22SJunchao Zhang un = xv(j - 1, i); 74c5566c22SJunchao Zhang us = xv(j + 1, i); 75c5566c22SJunchao Zhang 76c5566c22SJunchao Zhang /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 779371c9d4SSatish Balay if (i - 1 == 0) { 789371c9d4SSatish Balay c.x = (i - 1) * hx; 799371c9d4SSatish Balay c.y = j * hy; 803ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &uw)); 819371c9d4SSatish Balay } 829371c9d4SSatish Balay if (i + 1 == mx - 1) { 839371c9d4SSatish Balay c.x = (i + 1) * hx; 849371c9d4SSatish Balay c.y = j * hy; 853ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &ue)); 869371c9d4SSatish Balay } 879371c9d4SSatish Balay if (j - 1 == 0) { 889371c9d4SSatish Balay c.x = i * hx; 899371c9d4SSatish Balay c.y = (j - 1) * hy; 903ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &un)); 919371c9d4SSatish Balay } 929371c9d4SSatish Balay if (j + 1 == my - 1) { 939371c9d4SSatish Balay c.x = i * hx; 949371c9d4SSatish Balay c.y = (j + 1) * hy; 953ba16761SJacob Faibussowitsch static_cast<void>(MMSSolution1(user, &c, &us)); 969371c9d4SSatish Balay } 97c5566c22SJunchao Zhang 98c5566c22SJunchao Zhang uxx = (2.0 * u - uw - ue) * hydhx; 99c5566c22SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 100c5566c22SJunchao Zhang mms_forcing = 0; 1019371c9d4SSatish Balay c.x = i * hx; 1029371c9d4SSatish Balay c.y = j * hy; 1033ba16761SJacob Faibussowitsch static_cast<void>(MMSForcing1(user_param, &c, &mms_forcing)); 104c5566c22SJunchao Zhang fv(j, i) = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing); 105c5566c22SJunchao Zhang } 106c5566c22SJunchao Zhang })); 107c5566c22SJunchao Zhang 1083ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 1093ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(info->da, f, &fv)); 110c5566c22SJunchao Zhang 111c5566c22SJunchao Zhang PetscCall(PetscLogFlops(11.0 * info->ym * info->xm)); 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113c5566c22SJunchao Zhang } 114c5566c22SJunchao Zhang 115d71ae5a4SJacob Faibussowitsch PetscErrorCode FormObjectiveLocalVec(DMDALocalInfo *info, Vec x, PetscReal *obj, AppCtx *user) 116d71ae5a4SJacob Faibussowitsch { 117c5566c22SJunchao Zhang PetscInt xs = info->xs, ys = info->ys, xm = info->xm, ym = info->ym, mx = info->mx, my = info->my; 118c5566c22SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0; 119c5566c22SJunchao Zhang MPI_Comm comm; 120c5566c22SJunchao Zhang 121c5566c22SJunchao Zhang ConstPetscScalarKokkosOffsetView2D xv; 122c5566c22SJunchao Zhang 123c5566c22SJunchao Zhang PetscFunctionBeginUser; 124c5566c22SJunchao Zhang *obj = 0; 125c5566c22SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm)); 126c5566c22SJunchao Zhang lambda = user->param; 127c5566c22SJunchao Zhang hx = 1.0 / (PetscReal)(mx - 1); 128c5566c22SJunchao Zhang hy = 1.0 / (PetscReal)(my - 1); 129c5566c22SJunchao Zhang sc = hx * hy * lambda; 130c5566c22SJunchao Zhang hxdhy = hx / hy; 131c5566c22SJunchao Zhang hydhx = hy / hx; 132c5566c22SJunchao Zhang /* 133c5566c22SJunchao Zhang Compute function over the locally owned part of the grid 134c5566c22SJunchao Zhang */ 1353ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 136c5566c22SJunchao Zhang 1379371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_reduce( 1389371c9d4SSatish Balay "FormObjectiveLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>({ys, xs}, {ys + ym, xs + xm}), 1399371c9d4SSatish Balay KOKKOS_LAMBDA(PetscInt j, PetscInt i, PetscReal &update) { 140c5566c22SJunchao Zhang PetscScalar u, ue, uw, un, us, uxux, uyuy; 141c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 142c5566c22SJunchao Zhang update += PetscRealPart((hydhx + hxdhy) * xv(j, i) * xv(j, i)); 143c5566c22SJunchao Zhang } else { 144c5566c22SJunchao Zhang u = xv(j, i); 145c5566c22SJunchao Zhang uw = xv(j, i - 1); 146c5566c22SJunchao Zhang ue = xv(j, i + 1); 147c5566c22SJunchao Zhang un = xv(j - 1, i); 148c5566c22SJunchao Zhang us = xv(j + 1, i); 149c5566c22SJunchao Zhang 150c5566c22SJunchao Zhang if (i - 1 == 0) uw = 0.; 151c5566c22SJunchao Zhang if (i + 1 == mx - 1) ue = 0.; 152c5566c22SJunchao Zhang if (j - 1 == 0) un = 0.; 153c5566c22SJunchao Zhang if (j + 1 == my - 1) us = 0.; 154c5566c22SJunchao Zhang 155c5566c22SJunchao Zhang /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 156c5566c22SJunchao Zhang 157c5566c22SJunchao Zhang uxux = u * (2. * u - ue - uw) * hydhx; 158c5566c22SJunchao Zhang uyuy = u * (2. * u - un - us) * hxdhy; 159c5566c22SJunchao Zhang 160c5566c22SJunchao Zhang update += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u)); 161c5566c22SJunchao Zhang } 1629371c9d4SSatish Balay }, 1639371c9d4SSatish Balay lobj)); 164c5566c22SJunchao Zhang 1653ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 166c5566c22SJunchao Zhang PetscCall(PetscLogFlops(12.0 * info->ym * info->xm)); 1670bf52853SStefano Zampini *obj = lobj; 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 169c5566c22SJunchao Zhang } 170c5566c22SJunchao Zhang 171d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocalVec(DMDALocalInfo *info, Vec x, Mat jac, Mat jacpre, AppCtx *user) 172d71ae5a4SJacob Faibussowitsch { 173c5566c22SJunchao Zhang PetscInt i, j; 174c5566c22SJunchao Zhang PetscInt xs = info->xs, ys = info->ys, xm = info->xm, ym = info->ym, mx = info->mx, my = info->my; 175c5566c22SJunchao Zhang MatStencil col[5], row; 176c5566c22SJunchao Zhang PetscScalar lambda, hx, hy, hxdhy, hydhx, sc; 177c5566c22SJunchao Zhang DM coordDA; 178c5566c22SJunchao Zhang Vec coordinates; 179c5566c22SJunchao Zhang DMDACoor2d **coords; 180c5566c22SJunchao Zhang 181c5566c22SJunchao Zhang PetscFunctionBeginUser; 182c5566c22SJunchao Zhang lambda = user->param; 183c5566c22SJunchao Zhang /* Extract coordinates */ 184c5566c22SJunchao Zhang PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 185c5566c22SJunchao Zhang PetscCall(DMGetCoordinates(info->da, &coordinates)); 186c5566c22SJunchao Zhang 187c5566c22SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 188c5566c22SJunchao Zhang hx = xm > 1 ? PetscRealPart(coords[ys][xs + 1].x) - PetscRealPart(coords[ys][xs].x) : 1.0; 189c5566c22SJunchao Zhang hy = ym > 1 ? PetscRealPart(coords[ys + 1][xs].y) - PetscRealPart(coords[ys][xs].y) : 1.0; 190c5566c22SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 191c5566c22SJunchao Zhang 192c5566c22SJunchao Zhang hxdhy = hx / hy; 193c5566c22SJunchao Zhang hydhx = hy / hx; 194c5566c22SJunchao Zhang sc = hx * hy * lambda; 195c5566c22SJunchao Zhang 196c5566c22SJunchao Zhang /* ----------------------------------------- */ 197c5566c22SJunchao Zhang /* MatSetPreallocationCOO() */ 198c5566c22SJunchao Zhang /* ----------------------------------------- */ 1990aeb1f43SStefano Zampini if (!user->ncoo) { 200c5566c22SJunchao Zhang PetscCount ncoo = ((PetscCount)xm) * ((PetscCount)ym) * 5; 201c5566c22SJunchao Zhang PetscInt *coo_i, *coo_j, *ip, *jp; 202c5566c22SJunchao Zhang PetscCall(PetscMalloc2(ncoo, &coo_i, ncoo, &coo_j)); /* 5-point stencil such that each row has at most 5 nonzeros */ 203c5566c22SJunchao Zhang 204c5566c22SJunchao Zhang ip = coo_i; 205c5566c22SJunchao Zhang jp = coo_j; 206c5566c22SJunchao Zhang for (j = ys; j < ys + ym; j++) { 207c5566c22SJunchao Zhang for (i = xs; i < xs + xm; i++) { 2089371c9d4SSatish Balay row.j = j; 2099371c9d4SSatish Balay row.i = i; 210c5566c22SJunchao Zhang /* Initialize neighbors with negative indices */ 211c5566c22SJunchao Zhang col[0].j = col[1].j = col[3].j = col[4].j = -1; 212c5566c22SJunchao Zhang /* boundary points */ 213c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 214c5566c22SJunchao Zhang col[2].j = row.j; 215c5566c22SJunchao Zhang col[2].i = row.i; 216c5566c22SJunchao Zhang } else { 217c5566c22SJunchao Zhang /* interior grid points */ 218c5566c22SJunchao Zhang if (j - 1 != 0) { 219c5566c22SJunchao Zhang col[0].j = j - 1; 220c5566c22SJunchao Zhang col[0].i = i; 221c5566c22SJunchao Zhang } 222c5566c22SJunchao Zhang 223c5566c22SJunchao Zhang if (i - 1 != 0) { 224c5566c22SJunchao Zhang col[1].j = j; 225c5566c22SJunchao Zhang col[1].i = i - 1; 226c5566c22SJunchao Zhang } 227c5566c22SJunchao Zhang 228c5566c22SJunchao Zhang col[2].j = row.j; 229c5566c22SJunchao Zhang col[2].i = row.i; 230c5566c22SJunchao Zhang 231c5566c22SJunchao Zhang if (i + 1 != mx - 1) { 232c5566c22SJunchao Zhang col[3].j = j; 233c5566c22SJunchao Zhang col[3].i = i + 1; 234c5566c22SJunchao Zhang } 235c5566c22SJunchao Zhang 236c5566c22SJunchao Zhang if (j + 1 != mx - 1) { 237c5566c22SJunchao Zhang col[4].j = j + 1; 238c5566c22SJunchao Zhang col[4].i = i; 239c5566c22SJunchao Zhang } 240c5566c22SJunchao Zhang } 241c5566c22SJunchao Zhang PetscCall(DMDAMapMatStencilToGlobal(info->da, 5, col, jp)); 242c5566c22SJunchao Zhang for (PetscInt k = 0; k < 5; k++) ip[k] = jp[2]; 243c5566c22SJunchao Zhang ip += 5; 244c5566c22SJunchao Zhang jp += 5; 245c5566c22SJunchao Zhang } 246c5566c22SJunchao Zhang } 247c5566c22SJunchao Zhang PetscCall(MatSetPreallocationCOO(jacpre, ncoo, coo_i, coo_j)); 248c5566c22SJunchao Zhang PetscCall(PetscFree2(coo_i, coo_j)); 2490aeb1f43SStefano Zampini user->ncoo = ncoo; 2500aeb1f43SStefano Zampini } 251c5566c22SJunchao Zhang 252c5566c22SJunchao Zhang /* ----------------------------------------- */ 253c5566c22SJunchao Zhang /* MatSetValuesCOO() */ 254c5566c22SJunchao Zhang /* ----------------------------------------- */ 2550aeb1f43SStefano Zampini PetscScalarKokkosView coo_v("coo_v", user->ncoo); 256c5566c22SJunchao Zhang ConstPetscScalarKokkosOffsetView2D xv; 257c5566c22SJunchao Zhang 2583ba16761SJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(info->da, x, &xv)); 259c5566c22SJunchao Zhang 2609371c9d4SSatish Balay PetscCallCXX(Kokkos::parallel_for( 2619371c9d4SSatish Balay "FormFunctionLocalVec", MDRangePolicy<Rank<2, Iterate::Right, Iterate::Right>>({ys, xs}, {ys + ym, xs + xm}), KOKKOS_LAMBDA(PetscCount j, PetscCount i) { 262c5566c22SJunchao Zhang PetscInt p = ((j - ys) * xm + (i - xs)) * 5; 263c5566c22SJunchao Zhang /* boundary points */ 264c5566c22SJunchao Zhang if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 265c5566c22SJunchao Zhang coo_v(p + 2) = 2.0 * (hydhx + hxdhy); 266c5566c22SJunchao Zhang } else { 267c5566c22SJunchao Zhang /* interior grid points */ 268ad540459SPierre Jolivet if (j - 1 != 0) coo_v(p + 0) = -hxdhy; 269ad540459SPierre Jolivet if (i - 1 != 0) coo_v(p + 1) = -hydhx; 270c5566c22SJunchao Zhang 271c5566c22SJunchao Zhang coo_v(p + 2) = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(xv(j, i)); 272c5566c22SJunchao Zhang 273ad540459SPierre Jolivet if (i + 1 != mx - 1) coo_v(p + 3) = -hydhx; 274ad540459SPierre Jolivet if (j + 1 != mx - 1) coo_v(p + 4) = -hxdhy; 275c5566c22SJunchao Zhang } 276c5566c22SJunchao Zhang })); 277c5566c22SJunchao Zhang PetscCall(MatSetValuesCOO(jacpre, coo_v.data(), INSERT_VALUES)); 2783ba16761SJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(info->da, x, &xv)); 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 280c5566c22SJunchao Zhang } 2810aeb1f43SStefano Zampini 2820aeb1f43SStefano Zampini #else 2830aeb1f43SStefano Zampini #include "ex55.h" 2840aeb1f43SStefano Zampini 2850aeb1f43SStefano Zampini PetscErrorCode FormObjectiveLocalVec(DMDALocalInfo *info, Vec x, PetscReal *obj, AppCtx *user) 2860aeb1f43SStefano Zampini { 2870aeb1f43SStefano Zampini PetscFunctionBeginUser; 2880aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 2890aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2900aeb1f43SStefano Zampini } 2910aeb1f43SStefano Zampini 2920aeb1f43SStefano Zampini PetscErrorCode FormFunctionLocalVec(DMDALocalInfo *info, Vec x, Vec f, AppCtx *user) 2930aeb1f43SStefano Zampini { 2940aeb1f43SStefano Zampini PetscFunctionBeginUser; 2950aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 2960aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2970aeb1f43SStefano Zampini } 2980aeb1f43SStefano Zampini 2990aeb1f43SStefano Zampini PetscErrorCode FormJacobianLocalVec(DMDALocalInfo *info, Vec x, Mat jac, Mat jacpre, AppCtx *user) 3000aeb1f43SStefano Zampini { 3010aeb1f43SStefano Zampini PetscFunctionBeginUser; 3020aeb1f43SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Need to reconfigure with --download-kokkos-kernels"); 3030aeb1f43SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 3040aeb1f43SStefano Zampini } 3050aeb1f43SStefano Zampini #endif 30677e5a1f9SBarry Smith 30777e5a1f9SBarry Smith /*TEST 30877e5a1f9SBarry Smith 30977e5a1f9SBarry Smith build: 310*0338c944SBarry Smith TODO: 31177e5a1f9SBarry Smith 31277e5a1f9SBarry Smith TEST*/ 313