1c0558f20SBarry Smith static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel. Uses Kokkos\n\\n"; 2c0558f20SBarry Smith 341e469ddSJunchao Zhang #include <petscdmda_kokkos.hpp> 4c0558f20SBarry Smith #include <petscsnes.h> 5c0558f20SBarry Smith 6c0558f20SBarry Smith /* 7c0558f20SBarry Smith User-defined application context 8c0558f20SBarry Smith */ 9*9371c9d4SSatish Balay typedef struct { 10c0558f20SBarry Smith DM da; /* distributed array */ 11c0558f20SBarry Smith Vec F; /* right-hand-side of PDE */ 12c0558f20SBarry Smith PetscReal h; /* mesh spacing */ 13c0558f20SBarry Smith } ApplicationCtx; 14c0558f20SBarry Smith 15c0558f20SBarry Smith /* ------------------------------------------------------------------- */ 16c0558f20SBarry Smith /* 17c0558f20SBarry Smith FormInitialGuess - Computes initial guess. 18c0558f20SBarry Smith 19c0558f20SBarry Smith Input/Output Parameter: 20c0558f20SBarry Smith . x - the solution vector 21c0558f20SBarry Smith */ 22*9371c9d4SSatish Balay PetscErrorCode FormInitialGuess(Vec x) { 23c0558f20SBarry Smith PetscScalar pfive = .50; 24c0558f20SBarry Smith 25c0558f20SBarry Smith PetscFunctionBeginUser; 269566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 27c0558f20SBarry Smith PetscFunctionReturn(0); 28c0558f20SBarry Smith } 29c0558f20SBarry Smith 30c0558f20SBarry Smith /* ------------------------------------------------------------------- */ 31c0558f20SBarry Smith /* 3241e469ddSJunchao Zhang CpuFunction - Evaluates nonlinear function, F(x) on CPU 33c0558f20SBarry Smith 34c0558f20SBarry Smith Input Parameters: 35c0558f20SBarry Smith . snes - the SNES context 36c0558f20SBarry Smith . x - input vector 37c0558f20SBarry Smith . ctx - optional user-defined context, as set by SNESSetFunction() 38c0558f20SBarry Smith 39c0558f20SBarry Smith Output Parameter: 4041e469ddSJunchao Zhang . r - function vector 41c0558f20SBarry Smith 42c0558f20SBarry Smith Note: 43c0558f20SBarry Smith The user-defined context can contain any application-specific 44c0558f20SBarry Smith data needed for the function evaluation. 45c0558f20SBarry Smith */ 46*9371c9d4SSatish Balay PetscErrorCode CpuFunction(SNES snes, Vec x, Vec r, void *ctx) { 47c0558f20SBarry Smith ApplicationCtx *user = (ApplicationCtx *)ctx; 48c0558f20SBarry Smith DM da = user->da; 4941e469ddSJunchao Zhang PetscScalar *X, *R, *F, d; 50c0558f20SBarry Smith PetscInt i, M, xs, xm; 5141e469ddSJunchao Zhang Vec xl; 52c0558f20SBarry Smith 53c0558f20SBarry Smith PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &xl)); 559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl)); 569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, xl, &X)); 579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, r, &R)); 589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, user->F, &F)); 59c0558f20SBarry Smith 609566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 619566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 62c0558f20SBarry Smith 63c0558f20SBarry Smith if (xs == 0) { /* left boundary */ 6441e469ddSJunchao Zhang R[0] = X[0]; 6511486bccSBarry Smith xs++; 6611486bccSBarry Smith xm--; 67c0558f20SBarry Smith } 68c0558f20SBarry Smith if (xs + xm == M) { /* right boundary */ 6941e469ddSJunchao Zhang R[xs + xm - 1] = X[xs + xm - 1] - 1.0; 70c0558f20SBarry Smith xm--; 71c0558f20SBarry Smith } 72c0558f20SBarry Smith d = 1.0 / (user->h * user->h); 7341e469ddSJunchao Zhang for (i = xs; i < xs + xm; i++) R[i] = d * (X[i - 1] - 2.0 * X[i] + X[i + 1]) + X[i] * X[i] - F[i]; 74c0558f20SBarry Smith 759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, xl, &X)); 769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, r, &R)); 779566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, user->F, &F)); 789566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &xl)); 79c0558f20SBarry Smith PetscFunctionReturn(0); 80c0558f20SBarry Smith } 81c0558f20SBarry Smith 8241e469ddSJunchao Zhang using DefaultExecutionSpace = Kokkos::DefaultExecutionSpace; 8341e469ddSJunchao Zhang using DefaultMemorySpace = Kokkos::DefaultExecutionSpace::memory_space; 8441e469ddSJunchao Zhang using PetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<PetscScalar *, DefaultMemorySpace>; 8541e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<const PetscScalar *, DefaultMemorySpace>; 8641e469ddSJunchao Zhang 87*9371c9d4SSatish Balay PetscErrorCode KokkosFunction(SNES snes, Vec x, Vec r, void *ctx) { 8841e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx *)ctx; 8941e469ddSJunchao Zhang DM da = user->da; 9041e469ddSJunchao Zhang PetscScalar d; 9141e469ddSJunchao Zhang PetscInt M; 9241e469ddSJunchao Zhang Vec xl; 9341e469ddSJunchao Zhang PetscScalarKokkosOffsetView R; 9441e469ddSJunchao Zhang ConstPetscScalarKokkosOffsetView X, F; 9541e469ddSJunchao Zhang 9641e469ddSJunchao Zhang PetscFunctionBeginUser; 979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &xl)); 989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl)); 9941e469ddSJunchao Zhang d = 1.0 / (user->h * user->h); 1009566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da, xl, &X)); /* read only */ 1029566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(da, r, &R)); /* write only */ 1039566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da, user->F, &F)); /* read only */ 10411486bccSBarry Smith Kokkos::parallel_for( 10511486bccSBarry Smith Kokkos::RangePolicy<>(R.begin(0), R.end(0)), KOKKOS_LAMBDA(int i) { 10641e469ddSJunchao Zhang if (i == 0) R(0) = X(0); /* left boundary */ 10741e469ddSJunchao Zhang else if (i == M - 1) R(i) = X(i) - 1.0; /* right boundary */ 10841e469ddSJunchao Zhang else R(i) = d * (X(i - 1) - 2.0 * X(i) + X(i + 1)) + X(i) * X(i) - F(i); /* interior */ 10941e469ddSJunchao Zhang }); 1109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da, xl, &X)); 1119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(da, r, &R)); 1129566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da, user->F, &F)); 1139566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &xl)); 11441e469ddSJunchao Zhang PetscFunctionReturn(0); 11541e469ddSJunchao Zhang } 11641e469ddSJunchao Zhang 117*9371c9d4SSatish Balay PetscErrorCode StubFunction(SNES snes, Vec x, Vec r, void *ctx) { 11841e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx *)ctx; 11941e469ddSJunchao Zhang DM da = user->da; 12041e469ddSJunchao Zhang Vec rk; 12141e469ddSJunchao Zhang PetscReal norm = 0; 12241e469ddSJunchao Zhang 12341e469ddSJunchao Zhang PetscFunctionBeginUser; 1249566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &rk)); 1259566063dSJacob Faibussowitsch PetscCall(CpuFunction(snes, x, r, ctx)); 1269566063dSJacob Faibussowitsch PetscCall(KokkosFunction(snes, x, rk, ctx)); 1279566063dSJacob Faibussowitsch PetscCall(VecAXPY(rk, -1.0, r)); 1289566063dSJacob Faibussowitsch PetscCall(VecNorm(rk, NORM_2, &norm)); 1299566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &rk)); 13063a3b9bcSJacob Faibussowitsch PetscCheck(norm <= 1e-6, PETSC_COMM_SELF, PETSC_ERR_PLIB, "KokkosFunction() different from CpuFunction() with a diff norm = %g", (double)norm); 13141e469ddSJunchao Zhang PetscFunctionReturn(0); 13241e469ddSJunchao Zhang } 133c0558f20SBarry Smith /* ------------------------------------------------------------------- */ 134c0558f20SBarry Smith /* 135c0558f20SBarry Smith FormJacobian - Evaluates Jacobian matrix. 136c0558f20SBarry Smith 137c0558f20SBarry Smith Input Parameters: 138c0558f20SBarry Smith . snes - the SNES context 139c0558f20SBarry Smith . x - input vector 140c0558f20SBarry Smith . dummy - optional user-defined context (not used here) 141c0558f20SBarry Smith 142c0558f20SBarry Smith Output Parameters: 143c0558f20SBarry Smith . jac - Jacobian matrix 144c0558f20SBarry Smith . B - optionally different preconditioning matrix 145c0558f20SBarry Smith . flag - flag indicating matrix structure 146c0558f20SBarry Smith */ 147*9371c9d4SSatish Balay PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx) { 148c0558f20SBarry Smith ApplicationCtx *user = (ApplicationCtx *)ctx; 149c0558f20SBarry Smith PetscScalar *xx, d, A[3]; 150c0558f20SBarry Smith PetscInt i, j[3], M, xs, xm; 151c0558f20SBarry Smith DM da = user->da; 152c0558f20SBarry Smith 153c0558f20SBarry Smith PetscFunctionBeginUser; 154c0558f20SBarry Smith /* 155c0558f20SBarry Smith Get pointer to vector data 156c0558f20SBarry Smith */ 1579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, x, &xx)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 159c0558f20SBarry Smith 160c0558f20SBarry Smith /* 161c0558f20SBarry Smith Get range of locally owned matrix 162c0558f20SBarry Smith */ 1639566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 164c0558f20SBarry Smith 165c0558f20SBarry Smith /* 166c0558f20SBarry Smith Determine starting and ending local indices for interior grid points. 167c0558f20SBarry Smith Set Jacobian entries for boundary points. 168c0558f20SBarry Smith */ 169c0558f20SBarry Smith 170c0558f20SBarry Smith if (xs == 0) { /* left boundary */ 17111486bccSBarry Smith i = 0; 17211486bccSBarry Smith A[0] = 1.0; 173c0558f20SBarry Smith 1749566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 17511486bccSBarry Smith xs++; 17611486bccSBarry Smith xm--; 177c0558f20SBarry Smith } 178c0558f20SBarry Smith if (xs + xm == M) { /* right boundary */ 179c0558f20SBarry Smith i = M - 1; 180c0558f20SBarry Smith A[0] = 1.0; 1819566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 182c0558f20SBarry Smith xm--; 183c0558f20SBarry Smith } 184c0558f20SBarry Smith 185c0558f20SBarry Smith /* 186c0558f20SBarry Smith Interior grid points 187c0558f20SBarry Smith - Note that in this case we set all elements for a particular 188c0558f20SBarry Smith row at once. 189c0558f20SBarry Smith */ 190c0558f20SBarry Smith d = 1.0 / (user->h * user->h); 191c0558f20SBarry Smith for (i = xs; i < xs + xm; i++) { 19211486bccSBarry Smith j[0] = i - 1; 19311486bccSBarry Smith j[1] = i; 19411486bccSBarry Smith j[2] = i + 1; 19511486bccSBarry Smith A[0] = A[2] = d; 19611486bccSBarry Smith A[1] = -2.0 * d + 2.0 * xx[i]; 1979566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES)); 198c0558f20SBarry Smith } 199c0558f20SBarry Smith 200c0558f20SBarry Smith /* 201c0558f20SBarry Smith Assemble matrix, using the 2-step process: 202c0558f20SBarry Smith MatAssemblyBegin(), MatAssemblyEnd(). 203c0558f20SBarry Smith By placing code between these two statements, computations can be 204c0558f20SBarry Smith done while messages are in transition. 205c0558f20SBarry Smith 206c0558f20SBarry Smith Also, restore vector. 207c0558f20SBarry Smith */ 208c0558f20SBarry Smith 2099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, x, &xx)); 2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 212c0558f20SBarry Smith 213c0558f20SBarry Smith PetscFunctionReturn(0); 214c0558f20SBarry Smith } 215c0558f20SBarry Smith 216*9371c9d4SSatish Balay int main(int argc, char **argv) { 21741e469ddSJunchao Zhang SNES snes; /* SNES context */ 21841e469ddSJunchao Zhang Mat J; /* Jacobian matrix */ 21941e469ddSJunchao Zhang ApplicationCtx ctx; /* user-defined context */ 22041e469ddSJunchao Zhang Vec x, r, U, F; /* vectors */ 22141e469ddSJunchao Zhang PetscScalar none = -1.0; 22241e469ddSJunchao Zhang PetscInt its, N = 5, maxit, maxf; 22341e469ddSJunchao Zhang PetscReal abstol, rtol, stol, norm; 22441e469ddSJunchao Zhang PetscBool viewinitial = PETSC_FALSE; 22541e469ddSJunchao Zhang 226327415f7SBarry Smith PetscFunctionBeginUser; 2279566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL)); 22941e469ddSJunchao Zhang ctx.h = 1.0 / (N - 1); 2309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL)); 23141e469ddSJunchao Zhang 23241e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 23341e469ddSJunchao Zhang Create nonlinear solver context 23441e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2359566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 23641e469ddSJunchao Zhang 237c0558f20SBarry Smith /* 23841e469ddSJunchao Zhang Create distributed array (DMDA) to manage parallel grid and vectors 239c0558f20SBarry Smith */ 2409566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da)); 2419566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(ctx.da)); 2429566063dSJacob Faibussowitsch PetscCall(DMSetUp(ctx.da)); 24341e469ddSJunchao Zhang 24441e469ddSJunchao Zhang /* 24541e469ddSJunchao Zhang Extract global and local vectors from DMDA; then duplicate for remaining 24641e469ddSJunchao Zhang vectors that are the same types 24741e469ddSJunchao Zhang */ 2489566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(ctx.da, &x)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 2509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 25111486bccSBarry Smith PetscCall(VecDuplicate(x, &F)); 25211486bccSBarry Smith ctx.F = F; 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)F, "Forcing function")); 2549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 25641e469ddSJunchao Zhang 25741e469ddSJunchao Zhang /* 25841e469ddSJunchao Zhang Set function evaluation routine and vector. Whenever the nonlinear 25941e469ddSJunchao Zhang solver needs to compute the nonlinear function, it will call this 26041e469ddSJunchao Zhang routine. 26141e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 26241e469ddSJunchao Zhang context that provides application-specific data for the 26341e469ddSJunchao Zhang function evaluation routine. 26441e469ddSJunchao Zhang 265a5b23f4aSJose E. Roman At the beginning, one can use a stub function that checks the Kokkos version 26641e469ddSJunchao Zhang against the CPU version to quickly expose errors. 2679566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,StubFunction,&ctx)); 26841e469ddSJunchao Zhang */ 2699566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, KokkosFunction, &ctx)); 27041e469ddSJunchao Zhang 27141e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27241e469ddSJunchao Zhang Create matrix data structure; set Jacobian evaluation routine 27341e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2749566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(ctx.da, &J)); 27541e469ddSJunchao Zhang 27641e469ddSJunchao Zhang /* 27741e469ddSJunchao Zhang Set Jacobian matrix data structure and default Jacobian evaluation 27841e469ddSJunchao Zhang routine. Whenever the nonlinear solver needs to compute the 27941e469ddSJunchao Zhang Jacobian matrix, it will call this routine. 28041e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 28141e469ddSJunchao Zhang context that provides application-specific data for the 28241e469ddSJunchao Zhang Jacobian evaluation routine. 28341e469ddSJunchao Zhang */ 2849566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx)); 2859566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 2869566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf)); 28763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n", (double)abstol, (double)rtol, (double)stol, maxit, maxf)); 28841e469ddSJunchao Zhang 28941e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 29041e469ddSJunchao Zhang Initialize application: 29141e469ddSJunchao Zhang Store forcing function of PDE and exact solution 29241e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2939000cb80SJunchao Zhang { 2949000cb80SJunchao Zhang PetscScalarKokkosOffsetView FF, UU; 2959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, F, &FF)); 2969566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, U, &UU)); 29711486bccSBarry Smith Kokkos::parallel_for( 29811486bccSBarry Smith Kokkos::RangePolicy<>(FF.begin(0), FF.end(0)), KOKKOS_LAMBDA(int i) { 29941e469ddSJunchao Zhang PetscReal xp = i * ctx.h; 30041e469ddSJunchao Zhang FF(i) = 6.0 * xp + pow(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 30141e469ddSJunchao Zhang UU(i) = xp * xp * xp; 30241e469ddSJunchao Zhang }); 3039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, F, &FF)); 3049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, U, &UU)); 3059000cb80SJunchao Zhang } 30641e469ddSJunchao Zhang 30741e469ddSJunchao Zhang if (viewinitial) { 3089566063dSJacob Faibussowitsch PetscCall(VecView(U, NULL)); 3099566063dSJacob Faibussowitsch PetscCall(VecView(F, NULL)); 31041e469ddSJunchao Zhang } 31141e469ddSJunchao Zhang 31241e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 31341e469ddSJunchao Zhang Evaluate initial guess; then solve nonlinear system 31441e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 31541e469ddSJunchao Zhang 31641e469ddSJunchao Zhang /* 31741e469ddSJunchao Zhang Note: The user should initialize the vector, x, with the initial guess 31841e469ddSJunchao Zhang for the nonlinear solver prior to calling SNESSolve(). In particular, 31941e469ddSJunchao Zhang to employ an initial guess of zero, the user should explicitly set 32041e469ddSJunchao Zhang this vector to zero by calling VecSet(). 32141e469ddSJunchao Zhang */ 3229566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x)); 3239566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 3249566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 32563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 32641e469ddSJunchao Zhang 32741e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 32841e469ddSJunchao Zhang Check solution and clean up 32941e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 33041e469ddSJunchao Zhang /* 33141e469ddSJunchao Zhang Check the error 33241e469ddSJunchao Zhang */ 3339566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, none, U)); 3349566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 33563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its)); 33641e469ddSJunchao Zhang 33741e469ddSJunchao Zhang /* 33841e469ddSJunchao Zhang Free work space. All PETSc objects should be destroyed when they 33941e469ddSJunchao Zhang are no longer needed. 34041e469ddSJunchao Zhang */ 3419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 3459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 3469566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx.da)); 3489566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 349b122ec5aSJacob Faibussowitsch return 0; 35041e469ddSJunchao Zhang } 351c0558f20SBarry Smith 352c0558f20SBarry Smith /*TEST 353c0558f20SBarry Smith 354c0558f20SBarry Smith build: 35541e469ddSJunchao Zhang requires: kokkos_kernels 356c0558f20SBarry Smith 357c0558f20SBarry Smith test: 35887699089SJunchao Zhang requires: kokkos_kernels !complex !single 359c0558f20SBarry Smith nsize: 2 36041e469ddSJunchao Zhang args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor 361c0558f20SBarry Smith output_file: output/ex3k_1.out 362c0558f20SBarry Smith 363c0558f20SBarry Smith TEST*/ 364